Evaluation kernels#

Module containing accelerated (pyccelized) functions for evaluating metric coefficients corresponding to mappings (x, y, z) = F(eta_1, eta_2, eta_3).

struphy.geometry.evaluation_kernels.f(eta1: float, eta2: float, eta3: float, args: DomainArguments, f_out: float[:])[source]#

Point-wise evaluation of (x, y, z) = F(eta1, eta2, eta3).

Parameters:
  • eta1 (float) – Position on the unit cube.

  • eta2 (float) – Position on the unit cube.

  • eta3 (float) – Position on the unit cube.

  • args (DomainArguments) – Arguments for the mapping.

  • f_out (xp.array) – Output array of shape (3,).

struphy.geometry.evaluation_kernels.df(eta1: float, eta2: float, eta3: float, args: DomainArguments, df_out: float[:, :])[source]#

Point-wise evaluation of the Jacobian matrix DF = (dF_i/deta_j)_(i,j=1,2,3).

Parameters:
  • eta1 (float) – Position on the unit cube.

  • eta2 (float) – Position on the unit cube.

  • eta3 (float) – Position on the unit cube.

  • args (DomainArguments) – Arguments for the mapping.

  • df_out (xp.array) – Output array of shape (3, 3).

struphy.geometry.evaluation_kernels.det_df(eta1: float, eta2: float, eta3: float, args: DomainArguments, tmp1: float[:, :]) float[source]#

Point-wise evaluation of the Jacobian determinant det(dF) = dF/deta1.dot(dF/deta2 x dF/deta3).

Parameters:
  • eta1 (float) – Position on the unit cube.

  • eta2 (float) – Position on the unit cube.

  • eta3 (float) – Position on the unit cube.

  • args (DomainArguments) – Arguments for the mapping.

  • tmp1 (xp.array) – Temporary array of shape (3, 3).

struphy.geometry.evaluation_kernels.df_inv(eta1: float, eta2: float, eta3: float, args: DomainArguments, tmp1: float[:, :], avoid_round_off: bool, dfinv_out: float[:, :])[source]#

Point-wise evaluation of the inverse Jacobian matrix DF^(-1).

Parameters:
  • eta1 (float) – Position on the unit cube.

  • eta2 (float) – Position on the unit cube.

  • eta3 (float) – Position on the unit cube.

  • args (DomainArguments) – Arguments for the mapping.

  • tmp1 (xp.array) – Temporary array of shape (3, 3).

  • avoid_round_off (bool) – Whether to manually set exact zeros in arrays.

  • dfinv_out (xp.array) – Output array of shape (3, 3).

struphy.geometry.evaluation_kernels.g(eta1: float, eta2: float, eta3: float, args: DomainArguments, tmp1: float[:, :], tmp2: float[:, :], avoid_round_off: bool, g_out: float[:, :])[source]#

Point-wise evaluation of the metric tensor G = DF^T * DF.

Parameters:
  • eta1 (float) – Position on the unit cube.

  • eta2 (float) – Position on the unit cube.

  • eta3 (float) – Position on the unit cube.

  • args (DomainArguments) – Arguments for the mapping.

  • tmp1 (xp.array) – Temporary arrays of shape (3, 3).

  • tmp2 (xp.array) – Temporary arrays of shape (3, 3).

  • avoid_round_off (bool) – Whether to manually set exact zeros in arrays.

  • g_out (xp.array) – Output array of shape (3, 3).

struphy.geometry.evaluation_kernels.g_inv(eta1: float, eta2: float, eta3: float, args: DomainArguments, tmp1: float[:, :], tmp2: float[:, :], tmp3: float[:, :], avoid_round_off: bool, ginv_out: float[:, :])[source]#

Point-wise evaluation of the inverse metric tensor G^(-1) = DF^(-1) * DF^(-T).

