pyoculus 0.1.1
|
Public Member Functions | |
__init__ (self, ToroidalBfield problem, params=dict(), integrator=None, integrator_params=dict()) | |
Set up the class of the fixed point finder. | |
construct_qfms (self, plist, qlist, sguesslist=None, bounding_surfaces=[0.0, 1.0], rho_label="s", niter=5, verbose=True) | |
Construct the QFM surfaces for given list of p and q This subroutine will iteratively construct the QFM surfaces for the given list of p and q. | |
compute_tflux () | |
straighten_boundary (self, surfaces, tol=1e-9, niter=10) | |
Convert a boundary surface to have straight field line For a boundary surface with rho=constant, find the function \(\lambda(\vartheta, \zeta)\) with transformation. | |
action (self, int pp, int qq, sguess=0.5, root_method="hybr", tol=1e-8) | |
Construct a QFM surface based on Stuart's method, for a given orbit periodicity pp and qq. | |
action_gradient (self, xx, pp, qq, a, qN, Nfft) | |
Computes the action gradient, being used in root finding. | |
action_gradient_jacobi (self, xx, pp, qq, a, qN, Nfft) | |
Computes the jacobi matrix of action gradient, being used in root finding. | |
Public Member Functions inherited from pyoculus.solvers.base_solver.BaseSolver | |
is_successful (self) | |
Returns True if the computation is successfully completed. | |
Public Attributes | |
nfft_multiplier = params["nfft_multiplier"] | |
Nfp = problem.Nfp | |
sym = params["stellar_sym"] | |
Public Attributes inherited from pyoculus.solvers.base_solver.BaseSolver | |
bool | successful = False |
flagging if the computation is done and successful | |
Protected Member Functions | |
_straighten_boundary (self, rho=1, tol=1e-9, niter=10) | |
Convert a boundary surface to have straight field line, internal function For a boundary surface with rho=constant, find the function \(\lambda(\vartheta, \zeta)\) with transformation. | |
_unpack_dof (self, xx, qN) | |
Unpack the degrees of freedom into Fourier harmonics. | |
_pack_dof (self, nv, rcn, tsn, rsn, tcn) | |
Unpack the degrees of freedom into Fourier harmonics. | |
Protected Attributes | |
int | _MM = params["nfft_multiplier"] * 2 |
_pqNtor = params["pqNtor"] | |
_pqMpol = params["pqMpol"] | |
_dz = dz | |
_nlist = np.arange(0, qN + 1) | |
_zeta = np.arange(0, Nfft) * dz | |
_nzq = self._nlist[:, nax] * self._zeta[nax, :] / qq | |
_cnzq = np.cos(self._nzq) | |
_snzq = np.sin(self._nzq) | |
Protected Attributes inherited from pyoculus.solvers.base_solver.BaseSolver | |
_integrator_type = RKIntegrator | |
_params = dict(params) | |
_integrator = self._integrator_type(integrator_params) | |
_problem = problem | |
_integrator_params = dict(integrator_params) | |
pyoculus.solvers.qfm.QFM.__init__ | ( | self, | |
ToroidalBfield | problem, | ||
params = dict(), | |||
integrator = None, | |||
integrator_params = dict() ) |
Set up the class of the fixed point finder.
problem | must inherit pyoculus.problems.ToroidalBfield, the problem to solve |
params | dict, the parameters for the solver |
integrator | the integrator to use, must inherit \pyoculus.integrators.BaseIntegrator, if set to None by default using RKIntegrator (not used here) |
integrator_params | dict, the parmaters passed to the integrator (not used here) |
params['pqMpol']=8
– Fourier resolution multiplier for poloidal direction params['pqNtor']=4
– Fourier resolution multiplier for toroidal direction params['nfft_multiplier']=4
– the extended (multiplier) resolution for FFT params['stellar_sym']=True
– stellarator symmetry
Reimplemented from pyoculus.solvers.base_solver.BaseSolver.
|
protected |
Unpack the degrees of freedom into Fourier harmonics.
nv | |
rcn | |
tsn | |
rsn | |
tcn |
|
protected |
Convert a boundary surface to have straight field line, internal function For a boundary surface with rho=constant, find the function \(\lambda(\vartheta, \zeta)\) with transformation.
rho | the boundary surface coordinate |
tol | the tolerance to stop iteration |
niter | the number of iterations |
|
protected |
Unpack the degrees of freedom into Fourier harmonics.
xx | the packed degrees of freedom |
qN | Fourier resolution |
pyoculus.solvers.qfm.QFM.action | ( | self, | |
int | pp, | ||
int | qq, | ||
sguess = 0.5, | |||
root_method = "hybr", | |||
tol = 1e-8 ) |
Construct a QFM surface based on Stuart's method, for a given orbit periodicity pp and qq.
The surface in the old coordinate \((s, \theta, \zeta)\) will be expressed in Fourier harmonics of the reverse mapping
\[ s(\vartheta, \zeta) = \sum_{m,n} s_{c, m, n} \cos(m \vartheta - n \zeta) + s_{s, m, n} \sin(m \vartheta - n \zeta), \]
\[ \theta(\vartheta, \zeta) = \vartheta + \sum_{m,n} t_{c, m, n} \cos(m \vartheta - n \zeta) + t_{s, m, n} \sin(m \vartheta - n \zeta), \]
pp | the orbit peroidicity. The rotational transform of such a surface will be pp/qq. |
the orbit peroidicity. The rotational transform of such a surface will be pp/qq. | |
sguess | a guess of the s coordinate for the surfaces. |
root_method | root finding method being used in scipy.optimize.root |
tol | the tolarence of root finding |
pyoculus.solvers.qfm.QFM.action_gradient | ( | self, | |
xx, | |||
pp, | |||
qq, | |||
a, | |||
qN, | |||
Nfft ) |
Computes the action gradient, being used in root finding.
xx | the packed degrees of freedom. It should contain rcn, tsn, rsn, tcn, nv. |
pp | the poloidal periodicity of the island, should be an integer |
the toroidal periodicity of the island, should be an integer | |
a | the target area |
\[ B^\t / B^\z & = & f^c_0 + \sum_{n=1}^{qN} \left[ f^c_n \cos(n\z/q) + f^s_n \sin(n\z/q) \right], \label{eqn:f} \]
\[ B^\rho / B^\z + \bar \nu / \sqrt g B^\z & = & g^c_0 + \sum_{n=1}^{qN} \left[ g^c_n \cos(n\z/q) + g^s_n \sin(n\z/q) \right], \label{eqn:g} \]
\item The Fourier harmonics of $
pyoculus.solvers.qfm.QFM.action_gradient_jacobi | ( | self, | |
xx, | |||
pp, | |||
qq, | |||
a, | |||
qN, | |||
Nfft ) |
Computes the jacobi matrix of action gradient, being used in root finding.
xx | the packed degrees of freedom. It should contain rcn, tsn, rsn, tcn, nv. |
pp | the poloidal periodicity of the island, should be an integer |
the toroidal periodicity of the island, should be an integer | |
a | the target area |
pyoculus.solvers.qfm.QFM.compute_tflux | ( | ) |
pyoculus.solvers.qfm.QFM.construct_qfms | ( | self, | |
plist, | |||
qlist, | |||
sguesslist = None, | |||
bounding_surfaces = [0.0, 1.0], | |||
rho_label = "s", | |||
niter = 5, | |||
verbose = True ) |
Construct the QFM surfaces for given list of p and q This subroutine will iteratively construct the QFM surfaces for the given list of p and q.
If the system has boundaries with constants s that are flux surfaces, two bounding surfaces will be added. Otherwise, the outmost two QFM surfaces will be the new bounding surfaces.
plist | the list of p |
qlist | the list of q |
sguesslist | guess of the approximate s coordinate of each surface |
bounding_surfaces | instruct to add the bounding surfaces for these s coordinates, maximum two |
rho_label | what to use for the new radius label \(\rho\). Can be 's', 'tflux', 'sqrttflux'. |
niter | the number of iterations in constructing QFM based of straight field line coordinates |
verbose | if True, print the progress |
pyoculus.solvers.qfm.QFM.straighten_boundary | ( | self, | |
surfaces, | |||
tol = 1e-9, | |||
niter = 10 ) |
Convert a boundary surface to have straight field line For a boundary surface with rho=constant, find the function \(\lambda(\vartheta, \zeta)\) with transformation.
\[ \theta = \vartheta + \lambda(\vartheta, \zeta), \]
such that \(B^{\vartheta} / B^{\zeta} = 1 / q\) is a constant on the surface.
The transformation gives
\[ B^{\vartheta} = \frac{B^\theta - \lambda_\zeta B^\zeta}{1 + \lambda_\vartheta}, \]
and no change to \(B^\zeta\). Therefore, we get
\[ \frac{B^\theta}{B^\zeta}(\theta = \vartheta + \lambda(\vartheta, \zeta), \zeta) = \frac{1}{q} (1 + \lambda_\vartheta) + \lambda_\zeta \]
A fixed-point iteration method will be used. For a initial guess of \(\lambda \), we use it to compute the left-hand-side. An updated \(\lambda \) is then computed by solving the right-hand-side assuming the left-hand-side is given. This new \(\lambda \) is again substituted into the left hand side to generate a new right-hand-side. The process is repeated until converged (the difference between two iterations is lower than a threshold, or a given number of iteration is reached).
surfaces | class SurfacesToroidal, the surfaces for which the first and last surface are bounding surfaces |
tol | the tolerance to stop iteration |
niter | the number of iterations |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
pyoculus.solvers.qfm.QFM.nfft_multiplier = params["nfft_multiplier"] |
pyoculus.solvers.qfm.QFM.Nfp = problem.Nfp |
pyoculus.solvers.qfm.QFM.sym = params["stellar_sym"] |