2D AMRCLAW
Functions/Subroutines
flag2refine2.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine flag2refine2 (mx, my, mbc, mbuff, meqn, maux, xlower, ylower, dx, dy, t, level, tolsp, q, aux, amrflags, DONTFLAG, DOFLAG)
 User routine to control flagging of points for refinement. More...
 

Function/Subroutine Documentation

◆ flag2refine2()

subroutine flag2refine2 ( integer, intent(in)  mx,
integer, intent(in)  my,
integer, intent(in)  mbc,
integer, intent(in)  mbuff,
integer, intent(in)  meqn,
integer, intent(in)  maux,
real(kind=8), intent(in)  xlower,
real(kind=8), intent(in)  ylower,
real(kind=8), intent(in)  dx,
real(kind=8), intent(in)  dy,
real(kind=8), intent(in)  t,
integer, intent(in)  level,
real(kind=8), intent(in)  tolsp,
real(kind=8), dimension(meqn,1-mbc:mx+mbc,1-mbc:my+mbc), intent(in)  q,
real(kind=8), dimension(maux,1-mbc:mx+mbc,1-mbc:my+mbc), intent(in)  aux,
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 
)

User routine to control flagging of points for refinement.

Default version computes spatial difference dq in each direction and for each component of q and flags any point where this is greater than the tolerance tolsp. This is consistent with what the routine errsp did in earlier versions of amrclaw (4.2 and before).

This routine can be copied to an application directory and modified to implement some other desired refinement criterion.

Points may also be flagged for refining based on a Richardson estimate of the error, obtained by comparing solutions on the current grid and a coarsened grid. Points are flagged if the estimated error is larger than the parameter tol in amr2ez.data, provided flag_richardson is .true., otherwise the coarsening and Richardson estimation is not performed! Points are flagged via Richardson in a separate routine.

Once points are flagged via this routine and/or Richardson, the subroutine flagregions is applied to check each point against the min_level and max_level of refinement specified in any "region" set by the user. So flags set here might be over-ruled by region constraints.

output: amrflags

Parameters
mxnumber of cells in i direction
mynumber of cells in j direction
mbcwidth of ghost cell region
mbuffwidth of buffer region
meqnnumber of equations for the system
mauxnumber of auxiliary variables
xlowerx-coordinate of left physical boundary
ylowery-coordinate of lower physical boundary
dxspacing in i direction
dyspacing in j direction
tsimulation time on this grid
levelAMR level of this grid
tolsptolerance specified by user in input file amr.data, used in default version of this routine as a tolerance for spatial differences
qgrid values including ghost cells (bndry vals at specified time have already been set, so can use ghost cell values too)
auxauxiliary array on this grid patch
amrflagsarray to be flagged with either the value DONTFLAG or DOFLAG for each cell. It is enlarged from grid size to include buffer regions around the grid.
DONTFLAGvalue to be assigned to amrflags for cells that need no refinement
DOFLAGvalue to be assigned to amrflags for cells that do need refinement

Definition at line 55 of file flag2refine2.f90.

Referenced by flagger(), and spest2().

55 
56  use regions_module
57 
58  implicit none
59 
60  ! Subroutine arguments
61  integer, intent(in) :: mx,my,mbc,meqn,maux,level,mbuff
62  real(kind=8), intent(in) :: xlower,ylower,dx,dy,t,tolsp
63 
64  real(kind=8), intent(in) :: q(meqn,1-mbc:mx+mbc,1-mbc:my+mbc)
65  real(kind=8), intent(in) :: aux(maux,1-mbc:mx+mbc,1-mbc:my+mbc)
66 
67  ! Flagging
68  real(kind=8),intent(inout) :: amrflags(1-mbuff:mx+mbuff,1-mbuff:my+mbuff)
69  real(kind=8), intent(in) :: dontflag
70  real(kind=8), intent(in) :: doflag
71 
72  logical :: allowflag
73  external allowflag
74 
75  ! Locals
76  integer :: i,j,m
77  real(kind=8) :: x_c,y_c,x_low,y_low,x_hi,y_hi
78  real(kind=8) :: dqi(meqn), dqj(meqn), dq(meqn)
79 
80  ! Initialize flags
81  amrflags = dontflag
82 
83  ! Loop over interior points on this grid
84  ! (i,j) grid cell is [x_low,x_hi] x [y_low,y_hi], cell center at (x_c,y_c)
85  ! This information is not needed for the default flagging based on
86  ! undivided differences, but might be needed in a user's version.
87  ! Note that if you want to refine only in certain space-time regions,
88  ! it may be easiest to use the "regions" feature. The flags set here or
89  ! in the Richardson error estimator are potentially modified by the
90  ! min_level and max_level specified in any regions.
91 
92  y_loop: do j=1,my
93  y_low = ylower + (j - 1) * dy
94  y_c = ylower + (j - 0.5d0) * dy
95  y_hi = ylower + j * dy
96 
97  x_loop: do i = 1,mx
98  x_low = xlower + (i - 1) * dx
99  x_c = xlower + (i - 0.5d0) * dx
100  x_hi = xlower + i * dx
101 
102  ! -----------------------------------------------------------------
103  dq = 0.d0
104  dqi = abs(q(:,i+1,j) - q(:,i-1,j))
105  dqj = abs(q(:,i,j+1) - q(:,i,j-1))
106  dq = max(dq,dqi,dqj)
107 
108  ! default checks all components of undivided difference:
109  do m=1,meqn
110  if (dq(m) > tolsp) then
111  amrflags(i,j) = doflag
112  cycle x_loop
113  endif
114  enddo
115 
116  enddo x_loop
117  enddo y_loop
118 
real(kind=8) xlower
Definition: amr_module.f90:231
real(kind=8) tolsp
Definition: amr_module.f90:197
real(kind=8) ylower
Definition: amr_module.f90:231
logical function allowflag(x, y, t, level)
Indicate whether the grid point at (x,y,t) at this refinement level is allowed to be flagged for furt...
Definition: allowflag.f:20
Here is the caller graph for this function: