SPEC 3.20
Stepped Pressure Equilibrium Code
|
Functions/Subroutines | |
subroutine | bnorml (mn, Ntz, efmn, ofmn) |
Computes \({\bf B}_{Plasma} \cdot {\bf e}_\theta \times {\bf e}_\zeta \;\) on the computational boundary, \(\partial {\cal D}\). More... | |
subroutine | casing (teta, zeta, gBn, icasing) |
Constructs the field created by the plasma currents, at an arbitrary, external location using virtual casing. More... | |
subroutine | dvcfield (Ndim, tz, Nfun, vcintegrand) |
Differential virtual casing integrand. More... | |
subroutine bnorml | ( | integer, intent(in) | mn, |
integer, intent(in) | Ntz, | ||
real, dimension(1:mn), intent(out) | efmn, | ||
real, dimension(1:mn), intent(out) | ofmn | ||
) |
Computes \({\bf B}_{Plasma} \cdot {\bf e}_\theta \times {\bf e}_\zeta \;\) on the computational boundary, \(\partial {\cal D}\).
free-boundary constraint
compute the normal field on a regular grid on the computational boundary
[in] | mn | total number of Fourier harmonics |
[in] | Ntz | total number of grid points in \(\theta\) and \(zeta\) |
[out] | efmn | even Fourier coefficients |
[out] | ofmn | odd Fouier coefficients |
References allglobal::ate, allglobal::ato, allglobal::aze, allglobal::azo, bnorml(), casing(), allglobal::cfmn, allglobal::cpus, allglobal::dxyz, allglobal::globaljk, allglobal::gteta, allglobal::guvij, allglobal::gzeta, constants::half, inputlist::igeometry, allglobal::ijimag, allglobal::ijreal, allglobal::im, allglobal::in, allglobal::jiimag, allglobal::jireal, inputlist::lcheck, allglobal::lcoordinatesingularity, inputlist::lrad, fileunits::lunit, allglobal::mpi_comm_spec, allglobal::myid, allglobal::ncpu, allglobal::notstellsym, allglobal::nt, allglobal::nxyz, allglobal::nz, constants::one, fileunits::ounit, constants::pi, constants::pi2, allglobal::rij, allglobal::sfmn, allglobal::sg, numerical::small, constants::ten, allglobal::tetazeta, allglobal::tt, constants::two, inputlist::vcasingper, inputlist::vcasingtol, allglobal::virtualcasingfactor, inputlist::wmacros, constants::zero, and allglobal::zij.
Referenced by bnorml(), preset(), and spec().
subroutine casing | ( | real, intent(in) | teta, |
real, intent(in) | zeta, | ||
real, intent(out) | gBn, | ||
integer, intent(out) | icasing | ||
) |
Constructs the field created by the plasma currents, at an arbitrary, external location using virtual casing.
Compute the external magnetic field using virtual casing.
Theory and numerics
\begin{eqnarray} {\bf x}(\theta,\zeta) \equiv x(\theta,\zeta) {\bf i} + y(\theta,\zeta) {\bf j} + z(\theta,\zeta) {\bf k}, \end{eqnarray}
and the tangential field on this boundary,\begin{eqnarray} {\bf B}_s=B^\theta {\bf e}_\theta + B^\zeta {\bf e}_\zeta, \end{eqnarray}
where \(\theta\) and \(\zeta\) are arbitrary poloidal and toroidal angles, and \({\bf e}_\theta \equiv \partial {\bf x}/\partial \theta\), \({\bf e}_\zeta \equiv \partial {\bf x}/\partial \zeta\). This routine assumes that the plasma boundary is a flux surface, i.e. \({\bf B} \cdot {\bf e}_\theta \times {\bf e}_\zeta = 0\).\begin{eqnarray} {\bf j} = {\bf B}_s \times {\bf n}, \end{eqnarray}
where \({\bf n}\) is normal to the surface.\begin{eqnarray} {\bf B}({\bf \bar x}) = -\frac{1}{4\pi} \int_{\cal S} \frac{\left( {\bf B}_s \times d{\bf s} \right) \times {\bf \hat r}}{r^2}, \end{eqnarray}
where \(d{\bf s} \equiv {\bf e}_\theta \times {\bf e}_\zeta \; d\theta d\zeta\).\begin{eqnarray} {\bf J} & \equiv & {\bf B}_s \times d{\bf s} = \alpha \; {\bf e}_\theta - \beta \; {\bf e}_\zeta, \end{eqnarray}
where \(\alpha \equiv B_\zeta = B^\theta g_{\theta\zeta} + B^\zeta g_{\zeta\zeta}\) and \(\beta \equiv B_\theta = B^\theta g_{\theta\theta} + B^\zeta g_{\theta\zeta}\).\begin{eqnarray} j_x & = & \alpha \; x_\theta - \beta \; x_\zeta \\ j_y & = & \alpha \; y_\theta - \beta \; y_\zeta \\ j_z & = & \alpha \; z_\theta - \beta \; z_\zeta . \end{eqnarray}
\begin{eqnarray} {\bf j} & \equiv & \nabla \times {\bf B} = {\sqrt g}^{-1}(\partial_\theta B_\zeta -\partial_\zeta B_\theta) \; {\bf e}_s + {\sqrt g}^{-1}(\partial_\zeta B_s -\partial_s B_\zeta ) \; {\bf e}_\theta + {\sqrt g}^{-1}(\partial_s B_\theta-\partial_\theta B_s ) \; {\bf e}_\zeta, \end{eqnarray}
has no normal component to the surface, i.e. \({\bf j}\cdot \nabla s=0\), we obtain the condition \(\partial_\theta B_\zeta = \partial_\zeta B_\theta\), or \(\partial_\theta \alpha = \partial_\zeta \beta\). In axisymmetric configurations, where \(\partial_\zeta \beta=0\), we must have \(\partial_\theta \alpha=0\).\begin{eqnarray} {\bf r} \equiv r_x \; {\bf i} + r_y \; {\bf j} + r_z \; {\bf k} = (X-x) \; {\bf i} + (Y-y) \; {\bf j} + (Z-z) \; {\bf k}. \end{eqnarray}
\begin{eqnarray} B^x &=& \oint\!\!\!\oint \!d\theta d\zeta \,\,\, (j_y r_z - j_z r_y)/r^3,\\ B^y &=& \oint\!\!\!\oint \!d\theta d\zeta \,\,\, (j_z r_x - j_x r_z)/r^3,\\ B^z &=& \oint\!\!\!\oint \!d\theta d\zeta \,\,\, (j_x r_y - j_y r_x)/r^3 \end{eqnarray}
up to a scaling factorvirtualcasingfactor
\(=-1/4\pi\) that is taken into account at the end. \begin{eqnarray} \int_0^{2\pi} \int_0^{2\pi} \verb+vcintegrand+ \;\; d\theta d\zeta \end{eqnarray}
for a given \((X,Y,Z)\), wherevcintegrand
is given in Eqn. \((\ref{eq:integrand_casing})\). DCUHRE
, which uses an adaptive subdivision strategy and also computes absolute error estimates. The absolute and relative accuracy required are provided by the inputvar
vcasingtol
. The minimum number of function evaluations is provided by the inputvar
vcasingits
. Calculation of integrand
An adaptive integration is used to compute the integrals. Consequently, the magnetic field tangential to the plasma boundary is required at an arbitrary point. This is computed, as always, from \({\bf B} = \nabla \times {\bf A}\), and this provides \({\bf B} = B^\theta {\bf e}_\theta + B^\zeta {\bf e}_\zeta\). Recall that \(B^s=0\) by construction on the plasma boundary.
Then, the metric elements \(g_{\theta\theta}\), \(g_{\theta\zeta}\) and \(g_{\zeta\zeta}\) are computed. These are used to "lower" the components of the magnetic field, \({\bf B} = B_\theta \nabla \theta + B_\zeta \nabla \zeta\).
inputvar
vcasingeps
. Let the "distance" be \begin{eqnarray} D \equiv \sqrt{(X-x)^2 + (Y-y)^2 + (Z-Z)^2} + \epsilon^2. \end{eqnarray}
On taking the limit that \(\epsilon \rightarrow 0\), the virtual casing integrand is
\begin{eqnarray} \verb+vcintegrand+ \equiv ( B_x n_x + B_y n_y + B_z n_z ) ( 1 + 3 \epsilon^2 / D^2 ) / D^3, \label{eq:integrand_casing} \end{eqnarray}
where the normal vector is \({\bf n} \equiv n_x {\bf i} + n_y {\bf j} + n_z {\bf k}\). The normal vector, Nxyz
, to the computational boundary (which does not change) is computed in preset().
[in] | teta | \(\theta\) |
[in] | zeta | \(\zeta\) |
[out] | gBn | \( \sqrt g {\bf B} \cdot {\bf n}\) |
[out] | icasing | return flag from dcuhre() |
References casing(), allglobal::cpus, dvcfield(), allglobal::dxyz, allglobal::globaljk, allglobal::mpi_comm_spec, allglobal::myid, allglobal::ncpu, allglobal::nxyz, fileunits::ounit, constants::pi, constants::pi2, inputlist::vcasingits, inputlist::vcasingper, inputlist::vcasingtol, fileunits::vunit, inputlist::wmacros, and constants::zero.
Referenced by bnorml(), casing(), and dvcfield().
subroutine dvcfield | ( | integer, intent(in) | Ndim, |
real, dimension(1:ndim), intent(in) | tz, | ||
integer, intent(in) | Nfun, | ||
real, dimension(1:nfun), intent(out) | vcintegrand | ||
) |
Differential virtual casing integrand.
Differential virtual casing integrand
[in] | Ndim | number of parameters (==2) |
[in] | tz | \(\theta\) and \(\zeta\) |
[in] | Nfun | number of function values (==3) |
[out] | vcintegrand | cartesian components of magnetic field |
References allglobal::ate, allglobal::ato, allglobal::aze, allglobal::azo, casing(), allglobal::cpus, allglobal::dxyz, allglobal::first_free_bound, constants::four, allglobal::globaljk, constants::half, inputlist::igeometry, allglobal::im, allglobal::in, allglobal::irbc, allglobal::irbs, allglobal::izbc, allglobal::izbs, inputlist::lrad, allglobal::mn, allglobal::mpi_comm_spec, allglobal::myid, allglobal::ncpu, allglobal::notstellsym, inputlist::nvol, allglobal::nxyz, constants::one, fileunits::ounit, numerical::small, constants::three, allglobal::tt, inputlist::vcasingeps, fileunits::vunit, allglobal::yesstellsym, and constants::zero.
Referenced by casing().