gray/src/main.f90

139 lines
5.2 KiB
Fortran

program main_std
use const_and_precisions, only : wp_,one
use graycore, only : gray_main
use gray_params, only : read_params, antctrl_type,eqparam_type, &
prfparam_type,outparam_type,rtrparam_type,hcdparam_type
use beams, only : read_beam0, read_beam1, read_beam2
use equilibrium, only : read_equil_an,read_eqdsk,change_cocos,eq_scal, &
set_rhospl,setqphi_num,frhopolv
use coreprofiles, only : read_profiles_an,read_profiles,tene_scal
use reflections, only : range2rect
implicit none
type(antctrl_type) :: antp
type(eqparam_type) :: eqp
type(prfparam_type) :: prfp
type(outparam_type) :: outp
type(rtrparam_type) :: rtrp
type(hcdparam_type) :: hcdp
real(wp_), dimension(:), allocatable :: psrad, terad, derad, zfc
real(wp_), dimension(:), allocatable :: rv, zv, psinr, fpol, qpsi
real(wp_), dimension(:), allocatable :: rbnd, zbnd, rlim, zlim
real(wp_), dimension(:,:), allocatable :: psin
real(wp_) :: psia, rvac, rax, zax
integer :: iox0
real(wp_) :: p0mw, fghz, psipol0, chipol0
real(wp_) :: alpha0, beta0, x0, y0, z0, w1, w2, ri1, ri2, phiw, phir
real(wp_) :: pec,icd
integer :: ierr
real(wp_), dimension(:), allocatable :: xrad, rhot, dpdv, jcd
real(wp_) :: rwallm, rmxm, r0m, z0m, dzmx
! ======= read parameters BEGIN =======
call read_params('gray_params.data',rtrp,hcdp,antp,eqp,rwallm,prfp,outp)
! ======= read parameters END =======
! ======= read input data BEGIN =======
!------------ equilibrium ------------
if(eqp%iequil<2) then
call read_equil_an(eqp%filenm, rtrp%ipass, rv, zv, fpol, qpsi, rlim, zlim)
! psia sign set to give the correct sign to Iphi (COCOS=3: psia<0 for Iphi>0)
psia = sign(one,qpsi(2)*fpol(1))
else
call read_eqdsk(eqp%filenm, rv,zv,psin, psia, psinr,fpol,qpsi, rvac, &
rax,zax, rbnd,zbnd, rlim,zlim, eqp%ipsinorm,eqp%idesc,eqp%ifreefmt)
call change_cocos(psia, fpol, qpsi, eqp%icocos, 3)
end if
! re-scale B/I and/or force signs. If sgn=0 on input, set to fpol/-psia signs on output
call eq_scal(psia, fpol, eqp%sgni, eqp%sgnb, eqp%factb)
! ??? analytical only? change for numerical!
! qpsi(1) = sign(qpsi(1),qpsi(1)*qpsi(2)*psia*fpol(1))
! qpsi(2) = sign(qpsi(2),psia*fpol(1))
!------------- profiles -------------
if(prfp%iprof==0) then
call read_profiles_an(prfp%filenm, terad, derad, zfc)
else
call read_profiles(prfp%filenm, xrad, terad, derad, zfc)
allocate(psrad(size(xrad)))
if(prfp%irho==0) then ! xrad==rhot
allocate(rhot(size(psinr)))
call setqphi_num(psinr,qpsi,psia,rhot)
call set_rhospl(sqrt(psinr),rhot)
deallocate(rhot)
psrad=frhopolv(xrad)**2
else if(prfp%irho == 1) then ! xrad==rhop
psrad=xrad**2
else
psrad=xrad
end if
deallocate(xrad)
end if
! re-scale input data
call tene_scal(terad,derad,prfp%factte,prfp%factne,eqp%factb,prfp%iscal, &
prfp%iprof)
!------------- antenna --------------
! interpolate beam table if antctrl%ibeam>0
select case (antp%ibeam)
case (2)
! to be completed: now 1st beamd always selected, iox read from table
call read_beam2(antp%filenm,1,antp%alpha,antp%beta,fghz,antp%iox,x0,y0,z0, &
w1,w2,ri1,ri2,phiw,phir)
case (1)
call read_beam1(antp%filenm,antp%alpha,antp%beta,fghz,x0,y0,z0, &
w1,w2,ri1,ri2,phiw,phir)
case default
call read_beam0(antp%filenm,fghz,x0,y0,z0,w1,w2,ri1,ri2,phiw,phir)
end select
alpha0=antp%alpha
beta0=antp%beta
p0mw=antp%power
psipol0=antp%psi
chipol0=antp%chi
iox0=antp%iox
!--------------- wall ---------------
! set simple limiter if not read from EQDSK
! need to clean up...
r0m=sqrt(x0**2+y0**2)*0.01_wp_
dzmx=abs(rtrp%ipass)*rtrp%dst*rtrp%nstep*0.01_wp_
z0m=z0*0.01_wp_
if (.not.allocated(rlim).or.rtrp%ipass<0) then
if (allocated(rlim)) deallocate(rlim)
if (allocated(zlim)) deallocate(zlim)
allocate(rlim(5))
allocate(zlim(5))
if (rtrp%ipass<0) rtrp%ipass = -rtrp%ipass
if(eqp%iequil<2) then
rmxm=(rv(1)+rv(2))*0.01_wp_
else
rmxm=rv(size(rv))
end if
call range2rect(rwallm,max(r0m,rmxm),z0m-dzmx,z0m+dzmx,rlim,zlim)
end if
! ======= read input data END =======
! ========================= MAIN SUBROUTINE CALL =========================
allocate(dpdv(outp%nrho),jcd(outp%nrho))
call gray_main(rv,zv,psin,psia,psinr,fpol,qpsi,rvac,rax,zax,rbnd,zbnd,eqp, &
psrad,terad,derad,zfc,prfp, rlim,zlim, &
p0mw,fghz,alpha0,beta0,(/x0,y0,z0/),w1,w2,ri1,ri2,phiw,phir,iox0, &
psipol0,chipol0, dpdv,jcd,pec,icd, outp,rtrp,hcdp,ierr)
! ========================================================================
! ======= control prints BEGIN =======
if(ierr/=0) print*,' IERR = ', ierr
print*,' '
print*,'Pabs (MW), Icd (kA) = ', pec,icd*1.0e3_wp_
! ======= control prints END =======
! ======= free memory BEGIN =======
if(allocated(psrad)) deallocate(psrad)
if(allocated(terad)) deallocate(terad, derad, zfc)
if(allocated(rv)) deallocate(rv, zv, fpol, qpsi)
if(allocated(psin)) deallocate(psin, psinr)
if(allocated(rbnd)) deallocate(rbnd,zbnd)
if(allocated(rlim)) deallocate(rlim,zlim)
if(allocated(dpdv)) deallocate(dpdv, jcd)
! ======= free memory END ======
end program main_std