Parameters:
  • eta1 (float) – Position on the unit cube.

  • eta2 (float) – Position on the unit cube.

  • eta3 (float) – Position on the unit cube.

  • args (DomainArguments) – Arguments for the mapping.

  • tmp1 (xp.array) – Temporary arrays of shape (3, 3).

  • tmp2 (xp.array) – Temporary arrays of shape (3, 3).

  • tmp3 (xp.array) – Temporary arrays of shape (3, 3).

  • avoid_round_off (bool) – Whether to manually set exact zeros in arrays.

  • ginv_out (xp.array) – Output array of shape (3, 3).

struphy.geometry.evaluation_kernels.select_metric_coeff(eta1: float, eta2: float, eta3: float, kind_coeff: int, args: DomainArguments, tmp0: float[:], tmp1: float[:, :], tmp2: float[:, :], tmp3: float[:, :], avoid_round_off: bool, out: float[:, :])[source]#

Point-wise evaluation of metric coefficients.

Parameters:
  • eta1 (float) – Position on the unit cube.

  • eta2 (float) – Position on the unit cube.

  • eta3 (float) – Position on the unit cube.

  • kind_coeff (int) –

    Which metric coefficient to evaluate.
    • -1 : identity

    • 0 : mapping

    • 1 : Jacobian matrix

    • 2 : Jacobian determinant

    • 3 : inverse Jacobian matrix

    • 4 : metric tensor

    • 5 : inverse metric tensor

  • args (DomainArguments) – Arguments for the mapping.

  • tmp0 (xp.array) – Temporary array of shape (3,).

  • tmp1 (xp.array) – Temporary arrays of shape (3, 3).

  • tmp2 (xp.array) – Temporary arrays of shape (3, 3).

  • tmp3 (xp.array) – Temporary arrays of shape (3, 3).

  • avoid_round_off (bool) – Whether to manually set exact zeros in arrays.

  • out (xp.array) – Output array of shape (3, 3).

struphy.geometry.evaluation_kernels.kernel_evaluate(eta1: float[:, :, :], eta2: float[:, :, :], eta3: float[:, :, :], kind_coeff: int, args: DomainArguments, mat_f: float[:, :, :, :, :], is_sparse_meshgrid: bool, avoid_round_off: bool)[source]#

Evaluation of metric coefficients on a given 3d grid of evaluation points.

Parameters:

is_sparse_meshgrid (bool) – Whether the 3d evaluation points were obtained from a sparse meshgrid.

struphy.geometry.evaluation_kernels.kernel_evaluate_pic(markers: float[:, :], kind_coeff: int, args: DomainArguments, mat_f: float[:, :, :], remove_outside: bool, avoid_round_off: bool) int[source]#

Evaluation of metric coefficients for given markers.

Parameters:

remove_outside (bool) – Whether to remove values that originate from markers outside of [0, 1]^d.

Returns:

counter – How many markers have been treated (not been skipped).

Return type:

int

Mapping kernels#

struphy.geometry.mappings_kernels.spline_3d(eta1: float, eta2: float, eta3: float, p: int[:], ind1: int[:, :], ind2: int[:, :], ind3: int[:, :], args: DomainArguments, f_out: float[:])[source]#

Point-wise evaluation of a 3d spline map \(F = (F_n)_{(n=x,y,z)}\) with

\[F_n = \sum_{ijk} c^n_{ijk} N_i(\eta_1) N_j(\eta_2) N_k(\eta_3)\,,\]

where \(c^n_{ijk}\) are the control points of component \(n\).

struphy.geometry.mappings_kernels.spline_3d_df(eta1: float, eta2: float, eta3: float, p: int[:], ind1: int[:, :], ind2: int[:, :], ind3: int[:, :], args: DomainArguments, df_out: float[:, :])[source]#

Jacobian matrix for struphy.geometry.mappings_kernels.spline_3d().

struphy.geometry.mappings_kernels.spline_2d_straight(eta1: float, eta2: float, eta3: float, p: int[:], ind1: int[:, :], ind2: int[:, :], args: DomainArguments, lz: float, f_out: float[:])[source]#

Point-wise evaluation of a 2d spline map \(F = (F_n)_{(n=x,y,z)}\) with

