Loading...
Searching...
No Matches
08DEIM.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12License
13 This file is part of ITHACA-FV
14 ITHACA-FV is free software: you can redistribute it and/or modify
15 it under the terms of the GNU Lesser General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18 ITHACA-FV is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU Lesser General Public License for more details.
22 You should have received a copy of the GNU Lesser General Public License
23 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
24Description
25 Example of the reconstruction of a non-linear function using the DEIM
26SourceFiles
27 08DEIM.C
28\*---------------------------------------------------------------------------*/
29
30#include "fvCFD.H"
31#include "fvOptions.H"
32#include "simpleControl.H"
33#include "simpleControl.H"
34#include "fvMeshSubset.H"
35#include "ITHACAutilities.H"
36#include "Foam2Eigen.H"
37#include "ITHACAstream.H"
38#include "ITHACAPOD.H"
40#include <chrono>
41class DEIM_function : public HyperReduction<PtrList<volScalarField>&>
42{
43 public:
45 static volScalarField evaluate_expression(volScalarField& S, Eigen::MatrixXd mu)
46 {
47 volScalarField yPos = S.mesh().C().component(vector::Y).ref();
48 volScalarField xPos = S.mesh().C().component(vector::X).ref();
49
50 for (auto i = 0; i < S.size(); i++)
51 {
52 S[i] = std::exp( - 2 * std::pow(xPos[i] - mu(0) - 1,
53 2) - 2 * std::pow(yPos[i] - mu(1) - 0.5, 2));
54 }
55
56 return S;
57 }
58 Eigen::VectorXd onlineCoeffs(Eigen::MatrixXd mu)
59 {
60 theta.resize(nodePoints().size());
61 auto f = evaluate_expression(subField(), mu);
62
63 for (int i = 0; i < nodePoints().size(); i++)
64 {
65 // double on_coeff = f[localMagicPoints[i]];
66 theta(i) = f[localNodePoints[i]];
67 }
68
69 return theta;
70 }
71
72 Eigen::VectorXd theta;
73 PtrList<volScalarField> fields;
74 autoPtr<volScalarField> subField;
75};
76
77int main(int argc, char* argv[])
78{
79 // Read parameters from ITHACAdict file
80#include "setRootCase.H"
81 Foam::Time runTime(Foam::Time::controlDictName, args);
82 Foam::fvMesh mesh
83 (
84 Foam::IOobject
85 (
86 Foam::fvMesh::defaultRegion,
87 runTime.timeName(),
88 runTime,
89 Foam::IOobject::MUST_READ
90 )
91 );
93 int NDEIM = para->ITHACAdict->lookupOrDefault<int>("NDEIM", 15);
94 simpleControl simple(mesh);
95#include "createFields.H"
96 // List of volScalarField where the snapshots are stored
97 PtrList<volScalarField> Sp;
98 // Non linear function
99 volScalarField S
100 (
101 IOobject
102 (
103 "S",
104 runTime.timeName(),
105 mesh,
106 IOobject::NO_READ,
107 IOobject::AUTO_WRITE
108 ),
109 T.mesh(),
110 dimensionedScalar("zero", dimensionSet(0, 0, -1, 1, 0, 0, 0), 0)
111 );
112 // Parameters used to train the non-linear function
113 Eigen::MatrixXd pars = ITHACAutilities::rand(100, 2, -0.5, 0.5);
114
115 // To compare with the same parameters, save data and then load it
116 // cnpy::load(pars, "./random100.npy");
117 // cnpy::save(pars, "./random100.npy");
118
119 // Perform the offline phase
120 for (int i = 0; i < 100; i++)
121 {
123 Sp.append((S).clone());
124 ITHACAstream::exportSolution(S, name(i + 1), "./ITHACAoutput/Offline/");
125 }
126
127 Eigen::MatrixXd snapshotsModes;
128
129 // To use the DEIMmodes method from ITHACAPOD :
130 Eigen::VectorXd normalizingWeights = ITHACAutilities::getMassMatrixFV(Sp[0]).array();
131 PtrList<volScalarField> modes = ITHACAPOD::DEIMmodes(Sp, NDEIM, "GappyDEIM", S.name());
132 snapshotsModes = Foam2Eigen::PtrList2Eigen(modes);
133
134 // To use the SVD modes from the HyperReduction class :
135 // Eigen::VectorXd normalizingWeights;
136 // c.getModesSVD(c.snapshotsListTuple, snapshotsModes, normalizingWeights);
137
138 // Create DEIM object with given number of basis functions
139 Eigen::VectorXi initSeeds(0);
140 DEIM_function c(NDEIM, NDEIM, initSeeds, "DEIM", Sp);
141
142 // Compute the DEIM
143 c.offlineGappyDEIM(snapshotsModes, normalizingWeights);
144 // To access the same folder hierarchy as before :
145 // c.problemName = "Gaussian_function";
146 // c.offlineGappyDEIM(snapshotsModes, normalizingWeights, "ITHACAoutput/DEIM/"+c.problemName);
147
148 // Generate the submeshes with the depth of the layer
149 c.generateSubmesh(2, mesh);
150 c.subField = c.interpolateField<volScalarField>(Sp[0]);
151
152 // Define a new online parameter
153 Eigen::MatrixXd par_new(2, 1);
154 par_new(0, 0) = 0;
155 par_new(1, 0) = 0;
156 // Online evaluation of the non linear function
157 Eigen::VectorXd aprfield = c.MatrixOnline * c.onlineCoeffs(par_new);
158 // Transform to an OpenFOAM field and export
159 volScalarField S2("S_online", Foam2Eigen::Eigen2field(S, aprfield));
160 ITHACAstream::exportSolution(S2, name(1), "./ITHACAoutput/Online/");
161 // Evaluate the full order function and export it
163 ITHACAstream::exportSolution(S, name(1), "./ITHACAoutput/Online/");
164 // Compute the L2 error and print it
165 Info << ITHACAutilities::errorL2Rel(S2, S) << endl;
166 return 0;
167}
168
172
265
266
int main(int argc, char *argv[])
Definition 08DEIM.C:77
Header file of the Foam2Eigen class.
Header file of the ITHACAPOD class.
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the ITHACAutilities namespace.
static volScalarField evaluate_expression(volScalarField &S, Eigen::MatrixXd mu)
Definition 08DEIM.C:45
autoPtr< volScalarField > subField
Definition 08DEIM.C:74
Eigen::VectorXd theta
Definition 08DEIM.C:72
Eigen::VectorXd onlineCoeffs(Eigen::MatrixXd mu)
Definition 08DEIM.C:58
PtrList< volScalarField > fields
Definition 08DEIM.C:73
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
void generateSubmesh(label layers, const fvMesh &mesh)
Compute the submesh common to all fields in SnapshotsLists.
autoPtr< FieldType > interpolateField(const FieldType &field)
TODO.
void offlineGappyDEIM(Eigen::MatrixXd &snapshotsModes, Eigen::VectorXd &normalizingWeights)
Methods implemented: 'GappyDEIM' from "DEIM, Chaturantabut, Saifon, and Danny C. Sorensen....
Eigen::MatrixXd MatrixOnline
Online Matrix.
List< label > localNodePoints
Indices of the local node points in the subMesh.
HyperReduction(label n_modes, label n_nodes, Eigen::VectorXi initialSeeds, word problemName, SnapshotsLists &&...snapshotsLists)
Construct HyperReduction class, interpolation-based.
autoPtr< IOList< label > > nodePoints
Nodes in the case of the a nonlinear function.
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
IOdictionary * ITHACAdict
Dictionary for input objects from file.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
volScalarField & T
Definition createT.H:46
std::tuple< List< Eigen::SparseMatrix< double > >, List< Eigen::VectorXd > > DEIMmodes(List< Eigen::SparseMatrix< double > > &A, List< Eigen::VectorXd > &b, label nmodesA, label nmodesB, word MatrixName)
Get the DEIM modes for a generic non a parametrized matrix coming from a differential operator functi...
Definition ITHACAPOD.C:784
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
double errorL2Rel(GeometricField< T, fvPatchField, volMesh > &field1, GeometricField< T, fvPatchField, volMesh > &field2, List< label > *labels)
Computes the relative error between two geometric Fields in L2 norm.
Eigen::VectorXd getMassMatrixFV(GeometricField< Type, PatchField, GeoMesh > &snapshot)
Gets a vector containing the volumes of each cell of the mesh.
Eigen::MatrixXd rand(label rows, label cols, double min, double max)
Generates random matrix with random values in an interval.
simpleControl simple(mesh)
label i
Definition pEqn.H:46
volScalarField xPos
volScalarField yPos
volScalarField S(IOobject("S", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), T.mesh(), dimensionedScalar("zero", dimensionSet(0, 0, -1, 1, 0, 0, 0), 0))