Loading...
Searching...
No Matches
09DEIM_ROM.C

Introduction to tutorial 9

In this tutorial we implement a test where we use the Discrete Empirical Interpolation Method for a case where we have a non-linear dependency with respect to the input parameters.

The following image illustrates the computational domain which is the same as the previous example

The physical problem is given by a heat transfer problem which is described by the Poisson equation:

\[ \nabla \cdot (\nu \nabla T) = S \]

The parametric diffusivity is described by a parametric Gaussian function:

\[ \nu(\mathbf{x},\mathbf{\mu}) = e^{-2(x-\mu_x-1)^2 - 2(y-\mu_y-0.5)^2}, \]

The problem is then discretized as:

\[ A(\mu)T = b \]

In this case, even if the problem is linear, due to non-linearity with respect to the input parameter of the conductivity constant it is not possible to have an affine decomposition of the discretized differential operator.

We seek therefore an approximate affine expansion of the differential operator of this type:

\[ A(\mu) = \sum_{i = 1}^{N_D} \theta_i(\mu) A_i \]

using the Discrete Empirical Interpolation Method

The plain program

Here there's the plain code

/*---------------------------------------------------------------------------*\
██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
* In real Time Highly Advanced Computational Applications for Finite Volumes
* Copyright (C) 2017 by the ITHACA-FV authors
-------------------------------------------------------------------------------
License
This file is part of ITHACA-FV
ITHACA-FV is free software: you can redistribute it and/or modify
it under the terms of the GNU Lesser General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
ITHACA-FV is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License
along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
Description
Example of model reduction problem using the DEIM for a Heat Transfer Problem
SourceFiles
09DEIM_ROM.C
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "IOmanip.H"
#include "fvOptions.H"
#include "simpleControl.H"
#include "Foam2Eigen.H"
#include <chrono>
#include "fvMeshSubset.H"
#include "ITHACAstream.H"
#include "ITHACAPOD.H"
#include "EigenFunctions.H"
#include "DEIM.H"
#include <chrono>
#include <Eigen/SVD>
#include <Eigen/SparseLU>
class DEIM_function : public DEIM<fvScalarMatrix>
{
using DEIM::DEIM;
public:
static fvScalarMatrix evaluate_expression(volScalarField& T, Eigen::MatrixXd mu)
{
volScalarField yPos = T.mesh().C().component(vector::Y).ref();
volScalarField xPos = T.mesh().C().component(vector::X).ref();
volScalarField nu(T);
for (auto i = 0; i < nu.size(); i++)
{
nu[i] = std::exp( - 2 * std::pow(xPos[i] - mu(0) - 1,
2) - 2 * std::pow(yPos[i] - mu(1) - 0.5, 2)) + 1;
}
nu.correctBoundaryConditions();
fvScalarMatrix TiEqn22
(
fvm::laplacian(nu, T, "Gauss linear")
);
return TiEqn22;
}
Eigen::MatrixXd onlineCoeffsA(Eigen::MatrixXd mu)
{
Eigen::MatrixXd theta(magicPointsAcol().size(), 1);
fvScalarMatrix Aof = evaluate_expression(fieldA(), mu);
Eigen::SparseMatrix<double> Mr;
Eigen::VectorXd br;
for (int i = 0; i < magicPointsAcol().size(); i++)
{
int ind_row = localMagicPointsArow[i] + xyz_Arow()[i] *
fieldA().size();
int ind_col = localMagicPointsAcol[i] + xyz_Acol()[i] *
fieldA().size();
theta(i) = Mr.coeffRef(ind_row, ind_col);
}
return theta;
}
Eigen::MatrixXd onlineCoeffsB(Eigen::MatrixXd mu)
{
Eigen::MatrixXd theta(magicPointsB().size(), 1);
fvScalarMatrix Aof = evaluate_expression(fieldB(), mu);
Eigen::SparseMatrix<double> Mr;
Eigen::VectorXd br;
for (int i = 0; i < magicPointsB().size(); i++)
{
int ind_row = localMagicPointsB[i] + xyz_B()[i] * fieldB().size();
theta(i) = br(ind_row);
}
return theta;
}
PtrList<volScalarField> fieldsA;
autoPtr<volScalarField> fieldA;
autoPtr<volScalarField> fieldB;
PtrList<volScalarField> fieldsB;
};
{
public:
explicit DEIMLaplacian(int argc, char* argv[])
:
laplacianProblem(argc, argv),
nu(_nu()),
S(_S()),
T(_T())
{
fvMesh& mesh = _mesh();
ITHACAdict = new IOdictionary
(
IOobject
(
"ITHACAdict",
"./system",
IOobject::MUST_READ,
IOobject::NO_WRITE
)
);
NTmodes = readInt(ITHACAdict->lookup("N_modes_T"));
NmodesDEIMA = readInt(ITHACAdict->lookup("N_modes_DEIM_A"));
NmodesDEIMB = readInt(ITHACAdict->lookup("N_modes_DEIM_B"));
}
volScalarField& nu;
volScalarField& S;
volScalarField& T;
PtrList<fvScalarMatrix> Mlist;
Eigen::MatrixXd ModesTEig;
std::vector<Eigen::MatrixXd> ReducedMatricesA;
std::vector<Eigen::MatrixXd> ReducedVectorsB;
int NTmodes;
double time_full;
double time_rom;
void OfflineSolve(Eigen::MatrixXd par, word Folder)
{
if (offline)
{
ITHACAstream::read_fields(Tfield, T, "./ITHACAoutput/Offline/");
}
else
{
for (int i = 0; i < par.rows(); i++)
{
fvScalarMatrix Teqn = DEIMmatrice->evaluate_expression(T, par.row(i));
Teqn.solve();
Mlist.append((Teqn).clone());
ITHACAstream::exportSolution(T, name(i + 1), "./ITHACAoutput/" + Folder);
Tfield.append((T).clone());
}
}
};
void OnlineSolveFull(Eigen::MatrixXd par, word Folder)
{
auto t1 = std::chrono::high_resolution_clock::now();
auto t2 = std::chrono::high_resolution_clock::now();
auto time_span = std::chrono::duration_cast<std::chrono::duration<double>>
(t2 - t1);
time_full = 0;
for (int i = 0; i < par.rows(); i++)
{
fvScalarMatrix Teqn = DEIMmatrice->evaluate_expression(T, par.row(i));
t1 = std::chrono::high_resolution_clock::now();
Teqn.solve();
t2 = std::chrono::high_resolution_clock::now();
time_span = std::chrono::duration_cast<std::chrono::duration<double>>(t2 - t1);
time_full += time_span.count();
ITHACAstream::exportSolution(T, name(i + 1), "./ITHACAoutput/" + Folder);
Tfield.append((T).clone());
}
};
void PODDEIM()
{
}
void PODDEIM(int NmodesT, int NmodesDEIMA, int NmodesDEIMB)
{
fvMesh& mesh = const_cast<fvMesh&>(T.mesh());
// Differential Operator
DEIMmatrice->fieldA = autoPtr<volScalarField>(new volScalarField(
DEIMmatrice->fieldB = autoPtr<volScalarField>(new volScalarField(
// Source Terms
ModesTEig.conservativeResize(ModesTEig.rows(), NmodesT);
for (int i = 0; i < NmodesDEIMA; i++)
{
}
for (int i = 0; i < NmodesDEIMB; i++)
{
}
};
void OnlineSolve(Eigen::MatrixXd par_new, word Folder)
{
auto t1 = std::chrono::high_resolution_clock::now();
auto t2 = std::chrono::high_resolution_clock::now();
auto time_span = std::chrono::duration_cast<std::chrono::duration<double>>
(t2 - t1);
time_rom = 0;
for (int i = 0; i < par_new.rows(); i++)
{
// solve
t1 = std::chrono::high_resolution_clock::now();
Eigen::MatrixXd thetaonA = DEIMmatrice->onlineCoeffsA(par_new.row(i));
Eigen::MatrixXd thetaonB = DEIMmatrice->onlineCoeffsB(par_new.row(i));
Eigen::MatrixXd A = EigenFunctions::MVproduct(ReducedMatricesA, thetaonA);
Eigen::VectorXd B = EigenFunctions::MVproduct(ReducedVectorsB, thetaonB);
Eigen::VectorXd x = A.fullPivLu().solve(B);
Eigen::VectorXd full = ModesTEig * x;
t2 = std::chrono::high_resolution_clock::now();
time_span = std::chrono::duration_cast<std::chrono::duration<double>>(t2 - t1);
time_rom += time_span.count();
// Export
volScalarField Tred("Tred", T);
Tred = Foam2Eigen::Eigen2field(Tred, full);
ITHACAstream::exportSolution(Tred, name(i + 1), "./ITHACAoutput/" + Folder);
Tonline.append((Tred).clone());
}
}
};
int main(int argc, char* argv[])
{
// Construct the case
DEIMLaplacian example(argc, argv);
// Read some parameters from file
example._runTime());
// Create the offline parameters for the solve
example.mu = ITHACAutilities::rand(100, 2, -0.5, 0.5);
// Solve the offline problem to compute the snapshots for the projections
example.OfflineSolve(example.mu, "Offline");
// Compute the POD modes
ITHACAPOD::getModes(example.Tfield, example.Tmodes, example._T().name(),
example.podex, 0, 0, 20);
// Compute the offline part of the DEIM procedure
example.PODDEIM();
// Construct a new set of parameters
Eigen::MatrixXd par_new1 = ITHACAutilities::rand(100, 2, -0.5, 0.5);
// Solve the online problem with the new parameters
example.OnlineSolve(par_new1, "Online_red");
// Solve a new full problem with the new parameters (necessary to compute speed up and error)
DEIMLaplacian example_new(argc, argv);
example_new.OnlineSolveFull(par_new1, "Online_full");
// Output some infos
std::cout << std::endl << "The FOM Solve took: " << example_new.time_full <<
" seconds." << std::endl;
std::cout << std::endl << "The ROM Solve took: " << example.time_rom <<
" seconds." << std::endl;
std::cout << std::endl << "The Speed-up is: " << example_new.time_full /
example.time_rom << std::endl << std::endl;
Eigen::MatrixXd error = ITHACAutilities::errorL2Abs(example_new.Tfield,
example.Tonline);
std::cout << "The mean L2 error is: " << error.mean() << std::endl;
exit(0);
}
Header file of the EigenFunctions class.
Header file of the Foam2Eigen class.
Header file of the ITHACAPOD class.
Foam::fvMesh & mesh
Definition createMesh.H:47
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the ITHACAutilities namespace.
double time_rom
Definition 09DEIM_ROM.C:157
PtrList< fvScalarMatrix > Mlist
Definition 09DEIM_ROM.C:147
std::vector< Eigen::MatrixXd > ReducedMatricesA
Definition 09DEIM_ROM.C:149
std::vector< Eigen::MatrixXd > ReducedVectorsB
Definition 09DEIM_ROM.C:150
DEIMLaplacian(int argc, char *argv[])
Definition 09DEIM_ROM.C:117
double time_full
Definition 09DEIM_ROM.C:156
DEIM_function * DEIMmatrice
Definition 09DEIM_ROM.C:146
void OnlineSolve(Eigen::MatrixXd par_new, word Folder)
Definition 09DEIM_ROM.C:231
void OnlineSolveFull(Eigen::MatrixXd par, word Folder)
Definition 09DEIM_ROM.C:178
Eigen::MatrixXd ModesTEig
Definition 09DEIM_ROM.C:148
volScalarField & S
Definition 09DEIM_ROM.C:142
void OfflineSolve(Eigen::MatrixXd par, word Folder)
Definition 09DEIM_ROM.C:159
volScalarField & T
Definition 09DEIM_ROM.C:143
volScalarField & nu
Definition 09DEIM_ROM.C:141
Eigen::MatrixXd onlineCoeffsA(Eigen::MatrixXd mu)
Definition 09DEIM_ROM.C:71
PtrList< volScalarField > fieldsB
Definition 09DEIM_ROM.C:111
static volScalarField evaluate_expression(volScalarField &S, Eigen::MatrixXd mu)
Definition 08DEIM.C:45
Eigen::MatrixXd onlineCoeffsB(Eigen::MatrixXd mu)
Definition 09DEIM_ROM.C:91
Eigen::VectorXd theta
Definition 08DEIM.C:72
autoPtr< volScalarField > fieldB
Definition 09DEIM_ROM.C:110
autoPtr< volScalarField > fieldA
Definition 09DEIM_ROM.C:109
PtrList< volScalarField > fieldsA
Definition 09DEIM_ROM.C:108
Definition DEIM.H:41
autoPtr< IOList< label > > magicPointsAcol
Magic points in the case of the a matrix function (cols indices)
Definition DEIM.H:125
List< label > localMagicPointsAcol
Definition DEIM.H:161
S generateSubmeshMatrix(label layers, const fvMesh &mesh, S field, label secondTime=0)
Function to generate the submesh for the nonlinear matrix function case (RHS)
Definition DEIM.C:450
Eigen::MatrixXd MatrixOnlineB
Online Matrices.
Definition DEIM.H:110
List< label > localMagicPointsB
Definition DEIM.H:162
S generateSubmeshVector(label layers, const fvMesh &mesh, S field, label secondTime=0)
Function to generate the submesh for the nonlinear matrix function case (LHS)
Definition DEIM.C:543
autoPtr< IOList< label > > xyz_Acol
Definition DEIM.H:153
DEIM(PtrList< T > &SnapShotsMatrix, label MaxModes, word FunctionName, word FieldName)
Construct DEIM for non-linear function.
Definition DEIM.C:34
List< label > localMagicPointsArow
Definition DEIM.H:160
autoPtr< IOList< label > > xyz_Arow
Definition DEIM.H:152
List< Eigen::SparseMatrix< double > > MatrixOnlineA
Online Matrices.
Definition DEIM.H:109
autoPtr< IOList< label > > xyz_B
Definition DEIM.H:154
autoPtr< IOList< label > > magicPointsB
Magic points in the case of the a matrix function, right hand side.
Definition DEIM.H:127
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.
Definition Foam2Eigen.C:506
static Eigen::MatrixXd PtrList2Eigen(PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, label Nfields=-1)
Convert a PtrList of snapshots to Eigen matrix (only internal field)
Definition Foam2Eigen.C:638
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.
Class to implement a full order laplacian parametrized problem.
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
volScalarField & T
Definition createT.H:46
int main(int argc, char *argv[])
Header file of the laplacianProblem class.
volScalarField & A
dimensionedScalar & nu
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.
Definition ITHACAPOD.C:93
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.
double errorL2Abs(GeometricField< vector, fvPatchField, volMesh > &field1, GeometricField< vector, fvPatchField, volMesh > &field2, volScalarField &Volumes)
label i
Definition pEqn.H:46
volScalarField xPos
volScalarField yPos