2D AMRCLAW
Functions/Subroutines
flagregions2.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

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. More...
 

Function/Subroutine Documentation

◆ flagregions2()

subroutine flagregions2 ( integer, intent(in)  mx,
integer, intent(in)  my,
integer, intent(in)  mbuff,
real(kind=8), intent(in)  xlower,
real(kind=8), intent(in)  ylower,
real(kind=8), intent(in)  dx,
real(kind=8), intent(in)  dy,
integer, intent(in)  level,
real(kind=8), intent(in)  t,
real(kind=8), dimension(1-mbuff:mx+mbuff,1-mbuff:my+mbuff), intent(inout)  amrflags,
real(kind=8), intent(in)  DONTFLAG,
real(kind=8), intent(in)  DOFLAG 
)

Modify array of flagged points to respect minlevels and maxlevels specified by regions.

This subroutine processes ONE grid associated with this amrflags

Second version with outer loop on regions, should be faster.

On input, amrflags is already set by flagger routine using Richardson extrapolation and/or flag2refine routine, as requested. This routine may change flags only in cells that are (partially) covered by one or more regions.

If any part of a grid cell is covered by one or more regions, then refinement is required to at least the max of all region min_levels and is allowed to at most to the max of all region max_levels.

Note that buffering is done after this, so additional cells may be refined in areas covered by regions that do not allow it!

Definition at line 42 of file flagregions2.f90.

References regions_module::num_regions, and regions_module::regions.

Referenced by bufnst2().

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 
real(kind=8) xupper
Definition: amr_module.f90:231
type(region_type), dimension(:), allocatable regions
real(kind=8) xlower
Definition: amr_module.f90:231
real(kind=8) yupper
Definition: amr_module.f90:231
real(kind=8) ylower
Definition: amr_module.f90:231
Here is the caller graph for this function: