src/gray_params.f90: formatting

This commit is contained in:
Michele Guerini Rocco 2021-12-15 02:30:52 +01:00
parent 75c30ee834
commit 023facac6b
Signed by: rnhmjoj
GPG Key ID: BFBAF4C975F76450

View File

@ -40,7 +40,6 @@ module gray_params
integer :: ipec, nrho, istpr, istpl
end type outparam_type
integer, save :: iequil,iprof,ipol
integer, save :: iwarm,ilarm,imx,ieccd
integer, save :: igrad,idst,ipass
@ -48,11 +47,12 @@ module gray_params
integer, save :: ipec,nnd
contains
subroutine print_params(rtrparam,hcdparam,antctrl,eqparam,rwall, &
prfparam,outparam,strout)
prfparam,outparam,strout)
implicit none
! arguments
! arguments
type(rtrparam_type), intent(in) :: rtrparam
type(hcdparam_type), intent(in) :: hcdparam
type(antctrl_type), intent(in) :: antctrl
@ -61,29 +61,33 @@ contains
type(prfparam_type), intent(in) :: prfparam
type(outparam_type), intent(in) :: outparam
character(len=*), dimension(:), intent(out) :: strout ! min len=110, dimension(21)
! local variables
! local variables
character(len=8) :: rdat
character(len=10) :: rtim
#ifndef REVISION
character(len=*), parameter :: REVISION="unknown"
#endif
! date and time
call date_and_time(rdat,rtim)
! date and time
call date_and_time(rdat,rtim)
write(strout(1),'("# Run date/time: ",a4,2("/",a2),1x,2(a2,":"),a6)') &
rdat(1:4),rdat(5:6),rdat(7:8),rtim(1:2),rtim(3:4),rtim(5:10)
! SVN revision
write(strout(2),'("# GRAY SVN revision: ",a)') REVISION
! equilibrium input data
! Git revision
write(strout(2),'("# GRAY Git revision: ",a)') REVISION
! equilibrium input data
if (eqparam%iequil > 0) then
write(strout(3),'("# EQL input: ",a)') trim(eqparam%filenm)
!!!!!!! missing values
!!! missing values
write(strout(7),'("# EQL B0 R0 aminor Rax zax:",5(1x,e12.5))') &
0._wp_, 0._wp_, 0._wp_, 0._wp_, 0._wp_
else
write(strout(3),'("# EQL input: N/A (vacuum)")')
write(strout(7),'("# EQL B0 R0 aminor Rax zax: N/A (vacuum)")')
end if
write(strout(4),'("# EQL iequil sgnb sgni factb:",3(1x,i4),1x,e12.5)') &
eqparam%iequil, eqparam%sgnb, eqparam%sgni, eqparam%factb
if (eqparam%iequil > 1) then
@ -95,7 +99,8 @@ contains
write(strout(5),'("# EQL icocos ipsinorm idesc ifreefmt: N/A (analytical)")')
write(strout(6),'("# EQL ssplps ssplf ixp: N/A (analytical)")')
end if
! profiles input data
! profiles input data
if (eqparam%iequil > 0) then
write(strout(8),'("# PRF input: ",a)') trim(prfparam%filenm)
write(strout(9),'("# PRF iprof iscal factne factte:",2(1x,i4),2(1x,e12.5))') &
@ -106,7 +111,7 @@ contains
else
write(strout(10),'("# PRF irho psnbnd sspld: N/A (analytical)")')
end if
!!!!!!! missing values
!!! missing values
write(strout(11),'("# PRF Te0 ne0 Zeff0:",3(1x,e12.5))') &
0._wp_, 0._wp_, 0._wp_
else
@ -115,21 +120,24 @@ contains
write(strout(10),'("# PRF irho psnbnd sspld: N/A (vacuum)")')
write(strout(11),'("# PRF Te0 ne0 Zeff0: N/A (vacuum)")')
end if
! launch parameters
! launch parameters
write(strout(12),'("# ANT input: ",a)') trim(antctrl%filenm)
write(strout(13),'("# ANT ibeam iox psi chi:",2(1x,i4),2(1x,e12.5))') &
antctrl%ibeam, antctrl%iox, antctrl%psi, antctrl%chi
write(strout(14),'("# ANT alpha beta power:",3(1x,e12.5))') &
antctrl%alpha, antctrl%beta, antctrl%power
!!!!!!! missing values
!!! missing values
write(strout(15),'("# ANT x0 y0 z0:",3(1x,e12.5))') &
0._wp_, 0._wp_, 0._wp_
!!!!!!! missing values
!!! missing values
write(strout(16),'("# ANT wx wy Rcix Rciy psiw psir:",6(1x,e12.5))') &
0._wp_, 0._wp_, 0._wp_, 0._wp_, 0._wp_, 0._wp_
! wall parameters
! wall parameters
write(strout(17),'("# RFL rwall:",1x,e12.5)') rwall
! code parameters
! code parameters
write(strout(18),'("# COD igrad idst ipass ipol:",4(1x,i4))') &
rtrparam%igrad, rtrparam%idst, rtrparam%ipass, rtrparam%ipol
write(strout(19),'("# COD nrayr nrayth nstep rwmax dst:",3(1x,i4),2(1x,e12.5))') &
@ -139,12 +147,14 @@ contains
write(strout(21),'("# COD ipec nrho istpr istpl:",4(1x,i4))') &
outparam%ipec, outparam%nrho, outparam%istpr, outparam%istpl
end subroutine print_params
subroutine read_params(filenm,rtrparam,hcdparam,antctrl,eqparam,rwall, &
prfparam,outparam,unit)
use utils, only : get_free_unit
implicit none
! arguments
! arguments
character(len=*), intent(in) :: filenm
type(rtrparam_type), intent(out) :: rtrparam
type(hcdparam_type), intent(out) :: hcdparam
@ -154,7 +164,8 @@ contains
type(prfparam_type), intent(out) :: prfparam
type(outparam_type), intent(out) :: outparam
integer, intent(in), optional :: unit
! local variables
! local variables
integer :: u
if (present(unit)) then
@ -162,91 +173,103 @@ contains
else
u = get_free_unit()
end if
open(u,file=filenm,status= 'old',action='read')
open(u,file=filenm,status='old',action='read',iostat=iostat)
! ==========================================================================
! nrayr number of rays in radial direction
! nrayth number of rays in angular direction
! rwmax normalized maximum radius of beam power
! rwmax=1 -> last ray at radius = waist
! Raytracing
! ========================================================================
! nrayr :number of rays in radial direction
! nrayth :number of rays in angular direction
! rwmax :normalized maximum radius of beam power
! rwmax=1 -> :last ray at radius = waist
read(u,*) rtrparam%nrayr, rtrparam%nrayth, rtrparam%rwmax
! igrad=0 optical ray-tracing, initial conditions as for beam
! igrad=1 quasi-optical ray-tracing
! igrad=-1 ray-tracing, init. condit.
! from center of mirror and with angular spread
! ipass=1/2 1 or 2 passes into plasma
! ipol=0 compute mode polarization at antenna, ipol=1 use polariz angles
! igrad=0 :optical ray-tracing, initial conditions as for beam
! igrad=1 :quasi-optical ray-tracing
! igrad=-1 :ray-tracing, init. condit.
! from center of mirror and with angular spread
! ipass=1/2 :1 or 2 passes into plasma
! ipol=0 :compute mode polarization at antenna, ipol=1 use polariz angles
read(u,*) rtrparam%igrad, rtrparam%ipass, rtrparam%ipol
! dst integration step
! nstep maximum number of integration steps
! idst=0/1/2 0 integration in s, 1 integr. in ct, 2 integr. in Sr
! dst :integration step
! nstep :maximum number of integration steps
! idst=0/1/2 :0 integration in s, 1 integr. in ct, 2 integr. in Sr
read(u,*) rtrparam%dst, rtrparam%nstep, rtrparam%idst
! ==========================================================================
! iwarm=0 :no absorption and cd
! iwarm=1 :weakly relativistic absorption
! iwarm=2 :relativistic absorption, n<1 asymptotic expansion
! iwarm=3 :relativistic absorption, numerical integration
! ilarm :order of larmor expansion
! imx :max n of iterations in dispersion, imx<0 uses 1st
! iteration in case of failure after |imx| iterations
! Heating & Current drive
! ========================================================================
! iwarm=0 :no absorption and cd
! iwarm=1 :weakly relativistic absorption
! iwarm=2 :relativistic absorption, n<1 asymptotic expansion
! iwarm=3 :relativistic absorption, numerical integration
! ilarm :order of larmor expansion
! imx :max n of iterations in dispersion, imx<0 uses 1st
! iteration in case of failure after |imx| iterations
read(u,*) hcdparam%iwarm,hcdparam%ilarm,hcdparam%imx
! ieccd 0/1 NO/YES ECCD calculation ieccd>0 different CD models
! ieccd 0/1 NO/YES ECCD calculation ieccd>0 different CD models
read(u,*) hcdparam%ieccd
! ==========================================================================
! alpha0, beta0 (cartesian) launching angles
read(u,*) antctrl%alpha, antctrl%beta
! p0mw injected power (MW)
read(u,*) antctrl%power
! abs(iox)=1/2 OM/XM
! psipol0,chipol0 polarization angles at the antenna (if iox<0)
read(u,*) antctrl%iox, antctrl%psi, antctrl%chi
! ibeam=0 :read data for beam as above
! ibeam=1 :read data from file simple astigmatic beam
! ibeam=2 :read data from file general astigmatic beam
! Wave launcher
! ========================================================================
! alpha0, beta0 (cartesian) launching angles
read(u,*) antctrl%alpha, antctrl%beta
! p0mw injected power (MW)
read(u,*) antctrl%power
! abs(iox)=1/2 OM/XM
! psipol0,chipol0 polarization angles at the antenna (if iox<0)
read(u,*) antctrl%iox, antctrl%psi, antctrl%chi
! ibeam=0 :read data for beam as above
! ibeam=1 :read data from file simple astigmatic beam
! ibeam=2 :read data from file general astigmatic beam
read(u,*) antctrl%ibeam
read(u,*) antctrl%filenm
! ==========================================================================
! iequil=0 :vacuum
! iequil=1 :analytical equilibrium
! iequil=2 :read eqdsk
! Magnetic equilibrium
! ========================================================================
! iequil=0 :vacuum
! iequil=1 :analytical equilibrium
! iequil=2 :read eqdsk
read(u,*) eqparam%iequil
read(u,*) eqparam%filenm
! icocos :index for equilibrium from COCOS - O. Sauter Feb 2012
! ipsinorm :0 standard EQDSK format, 1 format Portone summer 2004
! icocos :index for equilibrium from COCOS - O. Sauter Feb 2012
! ipsinorm :0 standard EQDSK format, 1 format Portone summer 2004
read(u,*) eqparam%icocos, eqparam%ipsinorm, eqparam%idesc, eqparam%ifreefmt
! ixp=0,-1,+1 : no X point , bottom/up X point
! ssplps : spline parameter for psi interpolation
! ixp=0,-1,+1 : no X point , bottom/up X point
! ssplps : spline parameter for psi interpolation
read(u,*) eqparam%ixp, eqparam%ssplps !, eqparam%ssplf
eqparam%ssplf=0.01_wp_
! signum of toroidal B and I
! factb factor for magnetic field (only for numerical equil)
! scaling adopted: beta=const, qpsi=const, nustar=const
! signum of toroidal B and I
! factb factor for magnetic field (only for numerical equil)
! scaling adopted: beta=const, qpsi=const, nustar=const
read(u,*) eqparam%sgnb, eqparam%sgni, eqparam%factb
! ==========================================================================
! Wall
! ========================================================================
read(u,*) rwall
! ==========================================================================
! iprof=0 :analytical density and temp. profiles
! iprof>0 :numerical density and temp. profiles
! Profiles
! ========================================================================
! iprof=0 :analytical density and temp. profiles
! iprof>0 :numerical density and temp. profiles
read(u,*) prfparam%iprof, prfparam%irho ! irho=0,1,2 -> num profiles vs rhot,rhop,psin
read(u,*) prfparam%filenm
! psbnd value of psi ( > 1 ) of density boundary
! psbnd value of psi ( > 1 ) of density boundary
read(u,*) prfparam%psnbnd, prfparam%sspld
! prfparam%sspld=0.001_wp_
! iscal ne Te scaling 0: nustar=const, 1: n_greenw=const; 2 no rescaling
! factT factn factor for Te&ne scaling
! prfparam%sspld=0.001_wp_
! iscal :ne Te scaling 0: nustar=const, 1: n_greenw=const; 2 no rescaling
! factT factn :factor for Te&ne scaling
read(u,*) prfparam%factte, prfparam%factne, prfparam%iscal
! ==========================================================================
! ipec=0/1 :pec profiles grid in psi/rhop
! nrho :number of grid steps for pec profiles +1
! Output
! ========================================================================
! ipec=0/1 :pec profiles grid in psi/rhop
! nrho :number of grid steps for pec profiles +1
read(u,*) outparam%ipec, outparam%nrho
! istpr0 projection step = dsdt*istprj
! istpl0 plot step = dsdt*istpl
! istpr0 :projection step = dsdt*istprj
! istpl0 :plot step = dsdt*istpl
read(u,*) outparam%istpr, outparam%istpl
close(u)
end subroutine read_params
subroutine set_codepar(eqparam,prfparam,outparam,rtrparam,hcdparam)
implicit none
type(eqparam_type), intent(in) :: eqparam