libROM
v1.0
Data-driven physical simulation library
|
Enumerations | |
enum | SamplingType { deim , gnat , qdeim , sopt } |
enum class | NNLS_termination { LINF , L2 } |
Termination criterion for NNLS solver: LINF: L-infinity norm (maximum value) of the residual L2: L2 norm of the residual. | |
Functions | |
struct GreedyErrorIndicatorPoint | createGreedyErrorIndicatorPoint (Vector *point, Vector *localROM) |
Create a greedy error indicator point. More... | |
struct GreedyErrorIndicatorPoint | createGreedyErrorIndicatorPoint (Vector *point, std::shared_ptr< Vector > &localROM) |
Create a greedy error indicator point. More... | |
Vector | getNearestPoint (std::vector< Vector > paramPoints, Vector point) |
Given a vector, find the nearest point in a domain. More... | |
int | getNearestPointIndex (std::vector< Vector > paramPoints, Vector point) |
Given a vector, find the nearest point in a domain. More... | |
double | getNearestPoint (std::vector< double > paramPoints, double point) |
Given a double, find the nearest point in a domain. More... | |
int | getNearestPointIndex (std::vector< double > paramPoints, double point) |
Given a double, find the nearest point in a domain. More... | |
std::vector< double > | obtainRBFToTrainingPoints (std::vector< Vector * > parameter_points, std::string interp_method, std::string rbf, double epsilon, Vector *point) |
Compute the RBF from the parameter points with the unsampled parameter point. More... | |
double | rbfWeightedSum (std::vector< double > &rbf) |
Compute the sum of the RBF weights. More... | |
double | obtainRBF (std::string rbf, double epsilon, Vector *point1, Vector *point2) |
Compute the RBF between two points. More... | |
double | convertClosestRBFToEpsilon (std::vector< Vector * > parameter_points, std::string rbf, double closest_rbf_val) |
Convert closest RBF value to an epsilon value. More... | |
std::vector< Matrix * > | obtainRotationMatrices (std::vector< Vector * > parameter_points, std::vector< Matrix * > bases, int ref_point) |
Obtain the rotation matrices for all the parameter points using the basis of the reference point. More... | |
Vector * | obtainInterpolatedVector (std::vector< Vector * > data, Matrix *f_T, std::string interp_method, std::vector< double > &rbf) |
Matrix * | solveLinearSystem (std::vector< Vector * > parameter_points, std::vector< Vector * > data, std::string interp_method, std::string rbf, double &epsilon) |
template<class T > | |
void | getParametricDMD (T *¶metric_dmd, std::vector< Vector * > ¶meter_points, std::vector< T * > &dmds, Vector *desired_point, std::string rbf="G", std::string interp_method="LS", double closest_rbf_val=0.9, bool reorthogonalize_W=false) |
Constructor. More... | |
template<class T > | |
void | getParametricDMD (T *¶metric_dmd, std::vector< Vector * > ¶meter_points, std::vector< std::string > &dmd_paths, Vector *desired_point, std::string rbf="G", std::string interp_method="LS", double closest_rbf_val=0.9, bool reorthogonalize_W=false) |
Constructor. More... | |
template<class T > | |
void | getParametricDMDc (T *¶metric_dmdc, std::vector< Vector * > ¶meter_points, std::vector< T * > &dmdcs, std::vector< Matrix * > controls, Matrix *&controls_interpolated, Vector *desired_point, std::string rbf="G", std::string interp_method="LS", double closest_rbf_val=0.9, bool reorthogonalize_W=false) |
Constructor. More... | |
template<class T > | |
void | getParametricDMDc (T *¶metric_dmdc, std::vector< Vector * > ¶meter_points, std::vector< std::string > &dmdc_paths, std::vector< Matrix * > controls, Matrix *&controls_interpolated, Vector *desired_point, std::string rbf="G", std::string interp_method="LS", double closest_rbf_val=0.9, bool reorthogonalize_W=false) |
Constructor. More... | |
void | DEIM (const Matrix *f_basis, int num_f_basis_vectors_used, std::vector< int > &f_sampled_row, std::vector< int > &f_sampled_rows_per_proc, Matrix &f_basis_sampled_inv, int myid, int num_procs) |
Computes the DEIM algorithm on the given basis. More... | |
void | GNAT (const Matrix *f_basis, const int num_f_basis_vectors_used, std::vector< int > &f_sampled_row, std::vector< int > &f_sampled_rows_per_proc, Matrix &f_basis_sampled_inv, const int myid, const int num_procs, const int num_samples_req=-1, std::vector< int > *init_samples=NULL) |
Computes the GNAT algorithm on the given basis. More... | |
void | QDEIM (const Matrix *f_basis, int num_f_basis_vectors_used, std::vector< int > &f_sampled_row, std::vector< int > &f_sampled_rows_per_proc, Matrix &f_basis_sampled_inv, const int myid, const int num_procs, const int num_samples_req) |
Computes the QDEIM algorithm on the given basis. More... | |
void | S_OPT (const Matrix *f_basis, int num_f_basis_vectors_used, std::vector< int > &f_sampled_row, std::vector< int > &f_sampled_rows_per_proc, Matrix &f_basis_sampled_inv, const int myid, const int num_procs, const int num_samples_req=-1, std::vector< int > *init_samples=NULL, bool qr_factorize=false) |
Computes the S_OPT algorithm on the given basis. More... | |
void | SampleTemporalIndices (const Matrix *s_basis, const Matrix *t_basis, const int num_f_basis_vectors_used, int *t_samples, const int myid, const int num_procs, const int num_samples_req, const bool excludeFinalTime) |
void | SampleSpatialIndices (const Matrix *s_basis, const Matrix *t_basis, const int num_f_basis_vectors_used, const int num_t_samples, int *t_samples, int *s_sampled_row, int *s_sampled_rows_per_proc, Matrix &f_basis_sampled, const int myid, const int num_procs, const int num_samples_req) |
void | SpaceTimeSampling (const Matrix *s_basis, const Matrix *t_basis, const int num_f_basis_vectors_used, std::vector< int > &t_samples, int *f_sampled_row, int *f_sampled_rows_per_proc, Matrix &s_basis_sampled, const int myid, const int num_procs, const int num_t_samples_req, const int num_s_samples_req, const bool excludeFinalTime) |
void | GetSampledSpaceTimeBasis (std::vector< int > const &t_samples, const Matrix *t_basis, Matrix const &s_basis_sampled, Matrix &f_basis_sampled_inv) |
Matrix | outerProduct (const Vector &v, const Vector &w) |
Computes the outer product of two Vectors, v and w. More... | |
Matrix | DiagonalMatrixFactory (const Vector &v) |
Factory function to make a diagonal matrix with nonzero entries as in its Vector argument. The rows of this diagonal matrix are distributed in the same fashion as its Vector argument. More... | |
Matrix | IdentityMatrixFactory (const Vector &v) |
Factory function to make an identity matrix with rows distributed in the same fashion as its Vector argument. This function is provided for convenience due to the ubiquity of identity matrices (and operators) in mathematics. More... | |
struct EigenPair | SymmetricRightEigenSolve (Matrix *A) |
Computes the eigenvectors/eigenvalues of an NxN real symmetric matrix. More... | |
struct ComplexEigenPair | NonSymmetricRightEigenSolve (Matrix *A) |
Computes the eigenvectors/eigenvalues of an NxN real nonsymmetric matrix. More... | |
void | SerialSVD (Matrix *A, Matrix *U, Vector *S, Matrix *V) |
Computes the SVD of a undistributed matrix in serial. More... | |
struct SerialSVDDecomposition | SerialSVD (Matrix *A) |
Computes the SVD of a undistributed matrix in serial. More... | |
Matrix * | SpaceTimeProduct (const CAROM::Matrix *As, const CAROM::Matrix *At, const CAROM::Matrix *Bs, const CAROM::Matrix *Bt, const std::vector< double > *tscale, const bool At0at0, const bool Bt0at0, const bool lagB, const bool skip0) |
void | blacs_get_ (int *, int *, int *) |
void | blacs_pinfo_ (int *, int *) |
void | blacs_gridinit_ (int *, char *, int *, int *) |
void | blacs_gridinfo_ (int *, int *, int *, int *, int *) |
void | blacs_gridexit_ (int *icontxt) |
void | blacs_freebuff_ (int *icontxt, int *wait) |
void | blacs_exit_ (int *cont) |
void | descinit_ (int *DESC, int *M, int *N, int *MB, int *NB, int *IRSRC, int *ICSRC, int *ICTXT, int *LLD, int *INFO) |
void | pdormqr_ (char *SIDE, char *TRANS, int *M, int *N, int *K, double *A, int *IA, int *JA, int *DESCA, double *TAU, double *C, int *IC, int *JC, int *DESCC, double *WORK, int *LWORK, int *INFO) |
void | pdgeqrf_ (int *M, int *N, double *A, int *IA, int *JA, int *DESCA, double *TAU, double *WORK, int *LWORK, int *INFO) |
void | pdtrsm_ (char *SIDE, char *UPLO, char *TRANS, char *DIAG, int *M, int *N, double *ALPHA, double *A, int *IA, int *JA, int *DESCA, double *B, int *IB, int *JB, int *DESCB) |
void | pdgemr2d_ (int *m, int *n, double *A, int *IA, int *JA, int *descA, double *B, int *IB, int *JB, int *descB, int *gcontext) |
void | pdlacpy_ (char *UPLO, int *M, int *N, double *A, int *IA, int *JA, int *DESCA, double *B, int *IB, int *JB, int *DESCB) |
int | numroc_ (int *N, int *NB, int *IPROC, int *ISRCPROC, int *NPROCS) |
void | pdgemv_ (char *TRANS, int *M, int *N, double *ALPHA, double *A, int *IA, int *JA, int *DESCA, double *X, int *IX, int *JX, int *DESCX, int *INCX, double *BETA, double *Y, int *IY, int *JY, int *DESCY, int *INCY) |
int | getCenterPoint (std::vector< Vector * > &points, bool use_centroid) |
Get center point of a group of points. | |
int | getCenterPoint (std::vector< Vector > &points, bool use_centroid) |
Get center point of a group of points. | |
int | getClosestPoint (std::vector< Vector * > &points, Vector *test_point) |
Get closest point to a test point among a group of points. | |
int | getClosestPoint (std::vector< Vector > &points, Vector test_point) |
Get closest point to a test point among a group of points. | |
void | FindStencilElements (const vector< int > &sample_dofs_gid, set< int > &elements, ParFiniteElementSpace &fespace) |
void | FindSampledMeshEntities (const int type, const vector< int > &sample_dofs_gid, set< int > &entities, ParFiniteElementSpace &fespace) |
void | GetLocalSampleMeshElements (ParMesh &pmesh, ParFiniteElementSpace &fespace, const vector< int > &sample_dofs, const vector< int > &local_num_sample_dofs, set< int > &elems, bool getSampledEntities, set< int > &intElems, set< int > &faces, set< int > &edges, set< int > &vertices) |
void | SplitDofsIntoBlocks (const vector< int > &Ntrue, const vector< int > &dofs, const vector< int > &local_num_dofs, vector< vector< int >> &dofs_block, vector< vector< int >> &dofs_block_todofs, vector< vector< int >> &local_num_dofs_block) |
void | InsertElementDofs (ParFiniteElementSpace &fespace, const int elId, const int offset, set< int > &element_dofs) |
void | AugmentDofListWithOwnedDofs (vector< int > &mixedDofs, vector< ParFiniteElementSpace * > &fespace) |
void | BuildSampleMesh (ParMesh &pmesh, vector< ParFiniteElementSpace * > &fespace, const set< int > &elems, Mesh *&sample_mesh, vector< int > &stencil_dofs, vector< int > &elemLocalIndices, vector< map< int, int > > &elemLocalIndicesInverse) |
void | GetLocalDofsToLocalElementMap (ParFiniteElementSpace &fespace, const vector< int > &dofs, const vector< int > &localNumDofs, const set< int > &elems, vector< int > &dofToElem, vector< int > &dofToElemDof, const bool useTDof) |
void | Set_s2sp (const int myid, const int num_procs, vector< int > const &spNtrue, const int global_num_sample_dofs, const vector< int > &local_num_sample_dofs, const vector< vector< int > > &local_num_sample_dofs_sub, const vector< vector< int > > &localSampleDofsToElem_sub, const vector< vector< int > > &localSampleDofsToElemDof_sub, const vector< vector< int > > &sample_dofs_sub_to_sample_dofs, const vector< map< int, int > > &elemLocalIndicesInverse, vector< ParFiniteElementSpace * > &spfespace, vector< int > &s2sp) |
void | Finish_s2sp_augmented (const int rank, const int nprocs, vector< ParFiniteElementSpace * > &fespace, vector< vector< int >> &dofs_block, vector< vector< int > > &dofs_sub_to_sdofs, vector< vector< int > > &local_num_dofs_sub, const bool dofsTrue, vector< int > &s2sp_) |
void | ParaViewPrintAttributes (const string &fname, Mesh &mesh, int entity_dim, const Array< int > *el_number=nullptr, const Array< int > *vert_number=nullptr) |
void | GatherDistributedMatrixRows_aux (const CAROM::Matrix &B, const int rdim, const int os0, const int os1, const int ossp, ParFiniteElementSpace &fespace, const vector< int > &st2sp, const vector< int > &sprows, const vector< int > &all_sprows, CAROM::Matrix &Bsp) |
void | SampleVisualization (ParMesh *pmesh, set< int > const &elems, set< int > const &intElems, set< int > const &faces, set< int > const &edges, set< int > const &vertices, string const &filename, mfem::Vector *elemCount, double elementScaling) |
void | ComputeCtAB (const HypreParMatrix &A, const CAROM::Matrix &B, const CAROM::Matrix &C, CAROM::Matrix &CtAB) |
This function computes a reduced, non-distributed matrix C^t AB stored identically (redundantly) on every process. More... | |
void | ComputeCtAB_vec (const HypreParMatrix &A, const HypreParVector &B, const CAROM::Matrix &C, CAROM::Vector &CtAB_vec) |
This function computes a reduced, non-distributed vector C^t AB stored identically (redundantly) on every process. More... | |
void | verify_within_portion (const mfem::Vector &bb_min, const mfem::Vector &bb_max, const mfem::Vector &t, const double limit) |
Helper function to ensure that t is within a given percentage of the domain relative to the center of the mesh. Performs the check for each dimension of the mesh (works if mesh is 2D or 3D). More... | |
double | map_to_ref_mesh (const double &bb_min, const double &bb_max, const double &fraction) |
Maps a value from [-1, 1] to the corresponding mesh domain [bb_min, bb_max]. More... | |
double | map_from_ref_mesh (const double &bb_min, const double &bb_max, const double &value) |
Maps a value within mesh domain [bb_min, bb_max] to the corresponding value between [-1, 1]. More... | |
bool | fileExists (const std::string &name) |
int | split_dimension (const int dim, const MPI_Comm &comm=MPI_COMM_WORLD) |
Distribute the global size dim into MPI processes as equally as possible. More... | |
int | get_global_offsets (const int local_dim, std::vector< int > &offsets, const MPI_Comm &comm=MPI_COMM_WORLD) |
Save integer offsets for each MPI rank under MPI communicator comm, where each rank as the size of local_dim. the sum of local_dim over all MPI processes is returned as the total dimension. More... | |
bool | is_same (int x, const MPI_Comm &comm=MPI_COMM_WORLD) |
Check if an integer is equal over all MPI processes. More... | |
Copyright (c) 2013-2024, Lawrence Livermore National Security, LLC and other libROM project developers. See the top-level COPYRIGHT file for details.
SPDX-License-Identifier: (Apache-2.0 OR MIT)
Author: Milos Stojanovic Stojke (milsto) Adapted by Kevin Huynh from DifferentialEvolution.h from milsto's Github repo, https://github.com/milsto/differential-evolution, permission granted by milsto/differential-evolution's MIT license.
SPDX-License-Identifier: (MIT)
void CAROM::ComputeCtAB | ( | const HypreParMatrix & | A, |
const CAROM::Matrix & | B, | ||
const CAROM::Matrix & | C, | ||
CAROM::Matrix & | CtAB | ||
) |
This function computes a reduced, non-distributed matrix C^t AB stored identically (redundantly) on every process.
[in] | A | The distributed HypreParMatrix (an MFEM class) A in C^t AB. |
[in] | B | The distributed Matrix B in C^t AB. |
[in] | C | The distributed Matrix C in C^t AB. |
[out] | CtAB | The non-distributed Matrix C^t AB. |
Definition at line 15 of file Utilities.cpp.
void CAROM::ComputeCtAB_vec | ( | const HypreParMatrix & | A, |
const HypreParVector & | B, | ||
const CAROM::Matrix & | C, | ||
CAROM::Vector & | CtAB_vec | ||
) |
This function computes a reduced, non-distributed vector C^t AB stored identically (redundantly) on every process.
[in] | A | The distributed HypreParMatrix (an MFEM class) A in C^t AB. |
[in] | B | The distributed HypreParVector B in C^t AB. |
[in] | C | The distributed Matrix C in C^t AB. |
[out] | CtAB | The non-distributed Vector C^t AB. |
Definition at line 47 of file Utilities.cpp.
double CAROM::convertClosestRBFToEpsilon | ( | std::vector< Vector * > | parameter_points, |
std::string | rbf, | ||
double | closest_rbf_val | ||
) |
Convert closest RBF value to an epsilon value.
[in] | parameter_points | The parameter points. |
[in] | rbf | Which RBF to compute. |
[in] | closest_rbf_val | The RBF parameter determines the width of influence. Set the RBF value of the nearest two parameter points to a value between 0.0 to 1.0 |
Definition at line 176 of file Interpolator.cpp.
struct GreedyErrorIndicatorPoint CAROM::createGreedyErrorIndicatorPoint | ( | Vector * | point, |
std::shared_ptr< Vector > & | localROM | ||
) |
Create a greedy error indicator point.
[in] | point | The error indicator point. |
[in] | localROM | The point of the nearest ROM. |
Definition at line 1 of file GreedySampler.cpp.
struct GreedyErrorIndicatorPoint CAROM::createGreedyErrorIndicatorPoint | ( | Vector * | point, |
Vector * | localROM | ||
) |
Create a greedy error indicator point.
[in] | point | The error indicator point. |
[in] | localROM | The point of the nearest ROM. |
Definition at line 1 of file GreedySampler.cpp.
void CAROM::DEIM | ( | const Matrix * | f_basis, |
int | num_f_basis_vectors_used, | ||
std::vector< int > & | f_sampled_row, | ||
std::vector< int > & | f_sampled_rows_per_proc, | ||
Matrix & | f_basis_sampled_inv, | ||
int | myid, | ||
int | num_procs | ||
) |
Computes the DEIM algorithm on the given basis.
Implemented from Saifon Chaturantabut and Danny C. Sorensen "Nonlinear Model Reduction via Discrete Empirical Interpolation", SIAM J. Sci. Comput., 32(5), 2737–2764. (28 pages)
[in] | f_basis | The basis vectors for the RHS. |
[in] | num_f_basis_vectors_used | The number of basis vectors in f_basis to use in the algorithm. |
[out] | f_sampled_row | The local row ids of each sampled row. This will contain the sampled rows from all processors. Use f_sampled_rows_per_proc to index into this array for the sampled rows from a specific processor. |
[out] | f_sampled_rows_per_proc | The number of sampled rows for each processor. |
[out] | f_basis_sampled_inv | The inverse of the sampled basis of the RHS. |
[in] | myid | The rank of this process. |
[in] | num_procs | The total number of processes. |
Factory function to make a diagonal matrix with nonzero entries as in its Vector argument. The rows of this diagonal matrix are distributed in the same fashion as its Vector argument.
[in] | v | Vector of elements that will be on the diagonal of the output matrix. |
Get the number of rows on each process; these process row counts must be summed to get the number of columns on each process.
Compute the row starting index and total number of rows over all processes – which is also the total number of columns. Also compute the starting row index on each process.
Assume that Matrix rows are contiguous sets of row indices, e.g., if a four row matrix is partititioned over two processes, then process zero on the MPI_Comm owns rows 0 and 1, and process one on the MPI_Comm owns rows 2 and 3.
Create the diagonal matrix and assign process local entries in v to process local entries in the diagonal matrix.
Off-diagonal matrix entries are zero; diagonal matrix entries come from the input Vector.
Definition at line 2073 of file Matrix.cpp.
int CAROM::get_global_offsets | ( | const int | local_dim, |
std::vector< int > & | offsets, | ||
const MPI_Comm & | comm = MPI_COMM_WORLD |
||
) |
Save integer offsets for each MPI rank under MPI communicator comm, where each rank as the size of local_dim. the sum of local_dim over all MPI processes is returned as the total dimension.
[in] | local_dim | Input local dimension specified for each MPI rank. |
[out] | offsets | Resulting global integer offsets split by local_dim. |
[in] | comm | MPI communicator. default value MPI_COMM_WORLD. |
[out] | dim | (Returned value) Global dimension as the sum of all local_dim. |
Definition at line 41 of file mpi_utils.cpp.
double CAROM::getNearestPoint | ( | std::vector< double > | paramPoints, |
double | point | ||
) |
Given a double, find the nearest point in a domain.
[in] | paramPoints | The domain to search. |
[in] | point | The specified point. |
Definition at line 74 of file GreedySampler.cpp.
Given a vector, find the nearest point in a domain.
[in] | paramPoints | The domain to search. |
[in] | point | The specified point. |
Definition at line 44 of file GreedySampler.cpp.
int CAROM::getNearestPointIndex | ( | std::vector< double > | paramPoints, |
double | point | ||
) |
Given a double, find the nearest point in a domain.
[in] | paramPoints | The domain to search. |
[in] | point | The specified point. |
Definition at line 82 of file GreedySampler.cpp.
Given a vector, find the nearest point in a domain.
[in] | paramPoints | The domain to search. |
[in] | point | The specified point. |
Definition at line 52 of file GreedySampler.cpp.
void CAROM::getParametricDMD | ( | T *& | parametric_dmd, |
std::vector< Vector * > & | parameter_points, | ||
std::vector< std::string > & | dmd_paths, | ||
Vector * | desired_point, | ||
std::string | rbf = "G" , |
||
std::string | interp_method = "LS" , |
||
double | closest_rbf_val = 0.9 , |
||
bool | reorthogonalize_W = false |
||
) |
Constructor.
[in] | parameter_points | The parameter points. |
[in] | dmd_paths | The paths to the saved DMD objects associated with each parameter point. |
[in] | desired_point | The desired point at which to create a parametric DMD. |
[in] | rbf | The RBF type ("G" == gaussian, "IQ" == inverse quadratic, "IMQ" == inverse multiquadric) |
[in] | interp_method | The interpolation method type ("LS" == linear solve, "IDW" == inverse distance weighting, "LP" == lagrangian polynomials) |
[in] | closest_rbf_val | The RBF parameter determines the width of influence. Set the RBF value of the nearest two parameter points to a value between 0.0 to 1.0 |
[in] | reorthogonalize_W | Whether to reorthogonalize the interpolated W (basis) matrix. |
Definition at line 141 of file ParametricDMD.h.
void CAROM::getParametricDMD | ( | T *& | parametric_dmd, |
std::vector< Vector * > & | parameter_points, | ||
std::vector< T * > & | dmds, | ||
Vector * | desired_point, | ||
std::string | rbf = "G" , |
||
std::string | interp_method = "LS" , |
||
double | closest_rbf_val = 0.9 , |
||
bool | reorthogonalize_W = false |
||
) |
Constructor.
[in] | parametric_dmd | The interpolant DMD model at the desired point. |
[in] | parameter_points | The training parameter points. |
[in] | dmds | The DMD objects associated with each training parameter point. |
[in] | desired_point | The desired point at which to create a parametric DMD. |
[in] | rbf | The RBF type ("G" == gaussian, "IQ" == inverse quadratic, "IMQ" == inverse multiquadric) |
[in] | interp_method | The interpolation method type ("LS" == linear solve, "IDW" == inverse distance weighting, "LP" == lagrangian polynomials) |
[in] | closest_rbf_val | The RBF parameter determines the width of influence. Set the RBF value of the nearest two parameter points to a value between 0.0 to 1.0 |
[in] | reorthogonalize_W | Whether to reorthogonalize the interpolated W (basis) matrix. |
Definition at line 52 of file ParametricDMD.h.
void CAROM::getParametricDMDc | ( | T *& | parametric_dmdc, |
std::vector< Vector * > & | parameter_points, | ||
std::vector< std::string > & | dmdc_paths, | ||
std::vector< Matrix * > | controls, | ||
Matrix *& | controls_interpolated, | ||
Vector * | desired_point, | ||
std::string | rbf = "G" , |
||
std::string | interp_method = "LS" , |
||
double | closest_rbf_val = 0.9 , |
||
bool | reorthogonalize_W = false |
||
) |
Constructor.
[in] | parameter_points | The parameter points. |
[in] | dmdc_paths | The paths to the saved DMD objects associated with each parameter point. |
[in] | desired_point | The desired point at which to create a parametric DMD. |
[in] | controls | The matrices of controls from previous runs which we use to interpolate. |
[in] | controls_interpolated | The interpolated controls. |
[in] | rbf | The RBF type ("G" == gaussian, "IQ" == inverse quadratic, "IMQ" == inverse multiquadric) |
[in] | interp_method | The interpolation method type ("LS" == linear solve, "IDW" == inverse distance weighting, "LP" == lagrangian polynomials) |
[in] | closest_rbf_val | The RBF parameter determines the width of influence. Set the RBF value of the nearest two parameter points to a value between 0.0 to 1.0 |
[in] | reorthogonalize_W | Whether to reorthogonalize the interpolated W (basis) matrix. |
Definition at line 164 of file ParametricDMDc.h.
void CAROM::getParametricDMDc | ( | T *& | parametric_dmdc, |
std::vector< Vector * > & | parameter_points, | ||
std::vector< T * > & | dmdcs, | ||
std::vector< Matrix * > | controls, | ||
Matrix *& | controls_interpolated, | ||
Vector * | desired_point, | ||
std::string | rbf = "G" , |
||
std::string | interp_method = "LS" , |
||
double | closest_rbf_val = 0.9 , |
||
bool | reorthogonalize_W = false |
||
) |
Constructor.
[in] | parametric_dmdc | The interpolant DMDc model at the desired point. |
[in] | parameter_points | The training parameter points. |
[in] | dmdcs | The DMDc objects associated with each training parameter point. |
[in] | controls | The matrices of controls from previous runs which we use to interpolate. |
[in] | controls_interpolated | The interpolated controls. |
[in] | desired_point | The desired point at which to create a parametric DMDc. |
[in] | rbf | The RBF type ("G" == gaussian, "IQ" == inverse quadratic, "IMQ" == inverse multiquadric) |
[in] | interp_method | The interpolation method type ("LS" == linear solve, "IDW" == inverse distance weighting, "LP" == lagrangian polynomials) |
[in] | closest_rbf_val | The RBF parameter determines the width of influence. Set the RBF value of the nearest two parameter points to a value between 0.0 to 1.0 |
[in] | reorthogonalize_W | Whether to reorthogonalize the interpolated W (basis) matrix. |
Definition at line 56 of file ParametricDMDc.h.
void CAROM::GNAT | ( | const Matrix * | f_basis, |
const int | num_f_basis_vectors_used, | ||
std::vector< int > & | f_sampled_row, | ||
std::vector< int > & | f_sampled_rows_per_proc, | ||
Matrix & | f_basis_sampled_inv, | ||
const int | myid, | ||
const int | num_procs, | ||
const int | num_samples_req = -1 , |
||
std::vector< int > * | init_samples = NULL |
||
) |
Computes the GNAT algorithm on the given basis.
Implemented from Kevin Carlberg, Charbel Farhat, Julien Cortial, and David Amsallem "The GNAT method for nonlinear model reduction: Effective implementation and application to computational fluid dynamics and turbulent flows", Journal of Computational Physics Volume 242, 1 June 2013, Pages 623-647.
[in] | f_basis | The basis vectors for the RHS. |
[in] | num_f_basis_vectors_used | The number of basis vectors in f_basis to use in the algorithm. |
[out] | f_sampled_row | The local row ids of each sampled row. This will contain the sampled rows from all processors. Use f_sampled_rows_per_proc to index into this array for the sampled rows from a specific processor. |
[out] | f_sampled_rows_per_proc | The number of sampled rows for each processor. |
[out] | f_basis_sampled_inv | The inverse of the sampled basis of the RHS. |
[in] | myid | The rank of this process. |
[in] | num_procs | The total number of processes. |
[in] | num_samples_req | The minimum number of samples required. |
[in] | init_samples | Samples to initialize the GNAT algorithm. |
Factory function to make an identity matrix with rows distributed in the same fashion as its Vector argument. This function is provided for convenience due to the ubiquity of identity matrices (and operators) in mathematics.
[in] | v | Vector of elements that will be on the diagonal of the output matrix. |
Definition at line 2149 of file Matrix.cpp.
bool CAROM::is_same | ( | int | x, |
const MPI_Comm & | comm = MPI_COMM_WORLD |
||
) |
Check if an integer is equal over all MPI processes.
[in] | x | Input integer value to test equality. |
[in] | comm | MPI communicator. default value MPI_COMM_WORLD. |
Definition at line 72 of file mpi_utils.cpp.
double CAROM::map_from_ref_mesh | ( | const double & | bb_min, |
const double & | bb_max, | ||
const double & | value | ||
) |
Maps a value within mesh domain [bb_min, bb_max] to the corresponding value between [-1, 1].
bb_min | Minimum value of domain |
bb_max | Maximum value of domain |
value | Value between [bb_min, bb_max] to map |
value
mapped to [-1, 1] Definition at line 96 of file Utilities.cpp.
double CAROM::map_to_ref_mesh | ( | const double & | bb_min, |
const double & | bb_max, | ||
const double & | fraction | ||
) |
Maps a value from [-1, 1] to the corresponding mesh domain [bb_min, bb_max].
bb_min | Minimum value of domain |
bb_max | Maximum value of domain |
fraction | Value between [-1, 1] to map |
fraction
mapped to [bb_min, bb_max] Definition at line 88 of file Utilities.cpp.
struct ComplexEigenPair CAROM::NonSymmetricRightEigenSolve | ( | Matrix * | A | ) |
Computes the eigenvectors/eigenvalues of an NxN real nonsymmetric matrix.
[in] | A | The NxN real nonsymmetric matrix to be eigendecomposed. |
Definition at line 2149 of file Matrix.cpp.
Compute the RBF between two points.
[in] | rbf | Which RBF to compute. |
[in] | epsilon | The RBF parameter that determines the width of influence. |
[in] | point1 | The first point. |
[in] | point2 | The second point. |
Definition at line 149 of file Interpolator.cpp.
std::vector< double > CAROM::obtainRBFToTrainingPoints | ( | std::vector< Vector * > | parameter_points, |
std::string | interp_method, | ||
std::string | rbf, | ||
double | epsilon, | ||
Vector * | point | ||
) |
Compute the RBF from the parameter points with the unsampled parameter point.
[in] | parameter_points | The parameter points. |
[in] | interp_method | The interpolation method type ("LS" == linear solve, "IDW" == inverse distance weighting, "LP" == lagrangian polynomials) |
[in] | rbf | Which RBF to compute. |
[in] | epsilon | The RBF parameter that determines the width of influence. |
[in] | point | The unsampled parameter point. |
Definition at line 70 of file Interpolator.cpp.
std::vector< Matrix * > CAROM::obtainRotationMatrices | ( | std::vector< Vector * > | parameter_points, |
std::vector< Matrix * > | bases, | ||
int | ref_point | ||
) |
Obtain the rotation matrices for all the parameter points using the basis of the reference point.
[in] | parameter_points | The parameter points. |
[in] | bases | The spatial basis associated with each parameter point. |
[in] | ref_point | The index within the vector of parameter points to the reference point |
Definition at line 219 of file Interpolator.cpp.
Computes the outer product of two Vectors, v and w.
[in] | v | The first vector in the outer product |
[in] | w | The second vector in the outer product |
Definition at line 1973 of file Matrix.cpp.
void CAROM::QDEIM | ( | const Matrix * | f_basis, |
int | num_f_basis_vectors_used, | ||
std::vector< int > & | f_sampled_row, | ||
std::vector< int > & | f_sampled_rows_per_proc, | ||
Matrix & | f_basis_sampled_inv, | ||
const int | myid, | ||
const int | num_procs, | ||
const int | num_samples_req | ||
) |
Computes the QDEIM algorithm on the given basis.
Implemented from Zlatko Drmac, Serkan Gugercin, "A new selection operator for the discrete empirical interpolation method – Improved a priori error bound and extensions", SIAM Journal on Scientific Computing, Volume 38, Issue 2, pages A631-A648, 2016.
Oversampling uses GappyPOD+E from Peherstorfer, Drmac, Gugercin, "Stability of discrete empirical interpolation and gappy proper orthogonal decomposition with randomized and deterministic sampling points", preprint May 20, 2020.
[in] | f_basis | The basis vectors for the RHS. |
[in] | num_f_basis_vectors_used | The number of basis vectors in f_basis to use in the algorithm. |
[out] | f_sampled_row | The local row ids of each sampled row. This will contain the sampled rows from all processors. Use f_sampled_rows_per_proc to index into this array for the sampled rows from a specific processor. |
[out] | f_sampled_rows_per_proc | The number of sampled rows for each processor. |
[out] | f_basis_sampled_inv | The inverse of the sampled basis of the RHS. |
[in] | myid | The rank of this process. |
[in] | num_procs | The total number of processes. |
[in] | num_samples_req | The minimum number of samples required. |
double CAROM::rbfWeightedSum | ( | std::vector< double > & | rbf | ) |
Compute the sum of the RBF weights.
[in] | rbf | The vector holding the rbfs of the training points. |
Definition at line 139 of file Interpolator.cpp.
void CAROM::S_OPT | ( | const Matrix * | f_basis, |
int | num_f_basis_vectors_used, | ||
std::vector< int > & | f_sampled_row, | ||
std::vector< int > & | f_sampled_rows_per_proc, | ||
Matrix & | f_basis_sampled_inv, | ||
const int | myid, | ||
const int | num_procs, | ||
const int | num_samples_req = -1 , |
||
std::vector< int > * | init_samples = NULL , |
||
bool | qr_factorize = false |
||
) |
Computes the S_OPT algorithm on the given basis.
Implemented from Yeonjong Shin and Dongbin Xiu's "Nonadaptive Quasi-Optimal Points Selection for Least Squares Linear Regression", SIAM J. Sci. Comput., 38(1), A385-A411. (26 pages)
[in] | f_basis | The basis vectors for the RHS. |
[in] | num_f_basis_vectors_used | The number of basis vectors in f_basis to use in the algorithm. |
[out] | f_sampled_row | The local row ids of each sampled row. This will contain the sampled rows from all processors. Use f_sampled_rows_per_proc to index into this array for the sampled rows from a specific processor. |
[out] | f_sampled_rows_per_proc | The number of sampled rows for each processor. |
[out] | f_basis_sampled_inv | The inverse of the sampled basis of the RHS. |
[in] | myid | The rank of this process. |
[in] | num_procs | The total number of processes. |
[in] | num_samples_req | The minimum number of samples required. |
[in] | init_samples | Samples to initialize the S_OPT algorithm. |
[in] | qr_factorize | Whether to factorize the incoming matrix. If true and if the incoming matrix is a basis, the unnecessary computation of a QR factorization will be performed. |
struct SerialSVDDecomposition CAROM::SerialSVD | ( | Matrix * | A | ) |
Computes the SVD of a undistributed matrix in serial.
[in] | A | The MxN undistributed matrix to be eigendecomposed. |
Definition at line 2279 of file Matrix.cpp.
Computes the SVD of a undistributed matrix in serial.
[in] | A | The MxN undistributed matrix to be eigendecomposed. |
[in] | U | The matrix to write the left singular vectors into. |
[in] | S | The vector to write the singular values into. |
[in] | V | The matrix to write the right singular vectors into. |
Definition at line 2279 of file Matrix.cpp.
void CAROM::SpaceTimeSampling | ( | const Matrix * | s_basis, |
const Matrix * | t_basis, | ||
const int | num_f_basis_vectors_used, | ||
std::vector< int > & | t_samples, | ||
int * | f_sampled_row, | ||
int * | f_sampled_rows_per_proc, | ||
Matrix & | s_basis_sampled, | ||
const int | myid, | ||
const int | num_procs, | ||
const int | num_t_samples_req = -1 , |
||
const int | num_s_samples_req = -1 , |
||
const bool | excludeFinalTime = false |
||
) |
[in] | s_basis | The spatial basis vectors. |
[in] | t_basis | The temporal basis vectors. |
[in] | num_f_basis_vectors_used | The number of basis vectors in f_basis to use in the algorithm. |
[out] | t_samples | Array of temporal samples. |
[out] | f_sampled_row | The local row ids of each sampled row. This will contain the sampled rows from all processors. Use f_sampled_rows_per_proc to index into this array for the sampled rows from a specific processor. |
[out] | f_sampled_rows_per_proc | The number of sampled rows for each processor. |
[out] | f_basis_sampled_inv | The inverse of the sampled basis of the RHS. |
[in] | myid | The rank of this process. |
[in] | num_procs | The total number of processes. |
[in] | num_t_samples_req | The number of temporal samples to compute. |
[in] | num_s_samples_req | The number of spatial samples to compute. |
[in] | excludeFinalTime | Whether to exclude the final time index as a temporal sample. |
Definition at line 510 of file STSampling.cpp.
int CAROM::split_dimension | ( | const int | dim, |
const MPI_Comm & | comm = MPI_COMM_WORLD |
||
) |
Distribute the global size dim into MPI processes as equally as possible.
[in] | dim | Input global size. |
[in] | comm | MPI communicator. default value MPI_COMM_WORLD. |
[out] | local_dim | (Returned value) Local size assigned to the current MPI process. |
Definition at line 23 of file mpi_utils.cpp.
Computes the eigenvectors/eigenvalues of an NxN real symmetric matrix.
[in] | A | The NxN real symmetric matrix to be eigendecomposed. |
Definition at line 2149 of file Matrix.cpp.
void CAROM::verify_within_portion | ( | const mfem::Vector & | bb_min, |
const mfem::Vector & | bb_max, | ||
const mfem::Vector & | t, | ||
const double | limit | ||
) |
Helper function to ensure that t
is within a given percentage of the domain relative to the center of the mesh. Performs the check for each dimension of the mesh (works if mesh is 2D or 3D).
bb_min | Minimum corner of mesh bounding box. See mfem::Mesh::GetBoundingBox(). |
bb_max | Maximum corner of mesh bounding box. See mfem::Mesh::GetBoundingBox(). |
t | Point to check if within given limit percentage of mesh domain relative to mesh center. |
limit | Fractional percentage (from [0, 1]) of mesh domain to check bounds of t . |
Definition at line 65 of file Utilities.cpp.