33#include "simpleControl.H"
37#include "fvMeshSubset.H"
44#include <Eigen/SparseLU>
51 static fvScalarMatrix evaluate_expression(volScalarField& T, Eigen::MatrixXd mu)
53 volScalarField yPos = T.mesh().C().component(vector::Y).ref();
54 volScalarField xPos = T.mesh().C().component(vector::X).ref();
57 for (
auto i = 0; i < nu.size(); i++)
59 nu[i] = std::exp(- 2 * std::pow(xPos[i] - mu(0) - 1,
60 2) - 2 * std::pow(yPos[i] - mu(1) - 0.5, 2)) + 1;
63 nu.correctBoundaryConditions();
64 fvScalarMatrix TiEqn22
66 fvm::laplacian(nu, T,
"Gauss linear")
71 Eigen::MatrixXd onlineCoeffsA(Eigen::MatrixXd mu)
74 fvScalarMatrix Aof = evaluate_expression(fieldA(), mu);
75 Eigen::SparseMatrix<double> Mr;
85 theta(i) = Mr.coeffRef(ind_row, ind_col);
91 Eigen::MatrixXd onlineCoeffsB(Eigen::MatrixXd mu)
94 fvScalarMatrix Aof = evaluate_expression(fieldB(), mu);
95 Eigen::SparseMatrix<double> Mr;
102 theta(i) = br(ind_row);
108 PtrList<volScalarField> fieldsA;
109 autoPtr<volScalarField> fieldA;
110 autoPtr<volScalarField> fieldB;
111 PtrList<volScalarField> fieldsB;
114class DEIMLaplacian:
public laplacianProblem
117 explicit DEIMLaplacian(
int argc,
char* argv[])
119 laplacianProblem(argc, argv),
124 fvMesh& mesh =
_mesh();
136 NTmodes = readInt(
ITHACAdict->lookup(
"N_modes_T"));
137 NmodesDEIMA = readInt(
ITHACAdict->lookup(
"N_modes_DEIM_A"));
138 NmodesDEIMB = readInt(
ITHACAdict->lookup(
"N_modes_DEIM_B"));
147 PtrList<fvScalarMatrix> Mlist;
148 Eigen::MatrixXd ModesTEig;
149 std::vector<Eigen::MatrixXd> ReducedMatricesA;
150 std::vector<Eigen::MatrixXd> ReducedVectorsB;
159 void OfflineSolve(Eigen::MatrixXd par, word Folder)
167 for (
int i = 0; i < par.rows(); i++)
169 fvScalarMatrix Teqn = DEIMmatrice->evaluate_expression(T, par.row(i));
171 Mlist.append((Teqn).clone());
173 Tfield.append((T).clone());
178 void OnlineSolveFull(Eigen::MatrixXd par, word Folder)
180 auto t1 = std::chrono::high_resolution_clock::now();
181 auto t2 = std::chrono::high_resolution_clock::now();
182 auto time_span = std::chrono::duration_cast<std::chrono::duration<double >>
186 for (
int i = 0; i < par.rows(); i++)
188 fvScalarMatrix Teqn = DEIMmatrice->evaluate_expression(T, par.row(i));
189 t1 = std::chrono::high_resolution_clock::now();
191 t2 = std::chrono::high_resolution_clock::now();
192 time_span = std::chrono::duration_cast<std::chrono::duration<double >>
194 time_full += time_span.count();
196 Tfield.append((T).clone());
202 PODDEIM(NTmodes, NmodesDEIMA, NmodesDEIMB);
205 void PODDEIM(
int NmodesT,
int NmodesDEIMA,
int NmodesDEIMB)
207 DEIMmatrice =
new DEIM_function(Mlist, NmodesDEIMA, NmodesDEIMB,
"T_matrix");
208 fvMesh& mesh =
const_cast<fvMesh&
>(T.mesh());
210 DEIMmatrice->fieldA = autoPtr<volScalarField>(
new volScalarField(
211 DEIMmatrice->generateSubmeshMatrix(2, mesh, T)));
212 DEIMmatrice->fieldB = autoPtr<volScalarField>(
new volScalarField(
213 DEIMmatrice->generateSubmeshVector(2, mesh, T)));
216 ModesTEig.conservativeResize(ModesTEig.rows(), NmodesT);
217 ReducedMatricesA.resize(NmodesDEIMA);
218 ReducedVectorsB.resize(NmodesDEIMB);
220 for (
int i = 0; i < NmodesDEIMA; i++)
222 ReducedMatricesA[i] = ModesTEig.transpose() * DEIMmatrice->MatrixOnlineA[i] *
226 for (
int i = 0; i < NmodesDEIMB; i++)
228 ReducedVectorsB[i] = ModesTEig.transpose() * DEIMmatrice->MatrixOnlineB;
232 void OnlineSolve(Eigen::MatrixXd par_new, word Folder)
234 auto t1 = std::chrono::high_resolution_clock::now();
235 auto t2 = std::chrono::high_resolution_clock::now();
236 auto time_span = std::chrono::duration_cast<std::chrono::duration<double >>
240 for (
int i = 0; i < par_new.rows(); i++)
243 t1 = std::chrono::high_resolution_clock::now();
244 Eigen::MatrixXd thetaonA = DEIMmatrice->onlineCoeffsA(par_new.row(i));
245 Eigen::MatrixXd thetaonB = DEIMmatrice->onlineCoeffsB(par_new.row(i));
248 Eigen::VectorXd x = A.fullPivLu().solve(B);
249 Eigen::VectorXd full = ModesTEig * x;
250 t2 = std::chrono::high_resolution_clock::now();
251 time_span = std::chrono::duration_cast<std::chrono::duration<double >>
253 time_rom += time_span.count();
255 volScalarField Tred(
"Tred", T);
258 Tonline.append((Tred).clone());
263int main(
int argc,
char* argv[])
273 example.OfflineSolve(example.mu,
"Offline");
276 example.podex, 0, 0, 20);
282 example.OnlineSolve(par_new1,
"Online_red");
285 example_new.OnlineSolveFull(par_new1,
"Online_full");
287 Info << endl <<
"The FOM Solve took: " << example_new.time_full <<
289 Info << endl <<
"The ROM Solve took: " << example.time_rom <<
291 Info << endl <<
"The Speed-up is: " << example_new.time_full /
292 example.time_rom << endl << endl;
293 Eigen::MatrixXd error = ITHACAutilities::errorL2Abs(example_new.Tfield,
295 Info <<
"The mean L2 error is: " << error.mean() << endl;
Header file of the EigenFunctions class.
Header file of the Foam2Eigen class.
Header file of the ITHACAPOD class.
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the ITHACAutilities namespace.
autoPtr< IOList< label > > xyz_Arow
List< label > localMagicPointsAcol
autoPtr< IOList< label > > magicPointsAcol
List< label > localMagicPointsB
autoPtr< IOList< label > > xyz_B
DEIM(PtrList< T > &SnapShotsMatrix, label MaxModes, word FunctionName, word FieldName)
Construct DEIM for non-linear function.
List< label > localMagicPointsArow
autoPtr< IOList< label > > xyz_Acol
autoPtr< IOList< label > > magicPointsB
static GeometricField< scalar, PatchField, GeoMesh > Eigen2field(GeometricField< scalar, PatchField, GeoMesh > &field, Eigen::VectorXd &eigen_vector, bool correctBC=true)
Convert a vector in Eigen format into an OpenFOAM scalar GeometricField.
static Eigen::MatrixXd PtrList2Eigen(PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, label Nfields=-1)
Convert a PtrList of snapshots to Eigen matrix (only internal field).
static void fvMatrix2Eigen(fvMatrix< type_foam_matrix > foam_matrix, type_A &A, type_B &b)
Convert a FvMatrix OpenFOAM matrix (Linear System) into a Eigen Matrix A and a source vector b.
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
PtrList< volScalarField > Tonline
List of snapshots for the solution.
autoPtr< volScalarField > _T
Temperature field.
autoPtr< volScalarField > _nu
Diffusivity.
autoPtr< fvMesh > _mesh
Mesh.
volScalarModes Tmodes
List of POD modes.
autoPtr< volScalarField > _S
Source Term.
PtrList< volScalarField > Tfield
List of snapshots for the solution.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
IOdictionary * ITHACAdict
dictionary to store input output infos
Header file of the laplacianProblem class.
Eigen::SparseMatrix< T > MVproduct(List< Eigen::SparseMatrix< T > > &A, Eigen::DenseBase< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > > &C)
Sparse Matrix-Vector product between a list of sparse matrices and a vector of coefficients.
void getModes(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
Computes the bases or reads them for a field.
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
void read_fields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, word Name, fileName casename, int first_snap, int n_snap)
Function to read a list of fields from the name of the field and casename.
Eigen::MatrixXd rand(label rows, label cols, double min, double max)
Generates random matrix with random values in an interval.