2D AMRCLAW
flagregions.f90
Go to the documentation of this file.
1 ! ::::::::::::::::::::: flagregions ::::::::::::::::::::::::::::::::::
2 !
3 ! Modify array of flagged points to respect minlevels and maxlevels
4 ! specified by regions.
5 !
6 !
7 ! amrflags = array to be flagged with either the value
8 ! DONTFLAG (no refinement needed) or
9 ! DOFLAG (refinement desired)
10 !
11 ! ::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
12 subroutine flagregions(mx,my,mbuff,xlower,ylower,dx,dy,level,t, &
13  amrflags,dontflag,doflag)
14 
15  use regions_module
16 
17  implicit none
18 
19  ! Subroutine arguments
20  integer, intent(in) :: mx,my,level,mbuff
21  real(kind=8), intent(in) :: xlower,ylower,dx,dy,t
22 
23  ! Flagging
24  real(kind=8),intent(inout) :: amrflags(1-mbuff:mx+mbuff,1-mbuff:my+mbuff)
25  real(kind=8), intent(in) :: dontflag
26  real(kind=8), intent(in) :: doflag
27 
28  logical :: allowflag
29  external allowflag
30 
31  ! Locals
32  integer :: i,j,m,minlevel,maxlevel
33  real(kind=8) :: x_low,y_low,x_hi,y_hi
34 
35 
36  ! Loop over interior points on this grid
37  ! (i,j) grid cell is [x_low,x_hi] x [y_low,y_hi]
38  do j=1,my
39  y_low = ylower + (j - 1) * dy
40  y_hi = ylower + j * dy
41 
42  do i = 1,mx
43  x_low = xlower + (i - 1) * dx
44  x_hi = xlower + i * dx
45 
46  minlevel = 0
47  maxlevel = 0
48 
49  do m=1,num_regions
50  if (t >= regions(m)%t_low .and. t <= regions(m)%t_hi .and. &
51  x_hi > regions(m)%x_low .and. x_low < regions(m)%x_hi .and. &
52  y_hi > regions(m)%y_low .and. y_low < regions(m)%y_hi ) then
53  minlevel = max(minlevel, regions(m)%min_level)
54  maxlevel = max(maxlevel, regions(m)%max_level)
55  endif
56  enddo
57 
58  if (minlevel > maxlevel) then
59  write(6,*) '*** Error: this should never happen!'
60  write(6,*) '*** minlevel > maxlevel in flagregions'
61  stop
62  endif
63 
64  if (maxlevel /= 0) then
65  ! this point lies in at least one region, so may need
66  ! to modify the exisiting flag at this point...
67  if (level < minlevel) then
68  ! Require refinement of this cell:
69  amrflags(i,j) = doflag
70  else if (level >= maxlevel) then
71  ! Do not refine of this cell:
72  amrflags(i,j) = dontflag
73  ! else leave amrflags(i,j) alone.
74  endif
75  endif
76 
77  enddo
78  enddo
79 
80 end subroutine flagregions
subroutine flagregions(mx, my, mbuff, xlower, ylower, dx, dy, level, t, amrflags, DONTFLAG, DOFLAG)
Definition: flagregions.f90:12
logical function allowflag(x, y, t, level)
Definition: allowflag.f:3