pyoculus 0.1.1
Loading...
Searching...
No Matches
pyoculus.problems.interpolate_coordinates.SurfacesToroidal Class Reference

Toroidal surfaces object Define surfaces \(s(\vartheta, \zeta) \) and the angle transformations \( \theta (\vartheta, \zeta) \). More...

Classes

class  CoordsOutput
 Output class. More...
 

Public Member Functions

 __init__ (self, nsurfaces=2, mpol=10, ntor=10, Nfp=1, stellar_sym=True)
 Initialise the toroidal surface object.
 
 remove_surface (self, int i)
 Remove a surface.
 
 add_surface (self, float rho, scn, tsn, ssn=None, tcn=None)
 Adding a surface into the system with radial label rho.
 
 replace_surface (self, int idx, float rho=None, scn=None, tsn=None, ssn=None, tcn=None)
 Replacing a surface by the new one.
 
 get_coords (self, sarr=[0], tarr=[0], zarr=[0], derivative=0, input1D=True)
 Compute the coordinates and their derivatives given \(\rho, \vartheta, \zeta\).
 
 jacobi_transform (self, jacobi, CoordsOutput coords)
 Compute the derivatives wrt the new coordinate given the derivatives in the old coordinates.
 
 metric_transform (self, g, CoordsOutput coords)
 Compute the coordinate transformation for a metric tensor \(g_{ij}\) (lower!) (derivative computation under construction)
 
 contra_vector_transform (self, v, CoordsOutput coords, has_jacobian=False, derivative=False, dv=None)
 Compute the coordinate transformation for a contravariant vector \(v^i\) (upper!)
 
 construct_interpolant (self, rhosurfs=None, method="cubic_spline", **kwargs)
 Construct the interpolant between surfaces for the given method.
 
 plot (self, zeta=0, npoints=129, **kwargs)
 Plot the surfaces.
 
 read_surfaces_from_file (self, filename="data.npz", Nfp=None)
 Read the surfaces into a numpy array on disk.
 
 write_surfaces_to_file (self, filename="data.npz")
 Save the surfaces into a numpy array on disk.
 

Public Attributes

int mpol = mpol
 The poloidal mode number.
 
tuple ntor = ntor
 The toroidal mode number.
 
 Nfp = Nfp
 The field period.
 
