2D AMRCLAW
flagregions2.f90
Go to the documentation of this file.
1 ! ::::::::::::::::::::: flagregions ::::::::::::::::::::::::::::::::::
2 !
19 
20 !! amrflags = array to be flagged with either the value
21 !! DONTFLAG (no refinement needed) or
22 !! DOFLAG (refinement desired)
23 !!
24 !! \param mx number of cells in *i* direction
25 !! \param my number of cells in *j* direction
26 !! \param mbuff width of buffer region
27 !! \param maux number of auxiliary variables
28 !! \param xlower x-coordinate of left physical boundary
29 !! \param ylower y-coordinate of lower physical boundary
30 !! \param dx spacing in *i* direction
31 !! \param dy spacing in *j* direction
32 !! \param level AMR level of this grid
33 !! \param t simulation time on this grid
34 !! \param amrflags array to be flagged with either the value **DONTFLAG** or **DOFLAG** for each cell.
35 !! It is enlarged from grid size to include buffer regions around the grid.
36 !! \param DONTFLAG value to be assigned to amrflags for cells that need no refinement
37 !! \param DOFLAG value to be assigned to amrflags for cells that do need refinement
38 ! ::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
39 
40 subroutine flagregions2(mx,my,mbuff,xlower,ylower,dx,dy,level,t, &
41  amrflags,DONTFLAG,DOFLAG)
42 
43  use regions_module
44 
45  implicit none
46 
47  ! Subroutine arguments
48  integer, intent(in) :: mx,my,level,mbuff
49  real(kind=8), intent(in) :: xlower,ylower,dx,dy,t
50 
51  ! Flagging
52  real(kind=8),intent(inout) :: amrflags(1-mbuff:mx+mbuff,1-mbuff:my+mbuff)
53  real(kind=8), intent(in) :: DONTFLAG
54  real(kind=8), intent(in) :: DOFLAG
55 
56  ! Locals
57  integer :: i,j,m,i1,i2,j1,j2
58  real(kind=8) :: x_low,y_low,x_hi,y_hi, xupper,yupper
59  integer, allocatable :: minlevel(:,:), maxlevel(:,:)
60 
61  allocate(minlevel(mx,my), maxlevel(mx,my))
62 
63  minlevel = 0
64  maxlevel = 0
65 
66  xupper = xlower + mx*dx
67  yupper = ylower + my*dy
68 
69  rloop: do m=1,num_regions
70  if (t < regions(m)%t_low .or. t > regions(m)%t_hi) then
71  cycle rloop ! no intersection
72  endif
73 
74  if (xlower >= regions(m)%x_hi .or. xupper <= regions(m)%x_low) then
75  cycle rloop ! no intersection
76  else
77  i1 = max(floor((regions(m)%x_low - xlower) / dx) + 1, 1)
78  i2 = min(floor((regions(m)%x_hi -xlower) / dx) + 1, mx)
79  endif
80 
81  if (ylower >= regions(m)%y_hi .or. yupper <= regions(m)%y_low) then
82  cycle rloop ! no intersection
83  else
84  j1 = max(floor((regions(m)%y_low - ylower) / dy) + 1, 1)
85  j2 = min(floor((regions(m)%y_hi - ylower) / dy) + 1, my)
86  endif
87 
88  do j=j1,j2
89  do i=i1,i2
90  minlevel(i,j) = max(minlevel(i,j), regions(m)%min_level)
91  maxlevel(i,j) = max(maxlevel(i,j), regions(m)%max_level)
92  enddo
93  enddo
94  enddo rloop
95 
96  do j=1,my
97  do i=1,mx
98  if (minlevel(i,j) > maxlevel(i,j)) then
99  write(6,*) '*** Error: this should never happen!'
100  write(6,*) '*** minlevel > maxlevel in flagregions'
101  stop
102  endif
103 
104  if (maxlevel(i,j) /= 0) then
105  ! this point lies in at least one region, so may need
106  ! to modify the exisiting flag at this point...
107  if (level < minlevel(i,j)) then
108  ! Require refinement of this cell:
109  amrflags(i,j) = doflag
110  else if (level >= maxlevel(i,j)) then
111  ! Do not refine of this cell:
112  amrflags(i,j) = dontflag
113  ! else leave amrflags(i,j) alone.
114  endif
115  endif
116 
117  enddo
118  enddo
119 
120 end subroutine flagregions2
subroutine flagregions2(mx, my, mbuff, xlower, ylower, dx, dy, level, t, amrflags, DONTFLAG, DOFLAG)
Modify array of flagged points to respect minlevels and maxlevels specified by regions.
type(region_type), dimension(:), allocatable regions