MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_initialize_amr.t
Go to the documentation of this file.
2 
3  implicit none
4  private
5 
6  public :: initlevelone
7  public :: modify_ic
8  public :: initial_condition
9  {^nooned
11  }
12 
13 
14 contains
15 
16 
17  !> Generate and initialize all grids at the coarsest level (level one)
18  subroutine initlevelone
24 
25  integer :: iigrid, igrid{#ifdef evolvingboundary , morton_no}
26 
27  levmin=1
28  levmax=1
29 
30  call init_forest_root
31 
32  call getigrids
34 
35  ! fill solution space of all root grids
36  do iigrid=1,igridstail; igrid=igrids(iigrid);
37  call alloc_node(igrid)
38  ! in case gradient routine used in initial condition, ensure geometry known
39  call initial_condition(igrid)
40  end do
41  {#IFDEF EVOLVINGBOUNDARY
42  ! mark physical-boundary blocks on space-filling curve
43  do morton_no=morton_start(mype),morton_stop(mype)
44  igrid=sfc_to_igrid(morton_no)
45  if (phyboundblock(igrid)) sfc_phybound(morton_no)=1
46  end do
47  call mpi_allreduce(mpi_in_place,sfc_phybound,nleafs,mpi_integer,&
48  mpi_sum,icomm,ierrmpi)
49  }
50 
51  ! update ghost cells
52  call getbc(global_time,0.d0,ps,iwstart,nwgc)
53 
54  end subroutine initlevelone
55 
56  !> fill in initial condition
57  subroutine initial_condition(igrid)
58  ! Need only to set the mesh values (can leave ghost cells untouched)
61  use mod_comm_lib, only: mpistop
62 
63  integer, intent(in) :: igrid
64 
65  ! in case gradient routine used in initial condition, ensure geometry known
66  block=>ps(igrid)
67  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
68 
69  if (.not. associated(usr_init_one_grid)) then
70  call mpistop("usr_init_one_grid not defined")
71  else
72  call usr_init_one_grid(ixg^ll,ixm^ll,ps(igrid)%w,ps(igrid)%x)
73  end if
74 
75  end subroutine initial_condition
76 
77  !> modify initial condition
78  subroutine modify_ic
81  use mod_comm_lib, only: mpistop
82 
83  integer :: iigrid, igrid
84 
85  do iigrid=1,igridstail; igrid=igrids(iigrid);
86  block=>ps(igrid)
87  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
88 
89  if (.not. associated(usr_init_one_grid)) then
90  call mpistop("usr_init_one_grid not defined")
91  else
92  call usr_init_one_grid(ixg^ll,ixm^ll,ps(igrid)%w,ps(igrid)%x)
93  end if
94  end do
95 
96  end subroutine modify_ic
97 
98  {^nooned
99  !> improve initial condition after initialization
102  use mod_usr_methods
105  use mod_physics
107 
108  logical :: active
109 
110  if(associated(usr_improve_initial_condition)) then
112  else if(stagger_grid) then
113  if(associated(usr_init_vector_potential)) then
114  ! re-calculate magnetic field from the vector potential in a
115  ! completely divergency free way for AMR mesh in 3D
116  if(levmax>levmin.and.ndim==3) call recalculateb
117  end if
118  if(slab_uniform.and.associated(phys_clean_divb)) then
119  ! Project out the divB using multigrid poisson solver
120  ! if not initialised from vector potential
121  if(.not.use_multigrid) call mg_setup_multigrid()
122  call phys_clean_divb(global_time,0.d0,active)
123  call getbc(global_time,0.d0,ps,iwstart,nwgc)
124  end if
125  end if
126 
127  end subroutine improve_initial_condition
128 }
129 
130 end module mod_initialize_amr
subroutine, public alloc_node(igrid)
allocate arrays on igrid node
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
subroutine recalculateb
re-calculate the magnetic field from the vector potential in a completely divergency free way
subroutine, public init_forest_root
build root forest
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...
type(state), pointer block
Block pointer for using one block and its previous state.
double precision global_time
The global simulation time.
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical stagger_grid
True for using stagger grid.
logical, dimension(:), allocatable phyboundblock
True if a block has any physical boundary.
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
integer, dimension(:), allocatable, parameter d
integer ixm
the mesh range of a physical block without ghost cells
integer ierrmpi
A global MPI error return code.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
logical use_multigrid
Use multigrid (only available in 2D and 3D)
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
double precision, dimension(ndim) dxlevel
subroutine, public improve_initial_condition()
improve initial condition after initialization
subroutine, public initlevelone
Generate and initialize all grids at the coarsest level (level one)
subroutine, public modify_ic
modify initial condition
subroutine, public initial_condition(igrid)
fill in initial condition
Module to couple the octree-mg library to AMRVAC. This file uses the VACPP preprocessor,...
subroutine mg_setup_multigrid()
Setup multigrid for usage.
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
procedure(sub_clean_divb), pointer phys_clean_divb
Definition: mod_physics.t:86
Module with all the methods that users can customize in AMRVAC.
procedure(p_no_args), pointer usr_improve_initial_condition
procedure(init_one_grid), pointer usr_init_one_grid
Initialize earch grid block data.
procedure(init_vector_potential), pointer usr_init_vector_potential