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

MPI-parallel distributed array of per-cell Degrees-of-Freedom (conservative variable vectors). More...

#include <Euler.hpp>

Inheritance diagram for DNDS::Euler::ArrayDOFV< nVarsFixed >:
[legend]
Collaboration diagram for DNDS::Euler::ArrayDOFV< nVarsFixed >:
[legend]

Public Types

using t_self = ArrayDOFV< nVarsFixed >
 Self type alias for CRTP-style forwarding.
 
using t_base = CFV::tUDof< nVarsFixed >
 Base distributed DOF array type.
 
using t_element_mat = typename t_base::t_element_mat
 Per-element Eigen matrix/vector 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 DOF entries to a uniform scalar value.
 
void setConstant (const Eigen::Ref< const t_element_mat > &R)
 Set all DOF entries to copies of a given vector.
 
void operator+= (const t_self &R)
 Element-wise addition: this[i] += R[i] for every cell.
 
void operator-= (const 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= (const t_self &R)
 Deep copy assignment from another ArrayDOFV.
 
void addTo (const t_self &R, real r)
 Scaled addition: this[i] += r * R[i] for every cell (axpy-style).
 
void operator*= (const std::vector< real > &R)
 Per-cell scalar multiplication from a vector of cell-wise weights.
 
template<int nVarsFixed_T = nVarsFixed>
std::enable_if_t<!(nVarsFixed_T==1)> operator*= (const ArrayDOFV< 1 > &R)
 Component-wise multiplication by a scalar-per-cell ArrayDOFV<1>.
 
void operator+= (const Eigen::Vector< real, nVarsFixed > &R)
 Add a constant vector to every cell: this[i] += R for all cells.
 
void operator+= (real R)
 Add a scalar to every component of every cell.
 
void operator*= (const Eigen::Vector< real, nVarsFixed > &R)
 Component-wise multiplication by a constant vector.
 
void operator*= (const t_self &R)
 Element-wise (Hadamard) multiplication: this[i] = this[i] .* R[i].
 
void operator/= (const t_self &R)
 Element-wise division: this[i] = this[i] ./ R[i].
 
real norm2 ()
 Compute the global L2 norm across all MPI ranks.
 
Eigen::Vector< real, nVarsFixed > componentWiseNorm1 ()
 Compute the per-component global L1 norm.
 
Eigen::Vector< real, nVarsFixed > min ()
 Compute the global component-wise minimum across all MPI ranks.
 
real dot (const t_self &R)
 Compute the global dot product <this, R> with MPI reduction.
 
template<class TMultL , class TMultR >
real dot (const t_self &R, TMultL &&mL, TMultR &&mR)
 Compute a weighted global dot product with per-component multipliers.
 
template<class TR >
void setMaxWith (TR R)
 Clamp all components from below: this[i] = max(this[i], R).
 
- 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>
class DNDS::Euler::ArrayDOFV< nVarsFixed >

MPI-parallel distributed array of per-cell Degrees-of-Freedom (conservative variable vectors).

Wraps CFV::tUDof<nVarsFixed> (which is ArrayDof<nVarsFixed, 1>), adding element-wise arithmetic operators, norms, dot products, and MPI-collective reductions that operate across all ranks transparently.

Each cell stores a single column vector of nVarsFixed conservative variables (e.g. \([\rho,\; \rho u,\; \rho v,\; \rho w,\; E]\) for a 3D NS model). Instances of this class are used for the solution vector u, the solution increment uInc, the right-hand-side rhs, and other per-cell vector quantities throughout the EulerEvaluator.

The "father" sub-array holds locally-owned cells; the "son" sub-array holds ghost cells received from neighboring MPI ranks. Reduction operations (norm, dot, min) only iterate over father cells and then call MPI_Allreduce for global results.

Template Parameters
nVarsFixedCompile-time number of conservative variables per cell, or Eigen::Dynamic for runtime-sized vectors.
See also
ArrayRECV for reconstruction coefficient storage.
ArrayGRADV for gradient storage.

Definition at line 92 of file Euler.hpp.

Member Typedef Documentation

◆ t_base

template<int nVarsFixed>
using DNDS::Euler::ArrayDOFV< nVarsFixed >::t_base = CFV::tUDof<nVarsFixed>

Base distributed DOF array type.

Definition at line 96 of file Euler.hpp.

◆ t_element_mat

template<int nVarsFixed>
using DNDS::Euler::ArrayDOFV< nVarsFixed >::t_element_mat = typename t_base::t_element_mat

Per-element Eigen matrix/vector type.

Definition at line 97 of file Euler.hpp.

◆ t_self

template<int nVarsFixed>
using DNDS::Euler::ArrayDOFV< nVarsFixed >::t_self = ArrayDOFV<nVarsFixed>

Self type alias for CRTP-style forwarding.

Definition at line 95 of file Euler.hpp.

Member Function Documentation

◆ addTo()

template<int nVarsFixed>
void DNDS::Euler::ArrayDOFV< nVarsFixed >::addTo ( const t_self R,
real  r 
)
inline

Scaled addition: this[i] += r * R[i] for every cell (axpy-style).

Parameters
RSource array to add.
rScalar multiplier applied to R before addition.

Definition at line 153 of file Euler.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ componentWiseNorm1()

template<int nVarsFixed>
Eigen::Vector< real, nVarsFixed > DNDS::Euler::ArrayDOFV< nVarsFixed >::componentWiseNorm1 ( )
inline

Compute the per-component global L1 norm.

Returns a vector where each component j contains \(\sum_i |u_{i,j}|\) summed over all locally-owned cells across all ranks via MPI_Allreduce(SUM).

Returns
Vector of component-wise L1 norms (same on all ranks).

Definition at line 269 of file Euler.hpp.

Here is the call graph for this function:

◆ dot() [1/2]

template<int nVarsFixed>
real DNDS::Euler::ArrayDOFV< nVarsFixed >::dot ( const t_self R)
inline

Compute the global dot product <this, R> with MPI reduction.

Sums \(\sum_i \mathbf{u}_i \cdot \mathbf{R}_i\) over locally-owned cells, then calls MPI_Allreduce(SUM).

Parameters
RThe other ArrayDOFV operand.
Returns
Global dot product (same value on all ranks).

Definition at line 312 of file Euler.hpp.

Here is the call graph for this function:

◆ dot() [2/2]

template<int nVarsFixed>
template<class TMultL , class TMultR >
real DNDS::Euler::ArrayDOFV< nVarsFixed >::dot ( const t_self R,
TMultL &&  mL,
TMultR &&  mR 
)
inline

Compute a weighted global dot product with per-component multipliers.

Calculates \(\sum_i (u_i \circ m_L) \cdot (R_i \circ m_R)\) where \(\circ\) denotes component-wise (Hadamard) multiplication by the multiplier arrays. Performs MPI_Allreduce(SUM) for the global result.

This is useful for preconditioned inner products in iterative solvers where different variable components require different scaling.

Template Parameters
TMultLType of the left multiplier (typically Eigen array expression).
TMultRType of the right multiplier (typically Eigen array expression).
Parameters
RThe other ArrayDOFV operand.
mLPer-component multiplier applied to this entries.
mRPer-component multiplier applied to R entries.
Returns
Global weighted dot product (same value on all ranks).
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 338 of file Euler.hpp.

Here is the call graph for this function:

◆ min()

template<int nVarsFixed>
Eigen::Vector< real, nVarsFixed > DNDS::Euler::ArrayDOFV< nVarsFixed >::min ( )
inline

Compute the global component-wise minimum across all MPI ranks.

For each conservative variable component, finds the minimum value across all locally-owned cells on all ranks. Uses MPI_Allreduce(MIN).

Useful for checking physical validity (e.g. minimum density or pressure).

Returns
Vector of per-component global minimums (same on all ranks).
Note
Uses OpenMP custom reduction (EigenVecMin) when OMP is enabled.
Only iterates over father (locally-owned) cells, not ghost cells.

Definition at line 287 of file Euler.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ norm2()

template<int nVarsFixed>
real DNDS::Euler::ArrayDOFV< nVarsFixed >::norm2 ( )
inline

Compute the global L2 norm across all MPI ranks.

Calculates \(\sqrt{\sum_i \|\mathbf{u}_i\|^2}\) over all locally-owned cells, then performs an MPI_Allreduce(SUM) to obtain the global result.

Returns
Global L2 norm (same value on all ranks).

Definition at line 255 of file Euler.hpp.

Here is the call graph for this function:

◆ operator*=() [1/5]

template<int nVarsFixed>
template<int nVarsFixed_T = nVarsFixed>
std::enable_if_t<!(nVarsFixed_T==1)> DNDS::Euler::ArrayDOFV< nVarsFixed >::operator*= ( const ArrayDOFV< 1 > &  R)
inline

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

Each cell's vector is scaled by the single scalar stored in the corresponding cell of R. Disabled via SFINAE when nVarsFixed == 1 (to avoid ambiguity with the self-multiplication overload).

Template Parameters
nVarsFixed_TSFINAE helper; defaults to nVarsFixed.
Parameters
RScalar-per-cell array (ArrayDOFV<1>).

Definition at line 193 of file Euler.hpp.

Here is the call graph for this function:

◆ operator*=() [2/5]

template<int nVarsFixed>
void DNDS::Euler::ArrayDOFV< nVarsFixed >::operator*= ( const Eigen::Vector< real, nVarsFixed > &  R)
inline

Component-wise multiplication by a constant vector.

Multiplies each cell's DOF vector component-wise by R: this[i] = this[i].cwiseProduct(R).

Parameters
RPer-variable multiplier vector applied identically to every cell.

Definition at line 224 of file Euler.hpp.

Here is the call graph for this function:

◆ operator*=() [3/5]

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

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

Multiplies each cell's DOF vector by the corresponding scalar in R: this[i] *= R[i]. Useful for applying cell-volume or time-step scaling.

Parameters
RVector of per-cell scalar multipliers. Must have size >= number of locally-owned cells.
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 170 of file Euler.hpp.

Here is the call graph for this function:

◆ operator*=() [4/5]

template<int nVarsFixed>
void DNDS::Euler::ArrayDOFV< nVarsFixed >::operator*= ( const t_self R)
inline

Element-wise (Hadamard) multiplication: this[i] = this[i] .* R[i].

Parameters
RSource array; each cell's vector is multiplied component-wise.

Definition at line 233 of file Euler.hpp.

Here is the call graph for this function:

◆ operator*=() [5/5]

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

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

Parameters
RScalar multiplier applied to all components of all cells.

Definition at line 135 of file Euler.hpp.

Here is the call graph for this function:

◆ operator+=() [1/3]

template<int nVarsFixed>
void DNDS::Euler::ArrayDOFV< nVarsFixed >::operator+= ( const Eigen::Vector< real, nVarsFixed > &  R)
inline

Add a constant vector to every cell: this[i] += R for all cells.

Parameters
RVector added uniformly to each cell's DOF.

Definition at line 202 of file Euler.hpp.

Here is the call graph for this function:

◆ operator+=() [2/3]

template<int nVarsFixed>
void DNDS::Euler::ArrayDOFV< nVarsFixed >::operator+= ( const t_self R)
inline

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

Parameters
RSource array with the same distribution layout.

Definition at line 119 of file Euler.hpp.

Here is the call graph for this function:

◆ operator+=() [3/3]

template<int nVarsFixed>
void DNDS::Euler::ArrayDOFV< nVarsFixed >::operator+= ( real  R)
inline

Add a scalar to every component of every cell.

Parameters
RScalar value to add uniformly.

Definition at line 211 of file Euler.hpp.

Here is the call graph for this function:

◆ operator-=()

template<int nVarsFixed>
void DNDS::Euler::ArrayDOFV< nVarsFixed >::operator-= ( const t_self R)
inline

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

Parameters
RSource array with the same distribution layout.

Definition at line 127 of file Euler.hpp.

Here is the call graph for this function:

◆ operator/=()

template<int nVarsFixed>
void DNDS::Euler::ArrayDOFV< nVarsFixed >::operator/= ( const t_self R)
inline

Element-wise division: this[i] = this[i] ./ R[i].

Parameters
RDivisor array; each cell's vector is divided component-wise.

Definition at line 242 of file Euler.hpp.

Here is the call graph for this function:

◆ operator=()

template<int nVarsFixed>
void DNDS::Euler::ArrayDOFV< nVarsFixed >::operator= ( const t_self R)
inline

Deep copy assignment from another ArrayDOFV.

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

Definition at line 143 of file Euler.hpp.

Here is the call graph for this function:

◆ setConstant() [1/2]

template<int nVarsFixed>
void DNDS::Euler::ArrayDOFV< nVarsFixed >::setConstant ( const Eigen::Ref< const t_element_mat > &  R)
inline

Set all DOF entries to copies of a given vector.

Parameters
RReference vector; each cell's DOF vector is set to this value.

Definition at line 111 of file Euler.hpp.

Here is the call graph for this function:

◆ setConstant() [2/2]

template<int nVarsFixed>
void DNDS::Euler::ArrayDOFV< nVarsFixed >::setConstant ( real  R)
inline

Set all DOF entries to a uniform scalar value.

Parameters
RThe scalar value to assign to every component of every cell.

Definition at line 103 of file Euler.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ setMaxWith()

template<int nVarsFixed>
template<class TR >
void DNDS::Euler::ArrayDOFV< nVarsFixed >::setMaxWith ( TR  R)
inline

Clamp all components from below: this[i] = max(this[i], R).

Applies a component-wise maximum with the given value R, ensuring no component drops below R. Operates over all cells including ghosts.

Template Parameters
TRType of the lower bound (scalar or array expression broadcastable by Eigen's .max()).
Parameters
RLower bound 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 366 of file Euler.hpp.

Here is the call graph for this function:

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