\[ \begin{align}\begin{aligned}F_{x(y)} &= \sum_{ij} c^{x(y)}_{ij} N_i(\eta_1) N_j(\eta_2) \,,\\F_z &= L_z*\eta_3\,.\end{aligned}\end{align} \]

where \(c^{x(y)}_{ij}\) are the control points in the \(\eta_1-\eta_2\)-plane, independent of \(\eta_3\).

struphy.geometry.mappings_kernels.spline_2d_straight_df(eta1: float, eta2: float, p: int[:], ind1: int[:, :], ind2: int[:, :], args: DomainArguments, lz: float, df_out: float[:, :])[source]#

Jacobian matrix for struphy.geometry.mappings_kernels.spline_2d_straight().

struphy.geometry.mappings_kernels.spline_2d_torus(eta1: float, eta2: float, eta3: float, p: int[:], ind1: int[:, :], ind2: int[:, :], args: DomainArguments, tor_period: float, f_out: float[:])[source]#

Point-wise evaluation of a 2d spline map \(F = (F_n)_{(n=x,y,z)}\) with

\[ \begin{align}\begin{aligned}S_{R(z)}(\eta_1, \eta_2) &= \sum_{ij} c^{R(z)}_{ij} N_i(\eta_1) N_j(\eta_2) \,,\\F_x &= S_R(\eta_1, \eta_2) * \cos(2\pi\eta_3)\\F_y &= - S_R(\eta_1, \eta_2) * \sin(2\pi\eta_3)\\F_z &= S_z(\eta_1, \eta_2)\,.\end{aligned}\end{align} \]

where \(c^{R(z)}_{ij}\) are the control points in the \(\eta_1-\eta_2\)-plane, independent of \(\eta_3\).

struphy.geometry.mappings_kernels.spline_2d_torus_df(eta1: float, eta2: float, eta3: float, p: int[:], ind1: int[:, :], ind2: int[:, :], args: DomainArguments, tor_period: float, df_out: float[:, :])[source]#

Jacobian matrix for struphy.geometry.mappings_kernels.spline_2d_torus().

struphy.geometry.mappings_kernels.cuboid(eta1: float, eta2: float, eta3: float, l1: float, r1: float, l2: float, r2: float, l3: float, r3: float, f_out: float[:])[source]#

Point-wise evaluation of

\[ \begin{align}\begin{aligned}F_x &= l_1 + (r_1 - l_1)\,\eta_1\,,\\F_y &= l_2 + (r_2 - l_2)\,\eta_2\,,\\F_z &= l_3 + (r_3 - l_3)\,\eta_3\,.\end{aligned}\end{align} \]

Note

Example with paramters \(l_1=0\,,r_1=1\,,l_2=0\,,r_2=1\,,l_3=0\) and \(r_3=1\):

sections/pics/mappings/cuboid.png
Parameters:
  • eta1 (float) – Logical coordinate in [0, 1].

  • eta2 (float) – Logical coordinate in [0, 1].

  • eta3 (float) – Logical coordinate in [0, 1].

  • l1 (float) – Left domain boundary.

  • l2 (float) – Left domain boundary.

  • l3 (float) – Left domain boundary.

  • r1 (float) – Right domain boundary.

  • r2 (float) – Right domain boundary.

  • r3 (float) – Right domain boundary.

  • f_out (array[float]) – Output: (x, y, z) = F(eta1, eta2, eta3).

struphy.geometry.mappings_kernels.cuboid_df(l1: float, r1: float, l2: float, r2: float, l3: float, r3: float, df_out: float[:, :])[source]#

Jacobian matrix for struphy.geometry.mappings_kernels.cuboid().

struphy.geometry.mappings_kernels.orthogonal(eta1: float, eta2: float, eta3: float, lx: float, ly: float, alpha: float, lz: float, f_out: float[:])[source]#

Point-wise evaluation of

\[ \begin{align}\begin{aligned}\begin{split}F_x &= L_x\,\left[\,\eta_1 + \\alpha\sin(2\pi\,\eta_1)\,\\right]\,,\end{split}\\\begin{split}F_y &= L_y\,\left[\,\eta_2 + \\alpha\sin(2\pi\,\eta_2)\,\\right]\,,\end{split}\\F_z &= L_z\,\eta_3\,.\end{aligned}\end{align} \]

