src/coreprofiles.f90: use derived types

This commit is contained in:
Michele Guerini Rocco 2021-12-15 02:31:10 +01:00
parent 948a512254
commit f56e1cbc05
Signed by: rnhmjoj
GPG Key ID: BFBAF4C975F76450
2 changed files with 153 additions and 124 deletions

View File

@ -2,11 +2,11 @@ module coreprofiles
use const_and_precisions, only : wp_,zero,one use const_and_precisions, only : wp_,zero,one
implicit none implicit none
INTEGER, SAVE :: npp,nsfd integer, save :: npp, nsfd
REAL(wp_), SAVE :: psdbnd,psnpp,denpp,ddenpp,d2denpp real(wp_), save :: psdbnd, psnpp, denpp, ddenpp, d2denpp
REAL(wp_), DIMENSION(:), ALLOCATABLE, SAVE :: tfn,cfn,psrad real(wp_), dimension(:), allocatable, save :: tfn, cfn, psrad
REAL(wp_), DIMENSION(:,:), ALLOCATABLE, SAVE :: ct,cz real(wp_), dimension(:, :), allocatable, save :: ct, cz
REAL(wp_), SAVE :: dens0,aln1,aln2,te0,dte0,alt1,alt2,zeffan real(wp_), save :: dens0, aln1, aln2, te0, dte0, alt1, alt2, zeffan
contains contains
@ -125,35 +125,57 @@ contains
endif endif
end function fzeff end function fzeff
subroutine read_profiles(filenm,psin,te,ne,zeff,unit)
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 utils, only : get_free_unit
use gray_params, only : profiles_data
implicit none implicit none
! arguments
! subroutine arguments
character(len=*), intent(in) :: filenm character(len=*), intent(in) :: filenm
real(wp_), dimension(:), allocatable, intent(out) :: psin,te,ne,zeff type(profiles_data), intent(out) :: data
integer, optional, intent(in) :: unit integer, optional, intent(in) :: unit
! local variables ! local variables
integer :: u, i, n integer :: u, i, nrows
! 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)
if (present(unit)) then if (present(unit)) then
u = unit u = unit
else else
u = get_free_unit() u = get_free_unit()
end if end if
! Read number of rows and allocate the arrays
open(file=trim(filenm), status='old', action='read', unit=u) open(file=trim(filenm), status='old', action='read', unit=u)
read(u,*) n read(u, *) nrows
if(allocated(psin)) deallocate(psin) allocate(data%psrad(nrows), data%terad(nrows), &
if(allocated(te)) deallocate(te) data%derad(nrows), data%zfc(nrows))
if(allocated(ne)) deallocate(ne)
if(allocated(zeff)) deallocate(zeff) ! Read the table rows
allocate(psin(n),te(n),ne(n),zeff(n)) do i=1,nrows
do i=1,n read(u, *) data%psrad(i), data%terad(i), data%derad(i), data%zfc(i)
read(u,*) psin(i),te(i),ne(i),zeff(i)
end do end do
psin(1)=max(psin(1),zero)
close(u) close(u)
! ??
data%psrad(1) = max(data%psrad(1), zero)
end subroutine read_profiles end subroutine read_profiles
subroutine read_profiles_an(filenm,te,ne,zeff,unit) subroutine read_profiles_an(filenm,te,ne,zeff,unit)
use utils, only : get_free_unit use utils, only : get_free_unit
implicit none implicit none
@ -182,46 +204,65 @@ contains
close(u) close(u)
end subroutine read_profiles_an end subroutine read_profiles_an
subroutine tene_scal(te,ne,tfact,nfact,bfact,iscal,iprof)
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 implicit none
! arguments
real(wp_), dimension(:), intent(inout) :: te,ne ! subroutine arguments
real(wp_), intent(in) :: tfact,nfact,bfact type(profiles_parameters), intent(in) :: params
integer, intent(in) :: iscal,iprof real(wp_), intent(in) :: factb
type(profiles_data), intent(inout) :: data
! local variables ! local variables
real(wp_) :: aat, aan, ffact real(wp_) :: aat, aan, ffact
integer :: lastte,lastne integer :: last_te, last_ne
if (iscal==0) then if (params%iscal==0) then
aat = 2.0_wp_/3.0_wp_ aat = 2.0_wp_/3.0_wp_
aan = 4.0_wp_/3.0_wp_ aan = 4.0_wp_/3.0_wp_
else else
aat=1.0_wp_ aat = one
aan=1.0_wp_ aan = one
end if end if
if(iscal==2) then
ffact=1.0_wp_
else
ffact=bfact
end if
if (iprof==0) then
lastte=2
lastne=1
else
lastte=size(te)
lastne=size(ne)
end if
te(1:lastte)=te(1:lastte)*ffact**aat*tfact
ne(1:lastne)=ne(1:lastne)*ffact**aan*nfact
end subroutine tene_scal
subroutine set_prfspl(psin,te,ne,zeff,ssplne,psdbndmx) 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 simplespline, only : difcs
use dierckx, only : curfit, splev, splder use dierckx, only : curfit, splev, splder
use gray_params, only : profiles_parameters, profiles_data
implicit none implicit none
! arguments
real(wp_), dimension(:), intent(in) :: psin,te,ne,zeff ! subroutine arguments
real(wp_), intent(in) :: ssplne,psdbndmx type(profiles_parameters), intent(in) :: params
type(profiles_data), intent(inout) :: data
! local variables ! local variables
integer, parameter :: iopt=0, kspl=3 integer, parameter :: iopt=0, kspl=3
integer :: n, npest, lwrkf, ier integer :: n, npest, lwrkf, ier
@ -230,12 +271,12 @@ contains
integer, dimension(:), allocatable :: iwrkf integer, dimension(:), allocatable :: iwrkf
real(wp_), dimension(1) :: dedge,ddedge,d2dedge real(wp_), dimension(1) :: dedge,ddedge,d2dedge
n=size(psin) n=size(data%psrad)
npest=n+4 npest=n+4
lwrkf=n*4+npest*16 lwrkf=n*4+npest*16
allocate(wrkf(lwrkf),iwrkf(npest),wf(n)) allocate(wrkf(lwrkf),iwrkf(npest),wf(n))
ssplne_loc=ssplne ssplne_loc=params%sspld
! if necessary, reallocate spline arrays ! if necessary, reallocate spline arrays
if(.not.allocated(psrad)) then if(.not.allocated(psrad)) then
@ -255,37 +296,37 @@ contains
end if end if
end if end if
! spline approximation of temperature and Zeff ! spline approximation of temperature and data%zfc
call difcs(psin,te, n,iopt,ct,ier) call difcs(data%psrad,data%terad, n,iopt,ct,ier)
call difcs(psin,zeff,n,iopt,cz,ier) call difcs(data%psrad,data%zfc,n,iopt,cz,ier)
psrad=psin psrad=data%psrad
npp=n npp=n
! spline approximation of density ! spline approximation of density
xb=zero xb=zero
xe=psin(n) xe=data%psrad(n)
wf(:)=one wf(:)=one
call curfit(iopt,n,psin,ne,wf,xb,xe,kspl,ssplne_loc,npest, & call curfit(iopt,n,data%psrad,data%derad,wf,xb,xe,kspl,ssplne_loc,npest, &
nsfd,tfn,cfn,fp,wrkf,lwrkf,iwrkf,ier) nsfd,tfn,cfn,fp,wrkf,lwrkf,iwrkf,ier)
! if ier=-1 data are re-fitted using sspl=0 ! if ier=-1 data are re-fitted using sspl=0
if(ier==-1) then if(ier==-1) then
write(*,*) 'density curfit: ier=-1. Re-fitting with interpolating spline' write(*,*) 'density curfit: ier=-1. Re-fitting with interpolating spline'
ssplne_loc=0.0_wp_ ssplne_loc=0.0_wp_
call curfit(iopt,n,psin,ne,wf,xb,xe,kspl,ssplne_loc,npest, & call curfit(iopt,n,data%psrad,data%derad,wf,xb,xe,kspl,ssplne_loc,npest, &
nsfd,tfn,cfn,fp,wrkf,lwrkf,iwrkf,ier) nsfd,tfn,cfn,fp,wrkf,lwrkf,iwrkf,ier)
end if end if
! compute polinomial extrapolation matching the spline boundary up to the ! compute polinomial extrapolation matching the spline boundary up to the
! 2nd order derivative, extending the profile up to psi=psdbnd where ! 2nd order derivative, extending the profile up to psi=psdbnd where
! ne=ne'=ne''=0 ! data%derad=data%derad'=data%derad''=0
! spline value and derivatives at the edge ! spline value and derivatives at the edge
call splev(tfn,nsfd,cfn,kspl,psin(n:n),dedge(1:1),1,ier) call splev(tfn,nsfd,cfn,kspl,data%psrad(n:n),dedge(1:1),1,ier)
call splder(tfn,nsfd,cfn,kspl,1,psin(n:n),ddedge(1:1), 1,wrkf(1:nsfd),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,psin(n:n),d2dedge(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 ! determination of the boundary
psdbnd=psdbndmx psdbnd=params%psnbnd
psnpp=psin(n) psnpp=data%psrad(n)
denpp=dedge(1) denpp=dedge(1)
ddenpp=ddedge(1) ddenpp=ddedge(1)
d2denpp=d2dedge(1) d2denpp=d2dedge(1)
@ -302,12 +343,13 @@ contains
else if (xxp > psnpp) then else if (xxp > psnpp) then
psdbnd=min(psdbnd,xxp) psdbnd=min(psdbnd,xxp)
end if end if
write(*,*) 'density psdbnd=',psdbnd write(*,'(a,f5.3)') 'density psdbnd=', psdbnd
end if end if
deallocate(iwrkf,wrkf,wf) deallocate(iwrkf,wrkf,wf)
end subroutine set_prfspl end subroutine set_prfspl
subroutine unset_prfspl subroutine unset_prfspl
implicit none implicit none

View File

@ -21,7 +21,7 @@ program main
! Read the input data into set the global variables ! Read the input data into set the global variables
! of the respective module. Note: order matters. ! of the respective module. Note: order matters.
call init_equilibrium(params, data) call init_equilibrium(params, data)
call init_profiles(params, data) call init_profiles(params%profiles, params%equilibrium%factb, data%profiles)
call init_antenna(params%antenna) call init_antenna(params%antenna)
call init_misc(params, data) call init_misc(params, data)
@ -182,63 +182,50 @@ contains
end subroutine deinit_equilibrium end subroutine deinit_equilibrium
subroutine init_profiles(params, data) subroutine init_profiles(params, factb, data)
! Reads the plasma kinetic profiles file (containing the elecron ! Reads the plasma kinetic profiles file (containing the elecron
! temperature, density and plasma effective charge) and initialises ! temperature, density and plasma effective charge) and initialises
! the respective GRAY data structure. ! the respective GRAY data structure.
use gray_params, only : profiles_parameters, profiles_data use gray_params, only : profiles_parameters, profiles_data
use equilibrium, only : frhopolv use equilibrium, only : frhopolv
use coreprofiles, only : read_profiles_an, read_profiles, tene_scal, & use coreprofiles, only : read_profiles_an, read_profiles, &
set_prfan, set_prfspl scale_profiles, set_prfan, set_prfspl
implicit none implicit none
! subroutine arguments ! subroutine arguments
type(gray_parameters), intent(in) :: params type(profiles_parameters), intent(in) :: params
type(gray_data), intent(inout), target :: data real(wp_), intent(in) :: factb
type(profiles_data), intent(out) :: data
! local variables if(params%iprof == 0) then
type(profiles_parameters) :: profp
type(profiles_data), pointer :: profd
! Radial coordinate (depending on profp%irho: ρ_t, ρ_p, or ψ)
real(wp_), allocatable :: xrad(:)
profp = params%profiles
profd => data%profiles
if(params%profiles%iprof == 0) then
! Analytical profiles ! Analytical profiles
! TODO: rewrite using derived type ! TODO: rewrite using derived type
call read_profiles_an(profp%filenm, profd%terad, profd%derad, profd%zfc) call read_profiles_an(params%filenm, data%terad, data%derad, data%zfc)
else else
! Numerical profiles ! Numerical profiles
call read_profiles(profp%filenm, xrad, profd%terad, profd%derad, profd%zfc) call read_profiles(params%filenm, data)
allocate(profd%psrad(size(xrad))) ! Convert psrad to ψ
select case (params%irho)
select case (profp%irho) case (0) ! psrad is ρ_t
case (0) ! xrad is rhot data%psrad = frhopolv(data%psrad)**2
profd%psrad = frhopolv(xrad)**2 case (1) ! psrad is ρ_p
case (1) ! xrad is rhop data%psrad = data%psrad**2
profd%psrad = xrad**2 case default ! psrad is already ψ
case default ! xrad is psi
profd%psrad = xrad
end select end select
deallocate(xrad)
end if end if
! Rescale input data ! Rescale input data
! TODO: rewrite using derived type call scale_profiles(params, factb, data)
call tene_scal(profd%terad, profd%derad, profp%factte, profp%factne, &
params%equilibrium%factb, profp%iscal, profp%iprof)
! Set global variables ! Set global variables
if(params%iprof == 0) then
! Analytical profiles
! TODO: rewrite using derived type ! TODO: rewrite using derived type
if(params%profiles%iprof == 0) then call set_prfan(data%terad, data%derad, data%zfc)
call set_prfan(profd%terad, profd%derad, profd%zfc)
else else
call set_prfspl(profd%psrad, profd%terad, profd%derad, profd%zfc, & ! Numerical profiles
profp%sspld, profp%psnbnd) call set_prfspl(params, data)
end if end if
end subroutine init_profiles end subroutine init_profiles