collisions Module


Uses

  • module~~collisions~~UsesGraph module~collisions collisions iso_fortran_env iso_fortran_env module~collisions->iso_fortran_env module~phys_consts phys_consts module~collisions->module~phys_consts module~phys_consts->iso_fortran_env

Used by

  • module~~collisions~~UsedByGraph module~collisions collisions module~sources sources module~sources->module~collisions module~multifluid multifluid module~multifluid->module~collisions module~multifluid->module~sources module~potential_comm potential_comm module~potential_comm->module~collisions module~potential_root potential_root module~potential_root->module~potential_comm module~sources_mpi sources_mpi module~sources_mpi->module~sources program~gemini3d Gemini3D program~gemini3d->module~multifluid program~gemini3d->module~potential_comm module~potential_worker potential_worker module~potential_worker->module~potential_comm

Contents


Variables

TypeVisibility AttributesNameInitial
real(kind=wp), private, parameter:: Csn(lsp,ln) =reshape([-1.0_wp, 6.82e-10_wp, 6.64e-10_wp, -1.0_wp, 2.44e-10_wp, 4.34e-10_wp, 4.27e-10_wp, 0.69e-10_wp, 2.58e-10_wp, -1.0_wp, 4.49e-10_wp, 0.74e-10_wp, 2.31e-10_wp, 4.13e-10_wp, -1.0_wp, 0.65e-10_wp, 4.42e-10_wp, 7.47e-10_wp, 7.25e-10_wp, 1.45e-10_wp, -1.0_wp, 33.6e-10_wp, 32.0e-10_wp, -1.0_wp, -1.0_wp, -1.0_wp, -1.0_wp, -1.0_wp], shape(Csn), order=[2, 1])
real(kind=wp), private, parameter:: C2sn1(lsp,ln) =reshape([3.67e-11_wp, 0._wp, 0._wp, 4.63e-12_wp, 0._wp, 0._wp, 0._wp, 0._wp, 0._wp, 5.14e-11_wp, 0._wp, 0._wp, 0._wp, 0._wp, 2.59e-11_wp, 0._wp, 0._wp, 0._wp, 0._wp, 0._wp, 6.61e-11_wp, 0._wp, 0._wp, 2.65e-10_wp, -1._wp, -1._wp, -1._wp, -1._wp], shape(C2sn1), order=[2, 1])
real(kind=wp), private, parameter:: C2sn2(lsp,ln) =reshape([0.064_wp, 0._wp, 0._wp, -1._wp, 0._wp, 0._wp, 0._wp, 0._wp, 0._wp, 0.069_wp, 0._wp, 0._wp, 0._wp, 0._wp, 0.073_wp, 0._wp, 0._wp, 0._wp, 0._wp, 0._wp, 0.047_wp, 0._wp, 0._wp, 0.083_wp, -1._wp, -1._wp, -1._wp, -1._wp], shape(C2sn2), order=[2, 1])
real(kind=wp), private, parameter:: Csj(lsp,lsp) =reshape([0.22_wp, 0.26_wp, 0.25_wp, 0.26_wp, 0.22_wp, 0.077_wp, 1.87e-3_wp, 0.14_wp, 0.16_wp, 0.16_wp, 0.17_wp, 0.13_wp, 0.042_wp, 9.97e-4_wp, 0.15_wp, 0.17_wp, 0.17_wp, 0.18_wp, 0.14_wp, 0.045_wp, 1.07e-3_wp, 0.13_wp, 0.16_wp, 0.15_wp, 0.16_wp, 0.12_wp, 0.039_wp, 9.347e-4_wp, 0.25_wp, 0.28_wp, 0.28_wp, 0.28_wp, 0.24_wp, 0.088_wp, 2.136e-3_wp, 1.23_wp, 1.25_wp, 1.25_wp, 1.25_wp, 1.23_wp, 0.90_wp, 29.7e-3_wp, 54.5_wp, 54.5_wp, 54.5_wp, 54.5_wp, 54.5_wp, 54.5_wp, 38.537_wp], shape(Csj), order=[2, 1])

Subroutines

public subroutine maxwell_colln(isp, isp2, nn, Tn, Ts, nusn)

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: isp
integer, intent(in) :: isp2
real(kind=wp), intent(in), dimension(:,:,:,:):: nn
real(kind=wp), intent(in), dimension(:,:,:):: Tn
real(kind=wp), intent(in), dimension(-1:,-1:,-1:,:):: Ts
real(kind=wp), intent(out), dimension(1:size(Tn,1),1:size(Tn,2),1:size(Tn,3)):: nusn

public pure subroutine coulomb_colln(isp, isp2, ns, Ts, vs1, nusj, Phisj, Psisj)

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: isp
integer, intent(in) :: isp2
real(kind=wp), intent(in), dimension(-1:,-1:,-1:,:):: ns
real(kind=wp), intent(in), dimension(-1:,-1:,-1:,:):: Ts
real(kind=wp), intent(in), dimension(-1:,-1:,-1:,:):: vs1
real(kind=wp), intent(out), dimension(1:size(Ts,1)-4,1:size(Ts,2)-4,1:size(Ts,3)-4):: nusj
real(kind=wp), intent(out), dimension(1:size(Ts,1)-4,1:size(Ts,2)-4,1:size(Ts,3)-4):: Phisj
real(kind=wp), intent(out), dimension(1:size(Ts,1)-4,1:size(Ts,2)-4,1:size(Ts,3)-4):: Psisj

public subroutine thermal_conduct(isp, Ts, ns, nn, J1, lambda, beta)

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: isp
real(kind=wp), intent(in), dimension(-1:,-1:,-1:):: Ts
real(kind=wp), intent(in), dimension(-1:,-1:,-1:):: ns
real(kind=wp), intent(in), dimension(1:size(Ts,1)-4,1:size(Ts,2)-4,1:size(Ts,3)-4,ln):: nn
real(kind=wp), intent(in), dimension(1:size(Ts,1)-4,1:size(Ts,2)-4,1:size(Ts,3)-4):: J1
real(kind=wp), intent(out), dimension(1:size(Ts,1)-4,1:size(Ts,2)-4,1:size(Ts,3)-4):: lambda
real(kind=wp), intent(out), dimension(1:size(Ts,1)-4,1:size(Ts,2)-4,1:size(Ts,3)-4):: beta

public subroutine conductivities(nn, Tn, ns, Ts, vs1, B1, sig0, sigP, sigH, muP, muH, muPvn, muHvn)

cyclotron, abs() is sketch, needs to be checked. Basically a negative sign here is fine, while abs messes up direction of Hall current cyclotron, a negative sign from B1 here is fine for cartesian, but for dipole this should be the magnitude since the magnetic field is assumed to be along the x1-direction

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:,:,:,:):: nn
real(kind=wp), intent(in), dimension(:,:,:):: Tn
real(kind=wp), intent(in), dimension(-1:,-1:,-1:,:):: ns
real(kind=wp), intent(in), dimension(-1:,-1:,-1:,:):: Ts
real(kind=wp), intent(in), dimension(-1:,-1:,-1:,:):: vs1
real(kind=wp), intent(in), dimension(-1:,-1:,-1:):: B1
real(kind=wp), intent(out), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4):: sig0
real(kind=wp), intent(out), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4):: sigP
real(kind=wp), intent(out), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4):: sigH
real(kind=wp), intent(out), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4,lsp):: muP
real(kind=wp), intent(out), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4,lsp):: muH
real(kind=wp), intent(out), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4,lsp):: muPvn
real(kind=wp), intent(out), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4,lsp):: muHvn

public subroutine capacitance(ns, B1, flagcap, incap)

kludge the value to account for a magnetosheric contribution this is just a random guess that makes the KHI examples work well; a better value should be investigated

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(-1:,-1:,-1:,:):: ns
real(kind=wp), intent(in), dimension(-1:,-1:,-1:):: B1
integer, intent(in) :: flagcap
real(kind=wp), intent(out), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4):: incap