Note

Example with paramters \(L_x=1\,,L_y=1\,,\\alpha=0.1\) and \(L_z=1\):

sections/pics/mappings/orthogonal.png
Parameters:
  • eta1 (float) – Logical coordinate in [0, 1].

  • eta2 (float) – Logical coordinate in [0, 1].

  • eta3 (float) – Logical coordinate in [0, 1].

  • lx (float) – Length in x-direction.

  • ly (float) – Length in yy-direction.

  • alpha (float) – Distortion factor.

  • lz (float) – Length in third direction.

  • f_out (array[float]) – Output: (x, y, z) = F(eta1, eta2, eta3).

struphy.geometry.mappings_kernels.orthogonal_df(eta1: float, eta2: float, lx: float, ly: float, alpha: float, lz: float, df_out: float[:, :])[source]#

Jacobian matrix for struphy.geometry.mappings_kernels.orthogonal().

struphy.geometry.mappings_kernels.colella(eta1: float, eta2: float, eta3: float, lx: float, ly: float, alpha: float, lz: float, f_out: float[:])[source]#

Point-wise evaluation of

\[ \begin{align}\begin{aligned}\begin{split}F_x &= L_x\,\left[\,\eta_1 + \\alpha\sin(2\pi\,\eta_1)\sin(2\pi\,\eta_2)\,\\right]\,,\end{split}\\\begin{split}F_y &= L_y\,\left[\,\eta_2 + \\alpha\sin(2\pi\,\eta_2)\sin(2\pi\,\eta_1)\,\\right]\,,\end{split}\\F_z &= L_z\,\eta_3\,.\end{aligned}\end{align} \]

Note

Example with paramters \(L_x=1\,,L_y=1\,,\\alpha=0.1\) and \(L_z=1\):

sections/pics/mappings/colella.png
Parameters:
  • eta1 (float) – Logical coordinate in [0, 1].

  • eta2 (float) – Logical coordinate in [0, 1].

  • eta3 (float) – Logical coordinate in [0, 1].

  • lx (float) – Length in x-direction.

  • ly (float) – Length in y-direction.

  • alpha (float) – Distortion factor.

  • lz (float) – Length in z-direction.

  • f_out (array[float]) – Output: (x, y, z) = F(eta1, eta2, eta3).

struphy.geometry.mappings_kernels.colella_df(eta1: float, eta2: float, lx: float, ly: float, alpha: float, lz: float, df_out: float[:, :])[source]#

Jacobian matrix for struphy.geometry.mappings_kernels.colella().

struphy.geometry.mappings_kernels.hollow_cyl(eta1: float, eta2: float, eta3: float, a1: float, a2: float, lz: float, poc: float, f_out: float[:])[source]#

Point-wise evaluation of

\[ \begin{align}\begin{aligned}\begin{split}F_x &= \left[\,a_1 + (a_2-a_1)\,\eta_1\,\\right]\cos(2\pi\,\eta_2 / poc)\,,\end{split}\\\begin{split}F_y &= \left[\,a_1 + (a_2-a_1)\,\eta_1\,\\right]\sin(2\pi\,\eta_2 / poc)\,,\end{split}\\F_z &= L_z\,\eta_3\,.\end{aligned}\end{align} \]

Note

Example with paramters \(a_1=0.2\,,a_2=1\) and \(L_z=3\):

sections/pics/mappings/hollow_cylinder.png
Parameters:
  • eta1 (float) – Logical coordinate in [0, 1].

  • eta2 (float) – Logical coordinate in [0, 1].

  • eta3 (float) – Logical coordinate in [0, 1].

  • a1 (float) – Inner radius.

  • a2 (float) – Outer radius.

  • lz (float) – Length in third direction.

  • poc (int) – periodicity in second direction.

  • f_out (array[float]) – Output: (x, y, z) = F(eta1, eta2, eta3).

