Loading...
Searching...
No Matches
ITHACAPOD Namespace Reference

namespace for the computation of the POD, it exploits bot the method of snapshots and SVD. More...

Classes

struct  indexTri
struct  indexSquare
class  PODTemplate
class  PODTemplateH1

Functions

void weightedGramSchmidt (Eigen::MatrixXd &matrix, Eigen::VectorXd &weights)
 Weighted Gram-Schmidt orthogonalization.
template<class Type, template< class > class PatchField, class GeoMesh>
void getNestedSnapshotMatrix (PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, label Npar, label NnestedOut)
 Nested-POD approach.
template void getNestedSnapshotMatrix (PtrList< volScalarField > &snapshots, PtrList< volScalarField > &ModesGlobal, word fieldName, label Npar, label NnestedOut)
template void getNestedSnapshotMatrix (PtrList< volVectorField > &snapshots, PtrList< volVectorField > &ModesGlobal, word fieldName, label Npar, label NnestedOut)
template<class Type, template< class > class PatchField, class GeoMesh>
void getModes (PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex=0, bool sup=0, label nmodes=0, bool correctBC=true)
 Computes the bases or reads them for a field.
template void getModes (PtrList< volVectorField > &snapshots, PtrList< volVectorField > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template void getModes (PtrList< volScalarField > &snapshots, PtrList< volScalarField > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template void getModes (PtrList< surfaceScalarField > &snapshots, PtrList< surfaceScalarField > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template void getModes (PtrList< pointVectorField > &snapshots, PtrList< pointVectorField > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template<class Type, template< class > class PatchField, class GeoMesh>
void getModesMemoryEfficient (GeometricField< Type, PatchField, GeoMesh > &templateField, word snapshotsPath, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex=false, bool supex=false, bool sup=false, label nmodes=0, bool correctBC=true, autoPtr< GeometricField< Type, PatchField, GeoMesh > > meanField=NULL)
 Gets the modes in a memory-efficient manner.
template void getModesMemoryEfficient (GeometricField< scalar, fvPatchField, volMesh > &, word, PtrList< GeometricField< scalar, fvPatchField, volMesh > > &, word, bool, bool, bool, label, bool, autoPtr< GeometricField< scalar, fvPatchField, volMesh > >)
template void getModesMemoryEfficient (GeometricField< vector, fvPatchField, volMesh > &, word, PtrList< GeometricField< vector, fvPatchField, volMesh > > &, word, bool, bool, bool, label, bool, autoPtr< GeometricField< vector, fvPatchField, volMesh > >)
template<class Type, template< class > class PatchField, class GeoMesh>
void getMeanMemoryEfficient (GeometricField< Type, PatchField, GeoMesh > &templateField, word snapshotsPath, autoPtr< GeometricField< Type, PatchField, GeoMesh > > &meanField, bool meanex=false)
 Gets the mean field in a memory-efficient manner.
template void getMeanMemoryEfficient (GeometricField< scalar, fvPatchField, volMesh > &, word snapshotsPath, autoPtr< GeometricField< scalar, fvPatchField, volMesh > > &, bool)
template void getMeanMemoryEfficient (GeometricField< vector, fvPatchField, volMesh > &, word snapshotsPath, autoPtr< GeometricField< vector, fvPatchField, volMesh > > &, bool)
template<class Type, template< class > class PatchField, class GeoMesh>
void getWeightedModes (PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex=0, bool sup=0, label nmodes=0, bool correctBC=true)
 Computes the weighted bases (using the nested-pod approach) or read them for a vector field.
template void getWeightedModes (PtrList< volScalarField > &snapshots, PtrList< volScalarField > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template void getWeightedModes (PtrList< volVectorField > &snapshots, PtrList< volVectorField > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template<class Type, template< class > class PatchField, class GeoMesh>
void getModesSVD (PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex=0, bool sup=0, label nmodes=0, bool correctBC=true)
 Gets the bases for a scalar field using SVD instead of the method of snapshots.
template void getModesSVD (PtrList< volScalarField > &snapshots, PtrList< volScalarField > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template void getModesSVD (PtrList< volVectorField > &snapshots, PtrList< volVectorField > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template<>
Eigen::MatrixXd corMatrix (PtrList< volScalarField > &snapshots)
 Construct the Correlation Matrix for Scalar Field.
template<>
Eigen::MatrixXd corMatrix (PtrList< volVectorField > &snapshots)
 Construct the Correlation Matrix for Vector Field.
template<>
Eigen::MatrixXd corMatrix (List< Eigen::SparseMatrix< double > > &snapshots)
 Construct the Correlation Matrix for Vector Field.
template<>
Eigen::MatrixXd corMatrix (List< Eigen::VectorXd > &snapshots)
 Construct the Correlation Matrix for Vector Field.
template<class Type, template< class > class PatchField, class GeoMesh>
void exportBases (PtrList< GeometricField< Type, PatchField, GeoMesh > > &s, PtrList< GeometricField< Type, PatchField, GeoMesh > > &bases, word fieldName, bool sup)
 Export the Bases.
template void exportBases (PtrList< volVectorField > &s, PtrList< volVectorField > &bases, word fieldName, bool sup)
template void exportBases (PtrList< volScalarField > &s, PtrList< volScalarField > &bases, word fieldName, bool sup)
void exportEigenvalues (scalarField Eigenvalues, fileName name, bool sup=0)
 Exports the eigenvalues as a txt file.
void exportcumEigenvalues (scalarField cumEigenvalues, fileName name, bool sup=0)
 Export the cumulative eigenvalues as a txt file.
std::tuple< List< Eigen::SparseMatrix< double > >, List< Eigen::VectorXd > > DEIMmodes (List< Eigen::SparseMatrix< double > > &A, List< Eigen::VectorXd > &b, label nmodesA, label nmodesB, word MatrixName="Matrix")
 Get the DEIM modes for a generic non a parametrized matrix coming from a differential operator function.
void GrammSchmidt (Eigen::MatrixXd &Matrix)
 Performs GrammSchmidt orthonormalization on an Eigen Matrix.
template<class Type, template< class > class PatchField, class GeoMesh>
void getModes (PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, PtrList< volScalarField > &Volumes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template void getModes (PtrList< volScalarField > &snapshots, PtrList< volScalarField > &modes, PtrList< volScalarField > &Volumes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template void getModes (PtrList< volVectorField > &snapshots, PtrList< volVectorField > &modes, PtrList< volScalarField > &Volumes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template<typename type_matrix>
std::tuple< List< Eigen::SparseMatrix< double > >, List< Eigen::VectorXd > > DEIMmodes (PtrList< type_matrix > &MatrixList, label nmodesA, label nmodesB, word MatrixName="Matrix")
 Get the DEIM modes for a generic non-parametrized matrix coming from a differential operator function.
template std::tuple< List< Eigen::SparseMatrix< double > >, List< Eigen::VectorXd > > DEIMmodes (PtrList< fvScalarMatrix > &MatrixList, label nmodesA, label nmodesB, word MatrixName)
template std::tuple< List< Eigen::SparseMatrix< double > >, List< Eigen::VectorXd > > DEIMmodes (PtrList< fvVectorMatrix > &MatrixList, label nmodesA, label nmodesB, word MatrixName)
template<class Type, template< class > class PatchField, class GeoMesh>
PtrList< GeometricField< Type, PatchField, GeoMesh > > DEIMmodes (PtrList< GeometricField< Type, PatchField, GeoMesh > > &SnapShotsMatrix, label nmodes, word FunctionName, word FieldName)
 Get the DEIM modes for a generic non linear function.
template<class Field_type, class Field_type_2>
void getModes (PtrList< Field_type > &snapshots, PtrList< Field_type > &modes, PtrList< Field_type_2 > &fields2, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template void getModes (PtrList< surfaceScalarField > &snapshots, PtrList< surfaceScalarField > &modes, PtrList< volVectorField > &fields2, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template void getModes (PtrList< volScalarField > &snapshots, PtrList< volScalarField > &modes, PtrList< volVectorField > &fields2, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
template PtrList< volScalarField > DEIMmodes (PtrList< volScalarField > &SnapShotsMatrix, label nmodes, word FunctionName, word FieldName)
template PtrList< volVectorField > DEIMmodes (PtrList< volVectorField > &SnapShotsMatrix, label nmodes, word FunctionName, word FieldName)
template<>
scalar computeInnerProduct (const GeometricField< scalar, fvPatchField, volMesh > &field1, const GeometricField< scalar, fvPatchField, volMesh > &field2)
template<>
scalar computeInnerProduct (const GeometricField< vector, fvPatchField, volMesh > &field1, const GeometricField< vector, fvPatchField, volMesh > &field2)
template<>
scalar computeFrobeniusInnerProduct (const GeometricField< scalar, fvPatchField, volMesh > &field1, const GeometricField< scalar, fvPatchField, volMesh > &field2)
template<>
scalar computeFrobeniusInnerProduct (const GeometricField< vector, fvPatchField, volMesh > &field1, const GeometricField< vector, fvPatchField, volMesh > &field2)
template<class Field_type>
Eigen::MatrixXd corMatrix (Field_type &snapshots)
 Computes the correlation matrix given a vector field snapshot Matrix using different norms depending on the input snapshots.
template<class Type, template< class > class PatchField, class GeoMesh>
void getModes (PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, PtrList< volScalarField > &Volumes, word fieldName, bool podex, bool supex=0, bool sup=0, label nmodes=0)
 Gets the modes.
template<class Field_type, class Field_type_2>
void getModes (PtrList< Field_type > &snapshots, PtrList< Field_type > &modes, PtrList< Field_type_2 > &fields2, word fieldName, bool podex, bool supex=0, bool sup=0, label nmodes=0)
 Gets Modes based on eigenvalue decomposition of another Field.
template<class Type, template< class > class PatchField, class GeoMesh>
scalar computeInnerProduct (const GeometricField< Type, PatchField, GeoMesh > &field1, const GeometricField< Type, PatchField, GeoMesh > &field2)
 Helper function to compute inner product value for different field types.
template<class Type, template< class > class PatchField, class GeoMesh>
scalar computeFrobeniusInnerProduct (const GeometricField< Type, PatchField, GeoMesh > &field1, const GeometricField< Type, PatchField, GeoMesh > &field2)
 Helper function to compute Frobenius norm inner product.
void computeLift (PtrList< volTensorField > &Lfield, PtrList< volTensorField > &liftfield, PtrList< volTensorField > &omfield, Eigen::MatrixXi inletIndex)
void computeLift (PtrList< volVectorField > &Lfield, PtrList< volVectorField > &liftfield, PtrList< volVectorField > &omfield, Eigen::MatrixXi inletIndex)
void computeLift (PtrList< volScalarField > &Lfield, PtrList< volScalarField > &liftfield, PtrList< volScalarField > &omfield, Eigen::MatrixXi inletIndex)
void computeLift (Foam::PtrList< Foam::volScalarField > &Lfield, Foam::PtrList< Foam::volScalarField > &liftfield, Foam::PtrList< Foam::volScalarField > &omfield, Eigen::MatrixXi inletIndex)
void computeLift (Foam::PtrList< Foam::volVectorField > &Lfield, Foam::PtrList< Foam::volVectorField > &liftfield, Foam::PtrList< Foam::volVectorField > &omfield, Eigen::MatrixXi inletIndex)
void computeLift (Foam::PtrList< Foam::volTensorField > &Lfield, Foam::PtrList< Foam::volTensorField > &liftfield, Foam::PtrList< Foam::volTensorField > &omfield, Eigen::MatrixXi inletIndex)

Detailed Description

namespace for the computation of the POD, it exploits bot the method of snapshots and SVD.


Class Documentation

◆ ITHACAPOD::indexTri

struct ITHACAPOD::indexTri

Definition at line 25 of file PODTemplate.H.

Class Members
label index_end
label index_start

◆ ITHACAPOD::indexSquare

struct ITHACAPOD::indexSquare

Definition at line 33 of file PODTemplate.H.

Class Members
label index1_end
label index1_start
label index2_end
label index2_start

Function Documentation

◆ computeFrobeniusInnerProduct() [1/3]

template<>
scalar ITHACAPOD::computeFrobeniusInnerProduct ( const GeometricField< scalar, fvPatchField, volMesh > & field1,
const GeometricField< scalar, fvPatchField, volMesh > & field2 )

Definition at line 2117 of file ITHACAPOD.C.

◆ computeFrobeniusInnerProduct() [2/3]

template<class Type, template< class > class PatchField, class GeoMesh>
scalar ITHACAPOD::computeFrobeniusInnerProduct ( const GeometricField< Type, PatchField, GeoMesh > & field1,
const GeometricField< Type, PatchField, GeoMesh > & field2 )

Helper function to compute Frobenius norm inner product.

Parameters
[in]field1First field
[in]field2Second field
Returns
Inner product value as scalar

◆ computeFrobeniusInnerProduct() [3/3]

template<>
scalar ITHACAPOD::computeFrobeniusInnerProduct ( const GeometricField< vector, fvPatchField, volMesh > & field1,
const GeometricField< vector, fvPatchField, volMesh > & field2 )

Definition at line 2125 of file ITHACAPOD.C.

◆ computeInnerProduct() [1/3]

template<>
scalar ITHACAPOD::computeInnerProduct ( const GeometricField< scalar, fvPatchField, volMesh > & field1,
const GeometricField< scalar, fvPatchField, volMesh > & field2 )

Definition at line 2101 of file ITHACAPOD.C.

◆ computeInnerProduct() [2/3]

template<class Type, template< class > class PatchField, class GeoMesh>
scalar ITHACAPOD::computeInnerProduct ( const GeometricField< Type, PatchField, GeoMesh > & field1,
const GeometricField< Type, PatchField, GeoMesh > & field2 )

Helper function to compute inner product value for different field types.

Parameters
[in]field1First field
[in]field2Second field
Returns
Inner product value as scalar

◆ computeInnerProduct() [3/3]

template<>
scalar ITHACAPOD::computeInnerProduct ( const GeometricField< vector, fvPatchField, volMesh > & field1,
const GeometricField< vector, fvPatchField, volMesh > & field2 )

Definition at line 2109 of file ITHACAPOD.C.

◆ computeLift() [1/3]

void ITHACAPOD::computeLift ( PtrList< volScalarField > & Lfield,
PtrList< volScalarField > & liftfield,
PtrList< volScalarField > & omfield,
Eigen::MatrixXi inletIndex )

Definition at line 1218 of file PODTemplate.C.

◆ computeLift() [2/3]

void ITHACAPOD::computeLift ( PtrList< volTensorField > & Lfield,
PtrList< volTensorField > & liftfield,
PtrList< volTensorField > & omfield,
Eigen::MatrixXi inletIndex )

Definition at line 1142 of file PODTemplate.C.

◆ computeLift() [3/3]

void ITHACAPOD::computeLift ( PtrList< volVectorField > & Lfield,
PtrList< volVectorField > & liftfield,
PtrList< volVectorField > & omfield,
Eigen::MatrixXi inletIndex )

Definition at line 1180 of file PODTemplate.C.

◆ corMatrix() [1/5]

template<class Field_type>
Eigen::MatrixXd ITHACAPOD::corMatrix ( Field_type & snapshots)

Computes the correlation matrix given a vector field snapshot Matrix using different norms depending on the input snapshots.

Parameters
[in]snapshotsList of snapshots.
Template Parameters
Field_typeType of the input list you are passing, it can be: PtrList<volVectorField>, PtrList<volScalarField>, List<Eigen::SparseMatrix<double>>, List<Eigen::VectorXd>.
Returns
the Eigen::MatrixXd correlation matrix.

◆ corMatrix() [2/5]

template<>
Eigen::MatrixXd ITHACAPOD::corMatrix ( List< Eigen::SparseMatrix< double > > & snapshots)

Construct the Correlation Matrix for Vector Field.

Definition at line 1010 of file ITHACAPOD.C.

◆ corMatrix() [3/5]

template<>
Eigen::MatrixXd ITHACAPOD::corMatrix ( List< Eigen::VectorXd > & snapshots)

Construct the Correlation Matrix for Vector Field.

Definition at line 1045 of file ITHACAPOD.C.

◆ corMatrix() [4/5]

template<>
Eigen::MatrixXd ITHACAPOD::corMatrix ( PtrList< volScalarField > & snapshots)

Construct the Correlation Matrix for Scalar Field.

Definition at line 955 of file ITHACAPOD.C.

◆ corMatrix() [5/5]

template<>
Eigen::MatrixXd ITHACAPOD::corMatrix ( PtrList< volVectorField > & snapshots)

Construct the Correlation Matrix for Vector Field.

Definition at line 983 of file ITHACAPOD.C.

◆ DEIMmodes() [1/3]

std::tuple< List< Eigen::SparseMatrix< double > >, List< Eigen::VectorXd > > ITHACAPOD::DEIMmodes ( List< Eigen::SparseMatrix< double > > & A,
List< Eigen::VectorXd > & b,
label nmodesA,
label nmodesB,
word MatrixName = "Matrix" )

Get the DEIM modes for a generic non a parametrized matrix coming from a differential operator function.

Parameters
AList of matrices for the lhs
bList of vectors for the rhs
[in]nmodesAThe number of modes
[in]nmodesBThe nmodes b
[in]MatrixNameThe matrix name, used to save the matrix
[in]MatrixListThe matrix list as a list of pointers
Template Parameters
Thetype of matrix, can be fvScalarMatrix or fvVectorMatrix
Returns
It returns a tuple containing the list of modes for the matrix and the list of modes for the source term

Definition at line 1172 of file ITHACAPOD.C.

◆ DEIMmodes() [2/3]

template<class Type, template< class > class PatchField, class GeoMesh>
PtrList< GeometricField< Type, PatchField, GeoMesh > > ITHACAPOD::DEIMmodes ( PtrList< GeometricField< Type, PatchField, GeoMesh > > & SnapShotsMatrix,
label nmodes,
word FunctionName,
word FieldName )

Get the DEIM modes for a generic non linear function.

Parameters
[in]SnapShotsMatrixThe snapshots matrix
[in]nmodesThe number of modes
[in]FunctionNameThe function name
Template Parameters
Typetemplated object for the snapshots matrix type
Returns
Type the POD modes

Definition at line 1740 of file ITHACAPOD.C.

◆ DEIMmodes() [3/3]

template<typename type_matrix>
std::tuple< List< Eigen::SparseMatrix< double > >, List< Eigen::VectorXd > > ITHACAPOD::DEIMmodes ( PtrList< type_matrix > & MatrixList,
label nmodesA,
label nmodesB,
word MatrixName = "Matrix" )

Get the DEIM modes for a generic non-parametrized matrix coming from a differential operator function.

Parameters
[in]MatrixListThe matrix list as a list of pointers
[in]nmodesAThe number of modes for A
[in]nmodesBThe number of modes for B
[in]MatrixNameThe matrix name, used to save the matrix
Template Parameters
type_matrixThe type of matrix, can be fvScalarMatrix or fvVectorMatrix
Returns
It returns a tuple containing the list of modes for the matrix and the list of modes for the source term

Definition at line 1555 of file ITHACAPOD.C.

◆ exportBases()

template<class Type, template< class > class PatchField, class GeoMesh>
void ITHACAPOD::exportBases ( PtrList< GeometricField< Type, PatchField, GeoMesh > > & snapshots,
PtrList< GeometricField< Type, PatchField, GeoMesh > > & bases,
word fieldName,
bool sup = 0 )

Export the Bases.

Exports the basis for an OpenFOAM GeometricField into the ITHACAOutput/POD or ITHACAOutput/supremizer.

Parameters
snapshotsThe snapshots
[in]basesa PtrList of GeometricField where the bases are stored.
[in]fieldNameThe field name
[in]supa boolean variable 1 if you want to export the supremizer bases (in ITHACAOutput/supremizer) 0 elsewhere (Default is 0).
[in]sa PtrList of GeometricField where the snapshots associated with the bases are stored.
Template Parameters
Typevector or scalar.
PatchFieldfvPatchField or fvsPatchField.
GeoMeshvolMesh or surfaceMesh.

Definition at line 1074 of file ITHACAPOD.C.

◆ exportcumEigenvalues()

void ITHACAPOD::exportcumEigenvalues ( scalarField cumEigenvalues,
fileName name,
bool sup = 0 )

Export the cumulative eigenvalues as a txt file.

Parameters
[in]cumEigenvaluesa scalarField of cumEigenvalues.
[in]namea fileName to indicate the name of the txt file.
[in]supIf 1, it exports the cumEigenValues in "ITHACAOutput/supremizer". If 0, no actions (Default is 0).

Definition at line 1141 of file ITHACAPOD.C.

◆ exportEigenvalues()

void ITHACAPOD::exportEigenvalues ( scalarField Eigenvalues,
fileName name,
bool sup = 0 )

Exports the eigenvalues as a txt file.

Parameters
[in]Eigenvaluesa scalarField of eigenvalues.
[in]namea fileName to indicate the name of the txt file.
[in]supIf 1, it exports the eigenValues in "ITHACAOutput/supremizer". If 0, no actions (Default is 0).

Definition at line 1112 of file ITHACAPOD.C.

◆ getMeanMemoryEfficient()

template<class Type, template< class > class PatchField, class GeoMesh>
void ITHACAPOD::getMeanMemoryEfficient ( GeometricField< Type, PatchField, GeoMesh > & templateField,
word snapshotsPath,
autoPtr< GeometricField< Type, PatchField, GeoMesh > > & meanField,
bool meanex )

Gets the mean field in a memory-efficient manner.

Definition at line 660 of file ITHACAPOD.C.

◆ getModes() [1/5]

template<class Field_type, class Field_type_2>
void ITHACAPOD::getModes ( PtrList< Field_type > & snapshots,
PtrList< Field_type > & modes,
PtrList< Field_type_2 > & fields2,
word fieldName,
bool podex,
bool supex,
bool sup,
label nmodes,
bool correctBC )

Definition at line 1940 of file ITHACAPOD.C.

◆ getModes() [2/5]

template<class Field_type, class Field_type_2>
void ITHACAPOD::getModes ( PtrList< Field_type > & snapshots,
PtrList< Field_type > & modes,
PtrList< Field_type_2 > & fields2,
word fieldName,
bool podex,
bool supex = 0,
bool sup = 0,
label nmodes = 0 )

Gets Modes based on eigenvalue decomposition of another Field.

Parameters
snapshotsThe snapshots of which we want to create the modes
modesThe modes of the snapshots
fields2The fields that we want to use for the eigenvalue decomposition
[in]fieldNameThe field name
[in]podexboolean variable 1 if the modes have been already computed and stored (in this case the function is reading them) 0 elsewhere.
[in]supexboolean variable 1 if the supremizer modes have been already computed and stored (in this case the function is reading them) 0 elsewhere.
[in]supboolean variable 1 if you want to compute the supremizer modes 0 elsewhere.
[in]nmodeslabel variable to set the number of modes to be stored, if set to 0 the maximum number of modes will computed.
Template Parameters
Field_typeType of field snapshots
Field_type_2Type of the other fields

◆ getModes() [3/5]

template<class Type, template< class > class PatchField, class GeoMesh>
void ITHACAPOD::getModes ( PtrList< GeometricField< Type, PatchField, GeoMesh > > & snapshots,
PtrList< GeometricField< Type, PatchField, GeoMesh > > & modes,
PtrList< volScalarField > & Volumes,
word fieldName,
bool podex,
bool supex,
bool sup,
label nmodes,
bool correctBC )

Definition at line 1392 of file ITHACAPOD.C.

◆ getModes() [4/5]

template<class Type, template< class > class PatchField, class GeoMesh>
void ITHACAPOD::getModes ( PtrList< GeometricField< Type, PatchField, GeoMesh > > & snapshots,
PtrList< GeometricField< Type, PatchField, GeoMesh > > & modes,
PtrList< volScalarField > & Volumes,
word fieldName,
bool podex,
bool supex = 0,
bool sup = 0,
label nmodes = 0 )

Gets the modes.

Parameters
[in]snapshotsList of snapshots.
[out]modesA PtrList where modes are stored (it must be passed empty).
VolumesThe volumes
[in]fieldNameThe field name
[in]podexIf 1, the functions read the stored mode. If 0, the function computes the modes and stores them.
[in]supexIf 1, the functions read the stored supremizer mode. If 0, the function computes and stores them.
[in]supIf 1 it computes the supremizer modes.
[in]nmodesNumber of modes to be stored. If 0, the maximum number of modes will computed.
Template Parameters
Typevector or scalar.
PatchFieldfvPatchField or fvsPatchField.
GeoMeshvolMesh or surfaceMesh.

◆ getModes() [5/5]

template<class Type, template< class > class PatchField, class GeoMesh>
void ITHACAPOD::getModes ( PtrList< GeometricField< Type, PatchField, GeoMesh > > & snapshots,
PtrList< GeometricField< Type, PatchField, GeoMesh > > & modes,
word fieldName,
bool podex,
bool supex = 0,
bool sup = 0,
label nmodes = 0,
bool correctBC = true )

Computes the bases or reads them for a field.

Parameters
[in]snapshotsList of snapshots.
[out]modesA PtrList where modes are stored (it must be passed empty).
[in]fieldNameThe field name
[in]podexIf 1, the functions read the stored mode. If 0, the function computes the modes and stores them.
[in]supexIf 1, the functions read the stored supremizer mode. If 0, the function computes and stores them.
[in]supIf 1 it computes the supremizer modes.
[in]nmodesNumber of modes to be stored. If 0, the maximum number of modes will computed.
Template Parameters
Typevector or scalar.
PatchFieldfvPatchField or fvsPatchField.
GeoMeshvolMesh or surfaceMesh.

Definition at line 93 of file ITHACAPOD.C.

◆ getModesMemoryEfficient()

template<class Type, template< class > class PatchField, class GeoMesh>
void ITHACAPOD::getModesMemoryEfficient ( GeometricField< Type, PatchField, GeoMesh > & templateField,
word snapshotsPath,
PtrList< GeometricField< Type, PatchField, GeoMesh > > & modes,
word fieldName,
bool podex = false,
bool supex = false,
bool sup = false,
label nmodes = 0,
bool correctBC = true,
autoPtr< GeometricField< Type, PatchField, GeoMesh > > meanField = NULL )

Gets the modes in a memory-efficient manner.

Parameters
[in]templateFieldThe template field
[in]snapshotsPathThe path to the snapshots
[out]modesThe modes
[in]fieldNameThe field name
[in]podexIf 1, the functions read the stored mode. If 0, the function computes the modes and stores them.
[in]supexIf 1, the functions read the stored supremizer mode. If 0, the function computes and stores them.
[in]supIf 1, it computes the supremizer modes.
[in]nmodesNumber of modes to be stored. If 0, the maximum number of modes will computed.
[in]correctBCIf true, correct the boundary conditions.
Template Parameters
TypeThe type of the field
PatchFieldThe patch field type
GeoMeshThe mesh type

Definition at line 345 of file ITHACAPOD.C.

◆ getModesSVD()

template<class Type, template< class > class PatchField, class GeoMesh>
void ITHACAPOD::getModesSVD ( PtrList< GeometricField< Type, PatchField, GeoMesh > > & snapshots,
PtrList< GeometricField< Type, PatchField, GeoMesh > > & modes,
word fieldName,
bool podex,
bool supex = 0,
bool sup = 0,
label nmodes = 0,
bool correctBC = true )

Gets the bases for a scalar field using SVD instead of the method of snapshots.

Parameters
[in]snapshotsList of snapshots.
[out]modesA PtrList where modes are stored (it must be passed empty).
[in]fieldNameThe file name
[in]podexIf 1, the functions read the stored mode. If 0, the function computes the modes and stores them.
[in]supexIf 1, the functions read the stored supremizer mode. If 0, the function computes and stores them.
[in]supIf 1 it computes the supremizer modes.
[in]nmodesNumber of modes to be stored. If 0, the maximum number of modes will computed.
Template Parameters
Typevector or scalar.
PatchFieldfvPatchField or fvsPatchField.
GeoMeshvolMesh or surfaceMesh.

Definition at line 858 of file ITHACAPOD.C.

◆ getNestedSnapshotMatrix()

template<class Type, template< class > class PatchField, class GeoMesh>
void ITHACAPOD::getNestedSnapshotMatrix ( PtrList< GeometricField< Type, PatchField, GeoMesh > > & snapshots,
PtrList< GeometricField< Type, PatchField, GeoMesh > > & modes,
word fieldName,
label Npar,
label NnestedOut )

Nested-POD approach.

Computes the nested snapshot matrix and weighted bases for a vector field

Parameters
[in]snapshotsList of snapshots.
modesThe modes
[in]fieldNameThe field name
[in]NparNumber of parameters
[in]NnestedOutNumber of weighted modes to be stored for the assembly of the nested matrix
[out]UModesGlobalList of snapshots where modes are stored (it must be passed empty).
Template Parameters
Typevector or scalar.
PatchFieldfvPatchField or fvsPatchField.
GeoMeshvolMesh or surfaceMesh.

Definition at line 41 of file ITHACAPOD.C.

◆ getWeightedModes()

template<class Type, template< class > class PatchField, class GeoMesh>
void ITHACAPOD::getWeightedModes ( PtrList< GeometricField< Type, PatchField, GeoMesh > > & snapshots,
PtrList< GeometricField< Type, PatchField, GeoMesh > > & modes,
word fieldName,
bool podex,
bool supex = 0,
bool sup = 0,
label nmodes = 0,
bool correctBC = true )

Computes the weighted bases (using the nested-pod approach) or read them for a vector field.

Parameters
[in]snapshotsList of snapshots.
[out]modesA PtrList where modes are stored (it must be passed empty).
[in]fieldNameThe field name
[in]podexIf 1, the functions read the stored mode. If 0, the function computes the modes and stores them.
[in]supexIf 1, the functions read the stored supremizer mode. If 0, the function computes and stores them.
[in]supIf 1, it computes the supremizer modes.
[in]nmodesNumber of modes to be stored. If 0, the maximum number of modes will computed.
Template Parameters
Typevector or scalar.
PatchFieldfvPatchField or fvsPatchField.
GeoMeshvolMesh or surfaceMesh.

Definition at line 715 of file ITHACAPOD.C.

◆ GrammSchmidt()

void ITHACAPOD::GrammSchmidt ( Eigen::MatrixXd & Matrix)

Performs GrammSchmidt orthonormalization on an Eigen Matrix.

Parameters
[out]MatrixThe matrix

Definition at line 1370 of file ITHACAPOD.C.

◆ weightedGramSchmidt()

void ITHACAPOD::weightedGramSchmidt ( Eigen::MatrixXd & matrix,
Eigen::VectorXd & weights )

Weighted Gram-Schmidt orthogonalization.

Parameters
[in]matrixMatrix to be orthogonalized
[in]weightsVector of weights

Definition at line 190 of file incrementalPOD.H.