libROM  v1.0
Data-driven physical simulation library
Public Member Functions | List of all members
CAROM::SampleMeshManager Class Reference

#include <SampleMesh.hpp>

Public Member Functions

 SampleMeshManager (vector< ParFiniteElementSpace * > &fespace_, string visFileName="", double visScale=1.0)
 Constructor creating a SampleMeshManager with arbitrarily many ParFiniteElementSpaces on a full-order mesh. More...
 
void RegisterSampledVariable (const string variable, const int space, vector< int > const &sample_dofs_v, vector< int > const &num_sample_dofs_per_proc)
 Register a variable and set its sampled DOFs. More...
 
void ConstructSampleMesh ()
 Construct the sample mesh, after registering all sampled variables.
 
void GatherDistributedMatrixRows (const string variable, CAROM::Matrix const &B, const int rdim, CAROM::Matrix &Bsp) const
 Gather the rows of a distributed basis CAROM::Matrix corresponding to a sample mesh space. More...
 
ParFiniteElementSpace * GetSampleFESpace (const int space) const
 Returns a sample mesh space. More...
 
ParMesh * GetSampleMesh () const
 Returns the sample mesh. More...
 
int GetNumVarSamples (const string variable) const
 Returns the number of sampled DOFs on all processes for a variable. More...
 
void GetSampledValues (const string variable, mfem::Vector const &v, CAROM::Vector &s) const
 Sets a vector of sampled DOFs from a vector on the sample mesh space. More...
 
set< int > * GetSampleElements ()
 Returns a set of indices of local FOM mesh elements corresponding to sample elements. More...
 
void WriteVariableSampleMap (const string variable, string file_name) const
 Writes a variable sample DOF map to file, which can be read by SampleDOFSelector::ReadMapFromFile in order to use SampleDOFSelector::GetSampledValues when this SampleMeshManager object is not available. More...
 
 ~SampleMeshManager ()
 Destructor.
 

Detailed Description

Class SampleMeshManager constructs and manages a sample mesh and finite element spaces on the sample mesh. The spaces defined on the sample mesh are referred to as sample mesh spaces, and they generally contain the sampled DOFs along with other unsampled DOFs. This class can be used to work with vectors and basis matrices in the sample mesh spaces.

Definition at line 41 of file SampleMesh.hpp.

Constructor & Destructor Documentation

◆ SampleMeshManager()

CAROM::SampleMeshManager::SampleMeshManager ( vector< ParFiniteElementSpace * > &  fespace_,
string  visFileName = "",
double  visScale = 1.0 
)

Constructor creating a SampleMeshManager with arbitrarily many ParFiniteElementSpaces on a full-order mesh.

Precondition
fespace_.size() > 0
Parameters
[in]fespace_Arbitrarily many ParFiniteElementSpaces defined on the same full-order mesh. Each sampled variable must be defined on one of these spaces. All spaces must be specified in the constructor.
[in]visFileNameIf non-empty, this filename is used for VisIt and ParaView output of the sample mesh and the sampled DOFs on the full-order mesh.
[in]visScaleConstant value for indicating the sample mesh elements in visualization.

Definition at line 1069 of file SampleMesh.cpp.

Member Function Documentation

◆ GatherDistributedMatrixRows()

void CAROM::SampleMeshManager::GatherDistributedMatrixRows ( const string  variable,
CAROM::Matrix const &  B,
const int  rdim,
CAROM::Matrix Bsp 
) const

Gather the rows of a distributed basis CAROM::Matrix corresponding to a sample mesh space.

Parameters
[in]variableString containing the name of the variable, used for look-up.
[in]BDistributed matrix with columns containing full-order vectors.
[in]rdimNumber of columns of B to extract (starting with the first).
[out]BspMatrix with rows corresponding to sample mesh space DOFs, gathered only to the root process (MPI rank 0).

Definition at line 1510 of file SampleMesh.cpp.

◆ GetNumVarSamples()

int CAROM::SampleMeshManager::GetNumVarSamples ( const string  variable) const

Returns the number of sampled DOFs on all processes for a variable.

Parameters
[in]variableString containing the name of the variable, used for look-up.
Returns
Number of sampled DOFs on all processes for a variable.

Definition at line 1298 of file SampleMesh.cpp.

◆ GetSampledValues()

void CAROM::SampleMeshManager::GetSampledValues ( const string  variable,
mfem::Vector const &  v,
CAROM::Vector s 
) const

Sets a vector of sampled DOFs from a vector on the sample mesh space.

Parameters
[in]variableString containing the name of the variable, used for look-up.
[in]vVector on the sample mesh space.
[out]sVector of sampled DOFs on all processes.

Definition at line 1304 of file SampleMesh.cpp.

◆ GetSampleElements()

set<int>* CAROM::SampleMeshManager::GetSampleElements ( )
inline

Returns a set of indices of local FOM mesh elements corresponding to sample elements.

Returns
Pointer to a set of local FOM mesh element indices.

Definition at line 142 of file SampleMesh.hpp.

◆ GetSampleFESpace()

ParFiniteElementSpace* CAROM::SampleMeshManager::GetSampleFESpace ( const int  space) const
inline

Returns a sample mesh space.

Parameters
[in]spaceIndex of the finite element space (corresponds to fespace_ in constructor).
Returns
The sample mesh space with the given index.

Definition at line 100 of file SampleMesh.hpp.

◆ GetSampleMesh()

ParMesh* CAROM::SampleMeshManager::GetSampleMesh ( ) const
inline

Returns the sample mesh.

Returns
The sample mesh, defined only on the root process (MPI rank 0).

Definition at line 109 of file SampleMesh.hpp.

◆ RegisterSampledVariable()

void CAROM::SampleMeshManager::RegisterSampledVariable ( const string  variable,
const int  space,
vector< int > const &  sample_dofs_v,
vector< int > const &  num_sample_dofs_per_proc 
)

Register a variable and set its sampled DOFs.

Parameters
[in]variableString containing the name of the variable, used for look-up.
[in]spaceIndex of the space in fespace (see constructor) where the variable is defined.
[in]sample_dofs_vLocal sampled DOF indices on all MPI processes (cf. DEIM function).
[in]num_sample_dofs_per_procNumber of sampled DOFs on each MPI process (cf. DEIM function).

Definition at line 1098 of file SampleMesh.cpp.

◆ WriteVariableSampleMap()

void CAROM::SampleMeshManager::WriteVariableSampleMap ( const string  variable,
string  file_name 
) const

Writes a variable sample DOF map to file, which can be read by SampleDOFSelector::ReadMapFromFile in order to use SampleDOFSelector::GetSampledValues when this SampleMeshManager object is not available.

Parameters
[in]variableString containing the name of the variable, used for look-up.
[in]file_nameName of the output file.

Definition at line 1316 of file SampleMesh.cpp.


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