gray/src/magsurf_data.f90
Michele Guerini Rocco 91a2e6cf07
src: implement toggling of output units
This actually implements the --units option
2022-05-11 01:15:04 +02:00

592 lines
18 KiB
Fortran

module magsurf_data
use const_and_precisions, only : wp_
implicit none
integer, save :: npsi, npoints !# sup mag, # punti per sup
integer, save :: njpt, nlmt
real(wp_), save :: rarea
real(wp_), dimension(:), allocatable, save :: psicon,pstab,rhot_eq, &
rhotqv,bav,varea,vcurrp,vajphiav,qqv,ffc,vratja,vratjb
real(wp_), dimension(:), allocatable, save :: rpstab
real(wp_), dimension(:), allocatable, save :: vvol,rri,rbav,bmxpsi,bmnpsi
real(wp_), dimension(:), allocatable, save :: tjp,tlm,ch,ch01
real(wp_), dimension(:,:), allocatable, save :: rcon,zcon
real(wp_), dimension(:,:), allocatable, save :: cdadrhot,cdvdrhot
real(wp_), dimension(:,:), allocatable, save :: cvol,crri,crbav,cbmx,cbmn,carea,cfc
real(wp_), dimension(:,:), allocatable, save :: crhotq
real(wp_), dimension(:,:), allocatable, save :: cratja,cratjb,cratjpl
contains
subroutine alloc_cnt(ierr)
implicit none
integer, intent(out) :: ierr
if(npsi.le.0.or.npoints.le.0) then
ierr = -1
return
end if
call dealloc_cnt
allocate(psicon(npsi),rcon(npoints,npsi),zcon(npoints,npsi))
end subroutine alloc_cnt
subroutine dealloc_cnt
implicit none
if(allocated(psicon)) deallocate(psicon)
if(allocated(rcon)) deallocate(rcon)
if(allocated(zcon)) deallocate(zcon)
end subroutine dealloc_cnt
subroutine alloc_surfvec(ierr)
implicit none
integer, intent(out) :: ierr
if(npsi.le.0.or.npoints.le.0) then
ierr = -1
return
end if
call dealloc_surfvec
call alloc_cnt(ierr)
allocate(pstab(npsi), &
rhot_eq(npsi),rhotqv(npsi),bav(npsi),bmxpsi(npsi),bmnpsi(npsi),varea(npsi), &
vvol(npsi),vcurrp(npsi),vajphiav(npsi),qqv(npsi),ffc(npsi),vratja(npsi), &
vratjb(npsi),rpstab(npsi),rri(npsi),rbav(npsi),cdadrhot(npsi,4), &
cdvdrhot(npsi,4),cbmx(npsi,4),cbmn(npsi,4),crbav(npsi,4),cvol(npsi,4), &
crri(npsi,4),carea(npsi,4),cfc(npsi,4),crhotq(npsi,4),cratjpl(npsi,4), &
cratja(npsi,4),cratjb(npsi,4))
end subroutine alloc_surfvec
subroutine dealloc_surfvec
implicit none
call dealloc_cnt
if(allocated(pstab)) deallocate(pstab)
if(allocated(rhot_eq)) deallocate(rhot_eq)
if(allocated(rhotqv)) deallocate(rhotqv)
if(allocated(bav)) deallocate(bav)
if(allocated(bmxpsi)) deallocate(bmxpsi)
if(allocated(bmnpsi)) deallocate(bmnpsi)
if(allocated(varea)) deallocate(varea)
if(allocated(vvol)) deallocate(vvol)
if(allocated(vcurrp)) deallocate(vcurrp)
if(allocated(vajphiav)) deallocate(vajphiav)
if(allocated(qqv)) deallocate(qqv)
if(allocated(ffc)) deallocate(ffc)
if(allocated(vratja)) deallocate(vratja)
if(allocated(vratjb)) deallocate(vratjb)
if(allocated(rpstab)) deallocate(rpstab)
if(allocated(rri)) deallocate(rri)
if(allocated(rbav)) deallocate(rbav)
if(allocated(cdadrhot)) deallocate(cdadrhot)
if(allocated(cdvdrhot)) deallocate(cdvdrhot)
if(allocated(cbmx)) deallocate(cbmx)
if(allocated(cbmn)) deallocate(cbmn)
if(allocated(crbav)) deallocate(crbav)
if(allocated(cvol)) deallocate(cvol)
if(allocated(crri)) deallocate(crri)
if(allocated(carea)) deallocate(carea)
if(allocated(cfc)) deallocate(cfc)
if(allocated(crhotq)) deallocate(crhotq)
if(allocated(cratjpl)) deallocate(cratjpl)
if(allocated(cratja)) deallocate(cratja)
if(allocated(cratjb)) deallocate(cratjb)
if(allocated(tjp)) deallocate(tjp,tlm,ch)
end subroutine dealloc_surfvec
subroutine flux_average
use const_and_precisions, only : wp_,zero,one,pi,ccj=>mu0inv
use gray_params, only : iequil
use equilibrium, only : btrcen,btaxis,rmaxis,zmaxis,phitedge,zbsup,zbinf, &
equian,equinum_psi,bfield,frhotor,fq,tor_curr
use simplespline, only : difcs
use dierckx, only : regrid,coeff_parder
implicit none
! local constants
integer, parameter :: nnintp=101,ncnt=100,nlam=101,ksp=3, &
njest=nnintp+ksp+1,nlest=nlam+ksp+1, &
lwrk=4*(nnintp+nlam)+11*(njest+nlest)+njest*nnintp+nlest+54, &
kwrk=nnintp+nlam+njest+nlest+3
! local variables
integer :: ier,ierr,l,jp,inc,inc1,iopt,njp,nlm,ninpr
integer, dimension(kwrk) :: iwrk
real(wp_) :: lam,dlam,anorm,b2av,dvdpsi,dadpsi,ratio_cdator, &
ratio_cdbtor,ratio_pltor,fc,height,r2iav,currp, &
area,volume,ajphiav,bbav,bmmx,bmmn,btot0,bpoloid0,rpsim0,dla,dlb, &
dlp,drc,ph,area2,rzp,rz,rpsim,zpsim,btot,bpoloid,dlph,ajphi0, &
shlam,srl,rl2,rl0,rl,dhlam,dhlam0,ccfh,s,ajphi, &
bphi,brr,bzz,riav,fp,psinjp,rhopjp,rhotjp,qq,rup,rlw,zup,zlw
real(wp_), dimension(nnintp) :: dadrhotv,dvdrhotv,vratjpl
real(wp_), dimension(2*ncnt) :: dlpv
real(wp_), dimension(2*ncnt+1) :: bv,bpv
real(wp_), dimension(nlam) :: alam,weights
real(wp_), dimension(nnintp,nlam) :: fhlam
real(wp_), dimension(nnintp*nlam) :: ffhlam,dffhlam
real(wp_), dimension(lwrk) :: wrk
real(wp_), dimension(:), allocatable :: rctemp,zctemp
! common/external functions/variables
real(wp_) :: fpolv,ddpsidrr,ddpsidzz
npsi=nnintp
npoints = 2*ncnt+1
call alloc_surfvec(ierr)
if(allocated(tjp)) deallocate(tjp)
if(allocated(tlm)) deallocate(tlm)
if(allocated(ch)) deallocate(ch)
allocate(tjp(njest),tlm(nlest),ch((njest-4)*(nlest-4)), &
rctemp(npoints),zctemp(npoints),stat=ierr)
if (ierr.ne.0) return
! computation of flux surface averaged quantities
dlam=1.0_wp_/dble(nlam-1)
do l=1,nlam-1
alam(l)=dble(l-1)*dlam
fhlam(1,l)=sqrt(1.0_wp_-alam(l))
ffhlam(l)=fhlam(1,l)
dffhlam(l)=-0.5_wp_/sqrt(1.0_wp_-alam(l))
weights(l)=1.0_wp_
end do
weights(1)=0.5_wp_
weights(nlam)=0.5_wp_
alam(nlam)=1.0_wp_
fhlam(1,nlam)=0.0_wp_
ffhlam(nlam)=0.0_wp_
dffhlam(nlam)=-99999.0_wp_
jp=1
anorm=2.0_wp_*pi*rmaxis/abs(btaxis)
dvdpsi=2.0_wp_*pi*anorm
dadpsi=2.0_wp_*pi/abs(btaxis)
b2av=btaxis**2
ratio_cdator=abs(btaxis/btrcen)
ratio_cdbtor=1.0_wp_
ratio_pltor=1.0_wp_
fc=1.0_wp_
if(iequil < 2) then
call equian(rmaxis,zmaxis,ddpsidrr=ddpsidrr,ddpsidzz=ddpsidzz)
else
call equinum_psi(rmaxis,zmaxis,ddpsidrr=ddpsidrr,ddpsidzz=ddpsidzz)
end if
qq=abs(btaxis)/sqrt(ddpsidrr*ddpsidzz)
ajphiav=-ccj*(ddpsidrr+ddpsidzz)/rmaxis
psicon(1)=0.0_wp_
rcon(:,1)=rmaxis
zcon(:,1)=zmaxis
pstab(1)=0.0_wp_
rpstab(1)=0.0_wp_
vcurrp(1)=0.0_wp_
vajphiav(1)=ajphiav
bmxpsi(1)=abs(btaxis)
bmnpsi(1)=abs(btaxis)
bav(1)=abs(btaxis)
rbav(1)=1.0_wp_
rri(1)=rmaxis
varea(1)=0.0_wp_
vvol(1)=0.0_wp_
vratjpl(1)=ratio_pltor
vratja(1)=ratio_cdator
vratjb(1)=ratio_cdbtor
ffc(1)=fc
qqv(1)=qq
dadrhotv(1)=0.0_wp_
dvdrhotv(1)=0.0_wp_
rup=rmaxis
rlw=rmaxis
zup=zmaxis+(zbsup-zmaxis)/10.0_wp_
zlw=zmaxis-(zmaxis-zbinf)/10.0_wp_
do jp=2,npsi
height=dble(jp-1)/dble(npsi-1)
if(jp.eq.npsi) height=0.9999_wp_
rhopjp=height
psinjp=height*height
rhotjp=frhotor(rhopjp)
psicon(jp)=height
call contours_psi(psinjp,rctemp,zctemp,rup,zup,rlw,zlw)
rcon(:,jp) = rctemp
zcon(:,jp) = zctemp
r2iav=0.0_wp_
anorm=0.0_wp_
dadpsi=0.0_wp_
currp=0.0_wp_
b2av=0.0_wp_
area=0.0_wp_
volume=0.0_wp_
ajphiav=0.0_wp_
bbav=0.0_wp_
bmmx=-1.0e+30_wp_
bmmn=1.0e+30_wp_
call tor_curr(rctemp(1),zctemp(1),ajphi0)
call bfield(rctemp(1),zctemp(1),bphi,br=brr,bz=bzz)
fpolv=bphi*rctemp(1)
btot0=sqrt(bphi**2+brr**2+bzz**2)
bpoloid0=sqrt(brr**2+bzz**2)
bv(1)=btot0
bpv(1)=bpoloid0
rpsim0=rctemp(1)
do inc=1,npoints-1
inc1=inc+1
dla=sqrt((rctemp(inc)-rmaxis)**2+(zctemp(inc)-zmaxis)**2)
dlb=sqrt((rctemp(inc1)-rmaxis)**2+(zctemp(inc1)-zmaxis)**2)
dlp=sqrt((rctemp(inc1)-rctemp(inc))**2+(zctemp(inc1)-zctemp(inc))**2)
drc=(rctemp(inc1)-rctemp(inc))
! compute length, area and volume defined by psi=psinjp=height^2
ph=0.5_wp_*(dla+dlb+dlp)
area2=ph*(ph-dla)*(ph-dlb)*(ph-dlp)
area=area+sqrt(area2)
rzp=rctemp(inc1)*zctemp(inc1)
rz=rctemp(inc)*zctemp(inc)
volume=pi*(rzp+rz)*drc+volume
! compute line integrals on the contour psi=psinjp=height^2
rpsim=rctemp(inc1)
zpsim=zctemp(inc1)
call bfield(rpsim,zpsim,br=brr,bz=bzz)
call tor_curr(rpsim,zpsim,ajphi)
bphi=fpolv/rpsim
btot=sqrt(bphi**2+brr**2+bzz**2)
bpoloid=sqrt(brr**2+bzz**2)
dlpv(inc)=dlp
bv(inc1)=btot
bpv(inc1)=bpoloid
dlph=0.5_wp_*dlp
anorm=anorm+dlph*(1.0_wp_/bpoloid+1.0_wp_/bpoloid0)
dadpsi=dadpsi+dlph*(1.0_wp_/(bpoloid*rpsim)+1.0_wp_/(bpoloid0*rpsim0))
currp=currp+dlph*(bpoloid+bpoloid0)
b2av=b2av+dlph*(btot0**2/bpoloid0+btot**2/bpoloid)
bbav=bbav+dlph*(btot/bpoloid+btot0/bpoloid0)
r2iav=r2iav+dlph*(1.0_wp_/(bpoloid*rpsim**2)+1.0_wp_/(bpoloid0*rpsim0**2))
ajphiav=ajphiav+dlph*(ajphi0/(bpoloid0*rpsim0)+ajphi/(bpoloid*rpsim))
ajphi0=ajphi
rpsim0=rpsim
bpoloid0=bpoloid
btot0=btot
! computation maximum/minimum B values on given flux surface
if(btot.le.bmmn) bmmn=btot
if(btot.ge.bmmx) bmmx=btot
end do
! bav=<B> [T] , b2av=<B^2> [T^2] , rbav=<B>/b_min
! anorm = int d l_p/B_p = dV/dpsi/(2pi)
! r2iav=<1/R^2> [m^-2] ,
! riav=<1/R> [m^-1] = dA/dpsi/(dV/dpsi/(2pi)),
! rri = <B>/(|R B_tor|<1/R^2>) , used to compute I_tor [m^-1]
! currp = plasma current within psi=const
bbav=bbav/anorm
r2iav=r2iav/anorm
dvdpsi=2.0_wp_*pi*anorm
riav=dadpsi/anorm
b2av=b2av/anorm
vcurrp(jp)=ccj*currp
vajphiav(jp)=ajphiav/dadpsi
! area == varea, volume == vvol
! flux surface minor radius == (area/pi)^1/2
! ratio_cdator = Jcd_astra/J_phi Jcd_astra = <Jcd.B>/B0
! ratio_cdbtor = Jcd_jintrac/J_phi Jcd_jintrac = <Jcd.B>/<B>
! ratio_pltor = Jcd_||/J_phi Jcd_|| = <Jcd.b>
pstab(jp)=psinjp
rpstab(jp)=rhopjp
vvol(jp)=abs(volume)
varea(jp)=area
bav(jp)=bbav
rbav(jp)=bbav/bmmn
bmxpsi(jp)=bmmx
bmnpsi(jp)=bmmn
rri(jp)=bav(jp)/abs(fpolv*r2iav)
ratio_cdator=abs(b2av*riav/(fpolv*r2iav*btrcen))
ratio_cdbtor=abs(b2av*riav/(fpolv*r2iav*bbav))
ratio_pltor=abs(bbav*riav/(fpolv*r2iav))
vratjpl(jp)=ratio_pltor
vratja(jp)=ratio_cdator
vratjb(jp)=ratio_cdbtor
qq=abs(dvdpsi*fpolv*r2iav/(4.0_wp_*pi*pi))
qqv(jp)=qq
dadrhotv(jp)=phitedge*rhotjp/fq(psinjp)*dadpsi/pi
dvdrhotv(jp)=phitedge*rhotjp/fq(psinjp)*dvdpsi/pi
! computation of fraction of circulating/trapped fraction fc, ft
! and of function H(lambda,rhop)
! ffhlam = Bmn/Bmx/fc integral_lambda^1 dlam/<sqrt(1-lam*B(rhop)/Bmx)>
fc=0.0_wp_
shlam=0.0_wp_
do l=nlam,1,-1
lam=alam(l)
srl=0.0_wp_
rl2=1.0_wp_-lam*bv(1)/bmmx
rl0=0.0_wp_
if(rl2.gt.0) rl0=sqrt(rl2)
do inc=1,npoints-1
rl2=1.0_wp_-lam*bv(inc+1)/bmmx
rl=0.0_wp_
if(rl2.gt.0) rl=sqrt(rl2)
srl=srl+0.5_wp_*dlpv(inc)*(rl/bpv(inc+1)+rl0/bpv(inc))
rl0=rl
end do
srl=srl/anorm
dhlam=0.5_wp_/srl
fc=fc+lam/srl*weights(l)
if(l.eq.nlam) then
fhlam(jp,l)=0.0_wp_
ffhlam(nlam*(jp-1)+l)=0.0_wp_
dffhlam(nlam*(jp-1)+l)=-dhlam
dhlam0=dhlam
else
shlam=shlam+0.5_wp_*(dhlam+dhlam0)*dlam
fhlam(jp,l)=shlam
dffhlam(nlam*(jp-1)+l)=-dhlam
dhlam0=dhlam
end if
end do
fc=0.75_wp_*b2av/bmmx**2*fc*dlam
ffc(jp)=fc
ccfh=bmmn/bmmx/fc
do l=1,nlam
ffhlam(nlam*(jp-1)+l)=ccfh*fhlam(jp,l)
dffhlam(nlam*(jp-1)+l)=ccfh*dffhlam(nlam*(jp-1)+l)
end do
end do
rpstab(npsi)=1.0_wp_
pstab(npsi)=1.0_wp_
! spline coefficients of area,vol,rbav,rri,bmxpsi,bmnpsi,fc,dadrhot,dvdrhot,ratioJs
! used for computations of dP/dV and J_cd
iopt=0
call difcs(rpstab,vvol,npsi,iopt,cvol,ier)
iopt=0
call difcs(rpstab,rbav,npsi,iopt,crbav,ier)
iopt=0
call difcs(rpstab,rri,npsi,iopt,crri,ier)
iopt=0
call difcs(rpstab,bmxpsi,npsi,iopt,cbmx,ier)
iopt=0
call difcs(rpstab,bmnpsi,npsi,iopt,cbmn,ier)
iopt=0
call difcs(rpstab,vratja,npsi,iopt,cratja,ier)
iopt=0
call difcs(rpstab,vratjb,npsi,iopt,cratjb,ier)
iopt=0
call difcs(rpstab,vratjpl,npsi,iopt,cratjpl,ier)
iopt=0
call difcs(rpstab,varea,npsi,iopt,carea,ier)
iopt=0
call difcs(rpstab,ffc,npsi,iopt,cfc,ier)
iopt=0
call difcs(rpstab,dadrhotv,npsi,iopt,cdadrhot,ier)
iopt=0
call difcs(rpstab,dvdrhotv,npsi,iopt,cdvdrhot,ier)
! iopt=0
! call difcs(rpstab,qqv,npsi,iopt,cqq,ier)
! spline interpolation of H(lambda,rhop) and dH/dlambda
iopt=0
s=0.0_wp_
call regrid(iopt,npsi,rpstab,nlam,alam,ffhlam,zero,one,zero,one, &
ksp,ksp,s,njest,nlest,njp,tjp,nlm,tlm,ch,fp, &
wrk,lwrk,iwrk,kwrk,ier)
njpt=njp
nlmt=nlm
do jp=1,npsi
call print_fluxav(rpstab(jp),frhotor(rpstab(jp)),bav(jp),bmxpsi(jp), &
bmnpsi(jp),varea(jp),vvol(jp),vcurrp(jp),vajphiav(jp), &
ffc(jp),vratja(jp),vratjb(jp),qqv(jp))
end do
ninpr=(npsi-1)/10
do jp=ninpr+1,npsi,ninpr
call print_contour(psicon(jp),rcon(:,jp),zcon(:,jp))
end do
deallocate(rctemp, zctemp)
end subroutine flux_average
subroutine fluxval(rhop,area,vol,dervol,dadrhot,dvdrhot, &
rri,rbav,bmn,bmx,fc,ratja,ratjb,ratjpl)
use const_and_precisions, only : wp_
use utils, only : locate
use simplespline, only :spli,splid
implicit none
! arguments
real(wp_), intent(in) :: rhop
real(wp_), intent(out), optional :: vol,area,rri,rbav,dervol,bmn,bmx,fc, &
ratja,ratjb,ratjpl,dadrhot,dvdrhot
! local variables
integer :: ip
real(wp_) :: drh
call locate(rpstab,npsi,rhop,ip)
ip=min(max(1,ip),npsi-1)
drh=rhop-rpstab(ip)
if (present(area)) area=spli(carea,npsi,ip,drh)
if (present(vol)) vol=spli(cvol,npsi,ip,drh)
if (present(dervol)) dervol=splid(cvol,npsi,ip,drh)
if (present(dadrhot)) dadrhot=spli(cdadrhot,npsi,ip,drh)
if (present(dvdrhot)) dvdrhot=spli(cdvdrhot,npsi,ip,drh)
if (present(rri)) rri=spli(crri,npsi,ip,drh)
if (present(rbav)) rbav=spli(crbav,npsi,ip,drh)
if (present(bmn)) bmn=spli(cbmn,npsi,ip,drh)
if (present(bmx)) bmx=spli(cbmx,npsi,ip,drh)
if (present(fc)) fc=spli(cfc,npsi,ip,drh)
if (present(ratja)) ratja=spli(cratja,npsi,ip,drh)
if (present(ratjb)) ratjb=spli(cratjb,npsi,ip,drh)
if (present(ratjpl)) ratjpl=spli(cratjpl,npsi,ip,drh)
end subroutine fluxval
subroutine contours_psi(h,rcn,zcn,rup,zup,rlw,zlw)
use const_and_precisions, only : wp_,pi
use gray_params, only : iequil
use dierckx, only : profil,sproota
use equilibrium, only : rmaxis,zmaxis,aminor,frhotor,tr,nsr,tz,nsz,cceq, &
kspl,psiant,psinop,points_tgo
use limiter, only : rwallm
implicit none
! local constants
integer, parameter :: mest=4
! arguments
real(wp_), intent(in) :: h
real(wp_), dimension(:), intent(out) :: rcn,zcn
real(wp_), intent(inout) :: rup,zup,rlw,zlw
! local variables
integer :: npoints,np,info,ic,ier,iopt,m
real(wp_) :: ra,rb,za,zb,rn,th,zc,val
real(wp_), dimension(mest) :: zeroc
real(wp_), dimension(nsr) :: czc
npoints=size(rcn)
np=(npoints-1)/2
th=pi/dble(np)
if (iequil<2) then
rn=frhotor(sqrt(h))
do ic=1,npoints
zcn(ic)=zmaxis+aminor*rn*sin(th*(ic-1))
rcn(ic)=rmaxis+aminor*rn*cos(th*(ic-1))
end do
else
ra=rup
rb=rlw
za=zup
zb=zlw
call points_tgo(ra,za,rup,zup,h,info)
call points_tgo(rb,zb,rlw,zlw,h,info)
rcn(1)=rlw
zcn(1)=zlw
rcn(npoints)=rlw
zcn(npoints)=zlw
rcn(np+1)=rup
zcn(np+1)=zup
do ic=2,np
zc=zlw+(zup-zlw)*(1.0_wp_-cos(th*(ic-1)))/2.0_wp_
iopt=1
call profil(iopt,tr,nsr,tz,nsz,cceq,kspl,kspl,zc,nsr,czc,ier)
if(ier.gt.0) print*,' profil =',ier
val=h*psiant+psinop
call sproota(val,tr,nsr,czc,zeroc,mest,m,ier)
if (zeroc(1).gt.rwallm) then
rcn(ic)=zeroc(1)
rcn(npoints+1-ic)=zeroc(2)
else
rcn(ic)=zeroc(2)
rcn(npoints+1-ic)=zeroc(3)
end if
zcn(ic)=zc
zcn(npoints+1-ic)=zc
end do
end if
end subroutine contours_psi
subroutine print_contour(psin, rc, zc)
! Prints the flux surface contours (unit 71)
use const_and_precisions, only : wp_, comp_tiny
use units, only : ucnt, unit_active
implicit none
! subroutine arguments
real(wp_), intent(in) :: psin
real(wp_), dimension(:), intent(in) :: rc, zc
! local variables
integer :: npoints, ic
if (.not. unit_active(ucnt)) return
if (psin < comp_tiny) then
write(ucnt, *)' #i psin R z'
else
npoints = size(rc)
do ic=1,npoints
write(ucnt, '(i6,12(1x,e12.5))') ic, psin, rc(ic), zc(ic)
end do
write(ucnt, *)
write(ucnt, *)
end if
end subroutine print_contour
subroutine print_fluxav(psin, rhot, bav, bmx, bmn, area, vol, &
currp, ajphiav, ffc, ratja, ratjb, qq)
! Prints the flux-averaged quantities table (unit 56)
use const_and_precisions, only : wp_, comp_tiny
use units, only : uflx, unit_active
implicit none
! subroutine arguments
real(wp_), intent(in) :: psin, rhot, bav, bmx, bmn, area, vol, &
currp, ajphiav, ffc, ratja, ratjb, qq
if (.not. unit_active(uflx)) return
if (psin < comp_tiny) &
write(uflx, *)' #rhop rhot |<B>| |Bmx| |Bmn| Area Vol' // &
' |I_pl| <J_phi> fc ratJa ratJb qq'
write(uflx, '(20(1x,e12.5))') &
psin, rhot, bav, bmx, bmn, area, vol, currp, ajphiav, &
ffc, ratja, ratjb, qq
end subroutine print_fluxav
end module magsurf_data