struphy.geometry.mappings_kernels.hollow_cyl_df(eta1: float, eta2: float, a1: float, a2: float, lz: float, poc: float, df_out: float[:, :])[source]#

Jacobian matrix for struphy.geometry.mappings_kernels.hollow_cyl().

struphy.geometry.mappings_kernels.powered_ellipse(eta1: float, eta2: float, eta3: float, rx: float, ry: float, lz: float, s: float, f_out: float[:])[source]#

Point-wise evaluation of

\[ \begin{align}\begin{aligned}F_x &= r_x\,\eta_1^s\cos(2\pi\,\eta_2)\,,\\F_y &= r_y\,\eta_1^s\sin(2\pi\,\eta_2)\,,\\F_z &= L_z\,\eta_3\,.\end{aligned}\end{align} \]

Note

Example with paramters \(r_x=1\,,r_y=2,s=0.5\) and \(L_z=1\):

sections/pics/mappings/ellipse.png
Parameters:
  • eta1 (float) – Logical coordinate in [0, 1].

  • eta2 (float) – Logical coordinate in [0, 1].

  • eta3 (float) – Logical coordinate in [0, 1].

  • rx (float) – Axes lengths.

  • ry (float) – Axes lengths.

  • lz (float) – Length in third direction.

  • s (float) – Power of eta1

  • f_out (array[float]) – Output: (x, y, z) = F(eta1, eta2, eta3).

struphy.geometry.mappings_kernels.powered_ellipse_df(eta1: float, eta2: float, eta3: float, rx: float, ry: float, lz: float, s: float, df_out: float[:, :])[source]#

Jacobian matrix for struphy.geometry.mappings_kernels.powered_ellipse().

struphy.geometry.mappings_kernels.hollow_torus(eta1: float, eta2: float, eta3: float, a1: float, a2: float, r0: float, sfl: float, pol_period: float, tor_period: float, f_out: float[:])[source]#

Point-wise evaluation of

\[ \begin{align}\begin{aligned}\begin{split}F_x &= \lbrace\left[\,a_1 + (a_2-a_1)\,\eta_1\,\\right]\cos(\theta(\eta_1,\eta_2))+R_0\\rbrace\cos(2\pi\,\eta_3)\,,\end{split}\\\begin{split}F_y &= \lbrace\left[\,a_1 + (a_2-a_1)\,\eta_1\,\\right]\cos(\theta(\eta_1,\eta_2))+R_0\\rbrace\sin(2\pi\,\eta_3) \,,\end{split}\\\begin{split}F_z &= \,\,\,\left[\,a_1 + (a_2-a_1)\,\eta_1\,\\right]\sin(\theta(\eta_1,\eta_2)) \,,\end{split}\end{aligned}\end{align} \]

Note

Example with paramters \(a_1=0.2\,,a_2=1\) and \(R_0=3\):

sections/pics/mappings/hollow_torus.png
Parameters:
  • eta1 (float) – Logical coordinate in [0, 1].

  • eta2 (float) – Logical coordinate in [0, 1].

  • eta3 (float) – Logical coordinate in [0, 1].

  • a1 (float) – Inner radius.

  • a2 (float) – Outer radius.

  • r0 (float) – Major radius.

  • sfl (float) – Whether to use straight field line angular parametrization (yes: 1., no: 0.).

  • pol_period (float) – periodicity of theta used in the mapping: theta = 2*pi * eta2 / pol_period (if not sfl)

  • tor_period (int) – Toroidal periodicity built into the mapping: phi = 2*pi * eta3 / tor_period

  • f_out (array[float]) – Output: (x, y, z) = F(eta1, eta2, eta3).

struphy.geometry.mappings_kernels.hollow_torus_df(eta1: float, eta2: float, eta3: float, a1: float, a2: float, r0: float, sfl: float, pol_period: float, tor_period: float, df_out: float[:, :])[source]#

Jacobian matrix for struphy.geometry.mappings_kernels.hollow_torus().

struphy.geometry.mappings_kernels.shafranov_shift(eta1: float, eta2: float, eta3: float, rx: float, ry: float, lz: float, de: float, f_out: float[:])[source]#

Point-wise evaluation of

\[ \begin{align}\begin{aligned}F_x &= r_x\,\eta_1\cos(2\pi\,\eta_2)+(1-\eta_1^2)r_x\Delta\,,\\F_y &= r_y\,\eta_1\sin(2\pi\,\eta_2)\,,\\F_z &= L_z\,\eta_3\,.\end{aligned}\end{align} \]

Note

Example with paramters \(r_x=1\,,r_y=1\,,L_z=1\) and \(\Delta=0.2\):

sections/pics/mappings/shafranov_shift.png
Parameters:
  • eta1 (float) – Logical coordinate in [0, 1].

  • eta2 (float) – Logical coordinate in [0, 1].

  • eta3 (float) – Logical coordinate in [0, 1].

  • rx (float) – Axes lengths.

  • ry (float) – Axes lengths.

  • lz (float) – Length in third direction.

  • de (float) – Shift factor, should be in [0, 0.1].

  • f_out (array[float]) – Output: (x, y, z) = F(eta1, eta2, eta3).

struphy.geometry.mappings_kernels.shafranov_shift_df(eta1: float, eta2: float, eta3: float, rx: float, ry: float, lz: float, de: float, df_out: float[:, :])[source]#

Jacobian matrix for struphy.geometry.mappings_kernels.shafranov_shift().

struphy.geometry.mappings_kernels.shafranov_sqrt(eta1: float, eta2: float, eta3: float, rx: float, ry: float, lz: float, de: float, f_out: float[:])[source]#

Point-wise evaluation of

\[ \begin{align}\begin{aligned}F_x &= r_x\,\eta_1\cos(2\pi\,\eta_2)+(1-\sqrt \eta_1)r_x\Delta\,,\\F_y &= r_y\,\eta_1\sin(2\pi\,\eta_2)\,,\\F_z &= L_z\,\eta_3\,.\end{aligned}\end{align} \]

Note

No example plot yet.

Parameters:
  • eta1 (float) – Logical coordinate in [0, 1].

  • eta2 (float) – Logical coordinate in [0, 1].

  • eta3 (float) – Logical coordinate in [0, 1].

  • rx (float) – Axes lengths.

  • ry (float) – Axes lengths.

  • lz (float) – Length in third direction.

  • de (float) – Shift factor, should be in [0, 0.1].

  • f_out (array[float]) – Output: (x, y, z) = F(eta1, eta2, eta3).

struphy.geometry.mappings_kernels.shafranov_sqrt_df(eta1: float, eta2: float, eta3: float, rx: float, ry: float, lz: float, de: float, df_out: float[:, :])[source]#

Jacobian matrix for struphy.geometry.mappings_kernels.shafranov_sqrt().

struphy.geometry.mappings_kernels.shafranov_dshaped(eta1: float, eta2: float, eta3: float, r0: float, lz: float, dx: float, dy: float, dg: float, eg: float, kg: float, f_out: float[:])[source]#

Point-wise evaluation of

\[ \begin{align}\begin{aligned}\begin{split}x &= R_0\left[1 + (1 - \eta_1^2)\Delta_x + \eta_1\epsilon\cos(2\pi\,\eta_2 + \\arcsin(\delta)\eta_1\sin(2\pi\,\eta_2)) \\right]\,,\end{split}\\\begin{split}y &= R_0\left[ (1 - \eta_1^2)\Delta_y + \eta_1\epsilon\kappa\sin(2\pi\,\eta_2)\\right]\,,\end{split}\\z &= L_z\,\eta_3\,.\end{aligned}\end{align} \]

Note

Example with paramters \(R_0=3\,,L_z=1\,,\Delta_x=0.1\,,\Delta_y=0\,,\delta=0.2\,,\epsilon=1/3\) and \(\kappa=1.5\):

