MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_b0.t
Go to the documentation of this file.
1 module mod_b0
2 
3  implicit none
4  private
5 
6 
7  public :: set_b0_grid
8 
9 contains
10 
11 
12  subroutine set_b0_grid(igrid)
14 
15  integer, intent(in) :: igrid
16 
17  integer :: ixcog^l
18 
19  ixcogmin^d=1;
20  ixcogmax^d=(ixghi^d-2*nghostcells)/2+2*nghostcells;
21 
22  call set_b0_cell(ps(igrid)%B0(ixg^t,:,0),ps(igrid)%x,ixg^ll,ixg^ll)
23  if(b0fieldalloccoarse) call set_b0_cell(psc(igrid)%B0(ixcog^s,:,0),psc(igrid)%x,ixcog^l,ixcog^l)
24  call set_j0_cell(igrid,ps(igrid)%J0,ixg^ll,ixm^ll^ladd1)
25  call set_b0_face(igrid,ps(igrid)%x,ixg^ll,ixm^ll)
26 
27  end subroutine set_b0_grid
28 
29  subroutine set_b0_cell(wB0,x,ixI^L,ix^L)
30  use mod_usr_methods, only: usr_set_b0
32  use mod_geometry
33 
34  integer, intent(in):: ixI^L,ix^L
35  double precision, intent(inout) :: wB0(ixI^S,1:ndir)
36  double precision, intent(in) :: x(ixI^S,1:ndim)
37 
38  wb0(ix^s,1:ndir)=zero
39 
40  ! approximate cell-averaged B0 as cell-centered B0
41  select case (coordinate)
42  case (spherical)
43  {^nooned
44  if (dabs(bdip)>smalldouble) then
45  wb0(ix^s,1)=2.0d0*bdip*dcos(x(ix^s,2))/x(ix^s,1)**3
46  wb0(ix^s,2)=bdip*dsin(x(ix^s,2))/x(ix^s,1)**3
47  end if
48 
49  if (abs(bquad)>smalldouble) then
50  wb0(ix^s,1)=wb0(ix^s,1) &
51  +bquad*0.5d0*(1.0d0+3.0d0*dcos(2.0d0*x(ix^s,2)))/x(ix^s,1)**4
52  wb0(ix^s,2)=wb0(ix^s,2)+bquad*dsin(2.0d0*x(ix^s,2))/x(ix^s,1)**4
53  end if
54  if (abs(boct)>smalldouble) then
55  wb0(ix^s,1)=wb0(ix^s,1) &
56  +boct*(10.0d0*dcos(2.0d0*x(ix^s,2))-2.0d0) &
57  *dcos(x(ix^s,2))/x(ix^s,1)**5
58  wb0(ix^s,2)=wb0(ix^s,2) &
59  +boct*1.5d0*(3.0d0+5.0d0*dcos(2.0d0*x(ix^s,2))) &
60  *dsin(x(ix^s,2))/x(ix^s,1)**5
61  end if
62  }
63  end select
64 
65  if (associated(usr_set_b0)) call usr_set_b0(ixi^l,ix^l,x,wb0)
66 
67  end subroutine set_b0_cell
68 
69  subroutine set_j0_cell(igrid,wJ0,ixI^L,ix^L)
70  use mod_usr_methods, only: usr_set_j0
72  use mod_geometry
73 
74  integer, intent(in):: igrid,ixI^L,ix^L
75  double precision, intent(inout) :: wJ0(ixI^S,7-2*ndir:3)
76  integer :: idirmin0, idirmin
77 
78  if(associated(usr_set_j0)) then
79  call usr_set_j0(ixi^l,ix^l,ps(igrid)%x,wj0)
80  else
81  idirmin0 = 7-2*ndir
82  call curlvector(ps(igrid)%B0(ixi^s,:,0),ixi^l,ix^l,wj0,idirmin,idirmin0,ndir)
83  end if
84 
85  end subroutine set_j0_cell
86 
87  subroutine set_b0_face(igrid,x,ixI^L,ixO^L)
89 
90  integer, intent(in) :: igrid, ixI^L, ixO^L
91  double precision, intent(in) :: x(ixI^S,1:ndim)
92 
93  double precision :: delx(ixI^S,1:ndim)
94  double precision :: xC(ixI^S,1:ndim),xshift^D
95  integer :: idims, ixC^L, hxO^L, ix, idims2
96 
97  if(slab_uniform)then
98  ^d&delx(ixi^s,^d)=rnode(rpdx^d_,igrid)\
99  else
100  ! for all non-cartesian and stretched cartesian coordinates
101  delx(ixi^s,1:ndim)=ps(igrid)%dx(ixi^s,1:ndim)
102  endif
103 
104 
105  do idims=1,ndim
106  hxo^l=ixo^l-kr(idims,^d);
107  if(stagger_grid) then
108  ! ct needs all transverse cells
109  ixcmax^d=ixomax^d+nghostcells-nghostcells*kr(idims,^d); ixcmin^d=hxomin^d-nghostcells+nghostcells*kr(idims,^d);
110  else
111  ! ixC is centered index in the idims direction from ixOmin-1/2 to ixOmax+1/2
112  ixcmax^d=ixomax^d; ixcmin^d=hxomin^d;
113  end if
114  ! always xshift=0 or 1/2
115  xshift^d=half*(one-kr(^d,idims));
116  do idims2=1,ndim
117  select case(idims2)
118  {case(^d)
119  do ix = ixc^lim^d
120  ! xshift=half: this is the cell center coordinate
121  ! xshift=0: this is the cell edge i+1/2 coordinate
122  xc(ix^d%ixC^s,^d)=x(ix^d%ixC^s,^d)+(half-xshift^d)*delx(ix^d%ixC^s,^d)
123  end do\}
124  end select
125  end do
126  call set_b0_cell(ps(igrid)%B0(ixi^s,:,idims),xc,ixi^l,ixc^l)
127  end do
128 
129  end subroutine set_b0_face
130 
131 
132 end module mod_b0
133 
Definition: mod_b0.t:1
subroutine, public set_b0_grid(igrid)
Definition: mod_b0.t:13
subroutine set_b0_cell(wB0, x, ixIL, ixL)
Definition: mod_b0.t:30
subroutine set_b0_face(igrid, x, ixIL, ixOL)
Definition: mod_b0.t:88
subroutine set_j0_cell(igrid, wJ0, ixIL, ixL)
Definition: mod_b0.t:70
Module with geometry-related routines (e.g., divergence, curl)
Definition: mod_geometry.t:2
integer coordinate
Definition: mod_geometry.t:7
integer, parameter spherical
Definition: mod_geometry.t:11
subroutine curlvector(qvec, ixIL, ixOL, curlvec, idirmin, idirmin0, ndir0, fourthorder)
Calculate curl of a vector qvec within ixL Options to employ standard second order CD evaluations use...
Definition: mod_geometry.t:663
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer ixghi
Upper index of grid block arrays.
integer, dimension(3, 3) kr
Kronecker delta tensor.
logical stagger_grid
True for using stagger grid.
double precision bdip
amplitude of background dipolar, quadrupolar, octupolar, user's field
integer, dimension(:), allocatable, parameter d
integer ixm
the mesh range of a physical block without ghost cells
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer nghostcells
Number of ghost cells surrounding a grid.
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
Module with all the methods that users can customize in AMRVAC.
procedure(set_j0), pointer usr_set_j0
procedure(set_b0), pointer usr_set_b0