src/equilibrium.f90: use enums

This commit is contained in:
Michele Guerini Rocco 2024-01-30 10:14:24 +01:00
parent fac0c6ded8
commit 7c5b443847
Signed by: rnhmjoj
GPG Key ID: BFBAF4C975F76450
4 changed files with 380 additions and 329 deletions

View File

@ -8,6 +8,8 @@
module equilibrium
use const_and_precisions, only : wp_
use splines, only : spline_simple, spline_1d, spline_2d, linear_1d
use gray_params, only : EQ_VACUUM, EQ_ANALYTICAL, &
EQ_EQDSK_FULL, EQ_EQDSK_PARTIAL
implicit none
@ -360,9 +362,8 @@ contains
! 3. q = 1/2π Φ/ψ ~ Φ/rr/ψ < 0.
! 4. In general, sgn(q) = -sgn(I_p)sgn(B_φ).
if (iequil < 2) then
! Analytical model
select case(iequil)
case (EQ_ANALYTICAL) ! Analytical model
! Apply signs
if (params%sgni /= 0) then
model%q0 = sign(model%q0, -params%sgni*one)
@ -373,9 +374,8 @@ contains
end if
! Rescale
model%B0 = model%B0 * params%factb
else
! Numeric data
case (EQ_EQDSK_FULL, EQ_EQDSK_PARTIAL) ! Numeric data
! Apply signs
if (params%sgni /= 0) &
data%psia = sign(data%psia, -params%sgni*one)
@ -389,7 +389,7 @@ contains
! Note: In these cases the values sgni,sgnb from gray.ini are unused.
params%sgni = int(sign(one, -data%psia))
params%sgnb = int(sign(one, +data%fpol(size(data%fpol))))
end if
end select
end subroutine scale_equil
@ -398,7 +398,7 @@ contains
! in their respective global variables, see the top of this file.
use const_and_precisions, only : zero, one
use gray_params, only : equilibrium_parameters, equilibrium_data
use gray_params, only : iequil
use gray_params, only : iequil, X_AT_TOP, X_AT_BOTTOM, X_IS_MISSING
use utils, only : vmaxmin, vmaxmini, inside
use logger, only : log_info
@ -434,8 +434,10 @@ contains
! Spline interpolation of ψ(R, z)
if (iequil>2) then
! data valid only inside boundary (data%psin=0 outside), e.g. source==ESCO
select case (iequil)
case (EQ_EQDSK_PARTIAL)
! Data valid only inside boundary (data%psin=0 outside),
! presence of boundary anticipated here to filter invalid data
nbnd = min(size(data%rbnd), size(data%zbnd))
@ -503,8 +505,9 @@ contains
deallocate(rv1d, zv1d, wf, fvpsi)
! reset nrz to the total number of grid points for next allocations
nrz = nr*nz
else
! iequil==2: data are valid on the full R,z grid
case (EQ_EQDSK_FULL)
! Data are valid on the full R,z grid
! reshape 2D ψ array to 1D (transposed)
allocate(fvpsi(nrz))
@ -523,7 +526,7 @@ contains
err = 0
end if
deallocate(fvpsi)
end if
end select
if (err /= 0) then
err = 2
@ -592,11 +595,11 @@ contains
'r', rmaxis, 'z', zmaxis, 'ψ', psinoptmp
call log_info(msg, mod='equilibrium', proc='set_equil_spline')
! search for X-point if params%ixp /= 0
! Search for X-point
ixploc = params%ixp
if(ixploc/=0) then
if(ixploc<0) then
select case (ixploc)
case (X_AT_BOTTOM)
call points_ox(rbinf,zbinf,r1,z1,psinxptmp,info)
if(psinxptmp/=-1.0_wp_) then
write (msg, '("X-point found:", 3(x,a,"=",g0.3))') &
@ -611,7 +614,8 @@ contains
else
ixploc=0
end if
else
case (X_AT_TOP)
call points_ox(rbsup,zbsup,r1,z1,psinxptmp,info)
if(psinxptmp.ne.-1.0_wp_) then
write (msg, '("X-point found:", 3(x,a,"=",g0.3))') &
@ -626,18 +630,14 @@ contains
else
ixploc=0
end if
end if
end if
if (ixploc==0) then
case (X_IS_MISSING)
psinop=psinoptmp
psiant=one-psinop
! Find upper horizontal tangent point
call points_tgo(rmaxis,0.5_wp_*(zmaxis+zbsup),r1,z1,one,info)
zbsup=z1
rbsup=r1
! Find lower horizontal tangent point
call points_tgo(rmaxis,0.5_wp_*(zmaxis+zbinf),r1,z1,one,info)
zbinf=z1
@ -645,7 +645,7 @@ contains
write (msg, '("X-point not found in", 2(x,a,"∈[",g0.3,",",g0.3,"]"))') &
'r', rbinf, rbsup, 'z', zbinf, zbsup
call log_info(msg, mod='equilibrium', proc='set_equil_spline')
end if
end select
! Adjust all the B-spline coefficients
! Note: since ψ_n(R,z) = Σ_ij c_ij B_i(R)B_j(z), to correct ψ_n
@ -908,7 +908,16 @@ contains
real(wp_) :: dqdr, dqdz ! q
real(wp_) :: dphidr2, ddphidr2dr2 ! dΦ_n/d(r²), d²Φ_n/d(r²)²
if (iequil < 2) then
! Values for vacuum/outside the domain
if (present(psi_n)) psi_n = -1
if (present(dpsidr)) dpsidr = 0
if (present(dpsidz)) dpsidz = 0
if (present(ddpsidrr)) ddpsidrr = 0
if (present(ddpsidzz)) ddpsidzz = 0
if (present(ddpsidrz)) ddpsidrz = 0
select case (iequil)
case (EQ_ANALYTICAL)
! Analytical model
!
! The normalised poloidal flux ψ_n(R, z) is computed as follows:
@ -1010,7 +1019,8 @@ contains
if (present(ddpsidrr)) ddpsidrr = ddpsidphidr * dphidr + dpsidphi * ddphidrdr
if (present(ddpsidzz)) ddpsidzz = ddpsidphidz * dphidz + dpsidphi * ddphidzdz
if (present(ddpsidrz)) ddpsidrz = ddpsidphidr * dphidz + dpsidphi * ddphidrdz
else
case (EQ_EQDSK_FULL, EQ_EQDSK_PARTIAL)
! Numerical data
if (inside(psi_domain%R, psi_domain%z, R, z)) then
! Within the interpolation range
@ -1020,16 +1030,9 @@ contains
if (present(ddpsidrr)) ddpsidrr = psi_spline%deriv(R, z, 2, 0)
if (present(ddpsidzz)) ddpsidzz = psi_spline%deriv(R, z, 0, 2)
if (present(ddpsidrz)) ddpsidrz = psi_spline%deriv(R, z, 1, 1)
else
! Outside
if (present(psi_n)) psi_n = -1
if (present(dpsidr)) dpsidr = 0
if (present(dpsidz)) dpsidz = 0
if (present(ddpsidrr)) ddpsidrr = 0
if (present(ddpsidzz)) ddpsidzz = 0
if (present(ddpsidrz)) ddpsidrz = 0
end if
end if
end select
end subroutine pol_flux
@ -1043,12 +1046,19 @@ contains
real(wp_), intent(out) :: fpol ! poloidal current
real(wp_), intent(out), optional :: dfpol ! derivative
if (iequil < 2) then
select case (iequil)
case (EQ_VACUUM)
! Vacuum, no plasma
fpol = 0
if (present(dfpol)) dfpol = 0
case (EQ_ANALYTICAL)
! Analytical model
! F(ψ) = BR, a constant
fpol = model%B0 * model%R0
if (present(dfpol)) dfpol = 0
else
case (EQ_EQDSK_FULL, EQ_EQDSK_PARTIAL)
! Numerical data
if(psi_n <= 1 .and. psi_n >= 0) then
fpol = fpol_spline%eval(psi_n)
@ -1057,7 +1067,7 @@ contains
fpol = fpolas
if (present(dfpol)) dfpol = 0
end if
end if
end select
end subroutine pol_curr
@ -1069,7 +1079,9 @@ contains
real(wp_), intent(in) :: rho_p
real(wp_) :: frhotor
if (iequil < 2) then
select case (iequil)
case (EQ_ANALYTICAL)
! Analytical model
block
! The change of variable is obtained by integrating
@ -1089,10 +1101,12 @@ contains
frhotor = rho_p * sqrt((model%q0 + dq*rho_p**model%alpha) &
/ (model%q0 + dq))
end block
else
case (EQ_EQDSK_FULL, EQ_EQDSK_PARTIAL)
! Numerical data
frhotor = rhot_spline%eval(rho_p)
end if
end select
end function frhotor
@ -1105,7 +1119,12 @@ contains
real(wp_), intent(in) :: rho_t
real(wp_) :: frhopol
if (iequil < 2) then
select case (iequil)
case (EQ_VACUUM)
! Vacuum, no plasma
frhopol = 0
case (EQ_ANALYTICAL)
! Analytical model
block
! In general there is no closed form for ρ_p(ρ_t) in the
@ -1121,10 +1140,11 @@ contains
ldfjac=1, tol=comp_eps, info=info, wa=wa, lwa=7)
frhopol = rho_p(1)
end block
else
case (EQ_EQDSK_FULL, EQ_EQDSK_PARTIAL)
! Numerical data
frhopol = rhop_spline%eval(rho_t)
end if
end select
contains
@ -1162,7 +1182,12 @@ contains
real(wp_), intent(in) :: psin
real(wp_) :: fq
if (iequil < 2) then
select case(iequil)
case (EQ_VACUUM)
! Vacuum, q is undefined
fq = 0
case (EQ_ANALYTICAL)
! Analytical model
! The safety factor is a power law in ρ_p:
! q(ρ_p) = q + (q-q) ρ_p^α
@ -1171,10 +1196,11 @@ contains
rho_p = sqrt(psin)
fq = abs(model%q0 + (model%q1 - model%q0) * rho_p**model%alpha)
end block
else
case (EQ_EQDSK_FULL, EQ_EQDSK_PARTIAL)
! Numerical data
fq = q_spline%eval(psin)
end if
end select
end function fq
@ -1183,6 +1209,7 @@ contains
! (R, z) in cylindrical coordinates
!
! Note: all output arguments are optional.
use gray_params, only : iequil
! subroutine arguments
real(wp_), intent(in) :: R, z
@ -1191,6 +1218,14 @@ contains
! local variables
real(wp_) :: psi_n, fpol, dpsidr, dpsidz
if (iequil == EQ_VACUUM) then
! Vacuum, no plasma nor field
if (present(B_R)) B_R = 0
if (present(B_z)) B_z = 0
if (present(B_phi)) B_phi = 0
return
end if
call pol_flux(R, z, psi_n, dpsidr, dpsidz)
call pol_curr(psi_n, fpol)

View File

@ -12,7 +12,7 @@ contains
use coreprofiles, only : temp, fzeff
use dispersion, only : expinit
use gray_params, only : gray_parameters, gray_data, gray_results, &
print_parameters
print_parameters, EQ_VACUUM
use beams, only : xgygcoeff, launchangles2n
use beamdata, only : pweight, rayi2jk
use gray_errors, only : is_critical, print_err_raytracing, print_err_ecrh_cd
@ -100,12 +100,14 @@ contains
! Initialise the dispersion module
if(params%ecrh_cd%iwarm > 1) call expinit
if (params%equilibrium%iequil /= EQ_VACUUM) then
! Initialise the magsurf_data module
call flux_average ! requires frhotor for dadrhot,dvdrhot
! Initialise the output profiles
call pec_init(params%output%ipec, rhout)
nnd = size(rhop_tab) ! number of radial profile points
end if
call alloc_multipass(nnd, iwait, iroff, iop, iow, yynext, yypnext, yw0, ypw0, stnext, &
stv, p0ray, taus, tau1, etau1, cpls, cpl1, lgcpl1, jphi_beam, &
@ -138,11 +140,13 @@ contains
! print Btot=Bres
! print ne, Te, q, Jphi versus psi, rhop, rhot
if (params%equilibrium%iequil /= EQ_VACUUM) then
call print_bres(bres)
call print_prof(params%profiles)
call print_maps(bres, xgcn, &
norm2(params%antenna%pos(1:2)) * 0.01_wp_, &
sin(params%antenna%beta*degree))
end if
! ========= pre-proc prints END =========
! =========== main loop BEGIN ===========
@ -518,9 +522,11 @@ contains
icd_beam = sum(ccci(:,i))
call vmaxmin(tau0,params%raytracing%nray,taumn,taumx) ! taumn,taumx for print
if (params%equilibrium%iequil /= EQ_VACUUM) then
! compute power and current density profiles for all rays
call spec(psjki,ppabs,ccci,iiv,pabs_beam,icd_beam,dpdv_beam,jphi_beam,jcd_beam, &
pins_beam,currins_beam)
call spec(psjki,ppabs,ccci,iiv,pabs_beam,icd_beam,dpdv_beam, &
jphi_beam,jcd_beam,pins_beam,currins_beam)
end if
pabs_pass(iox) = pabs_pass(iox) + pabs_beam ! 0D results for current pass, sum on O/X mode beams
icd_pass(iox) = icd_pass(iox) + icd_beam
@ -563,17 +569,20 @@ contains
end if
end if
call print_pec(rhop_tab,rhot_tab,jphi_beam,jcd_beam,dpdv_beam,currins_beam, &
pins_beam,ip) ! *print power and current density profiles for current beam
if (params%equilibrium%iequil /= EQ_VACUUM) then
call print_pec(rhop_tab,rhot_tab,jphi_beam,jcd_beam, &
dpdv_beam,currins_beam,pins_beam,ip) ! *print power and current density profiles for current beam
call postproc_profiles(pabs_beam,icd_beam,rhot_tab,dpdv_beam,jphi_beam, &
rhotpav,drhotpav,rhotjava,drhotjava,dpdvp,jphip,rhotp,drhotp,rhotj, &
drhotj,dpdvmx,jphimx,ratjamx,ratjbmx) ! *compute profiles width for current beam
call postproc_profiles(pabs_beam,icd_beam,rhot_tab,dpdv_beam, &
jphi_beam, rhotpav,drhotpav,rhotjava,drhotjava,dpdvp,jphip, &
rhotp,drhotp,rhotj,drhotj,dpdvmx,jphimx,ratjamx,ratjbmx) ! *compute profiles width for current beam
call print_finals(pabs_beam,icd_beam,dpdvp,jphip,rhotpav, &
rhotjava,drhotpav,drhotjava,dpdvmx,jphimx,rhotp,rhotj, &
drhotp,drhotj,ratjamx,ratjbmx,stv(1),psipv(index_rt), &
chipv(index_rt),index_rt,sum(p0ray),cpl_beam1,cpl_beam2) ! *print 0D results for current beam
end if
call print_finals(pabs_beam,icd_beam,dpdvp,jphip,rhotpav,rhotjava, &
drhotpav,drhotjava,dpdvmx,jphimx,rhotp,rhotj,drhotp,drhotj,ratjamx, &
ratjbmx,stv(1),psipv(index_rt),chipv(index_rt),index_rt,sum(p0ray), &
cpl_beam1,cpl_beam2) ! *print 0D results for current beam
! ============ post-proc END ============
end do beam_loop
@ -1872,6 +1881,8 @@ contains
bv(1)=br*csphi-bphi*snphi
bv(2)=br*snphi+bphi*csphi
bv(3)=bz
if (norm2(bv) > 0) then
call pol_limit(anv,bv,bres,sox,ext0(jk),eyt0(jk))
if (jk == 1) then
@ -1880,6 +1891,11 @@ contains
psipol0=psipol0/degree ! convert from rad to degree
chipol0=chipol0/degree
end if
else
! X/O mode are undefined if B=0
psipol0 = 0
chipol0 = 0
end if
else
call stokes_ell(chipol0*degree,psipol0*degree,qq,uu,vv)
if(qq**2 < one) then

View File

@ -242,7 +242,7 @@ contains
! Set global variables
select case (params%equilibrium%iequil)
case (EQ_VACUUM, EQ_ANALYTICAL)
case (EQ_ANALYTICAL)
call set_equil_an
case (EQ_EQDSK_FULL, EQ_EQDSK_PARTIAL)
@ -387,7 +387,7 @@ contains
EQ_EQDSK_FULL, EQ_EQDSK_PARTIAL
use utils, only : range2rect
use limiter, only : limiter_set_globals=>set_globals
use const_and_precisions, only : cm
use const_and_precisions, only : cm, comp_huge
! subroutine arguments
type(gray_parameters), intent(inout) :: params
@ -415,7 +415,11 @@ contains
! Max radius, either due to the plasma extent or equilibrium grid
select case (params%equilibrium%iequil)
case (EQ_VACUUM, EQ_ANALYTICAL)
case (EQ_VACUUM)
! Use a very large R, ~ unbounded
R_max = comp_huge
case (EQ_ANALYTICAL)
! use R+a
block
use equilibrium, only : model

View File

@ -1,4 +0,0 @@
0.0 0.0 0.0 : rr0m,zr0m,rpam ! rhot[m] = min(sqrt((r-rr0m)**2+(z-zr0m)**2), rpam); tor flux phi = pi*b0*rhot**2
0.0 : b0 ! Bphi[T] @ rr0m[m]
1 3 2 : q0, qa, alq ! q = q0 + (qa-q0)*sqrt(psin)**alq
0 : nlim ! number of points in first wall (rlim,zlim) polygon