MPI-AMRVAC  2.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_initialize.t
Go to the documentation of this file.
1 !> This module handles the initialization of various components of amrvac
3 
4  implicit none
5  private
6 
7  logical :: initialized_already = .false.
8 
9  ! Public methods
10  public :: initialize_amrvac
11 
12 contains
13 
14  !> Initialize amrvac: read par files and initialize variables
15  subroutine initialize_amrvac()
19  use mod_bc_data, only: bc_data_init
20 
21  if (initialized_already) return
22 
23  ! Check whether the user has loaded a physics module
24  call phys_check()
25 
26  ! Read input files
27  call read_par_files()
28  call initialize_vars()
29  call init_comm_types()
30 
31  ! Possibly load boundary condition data
32  call bc_data_init()
33 
34  if(associated(usr_set_parameters)) call usr_set_parameters()
35 
36  call phys_check_params()
37 
38  initialized_already = .true.
39  end subroutine initialize_amrvac
40 
41  !> Initialize (and allocate) simulation and grid variables
42  !> @todo Explain which ones are not initialized here
43  subroutine initialize_vars
47  use mod_fix_conserve, only: pflux
49  use mod_geometry
50 
51  integer :: igrid, level, ipe, ig^D
52  logical :: ok
53 
54  allocate(ps(max_blocks))
55  allocate(ps1(max_blocks))
56  allocate(ps2(max_blocks))
57  allocate(ps3(max_blocks))
58  allocate(ps4(max_blocks))
59  allocate(pso(max_blocks))
60  allocate(psc(max_blocks))
61  allocate(ps_sub(max_blocks))
62  allocate(neighbor(2,-1:1^d&,max_blocks),neighbor_child(2,0:3^d&,max_blocks))
63  allocate(neighbor_type(-1:1^d&,max_blocks),neighbor_active(-1:1^d&,max_blocks))
64  allocate(neighbor_pole(-1:1^d&,max_blocks))
65  allocate(igrids(max_blocks),igrids_active(max_blocks),igrids_passive(max_blocks))
68  allocate(pflux(2,^nd,max_blocks))
69  ! allocate mesh for particles
70  if(use_particles) allocate(gridvars(max_blocks))
71  if(stagger_grid) then
72  allocate(pface(2,^nd,max_blocks),fine_neighbors(2^d&,^nd,max_blocks))
73  allocate(old_neighbor(2,-1:1^d,max_blocks))
74  end if
75 
76  it=it_init
78 
79  dt=zero
80 
81  ! set all dt to zero
82  dt_grid(1:max_blocks)=zero
83 
84  ! no poles initially
85  neighbor_pole=0
86 
87  ! check resolution
88  if ({mod(ixghi^d,2)/=0|.or.}) then
89  call mpistop("mesh widths must give even number grid points")
90  end if
91  ixm^ll=ixg^ll^lsubnghostcells;
92 
93  if (nbufferx^d>(ixmhi^d-ixmlo^d+1)|.or.) then
94  write(unitterm,*) "nbufferx^D bigger than mesh size makes no sense."
95  write(unitterm,*) "Decrease nbufferx or increase mesh size"
96  call mpistop("")
97  end if
98 
99  ! initialize dx arrays on finer (>1) levels
100  do level=2,refine_max_level
101  {dx(^d,level) = dx(^d,level-1) * half\} ! refine ratio 2
102  end do
103 
104  ! domain decomposition
105  ! physical extent of a grid block at level 1, per dimension
106  ^d&dg^d(1)=dx(^d,1)*dble(block_nx^d)\
107  ! number of grid blocks at level 1 in simulation domain, per dimension
108  ^d&ng^d(1)=nint((xprobmax^d-xprobmin^d)/dg^d(1))\
109  ! total number of grid blocks at level 1
110  nglev1={ng^d(1)*}
111 
112  do level=2,refine_max_level
113  dg^d(level)=half*dg^d(level-1);
114  ng^d(level)=ng^d(level-1)*2;
115  end do
116 
117  ! check that specified stepsize correctly divides domain
118  ok=({(abs(dble(ng^d(1))*dg^d(1)-(xprobmax^d-xprobmin^d))<=smalldouble)|.and.})
119  if (.not.ok) then
120  write(unitterm,*)"domain cannot be divided by meshes of given gridsize"
121  call mpistop("domain cannot be divided by meshes of given gridsize")
122  end if
123 
124 
125  poleb=.false.
126  if (.not.slab) call set_pole
127 
128  ! initialize connectivity data
129  igridstail=0
130 
131  ! allocate memory for forest data structures
133  do level=1,refine_max_level
134  nullify(level_head(level)%node,level_tail(level)%node)
135  end do
136 
137  allocate(igrid_to_node(max_blocks,0:npe-1))
138  do ipe=0,npe-1
139  do igrid=1,max_blocks
140  nullify(igrid_to_node(igrid,ipe)%node)
141  end do
142  end do
143 
144  allocate(sfc(1:3,max_blocks*npe))
145 
146  allocate(igrid_to_sfc(max_blocks))
147 
148  sfc=0
149  allocate(morton_start(0:npe-1),morton_stop(0:npe-1))
150  allocate(morton_sub_start(0:npe-1),morton_sub_stop(0:npe-1))
151 
152  allocate(nleafs_level(1:nlevelshi))
153 
154  allocate(coarsen(max_blocks,0:npe-1),refine(max_blocks,0:npe-1))
155  coarsen=.false.
156  refine=.false.
157  if (nbufferx^d/=0|.or.) then
158  allocate(buffer(max_blocks,0:npe-1))
159  buffer=.false.
160  end if
161  allocate(igrid_inuse(max_blocks,0:npe-1))
162  igrid_inuse=.false.
163 
164  allocate(tree_root(1:ng^d(1)))
165  {do ig^db=1,ng^db(1)\}
166  nullify(tree_root(ig^d)%node)
167  {end do\}
168 
169  ! define index ranges and MPI send/receive derived datatype for ghost-cell swap
170  call init_bc()
177  call create_bc_mpi_datatype(1,nwflux+nwaux)
178 
179  end subroutine initialize_vars
180 
181 
182 end module mod_initialize
integer, dimension(-1:1^d &, 0:3^d &), target type_recv_p_f
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, parameter rnodehi
grid location info (corner coordinates and grid spacing)
Module for reading input and writing output.
integer, dimension(-1:2^d &,-1:1^d &), target type_recv_srl_f
update ghost cells of all blocks including physical boundaries
integer, dimension(:^d &,:^d &), pointer type_recv_srl
This module handles the initialization of various components of amrvac.
Definition: mod_initialize.t:2
subroutine create_bc_mpi_datatype(nwstart, nwbc)
integer, dimension(:^d &,:^d &), pointer type_send_r
integer max_blocks
The maximum number of grid blocks in a processor.
subroutine set_pole
Definition: mod_geometry.t:95
type(tree_node_ptr), dimension(:,:), allocatable, save igrid_to_node
Array to go from an [igrid, ipe] index to a node pointer.
Definition: mod_forest.t:32
subroutine initialize_vars
Initialize (and allocate) simulation and grid variables.
double precision, dimension(:), allocatable dg
extent of grid blocks in domain per dimension, in array over levels
integer, dimension(:,:), allocatable node_sub
double precision, dimension(:,:), allocatable rnode_sub
Module with geometry-related routines (e.g., divergence, curl)
Definition: mod_geometry.t:2
integer nglev1
Definition: mod_forest.t:78
integer it_init
initial iteration count
subroutine, public initialize_amrvac()
Initialize amrvac: read par files and initialize variables.
integer, dimension(-1:2^d &,-1:1^d &), target type_send_srl_f
Module with basic grid data structures.
Definition: mod_forest.t:2
Module for flux conservation near refinement boundaries.
type(fake_neighbors), dimension(:^d &,:,:), allocatable, public fine_neighbors
Definition: mod_amr_fct.t:16
integer, dimension(:^d &,:^d &), pointer type_recv_r
integer, dimension(:^d &,:^d &), pointer type_recv_p
integer npe
The number of MPI tasks.
integer, dimension(-1:1^d &, 0:3^d &), target type_recv_r_f
integer, dimension(:,:^d &,:), allocatable, public old_neighbor
Definition: mod_amr_fct.t:18
integer, dimension(:,:), allocatable, save sfc
Array to go from a Morton number to an igrid and processor index. Sfc(1:3, MN) contains [igrid...
Definition: mod_forest.t:43
integer, dimension(:), allocatable ng
number of grid blocks in domain per dimension, in array over levels
Module to set boundary conditions from user data.
Definition: mod_bc_data.t:2
logical, dimension(:,:), allocatable, save igrid_inuse
Definition: mod_forest.t:70
Module with all the methods that users can customize in AMRVAC.
integer, parameter nodehi
grid hierarchy info (level and grid indices)
integer ixghi
Upper index of grid block arrays.
logical, dimension(2, ndim) poleb
Indicates whether there is a pole at a boundary.
logical stagger_grid
True for using stagger grid.
type(tree_node_ptr), dimension(:^d &), allocatable, save tree_root
Pointers to the coarse grid.
Definition: mod_forest.t:29
subroutine, public bc_data_init()
Definition: mod_bc_data.t:37
double precision, dimension(:), allocatable dt_grid
procedure(p_no_args), pointer usr_set_parameters
Initialize the user&#39;s settings (after initializing amrvac)
type(tree_node_ptr), dimension(:), allocatable, save level_head
The head pointer of the linked list per refinement level.
Definition: mod_forest.t:35
subroutine phys_check()
Definition: mod_physics.t:296
integer, dimension(:), allocatable, save morton_start
First Morton number per processor.
Definition: mod_forest.t:62
integer, parameter unitterm
Unit for standard output.
integer nbufferx
Number of cells as buffer zone.
double precision, dimension(:,:), allocatable dx
logical, dimension(:,:), allocatable, save coarsen
AMR flags and grids-in-use identifier per processor (igrid,ipe)
Definition: mod_forest.t:70
integer, dimension(:), allocatable, save morton_stop
Last Morton number per processor.
Definition: mod_forest.t:65
logical, dimension(:,:), allocatable, save refine
Definition: mod_forest.t:70
double precision global_time
The global simulation time.
integer block_nx
number of cells for each dimension in grid block excluding ghostcells
subroutine read_par_files()
Read in the user-supplied parameter-file.
subroutine init_comm_types
Create and store the MPI types that will be used for parallel communication.
Definition: comm_lib.t:48
logical use_particles
Use particles module or not.
integer, dimension(:), allocatable, save nleafs_level
How many leaves are present per refinement level.
Definition: mod_forest.t:81
subroutine mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: comm_lib.t:198
integer, dimension(-1:1^d &, 0:3^d &), target type_send_p_f
logical, dimension(:,:), allocatable, save buffer
Definition: mod_forest.t:70
procedure(sub_check_params), pointer phys_check_params
Definition: mod_physics.t:59
type(tree_node_ptr), dimension(:), allocatable, save level_tail
The tail pointer of the linked list per refinement level.
Definition: mod_forest.t:38
integer, dimension(:^d &,:^d &), pointer type_send_p
logical slab
Cartesian geometry or not.
type(fluxalloc), dimension(:,:,:), allocatable, public pflux
store flux to fix conservation
integer, parameter nlevelshi
The maximum number of levels in the grid refinement.
type(facealloc), dimension(:,:,:), allocatable, public pface
Definition: mod_amr_fct.t:14
integer refine_max_level
Maximal number of AMR levels.
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
integer, dimension(:,:), allocatable node
integer, dimension(-1:1^d &,-1:1^d &), target type_send_r_f
integer, dimension(:^d &,:^d &), pointer type_send_srl
integer, dimension(:), allocatable, save igrid_to_sfc
Go from a grid index to Morton number (for a single processor)
Definition: mod_forest.t:56
integer it
Number of time steps taken.
logical, dimension(:), allocatable phyboundblock
True if a block has any physical boundary.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer, dimension(:), allocatable, save morton_sub_start
Definition: mod_forest.t:67
double precision time_init
Start time for the simulation.
integer, dimension(:), allocatable, save morton_sub_stop
Definition: mod_forest.t:67