MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
mod_amr_grid.t
Go to the documentation of this file.
2 
3  implicit none
4  private
5 
6  public :: settree
7  public :: resettree
8  public :: resettree_convert
9 
10 contains
11 
12  !> Build up AMR
13  subroutine settree
17  use mod_errest, only: errest
18 
19  ! create and initialize grids on all levels > 1. On entry, all
20  ! level=1 grids have been formed and initialized. This subroutine
21  ! creates and initializes new level grids
22 
23  integer :: igrid, iigrid, levnew
24 
25  ! when only one level allowed, there is nothing to do anymore
26  if (refine_max_level == 1) return
27 
28  do levnew=2,refine_max_level
29 
30  call errest
31 
33 
34  if (.not.reset_grid) then
35  ! if no finer level grids created: exit
36  if (levmax/=levnew) exit
37  end if
38  end do
39 
40  ! set up boundary flux conservation arrays
41  if(levmax>levmin) call allocatebflux
42 
43  end subroutine settree
44 
45  !> reset AMR and (de)allocate boundary flux storage at level changes
46  subroutine resettree
49  use mod_amr_fct
51  use mod_errest, only: errest
52  use mod_particles
53 
56 
57  ! deallocate old grids for particles
58  if(use_particles) then
59  call finish_gridvars()
60  end if
61 
62  call errest
63 
65 
66  ! set up boundary flux conservation arrays
67  if(levmax>levmin) call allocatebflux
68 
69  ! generate new grids for particles
70  if(use_particles) then
71  call init_gridvars()
72  end if
73 
74  end subroutine resettree
75 
76  !> Force the tree to desired level(s) from level_io(_min/_max)
77  !> used for conversion to vtk output
78  subroutine resettree_convert
83 
84  integer :: igrid,iigrid, my_levmin, my_levmax
85 
86  if(level_io > 0) then
87  my_levmin = level_io
88  my_levmax = level_io
89  else
90  my_levmin = max(1,level_io_min)
91  my_levmax = min(refine_max_level,level_io_max)
92  end if
93 
94  do while(levmin<my_levmin.or.levmax>my_levmax)
95  call getbc(global_time,0.d0,ps,iwstart,nwgc)
96  do iigrid=1,igridstail; igrid=igrids(iigrid);
97  call forcedrefine_grid_io(igrid,ps(igrid)%w)
98  end do
99 
100  call amr_coarsen_refine
101  end do
102 
103  end subroutine resettree_convert
104 
105 end module mod_amr_grid
subroutine, public deallocatebfaces
Definition: mod_amr_fct.t:653
subroutine, public resettree
reset AMR and (de)allocate boundary flux storage at level changes
Definition: mod_amr_grid.t:47
subroutine, public resettree_convert
Force the tree to desired level(s) from level_io(_min/_max) used for conversion to vtk output.
Definition: mod_amr_grid.t:79
subroutine, public settree
Build up AMR.
Definition: mod_amr_grid.t:14
Module to coarsen and refine grids for AMR.
subroutine, public amr_coarsen_refine
coarsen and refine blocks to update AMR grid
subroutine, public errest
Do all local error estimation which determines (de)refinement.
Definition: mod_errest.t:12
subroutine, public forcedrefine_grid_io(igrid, w)
Definition: mod_errest.t:364
Module for flux conservation near refinement boundaries.
subroutine, public deallocatebflux
subroutine, public allocatebflux
update ghost cells of all blocks including physical boundaries
subroutine getbc(time, qdt, psb, nwstart, nwbc, req_diag)
do update ghost cells of all blocks including physical boundaries
This module contains definitions of global parameters and variables and some generic functions/subrou...
double precision global_time
The global simulation time.
logical stagger_grid
True for using stagger grid.
logical use_particles
Use particles module or not.
logical reset_grid
If true, rebuild the AMR grid upon restarting.
integer refine_max_level
Maximal number of AMR levels.
Module containing all the particle routines.
Definition: mod_particles.t:2