MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_functions_bfield.t
Go to the documentation of this file.
2
3 implicit none
4 private
5
6 public :: get_divb
7
8 !> Indices of the magnetic field
9 integer, allocatable, public :: mag(:)
10
11contains
12
13 !> Calculate div B within ixO
14 subroutine get_divb(w,ixI^L,ixO^L,divb,nth_in)
16 use mod_geometry
17
18 integer, intent(in) :: ixi^l, ixo^l
19 double precision, intent(in) :: w(ixi^s,1:nw)
20 double precision, intent(inout) :: divb(ixi^s)
21 integer, intent(in), optional :: nth_in
22 integer :: ixc^l, idir, nth
23
24 if(present(nth_in)) then
25 nth=nth_in
26 else
27 nth=1
28 endif
29
30 if(stagger_grid) then
31 divb(ixo^s)=0.d0
32 do idir=1,ndim
33 ixc^l=ixo^l-kr(idir,^d);
34 divb(ixo^s)=divb(ixo^s)+block%ws(ixo^s,idir)*block%surfaceC(ixo^s,idir)-&
35 block%ws(ixc^s,idir)*block%surfaceC(ixc^s,idir)
36 end do
37 divb(ixo^s)=divb(ixo^s)/block%dvolume(ixo^s)
38 else
39 select case(typediv)
40 case("central")
41 call divvector(w(ixi^s,mag(1:ndir)),ixi^l,ixo^l,divb,nth)
42 case("limited")
43 call divvectors(w(ixi^s,mag(1:ndir)),ixi^l,ixo^l,divb)
44 end select
45 end if
46 end subroutine get_divb
47
48
49end module mod_functions_bfield
subroutine, public get_divb(w, ixil, ixol, divb, nth_in)
Calculate div B within ixO.
integer, dimension(:), allocatable, public mag
Indices of the magnetic field.
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
subroutine divvector(qvec, ixil, ixol, divq, nth_in)
subroutine divvectors(qvec, ixil, ixol, divq)
Calculate divergence of a vector qvec within ixL using limited extrapolation to cell edges.
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.
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical stagger_grid
True for using stagger grid.
character(len=std_len) typediv
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.