gray/src/gray_plasma.f90
Michele Guerini Rocco 166086d369
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-11-03 09:18:33 +01:00

530 lines
18 KiB
Fortran
Raw Blame History

This file contains ambiguous Unicode characters

This file contains Unicode characters that might be confused with other characters. If you think that this is intentional, you can safely ignore this warning. Use the Escape button to reveal them.

! This module handles the loading, interpolation and evaluation of the
! plasma profiles (density, temperature, effective charge)
!
! The type `abstract_plasma` provides an abstract interface for
! accessing the profiles data; these implementations are available:
! - `numeric_plasma`, numerical data interpolated using splines
! - `analytic_plasma`, simple analytic model
!
module gray_plasma
use const_and_precisions, only : wp_, zero
use splines, only : spline_simple, spline_1d
implicit none
! Generic plasma interface
type, abstract :: abstract_plasma
contains
procedure(density_sub), deferred :: density
procedure(temp_fun), deferred :: temp
procedure(zeff_fun), deferred :: zeff
end type
abstract interface
subroutine density_sub(self, psin, dens, ddens)
! Computes the density its first derivative as a function of
! normalised poloidal flux.
!
! Note: density has units of 10¹⁹ m⁻³.
import :: abstract_plasma, wp_
class(abstract_plasma), intent(in) :: self
real(wp_), intent(in) :: psin ! normalised poloidal flux
real(wp_), intent(out) :: dens, ddens ! density and first derivative
end subroutine density_sub
function temp_fun(self, psin) result(temp)
! Computes the temperature as a function of the
! normalised poloidal flux.
!
! Note: temperature has units of keV.
import :: abstract_plasma, wp_
class(abstract_plasma), intent(in) :: self
real(wp_), intent(in) :: psin
real(wp_) :: temp
end function temp_fun
function zeff_fun(self, psin) result(zeff)
! Computes the effective charge Z_eff as a
! function of the normalised poloidal flux.
import :: abstract_plasma, wp_
class(abstract_plasma), intent(in) :: self
real(wp_), intent(in) :: psin
real(wp_) :: zeff
end function zeff_fun
end interface
! Analytical plasma description
type, extends(abstract_plasma) :: analytic_plasma
private
real(wp_) :: dens0 ! Density scaling factor
real(wp_) :: n1, n2 ! Density exponents
real(wp_) :: te0, te1 ! Temperature at ψ=0, ψ=1
real(wp_) :: t1, t2 ! Temperature exponents
real(wp_) :: z ! Effective charge
contains
procedure :: density => analytic_density
procedure :: temp => analytic_temp
procedure :: zeff => analytic_zeff
end type
! Parameters of the C² polynomial tail of the density spline
type density_tail
real(wp_) :: start ! ψ₀, start of the tail
real(wp_) :: end ! ψ₁, end of the end
real(wp_) :: value ! s(ψ₀), value at the start
real(wp_) :: deriv1 ! s'(ψ₀), first derivative at the start
real(wp_) :: deriv2 ! s"(ψ₀), second derivative at the start
end type
! Numerical plasma description
type, extends(abstract_plasma) :: numeric_plasma
private
type(spline_1d) :: dens_spline
type(spline_simple) :: temp_spline
type(spline_simple) :: zeff_spline
type(density_tail) :: tail
contains
procedure :: init => init_splines
procedure :: density => numeric_density
procedure :: temp => numeric_temp
procedure :: zeff => numeric_zeff
end type
private
public abstract_plasma ! The abstract plasma object
public analytic_plasma, numeric_plasma ! Implementations
public load_plasma ! To load plasma from file
contains
subroutine analytic_density(self, psin, dens, ddens)
! subroutine arguments
class(analytic_plasma), intent(in) :: self
real(wp_), intent(in) :: psin ! normalised poloidal flux
real(wp_), intent(out) :: dens, ddens ! density and first derivative
if (psin > 1) then
! Outside: the density and every derivative is identically zero
dens = 0
ddens = 0
else
! Inside: the density is given by n(ψ) = dens0⋅(1 - ψ^n1)^n2
dens = self%dens0 * (1 - psin**self%n1)**self%n2
ddens = -self%dens0 * self%n1*self%n2 * psin**(self%n1 - 1) &
* (1 - psin**self%n1)**(self%n2 - 1)
end if
end subroutine analytic_density
function analytic_temp(self, psin) result(temp)
! function arguments
class(analytic_plasma), intent(in) :: self
real(wp_), intent(in) :: psin
real(wp_) :: temp
! local variables
real(wp_) :: proft
if (psin >= 1 .or. psin < 0) then
! Outside: temperature is undefined, return 0
temp = 0
else
! Inside: the temperature is given by
!
! T(ψ) = (te0 - te1)⋅(1 - ψ^t1)^t2 + te1
!
proft = (1 - psin**self%t1)**self%t2
temp = (self%te0 - self%te1)*proft + self%te1
end if
end function analytic_temp
function analytic_zeff(self, psin) result(zeff)
! function arguments
class(analytic_plasma), intent(in) :: self
real(wp_), intent(in) :: psin
real(wp_) :: zeff
if (psin >= 1 .or. psin < 0) then
! Outside: Zeff is undefined, return 0
zeff = 0
else
! Inside: just a constant
zeff = self%z
end if
end function analytic_zeff
subroutine numeric_density(self, psin, dens, ddens)
use logger, only : log_error
! subroutine arguments
class(numeric_plasma), intent(in) :: self
real(wp_), intent(in) :: psin
real(wp_), intent(out) :: dens, ddens
! local variables
character(256) :: msg ! for log messages formatting
! Initialise both to zero
dens = 0
ddens = 0
! Bail out when ψ is not available
if (psin < 0) return
! Past the tail end both density and its
! derivatives are identically zero
if (psin >= self%tail%end) return
if (psin < self%tail%start) then
! Use the interpolating spline when in range
! Evaluate the spline
dens = self%dens_spline%eval(psin)
ddens = self%dens_spline%deriv(psin)
! Evaluate the spline 1st derivative
if (abs(dens) < 1.0e-10_wp_) dens = zero
else
! Use a C² polynomial extension outside (ψ > ψ₀)
! The tail consists of the product p(ψ)⋅h(t), where:
!
! - p(ψ) is the 2nd order Taylor polynomial of the spline,
! centered at ψ₀. See set_profiles_spline for details.
!
! - h(t) is a "smoothing" polynomial in the variable
! t = (ψ - ψ₀)/(ψ₁ - ψ₀), defined as:
!
! h(t) = (1 - t)³⋅(1 + 3t + 6t²)
!
! with the following properties:
!
! h(0) = 1 h'(0)=0 h"(0)=0
! h(1) = 0 h'(1)=0 h"(1)=0
block
real(wp_) :: dpsi, t, p, dp, h, dh
dpsi = psin - self%tail%start ! Δψ = (ψ - ψ₀)
! Taylor polynomial p(ψ) and its derivative
p = self%tail%value + dpsi*self%tail%deriv1 + dpsi**2*self%tail%deriv2/2
dp = self%tail%deriv1 + dpsi*self%tail%deriv2
! Smoothing polynomial h(t) and its derivative
t = dpsi/(self%tail%end - self%tail%start)
h = (1 - t)**3 * (1 + 3*t + 6*t**2)
dh = -30*(1 - t)**2 * t**2 / (self%tail%end - self%tail%start)
dens = p*h
ddens = dp*h + p*dh
end block
end if
if (dens < 0) then
write (msg, '("negative density:", 2(x,a,"=",g0.3))') &
'ne', dens, 'ψ', psin
call log_error(msg, mod='coreprofiles', proc='density')
end if
end subroutine numeric_density
function numeric_temp(self, psin) result(temp)
! function arguments
class(numeric_plasma), intent(in) :: self
real(wp_), intent(in) :: psin
real(wp_) :: temp
temp = 0
! Outside: temperature is undefined, return 0
if (psin >= 1 .or. psin < 0) return
! Inside: use the interpolated numerical data
temp = self%temp_spline%eval(psin)
end function numeric_temp
function numeric_zeff(self, psin) result(zeff)
! function arguments
class(numeric_plasma), intent(in) :: self
real(wp_), intent(in) :: psin
real(wp_) :: zeff
zeff = 1
! Outside: Zeff is undefined, return 1
if (psin >= 1 .or. psin < 0) return
! Inside: use the interpolated numerical data
zeff = self%zeff_spline%eval(psin)
end function numeric_zeff
subroutine load_plasma(params, equil, plasma, err)
! Loads a generic plasma description from file (params%filenm)
use gray_params, only : gray_parameters
use gray_params, only : PROF_ANALYTIC, PROF_NUMERIC
use gray_params, only : RHO_TOR, RHO_POL, RHO_PSI
use gray_params, only : SCALE_COLLISION, SCALE_GREENWALD, SCALE_OFF
use gray_equil, only : abstract_equil
use logger, only : log_error, log_debug
! subroutine arguments
type(gray_parameters), intent(inout) :: params
class(abstract_equil), intent(in) :: equil
class(abstract_plasma), allocatable, intent(out) :: plasma
integer, intent(out) :: err
! local variables
integer :: u ! unit number
real(wp_) :: at, an, fact ! rescaling exponents and factor
type(numeric_plasma) :: np
type(analytic_plasma) :: ap
open(file=params%profiles%filenm, status='old', &
action='read', newunit=u, iostat=err)
if (err /= 0) then
call log_error( &
'opening profiles file ('//trim(params%profiles%filenm)//') failed!', &
mod='coreprofiles', proc='load_plasma')
err = 1
return
end if
! Set appropriate rescaling factors
select case(params%profiles%iscal)
case (SCALE_COLLISION)
call log_debug('rescaling preserving ν*', mod='coreprofiles', proc='load_plasma')
at = 2/3.0_wp_
an = 4/3.0_wp_
fact = params%equilibrium%factb
case (SCALE_GREENWALD)
call log_debug('rescaling preserving n_G', mod='coreprofiles', proc='load_plasma')
at = 1
an = 1
fact = params%equilibrium%factb
case (SCALE_OFF)
at = 1
an = 1
fact = 1
end select
select case(params%profiles%iprof)
case (PROF_ANALYTIC)
! The analytical file format is:
!
! 1 dens0 n1 n2
! 2 te0 te1 t1 t2
! 3 zeff
call log_debug('loading analytical file', mod='coreprofiles', proc='load_plasma')
! Read model parameters
read (u, *) ap%dens0, ap%n1, ap%n2
read (u, *) ap%te0, ap%te1, ap%t1, ap%t2
read (u, *) ap%z
! Apply rescaling factors
ap%dens0 = ap%dens0 * fact**an * params%profiles%factne
ap%te0 = ap%te0 * fact**at * params%profiles%factte
ap%te1 = ap%te1 * fact**at * params%profiles%factte
allocate(plasma, source=ap)
case (PROF_NUMERIC)
! The numerical file format is:
!
! 1 nrows
! 2 ψ_n T_e n_e Z_eff
! 3 ...
call log_debug('loading numerical file', mod='coreprofiles', proc='load_plasma')
block
integer :: i, nrows
real(wp_), allocatable, dimension(:) :: psi, temp, dens, zeff
! Read header and allocate arrays
read(u, *) nrows
allocate(psi(nrows), temp(nrows), dens(nrows), zeff(nrows))
! Read the table rows
do i=1, nrows
read(u, *) psi(i), temp(i), dens(i), zeff(i)
end do
! ??
psi(1) = max(psi(1), zero)
! Apply rescaling factors
dens = dens * fact**an * params%profiles%factne
temp = temp * fact**at * params%profiles%factte
! Convert psi to ψ_n (normalised poloidal flux)
select case (params%profiles%irho)
case (RHO_TOR) ! psi is ρ_t = √Φ_n (toroidal flux)
psi = [(equil%tor2pol(psi(i))**2, i=1,nrows)]
case (RHO_POL) ! psi is ρ_p = √ψ_n (poloidal flux)
psi = psi**2
case (RHO_PSI) ! psi is already ψ_n
end select
! Interpolate
call np%init(params, equil, psi, temp, dens, zeff, err)
end block
allocate(plasma, source=np)
end select
close(u)
end subroutine load_plasma
subroutine init_splines(self, params, equil, psi, temp, dens, zeff, err)
! Computes splines for the plasma profiles data and stores them
! in their respective global variables, see the top of this file.
!
! `err` is 1 if I/O errors occured, 2 if other initialisation failed.
use gray_params, only : gray_parameters
use gray_equil, only : abstract_equil
use logger, only : log_debug, log_info, log_warning, log_error
! subroutine arguments
class(numeric_plasma), intent(out) :: self
type(gray_parameters), intent(inout) :: params
class(abstract_equil), intent(in) :: equil
real(wp_), dimension(:), intent(in) :: psi, temp, dens, zeff
integer, intent(out) :: err
! local variables
integer :: n
! for log messages formatting
character(256) :: msg
n = size(psi)
! Spline interpolation of temperature and effective charge
call self%temp_spline%init(psi, temp, n)
call self%zeff_spline%init(psi, zeff, n)
! Spline interpolation of density (smooth spline)
call self%dens_spline%init(psi, dens, n, range=[zero, psi(n)], &
tension=params%profiles%sspld, err=err)
! Tension is very small, re-fit with an interpolating spline
if (err == -1) then
call log_warning('re-fitting with density curve with zero tension', &
mod='coreprofiles', proc='init_splines')
call self%dens_spline%init(psi, dens, n, range=[zero, psi(n)], tension=zero)
else if (err > 0) then
write (msg, '(a, g0)') 'density fit failed with error ', err
call log_error(msg, mod='coreprofiles', proc='init_splines')
err = 2
return
else
err = 0
end if
! Computation of the polynomial tail parameters
!
! Note: The density is the only quantity that needs to be evaluated
! at the edge. The spline thus has to be extended to transition
! smoothly from the last profile point to 0 outside the plasma.
block
real(wp_) :: s0, s1, s2 ! spline, 1st, 2nd derivative
real(wp_) :: delta4 ! discriminant Δ/4 of q(x)
real(wp_) :: x0, x1 ! vertex of q(x), solution
! Compute the coefficients of a 2nd order Taylor polinomial to
! extend the spline beyond the last point:
!
! p(ψ) = s(ψ₀) + (ψ - ψ₀)s'(ψ₀) + ½(ψ - ψ₀)²s"(ψ₀)
!
! where s(ψ) is the spline and ψ₀ the last point.
!
s0 = self%dens_spline%eval(psi(n))
s1 = self%dens_spline%deriv(psi(n), order=1)
s2 = self%dens_spline%deriv(psi(n), order=2)
! Determine where to end the tail (to ensure the density remains
! positive) from the zeros of the Taylor polynomial p(ψ)
!
! Define x=(ψ - ψ₀), then p(ψ)=0 is rewritten as
!
! q(x) = x² + 2s'/s" x + 2s/s" = 0
!
! The discriminant is Δ/4 = (s'/s")² - 2(s/s") and
! the solutions are x = x₀ ± √(Δ/4), with x₀ = -s'/s".
!
x0 = -s1 / s2 ! vertex of parabola y=q(x)
delta4 = (s1 / s2)**2 - 2*s0/s2 ! Δ/4 of q(x)
if (delta4 > 0) then
! Pick the smallest positive zero (implying >ψ₀)
x1 = x0 + sign(sqrt(delta4), sqrt(delta4) - x0)
else
! There are no zeros, use the parabola vertex
x1 = x0
call log_debug('spline extension has no zeros', &
mod='coreprofiles', proc='init_splines')
end if
! Store the tail parameters
self%tail%start = psi(n)
self%tail%end = self%tail%start + x1
self%tail%value = s0
self%tail%deriv1 = s1
self%tail%deriv2 = s2
end block
! Make sure the wave launcher does not fall inside the tail
! Note: if it does, the initial wave conditions become
! invalid as they are given assuming a vacuum (N=1)
block
use const_and_precisions, only : cm
real(wp_) :: R, Z, psi
! Convert from cartesian to cylindrical coordinates
associate (launch_pos => params%antenna%pos)
R = hypot(launch_pos(1), launch_pos(2)) * cm ! R = √(x²+y²)
z = launch_pos(3) * cm
end associate
! Get the poloidal flux at the launcher
! Note: this returns -1 when the data is not available
call equil%pol_flux(R, z, psi)
if (psi > self%tail%start .and. psi < self%tail%end) then
! Fall back to the midpoint of ψ₀ and the launcher ψ
call log_warning('downscaled tail to not reach the wave launcher', &
mod='coreprofiles', proc='init_splines')
write (msg, '(a, g0.3, a, g0.3, " → ", g0.3)') &
'launcher: ψ=', psi, ' boundary: ψ=', self%tail%end, (self%tail%start + psi)/2
call log_debug(msg, mod='coreprofiles', proc='init_splines')
self%tail%end = (self%tail%start + psi)/2
end if
if (psi > 0 .and. psi < self%tail%start) then
! This must be a user error, stop here
write (msg, '(a, a, g0.3, a, g0.3)') &
'wave launcher is inside the plasma! ', &
'launcher: ψ=', psi, ' boundary: ψ=', self%tail%end
call log_error(msg, mod='coreprofiles', proc='init_splines')
err = 2
return
end if
end block
! Set the density boundary ψ
! Note: this is used to detect entrance in the plasma
params%profiles%psnbnd = self%tail%end
write (msg, '(a,g0.4)') 'density boundary: ψ=', self%tail%end
call log_info(msg, mod='coreprofiles', proc='init_splines')
end subroutine init_splines
end module gray_plasma