mirror of
https://github.com/JHUAPL/kaiju.git
synced 2026-01-08 22:58:05 -05:00
74 lines
1.9 KiB
Plaintext
74 lines
1.9 KiB
Plaintext
module testmhd2mix
|
|
use testHelper
|
|
use voltapp
|
|
use gamapp
|
|
use uservoltic
|
|
|
|
implicit none
|
|
|
|
contains
|
|
|
|
@before
|
|
subroutine firstSerial()
|
|
end subroutine firstSerial
|
|
|
|
@after
|
|
subroutine lastSerial()
|
|
end subroutine lastSerial
|
|
|
|
@test
|
|
subroutine dipoleTest()
|
|
! testing that a dipole field from gamera results in 0 FAC in remix
|
|
type(voltApp_T) :: voltronApp
|
|
|
|
procedure(VectorField_T), pointer :: Axyz
|
|
character(len=strLen) :: caseInput = 'cmiD.xml'
|
|
real(rp) testValue
|
|
type(XML_Input_T) :: xmlInp
|
|
|
|
voltronApp%vOptions%gamUserInitFunc => initUser
|
|
call initVoltron(voltronApp, caseInput)
|
|
|
|
associate(gameraApp=>voltronApp%gApp)
|
|
|
|
! create a dipole field in Gamera. This code is copied from prob.F90
|
|
Axyz => VP_Dipole
|
|
call VectorPot2Flux(gameraApp%Model,gameraApp%Grid,gameraApp%State,Axyz)
|
|
|
|
call bFlux2Fld(gameraApp%Model, gameraApp%Grid, gameraApp%State%magFlux, gameraApp%State%Bxyz)
|
|
|
|
call convertGameraToRemix(voltronApp%mhd2mix, gameraApp, voltronApp%remixApp)
|
|
|
|
call mapGameraToRemix(voltronApp%mhd2mix, voltronApp%remixApp)
|
|
|
|
! verify that remix input FAC is all zeroes
|
|
testValue = sum(abs(voltronApp%remixApp%ion(NORTH)%St%Vars(:,:,FAC))) + sum(abs(voltronApp%remixApp%ion(SOUTH)%St%Vars(:,:,FAC)))
|
|
|
|
@assertLessThanOrEqual(testValue, 1e-1_rp, 'Remix did not get all 0 FAC from an input dipole field')
|
|
|
|
end associate
|
|
|
|
contains
|
|
|
|
! modified to hard-code M2 constant to 1.0 (therefore removed)
|
|
subroutine VP_Dipole(x,y,z,Ax,Ay,Az)
|
|
real(rp), intent(in) :: x,y,z
|
|
real(rp), intent(out) :: Ax,Ay,Az
|
|
|
|
real(rp), dimension(NDIM) :: A,m,r,rhat
|
|
m = [0,0,1]
|
|
r = [x,y,z]
|
|
rhat = r/norm2(r)
|
|
|
|
A = cross(m,rhat)/(dot_product(r,r))
|
|
Ax = A(XDIR)
|
|
Ay = A(YDIR)
|
|
Az = A(ZDIR)
|
|
|
|
end subroutine VP_Dipole
|
|
|
|
end subroutine dipoleTest
|
|
|
|
end module testmhd2mix
|
|
|