sections/pics/mappings/shafranov_dshaped.png
Parameters:
  • eta1 (float) – Logical coordinate in [0, 1].

  • eta2 (float) – Logical coordinate in [0, 1].

  • eta3 (float) – Logical coordinate in [0, 1].

  • r0 (float) – Base radius.

  • lz (float) – Length in third direction.

  • dx (float) – Shafranov shift in x-direction.

  • dy (float) – Shafranov shift in y-direction.

  • dg (float) – Delta = sin(alpha): Triangularity, shift of high point.

  • eg (float) – Epsilon: Inverse aspect ratio a/r0.

  • kg (float) – Kappa: Ellipticity (elongation).

  • f_out (array[float]) – Output: (x, y, z) = F(eta1, eta2, eta3).

struphy.geometry.mappings_kernels.shafranov_dshaped_df(eta1: float, eta2: float, eta3: float, r0: float, lz: float, dx: float, dy: float, dg: float, eg: float, kg: float, df_out: float[:, :])[source]#

Jacobian matrix for struphy.geometry.mappings_kernels.shafranov_dshaped().

Transform kernels#

  1. Basic pull-back (physical –> logical) operations between scalar fields, vector fields and differential p-forms:

  • 0-form : a^0 = a

  • 1-form : (a^1_1, a^1_2, a^1_3) = DF^T (ax, ay, az)

  • 2-form : (a^2_1, a^2_2, a^2_3) = |det(DF)| DF^(-1) (ax, ay, az)

  • 3-form : a^3 = |det(DF)| a

  • vector : (a_1 , a_2 , a_3 ) = DF^(-1) (ax, ay, az)

  1. Basic push-forward (logical –> physical) operations between scalar fields, vector fields and differential p-forms:

  • 0-form : a = a^0

  • 1-form : (ax, ay, az) = DF^(-T) (a^1_1, a^1_2, a^1_3)

  • 2-form : (ax, ay, az) = 1/|det(DF)| DF (a^2_1, a^2_2, a^2_3)

  • 3-form : a = 1/|det(DF)| a^3

  • vector : (ax, ay, az) = DF (a_1 , a_2 , a_3 )

  1. Basic transformations between scalar fields, vector fields and differential p-forms:

  • 0-form –> 3-form : a^3 = a^0 * |det(DF)|

  • 3-form –> 0-form : a^0 = a^3 / |det(DF)|

  • 1-form –> 2-form : (a^2_1, a^2_2, a^2_3) = G^(-1) * (a^1_1, a^1_2, a^1_3) * |det(DF)|

  • 2-form –> 1-form : (a^1_1, a^1_2, a^1_3) = G * (a^2_1, a^2_2, a^2_3) / |det(DF)|

  • norm vector –> vector(a_1, a_2, a_3) = (a^*_1/sqrt(sum(DF[:,0]^2)),

    a^*_2/sqrt(sum(DF[:,1]^2)), a^*_3/sqrt(sum(DF[:,2]^2)))

  • norm vector –> 1-form(a^1_1, a^1_2, a^1_3) = G * (a^*_1/sqrt(sum(DF[:,0]^2)),

    a^*_2/sqrt(sum(DF[:,1]^2)), a^*_3/sqrt(sum(DF[:,2]^2)))

  • norm vector –> 2-form(a^2_1, a^2_2, a^2_3) = (a^*_1/sqrt(sum(DF[:,0]^2)),

    a^*_2/sqrt(sum(DF[:,1]^2)), a^*_3/sqrt(sum(DF[:,2]^2))) * |det(DF)|

  • vector –> 1-form : (a^1_1, a^1_2, a^1_3) = G * (a_1, a_2, a_3)

  • vector –> 2-form : (a^2_1, a^2_2, a^2_3) = (a_1, a_2, a_3) * |det(DF)|

  • 1-form –> vector : (a_1, a_2, a_3) = G^(-1) * (a^1_1, a^1_2, a^1_3)

  • 2-form –> vector : (a_1, a_2, a_3) = (a^2_1, a^2_2, a^2_3) / |det(DF)|

struphy.geometry.transform_kernels.pull(a: float[:], eta1: float, eta2: float, eta3: float, kind_fun: int, args_domain: DomainArguments, out: float[:])[source]#

