SPEC 3.20
Stepped Pressure Equilibrium Code
Integrals

Functions/Subroutines

subroutine df00ab (pNN, xi, Fxi, DFxi, Ldfjac, iflag)
 Evaluates volume integrals, and their derivatives w.r.t. interface geometry, using "packed" format. More...
 
subroutine ma00aa (lquad, mn, lvol, lrad)
 Calculates volume integrals of Chebyshev polynomials and metric element products. More...
 
subroutine spsint (lquad, mn, lvol, lrad)
 Calculates volume integrals of Chebyshev-polynomials and metric elements for preconditioner. More...
 

Detailed Description

Function/Subroutine Documentation

◆ df00ab()

subroutine df00ab ( integer, intent(in)  pNN,
real, dimension(0:pnn-1), intent(in)  xi,
real, dimension(0:pnn-1), intent(out)  Fxi,
real, dimension(0:ldfjac-1,0:pnn-1), intent(out)  DFxi,
integer, intent(in)  Ldfjac,
integer, intent(in), value  iflag 
)

Evaluates volume integrals, and their derivatives w.r.t. interface geometry, using "packed" format.

Parameters
[in]pNN
[in]xi
[out]Fxi
[out]DFxi
[in]Ldfjac
[in]iflag

References allglobal::cpus, df00ab(), allglobal::dma, allglobal::dmd, constants::half, inputlist::helicity, allglobal::ivol, allglobal::mbpsi, allglobal::mpi_comm_spec, allglobal::myid, inputlist::nvol, constants::one, fileunits::ounit, numerical::small, constants::two, and constants::zero.

Referenced by df00ab(), and ma02aa().

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

◆ ma00aa()

subroutine ma00aa ( integer, intent(in)  lquad,
integer, intent(in)  mn,
integer, intent(in)  lvol,
integer, intent(in)  lrad 
)

Calculates volume integrals of Chebyshev polynomials and metric element products.

