SPEC 3.20
Stepped Pressure Equilibrium Code
Output file(s)

Modules

module  sphdf5
 writing the HDF5 output file
 

Functions/Subroutines

subroutine ra00aa (writeorread)
 Writes vector potential to .ext.sp.A . More...
 
subroutine sphdf5::init_outfile
 Initialize the interface to the HDF5 library and open the output file.
 
subroutine sphdf5::mirror_input_to_outfile
 Mirror input variables into output file. More...
 
subroutine sphdf5::init_convergence_output
 Prepare convergence evolution output. More...
 
subroutine sphdf5::write_convergence_output (nDcalls, ForceErr)
 Write convergence output (evolution of interface geometry, force, etc).
 
subroutine sphdf5::write_grid
 Write the magnetic field on a grid. More...
 
subroutine sphdf5::init_flt_output (numTrajTotal)
 Initialize field line tracing output group and create array datasets. More...
 
subroutine sphdf5::write_poincare (offset, data, success)
 Write a hyperslab of Poincare data corresponding to the output of one parallel worker. More...
 
subroutine sphdf5::write_transform (offset, length, lvol, diotadxup, fiota)
 Write the rotational transform output from field line following. More...
 
subroutine sphdf5::finalize_flt_output
 Finalize Poincare output. More...
 
subroutine sphdf5::write_vector_potential (sumLrad, allAte, allAze, allAto, allAzo)
 Write the magnetic vector potential Fourier harmonics to the output file group /vector_potential . More...
 
subroutine sphdf5::hdfint
 Write the final state of the equilibrium to the output file. More...
 
subroutine sphdf5::finish_outfile
 Close all open HDF5 objects (we know of) and list any remaining still-open objects.
 

Detailed Description

Function/Subroutine Documentation

◆ ra00aa()

subroutine ra00aa ( character, intent(in)  writeorread)

Writes vector potential to .ext.sp.A .

representation of vector potential

  • The components of the vector potential, \({\bf A}=A_\theta \nabla + A_\zeta \nabla \zeta\), are

    \begin{eqnarray} A_\theta(s,\theta,\zeta) &=& \sum_{i,l} {\color{red} A_{\theta,e,i,l}} \; {\overline T}_{l,i}(s) \cos\alpha_i + \sum_{i,l} {\color{Orange} A_{\theta,o,i,l}} \; {\overline T}_{l,i}(s) \sin\alpha_i, \label{eq:At_ra00aa} \\ A_\zeta( s,\theta,\zeta) &=& \sum_{i,l} {\color{blue} A_{\zeta, e,i,l}} \; {\overline T}_{l,i}(s) \cos\alpha_i + \sum_{i,l} {\color{Cerulean}A_{\zeta ,o,i,l}} \; {\overline T}_{l,i}(s) \sin\alpha_i, \label{eq:Az_ra00aa} \end{eqnarray}

    where \({\overline T}_{l,i}(s) \equiv \bar s^{m_i/2} \, T_l(s)\), \(T_l(s)\) is the Chebyshev polynomial, and \(\alpha_j \equiv m_j\theta-n_j\zeta\). The regularity factor, \(\bar s^{m_i/2}\), where \(\bar s \equiv (1+s)/2\), is only included if there is a coordinate singularity in the domain (i.e. only in the innermost volume, and only in cylindrical and toroidal geometry.)