Pull-back of a Cartesian scalar/vector field to a differential p-form.

Parameters:
  • a (float[:]) – Value of scalar field a[0] or values of Cartesian components of vector field (a[0], a[1], a[2]).

  • eta1 (float) – Logical evaluation points.

  • eta2 (float) – Logical evaluation points.

  • eta3 (float) – Logical evaluation points.

  • kind_fun (int) – Which pull-back to be performed.

  • args_domain (DomainArguments) – Domain info.

  • out (float[:]) – Output values.

struphy.geometry.transform_kernels.push(a: float[:], eta1: float, eta2: float, eta3: float, kind_fun: int, args_domain: DomainArguments, out: float[:])[source]#

Pushforward of a differential p-forms to a Cartesian scalar/vector field.

Parameters:
  • a (float[:]) – Value of scalar p-form a[0] or values of components of vector valued p-form (a[0], a[1], a[2]).

  • eta1 (float) – Logical evaluation points.

  • eta2 (float) – Logical evaluation points.

  • eta3 (float) – Logical evaluation points.

  • kind_fun (int) – Which pushforward to be performed.

  • args_domain (DomainArguments) – Domain info.

  • out (float[:]) – Output values.

struphy.geometry.transform_kernels.tran(a: float[:], eta1: float, eta2: float, eta3: float, kind_fun: int, args_domain: DomainArguments, out: float[:])[source]#

Transformations between differential p-forms and/or vector fields.

Parameters:
  • a (float[:]) – Value of scalar function a[0] or values of components of vector valued functions (a[0], a[1], a[2]).

  • eta1 (float) – Logical evaluation points.

  • eta2 (float) – Logical evaluation points.

  • eta3 (float) – Logical evaluation points.

  • kind_fun (int) – Which transformation to be performed.

  • args_domain (DomainArguments) – Domain info.

  • out (float[:]) – Output values.

struphy.geometry.transform_kernels.kernel_pullpush(a: float[:, :, :, :], eta1: float[:, :, :], eta2: float[:, :, :], eta3: float[:, :, :], kind_transform: int, kind_fun: int, args_domain: DomainArguments, is_sparse_meshgrid: bool, out: float[:, :, :, :])[source]#

Pull-backs, pushforwards and transformations on a given 3d grid of evaluation points.

Parameters:
  • a (float[:,:,:,:]) – 3d values of scalar function a[0, i, j, k] or 3d values of components of vector valued function a[:, i, j, k].

  • eta1 (float[:,:,:]) – 3d evaluation point sets.

  • eta2 (float[:,:,:]) – 3d evaluation point sets.

  • eta3 (float[:,:,:]) – 3d evaluation point sets.

  • kind_transform (int) – Which general transformation to be performed (pull, push or tran).

  • kind_fun (int) – Which detailed transformation to be performed.

  • args_domain (DomainArguments) – Domain info.

  • is_sparse_meshgrid (bool) – Whether the evaluation points were obtained from a sparse meshgrid.

  • out (float[:,:,:,:]) – Output values.

struphy.geometry.transform_kernels.kernel_pullpush_pic(a: float[:, :], markers: float[:, :], kind_transform: int, kind_fun: int, args_domain: DomainArguments, out: float[:, :], remove_outside: bool) int[source]#

Pull-backs, pushforwards and transformations for given markers.

Parameters:
  • a (float[:,:]) – Values of scalar function a[0, ip] or values of components of a vector valued function (a[0, ip], a[1, ip], a[2, ip]).

  • markers (float[:,:]) – Evaluation points in marker format (eta1 = markers[:, 0], eta2 = markers[:, 1], eta3 = markers[:, 2]).

  • kind_transform (int) – Which general transformation to be performed (pull, push or tran).

  • kind_fun (int) – Which detailed transformation to be performed.

  • args_domain (DomainArguments) – Domain info.

  • out (float[:,:]) – Output values.

  • remove_outside (bool) – Whether to remove values that originate from markers outside of [0, 1]^d.