Chebyshev-metric information

  • The following quantities are calculated:

    \begin{eqnarray} \verb+DToocc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j} \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \cos\alpha_j \\ \verb+DToocs(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j} \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \sin\alpha_j \\ \verb+DToosc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j} \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \cos\alpha_j \\ \verb+DTooss(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j} \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \sin\alpha_j \end{eqnarray}

    \begin{eqnarray} \verb+TTsscc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i} \; {\overline T}_{p,j} \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \cos\alpha_j \; \bar g_{ss} \\ \verb+TTsscs(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i} \; {\overline T}_{p,j} \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \sin\alpha_j \; \bar g_{ss} \\ \verb+TTsssc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i} \; {\overline T}_{p,j} \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \cos\alpha_j \; \bar g_{ss} \\ \verb+TTssss(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i} \; {\overline T}_{p,j} \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \sin\alpha_j \; \bar g_{ss} \end{eqnarray}

    \begin{eqnarray} \verb+TDstcc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i} \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \cos\alpha_j \; \bar g_{s\theta} \\ \verb+TDstcs(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i} \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \sin\alpha_j \; \bar g_{s\theta} \\ \verb+TDstsc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i} \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \cos\alpha_j \; \bar g_{s\theta} \\ \verb+TDstss(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i} \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \sin\alpha_j \; \bar g_{s\theta} \end{eqnarray}

    \begin{eqnarray} \verb+TDstcc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i} \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \cos\alpha_j \; \bar g_{s\zeta} \\ \verb+TDstcs(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i} \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \sin\alpha_j \; \bar g_{s\zeta} \\ \verb+TDstsc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i} \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \cos\alpha_j \; \bar g_{s\zeta} \\ \verb+TDstss(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i} \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \sin\alpha_j \; \bar g_{s\zeta} \end{eqnarray}

    \begin{eqnarray} \verb+DDstcc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \cos\alpha_j \; \bar g_{\theta\theta} \\ \verb+DDstcs(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \sin\alpha_j \; \bar g_{\theta\theta} \\ \verb+DDstsc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \cos\alpha_j \; \bar g_{\theta\theta} \\ \verb+DDstss(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \sin\alpha_j \; \bar g_{\theta\theta} \end{eqnarray}

    \begin{eqnarray} \verb+DDstcc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \cos\alpha_j \; \bar g_{\theta\zeta} \\ \verb+DDstcs(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \sin\alpha_j \; \bar g_{\theta\zeta} \\ \verb+DDstsc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \cos\alpha_j \; \bar g_{\theta\zeta} \\ \verb+DDstss(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \sin\alpha_j \; \bar g_{\theta\zeta} \end{eqnarray}

    \begin{eqnarray} \verb+DDstcc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \cos\alpha_j \; \bar g_{\zeta\zeta} \\ \verb+DDstcs(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \cos\alpha_i \sin\alpha_j \; \bar g_{\zeta\zeta} \\ \verb+DDstsc(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \cos\alpha_j \; \bar g_{\zeta\zeta} \\ \verb+DDstss(l,p,i,j)+ & \equiv & \int ds \; {\overline T}_{l,i}^\prime \; {\overline T}_{p,j}^\prime \; \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \sin\alpha_i \sin\alpha_j \; \bar g_{\zeta\zeta} \end{eqnarray}

    where \({\overline T}_{l,i}\equiv T_l \, \bar s^{m_i/2}\) if the domain includes the coordinate singularity, and \({\overline T}_{l,i}\equiv T_l\) if not; and \(\bar g_{\mu\nu} \equiv g_{\mu\nu} / \sqrt g\).

  • The double-angle formulae are used to reduce the above expressions to the Fourier harmonics of \(\bar g_{\mu\nu}\): see kija and kijs, which are defined in preset.f90 .

Parameters
[in]lquaddegree of quadrature
[in]mnnumber of Fourier harmonics
[in]lvolindex of nested volume
[in]lradorder of Chebychev polynomials

References compute_guvijsave(), allglobal::cpus, allglobal::dbdx, allglobal::ddttcc, allglobal::ddttcs, allglobal::ddttsc, allglobal::ddttss, allglobal::ddtzcc, allglobal::ddtzcs, allglobal::ddtzsc, allglobal::ddtzss, allglobal::ddzzcc, allglobal::ddzzcs, allglobal::ddzzsc, allglobal::ddzzss, allglobal::dtoocc, allglobal::dtoocs, allglobal::dtoosc, allglobal::dtooss, allglobal::gaussianabscissae, allglobal::gaussianweight, get_cheby(), get_zernike(), allglobal::goomne, allglobal::goomno, allglobal::gssmne, allglobal::gssmno, allglobal::gstmne, allglobal::gstmno, allglobal::gszmne, allglobal::gszmno, allglobal::gttmne, allglobal::gttmno, allglobal::gtzmne, allglobal::gtzmno, allglobal::gzzmne, allglobal::gzzmno, constants::half, allglobal::im, allglobal::in, allglobal::ki, allglobal::kija, allglobal::kijs, allglobal::lcoordinatesingularity, allglobal::lsavedguvij, ma00aa(), metrix(), allglobal::mne, allglobal::mpi_comm_spec, inputlist::mpol, allglobal::myid, allglobal::ncpu, allglobal::notstellsym, constants::one, fileunits::ounit, constants::pi, constants::pi2, allglobal::regumm, allglobal::tdstcc, allglobal::tdstcs, allglobal::tdstsc, allglobal::tdstss, allglobal::tdszcc, allglobal::tdszcs, allglobal::tdszsc, allglobal::tdszss, allglobal::ttsscc, allglobal::ttsscs, allglobal::ttsssc, allglobal::ttssss, constants::two, volume(), inputlist::wmacros, allglobal::yesstellsym, and constants::zero.

Referenced by dfp100(), get_lu_beltrami_matrices(), and ma00aa().

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

◆ spsint()

subroutine spsint ( integer, intent(in)  lquad,
integer, intent(in)  mn,
integer, intent(in)  lvol,
integer, intent(in)  lrad 
)

Calculates volume integrals of Chebyshev-polynomials and metric elements for preconditioner.

Computes the integrals needed for spsmat.f90. Same as ma00aa.f90, but only compute the relevant terms that are non-zero.

Parameters
lquad
mn
lvol
lrad

References allglobal::cpus, allglobal::ddttcc, allglobal::ddttcs, allglobal::ddttsc, allglobal::ddttss, allglobal::ddtzcc, allglobal::ddtzcs, allglobal::ddtzsc, allglobal::ddtzss, allglobal::ddzzcc, allglobal::ddzzcs, allglobal::ddzzsc, allglobal::ddzzss, allglobal::dtoocc, allglobal::dtoocs, allglobal::dtoosc, allglobal::dtooss, allglobal::gaussianabscissae, allglobal::gaussianweight, get_cheby(), get_zernike(), allglobal::guvijsave, constants::half, allglobal::im, allglobal::in, allglobal::ki, allglobal::kija, allglobal::kijs, allglobal::lcoordinatesingularity, allglobal::mne, allglobal::mpi_comm_spec, inputlist::mpol, allglobal::myid, allglobal::ncpu, allglobal::notstellsym, allglobal::ntz, constants::one, fileunits::ounit, constants::pi, constants::pi2, allglobal::regumm, numerical::small, spsint(), numerical::sqrtmachprec, allglobal::tdstcc, allglobal::tdstcs, allglobal::tdstsc, allglobal::tdstss, allglobal::tdszcc, allglobal::tdszcs, allglobal::tdszsc, allglobal::tdszss, allglobal::ttsscc, allglobal::ttsscs, allglobal::ttsssc, allglobal::ttssss, constants::two, numerical::vsmall, inputlist::wmacros, allglobal::yesstellsym, and constants::zero.

Referenced by dfp100(), and spsint().

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