398 lines
11 KiB
Fortran
398 lines
11 KiB
Fortran
module coreprofiles
|
|
use const_and_precisions, only : wp_,zero,one
|
|
implicit none
|
|
|
|
integer, save :: npp, nsfd
|
|
real(wp_), save :: psdbnd, psnpp, denpp, ddenpp, d2denpp
|
|
real(wp_), dimension(:), allocatable, save :: tfn, cfn, psrad
|
|
real(wp_), dimension(:, :), allocatable, save :: ct, cz
|
|
real(wp_), save :: dens0, aln1, aln2, te0, dte0, alt1, alt2, zeffan
|
|
|
|
contains
|
|
|
|
subroutine density(psin,dens,ddens)
|
|
use gray_params, only : iprof
|
|
use dierckx, only : splev,splder
|
|
use logger, only : log_error
|
|
|
|
implicit none
|
|
|
|
! subroutine arguments
|
|
real(wp_), intent(in) :: psin
|
|
real(wp_), intent(out) :: dens,ddens
|
|
|
|
! local variables
|
|
integer :: ier,nu
|
|
real(wp_) :: profd,dprofd,dpsib,tt,fp,dfp,fh,dfh
|
|
real(wp_), dimension(1) :: xxs,ffs
|
|
real(wp_), dimension(npp+4) :: wrkfd
|
|
character(256) :: msg
|
|
|
|
!
|
|
! Computation of density [10¹⁹ m⁻³] and derivative wrt ψ
|
|
!
|
|
dens=zero
|
|
ddens=zero
|
|
if((psin >= psdbnd) .or. (psin < zero)) return
|
|
|
|
if(iprof == 0) then
|
|
if(psin > one) return
|
|
profd=(one-psin**aln1)**aln2
|
|
dens=dens0*profd
|
|
dprofd=-aln1*aln2*psin**(aln1-one) &
|
|
*(one-psin**aln1)**(aln2-one)
|
|
ddens=dens0*dprofd
|
|
else
|
|
if(psin > psnpp) then
|
|
|
|
! Smooth interpolation for psnpp < psi < psdbnd
|
|
! dens = fp * fh
|
|
! fp: parabola matched at psi=psnpp with given profile density
|
|
! fh=(1-t)^3(1+3t+6t^2) is a smoothing function:
|
|
! fh(0)=1, fh(1)=0 and zero first and second deriv at t=0,1
|
|
!
|
|
dpsib=psin-psnpp
|
|
fp=denpp+dpsib*ddenpp+0.5_wp_*dpsib**2*d2denpp
|
|
dfp=ddenpp+dpsib*d2denpp
|
|
tt=dpsib/(psdbnd-psnpp)
|
|
fh=(one-tt)**3*(one+3.0_wp_*tt+6.0_wp_*tt**2)
|
|
dfh=-30.0_wp_*(one-tt)**2*tt**2/(psdbnd-psnpp)
|
|
dens=fp*fh
|
|
ddens=dfp*fh+fp*dfh
|
|
else
|
|
xxs(1)=psin
|
|
ier=0
|
|
call splev(tfn,nsfd,cfn,3,xxs,ffs,1,ier)
|
|
dens=ffs(1)
|
|
nu=1
|
|
ier=0
|
|
call splder(tfn,nsfd,cfn,3,nu,xxs,ffs,1,wrkfd,ier)
|
|
ddens=ffs(1)
|
|
if(abs(dens) < 1.0e-10_wp_) dens=zero
|
|
end if
|
|
if(dens < zero) then
|
|
write (msg, '("negative density:", 2(x,a,"=",g0.3))') &
|
|
'ne', dens, 'ψ', psin
|
|
call log_error(msg, mod='coreprofiles', proc='density')
|
|
end if
|
|
end if
|
|
end subroutine density
|
|
|
|
|
|
function temp(psin)
|
|
use const_and_precisions, only : wp_,zero,one
|
|
use gray_params, only : iprof
|
|
use utils, only : locate
|
|
use simplespline, only :spli
|
|
implicit none
|
|
! arguments
|
|
real(wp_), intent(in) :: psin
|
|
real(wp_) :: temp
|
|
! local variables
|
|
integer :: k
|
|
real(wp_) :: proft,dps
|
|
|
|
temp=zero
|
|
if((psin >= one).or.(psin < zero)) return
|
|
if(iprof == 0) then
|
|
proft=(1.0_wp_-psin**alt1)**alt2
|
|
temp=(te0-dte0)*proft+dte0
|
|
else
|
|
call locate(psrad,npp,psin,k)
|
|
k=max(1,min(k,npp-1))
|
|
dps=psin-psrad(k)
|
|
temp=spli(ct,npp,k,dps)
|
|
endif
|
|
end function temp
|
|
|
|
function fzeff(psin)
|
|
use const_and_precisions, only : wp_,zero,one
|
|
use gray_params, only : iprof
|
|
use utils, only : locate
|
|
use simplespline, only :spli
|
|
implicit none
|
|
! arguments
|
|
real(wp_), intent(in) :: psin
|
|
real(wp_) :: fzeff
|
|
! local variables
|
|
integer :: k
|
|
real(wp_) :: dps
|
|
|
|
fzeff=one
|
|
if((psin >= one).or.(psin < zero)) return
|
|
if(iprof == 0) then
|
|
fzeff=zeffan
|
|
else
|
|
call locate(psrad,npp,psin,k)
|
|
k=max(1,min(k,npp-1))
|
|
dps=psin-psrad(k)
|
|
fzeff=spli(cz,npp,k,dps)
|
|
endif
|
|
end function fzeff
|
|
|
|
|
|
subroutine read_profiles(filenm, data, unit)
|
|
! Reads the radial plasma profiles from `file` and store them
|
|
! into `data`. If given, the file is opened in the `unit` number.
|
|
! Format notes:
|
|
! 1. The file is formatted as a table with the following columns:
|
|
! radial coordinate, temperature, density, effective charge.
|
|
! 2. The first line is a header specifying the number of rows.
|
|
use utils, only : get_free_unit
|
|
use gray_params, only : profiles_data
|
|
use logger, only : log_error
|
|
|
|
implicit none
|
|
|
|
! subroutine arguments
|
|
character(len=*), intent(in) :: filenm
|
|
type(profiles_data), intent(out) :: data
|
|
integer, optional, intent(in) :: unit
|
|
|
|
! local variables
|
|
integer :: u, i, nrows
|
|
integer :: err
|
|
|
|
! Free the arrays when already allocated
|
|
if(allocated(data%psrad)) deallocate(data%psrad)
|
|
if(allocated(data%terad)) deallocate(data%terad)
|
|
if(allocated(data%derad)) deallocate(data%derad)
|
|
if(allocated(data%zfc)) deallocate(data%zfc)
|
|
|
|
u = get_free_unit(unit)
|
|
|
|
! Read number of rows and allocate the arrays
|
|
open(file=filenm, status='old', action='read', unit=u, iostat=err)
|
|
if (err /= 0) then
|
|
call log_error('opening profiles file ('//trim(filenm)//') failed!', &
|
|
mod='coreprofiles', proc="read_profiles")
|
|
call exit(1)
|
|
end if
|
|
|
|
read(u, *) nrows
|
|
allocate(data%psrad(nrows), data%terad(nrows), &
|
|
data%derad(nrows), data%zfc(nrows))
|
|
|
|
! Read the table rows
|
|
do i=1,nrows
|
|
read(u, *) data%psrad(i), data%terad(i), data%derad(i), data%zfc(i)
|
|
end do
|
|
close(u)
|
|
|
|
! ??
|
|
data%psrad(1) = max(data%psrad(1), zero)
|
|
|
|
end subroutine read_profiles
|
|
|
|
|
|
subroutine read_profiles_an(filenm,te,ne,zeff,unit)
|
|
use utils, only : get_free_unit
|
|
use logger, only : log_error
|
|
|
|
implicit none
|
|
! arguments
|
|
character(len=*), intent(in) :: filenm
|
|
real(wp_), dimension(:), allocatable, intent(out) :: te,ne,zeff
|
|
integer, optional, intent(in) :: unit
|
|
! local variables
|
|
integer :: u
|
|
integer :: err
|
|
|
|
u = get_free_unit(unit)
|
|
|
|
if(allocated(te)) deallocate(te)
|
|
if(allocated(ne)) deallocate(ne)
|
|
if(allocated(zeff)) deallocate(zeff)
|
|
allocate(te(4),ne(3),zeff(1))
|
|
|
|
open(file=filenm, status='old', action='read', unit=u, iostat=err)
|
|
if (err /= 0) then
|
|
call log_error('opening profiles file ('//trim(filenm)//') failed!', &
|
|
mod='coreprofiles', proc='read_profiles_an')
|
|
call exit(1)
|
|
end if
|
|
|
|
read(u,*) ne(1:3) ! dens0,aln1,aln2
|
|
read(u,*) te(1:4) ! te0,dte0,alt1,alt2
|
|
read(u,*) zeff(1) ! zeffan
|
|
close(u)
|
|
end subroutine read_profiles_an
|
|
|
|
|
|
subroutine scale_profiles(params, factb, data)
|
|
! Rescales the temperature and density profiles by `params%factte`
|
|
! and `params%factne`, respectively, according to the model
|
|
! specified by `params%iscal`.
|
|
! See the GRAY user manual for the explanation.
|
|
use gray_params, only : profiles_parameters, profiles_data
|
|
|
|
implicit none
|
|
|
|
! subroutine arguments
|
|
type(profiles_parameters), intent(in) :: params
|
|
real(wp_), intent(in) :: factb
|
|
type(profiles_data), intent(inout) :: data
|
|
|
|
! local variables
|
|
real(wp_) :: aat, aan, ffact
|
|
integer :: last_te, last_ne
|
|
|
|
if (params%iscal==0) then
|
|
aat = 2.0_wp_/3.0_wp_
|
|
aan = 4.0_wp_/3.0_wp_
|
|
else
|
|
aat = one
|
|
aan = one
|
|
end if
|
|
|
|
if(params%iscal==2) then
|
|
ffact = one
|
|
else
|
|
ffact = factb
|
|
end if
|
|
|
|
if(params%iprof==0) then
|
|
last_te = 2
|
|
last_ne = 1
|
|
else
|
|
last_te = size(data%terad)
|
|
last_ne = size(data%derad)
|
|
end if
|
|
|
|
data%terad(1:last_te) = data%terad(1:last_te) * ffact**aat * params%factte
|
|
data%derad(1:last_ne) = data%derad(1:last_ne) * ffact**aan * params%factne
|
|
end subroutine scale_profiles
|
|
|
|
|
|
subroutine set_prfspl(params, data)
|
|
! Computes splines for the plasma profiles data and stores them
|
|
! in their respective global variables, see the top of this file.
|
|
use simplespline, only : difcs
|
|
use dierckx, only : curfit, splev, splder
|
|
use gray_params, only : profiles_parameters, profiles_data
|
|
use logger, only : log_info, log_warning
|
|
|
|
implicit none
|
|
|
|
! subroutine arguments
|
|
type(profiles_parameters), intent(in) :: params
|
|
type(profiles_data), intent(inout) :: data
|
|
|
|
! local variables
|
|
integer, parameter :: iopt=0, kspl=3
|
|
integer :: n, npest, lwrkf, ier
|
|
real(wp_) :: xb, xe, fp, xnv, xxp,xxm,delta2,ssplne_loc
|
|
real(wp_), dimension(:), allocatable :: wf, wrkf
|
|
integer, dimension(:), allocatable :: iwrkf
|
|
real(wp_), dimension(1) :: dedge,ddedge,d2dedge
|
|
character(256) :: msg ! for log messages formatting
|
|
|
|
n=size(data%psrad)
|
|
npest=n+4
|
|
lwrkf=n*4+npest*16
|
|
allocate(wrkf(lwrkf),iwrkf(npest),wf(n))
|
|
|
|
ssplne_loc=params%sspld
|
|
|
|
! if necessary, reallocate spline arrays
|
|
if(.not.allocated(psrad)) then
|
|
allocate(psrad(n),ct(n,4),cz(n,4))
|
|
else
|
|
if(size(psrad)<n) then
|
|
deallocate(psrad,ct,cz)
|
|
allocate(psrad(n),ct(n,4),cz(n,4))
|
|
end if
|
|
end if
|
|
if(.not.allocated(cfn)) then
|
|
allocate(tfn(npest),cfn(npest))
|
|
else
|
|
if(size(cfn)<npest) then
|
|
deallocate(tfn,cfn)
|
|
allocate(tfn(npest),cfn(npest))
|
|
end if
|
|
end if
|
|
|
|
! spline approximation of temperature and data%zfc
|
|
call difcs(data%psrad,data%terad, n,iopt,ct,ier)
|
|
call difcs(data%psrad,data%zfc,n,iopt,cz,ier)
|
|
psrad=data%psrad
|
|
npp=n
|
|
|
|
! spline approximation of density
|
|
xb=zero
|
|
xe=data%psrad(n)
|
|
wf(:)=one
|
|
call curfit(iopt,n,data%psrad,data%derad,wf,xb,xe,kspl,ssplne_loc,npest, &
|
|
nsfd,tfn,cfn,fp,wrkf,lwrkf,iwrkf,ier)
|
|
! if ier=-1 data are re-fitted using sspl=0
|
|
if(ier == -1) then
|
|
call log_warning('curfit failed with error -1: re-fitting with '// &
|
|
's=0', mod='coreprofiles', proc='density')
|
|
ssplne_loc=0.0_wp_
|
|
call curfit(iopt,n,data%psrad,data%derad,wf,xb,xe,kspl,ssplne_loc,npest, &
|
|
nsfd,tfn,cfn,fp,wrkf,lwrkf,iwrkf,ier)
|
|
end if
|
|
|
|
! compute polinomial extrapolation matching the spline boundary up to the
|
|
! 2nd order derivative, extending the profile up to psi=psdbnd where
|
|
! data%derad=data%derad'=data%derad''=0
|
|
! spline value and derivatives at the edge
|
|
call splev(tfn,nsfd,cfn,kspl,data%psrad(n:n),dedge(1:1),1,ier)
|
|
call splder(tfn,nsfd,cfn,kspl,1,data%psrad(n:n),ddedge(1:1), 1,wrkf(1:nsfd),ier)
|
|
call splder(tfn,nsfd,cfn,kspl,2,data%psrad(n:n),d2dedge(1:1),1,wrkf(1:nsfd),ier)
|
|
! determination of the boundary
|
|
psdbnd=params%psnbnd
|
|
|
|
psnpp=data%psrad(n)
|
|
denpp=dedge(1)
|
|
ddenpp=ddedge(1)
|
|
d2denpp=d2dedge(1)
|
|
|
|
delta2=(ddenpp/d2denpp)**2-2.0_wp_*denpp/d2denpp
|
|
xnv=psnpp-ddenpp/d2denpp
|
|
if(delta2 < zero) then
|
|
! if(xnv > psnpp) psdbnd=min(psdbnd,xnv)
|
|
else
|
|
xxm=xnv-sqrt(delta2)
|
|
xxp=xnv+sqrt(delta2)
|
|
if(xxm > psnpp) then
|
|
psdbnd=min(psdbnd,xxm)
|
|
else if (xxp > psnpp) then
|
|
psdbnd=min(psdbnd,xxp)
|
|
end if
|
|
write (msg, '(a,g0.3)') 'density boundary: ψ=', psdbnd
|
|
call log_info(msg, mod="coreprofiles", proc="set_prfspl")
|
|
end if
|
|
|
|
deallocate(iwrkf,wrkf,wf)
|
|
end subroutine set_prfspl
|
|
|
|
|
|
subroutine unset_prfspl
|
|
implicit none
|
|
|
|
if(allocated(psrad)) deallocate(psrad)
|
|
if(allocated(ct)) deallocate(ct)
|
|
if(allocated(cz)) deallocate(cz)
|
|
if(allocated(tfn)) deallocate(tfn)
|
|
if(allocated(cfn)) deallocate(cfn)
|
|
end subroutine unset_prfspl
|
|
|
|
subroutine set_prfan(te,ne,zeff)
|
|
implicit none
|
|
REAL(wp_), dimension(:), intent(in) :: te,ne,zeff
|
|
|
|
te0=te(1)
|
|
dte0=te(2)
|
|
alt1=te(3)
|
|
alt2=te(4)
|
|
dens0=ne(1)
|
|
aln1=ne(2)
|
|
aln2=ne(3)
|
|
zeffan=zeff(1)
|
|
|
|
psdbnd=1.0_wp_
|
|
end subroutine set_prfan
|
|
|
|
end module coreprofiles
|