Loading...
Searching...
No Matches
ReducedSteadyNS.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-------------------------------------------------------------------------------
12
13License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26 You should have received a copy of the GNU Lesser General Public License
27 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
28
29\*---------------------------------------------------------------------------*/
30
33
34#include "ReducedSteadyNS.H"
35
36// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
37
38// Constructor
43
45 :
46 problem(&FOMproblem)
47{
48 N_BC = problem->inletIndex.rows();
49 Nphi_u = problem->B_matrix.rows();
50 Nphi_p = problem->K_matrix.cols();
51
52 for (int k = 0; k < problem->liftfield.size(); k++)
53 {
54 Umodes.append((problem->liftfield[k]).clone());
55 }
56
57 for (int k = 0; k < problem->NUmodes; k++)
58 {
59 Umodes.append((problem->Umodes[k]).clone());
60 }
61
62 for (int k = 0; k < problem->NSUPmodes; k++)
63 {
64 Umodes.append((problem->supmodes[k]).clone());
65 }
66
68}
69
70int newton_steadyNS::operator()(const Eigen::VectorXd& x,
71 Eigen::VectorXd& fvec) const
72{
73 Eigen::VectorXd a_tmp(Nphi_u);
74 Eigen::VectorXd b_tmp(Nphi_p);
75 a_tmp = x.head(Nphi_u);
76 b_tmp = x.tail(Nphi_p);
77 // Convective term
78 Eigen::MatrixXd cc(1, 1);
79 // Mom Term
80 Eigen::VectorXd M1 = problem->B_matrix * a_tmp * nu;
81 // Gradient of pressure
82 Eigen::VectorXd M2 = problem->K_matrix * b_tmp;
83 // Pressure Term
84 Eigen::VectorXd M3 = problem->P_matrix * a_tmp;
85 // Penalty term
86 Eigen::MatrixXd penaltyU = Eigen::MatrixXd::Zero(Nphi_u, N_BC);
87
88 // Term for penalty method
89 if (problem->bcMethod == "penalty")
90 {
91 for (int l = 0; l < N_BC; l++)
92 {
93 penaltyU.col(l) = BC(l) * problem->bcVelVec[l] - problem->bcVelMat[l] *
94 a_tmp;
95 }
96 }
97
98 for (int i = 0; i < Nphi_u; i++)
99 {
100 cc = a_tmp.transpose() * Eigen::SliceFromTensor(problem->C_tensor, 0,
101 i) * a_tmp;
102 fvec(i) = M1(i) - cc(0, 0) - M2(i);
103
104 if (problem->bcMethod == "penalty")
105 {
106 fvec(i) += (penaltyU.row(i) * tauU)(0, 0);
107 }
108 }
109
110 for (int j = 0; j < Nphi_p; j++)
111 {
112 int k = j + Nphi_u;
113 fvec(k) = M3(j);
114 }
115
116 if (problem->bcMethod == "lift")
117 {
118 for (int j = 0; j < N_BC; j++)
119 {
120 fvec(j) = x(j) - BC(j);
121 }
122 }
123
124 return 0;
125}
126
127
128int newton_steadyNS::df(const Eigen::VectorXd& x,
129 Eigen::MatrixXd& fjac) const
130{
131 Eigen::NumericalDiff<newton_steadyNS> numDiff(*this);
132 numDiff.df(x, fjac);
133 return 0;
134}
135
136
137// * * * * * * * * * * * * * * * Solve Functions * * * * * * * * * * * * * //
138
139void reducedSteadyNS::solveOnline_PPE(Eigen::MatrixXd vel_now)
140{
141 Info << "This function is still not implemented for the stationary case" <<
142 endl;
143 exit(0);
144}
145
146void reducedSteadyNS::solveOnline_sup(Eigen::MatrixXd vel)
147{
148 if (problem->bcMethod == "lift")
149 {
151 }
152 else if (problem->bcMethod == "penalty")
153 {
154 vel_now = vel;
155 }
156 else
157 {
158 M_Assert(false,
159 "The BC method must be set to lift or penalty in ITHACAdict");
160 }
161
162 y.resize(Nphi_u + Nphi_p, 1);
163 y.setZero();
164
165 // Change initial condition for the lifting function
166 if (problem->bcMethod == "lift")
167 {
168 for (int j = 0; j < N_BC; j++)
169 {
170 y(j) = vel_now(j, 0);
171 }
172 }
173
177 Eigen::HybridNonLinearSolver<newton_steadyNS> hnls(newton_object);
178 newton_object.BC.resize(N_BC);
180
181 for (int j = 0; j < N_BC; j++)
182 {
183 newton_object.BC(j) = vel_now(j, 0);
184 }
185
187 hnls.solve(y);
188 Eigen::VectorXd res(y);
189 newton_object.operator()(y, res);
190 Info << "################## Online solve N° " << count_online_solve <<
191 " ##################" << endl;
192
193 if (Pstream::master())
194 {
195 std::cout << "Solving for the parameter: " << vel_now << std::endl;
196 }
197
198 if (res.norm() < 1e-5 && Pstream::master())
199 {
200 std::cout << green << "|F(x)| = " << res.norm() << " - Minimun reached in " <<
201 hnls.iter << " iterations " << def << std::endl << std::endl;
202 }
203 else if (Pstream::master())
204 {
205 std::cout << red << "|F(x)| = " << res.norm() << " - Minimun reached in " <<
206 hnls.iter << " iterations " << def << std::endl << std::endl;
207 }
208
210}
211
212
213// * * * * * * * * * * * * * * * Jacobian Evaluation * * * * * * * * * * * * * //
214
215void reducedSteadyNS::reconstruct_PPE(fileName folder, int printevery)
216{
217 mkDir(folder);
219 int counter = 0;
220 int nextwrite = 0;
221
222 for (int i = 0; i < online_solution.size(); i++)
223 {
224 if (counter == nextwrite)
225 {
226 volVectorField U_rec("U_rec", Umodes[0] * 0);
227
228 for (int j = 0; j < Nphi_u; j++)
229 {
230 U_rec += Umodes[j] * online_solution[i](j + 1, 0);
231 }
232
233 ITHACAstream::exportSolution(U_rec, name(online_solution[i](0, 0)), folder);
234 volScalarField P_rec("P_rec", problem->Pmodes[0] * 0);
235
236 for (int j = 0; j < Nphi_p; j++)
237 {
238 P_rec += problem->Pmodes[j] * online_solution[i](j + Nphi_u + 1, 0);
239 }
240
241 ITHACAstream::exportSolution(P_rec, name(online_solution[i](0, 0)), folder);
242 nextwrite += printevery;
243 }
244
245 counter++;
246 }
247}
248
249void reducedSteadyNS::reconstruct(bool exportFields, fileName folder,
250 int printevery)
251{
252 if (exportFields)
253 {
254 mkDir(folder);
256 }
257
258 int counter = 0;
259 int nextwrite = 0;
260 List <Eigen::MatrixXd> CoeffU;
261 List <Eigen::MatrixXd> CoeffP;
262 CoeffU.resize(0);
263 CoeffP.resize(0);
264
265 for (int i = 0; i < online_solution.size(); i++)
266 {
267 if (counter == nextwrite)
268 {
269 Eigen::MatrixXd currentUCoeff;
270 Eigen::MatrixXd currentPCoeff;
271 currentUCoeff = online_solution[i].block(1, 0, Nphi_u, 1);
272 currentPCoeff = online_solution[i].bottomRows(Nphi_p);
273 CoeffU.append(currentUCoeff);
274 CoeffP.append(currentPCoeff);
275 nextwrite += printevery;
276 }
277
278 counter++;
279 }
280
281 volVectorField uRec("uRec", Umodes[0]);
282 volScalarField pRec("pRec", problem->Pmodes[0]);
283 uRecFields = problem->L_U_SUPmodes.reconstruct(uRec, CoeffU, "uRec");
284 pRecFields = problem->Pmodes.reconstruct(pRec, CoeffP, "pRec");
285
286 if (exportFields)
287 {
289 "uRec");
291 "pRec");
292 }
293}
294
296{
297 double a;
298 Eigen::VectorXd sup(Nphi_u);
299 Eigen::VectorXd inf(Nphi_p);
300
301 for (int i = 0; i < Nphi_p; i++)
302 {
303 for (int j = 0; j < Nphi_u; j++)
304 {
305 sup(j) = fvc::domainIntegrate(fvc::div(Umodes[j]) * Pmodes[i]).value() /
307 }
308
309 inf(i) = sup.maxCoeff();
310 }
311
312 a = inf.minCoeff();
313 return a;
314}
315
316
318 fileName folder)
319{
320 mkDir(folder);
321 system("ln -s ../../constant " + folder + "/constant");
322 system("ln -s ../../0 " + folder + "/0");
323 system("ln -s ../../system " + folder + "/system");
324 int NUmodes = problem.NUmodes;
325 int NSUPmodes = problem.NSUPmodes;
326 int NPmodes = problem.NPmodes;
327 int liftfieldSize = problem.liftfield.size();
328 int totalSize = NUmodes + NSUPmodes + liftfieldSize;
329 Eigen::VectorXd cl(online_solution.size());
330 Eigen::VectorXd cd(online_solution.size());
331 //Read FORCESdict
332 IOdictionary FORCESdict
333 (
334 IOobject
335 (
336 "FORCESdict",
337 "./system",
338 Umodes[0].mesh(),
339 IOobject::MUST_READ,
340 IOobject::NO_WRITE
341 )
342 );
343 fTau.setZero(online_solution.size(), 3);
344 fN.setZero(online_solution.size(), 3);
345
346 for (int i = 0; i < online_solution.size(); i++)
347 {
348 for (int j = 0; j < totalSize; j++)
349 {
350 fTau.row(i) += problem.tauMatrix.row(j) * online_solution[i](j + 1, 0);
351 }
352
353 for (int j = 0; j < NPmodes; j++)
354 {
355 fN.row(i) += problem.nMatrix.row(j) * online_solution[i](j + Nphi_u + 1, 0);
356 }
357 }
358
359 // Export the matrices
360 if (para->exportPython)
361 {
362 ITHACAstream::exportMatrix(fTau, "fTau", "python", folder);
363 ITHACAstream::exportMatrix(fN, "fN", "python", folder);
364 }
365
366 if (para->exportMatlab)
367 {
368 ITHACAstream::exportMatrix(fTau, "fTau", "matlab", folder);
369 ITHACAstream::exportMatrix(fN, "fN", "matlab", folder);
370 }
371
372 if (para->exportTxt)
373 {
374 ITHACAstream::exportMatrix(fTau, "fTau", "eigen", folder);
375 ITHACAstream::exportMatrix(fN, "fN", "eigen", folder);
376 }
377}
378
379Eigen::MatrixXd reducedSteadyNS::setOnlineVelocity(Eigen::MatrixXd vel)
380{
381 assert(problem->inletIndex.rows() == vel.rows()
382 && "Imposed boundary conditions dimensions do not match given values matrix dimensions");
383 Eigen::MatrixXd vel_scal;
384 vel_scal.resize(vel.rows(), vel.cols());
385
386 for (int k = 0; k < problem->inletIndex.rows(); k++)
387 {
388 int p = problem->inletIndex(k, 0);
389 int l = problem->inletIndex(k, 1);
390 scalar area = gSum(problem->liftfield[0].mesh().magSf().boundaryField()[p]);
391 scalar u_lf = gSum(problem->liftfield[k].mesh().magSf().boundaryField()[p] *
392 problem->liftfield[k].boundaryField()[p]).component(l) / area;
393 vel_scal(k, 0) = vel(k, 0) / u_lf;
394 }
395
396 return vel_scal;
397}
Foam::fvMesh & mesh
Definition createMesh.H:47
#define M_Assert(Expr, Msg)
Header file of the reducedSteadyNS class.
Class to change color to the output stream.
Definition colormod.H:24
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
GeometricField< Type, PatchField, GeoMesh > reconstruct(GeometricField< Type, PatchField, GeoMesh > &inputField, Eigen::MatrixXd Coeff, word Name)
Function to reconstruct the solution starting from the coefficients, in this case the field is passed...
Definition Modes.C:275
Eigen::VectorXd y
Vector to store the solution during the Newton procedure.
ITHACAparameters * para
parameters to be read from the ITHACAdict file
Eigen::MatrixXd setOnlineVelocity(Eigen::MatrixXd vel)
Sets the online velocity.
void solveOnline_sup(Eigen::MatrixXd vel_now)
Method to perform an online solve using a supremizer stabilisation method.
Eigen::MatrixXd fTau
Reduced matrix for tangent forces.
PtrList< volScalarField > Pmodes
List of pointers to store the modes for pressure.
List< Eigen::MatrixXd > online_solution
List of Eigen matrices to store the online solution.
int Nphi_p
Number of pressure modes.
scalar nu
Reduced viscosity in case of parametrized viscosity.
Eigen::MatrixXd tauU
Penalty Factor.
void reconstructLiftAndDrag(steadyNS &problem, fileName folder)
Method to compute the reduced order forces.
void solveOnline_PPE(Eigen::MatrixXd vel_now)
Method to perform an online solve using a PPE stabilisation method.
PtrList< volScalarField > pRecFields
Reconstructed pressure fields list.
steadyNS * problem
Pointer to the FOM problem.
int count_online_solve
Counter to count the online solutions.
reducedSteadyNS()
Construct Null.
Eigen::MatrixXd vel_now
Online inlet velocity vector.
Eigen::MatrixXd fN
Reduced matrix for normal forces.
PtrList< volVectorField > Umodes
List of pointers to store the modes for velocity.
PtrList< volVectorField > uRecFields
Recontructed velocity fields list.
int N_BC
Number of parametrized boundary conditions.
void reconstruct_PPE(fileName folder="./ITHACAoutput/online_rec", int printevery=1)
Method to reconstruct a solution from an online solve with a PPE stabilisation technique.
double inf_sup_constant()
Method to evaluate the online inf-sup constant.
int Nphi_u
Number of velocity modes.
void reconstruct(bool exportFields=false, fileName folder="./ITHACAoutput/online_rec", int printevery=1)
Method to reconstruct the solutions from an online solve.
newton_steadyNS newton_object
Newton object used to solve the non linear problem.
Eigen::MatrixXi inletIndex
Matrix that contains informations about the inlet boundaries.
Implementation of a parametrized full order steady NS problem and preparation of the the reduced ma...
Definition steadyNS.H:70
List< Eigen::MatrixXd > bcVelVec
Boundary term for penalty method - vector.
Definition steadyNS.H:218
label NPmodes
Number of pressure modes used for the projection.
Definition steadyNS.H:143
Eigen::MatrixXd nMatrix
Pressure forces.
Definition steadyNS.H:215
Eigen::MatrixXd tauMatrix
Viscous forces.
Definition steadyNS.H:212
Eigen::Tensor< double, 3 > C_tensor
Diffusion term.
Definition steadyNS.H:167
volVectorModes supmodes
List of pointers used to form the supremizer modes.
Definition steadyNS.H:104
volVectorModes Umodes
List of pointers used to form the velocity modes.
Definition steadyNS.H:101
label NUmodes
Number of velocity modes used for the projection.
Definition steadyNS.H:140
PtrList< volVectorField > liftfield
List of pointers used to form the list of lifting functions.
Definition steadyNS.H:110
volVectorModes L_U_SUPmodes
List of pointers containing the total number of lift, supremizer and velocity modes.
Definition steadyNS.H:116
Eigen::MatrixXd B_matrix
Diffusion term.
Definition steadyNS.H:157
label NSUPmodes
Number of supremizer modes used for the projection.
Definition steadyNS.H:146
Eigen::MatrixXd K_matrix
Gradient of pressure matrix.
Definition steadyNS.H:163
List< Eigen::MatrixXd > bcVelMat
Boundary term for penalty method - matrix.
Definition steadyNS.H:221
Eigen::MatrixXd P_matrix
Div of velocity.
Definition steadyNS.H:170
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Definition steadyNS.H:98
word bcMethod
Boundary Method.
Definition steadyNS.H:317
@ FG_GREEN
Definition colormod.H:14
@ FG_DEFAULT
Definition colormod.H:16
@ FG_RED
Definition colormod.H:13
Matrix< VectorType, Dynamic, Dynamic > SliceFromTensor(Eigen::Tensor< VectorType, 3 > &tensor, label dim, label index1)
void exportFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &field, word folder, word fieldname)
Function to export a scalar of vector 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 exportMatrix(Eigen::Matrix< T, -1, dim > &matrix, word Name, word type, word folder)
Export the reduced matrices in numpy (type=python), matlab (type=matlab) and txt (type=eigen) format ...
double H1Seminorm(GeometricField< scalar, fvPatchField, volMesh > &field)
double L2Norm(GeometricField< scalar, fvPatchField, volMesh > &field)
Definition ITHACAerror.C:41
void createSymLink(word folder)
Creates symbolic links to 0, system and constant.
volScalarField & p
label i
Definition pEqn.H:46
Structure to implement a newton object for a stationary NS problem.
Eigen::MatrixXd tauU
Eigen::VectorXd BC
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const