merge graysum and gray
Note that graysum, now as sum_mode, still doesn't support multipass.
This commit is contained in:
parent
37ee881024
commit
693af2a763
6
Makefile
6
Makefile
@ -38,7 +38,7 @@ MODULES = $(filter-out $(OBJDIR)/main%,$(OBJECTS))
|
|||||||
# Build outputs
|
# Build outputs
|
||||||
GRAY = $(BINDIR)/gray
|
GRAY = $(BINDIR)/gray
|
||||||
LIBGRAY = $(LIBDIR)/libgray
|
LIBGRAY = $(LIBDIR)/libgray
|
||||||
BINARIES = $(GRAY) $(GRAY)sum
|
BINARIES = $(GRAY)
|
||||||
LIBRARIES = $(LIBGRAY).so
|
LIBRARIES = $(LIBGRAY).so
|
||||||
|
|
||||||
##
|
##
|
||||||
@ -89,10 +89,6 @@ install: $(BINARIES) $(LIBRARIES) doc/gray.1
|
|||||||
$(GRAY): $(shell ./depend $(OBJDIR)/main.o) | $(BINDIR)
|
$(GRAY): $(shell ./depend $(OBJDIR)/main.o) | $(BINDIR)
|
||||||
$(LD) $(LDFLAGS) -o $@ $^
|
$(LD) $(LDFLAGS) -o $@ $^
|
||||||
|
|
||||||
# sum variant
|
|
||||||
$(GRAY)sum: $(shell ./depend $(OBJDIR)/main-sum.o) | $(BINDIR)
|
|
||||||
$(LD) $(LDFLAGS) -o $@ $^
|
|
||||||
|
|
||||||
# GRAY shared library
|
# GRAY shared library
|
||||||
$(LIBGRAY).so: $(MODULES) | $(LIBDIR)
|
$(LIBGRAY).so: $(MODULES) | $(LIBDIR)
|
||||||
$(LD) -shared $(LDFLAGS) -o $@ $^
|
$(LD) -shared $(LDFLAGS) -o $@ $^
|
||||||
|
165
src/graycore.f90
165
src/graycore.f90
@ -603,6 +603,171 @@ contains
|
|||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
subroutine sum_profiles(rv,zv,psin,psia,psinr,fpol,qpsi,rvac,rax,zax,rbnd,zbnd, &
|
||||||
|
eqp,psrad,terad,derad,zfc,prfp, rlim,zlim, &
|
||||||
|
p0,fghz,alpha0,beta0,xv0,w1,w2,ri1,ri2,phiw,phir,iox0, &
|
||||||
|
psipol0,chipol0, jphi,jcd,dpdv,currins,pins,pabs,icd, &
|
||||||
|
jphip,dpdvp, &
|
||||||
|
rhotj,rhotjava,rhotp,rhotpav,drhotjava,drhotpav,ratjamx,ratjbmx, &
|
||||||
|
outp,rtrp,hcdp,ierr, rhout)
|
||||||
|
use const_and_precisions, only : zero, one, degree
|
||||||
|
use coreprofiles, only : set_prfan, set_prfspl, temp, fzeff, unset_prfspl
|
||||||
|
use dispersion, only : expinit
|
||||||
|
use gray_params, only : eqparam_type, prfparam_type, outparam_type, &
|
||||||
|
rtrparam_type, hcdparam_type, antctrl_type, set_codepar, print_params, &
|
||||||
|
iequil, iprof, iwarm, ipec, istpr0, igrad, headw, headl
|
||||||
|
use beams, only : read_beam0, read_beam1, launchangles2n, xgygcoeff
|
||||||
|
use beamdata, only : pweight, rayi2jk
|
||||||
|
use equilibrium, only : set_equian, set_eqspl, setqphi_num, set_rhospl, &
|
||||||
|
zbinf, zbsup, unset_eqspl, unset_rhospl, unset_q
|
||||||
|
use errcodes, only : check_err, print_errn, print_errhcd
|
||||||
|
use magsurf_data, only : flux_average, dealloc_surfvec
|
||||||
|
use beamdata, only : init_btr, dealloc_beam, nray, nstep, dst
|
||||||
|
use pec, only : pec_init, spec, postproc_profiles, dealloc_pec, &
|
||||||
|
rhop_tab, rhot_tab
|
||||||
|
use limiter, only : set_lim, unset_lim
|
||||||
|
use utils, only : vmaxmin
|
||||||
|
implicit none
|
||||||
|
! arguments
|
||||||
|
type(eqparam_type), intent(in) :: eqp
|
||||||
|
type(prfparam_type), intent(in) :: prfp
|
||||||
|
type(outparam_type), intent(in) :: outp
|
||||||
|
type(rtrparam_type), intent(in) :: rtrp
|
||||||
|
type(hcdparam_type), intent(in) :: hcdp
|
||||||
|
|
||||||
|
real(wp_), dimension(:), intent(in) :: psrad, terad, derad, zfc
|
||||||
|
real(wp_), dimension(:), intent(in) :: rv, zv, psinr, fpol, qpsi
|
||||||
|
real(wp_), dimension(:), intent(in) :: rbnd, zbnd, rlim, zlim
|
||||||
|
real(wp_), dimension(:,:), intent(in) :: psin
|
||||||
|
real(wp_), intent(in) :: psia, rvac, rax, zax
|
||||||
|
integer, intent(in) :: iox0
|
||||||
|
real(wp_), intent(in) :: p0, fghz, psipol0, chipol0
|
||||||
|
real(wp_), intent(in) :: alpha0,beta0, w1,w2, ri1,ri2, phiw,phir
|
||||||
|
real(wp_), dimension(3), intent(in) :: xv0
|
||||||
|
|
||||||
|
real(wp_), intent(in) :: pabs,icd
|
||||||
|
real(wp_), dimension(:), intent(in) :: jphi,jcd,dpdv,currins,pins
|
||||||
|
real(wp_), intent(out) :: jphip,dpdvp, &
|
||||||
|
rhotj,rhotjava,rhotp,rhotpav,drhotjava,drhotpav,ratjamx,ratjbmx
|
||||||
|
|
||||||
|
real(wp_), dimension(:), intent(in), optional :: rhout
|
||||||
|
|
||||||
|
integer, intent(out) :: ierr
|
||||||
|
! local variables
|
||||||
|
real(wp_), parameter :: taucr = 12._wp_
|
||||||
|
|
||||||
|
real(wp_) :: sox,ak0,bres,xgcn,xg,yg,zzm,alpha,didp,anpl,anpr,anprim,anprre
|
||||||
|
real(wp_) :: chipol,psipol,btot,psinv,dens,tekev,dersdst,derdnm,st,st0
|
||||||
|
real(wp_) :: tau,pow,dids,ddr,ddi,taumn,taumx
|
||||||
|
real(wp_) :: drhotp,drhotj,dpdvmx,jphimx
|
||||||
|
|
||||||
|
real(wp_), dimension(3) :: xv,anv0,anv,bv
|
||||||
|
real(wp_), dimension(:,:), pointer :: yw=>null(),ypw=>null(),gri=>null()
|
||||||
|
real(wp_), dimension(:,:,:), pointer :: xc=>null(),du1=>null(),ggri=>null()
|
||||||
|
integer :: i,jk,iox,nharm,nhf,nnd,iokhawa,istop,ierrn,ierrhcd,index_rt=1
|
||||||
|
logical :: ins_pl, somein, allout
|
||||||
|
|
||||||
|
real(wp_), dimension(:,:), pointer :: psjki=>null(),ppabs=>null(),ccci=>null()
|
||||||
|
real(wp_), dimension(:), pointer :: tau0=>null(),alphaabs0=>null(),dids0=>null(), &
|
||||||
|
ccci0=>null()
|
||||||
|
real(wp_), dimension(:), pointer :: p0jk=>null()
|
||||||
|
complex(wp_), dimension(:), pointer :: ext=>null(), eyt=>null()
|
||||||
|
integer, dimension(:), pointer :: iiv=>null()
|
||||||
|
|
||||||
|
! parameters log in file headers
|
||||||
|
character(len=headw), dimension(headl) :: strheader
|
||||||
|
type(antctrl_type) :: antp
|
||||||
|
real(wp_) :: rwall
|
||||||
|
|
||||||
|
! ======== set environment BEGIN ========
|
||||||
|
call set_codepar(eqp,prfp,outp,rtrp,hcdp)
|
||||||
|
|
||||||
|
call set_lim(rlim,zlim)
|
||||||
|
|
||||||
|
if(iequil<2) then
|
||||||
|
call set_equian(rv(1),zv(1),rv(2), fpol(1)/rv(1), qpsi(1),qpsi(2),qpsi(3))
|
||||||
|
else
|
||||||
|
call set_eqspl(rv,zv,psin, psia, psinr,fpol, qpsi, eqp%ssplps,eqp%ssplf, &
|
||||||
|
rvac, rax,zax, rbnd,zbnd, eqp%ixp)
|
||||||
|
! qpsi used for rho_pol/rho_tor mapping (initializes fq,frhotor,frhopol)
|
||||||
|
end if
|
||||||
|
! compute flux surface averaged quantities
|
||||||
|
call flux_average ! requires frhotor for dadrhot,dvdrhot
|
||||||
|
|
||||||
|
if(iprof==0) then
|
||||||
|
call set_prfan(terad,derad,zfc)
|
||||||
|
else
|
||||||
|
call set_prfspl(psrad, terad, derad, zfc, prfp%sspld, prfp%psnbnd)
|
||||||
|
end if
|
||||||
|
|
||||||
|
call xgygcoeff(fghz,ak0,bres,xgcn)
|
||||||
|
call launchangles2n(alpha0,beta0,xv0,anv0)
|
||||||
|
call init_btr(rtrp,yw,ypw,xc,du1,gri,ggri,psjki,ppabs,ccci, &
|
||||||
|
tau0,alphaabs0,dids0,ccci0,p0jk,ext,eyt,iiv)
|
||||||
|
|
||||||
|
if(iwarm > 1) call expinit
|
||||||
|
|
||||||
|
! ======= set environment END ======
|
||||||
|
|
||||||
|
! ======= pre-proc prints BEGIN ======
|
||||||
|
antp%alpha=alpha0
|
||||||
|
antp%beta=beta0
|
||||||
|
antp%power=p0
|
||||||
|
antp%psi=psipol0
|
||||||
|
antp%chi=chipol0
|
||||||
|
antp%iox=iox0
|
||||||
|
!!!!! missing values
|
||||||
|
antp%ibeam=0
|
||||||
|
antp%filenm=''
|
||||||
|
rwall=0._wp_
|
||||||
|
call print_params(rtrp,hcdp,antp,eqp,rwall,prfp,outp,strheader)
|
||||||
|
call print_headers(strheader, 0)
|
||||||
|
! print psi surface for q=1.5 and q=2 on file and psi,rhot,rhop on stdout
|
||||||
|
call print_surfq((/1.5_wp_,2.0_wp_/))
|
||||||
|
! print
|
||||||
|
print*,' '
|
||||||
|
print'(a,2f8.3)','alpha0, beta0 = ',alpha0,beta0
|
||||||
|
print'(a,4f8.3)','x00, y00, z00 = ',xv0
|
||||||
|
! print Btot=Bres
|
||||||
|
! print ne, Te, q, Jphi versus psi, rhop, rhot
|
||||||
|
call print_bres(bres)
|
||||||
|
call print_prof
|
||||||
|
call print_maps(bres,xgcn,0.01_wp_*sqrt(xv0(1)**2+xv0(2)**2), &
|
||||||
|
sin(beta0*degree))
|
||||||
|
! ======= pre-proc prints END ======
|
||||||
|
|
||||||
|
! ======= post-proc BEGIN ======
|
||||||
|
|
||||||
|
! compute power and current density profiles for all rays
|
||||||
|
call pec_init(ipec,rhout)
|
||||||
|
nnd=size(rhop_tab)
|
||||||
|
! print power and current density profiles
|
||||||
|
call print_pec(rhop_tab,rhot_tab,jphi,jcd,dpdv,currins,pins,index_rt)
|
||||||
|
! compute profiles width
|
||||||
|
call postproc_profiles(pabs,icd,rhot_tab,dpdv,jphi, &
|
||||||
|
rhotpav,drhotpav,rhotjava,drhotjava,dpdvp,jphip, &
|
||||||
|
rhotp,drhotp,rhotj,drhotj,dpdvmx,jphimx,ratjamx,ratjbmx)
|
||||||
|
! print 0D results
|
||||||
|
call print_finals(pabs,icd,dpdvp,jphip,rhotpav,rhotjava,drhotpav, &
|
||||||
|
drhotjava,dpdvmx,jphimx,rhotp,rhotj,drhotp,drhotj,ratjamx,ratjbmx, &
|
||||||
|
st,psipol,chipol,index_rt,p0,zero,zero) ! cpl1=cpl2=0 for a single pass
|
||||||
|
|
||||||
|
! ======= post-proc END ======
|
||||||
|
|
||||||
|
! ======= free memory BEGIN ======
|
||||||
|
call unset_eqspl
|
||||||
|
call unset_q
|
||||||
|
call unset_rhospl
|
||||||
|
call unset_prfspl
|
||||||
|
call unset_lim
|
||||||
|
call dealloc_surfvec
|
||||||
|
call dealloc_beam(yw,ypw,xc,du1,gri,ggri,psjki,ppabs,ccci, &
|
||||||
|
tau0,alphaabs0,dids0,ccci0,p0jk,ext,eyt,iiv)
|
||||||
|
call dealloc_pec
|
||||||
|
! ======= free memory END ======
|
||||||
|
end subroutine sum_profiles
|
||||||
|
|
||||||
|
|
||||||
subroutine vectinit(psjki,ppabs,ccci,tau0,alphaabs0,dids0,ccci0,iiv)
|
subroutine vectinit(psjki,ppabs,ccci,tau0,alphaabs0,dids0,ccci0,iiv)
|
||||||
use const_and_precisions, only : wp_, zero
|
use const_and_precisions, only : wp_, zero
|
||||||
implicit none
|
implicit none
|
||||||
|
195
src/main-sum.f90
195
src/main-sum.f90
@ -1,195 +0,0 @@
|
|||||||
program main_sum
|
|
||||||
use const_and_precisions, only : wp_,one
|
|
||||||
use sumcore, only : sum_profiles
|
|
||||||
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
|
|
||||||
real(wp_), dimension(:), allocatable :: jphi,jcd,dpdv,currins,pins,rtin,rpin
|
|
||||||
real(wp_) :: rwallm, rmxm, r0m, z0m, dzmx
|
|
||||||
|
|
||||||
integer :: i,j,k,n,ngam,irt
|
|
||||||
character(len=255) :: fn,dumstr
|
|
||||||
real(wp_), dimension(5) :: f48v
|
|
||||||
real(wp_) :: gam,alp,bet, jphip,dpdvp, &
|
|
||||||
rhotj,rhotjava,rhotp,rhotpav,drhotjava,drhotpav,ratjamx,ratjbmx
|
|
||||||
|
|
||||||
|
|
||||||
! ======= 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(jphi(outp%nrho),jcd(outp%nrho),dpdv(outp%nrho), &
|
|
||||||
currins(outp%nrho),pins(outp%nrho),rtin(outp%nrho),rpin(outp%nrho))
|
|
||||||
|
|
||||||
open(100,file='filelist.txt',action='read',status='old')
|
|
||||||
read(100,*) n,ngam
|
|
||||||
do i=1,n
|
|
||||||
read(100,*) fn
|
|
||||||
open(100+i,file=fn,action='read',status='old')
|
|
||||||
do j=1,22
|
|
||||||
read(100+i,*) dumstr
|
|
||||||
end do
|
|
||||||
end do
|
|
||||||
close(100)
|
|
||||||
|
|
||||||
open(100+n+1,file='f48sum.txt',action='write',status='unknown')
|
|
||||||
open(100+n+2,file='f7sum.txt',action='write',status='unknown')
|
|
||||||
|
|
||||||
do k=1,ngam
|
|
||||||
jphi=0.0_wp_
|
|
||||||
jcd=0.0_wp_
|
|
||||||
dpdv=0.0_wp_
|
|
||||||
currins=0.0_wp_
|
|
||||||
pins=0.0_wp_
|
|
||||||
do j=1,outp%nrho
|
|
||||||
do i=1,n
|
|
||||||
read(100+i,*) gam,alp,bet,rpin(j),rtin(j),f48v(1:5),irt
|
|
||||||
jphi(j)=jphi(j)+f48v(1)
|
|
||||||
jcd(j)=jcd(j)+f48v(2)
|
|
||||||
dpdv(j)=dpdv(j)+f48v(3)
|
|
||||||
currins(j)=currins(j)+f48v(4)
|
|
||||||
pins(j)=pins(j)+f48v(5)
|
|
||||||
end do
|
|
||||||
write(100+n+1,'(10(1x,e16.8e3),i5)') gam,alp,bet,rpin(j),rtin(j), &
|
|
||||||
jphi(j),jcd(j),dpdv(j),currins(j),pins(j),irt
|
|
||||||
end do
|
|
||||||
pec=pins(outp%nrho)
|
|
||||||
icd=currins(outp%nrho)
|
|
||||||
write(100+n+1,*)
|
|
||||||
call sum_profiles(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, jphi,jcd,dpdv,currins,pins,pec,icd, &
|
|
||||||
jphip,dpdvp, &
|
|
||||||
rhotj,rhotjava,rhotp,rhotpav,drhotjava,drhotpav,ratjamx,ratjbmx, &
|
|
||||||
outp,rtrp,hcdp,ierr)
|
|
||||||
write(100+n+2,'(15(1x,e12.5),i5,4(1x,e12.5))') gam,alp,bet,icd,pec, &
|
|
||||||
jphip,dpdvp, &
|
|
||||||
rhotj,rhotjava,rhotp,rhotpav,drhotjava,drhotpav,ratjamx,ratjbmx
|
|
||||||
end do
|
|
||||||
do i=1,n+2
|
|
||||||
close(100+i)
|
|
||||||
end do
|
|
||||||
! ========================================================================
|
|
||||||
|
|
||||||
! ======= control prints BEGIN =======
|
|
||||||
if(ierr/=0) print*,' IERR = ', ierr
|
|
||||||
print*,' '
|
|
||||||
print*,'Pabs (MW) = ', pec
|
|
||||||
print*,'Icd (kA) = ', 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(jphi,jcd,dpdv,currins,pins,rtin,rpin)
|
|
||||||
! ======= free memory END ======
|
|
||||||
end program main_sum
|
|
124
src/main.f90
124
src/main.f90
@ -1,13 +1,15 @@
|
|||||||
program main_std
|
program main_std
|
||||||
use const_and_precisions, only : wp_,one
|
use const_and_precisions, only : wp_,one
|
||||||
use graycore, only : gray_main
|
use graycore, only : gray_main, sum_profiles
|
||||||
use gray_params, only : read_params, antctrl_type, eqparam_type, &
|
use gray_params, only : read_params, antctrl_type, eqparam_type, &
|
||||||
prfparam_type,outparam_type,rtrparam_type,hcdparam_type
|
prfparam_type, outparam_type, rtrparam_type, &
|
||||||
|
hcdparam_type
|
||||||
use beams, only : read_beam0, read_beam1, read_beam2
|
use beams, only : read_beam0, read_beam1, read_beam2
|
||||||
use equilibrium, only : read_equil_an,read_eqdsk,change_cocos,eq_scal, &
|
use equilibrium, only : read_equil_an, read_eqdsk, change_cocos, &
|
||||||
set_rhospl,setqphi_num,frhopolv
|
eq_scal, set_rhospl, setqphi_num, frhopolv
|
||||||
use coreprofiles, only : read_profiles_an, read_profiles, tene_scal
|
use coreprofiles, only : read_profiles_an, read_profiles, tene_scal
|
||||||
use reflections, only : range2rect
|
use reflections, only : range2rect
|
||||||
|
|
||||||
implicit none
|
implicit none
|
||||||
type(antctrl_type) :: antp
|
type(antctrl_type) :: antp
|
||||||
type(eqparam_type) :: eqp
|
type(eqparam_type) :: eqp
|
||||||
@ -21,9 +23,8 @@ program main_std
|
|||||||
real(wp_), dimension(:), allocatable :: rbnd, zbnd, rlim, zlim
|
real(wp_), dimension(:), allocatable :: rbnd, zbnd, rlim, zlim
|
||||||
real(wp_), dimension(:,:), allocatable :: psin
|
real(wp_), dimension(:,:), allocatable :: psin
|
||||||
real(wp_) :: psia, rvac, rax, zax
|
real(wp_) :: psia, rvac, rax, zax
|
||||||
integer :: iox0
|
real(wp_) :: fghz
|
||||||
real(wp_) :: p0mw, fghz, psipol0, chipol0
|
real(wp_) :: x0, y0, z0, w1, w2, ri1, ri2, phiw, phir
|
||||||
real(wp_) :: alpha0, beta0, x0, y0, z0, w1, w2, ri1, ri2, phiw, phir
|
|
||||||
|
|
||||||
real(wp_) :: pec,icd
|
real(wp_) :: pec,icd
|
||||||
|
|
||||||
@ -31,11 +32,21 @@ program main_std
|
|||||||
real(wp_), dimension(:), allocatable :: xrad, rhot, dpdv, jcd
|
real(wp_), dimension(:), allocatable :: xrad, rhot, dpdv, jcd
|
||||||
real(wp_) :: rwallm, rmxm, r0m, z0m, dzmx
|
real(wp_) :: rwallm, rmxm, r0m, z0m, dzmx
|
||||||
|
|
||||||
! ======= read parameters BEGIN =======
|
logical :: sum_mode = .FALSE.
|
||||||
|
|
||||||
|
! ------- sum mode variables -------
|
||||||
|
real(wp_), dimension(:), allocatable :: jphi, currins, pins, rtin, rpin
|
||||||
|
integer :: i,j,k,n,ngam,irt
|
||||||
|
character(len=255) :: fn,dumstr
|
||||||
|
real(wp_), dimension(5) :: f48v
|
||||||
|
real(wp_) :: gam,alp,bet, jphip,dpdvp, &
|
||||||
|
rhotj,rhotjava,rhotp,rhotpav,drhotjava,drhotpav,ratjamx,ratjbmx
|
||||||
|
! ----------------------------------
|
||||||
|
|
||||||
call read_params('gray_params.data', rtrp, hcdp, antp, eqp, rwallm, prfp, outp)
|
call read_params('gray_params.data', rtrp, hcdp, antp, eqp, rwallm, prfp, outp)
|
||||||
! ======= read parameters END =======
|
|
||||||
|
|
||||||
! ======= read input data BEGIN =======
|
! ======= read input data BEGIN =======
|
||||||
|
|
||||||
!------------ equilibrium ------------
|
!------------ equilibrium ------------
|
||||||
if(eqp%iequil<2) then
|
if(eqp%iequil<2) then
|
||||||
call read_equil_an(eqp%filenm, rtrp%ipass, rv, zv, fpol, qpsi, rlim, zlim)
|
call read_equil_an(eqp%filenm, rtrp%ipass, rv, zv, fpol, qpsi, rlim, zlim)
|
||||||
@ -43,7 +54,8 @@ program main_std
|
|||||||
psia = sign(one,qpsi(2)*fpol(1))
|
psia = sign(one,qpsi(2)*fpol(1))
|
||||||
else
|
else
|
||||||
call read_eqdsk(eqp%filenm, rv,zv,psin, psia, psinr,fpol,qpsi, rvac, &
|
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)
|
rax,zax, rbnd,zbnd, rlim,zlim, &
|
||||||
|
eqp%ipsinorm,eqp%idesc,eqp%ifreefmt)
|
||||||
call change_cocos(psia, fpol, qpsi, eqp%icocos, 3)
|
call change_cocos(psia, fpol, qpsi, eqp%icocos, 3)
|
||||||
end if
|
end if
|
||||||
! re-scale B/I and/or force signs. If sgn=0 on input, set to fpol/-psia signs on output
|
! re-scale B/I and/or force signs. If sgn=0 on input, set to fpol/-psia signs on output
|
||||||
@ -51,6 +63,8 @@ program main_std
|
|||||||
! ??? analytical only? change for numerical!
|
! ??? analytical only? change for numerical!
|
||||||
! qpsi(1) = sign(qpsi(1),qpsi(1)*qpsi(2)*psia*fpol(1))
|
! qpsi(1) = sign(qpsi(1),qpsi(1)*qpsi(2)*psia*fpol(1))
|
||||||
! qpsi(2) = sign(qpsi(2),psia*fpol(1))
|
! qpsi(2) = sign(qpsi(2),psia*fpol(1))
|
||||||
|
! ----------------------------------
|
||||||
|
|
||||||
!------------- profiles -------------
|
!------------- profiles -------------
|
||||||
if(prfp%iprof == 0) then
|
if(prfp%iprof == 0) then
|
||||||
call read_profiles_an(prfp%filenm, terad, derad, zfc)
|
call read_profiles_an(prfp%filenm, terad, derad, zfc)
|
||||||
@ -71,27 +85,29 @@ program main_std
|
|||||||
deallocate(xrad)
|
deallocate(xrad)
|
||||||
end if
|
end if
|
||||||
! re-scale input data
|
! re-scale input data
|
||||||
call tene_scal(terad,derad,prfp%factte,prfp%factne,eqp%factb,prfp%iscal, &
|
call tene_scal(terad, derad, prfp%factte, prfp%factne, &
|
||||||
prfp%iprof)
|
eqp%factb, prfp%iscal, prfp%iprof)
|
||||||
|
! ----------------------------------
|
||||||
|
|
||||||
!------------- antenna --------------
|
!------------- antenna --------------
|
||||||
! interpolate beam table if antctrl%ibeam>0
|
! interpolate beam table if antctrl%ibeam>0
|
||||||
select case (antp%ibeam)
|
select case (antp%ibeam)
|
||||||
case (2)
|
case (2)
|
||||||
! to be completed: now 1st beamd always selected, iox read from table
|
! 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, &
|
call read_beam2(antp%filenm, 1, antp%alpha, antp%beta, &
|
||||||
|
fghz, antp%iox, x0, y0, z0, &
|
||||||
w1, w2, ri1, ri2, phiw, phir)
|
w1, w2, ri1, ri2, phiw, phir)
|
||||||
case (1)
|
case (1)
|
||||||
call read_beam1(antp%filenm,antp%alpha,antp%beta,fghz,x0,y0,z0, &
|
call read_beam1(antp%filenm, antp%alpha, antp%beta, &
|
||||||
|
fghz, x0, y0, z0, &
|
||||||
w1, w2, ri1, ri2, phiw, phir)
|
w1, w2, ri1, ri2, phiw, phir)
|
||||||
case default
|
case default
|
||||||
call read_beam0(antp%filenm,fghz,x0,y0,z0,w1,w2,ri1,ri2,phiw,phir)
|
call read_beam0(antp%filenm, &
|
||||||
|
fghz, x0, y0, z0, &
|
||||||
|
w1, w2, ri1, ri2, phiw, phir)
|
||||||
end select
|
end select
|
||||||
alpha0=antp%alpha
|
! ----------------------------------
|
||||||
beta0=antp%beta
|
|
||||||
p0mw=antp%power
|
|
||||||
psipol0=antp%psi
|
|
||||||
chipol0=antp%chi
|
|
||||||
iox0=antp%iox
|
|
||||||
!--------------- wall ---------------
|
!--------------- wall ---------------
|
||||||
! set simple limiter if not read from EQDSK
|
! set simple limiter if not read from EQDSK
|
||||||
! need to clean up...
|
! need to clean up...
|
||||||
@ -111,20 +127,79 @@ program main_std
|
|||||||
end if
|
end if
|
||||||
call range2rect(rwallm,max(r0m,rmxm),z0m-dzmx,z0m+dzmx,rlim,zlim)
|
call range2rect(rwallm,max(r0m,rmxm),z0m-dzmx,z0m+dzmx,rlim,zlim)
|
||||||
end if
|
end if
|
||||||
|
! ----------------------------------
|
||||||
|
|
||||||
! ======= read input data END =======
|
! ======= read input data END =======
|
||||||
|
|
||||||
! ========================= MAIN SUBROUTINE CALL =========================
|
! ========================= MAIN SUBROUTINE CALL =========================
|
||||||
allocate(dpdv(outp%nrho), jcd(outp%nrho))
|
allocate(dpdv(outp%nrho), jcd(outp%nrho))
|
||||||
|
if (sum_mode) then
|
||||||
|
allocate(jphi(outp%nrho), currins(outp%nrho), &
|
||||||
|
pins(outp%nrho), rtin(outp%nrho), rpin(outp%nrho))
|
||||||
|
|
||||||
|
open(100,file='filelist.txt',action='read',status='old')
|
||||||
|
read(100,*) n,ngam
|
||||||
|
do i=1,n
|
||||||
|
read(100,*) fn
|
||||||
|
open(100+i,file=fn,action='read',status='old')
|
||||||
|
do j=1,22
|
||||||
|
read(100+i,*) dumstr
|
||||||
|
end do
|
||||||
|
end do
|
||||||
|
close(100)
|
||||||
|
|
||||||
|
open(100+n+1,file='f48sum.txt',action='write',status='unknown')
|
||||||
|
open(100+n+2,file='f7sum.txt',action='write',status='unknown')
|
||||||
|
|
||||||
|
do k=1,ngam
|
||||||
|
jphi=0.0_wp_
|
||||||
|
jcd=0.0_wp_
|
||||||
|
dpdv=0.0_wp_
|
||||||
|
currins=0.0_wp_
|
||||||
|
pins=0.0_wp_
|
||||||
|
do j=1,outp%nrho
|
||||||
|
do i=1,n
|
||||||
|
read(100+i,*) gam,alp,bet,rpin(j),rtin(j),f48v(1:5),irt
|
||||||
|
jphi(j)=jphi(j)+f48v(1)
|
||||||
|
jcd(j)=jcd(j)+f48v(2)
|
||||||
|
dpdv(j)=dpdv(j)+f48v(3)
|
||||||
|
currins(j)=currins(j)+f48v(4)
|
||||||
|
pins(j)=pins(j)+f48v(5)
|
||||||
|
end do
|
||||||
|
write(100+n+1,'(10(1x,e16.8e3),i5)') gam,alp,bet,rpin(j),rtin(j), &
|
||||||
|
jphi(j),jcd(j),dpdv(j),currins(j),pins(j),irt
|
||||||
|
end do
|
||||||
|
pec=pins(outp%nrho)
|
||||||
|
icd=currins(outp%nrho)
|
||||||
|
write(100+n+1,*)
|
||||||
|
call sum_profiles(rv,zv,psin,psia,psinr,fpol,qpsi,rvac,rax,zax,rbnd,zbnd,eqp, &
|
||||||
|
psrad,terad,derad,zfc,prfp, rlim,zlim, &
|
||||||
|
antp%beta,fghz,antp%alpha,antp%beta, &
|
||||||
|
(/x0,y0,z0/),w1,w2,ri1,ri2,phiw,phir,antp%iox, &
|
||||||
|
antp%psi, antp%chi, jphi,jcd,dpdv,currins,pins,pec,icd, &
|
||||||
|
jphip,dpdvp, rhotj,rhotjava,rhotp,rhotpav, &
|
||||||
|
drhotjava,drhotpav, ratjamx,ratjbmx, outp,rtrp,hcdp,ierr)
|
||||||
|
write(100+n+2,'(15(1x,e12.5),i5,4(1x,e12.5))') gam,alp,bet,icd,pec, &
|
||||||
|
jphip,dpdvp, &
|
||||||
|
rhotj,rhotjava,rhotp,rhotpav,drhotjava,drhotpav,ratjamx,ratjbmx
|
||||||
|
end do
|
||||||
|
do i=1,n+2
|
||||||
|
close(100+i)
|
||||||
|
end do
|
||||||
|
else
|
||||||
call gray_main(rv,zv,psin,psia,psinr,fpol,qpsi,rvac,rax,zax,rbnd,zbnd,eqp, &
|
call gray_main(rv,zv,psin,psia,psinr,fpol,qpsi,rvac,rax,zax,rbnd,zbnd,eqp, &
|
||||||
psrad,terad,derad,zfc,prfp, rlim,zlim, &
|
psrad,terad,derad,zfc,prfp, rlim,zlim, &
|
||||||
p0mw,fghz,alpha0,beta0,(/x0,y0,z0/),w1,w2,ri1,ri2,phiw,phir,iox0, &
|
antp%power,fghz,antp%alpha,antp%beta, &
|
||||||
psipol0,chipol0, dpdv,jcd,pec,icd, outp,rtrp,hcdp,ierr)
|
(/x0,y0,z0/),w1,w2,ri1,ri2,phiw,phir,antp%iox, &
|
||||||
|
antp%psi,antp%chi, dpdv,jcd,pec,icd, outp,rtrp,hcdp,ierr)
|
||||||
|
end if
|
||||||
! ========================================================================
|
! ========================================================================
|
||||||
|
|
||||||
! ======= control prints BEGIN =======
|
! ======= control prints BEGIN =======
|
||||||
if(ierr/=0) print*,' IERR = ', ierr
|
if(ierr/=0) print*,' IERR = ', ierr
|
||||||
print*,' '
|
print*,' '
|
||||||
print*,'Pabs (MW), Icd (kA) = ', pec,icd*1.0e3_wp_
|
print*,'Pabs (MW) = ', pec
|
||||||
|
print*,'Icd (kA) = ', icd*1.0e3_wp_
|
||||||
! ======= control prints END =======
|
! ======= control prints END =======
|
||||||
|
|
||||||
! ======= free memory BEGIN =======
|
! ======= free memory BEGIN =======
|
||||||
@ -135,5 +210,6 @@ program main_std
|
|||||||
if(allocated(rbnd)) deallocate(rbnd, zbnd)
|
if(allocated(rbnd)) deallocate(rbnd, zbnd)
|
||||||
if(allocated(rlim)) deallocate(rlim, zlim)
|
if(allocated(rlim)) deallocate(rlim, zlim)
|
||||||
if(allocated(dpdv)) deallocate(dpdv, jcd)
|
if(allocated(dpdv)) deallocate(dpdv, jcd)
|
||||||
|
if(sum_mode) deallocate(jphi, currins, pins, rtin, rpin)
|
||||||
! ======= free memory END ======
|
! ======= free memory END ======
|
||||||
end program main_std
|
end program main_std
|
1907
src/sumcore.f90
1907
src/sumcore.f90
File diff suppressed because it is too large
Load Diff
Loading…
Reference in New Issue
Block a user