[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
assemble_ECSW_training_data_base.h
1 #ifndef __ASSEMBLE_ECSW_BASE__
2 #define __ASSEMBLE_ECSW_BASE__
3 
4 #include <eigen/Eigen/Dense>
5 #include <Epetra_MpiComm.h>
6 #include <Epetra_SerialComm.h>
7 #include <Epetra_CrsMatrix.h>
8 #include <Epetra_Map.h>
9 #include <Epetra_Vector.h>
10 #include <EpetraExt_MatrixMatrix.h>
11 #include "dg/dg_base.hpp"
12 #include "pod_basis_base.h"
13 #include "parameters/all_parameters.h"
14 #include "multi_core_helper_functions.h"
15 
16 namespace PHiLiP {
17 namespace HyperReduction {
18 using Eigen::MatrixXd;
19 using Eigen::RowVectorXd;
20 using Eigen::VectorXd;
21 
24 
25 /*
26 Reference for the ECSW training data assembly, find details in section 4.1 and Equation (17):
27 "Mesh sampling and weighting for the hyperreduction of nonlinear Petrov–Galerkin reduced-order models with local reduced-order bases"
28 Sebastian Grimberg, Charbel Farhat, Radek Tezaur, Charbel Bou-Mosleh
29 International Journal for Numerical Methods in Engineering, 2020
30 https://onlinelibrary.wiley.com/doi/10.1002/nme.6603
31 NOTE: The above presents the approach for training with the residual data
32 */
33 
34 template <int dim, int nstate>
36 {
37 public:
40  const PHiLiP::Parameters::AllParameters *const parameters_input,
41  const dealii::ParameterHandler &parameter_handler_input,
42  std::shared_ptr<DGBase<dim,double>> &dg_input,
44  MatrixXd snapshot_parameters_input,
46  Epetra_MpiComm &Comm);
47 
49  virtual ~AssembleECSWBase () {};
50 
52 
54  const dealii::ParameterHandler &parameter_handler;
55 
57  std::shared_ptr<DGBase<dim,double>> dg;
58 
60  std::shared_ptr<ProperOrthogonalDecomposition::PODBase<dim>> pod;
61 
63  mutable MatrixXd snapshot_parameters;
64 
65  const MPI_Comm mpi_communicator;
66  const int mpi_rank;
67 
69 
71  dealii::ConditionalOStream pcout;
72 
75 
77  Epetra_MpiComm Comm_;
78 
80  std::shared_ptr<dealii::TrilinosWrappers::SparseMatrix> A_T;
81 
83  dealii::LinearAlgebra::ReadWriteVector<double> b;
84 
86  mutable std::vector<dealii::LinearAlgebra::distributed::Vector<double>> fom_locations;
87 
89  std::shared_ptr<Epetra_CrsMatrix> local_generate_test_basis(Epetra_CrsMatrix &system_matrix, const Epetra_CrsMatrix &pod_basis);
90 
92  Parameters::AllParameters reinit_params(const RowVectorXd& parameter) const;
93 
95  void update_snapshots(dealii::LinearAlgebra::distributed::Vector<double> fom_solution);
96 
98  void update_POD_snaps(std::shared_ptr<ProperOrthogonalDecomposition::PODBase<dim>> pod_update, MatrixXd snapshot_parameters_update);
99 
101  virtual void build_problem() = 0;
102 };
103 
104 } // HyperReduction namespace
105 } // PHiLiP namespace
106 
107 #endif
const dealii::ParameterHandler & parameter_handler
Parameter handler for storing the .prm file being ran.
std::vector< dealii::LinearAlgebra::distributed::Vector< double > > fom_locations
Vector of parameter-ROMTestLocation pairs.
MatrixXd snapshot_parameters
Matrix of snapshot parameters.
std::shared_ptr< dealii::TrilinosWrappers::SparseMatrix > A_T
Matrix for the NNLS Problem.
void update_snapshots(dealii::LinearAlgebra::distributed::Vector< double > fom_solution)
Update POD and Snapshot Parameters.
Files for the baseline physics.
Definition: ADTypes.hpp:10
virtual void build_problem()=0
Fill entries of A and b.
Main parameter class that contains the various other sub-parameter classes.
dealii::LinearAlgebra::ReadWriteVector< double > b
RHS Vector for the NNLS Problem.
std::shared_ptr< ProperOrthogonalDecomposition::PODBase< dim > > pod
POD.
Parameters::ODESolverParam::ODESolverEnum ode_solver_type
ODE Solve Type/ Projection Type (galerkin or petrov-galerkin)
Epetra_MpiComm Comm_
Epetra Communicator Object with MPI.
std::shared_ptr< Epetra_CrsMatrix > local_generate_test_basis(Epetra_CrsMatrix &system_matrix, const Epetra_CrsMatrix &pod_basis)
Generate Test Basis from the pod and snapshot info depending on the ode_solve_type (copied from the O...
AssembleECSWBase(const PHiLiP::Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input, std::shared_ptr< DGBase< dim, double >> &dg_input, std::shared_ptr< ProperOrthogonalDecomposition::PODBase< dim >> pod, MatrixXd snapshot_parameters_input, Parameters::ODESolverParam::ODESolverEnum ode_solver_type, Epetra_MpiComm &Comm)
Constructor.
const Parameters::AllParameters *const all_parameters
Pointer to all parameters.
std::shared_ptr< DGBase< dim, double > > dg
dg
void update_POD_snaps(std::shared_ptr< ProperOrthogonalDecomposition::PODBase< dim >> pod_update, MatrixXd snapshot_parameters_update)
Update POD and Snapshot Parameters.
dealii::ConditionalOStream pcout
ConditionalOStream.
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
Parameters::AllParameters reinit_params(const RowVectorXd &parameter) const
Reinitialize parameters.