file format

  • The format of the files containing the vector potential is as follows:
    open(aunit, file="."//trim(ext)//".sp.A", status="replace", form="unformatted" )
    write(aunit) mvol, mpol, ntor, mn, nfp ! integers;
    write(aunit) im(1:mn) ! integers; poloidal modes;
    write(aunit) in(1:mn) ! integers; toroidal modes;
    do vvol = 1, mvol ! integers; loop over volumes;
    write(aunit) lrad(vvol) ! integers; the radial resolution in each volume may be different;
    do ii = 1, mn
    write(aunit) ate(vvol,ii)%s(0:lrad(vvol)) ! reals;
    write(aunit) aze(vvol,ii)%s(0:lrad(vvol)) ! reals;
    write(aunit) ato(vvol,ii)%s(0:lrad(vvol)) ! reals;
    write(aunit) azo(vvol,ii)%s(0:lrad(vvol)) ! reals;
    enddo ! end of do ii;
    enddo ! end of do vvol;
    close(aunit)
Parameters
[in]writeorread'W' to write the vector potential; 'R' to read it

References allglobal::ate, allglobal::ato, fileunits::aunit, allglobal::aze, allglobal::azo, allglobal::cpus, allglobal::im, allglobal::in, inputlist::lrad, allglobal::mn, allglobal::mpi_comm_spec, inputlist::mpol, allglobal::myid, allglobal::ncpu, inputlist::nfp, inputlist::ntor, fileunits::ounit, ra00aa(), inputlist::wmacros, sphdf5::write_vector_potential(), and constants::zero.

Referenced by preset(), ra00aa(), and spec().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ mirror_input_to_outfile()

subroutine sphdf5::mirror_input_to_outfile

Mirror input variables into output file.

The goal of this routine is to have an exact copy of the input file contents that were used to parameterize a given SPEC run. This also serves to check after the run if SPEC correctly understood the text-based input file.

References inputlist::absacc, inputlist::absreq, inputlist::adiabatic, inputlist::bnc, inputlist::bns, inputlist::bnsblend, inputlist::bnstol, inputlist::c05factor, inputlist::c05xmax, inputlist::c05xtol, inputlist::curpol, inputlist::curtor, inputlist::dpp, inputlist::dqq, inputlist::epsgmres, inputlist::epsilon, inputlist::epsilu, inputlist::epsr, inputlist::escale, sphdf5::file_id, inputlist::forcetol, inputlist::fudge, inputlist::gamma, inputlist::gbnbld, inputlist::gbntol, inputlist::helicity, inputlist::igeometry, inputlist::imethod, inputlist::impol, allglobal::in, inputlist::intor, inputlist::iorder, inputlist::iota, inputlist::iotatol, inputlist::iprecon, inputlist::istellsym, inputlist::isurf, inputlist::ivolume, inputlist::ladiabatic, inputlist::lbeltrami, inputlist::lcheck, inputlist::lconstraint, inputlist::lerrortype, inputlist::lextrap, inputlist::lfindzero, inputlist::lfreebound, inputlist::lgmresprec, inputlist::lhevalues, inputlist::lhevectors, inputlist::lhmatrix, inputlist::linitgues, inputlist::linitialize, inputlist::lmatsolver, inputlist::lp, inputlist::lperturbed, inputlist::lposdef, inputlist::lq, inputlist::lrad, inputlist::lreadgf, inputlist::lreflect, inputlist::lrzaxis, inputlist::lsparse, inputlist::lsvdiota, inputlist::ltiming, inputlist::ltransform, inputlist::lzerovac, inputlist::maxrndgues, inputlist::mcasingcal, inputlist::mfreeits, inputlist::mpol, inputlist::mregular, inputlist::mu, inputlist::mupfits, inputlist::mupftol, allglobal::myid, inputlist::ndiscrete, inputlist::nfp, inputlist::ngrid, inputlist::nppts, inputlist::nptrj, inputlist::nquad, inputlist::ntor, inputlist::ntoraxis, inputlist::nvol, inputlist::odetol, inputlist::oita, inputlist::opsilon, inputlist::pcondense, inputlist::pflux, inputlist::phiedge, inputlist::pl, inputlist::ppts, inputlist::pr, inputlist::pressure, inputlist::pscale, inputlist::ql, inputlist::qr, inputlist::rac, inputlist::ras, inputlist::rbc, inputlist::rbs, inputlist::relreq, inputlist::rp, inputlist::rpol, inputlist::rq, inputlist::rtor, inputlist::rwc, inputlist::rws, inputlist::scaling, inputlist::tflux, inputlist::upsilon, inputlist::vcasingeps, inputlist::vcasingits, inputlist::vcasingper, inputlist::vcasingtol, inputlist::vnc, inputlist::vns, volume(), inputlist::wpoloidal, inputlist::zac, inputlist::zas, inputlist::zbc, inputlist::zbs, inputlist::zwc, and inputlist::zws.

Referenced by xspech().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ init_convergence_output()

subroutine sphdf5::init_convergence_output

Prepare convergence evolution output.

  • The group iterations is created in the output file. This group contains the interface geometry at each iteration, which is useful for constructing movies illustrating the convergence. The data structure in use is an unlimited array of the following compound datatype:
    DATATYPE H5T_COMPOUND {
    H5T_NATIVE_INTEGER "nDcalls";
    H5T_NATIVE_DOUBLE "Energy";
    H5T_NATIVE_DOUBLE "ForceErr";
    H5T_ARRAY { [Mvol+1][mn] H5T_NATIVE_DOUBLE } "iRbc";
    H5T_ARRAY { [Mvol+1][mn] H5T_NATIVE_DOUBLE } "iZbs";
    H5T_ARRAY { [Mvol+1][mn] H5T_NATIVE_DOUBLE } "iRbs";
    H5T_ARRAY { [Mvol+1][mn] H5T_NATIVE_DOUBLE } "iZbc";
    }

References sphdf5::dt_energy_id, sphdf5::dt_forceerr_id, sphdf5::dt_irbc_id, sphdf5::dt_irbs_id, sphdf5::dt_izbc_id, sphdf5::dt_izbs_id, sphdf5::dt_ndcalls_id, sphdf5::file_id, sphdf5::hdfier, sphdf5::iteration_dset_id, sphdf5::memspace, allglobal::mn, allglobal::myid, and sphdf5::plist_id.

Referenced by xspech().

Here is the caller graph for this function:

◆ write_grid()

subroutine sphdf5::write_grid

Write the magnetic field on a grid.

The magnetic field is evaluated on a regular grid in \((s, \theta, \zeta)\) and the corresponding cylindrical coordinates \((R,Z)\) as well as the cylindrical components of the magnetic field \((B^R, B^\varphi, B^Z)\) are written out.

References bfield(), coords(), sphdf5::file_id, allglobal::gbzeta, inputlist::igeometry, allglobal::ijimag, allglobal::ijreal, allglobal::ivol, allglobal::jireal, allglobal::lcoordinatesingularity, allglobal::lplasmaregion, inputlist::lrad, allglobal::lvacuumregion, allglobal::mn, allglobal::myid, inputlist::ngrid, allglobal::node, allglobal::nt, allglobal::ntz, inputlist::nvol, allglobal::nz, constants::one, constants::pi2, allglobal::rij, inputlist::rpol, inputlist::rtor, allglobal::sg, constants::two, constants::zero, and allglobal::zij.

Referenced by xspech().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ init_flt_output()

subroutine sphdf5::init_flt_output ( integer, intent(in)  numTrajTotal)

Initialize field line tracing output group and create array datasets.

The field-line tracing diagnostic is parallelized over volumes, where all threads/ranks produce individual output. This is gathered in the output file, stacked over the radial dimension. The success flag signals if the integrator was successful in following the fieldline for the derired number of toroidal periods.

Parameters
[in]numTrajTotaltotal number of Poincare trajectories

References sphdf5::dset_id_diotadxup, sphdf5::dset_id_fiota, sphdf5::dset_id_r, sphdf5::dset_id_s, sphdf5::dset_id_success, sphdf5::dset_id_t, sphdf5::dset_id_z, sphdf5::file_id, sphdf5::filespace_diotadxup, sphdf5::filespace_fiota, sphdf5::filespace_r, sphdf5::filespace_s, sphdf5::filespace_success, sphdf5::filespace_t, sphdf5::filespace_z, sphdf5::grppoincare, sphdf5::grptransform, sphdf5::hdfier, allglobal::lmns, sphdf5::memspace_diotadxup, sphdf5::memspace_r, sphdf5::memspace_s, sphdf5::memspace_success, sphdf5::memspace_t, sphdf5::memspace_z, allglobal::myid, inputlist::nppts, allglobal::nz, sphdf5::rankp, and sphdf5::rankt.

Referenced by pp00aa().

Here is the caller graph for this function:

◆ write_poincare()

subroutine sphdf5::write_poincare ( integer, intent(in)  offset,
real, dimension(:,:,:), intent(in)  data,
integer, dimension(:), intent(in)  success 
)

Write a hyperslab of Poincare data corresponding to the output of one parallel worker.

Parameters
offsetradial offset at which the data belongs
dataoutput from field-line tracing
successflags to indicate if integrator was successful

References sphdf5::dset_id_r, sphdf5::dset_id_s, sphdf5::dset_id_success, sphdf5::dset_id_t, sphdf5::dset_id_z, sphdf5::filespace_r, sphdf5::filespace_s, sphdf5::filespace_success, sphdf5::filespace_t, sphdf5::filespace_z, sphdf5::hdfier, sphdf5::memspace_r, sphdf5::memspace_s, sphdf5::memspace_success, sphdf5::memspace_t, sphdf5::memspace_z, allglobal::myid, inputlist::nppts, and allglobal::nz.

Referenced by pp00aa().

Here is the caller graph for this function:

◆ write_transform()

subroutine sphdf5::write_transform ( integer, intent(in)  offset,
integer, intent(in)  length,
integer, intent(in)  lvol,
real, dimension(:), intent(in)  diotadxup,
real, dimension(:,:), intent(in)  fiota 
)

Write the rotational transform output from field line following.

Parameters
offsetradial offset at which the data belongs
lengthlength of dataset to write
lvolnested volume index
diotadxupderivative of rotational transform (?)
fiotarotational transform

References sphdf5::dset_id_diotadxup, sphdf5::dset_id_fiota, sphdf5::filespace_diotadxup, sphdf5::filespace_fiota, sphdf5::hdfier, sphdf5::memspace_diotadxup, sphdf5::memspace_fiota, and sphdf5::rankt.

Referenced by pp00aa().

Here is the caller graph for this function:

◆ finalize_flt_output()

subroutine sphdf5::finalize_flt_output

Finalize Poincare output.

This closes the still-open datasets related to field-line tracing, which had to be kept open during the tracing to be able to write the outputs directly when a given worker thread is finished.

References sphdf5::dset_id_diotadxup, sphdf5::dset_id_fiota, sphdf5::dset_id_r, sphdf5::dset_id_s, sphdf5::dset_id_success, sphdf5::dset_id_t, sphdf5::dset_id_z, sphdf5::filespace_diotadxup, sphdf5::filespace_fiota, sphdf5::filespace_r, sphdf5::filespace_s, sphdf5::filespace_success, sphdf5::filespace_t, sphdf5::filespace_z, sphdf5::grppoincare, sphdf5::grptransform, sphdf5::hdfier, sphdf5::memspace_diotadxup, sphdf5::memspace_r, sphdf5::memspace_s, sphdf5::memspace_success, sphdf5::memspace_t, and sphdf5::memspace_z.

Referenced by pp00aa().

Here is the caller graph for this function:

◆ write_vector_potential()

subroutine sphdf5::write_vector_potential ( integer, intent(in)  sumLrad,
real, dimension(:,:), intent(in)  allAte,
real, dimension(:,:), intent(in)  allAze,
real, dimension(:,:), intent(in)  allAto,
real, dimension(:,:), intent(in)  allAzo 
)

Write the magnetic vector potential Fourier harmonics to the output file group /vector_potential .

The data is stacked in the radial direction over Lrad , since Lrad can be different in each volume, but HDF5 only supports rectangular arrays. So, one needs to split the sumLrad dimension into chunks given by the input Lrad array.

Parameters
sumLradtotal sum over Lrad in all nested volumes
allAte\(A^{\theta}_\mathrm{even}\) for all nested volumes
allAze\(A^{\zeta}_\mathrm{even}\) for all nested volumes
allAto\(A^{\theta}_\mathrm{odd}\) for all nested volumes
allAzo\(A^{\zeta}_\mathrm{odd}\) for all nested volumes

References allglobal::ate, allglobal::ato, allglobal::aze, allglobal::azo, sphdf5::file_id, allglobal::mn, and allglobal::myid.

Referenced by ra00aa().

Here is the caller graph for this function:

◆ hdfint()

subroutine sphdf5::hdfint

Write the final state of the equilibrium to the output file.

  • In addition to the input variables, which are described in global(), the following quantities are written to ext.sp.h5 :

  • All quantities marked as real should be treated as double precision.

References inputlist::adiabatic, allglobal::beltramierror, inputlist::bnc, inputlist::bns, allglobal::bnserr, allglobal::bsupumn, allglobal::bsupvmn, allglobal::btemn, allglobal::btomn, allglobal::bzemn, allglobal::bzomn, allglobal::cpus, allglobal::drbc, allglobal::drbs, allglobal::dvolume, allglobal::dzbc, allglobal::dzbs, sphdf5::file_id, allglobal::forceerr, inputlist::helicity, allglobal::ibnc, allglobal::ibns, allglobal::im, allglobal::ims, allglobal::in, allglobal::ins, allglobal::irbc, allglobal::irbs, allglobal::ivnc, allglobal::ivns, inputlist::ivolume, allglobal::izbc, allglobal::izbs, inputlist::lcheck, allglobal::lmns, inputlist::lperturbed, inputlist::lrad, allglobal::mn, allglobal::mns, inputlist::mu, allglobal::myid, allglobal::ncpu, inputlist::nvol, fileunits::ounit, inputlist::pflux, inputlist::rbc, inputlist::rbs, inputlist::tflux, allglobal::tt, inputlist::vnc, inputlist::vns, volume(), allglobal::vvolume, inputlist::zbc, and inputlist::zbs.

Referenced by xspech().

Here is the call graph for this function:
Here is the caller graph for this function: