DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
DNDS::Euler::ArrayGRADV< nVarsFixed, gDim > Class Template Reference

MPI-parallel distributed array of per-cell gradient matrices. More...

#include <Euler.hpp>

Inheritance diagram for DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >:
[legend]
Collaboration diagram for DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >:
[legend]

Public Types

using t_self = ArrayGRADV< nVarsFixed, gDim >
 Self type alias.
 
using t_base = CFV::tUGrad< nVarsFixed, gDim >
 Base distributed gradient array type.
 
- Public Types inherited from DNDS::ArrayDof< n_m, n_n >
using t_base = ArrayEigenMatrixPair< n_m, n_n >
 
template<DeviceBackend B>
using t_deviceView = ArrayDofDeviceView< B, n_m, n_n >
 Mutable device view alias.
 
template<DeviceBackend B>
using t_deviceViewConst = ArrayDofDeviceViewConst< B, n_m, n_n >
 Const device view alias.
 
using t_self = ArrayDof< n_m, n_n >
 
template<DeviceBackend B>
using t_ops = ArrayDofOp< B, n_m, n_n >
 Static dispatcher alias selecting host / CUDA implementation.
 
using t_element_mat = Eigen::Matrix< real, RowSize_To_EigenSize(n_m), RowSize_To_EigenSize(n_n)>
 Shape of one DOF row as an Eigen matrix.
 
- Public Types inherited from DNDS::ArrayPair< TArray >
using t_self = ArrayPair< TArray >
 
using t_arr = TArray
 
using TTrans = typename ArrayTransformerType< TArray >::Type
 
template<DeviceBackend B>
using t_deviceView = ArrayPairDeviceView< B, TArray >
 Device-view template alias: t_deviceView<DeviceBackend::CUDA> gives the mutable CUDA view type for this pair.
 
template<DeviceBackend B>
using t_deviceViewConst = ArrayPairDeviceViewConst< B, TArray >
 Const-device-view template alias.
 

Public Member Functions

void setConstant (real R)
 Set all gradient entries to a uniform scalar value.
 
template<class TR >
void setConstant (const TR &R)
 Set all cells' gradient matrices to copies of a given matrix.
 
void operator+= (t_self &R)
 Element-wise addition: this[i] += R[i] for every cell.
 
void operator-= (t_self &R)
 Element-wise subtraction: this[i] -= R[i] for every cell.
 
void operator*= (real R)
 Uniform scalar multiplication: this[i] *= R for every cell.
 
void operator*= (std::vector< real > &R)
 Per-cell scalar multiplication from a vector of cell-wise weights.
 
void operator*= (ArrayDOFV< 1 > &R)
 Component-wise multiplication by a scalar-per-cell ArrayDOFV<1>.
 
void operator*= (const Eigen::Array< real, 1, nVarsFixed > &R)
 Column-wise scaling by a per-variable multiplier row-vector.
 
void operator= (t_self &R)
 Deep copy assignment from another ArrayGRADV.
 
- Public Member Functions inherited from DNDS::ArrayDof< n_m, n_n >
template<DeviceBackend B>
t_deviceView< B > deviceView ()
 Build a mutable device view (wraps the base-class implementation).
 
template<DeviceBackend B>
t_deviceViewConst< B > deviceView () const
 Build a const device view.
 
void clone (const t_self &R)
 Deep copy from another ArrayDof. Delegates to the base-class clone.
 
void setConstant (real R)
 Set every entry of every (father+son) row to the scalar R.
 
void setConstant (const Eigen::Ref< const t_element_mat > &R)
 Set every row to the matrix R (must have shape n_m x n_n).
 
void operator+= (const t_self &R)
 In-place element-wise add: this += R.
 
void operator+= (real R)
 Add the scalar R to every entry.
 
void operator+= (const Eigen::Ref< const t_element_mat > &R)
 Add a per-row matrix R (same to every row).
 
void operator-= (const t_self &R)
 In-place element-wise subtract: this -= R.
 
void operator*= (real R)
 Scalar multiply in place.
 
template<int n_m_T = n_m>
std::enable_if_t<!(n_m_T==1 &&n_n==1)> operator*= (const ArrayDof< 1, 1 > &R)
 Scale each row by a corresponding scalar stored in R (a 1x1 ArrayDof).
 
void operator*= (const Eigen::Ref< const t_element_mat > &R)
 In-place multiplication by a small fixed matrix (same applied to every row).
 
void operator*= (const t_self &R)
 Element-wise multiply: this *= R (Hadamard).
 
void operator/= (const t_self &R)
 Element-wise divide: this /= R.
 
void operator= (const t_self &R)
 Value-copy assignment from another ArrayDof of identical layout.
 
void addTo (const t_self &R, real r)
 AXPY: this += r * R. One of the hot-path solver primitives.
 
real norm2 ()
 Global L2 norm (MPI-collective). sqrt(sum_i sum_j x_ij^2).
 
real norm2 (const t_self &R)
 Global L2 distance between this and R (collective).
 
real min ()
 Global minimum across all entries (collective).
 
real max ()
 Global maximum across all entries (collective).
 
real sum ()
 Global sum of all entries (collective).
 
t_element_mat componentWiseNorm1 ()
 Per-component global L1 norm, returned as an n_m x n_n matrix (collective).
 
t_element_mat componentWiseNorm1 (const t_self &R)
 Per-component global L1 distance between this and R (collective).
 
real dot (const t_self &R)
 Global inner product: sum_i sum_j x_ij * R_ij (collective).
 
void to_device (DeviceBackend backend)
 Mirror both father and son to the given device backend.
 
void to_host ()
 Bring both father and son mirrors back to host memory.
 
- Public Member Functions inherited from DNDS::ArrayPair< TArray >
void clone (const t_self &R)
 Deep-copy: allocate new father / son and copy their data; rebind trans.
 
decltype(father->operator[](index(0))) operator[] (index i) const
 Read-only row-pointer access in the combined address space.
 
decltype(father->operator[](index(0))) operator[] (index i)
 Mutable row-pointer access in the combined address space.
 
template<class... TOthers>
decltype(autooperator() (index i, TOthers... aOthers)
 N-ary element access in the combined space (mutable). Arguments after the row index are forwarded to the underlying operator().
 
template<class... TOthers>
decltype(autooperator() (index i, TOthers... aOthers) const
 N-ary element access (const).
 
template<class TF >
auto runFunctionAppendedIndex (index i, TF &&F)
 Invoke F(array, localIndex) on either father or son depending on which range i falls into.
 
auto RowSize () const
 Uniform row width (delegates to father).
 
auto RowSize (index i) const
 Per-row width in the combined address space.
 
void ResizeRow (index i, rowsize rs)
 Resize a single row in the combined address space.
 
template<class... TOthers>
void ResizeRow (index i, TOthers... aOthers)
 Variadic ResizeRow overload that forwards extra args.
 
index Size () const
 Combined row count (father->Size() + son->Size()).
 
void TransAttach ()
 Bind the transformer to the current father / son pointers.
 
template<typename... Args>
void InitPair (const std::string &name, Args &&...args)
 Allocate both father and son arrays, forwarding all args to TArray constructor.
 
template<class TPrimaryPair >
void BorrowAndPull (TPrimaryPair &primary)
 Attach, borrow ghost indexing from a primary pair, create MPI types, and pull once.
 
template<class TPrimaryPair >
void BorrowSetup (TPrimaryPair &primary)
 Attach, borrow ghost indexing from a primary pair, and create MPI types (no pull).
 
void CompressBoth ()
 Compress both father and son CSR arrays (no-op for non-CSR layouts).
 
void CopyFather (t_self &R)
 Copy only the father's data from another pair (shallow).
 
void SwapDataFatherSon (t_self &R)
 Swap both father and son data with another pair of the same type.
 
std::size_t hash ()
 Combined hash across ranks. Used for determinism / equality checks in tests.
 
void WriteSerialize (Serializer::SerializerBaseSSP serializerP, const std::string &name, bool includePIG=true, bool includeSon=true)
 Writes the ArrayPair (father, optional son, optional ghost mapping).
 
void WriteSerialize (Serializer::SerializerBaseSSP serializerP, const std::string &name, const std::vector< index > &origIndex, bool includePIG=true, bool includeSon=true)
 Writes the ArrayPair with an origIndex companion dataset for redistribution support.
 
void ReadSerialize (Serializer::SerializerBaseSSP serializerP, const std::string &name, bool includePIG=true, bool includeSon=true)
 Reads an ArrayPair written by WriteSerialize (same partition count).
 
void ReadSerializeRedistributed (Serializer::SerializerBaseSSP serializerP, const std::string &name, const std::vector< index > &newOrigIndex)
 Reads ArrayPair data from HDF5 with redistribution support.
 
template<DeviceBackend B>
auto deviceView ()
 Produce a mutable device view; both father and son must be allocated.
 
template<DeviceBackend B>
auto deviceView () const
 Produce a const device view.
 
void to_device (DeviceBackend backend)
 Mirror both father and son to the given device backend.
 
void to_host ()
 Bring both father and son mirrors back to host memory.
 

Additional Inherited Members

- Static Public Member Functions inherited from DNDS::ArrayPair< TArray >
static constexpr bool IsCSR ()
 Whether the underlying array uses CSR storage.
 
- Public Attributes inherited from DNDS::ArrayPair< TArray >
ssp< TArray > father
 Owned-side array (must be resized before ghost setup).
 
ssp< TArray > son
 Ghost-side array (sized automatically by createMPITypes / BorrowAndPull).
 
TTrans trans
 Ghost-communication engine bound to father and son.
 

Detailed Description

template<int nVarsFixed, int gDim>
class DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >

MPI-parallel distributed array of per-cell gradient matrices.

Wraps CFV::tUGrad<nVarsFixed, gDim> (which is ArrayDof<gDim, nVarsFixed>), adding element-wise arithmetic operators for gradient manipulation.

Each cell stores a (gDim x nVars) matrix where row d contains the partial derivative of all conservative variables with respect to spatial coordinate d:

\[ \mathbf{G}_i = \begin{bmatrix} \partial u_1/\partial x_1 & \cdots & \partial u_{nVars}/\partial x_1 \\ \vdots & \ddots & \vdots \\ \partial u_1/\partial x_{gDim} & \cdots & \partial u_{nVars}/\partial x_{gDim} \end{bmatrix} \]

Template Parameters
nVarsFixedCompile-time number of conservative variables (matrix columns), or Eigen::Dynamic for runtime-sized.
gDimGeometry (mesh) spatial dimension (2 or 3), determining the number of gradient rows.
See also
ArrayDOFV for per-cell DOF vectors.
ArrayRECV for reconstruction coefficient storage.

Definition at line 642 of file Euler.hpp.

Member Typedef Documentation

◆ t_base

template<int nVarsFixed, int gDim>
using DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >::t_base = CFV::tUGrad<nVarsFixed, gDim>

Base distributed gradient array type.

Definition at line 646 of file Euler.hpp.

◆ t_self

template<int nVarsFixed, int gDim>
using DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >::t_self = ArrayGRADV<nVarsFixed, gDim>

Self type alias.

Definition at line 645 of file Euler.hpp.

Member Function Documentation

◆ operator*=() [1/4]

template<int nVarsFixed, int gDim>
void DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >::operator*= ( ArrayDOFV< 1 > &  R)
inline

Component-wise multiplication by a scalar-per-cell ArrayDOFV<1>.

Each cell's gradient matrix is scaled by the single scalar stored in the corresponding cell of R.

Parameters
RScalar-per-cell array (ArrayDOFV<1>).

Definition at line 729 of file Euler.hpp.

Here is the call graph for this function:

◆ operator*=() [2/4]

template<int nVarsFixed, int gDim>
void DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >::operator*= ( const Eigen::Array< real, 1, nVarsFixed > &  R)
inline

Column-wise scaling by a per-variable multiplier row-vector.

Each column (variable) of every cell's gradient matrix is multiplied by the corresponding element of R: this[i].col(j) *= R(j).

Parameters
RRow array of per-variable multipliers (1 x nVarsFixed).
Note
Only operates on host memory. Asserts if device storage is active.
OpenMP parallelized when DNDS_DIST_MT_USE_OMP is defined.

Definition at line 744 of file Euler.hpp.

Here is the call graph for this function:

◆ operator*=() [3/4]

template<int nVarsFixed, int gDim>
void DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >::operator*= ( real  R)
inline

Uniform scalar multiplication: this[i] *= R for every cell.

Parameters
RScalar multiplier applied to all entries of all cells.

Definition at line 697 of file Euler.hpp.

Here is the call graph for this function:

◆ operator*=() [4/4]

template<int nVarsFixed, int gDim>
void DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >::operator*= ( std::vector< real > &  R)
inline

Per-cell scalar multiplication from a vector of cell-wise weights.

Multiplies each cell's gradient matrix by the corresponding scalar in R: this[i] *= R[i]. Only operates on father (locally-owned) cells.

Parameters
RVector of per-cell scalar multipliers. Must have size >= number of locally-owned cells.
Note
OpenMP parallelized when DNDS_DIST_MT_USE_OMP is defined.

Definition at line 712 of file Euler.hpp.

Here is the call graph for this function:

◆ operator+=()

template<int nVarsFixed, int gDim>
void DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >::operator+= ( t_self R)
inline

Element-wise addition: this[i] += R[i] for every cell.

Parameters
RSource gradient array with the same distribution layout.

Definition at line 681 of file Euler.hpp.

Here is the call graph for this function:

◆ operator-=()

template<int nVarsFixed, int gDim>
void DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >::operator-= ( t_self R)
inline

Element-wise subtraction: this[i] -= R[i] for every cell.

Parameters
RSource gradient array with the same distribution layout.

Definition at line 689 of file Euler.hpp.

Here is the call graph for this function:

◆ operator=()

template<int nVarsFixed, int gDim>
void DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >::operator= ( t_self R)
inline

Deep copy assignment from another ArrayGRADV.

Parameters
RSource array to copy from. Must have the same distribution layout.

Definition at line 758 of file Euler.hpp.

Here is the call graph for this function:

◆ setConstant() [1/2]

template<int nVarsFixed, int gDim>
template<class TR >
void DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >::setConstant ( const TR R)
inline

Set all cells' gradient matrices to copies of a given matrix.

Template Parameters
TRMatrix type compatible with per-cell gradient storage.
Parameters
RReference matrix; each cell's gradient is set to this value.
Note
Only operates on host memory. Asserts if device storage is active.
OpenMP parallelized when DNDS_DIST_MT_USE_OMP is defined.

Definition at line 666 of file Euler.hpp.

Here is the call graph for this function:

◆ setConstant() [2/2]

template<int nVarsFixed, int gDim>
void DNDS::Euler::ArrayGRADV< nVarsFixed, gDim >::setConstant ( real  R)
inline

Set all gradient entries to a uniform scalar value.

Parameters
RScalar value to fill every matrix entry with.

Definition at line 652 of file Euler.hpp.

Here is the call graph for this function:

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