Gyselalib++
AdvectionPseudoCartesianDomain< Mapping > Class Template Reference

Define the pseudo-Cartesian domain for the advection. More...

Inheritance diagram for AdvectionPseudoCartesianDomain< Mapping >:
AdvectionDomain< Mapping >

Public Types

using Matrix2x2 = std::array< std::array< double, 2 >, 2 >
 Define a 2x2 matrix with an 2D array of an 2D array.
 
using X_adv = DimX_pC
 The first dimension in the advection domain.
 
using Y_adv = DimY_pC
 The second dimension in the advection domain.
 
using CoordXY_adv = Coord< X_adv, Y_adv >
 The coordinate type associated to the dimensions in the advection domain.
 

Public Member Functions

 AdvectionPseudoCartesianDomain (Mapping const &mapping, double epsilon=1e-12)
 Instantiate an AdvectionPseudoCartesianDomain advection domain. More...
 
void advect_feet (host_t< FieldRTheta< CoordRTheta >> feet_coords_rp, host_t< DConstVectorFieldRTheta< X_adv, Y_adv >> const &advection_field, double const dt) const
 Advect the characteristic feet. More...
 

Detailed Description

template<class Mapping>
class AdvectionPseudoCartesianDomain< Mapping >

Define the pseudo-Cartesian domain for the advection.

The pseudo-Cartesian domain is recommended for not analytically invertible mappings.

We introduce a circular mapping \( \mathcal{G}\) from the logical domain to the pseudo-Cartesian domain. The circular mapping is analytically invertible.

The advection field is defined in the physical domain, so we need to define it in the pseudo-Cartesian domain before advecting (AdvectionPseudoCartesianDomain::compute_advection_field). To do so, we use the Jacobian matrix \( J_{\mathcal{F}}J_{\mathcal{G}}^{-1} \). This matrix is invertible :

\( A_{\text{cart}} = (J_{\mathcal{F}}J_{\mathcal{G}}^{-1})^{-1} A \).

The main difficulty is the center point. So for \( r \in [0, \varepsilon] \), we linearize the Jacobian matrix:

\( (J_{\mathcal{F}}J_{\mathcal{G}}^{-1})^{-1} (r, \theta) = \left(1 - \frac{r}{\varepsilon} \right) (J_{\mathcal{F}}J_{\mathcal{G}}^{-1})^{-1} (0, \theta) + \frac{r}{\varepsilon} (J_{\mathcal{F}}J_{\mathcal{G}}^{-1})^{-1} (\varepsilon, \theta) \).

The \( (J_{\mathcal{F}}J_{\mathcal{G}}^{-1})^{-1} (0, \theta) \) matrices are implemented in the Curvilinear2DToCartesian::to_pseudo_cartesian_jacobian_center_matrix child classes.

The method to advect is here given by:

  • compute \( (x_{\text{cart}},y_{\text{cart}})_{ij} = \mathcal{G} (r,\theta)_{ij} \),
  • advect
    • \( \text{feet}_{x_{\text{cart}}, ij} = x_{\text{cart}, ij} - dt A_{\text{cart}, x}(x_{\text{cart}, ij}, y_{\text{cart}, ij}) \),
    • \( \text{feet}_{y_{\text{cart}}, ij} = y_{\text{cart}, ij} - dt A_{\text{cart}, y}(x_{\text{cart}, ij}, y_{\text{cart}, ij}) \),
  • compute \( (\text{feet}_{r, ij}, \text{feet}_{\theta, ij}) = \mathcal{G}^{-1} (\text{feet}_{x_{\text{cart}}, ij}, \text{feet}_{y_{\text{cart}}, ij}) \).

More details can be found in Edoardo Zoni's article (https://doi.org/10.1016/j.jcp.2019.108889).

See also
Curvilinear2DToCartesian
DiscreteToCartesian

Constructor & Destructor Documentation

◆ AdvectionPseudoCartesianDomain()

template<class Mapping >
AdvectionPseudoCartesianDomain< Mapping >::AdvectionPseudoCartesianDomain ( Mapping const &  mapping,
double  epsilon = 1e-12 
)
inline

Instantiate an AdvectionPseudoCartesianDomain advection domain.

Parameters
[in]mappingThe mapping from the logical domain to the physical domain.
[in]epsilon\( \varepsilon \) parameter used for the linearization of the advection field around the central point.

Member Function Documentation

◆ advect_feet()

template<class Mapping >
void AdvectionPseudoCartesianDomain< Mapping >::advect_feet ( host_t< FieldRTheta< CoordRTheta >>  feet_coords_rp,
host_t< DConstVectorFieldRTheta< X_adv, Y_adv >> const &  advection_field,
double const  dt 
) const
inline

Advect the characteristic feet.

In the Backward Semi-Lagrangian method, the advection of a function uses the conservation along the characteristic property. So, we firstly compute the characteristic feet and then interpolate the function at this characteristic feet.

The function implemented here deals with the computation of the characteristic feet. The IFootFinder class uses a time integration method to solve the characteristic equation. The BslAdvectionRTheta class calls advect_feet to compute the characteristic feet and interpolate the function we want to advect.

The advect_feet implemented here computes only

  • \( (\text{feet}_r, \text{feet}_\theta) = \mathcal{G}^{-1} (\text{feet}_x, \text{feet}_y) \)

with

  • \( \text{feet}_x = x_{ij} - dt A_x(x_{ij}, y_{ij}) \),
  • \( \text{feet}_y = y_{ij} - dt A_y(x_{ij}, y_{ij}) \),

and

  • \( (x,y)_{ij} = \mathcal{G}(r_{ij}, \theta_{ij})\), with \(\{(r, \theta)_{ij}\}_{ij} \) the logical mesh points,
  • \( A \) the advection field in the advection domain,
  • \( \mathcal{G} \) the mapping from the logical domain to the advection domain.
Parameters
[in,out]feet_coords_rpThe computed characteristic feet in the logical domain. On input: the points we want to advect. On output: the characteristic feet.
[in]advection_fieldThe advection field defined in the advection domain.
[in]dtThe time step.

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