bool sym = stellar_sym
 Stellarator symmetry.
 
 nsurfaces = nsurfaces
 Number of surfaces.
 
 rhosurfs = np.linspace(0, 1, nsurfaces)
 The rho coordinate for each surface.
 
 scn = np.zeros([nsurfaces, mpol + 1, 2 * ntor + 1])
 The cosine coefficients of \(s\), dimension (#interfaces, mpol, 2*ntor+1)
 
 tsn = self.scn.copy()
 The sine coefficients of \(\vartheta\), dimension (#interfaces, mpol, 2*ntor+1)
 
 ssn = self.tsn.copy()
 The sine coefficients of \(s\), dimension (#interfaces, mpol, 2*ntor+1)
 
 tcn = self.tsn.copy()
 The cosine coefficients of \(\vartheta\), dimension (#interfaces, mpol, 2*ntor+1)
 
str _tcn_int = "pchip":
 
int tsn = 4:
 

Protected Attributes

 _mlist = np.arange(0, self.mpol + 1)
 
tuple _nlist
 
 _scn_int = CubicSpline(rhosurfs, self.scn, axis=0, **kwargs)
 The radial interpolant for s cosine components.
 
 _tsn_int = CubicSpline(rhosurfs, self.tsn, axis=0, **kwargs)
 The radial interpolant for theta sine components.
 
 _ssn_int = CubicSpline(rhosurfs, self.ssn, axis=0, **kwargs)
 The radial interpolant for s sine components.
 
 _tcn_int = CubicSpline(rhosurfs, self.tcn, axis=0, **kwargs)
 The radial interpolant for theta cosine components.
 
 _scn_d1 = self._scn_int.derivative(1)
 The radial interpolant for s cosine components, first derivative.
 
 _scn_d2 = self._scn_int.derivative(2)
 The radial interpolant for s cosine components, second derivative.
 
 _tsn_d1 = self._tsn_int.derivative(1)
 The radial interpolant for theta sine components, first derivative.
 
 _tsn_d2 = self._tsn_int.derivative(2)
 The radial interpolant for theta sine components, second derivative.
 
 _ssn_d1 = self._ssn_int.derivative(1)
 The radial interpolant for s sine components,first derivative.
 
 _ssn_d2 = self._ssn_int.derivative(2)
 The radial interpolant for s sine components, second derivative.
 
 _tcn_d1 = self._tcn_int.derivative(1)
 The radial interpolant for theta cosine components, first derivative.
 
 _tcn_d2 = self._tcn_int.derivative(2)
 The radial interpolant for theta cosine components, second derivative.
 

Detailed Description

Toroidal surfaces object Define surfaces \(s(\vartheta, \zeta) \) and the angle transformations \( \theta (\vartheta, \zeta) \).

\(s, \vartheta\) are written in terms of Fourier harmonics

\[ s(\vartheta, \zeta) = \sum_{m=0}^{\text{mpol}}\sum_{n=-\text{ntor}}^{\text{ntor}} s_{c,m,n} \cos(m \vartheta - n \text{N} \zeta) + s_{s,m,n} \sin(m \vartheta - n \text{N} \zeta), \]

where \(N = \text{Nfp}\).

\[ \theta(\vartheta, \zeta) = \vartheta + \sum_{m=0}^{\text{mpol}}\sum_{n=-\text{ntor}}^{\text{ntor}} \theta_{c,m,n} \cos(m \vartheta - n \text{N} \zeta) + \theta_{s,m,n} \sin(m \vartheta - n \text{N} \zeta). \]

If stellarator symmetry is assumed, then \(s\) will only have cosine components and \(\theta \) will only have sine components.

Each surface is assigned a new radial label \(\rho_i\). Now a new set of coordinates \((\rho, \vartheta, \zeta)\) will be constructed given the relationship \(s(\rho, \vartheta, \zeta), \theta(\rho, \vartheta, \zeta)\). If \(\rho\) is not on one of the known surfaces, an interpolation will be performed on all the Fourier harmonics radially using the method of the user's choice.

Constructor & Destructor Documentation

◆ __init__()

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.__init__ ( self,
nsurfaces = 2,
mpol = 10,
ntor = 10,
Nfp = 1,
stellar_sym = True )

Initialise the toroidal surface object.

Parameters
mpolthe poloidal mode number
ntorthe toroidal mode number
Nfpthe toroidal periodicity
stellar_symassuming stellarator symmetry or not

Member Function Documentation

◆ add_surface()

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.add_surface ( self,
float rho,
scn,
tsn,
ssn = None,
tcn = None )

Adding a surface into the system with radial label rho.

Parameters
rhothe new coordinate \(\rho\) for this new surface
scnthe cosine components of \(s(\rho, \vartheta, \zeta)\)
tsnthe sine components of \(\theta(\rho, \vartheta, \zeta)\)
ssnthe sine componets of \(s(\rho, \vartheta, \zeta)\)
tcnthe cosine components of \(\theta(\rho, \vartheta, \zeta)\)

◆ construct_interpolant()

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.construct_interpolant ( self,
rhosurfs = None,
method = "cubic_spline",
** kwargs )

Construct the interpolant between surfaces for the given method.

Parameters
rhosurfsThe \(\rho\) coordinates for each surface. If input is None, the default or internally stored will be used.
methodThe method being used, can be one of ['cubic_spline', 'cubic_hermite', 'pchip']. The corresponding scipy class is being used.
**kwargsOther parameters passing onto the interpolant constructor

◆ contra_vector_transform()

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.contra_vector_transform ( self,
v,
CoordsOutput coords,
has_jacobian = False,
derivative = False,
dv = None )

Compute the coordinate transformation for a contravariant vector \(v^i\) (upper!)

Parameters
vthe vector \(v^i\) to be transformed, should have the dimension (3,...)
coordsthe output of get_coords, should match the dimension of v
has_jacobianif the given vector already contains the jacobian
derivativeif the derivatives are needed or not. If True, dv is required.
dvthe derivative of v wrt \((s,\theta,\zeta)\), have the dimension (3 derivatives, 3 component,...)
Returns
voutput, dvoutput, the transformed metric and the derivative wrt \((\rho,\vartheta,\zeta)\)

The vector transformation is given by

\[ \left(\begin{array}{c} v^{\rho} \\ v^{\vartheta} \\ v^{\zeta} \end{array}\right) = \mathbf{J}^{-1} \left(\begin{array}{c} v^{{s}} \\ v^{\theta} \\ v^{\zeta} \end{array}\right), \]

where

\[ \mathbf{J} = \frac{\partial(s,\theta,\zeta)}{\partial(\rho,\vartheta,\zeta)} \]

is the Jacobi matrix. For computing the result we use matrix solve instead of matrix inversion.

When the derivative is needed, we note that

\[ \frac{\partial }{\partial \rho}\mathbf{J} \left(\begin{array}{c} v^{\rho} \\ v^{\vartheta} \\ v^{\zeta} \end{array}\right) +\mathbf{J} \frac{\partial }{\partial \rho}\left(\begin{array}{c} v^{\rho} \\ v^{\vartheta} \\ v^{\zeta} \end{array}\right) = \frac{\partial }{\partial \rho}\left(\begin{array}{c} v^{{s}} \\ v^{\theta} \\ v^{\zeta} \end{array}\right), \]

therefore

\[ \frac{\partial}{\partial \rho}\left(\begin{array}{c} v^{\rho} \\ v^{\vartheta} \\ v^{\zeta} \end{array}\right) = \mathbf{J}^{-1} \left[\frac{\partial }{\partial \rho}\left(\begin{array}{c} v^{{s}} \\ v^{\theta} \\ v^{\zeta} \end{array}\right) -\frac{\partial \mathbf{J}}{\partial\rho} \left(\begin{array}{c} v^{\rho} \\ v^{\vartheta} \\ v^{\zeta} \end{array}\right)\right]. \]

Similarly we can compute the \(\vartheta\) and \(\zeta\) derivatives.

Now,

\[ \frac{\partial (v^{{s}}, v^{\theta}, v^{\zeta})}{\partial(\rho, \vartheta, \zeta)} = \frac{\partial (v^{{s}}, v^{\theta}, v^{\zeta})}{\partial(s, \theta, \zeta)} \times \frac{\partial(s, \theta, \zeta)}{\partial(\rho, \vartheta, \zeta)} \]

Summing up

\[ \frac{\partial ( v^{\rho}, v^{\vartheta}, v^{\zeta})}{\partial(\rho, \vartheta, \zeta)} = \mathbf{J}^{-1} \cdot \frac{\partial (v^{{s}}, v^{\theta}, v^{\zeta})}{\partial(s, \theta, \zeta)} \cdot \mathbf{J} - \mathbf{J}^{-1} \cdot (\mathbf{J}_\rho, \mathbf{J}_\vartheta, \mathbf{J}_\zeta) \cdot \left(\begin{array}{c} v^{\rho} \\ v^{\vartheta} \\ v^{\zeta} \end{array}\right) \]

◆ get_coords()

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.get_coords ( self,
sarr = [0],
tarr = [0],
zarr = [0],
derivative = 0,
input1D = True )

Compute the coordinates and their derivatives given \(\rho, \vartheta, \zeta\).

Parameters
sarr\(\rho\), an array
tarr\(\vartheta\), an array
zarr\(\zeta\), an array
derivative0:only return \(s, \theta\), 1: also return the derivatives, 2: also return the second derivatives
input1Dif False, create a meshgrid with sarr, tarr and zarr, if True, treat them as a list of points
Returns
coords coords.s, coords.t contains the coordinates \(s,\theta\). coords.ds and dt contains the first derivatives. coords.dds and ddt contains second derivates. coords.jacobi contains the jacobi matrix \(J = \partial(s,\theta,\zeta)/\partial(\rho,\vartheta,\zeta)\), coords.djacobi contains the derivative.

◆ jacobi_transform()

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.jacobi_transform ( self,
jacobi,
CoordsOutput coords )

Compute the derivatives wrt the new coordinate given the derivatives in the old coordinates.

Parameters
jacobithe Jacobi matrix \(\partial(f_1, \cdots, f_N) / \partial (s,\theta,\zeta)\)
coordsthe output of get_coords, should match the dimension of J
Returns
jacobi_output the transformed Jacobi matrix and the derivatives wrt \((\rho,\vartheta,\zeta)\)

Given \(N\) variables \(f_1, \cdots, f_N \), the jacobi matrix (derivatives wrt the new coordinate is given by the chain rule

\[ \frac{\partial(f_1, \cdots, f_N)}{\partial(\rho, \vartheta, \zeta)} = \frac{\partial(f_1, \cdots, f_N)}{\partial(s, \theta, \zeta)} \cdot \frac{\partial(s, \theta, \zeta)}{\partial(\rho, \vartheta, \zeta)} \]

◆ metric_transform()

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.metric_transform ( self,
g,
CoordsOutput coords )

Compute the coordinate transformation for a metric tensor \(g_{ij}\) (lower!) (derivative computation under construction)

Parameters
gthe metric \(g_{ij}\) of the old coordinate \((s,\theta,\zeta)\) to be transformed, should have the dimension (3,3,...)
coordsthe output of get_coords, should match the dimension of g
Returns
goutput, [dgoutput], the transformed metric and the derivative wrt \((\rho,\vartheta,\zeta)\)

Let the input metric for \((s, \theta, \zeta)\) being \(\mathbf{g}\), the metric for \((\rho, \vartheta, \zeta)\) is given by

\[ \mathbf{G} = \mathbf{J}^T \mathbf{g} \mathbf{J}, \]

where

\[ \mathbf{J} = \frac{\partial(s,\theta,\zeta)}{\partial(\rho,\vartheta,\zeta)} \]

is the Jacobi matrix.

◆ plot()

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.plot ( self,
zeta = 0,
npoints = 129,
** kwargs )

Plot the surfaces.

Parameters
zetathe toroidal plane
npointsthe number of points in the plot
**kwargsall other parameters going into plt.plot

◆ read_surfaces_from_file()

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.read_surfaces_from_file ( self,
filename = "data.npz",
Nfp = None )

Read the surfaces into a numpy array on disk.

Parameters
filenamethe filename to load from
Nfpthe toroidal periodicity, if known.

◆ remove_surface()

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.remove_surface ( self,
int i )

Remove a surface.

Parameters
ithe id of the surface

◆ replace_surface()

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.replace_surface ( self,
int idx,
float rho = None,
scn = None,
tsn = None,
ssn = None,
tcn = None )

Replacing a surface by the new one.

Parameters
idxthe index of the surface to be replaced
rhothe new coordinate \(\rho\) for this new surface. If this is None then keep same.
scnthe cosine components of \(s(\rho, \vartheta, \zeta)\). If this is None then keep same.
tsnthe sine components of \(\theta(\rho, \vartheta, \zeta)\). If this is None then keep same.
ssnthe sine componets of \(s(\rho, \vartheta, \zeta)\). If this is None then keep same.
tcnthe cosine components of \(\theta(\rho, \vartheta, \zeta)\). If this is None then keep same.

◆ write_surfaces_to_file()

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.write_surfaces_to_file ( self,
filename = "data.npz" )

Save the surfaces into a numpy array on disk.

Parameters
filenamethe filename to save to

Member Data Documentation

◆ _mlist

pyoculus.problems.interpolate_coordinates.SurfacesToroidal._mlist = np.arange(0, self.mpol + 1)
protected

◆ _nlist

tuple pyoculus.problems.interpolate_coordinates.SurfacesToroidal._nlist
protected
Initial value:
= (
np.concatenate([np.arange(0, self.ntor + 1), np.arange(-self.ntor, 0)])
* self.Nfp
)

◆ _scn_d1

pyoculus.problems.interpolate_coordinates.SurfacesToroidal._scn_d1 = self._scn_int.derivative(1)
protected

The radial interpolant for s cosine components, first derivative.

◆ _scn_d2

pyoculus.problems.interpolate_coordinates.SurfacesToroidal._scn_d2 = self._scn_int.derivative(2)
protected

The radial interpolant for s cosine components, second derivative.

◆ _scn_int

pyoculus.problems.interpolate_coordinates.SurfacesToroidal._scn_int = CubicSpline(rhosurfs, self.scn, axis=0, **kwargs)
protected

The radial interpolant for s cosine components.

◆ _ssn_d1

pyoculus.problems.interpolate_coordinates.SurfacesToroidal._ssn_d1 = self._ssn_int.derivative(1)
protected

The radial interpolant for s sine components,first derivative.

◆ _ssn_d2

pyoculus.problems.interpolate_coordinates.SurfacesToroidal._ssn_d2 = self._ssn_int.derivative(2)
protected

The radial interpolant for s sine components, second derivative.

◆ _ssn_int

pyoculus.problems.interpolate_coordinates.SurfacesToroidal._ssn_int = CubicSpline(rhosurfs, self.ssn, axis=0, **kwargs)
protected

The radial interpolant for s sine components.

◆ _tcn_d1

pyoculus.problems.interpolate_coordinates.SurfacesToroidal._tcn_d1 = self._tcn_int.derivative(1)
protected

The radial interpolant for theta cosine components, first derivative.

◆ _tcn_d2

pyoculus.problems.interpolate_coordinates.SurfacesToroidal._tcn_d2 = self._tcn_int.derivative(2)
protected

The radial interpolant for theta cosine components, second derivative.

◆ _tcn_int [1/2]

pyoculus.problems.interpolate_coordinates.SurfacesToroidal._tcn_int = CubicSpline(rhosurfs, self.tcn, axis=0, **kwargs)
protected

The radial interpolant for theta cosine components.

◆ _tcn_int [2/2]

str pyoculus.problems.interpolate_coordinates.SurfacesToroidal._tcn_int = "pchip":

◆ _tsn_d1

pyoculus.problems.interpolate_coordinates.SurfacesToroidal._tsn_d1 = self._tsn_int.derivative(1)
protected

The radial interpolant for theta sine components, first derivative.

◆ _tsn_d2

pyoculus.problems.interpolate_coordinates.SurfacesToroidal._tsn_d2 = self._tsn_int.derivative(2)
protected

The radial interpolant for theta sine components, second derivative.

◆ _tsn_int

pyoculus.problems.interpolate_coordinates.SurfacesToroidal._tsn_int = CubicSpline(rhosurfs, self.tsn, axis=0, **kwargs)
protected

The radial interpolant for theta sine components.

◆ mpol

int pyoculus.problems.interpolate_coordinates.SurfacesToroidal.mpol = mpol

The poloidal mode number.

◆ Nfp

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.Nfp = Nfp

The field period.

◆ nsurfaces

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.nsurfaces = nsurfaces

Number of surfaces.

◆ ntor

tuple pyoculus.problems.interpolate_coordinates.SurfacesToroidal.ntor = ntor

The toroidal mode number.

◆ rhosurfs

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.rhosurfs = np.linspace(0, 1, nsurfaces)

The rho coordinate for each surface.

◆ scn

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.scn = np.zeros([nsurfaces, mpol + 1, 2 * ntor + 1])

The cosine coefficients of \(s\), dimension (#interfaces, mpol, 2*ntor+1)

◆ ssn

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.ssn = self.tsn.copy()

The sine coefficients of \(s\), dimension (#interfaces, mpol, 2*ntor+1)

◆ sym

bool pyoculus.problems.interpolate_coordinates.SurfacesToroidal.sym = stellar_sym

Stellarator symmetry.

◆ tcn

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.tcn = self.tsn.copy()

The cosine coefficients of \(\vartheta\), dimension (#interfaces, mpol, 2*ntor+1)

◆ tsn [1/2]

pyoculus.problems.interpolate_coordinates.SurfacesToroidal.tsn = self.scn.copy()

The sine coefficients of \(\vartheta\), dimension (#interfaces, mpol, 2*ntor+1)

◆ tsn [2/2]

int pyoculus.problems.interpolate_coordinates.SurfacesToroidal.tsn = 4:

The documentation for this class was generated from the following file: