2015-11-18 17:34:33 +01:00
|
|
|
|
module pec
|
2024-07-07 13:18:55 +02:00
|
|
|
|
use const_and_precisions, only : wp_
|
2024-01-27 12:09:56 +01:00
|
|
|
|
|
2015-11-18 17:34:33 +01:00
|
|
|
|
implicit none
|
2024-01-27 12:09:56 +01:00
|
|
|
|
|
2024-07-07 13:18:55 +02:00
|
|
|
|
real(wp_), dimension(:), allocatable, save :: rhop_tab,rhot_tab
|
2015-11-18 17:34:33 +01:00
|
|
|
|
real(wp_), dimension(:), allocatable, save :: rtabpsi1
|
|
|
|
|
real(wp_), dimension(:), allocatable, save :: dvol,darea
|
|
|
|
|
real(wp_), dimension(:), allocatable, save :: ratjav,ratjbv,ratjplv
|
|
|
|
|
|
|
|
|
|
contains
|
|
|
|
|
|
replace equilibrium module with an object
Similarly to eb648039 this change replaces the `equilibrium` module with
a new `gray_equil` module providing the same functionality without using
global variables.
- `read_eqdsk`, `read_equil_an` are replaced by a single `load_equil`
routine that handles all equilibrium kind (analytical, numerical,
and vacuum).
- `scale_equil` is merged into `load_equil`, which besides reading
the equilibrium from file peforms the rescaling and interpolation based
on the `gray_parameters` settings and the equilibrium kind.
To operate on G-EQDSK data specifically, the `change_cocors` and
`scale_eqdsk` are still available. The numeric equilibrium must then
be initialised manually by calling equil%init().
- `set_equil_spline`, `set_equil_an`, `unset_equil_spline`
are completely removed as the module no longer has any internal state.
- `fq` is replaced by `equil%safety`; `bfield` by `equil%b_field`;
`frhotor`, `frhopol` by `equil%pol2tor` and `equil%pol2tor`;
and the remaining subroutines by other methods of `abstract_equil`
retaining the old name.
- the `contours_psi` subroutine is replaced by `equil%flux_contour`,
with a slightly changed invocation but same functionality.
- the `gray_data` type is no longer required ans has been removed: all
the core subroutines now access the input data only though either
`abstract_equil`, `abstract_plasma` or the `limiter` contour.
2024-08-29 17:16:33 +02:00
|
|
|
|
subroutine pec_init(params, equil, rt_in)
|
2024-07-07 13:18:55 +02:00
|
|
|
|
use gray_params, only : output_parameters, RHO_POL
|
replace equilibrium module with an object
Similarly to eb648039 this change replaces the `equilibrium` module with
a new `gray_equil` module providing the same functionality without using
global variables.
- `read_eqdsk`, `read_equil_an` are replaced by a single `load_equil`
routine that handles all equilibrium kind (analytical, numerical,
and vacuum).
- `scale_equil` is merged into `load_equil`, which besides reading
the equilibrium from file peforms the rescaling and interpolation based
on the `gray_parameters` settings and the equilibrium kind.
To operate on G-EQDSK data specifically, the `change_cocors` and
`scale_eqdsk` are still available. The numeric equilibrium must then
be initialised manually by calling equil%init().
- `set_equil_spline`, `set_equil_an`, `unset_equil_spline`
are completely removed as the module no longer has any internal state.
- `fq` is replaced by `equil%safety`; `bfield` by `equil%b_field`;
`frhotor`, `frhopol` by `equil%pol2tor` and `equil%pol2tor`;
and the remaining subroutines by other methods of `abstract_equil`
retaining the old name.
- the `contours_psi` subroutine is replaced by `equil%flux_contour`,
with a slightly changed invocation but same functionality.
- the `gray_data` type is no longer required ans has been removed: all
the core subroutines now access the input data only though either
`abstract_equil`, `abstract_plasma` or the `limiter` contour.
2024-08-29 17:16:33 +02:00
|
|
|
|
use gray_equil, only : abstract_equil
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
|
type(output_parameters), intent(in) :: params
|
replace equilibrium module with an object
Similarly to eb648039 this change replaces the `equilibrium` module with
a new `gray_equil` module providing the same functionality without using
global variables.
- `read_eqdsk`, `read_equil_an` are replaced by a single `load_equil`
routine that handles all equilibrium kind (analytical, numerical,
and vacuum).
- `scale_equil` is merged into `load_equil`, which besides reading
the equilibrium from file peforms the rescaling and interpolation based
on the `gray_parameters` settings and the equilibrium kind.
To operate on G-EQDSK data specifically, the `change_cocors` and
`scale_eqdsk` are still available. The numeric equilibrium must then
be initialised manually by calling equil%init().
- `set_equil_spline`, `set_equil_an`, `unset_equil_spline`
are completely removed as the module no longer has any internal state.
- `fq` is replaced by `equil%safety`; `bfield` by `equil%b_field`;
`frhotor`, `frhopol` by `equil%pol2tor` and `equil%pol2tor`;
and the remaining subroutines by other methods of `abstract_equil`
retaining the old name.
- the `contours_psi` subroutine is replaced by `equil%flux_contour`,
with a slightly changed invocation but same functionality.
- the `gray_data` type is no longer required ans has been removed: all
the core subroutines now access the input data only though either
`abstract_equil`, `abstract_plasma` or the `limiter` contour.
2024-08-29 17:16:33 +02:00
|
|
|
|
class(abstract_equil), intent(in) :: equil
|
2024-07-07 13:18:55 +02:00
|
|
|
|
real(wp_), intent(in), optional :: rt_in(params%nrho)
|
|
|
|
|
|
|
|
|
|
! local variables
|
2015-11-18 17:34:33 +01:00
|
|
|
|
integer :: it
|
2024-07-07 13:18:55 +02:00
|
|
|
|
real(wp_) :: drt, rt, rt1, rhop1
|
|
|
|
|
real(wp_) :: ratjai, ratjbi, ratjpli
|
|
|
|
|
real(wp_) :: voli0, voli1, areai0, areai1
|
|
|
|
|
|
|
|
|
|
associate(n => params%nrho)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
|
2024-07-07 13:18:55 +02:00
|
|
|
|
! rt_in present: read input grid
|
|
|
|
|
! else: build equidistant grid dimension n
|
2015-11-18 17:34:33 +01:00
|
|
|
|
|
2024-07-07 13:18:55 +02:00
|
|
|
|
! ipec=0 ρ_t grid
|
|
|
|
|
! ipec=1 ρ_p grid
|
2015-11-18 17:34:33 +01:00
|
|
|
|
call dealloc_pec
|
2024-07-07 13:18:55 +02:00
|
|
|
|
allocate(rhop_tab(n), rhot_tab(n), rtabpsi1(0:n))
|
|
|
|
|
allocate(dvol(n), darea(n), ratjav(n), ratjbv(n), ratjplv(n))
|
|
|
|
|
|
|
|
|
|
voli0 = 0
|
|
|
|
|
areai0 = 0
|
|
|
|
|
rtabpsi1(0) = 0
|
|
|
|
|
|
|
|
|
|
do it = 1, n
|
|
|
|
|
if(present(rt_in)) then
|
|
|
|
|
! read radial grid from input
|
2015-11-18 17:34:33 +01:00
|
|
|
|
rt = rt_in(it)
|
2024-07-07 13:18:55 +02:00
|
|
|
|
if (it < n) then
|
2015-11-18 17:34:33 +01:00
|
|
|
|
drt = rt_in(it+1)-rt
|
|
|
|
|
end if
|
|
|
|
|
else
|
2024-07-07 13:18:55 +02:00
|
|
|
|
! build equidistant radial grid
|
|
|
|
|
drt = 1/dble(n - 1)
|
|
|
|
|
rt = (it - 1)*drt
|
2015-11-18 17:34:33 +01:00
|
|
|
|
end if
|
2024-07-07 13:18:55 +02:00
|
|
|
|
! radial coordinate of i-(i+1) interval mid point
|
|
|
|
|
if (it < n) then
|
|
|
|
|
rt1 = rt + drt/2
|
2015-11-18 17:34:33 +01:00
|
|
|
|
else
|
2024-07-07 13:18:55 +02:00
|
|
|
|
rt1 = 1
|
2015-11-18 17:34:33 +01:00
|
|
|
|
end if
|
2024-07-07 13:18:55 +02:00
|
|
|
|
if (params%ipec == RHO_POL) then
|
2015-11-18 17:34:33 +01:00
|
|
|
|
rhop_tab(it) = rt
|
replace equilibrium module with an object
Similarly to eb648039 this change replaces the `equilibrium` module with
a new `gray_equil` module providing the same functionality without using
global variables.
- `read_eqdsk`, `read_equil_an` are replaced by a single `load_equil`
routine that handles all equilibrium kind (analytical, numerical,
and vacuum).
- `scale_equil` is merged into `load_equil`, which besides reading
the equilibrium from file peforms the rescaling and interpolation based
on the `gray_parameters` settings and the equilibrium kind.
To operate on G-EQDSK data specifically, the `change_cocors` and
`scale_eqdsk` are still available. The numeric equilibrium must then
be initialised manually by calling equil%init().
- `set_equil_spline`, `set_equil_an`, `unset_equil_spline`
are completely removed as the module no longer has any internal state.
- `fq` is replaced by `equil%safety`; `bfield` by `equil%b_field`;
`frhotor`, `frhopol` by `equil%pol2tor` and `equil%pol2tor`;
and the remaining subroutines by other methods of `abstract_equil`
retaining the old name.
- the `contours_psi` subroutine is replaced by `equil%flux_contour`,
with a slightly changed invocation but same functionality.
- the `gray_data` type is no longer required ans has been removed: all
the core subroutines now access the input data only though either
`abstract_equil`, `abstract_plasma` or the `limiter` contour.
2024-08-29 17:16:33 +02:00
|
|
|
|
rhot_tab(it) = equil%pol2tor(rt)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
rhop1 = rt1
|
|
|
|
|
else
|
|
|
|
|
rhot_tab(it) = rt
|
replace equilibrium module with an object
Similarly to eb648039 this change replaces the `equilibrium` module with
a new `gray_equil` module providing the same functionality without using
global variables.
- `read_eqdsk`, `read_equil_an` are replaced by a single `load_equil`
routine that handles all equilibrium kind (analytical, numerical,
and vacuum).
- `scale_equil` is merged into `load_equil`, which besides reading
the equilibrium from file peforms the rescaling and interpolation based
on the `gray_parameters` settings and the equilibrium kind.
To operate on G-EQDSK data specifically, the `change_cocors` and
`scale_eqdsk` are still available. The numeric equilibrium must then
be initialised manually by calling equil%init().
- `set_equil_spline`, `set_equil_an`, `unset_equil_spline`
are completely removed as the module no longer has any internal state.
- `fq` is replaced by `equil%safety`; `bfield` by `equil%b_field`;
`frhotor`, `frhopol` by `equil%pol2tor` and `equil%pol2tor`;
and the remaining subroutines by other methods of `abstract_equil`
retaining the old name.
- the `contours_psi` subroutine is replaced by `equil%flux_contour`,
with a slightly changed invocation but same functionality.
- the `gray_data` type is no longer required ans has been removed: all
the core subroutines now access the input data only though either
`abstract_equil`, `abstract_plasma` or the `limiter` contour.
2024-08-29 17:16:33 +02:00
|
|
|
|
rhop_tab(it) = equil%tor2pol(rt)
|
|
|
|
|
rhop1 = equil%tor2pol(rt1)
|
2024-07-07 13:18:55 +02:00
|
|
|
|
end if
|
|
|
|
|
|
|
|
|
|
! psi grid at mid points, size n+1, for use in pec_tab
|
2015-11-18 17:34:33 +01:00
|
|
|
|
rtabpsi1(it) = rhop1**2
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
2024-10-23 15:15:37 +02:00
|
|
|
|
call equil%flux_average(rhop1, area=areai1, volume=voli1)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
dvol(it) = abs(voli1 - voli0)
|
|
|
|
|
darea(it) = abs(areai1 - areai0)
|
|
|
|
|
voli0 = voli1
|
|
|
|
|
areai0 = areai1
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
2024-10-23 15:15:37 +02:00
|
|
|
|
call equil%flux_average(rhop_tab(it), ratio_astra_tor=ratjai, &
|
|
|
|
|
ratio_jintrac_tor=ratjbi, ratio_paral_tor=ratjpli)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
ratjav(it) = ratjai
|
|
|
|
|
ratjbv(it) = ratjbi
|
|
|
|
|
ratjplv(it) = ratjpli
|
|
|
|
|
end do
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
|
|
|
|
end associate
|
2015-11-18 17:34:33 +01:00
|
|
|
|
end subroutine pec_init
|
|
|
|
|
|
|
|
|
|
|
2024-07-07 13:18:55 +02:00
|
|
|
|
subroutine spec(psjki, ppabs, ccci, iiv, pabs, &
|
|
|
|
|
currt, dpdv, ajphiv, ajcd, pins, currins)
|
|
|
|
|
! subroutine arguments
|
|
|
|
|
real(wp_), dimension(:, :), intent(in) :: psjki, ppabs, ccci
|
|
|
|
|
integer, intent(in) :: iiv(:)
|
|
|
|
|
real(wp_), intent(in) :: pabs, currt
|
|
|
|
|
real(wp_), dimension(:), intent(out) :: dpdv, ajphiv, ajcd, pins, currins
|
|
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
|
integer :: i, ii, jk
|
|
|
|
|
real(wp_) :: spds, sccs, facpds, facjs
|
|
|
|
|
real(wp_), dimension(size(psjki, dim=2)) :: xxi, ypt, yamp
|
|
|
|
|
real(wp_), dimension(size(dpdv)) :: wdpdv, wajphiv
|
|
|
|
|
|
|
|
|
|
! calculation of dP and dI over radial grid
|
|
|
|
|
dpdv = 0
|
|
|
|
|
ajphiv = 0
|
|
|
|
|
do jk = 1, size(iiv)
|
|
|
|
|
ii = iiv(jk)
|
|
|
|
|
if (ii < size(psjki, dim=2)) then
|
|
|
|
|
if (psjki(jk,ii+1) /= 0) ii = ii + 1 !!! CHECK
|
2015-11-18 17:34:33 +01:00
|
|
|
|
end if
|
2024-07-07 13:18:55 +02:00
|
|
|
|
xxi = 0
|
|
|
|
|
ypt = 0
|
|
|
|
|
yamp = 0
|
|
|
|
|
do i = 1, ii
|
2015-11-18 17:34:33 +01:00
|
|
|
|
xxi(i)=abs(psjki(jk,i))
|
2024-07-07 13:18:55 +02:00
|
|
|
|
if(xxi(i) <= 1) then
|
|
|
|
|
ypt(i) = ppabs(jk,i)
|
|
|
|
|
yamp(i) = ccci(jk,i)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
end if
|
|
|
|
|
end do
|
2024-07-07 13:18:55 +02:00
|
|
|
|
call pec_tab(xxi, ypt, yamp, ii, rtabpsi1, wdpdv, wajphiv)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
dpdv = dpdv + wdpdv
|
|
|
|
|
ajphiv = ajphiv + wajphiv
|
|
|
|
|
end do
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
|
|
|
|
! here dpdv is still dP (not divided yet by dV)
|
|
|
|
|
! here ajphiv is still dI (not divided yet by dA)
|
|
|
|
|
spds = 0
|
|
|
|
|
sccs = 0
|
|
|
|
|
do i = 1, size(dpdv)
|
|
|
|
|
spds = spds + dpdv(i)
|
|
|
|
|
sccs = sccs + ajphiv(i)
|
|
|
|
|
pins(i) = spds
|
|
|
|
|
currins(i) = sccs
|
2015-11-18 17:34:33 +01:00
|
|
|
|
end do
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
|
|
|
|
facpds = 1
|
|
|
|
|
facjs = 1
|
|
|
|
|
if(spds > 0) facpds = pabs / spds
|
|
|
|
|
if(sccs /= 0) facjs = currt / sccs
|
|
|
|
|
|
|
|
|
|
dpdv = facpds * (dpdv / dvol)
|
|
|
|
|
ajphiv = facjs * (ajphiv / darea)
|
|
|
|
|
ajcd = ajphiv * ratjbv
|
|
|
|
|
pins = facpds * pins
|
|
|
|
|
currins = facjs * currins
|
|
|
|
|
|
|
|
|
|
! now dpdv is dP/dV [MW/m^3]
|
|
|
|
|
! now ajphiv is J_phi=dI/dA [MA/m^2]
|
|
|
|
|
end subroutine spec
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
subroutine pec_tab(xxi, ypt, yamp, ii, xtab1, wdpdv, wajphiv)
|
|
|
|
|
! Power and current projected on ψ grid - mid points
|
2024-09-11 11:52:42 +02:00
|
|
|
|
use utils, only : locate, linear_interp
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
|
real(wp_), dimension(:), intent(in) :: xxi, ypt, yamp
|
|
|
|
|
integer, intent(in) :: ii
|
|
|
|
|
real(wp_), dimension(0:), intent(in) :: xtab1
|
|
|
|
|
real(wp_), dimension(ubound(xtab1, 1)), intent(out) :: wdpdv, wajphiv
|
|
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
|
integer :: isev(21)
|
|
|
|
|
real(wp_) :: ppa1, ppa2, cci1, cci2, dppa, didst, rt1
|
|
|
|
|
integer :: i, is, ise0, idecr, iise0, iise, iis, iis1
|
|
|
|
|
integer :: ind1, ind2, iind, ind, indi, itb1
|
|
|
|
|
|
2015-11-18 17:34:33 +01:00
|
|
|
|
isev = 0
|
|
|
|
|
ise0 = 0
|
|
|
|
|
idecr = -1
|
|
|
|
|
is = 1
|
2024-07-07 13:18:55 +02:00
|
|
|
|
wdpdv = 0
|
|
|
|
|
wajphiv = 0
|
2015-11-18 17:34:33 +01:00
|
|
|
|
do i=1,ii
|
|
|
|
|
if(ise0 == 0) then
|
2024-07-07 13:18:55 +02:00
|
|
|
|
if(xxi(i) < 1) then
|
2015-11-18 17:34:33 +01:00
|
|
|
|
ise0 = i
|
2019-03-26 15:21:22 +01:00
|
|
|
|
isev(is) = max(i - 1, 1)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
is = is + 1
|
|
|
|
|
end if
|
|
|
|
|
else
|
|
|
|
|
if (idecr == -1) then
|
|
|
|
|
if(xxi(i) > xxi(i-1)) then
|
|
|
|
|
isev(is) = i - 1
|
|
|
|
|
is = is + 1
|
|
|
|
|
idecr = 1
|
|
|
|
|
end if
|
|
|
|
|
else
|
2024-07-07 13:18:55 +02:00
|
|
|
|
if(xxi(i) > 1) exit
|
2015-11-18 17:34:33 +01:00
|
|
|
|
if(xxi(i) < xxi(i-1)) then
|
|
|
|
|
isev(is) = i - 1
|
|
|
|
|
is = is + 1
|
|
|
|
|
idecr = -1
|
|
|
|
|
end if
|
|
|
|
|
end if
|
|
|
|
|
end if
|
|
|
|
|
end do
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
2015-11-18 17:34:33 +01:00
|
|
|
|
isev(is) = i-1
|
2024-07-07 13:18:55 +02:00
|
|
|
|
ppa1 = 0
|
|
|
|
|
cci1 = 0
|
|
|
|
|
|
|
|
|
|
do iis = 1, is-1
|
2015-11-18 17:34:33 +01:00
|
|
|
|
iis1 = iis + 1
|
|
|
|
|
iise0 = isev(iis)
|
|
|
|
|
iise = isev(iis1)
|
|
|
|
|
if (mod(iis,2) /= 0) then
|
|
|
|
|
idecr = -1
|
2024-07-07 13:18:55 +02:00
|
|
|
|
ind1 = size(wdpdv)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
ind2 = 2
|
|
|
|
|
iind = -1
|
|
|
|
|
else
|
|
|
|
|
idecr = 1
|
|
|
|
|
ind1 = 1
|
2024-07-07 13:18:55 +02:00
|
|
|
|
ind2 = size(wdpdv)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
iind = 1
|
|
|
|
|
end if
|
2024-07-07 13:18:55 +02:00
|
|
|
|
do ind = ind1, ind2, iind
|
|
|
|
|
indi = ind
|
2015-11-18 17:34:33 +01:00
|
|
|
|
if (idecr == -1) indi = ind - 1
|
|
|
|
|
rt1 = xtab1(indi)
|
2024-09-11 11:52:42 +02:00
|
|
|
|
call locate(xxi(iise0:iise), rt1, itb1)
|
|
|
|
|
if (itb1 >= 1 .and. itb1 <= iise) then
|
|
|
|
|
itb1 = itb1 + (iise0 - 1)
|
2019-03-26 15:21:22 +01:00
|
|
|
|
if(itb1 == iise) then
|
|
|
|
|
ppa2=ypt(itb1)
|
|
|
|
|
cci2=yamp(itb1)
|
|
|
|
|
else
|
2024-09-11 11:52:42 +02:00
|
|
|
|
ppa2 = linear_interp(xxi(itb1), ypt(itb1), xxi(itb1+1), ypt(itb1+1), rt1)
|
|
|
|
|
cci2 = linear_interp(xxi(itb1), yamp(itb1), xxi(itb1+1), yamp(itb1+1), rt1)
|
2019-03-26 15:21:22 +01:00
|
|
|
|
end if
|
2015-11-18 17:34:33 +01:00
|
|
|
|
dppa = ppa2 - ppa1
|
|
|
|
|
didst = cci2 - cci1
|
|
|
|
|
wdpdv(ind) = wdpdv(ind) + dppa
|
|
|
|
|
wajphiv(ind) = wajphiv(ind) + didst
|
|
|
|
|
ppa1 = ppa2
|
|
|
|
|
cci1 = cci2
|
|
|
|
|
end if
|
|
|
|
|
end do
|
|
|
|
|
end do
|
|
|
|
|
end subroutine pec_tab
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
|
|
|
|
|
replace equilibrium module with an object
Similarly to eb648039 this change replaces the `equilibrium` module with
a new `gray_equil` module providing the same functionality without using
global variables.
- `read_eqdsk`, `read_equil_an` are replaced by a single `load_equil`
routine that handles all equilibrium kind (analytical, numerical,
and vacuum).
- `scale_equil` is merged into `load_equil`, which besides reading
the equilibrium from file peforms the rescaling and interpolation based
on the `gray_parameters` settings and the equilibrium kind.
To operate on G-EQDSK data specifically, the `change_cocors` and
`scale_eqdsk` are still available. The numeric equilibrium must then
be initialised manually by calling equil%init().
- `set_equil_spline`, `set_equil_an`, `unset_equil_spline`
are completely removed as the module no longer has any internal state.
- `fq` is replaced by `equil%safety`; `bfield` by `equil%b_field`;
`frhotor`, `frhopol` by `equil%pol2tor` and `equil%pol2tor`;
and the remaining subroutines by other methods of `abstract_equil`
retaining the old name.
- the `contours_psi` subroutine is replaced by `equil%flux_contour`,
with a slightly changed invocation but same functionality.
- the `gray_data` type is no longer required ans has been removed: all
the core subroutines now access the input data only though either
`abstract_equil`, `abstract_plasma` or the `limiter` contour.
2024-08-29 17:16:33 +02:00
|
|
|
|
subroutine postproc_profiles(equil, pabs, currt, rhot_tab, dpdv, ajphiv, &
|
2015-11-23 18:55:27 +01:00
|
|
|
|
rhotpav,drhotpav,rhotjava,drhotjava,dpdvp, ajphip, &
|
|
|
|
|
rhotp, drhotp, rhotjfi, drhotjfi, dpdvmx,ajmxfi, ratjamx, ratjbmx)
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
|
|
|
|
! Radial average values over power and current density profile
|
2019-03-26 15:21:22 +01:00
|
|
|
|
use const_and_precisions, only : pi, comp_eps
|
replace equilibrium module with an object
Similarly to eb648039 this change replaces the `equilibrium` module with
a new `gray_equil` module providing the same functionality without using
global variables.
- `read_eqdsk`, `read_equil_an` are replaced by a single `load_equil`
routine that handles all equilibrium kind (analytical, numerical,
and vacuum).
- `scale_equil` is merged into `load_equil`, which besides reading
the equilibrium from file peforms the rescaling and interpolation based
on the `gray_parameters` settings and the equilibrium kind.
To operate on G-EQDSK data specifically, the `change_cocors` and
`scale_eqdsk` are still available. The numeric equilibrium must then
be initialised manually by calling equil%init().
- `set_equil_spline`, `set_equil_an`, `unset_equil_spline`
are completely removed as the module no longer has any internal state.
- `fq` is replaced by `equil%safety`; `bfield` by `equil%b_field`;
`frhotor`, `frhopol` by `equil%pol2tor` and `equil%pol2tor`;
and the remaining subroutines by other methods of `abstract_equil`
retaining the old name.
- the `contours_psi` subroutine is replaced by `equil%flux_contour`,
with a slightly changed invocation but same functionality.
- the `gray_data` type is no longer required ans has been removed: all
the core subroutines now access the input data only though either
`abstract_equil`, `abstract_plasma` or the `limiter` contour.
2024-08-29 17:16:33 +02:00
|
|
|
|
use gray_equil, only : abstract_equil
|
|
|
|
|
|
|
|
|
|
class(abstract_equil), intent(in) :: equil
|
|
|
|
|
real(wp_), intent(in) :: pabs,currt
|
|
|
|
|
real(wp_), intent(in) :: rhot_tab(:)
|
|
|
|
|
real(wp_), intent(in) :: dpdv(:), ajphiv(:)
|
|
|
|
|
real(wp_), intent(out) :: rhotpav, drhotpav, dpdvp
|
|
|
|
|
real(wp_), intent(out) :: rhotjava, drhotjava, ajphip
|
|
|
|
|
real(wp_), intent(out) :: rhotp, drhotp, dpdvmx
|
|
|
|
|
real(wp_), intent(out) :: rhotjfi, drhotjfi, ajmxfi
|
|
|
|
|
real(wp_), intent(out) :: ratjamx, ratjbmx
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
|
|
|
|
real(wp_) :: sccsa, ratjplmx, rhopjava, rhoppav
|
|
|
|
|
real(wp_) :: rhotjav, rhot2pav, rhot2java, dvdrhotav, dadrhotava
|
|
|
|
|
|
|
|
|
|
rhotpav = 0
|
|
|
|
|
rhot2pav = 0
|
|
|
|
|
rhotjav = 0
|
|
|
|
|
rhotjava = 0
|
|
|
|
|
rhot2java = 0
|
|
|
|
|
|
|
|
|
|
if (pabs > 0) then
|
2015-11-18 17:34:33 +01:00
|
|
|
|
rhotpav = sum(rhot_tab *dpdv*dvol)/pabs
|
|
|
|
|
rhot2pav = sum(rhot_tab**2*dpdv*dvol)/pabs
|
|
|
|
|
end if
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
|
|
|
|
if (abs(currt) > 0) then
|
2015-11-18 17:34:33 +01:00
|
|
|
|
rhotjav = sum(rhot_tab*ajphiv*darea)/currt
|
|
|
|
|
end if
|
|
|
|
|
sccsa = sum(abs(ajphiv)*darea)
|
2024-07-07 13:18:55 +02:00
|
|
|
|
if (sccsa > 0) then
|
2015-11-18 17:34:33 +01:00
|
|
|
|
rhotjava = sum(rhot_tab *abs(ajphiv)*darea)/sccsa
|
|
|
|
|
rhot2java = sum(rhot_tab**2*abs(ajphiv)*darea)/sccsa
|
|
|
|
|
end if
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
|
|
|
|
! factor sqrt(8) = 2sqrt(2) to match full width of gaussian profile
|
2019-03-26 15:21:22 +01:00
|
|
|
|
drhotpav = sqrt(8._wp_*(max(rhot2pav -rhotpav**2,comp_eps)))
|
|
|
|
|
drhotjava = sqrt(8._wp_*(max(rhot2java-rhotjava**2,comp_eps)))
|
2015-11-18 17:34:33 +01:00
|
|
|
|
|
replace equilibrium module with an object
Similarly to eb648039 this change replaces the `equilibrium` module with
a new `gray_equil` module providing the same functionality without using
global variables.
- `read_eqdsk`, `read_equil_an` are replaced by a single `load_equil`
routine that handles all equilibrium kind (analytical, numerical,
and vacuum).
- `scale_equil` is merged into `load_equil`, which besides reading
the equilibrium from file peforms the rescaling and interpolation based
on the `gray_parameters` settings and the equilibrium kind.
To operate on G-EQDSK data specifically, the `change_cocors` and
`scale_eqdsk` are still available. The numeric equilibrium must then
be initialised manually by calling equil%init().
- `set_equil_spline`, `set_equil_an`, `unset_equil_spline`
are completely removed as the module no longer has any internal state.
- `fq` is replaced by `equil%safety`; `bfield` by `equil%b_field`;
`frhotor`, `frhopol` by `equil%pol2tor` and `equil%pol2tor`;
and the remaining subroutines by other methods of `abstract_equil`
retaining the old name.
- the `contours_psi` subroutine is replaced by `equil%flux_contour`,
with a slightly changed invocation but same functionality.
- the `gray_data` type is no longer required ans has been removed: all
the core subroutines now access the input data only though either
`abstract_equil`, `abstract_plasma` or the `limiter` contour.
2024-08-29 17:16:33 +02:00
|
|
|
|
rhoppav = equil%tor2pol(rhotpav)
|
|
|
|
|
rhopjava = equil%tor2pol(rhotjava)
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
|
|
|
|
if (pabs > 0) then
|
2024-10-23 15:15:37 +02:00
|
|
|
|
call equil%flux_average(rhoppav, dVdrho_t=dvdrhotav)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
dpdvp = pabs*2.0_wp_/(sqrt(pi)*drhotpav*dvdrhotav)
|
2024-07-07 13:18:55 +02:00
|
|
|
|
call profwidth(rhot_tab,dpdv,rhotp,dpdvmx,drhotp)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
else
|
2024-07-07 13:18:55 +02:00
|
|
|
|
dpdvp = 0
|
|
|
|
|
rhotp = 0
|
|
|
|
|
dpdvmx = 0
|
|
|
|
|
drhotp = 0
|
2015-11-18 17:34:33 +01:00
|
|
|
|
end if
|
|
|
|
|
|
2024-07-07 13:18:55 +02:00
|
|
|
|
if (sccsa > 0) then
|
2024-10-23 15:15:37 +02:00
|
|
|
|
call equil%flux_average(rhopjava, dAdrho_t=dadrhotava, ratio_astra_tor=ratjamx, &
|
|
|
|
|
ratio_jintrac_tor=ratjbmx, ratio_paral_tor=ratjplmx)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
ajphip = currt*2.0_wp_/(sqrt(pi)*drhotjava*dadrhotava)
|
2024-07-07 13:18:55 +02:00
|
|
|
|
call profwidth(rhot_tab,ajphiv,rhotjfi,ajmxfi,drhotjfi)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
else
|
2024-07-07 13:18:55 +02:00
|
|
|
|
ajphip = 0
|
|
|
|
|
rhotjfi = 0
|
|
|
|
|
ajmxfi = 0
|
|
|
|
|
drhotjfi = 0
|
|
|
|
|
ratjamx = 0
|
|
|
|
|
ratjbmx = 0
|
2015-11-18 17:34:33 +01:00
|
|
|
|
end if
|
2024-07-07 13:18:55 +02:00
|
|
|
|
end subroutine postproc_profiles
|
2015-11-18 17:34:33 +01:00
|
|
|
|
|
|
|
|
|
|
2024-07-07 13:18:55 +02:00
|
|
|
|
subroutine profwidth(xx,yy,xpk,ypk,dxxe)
|
|
|
|
|
use const_and_precisions, only : emn1
|
2024-09-11 11:52:42 +02:00
|
|
|
|
use utils, only : locate, linear_interp, vmaxmin
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
|
real(wp_), intent(in) :: xx(:), yy(:)
|
|
|
|
|
real(wp_), intent(out) :: xpk, ypk, dxxe
|
|
|
|
|
|
|
|
|
|
! local variables
|
2015-11-18 17:34:33 +01:00
|
|
|
|
integer :: imn,imx,ipk,ie
|
|
|
|
|
real(wp_) :: xmn,xmx,ymn,ymx,xpkp,xpkm,yye,rte1,rte2
|
|
|
|
|
real(wp_) :: ypkp,ypkm
|
|
|
|
|
|
2024-09-11 11:52:42 +02:00
|
|
|
|
call vmaxmin(yy, ymn, ymx, imn, imx)
|
2024-07-07 13:18:55 +02:00
|
|
|
|
ypk = 0
|
2015-11-18 17:34:33 +01:00
|
|
|
|
xmx = xx(imx)
|
|
|
|
|
xmn = xx(imn)
|
|
|
|
|
if (abs(ymx) > abs(ymn)) then
|
|
|
|
|
ipk = imx
|
|
|
|
|
ypkp = ymx
|
|
|
|
|
xpkp = xmx
|
|
|
|
|
if(abs(ymn/ymx) < 1.0e-2_wp_) ymn = 0.0_wp_
|
|
|
|
|
ypkm = ymn
|
|
|
|
|
xpkm = xmn
|
|
|
|
|
else
|
|
|
|
|
ipk = imn
|
|
|
|
|
ypkp = ymn
|
|
|
|
|
xpkp = xmn
|
|
|
|
|
if(abs(ymx/ymn) < 1.0e-2_wp_) ymx = 0.0_wp_
|
|
|
|
|
ypkm = ymx
|
|
|
|
|
xpkm = xmx
|
|
|
|
|
end if
|
2024-07-07 13:18:55 +02:00
|
|
|
|
if(xpkp > 0) then
|
2015-11-18 17:34:33 +01:00
|
|
|
|
xpk = xpkp
|
|
|
|
|
ypk = ypkp
|
|
|
|
|
yye = ypk*emn1
|
2024-09-11 11:52:42 +02:00
|
|
|
|
call locate(yy(1:ipk), yye, ie)
|
2024-07-07 13:18:55 +02:00
|
|
|
|
if(ie > 0 .and. ie < size(yy)) then
|
2024-09-11 11:52:42 +02:00
|
|
|
|
rte1 = linear_interp(yy(ie), xx(ie), yy(ie+1), xx(ie+1), yye)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
else
|
2024-07-07 13:18:55 +02:00
|
|
|
|
rte1 = 0
|
2015-11-18 17:34:33 +01:00
|
|
|
|
end if
|
2024-09-11 11:52:42 +02:00
|
|
|
|
call locate(yy(ipk:size(yy)), yye, ie)
|
2024-07-07 13:18:55 +02:00
|
|
|
|
if(ie > 0 .and. ie < size(yy)) then
|
2024-09-11 11:52:42 +02:00
|
|
|
|
ie = ie + (ipk - 1)
|
|
|
|
|
rte2 = linear_interp(yy(ie), xx(ie), yy(ie+1), xx(ie+1), yye)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
else
|
2024-07-07 13:18:55 +02:00
|
|
|
|
rte2 = 0
|
2015-11-18 17:34:33 +01:00
|
|
|
|
end if
|
|
|
|
|
else
|
|
|
|
|
ipk=2
|
|
|
|
|
xpk=xx(2)
|
|
|
|
|
ypk=yy(2)
|
|
|
|
|
rte1=0.0_wp_
|
|
|
|
|
yye=ypk*emn1
|
2024-09-11 11:52:42 +02:00
|
|
|
|
call locate(yy, yye, ie)
|
2024-07-07 13:18:55 +02:00
|
|
|
|
if(ie > 0 .and. ie < size(yy)) then
|
2024-09-11 11:52:42 +02:00
|
|
|
|
rte2 = linear_interp(yy(ie), xx(ie), yy(ie+1), xx(ie+1), yye)
|
2015-11-18 17:34:33 +01:00
|
|
|
|
else
|
2024-07-07 13:18:55 +02:00
|
|
|
|
rte2 = 0
|
2015-11-18 17:34:33 +01:00
|
|
|
|
end if
|
|
|
|
|
end if
|
|
|
|
|
dxxe = rte2 - rte1
|
2024-07-07 13:18:55 +02:00
|
|
|
|
if(ymx /= 0 .and. ymn /= 0) dxxe = -dxxe
|
2015-11-18 17:34:33 +01:00
|
|
|
|
end subroutine profwidth
|
|
|
|
|
|
2024-07-07 13:18:55 +02:00
|
|
|
|
|
2015-11-18 17:34:33 +01:00
|
|
|
|
subroutine dealloc_pec
|
|
|
|
|
if (allocated(rhop_tab)) deallocate(rhop_tab)
|
|
|
|
|
if (allocated(rhot_tab)) deallocate(rhot_tab)
|
|
|
|
|
if (allocated(rtabpsi1)) deallocate(rtabpsi1)
|
|
|
|
|
if (allocated(dvol)) deallocate(dvol)
|
|
|
|
|
if (allocated(darea)) deallocate(darea)
|
|
|
|
|
if (allocated(ratjav)) deallocate(ratjav)
|
|
|
|
|
if (allocated(ratjbv)) deallocate(ratjbv)
|
|
|
|
|
if (allocated(ratjplv)) deallocate(ratjplv)
|
|
|
|
|
end subroutine dealloc_pec
|
|
|
|
|
|
|
|
|
|
end module pec
|