Loading...
Searching...
No Matches
DEIM.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
31#include "DEIM.H"
32// Template function constructor
33template<typename T>
34DEIM<T>::DEIM (PtrList<T>& s, label MaxModes, word FunctionName, word FieldName)
35 :
36 SnapShotsMatrix(s),
37 MaxModes(MaxModes),
38 FunctionName(FunctionName)
39{
41 Folder = "ITHACAoutput/DEIM/" + FunctionName;
42 magicPoints = autoPtr<IOList<label>>
43 (
44 new IOList<label>
45 (
46 IOobject
47 (
48 "magicPoints",
49 para->runTime.time().constant(),
50 "../" + Folder,
51 para->mesh,
52 IOobject::READ_IF_PRESENT,
53 IOobject::NO_WRITE
54 )
55 )
56 );
57 xyz = autoPtr<IOList<label>>
58 (
59 new IOList<label>
60 (
61 IOobject
62 (
63 "xyz",
64 para->runTime.time().constant(),
65 "../" + Folder,
66 para->mesh,
67 IOobject::READ_IF_PRESENT,
68 IOobject::NO_WRITE
69 )
70 )
71 );
73 FieldName);
74
75 if (!(magicPoints().headerOk() && xyz().headerOk()))
76 {
77 Eigen::MatrixXd A;
78 Eigen::VectorXd b;
79 Eigen::VectorXd c;
80 Eigen::VectorXd r;
81 Eigen::VectorXd rho(1);
83 Ncells = modes[0].size();
84 label ind_max, c1, xyz_in;
85 double max = MatrixModes.cwiseAbs().col(0).maxCoeff(&ind_max, &c1);
86 check3DIndices(ind_max, xyz_in);
87 rho(0) = max;
88 magicPoints().append(ind_max);
89 xyz().append(xyz_in);
90 U = MatrixModes.col(0);
91 P.resize(MatrixModes.rows(), 1);
92 P.insert(ind_max, 0) = 1;
93
94 for (label i = 1; i < MaxModes; i++)
95 {
96 A = P.transpose() * U;
97 b = P.transpose() * MatrixModes.col(i);
98 c = A.fullPivLu().solve(b);
99 r = MatrixModes.col(i) - U * c;
100 max = r.cwiseAbs().maxCoeff(&ind_max, &c1);
101 P.conservativeResize(MatrixModes.rows(), i + 1);
102 P.insert(ind_max, i) = 1;
103 U.conservativeResize(MatrixModes.rows(), i + 1);
104 U.col(i) = MatrixModes.col(i);
105 rho.conservativeResize(i + 1);
106 rho(i) = max;
107 check3DIndices(ind_max, xyz_in);
108 magicPoints().append(ind_max);
109 xyz().append(xyz_in);
110 }
111
112 MatrixOnline = U * ((P.transpose() * U).fullPivLu().inverse());
113 mkDir(Folder);
114 cnpy::save(MatrixOnline, Folder + "/MatrixOnline.npy");
115 magicPoints().write();
116 xyz().write();
117 }
118 else
119 {
120 cnpy::load(MatrixOnline, Folder + "/MatrixOnline.npy");
121 }
122}
123
124// constructor for matrix DEIM
125template<typename T>
126DEIM<T>::DEIM (PtrList<T>& s, label MaxModesA, label MaxModesB, word MatrixName)
127 :
128 SnapShotsMatrix(s),
129 MaxModesA(MaxModesA),
130 MaxModesB(MaxModesB),
131 MatrixName(MatrixName),
132 runSubMesh(false),
133 runSubMeshA(false),
134 runSubMeshB(false)
135
136{
138 FolderM = "ITHACAoutput/DEIM/" + MatrixName;
139 magicPointsArow = autoPtr<IOList<label>>
140 (
141 new IOList<label>
142 (
143 IOobject
144 (
145 "magicPointsArow",
146 para->runTime.time().constant(),
147 "../" + FolderM,
148 para->mesh,
149 IOobject::READ_IF_PRESENT,
150 IOobject::NO_WRITE
151 )
152 )
153 );
154 magicPointsAcol = autoPtr<IOList<label>>
155 (
156 new IOList<label>
157 (
158 IOobject
159 (
160 "magicPointsAcol",
161 para->runTime.time().constant(),
162 "../" + FolderM,
163 para->mesh,
164 IOobject::READ_IF_PRESENT,
165 IOobject::NO_WRITE
166 )
167 )
168 );
169 magicPointsB = autoPtr<IOList<label>>
170 (
171 new IOList<label>
172 (
173 IOobject
174 (
175 "magicPointsB",
176 para->runTime.time().constant(),
177 "../" + FolderM,
178 para->mesh,
179 IOobject::READ_IF_PRESENT,
180 IOobject::NO_WRITE
181 )
182 )
183 );
184 xyz_Arow = autoPtr<IOList<label>>
185 (
186 new IOList<label>
187 (
188 IOobject
189 (
190 "xyz_Arow",
191 para->runTime.time().constant(),
192 "../" + FolderM,
193 para->mesh,
194 IOobject::READ_IF_PRESENT,
195 IOobject::NO_WRITE
196 )
197 )
198 );
199 xyz_Acol = autoPtr<IOList<label>>
200 (
201 new IOList<label>
202 (
203 IOobject
204 (
205 "xyz_Acol",
206 para->runTime.time().constant(),
207 "../" + FolderM,
208 para->mesh,
209 IOobject::READ_IF_PRESENT,
210 IOobject::NO_WRITE
211 )
212 )
213 );
214 xyz_B = autoPtr<IOList<label>>
216 new IOList<label>
217 (
218 IOobject
219 (
220 "xyz_B",
221 para->runTime.time().constant(),
222 "../" + FolderM,
223 para->mesh,
224 IOobject::READ_IF_PRESENT,
225 IOobject::NO_WRITE
226 )
227 )
228 );
229
230 if (!(magicPointsArow().headerOk() && magicPointsAcol().headerOk() &&
231 magicPointsB().headerOk() && xyz_Arow().headerOk() &&
232 xyz_Acol().headerOk() && xyz_B().headerOk()))
233 {
234 Eigen::MatrixXd AA;
235 Eigen::VectorXd bA;
236 Eigen::MatrixXd cA;
237 Eigen::SparseMatrix<double> rA;
238 Eigen::VectorXd rhoA(1);
240 MatrixName);
241 Ncells = getNcells(std::get<1>(Matrix_Modes)[0].rows());
242 label ind_rowA, ind_colA, xyz_rowA, xyz_colA;
243 ind_rowA = ind_colA = xyz_rowA = xyz_colA = 0;
244 double maxA = EigenFunctions::max(std::get<0>(Matrix_Modes)[0], ind_rowA,
245 ind_colA);
246 label ind_rowAOF = ind_rowA;
247 label ind_colAOF = ind_colA;
248 check3DIndices(ind_rowAOF, ind_colAOF, xyz_rowA, xyz_colA);
249 xyz_Arow().append(xyz_rowA);
250 xyz_Acol().append(xyz_colA);
251 rhoA(0) = maxA;
252 magicPointsArow().append(ind_rowAOF);
253 magicPointsAcol().append(ind_colAOF);
254 UA.append(std::get<0>(Matrix_Modes)[0]);
255 Eigen::SparseMatrix<double> Pnow(std::get<0>(Matrix_Modes)[0].rows(),
256 std::get<0>(Matrix_Modes)[0].cols());
257 Pnow.insert(ind_rowA, ind_colA) = 1;
258 PA.append(Pnow);
259
260 for (label i = 1; i < MaxModesA; i++)
261 {
263 bA = EigenFunctions::innerProduct(PA, std::get<0>(Matrix_Modes)[i]);
264 cA = AA.fullPivLu().solve(bA);
265 rA = std::get<0>(Matrix_Modes)[i] - EigenFunctions::MVproduct(UA, cA);
266 double maxA = EigenFunctions::max(rA, ind_rowA, ind_colA);
267 rhoA.conservativeResize(i + 1);
268 rhoA(i) = maxA;
269 label ind_rowAOF = ind_rowA;
270 label ind_colAOF = ind_colA;
271 check3DIndices(ind_rowAOF, ind_colAOF, xyz_rowA, xyz_colA);
272 xyz_Arow().append(xyz_rowA);
273 xyz_Acol().append(xyz_colA);
274 magicPointsArow().append(ind_rowAOF);
275 magicPointsAcol().append(ind_colAOF);
276 UA.append(std::get<0>(Matrix_Modes)[i]);
277 Eigen::SparseMatrix<double> Pnow(std::get<0>(Matrix_Modes)[0].rows(),
278 std::get<0>(Matrix_Modes)[0].cols());
279 Pnow.insert(ind_rowA, ind_colA) = 1;
280 PA.append(Pnow);
281 }
283 Eigen::MatrixXd Aaux = EigenFunctions::innerProduct(PA,
284 UA).fullPivLu().inverse();
286 Eigen::MatrixXd AB;
287 Eigen::VectorXd bB;
288 Eigen::VectorXd cB;
289 Eigen::VectorXd rB;
290 Eigen::VectorXd rhoB(1);
291 label ind_rowB, xyz_rowB, c1;
292 double maxB = std::get<1>(Matrix_Modes)[0].cwiseAbs().maxCoeff(&ind_rowB,
293 &c1);
294 label ind_rowBOF = ind_rowB;
295 check3DIndices(ind_rowBOF, xyz_rowB);
296 rhoB(0) = maxB;
297 xyz_B().append(xyz_rowB);
298 magicPointsB().append(ind_rowBOF);
299 UB = std::get<1>(Matrix_Modes)[0];
300 PB.resize(UB.rows(), 1);
301 PB.insert(ind_rowB, 0) = 1;
303 for (label i = 1; i < MaxModesB; i++)
304 {
305 AB = PB.transpose() * UB;
306 bB = PB.transpose() * std::get<1>(Matrix_Modes)[i];
307 cB = AB.fullPivLu().solve(bB);
308 rB = std::get<1>(Matrix_Modes)[i] - UB * cB;
309 maxB = rB.cwiseAbs().maxCoeff(&ind_rowB, &c1);
310 ind_rowBOF = ind_rowB;
311 check3DIndices(ind_rowBOF, xyz_rowB);
312 xyz_B().append(xyz_rowB);
313 PB.conservativeResize((std::get<1>(Matrix_Modes)[i]).size(), i + 1);
314 PB.insert(ind_rowB, i) = 1;
315 UB.conservativeResize((std::get<1>(Matrix_Modes)[i]).size(), i + 1);
316 UB.col(i) = std::get<1>(Matrix_Modes)[i];
317 rhoB.conservativeResize(i + 1);
318 rhoB(i) = maxB;
319 magicPointsB().append(ind_rowBOF);
320 }
321
322 if (MaxModesB == 1 && std::get<1>(Matrix_Modes)[0].norm() < 1e-8)
323 {
324 MatrixOnlineB = Eigen::MatrixXd::Zero(std::get<1>(Matrix_Modes)[0].rows(), 1);
325 }
326 else if (MaxModesB != 1)
327 {
328 MatrixOnlineB = UB * ((PB.transpose() * UB).fullPivLu().inverse());
329 }
330 else
331 {
332 Eigen::MatrixXd aux = PB.transpose() * UB;
333 MatrixOnlineB = UB * 1 / aux(0, 0);
334 }
335
336 mkDir(FolderM + "/lhs");
337 mkDir(FolderM + "/rhs");
338 ITHACAstream::save(MatrixOnlineA, FolderM, "/lhs/MatrixOnlineA");
339 cnpy::save(MatrixOnlineB, FolderM + "/rhs/MatrixOnlineB.npy");
340 magicPointsArow().write();
341 magicPointsAcol().write();
342 magicPointsB().write();
343 xyz_Arow().write();
344 xyz_Acol().write();
345 xyz_B().write();
346 }
347 else
348 {
349 ITHACAstream::load(MatrixOnlineA, FolderM, "/lhs/MatrixOnlineA");
350 cnpy::load(MatrixOnlineB, FolderM + "/rhs/MatrixOnlineB.npy");
351 }
352}
353
354
355template<typename T>
356template<typename S>
357S DEIM<T>::generateSubmesh(label layers, const fvMesh& mesh, S field,
358 label secondTime)
359{
361 totalMagicPoints = autoPtr<IOList<labelList>>
362 (
363 new IOList<labelList>
364 (
365 IOobject
366 (
367 "totalMagicPoints",
368 para->runTime.time().constant(),
369 "../" + Folder,
370 para->mesh,
371 IOobject::READ_IF_PRESENT,
372 IOobject::NO_WRITE
373 )
374 )
375 );
376 uniqueMagicPoints = autoPtr<IOList<label>>
377 (
378 new IOList<label>
379 (
380 IOobject
381 (
382 "uniqueMagicPoints",
383 para->runTime.time().constant(),
384 "../" + Folder,
385 para->mesh,
386 IOobject::READ_IF_PRESENT,
387 IOobject::NO_WRITE
388 )
389 )
390 );
391 volScalarField Indici
392 (
393 IOobject
394 (
395 FunctionName + "_indices",
396 mesh.time().timeName(),
397 mesh,
398 IOobject::NO_READ,
399 IOobject::NO_WRITE
400 ),
401 mesh,
402 dimensionedScalar(FunctionName + "_indices", dimensionSet(0, 0, 0, 0, 0, 0, 0),
403 0.0)
404 );
405 submesh = autoPtr<fvMeshSubset>(new fvMeshSubset(mesh));
406
407 if (!totalMagicPoints().headerOk())
408 {
409 List<label> indices;
410
411 for (label i = 0; i < magicPoints().size(); i++)
412 {
413 indices = ITHACAutilities::getIndices(mesh, magicPoints()[i], layers);
414 totalMagicPoints().append(indices);
415 }
416
417 uniqueMagicPoints() = ITHACAutilities::combineList(totalMagicPoints());
418 }
419
420#if OPENFOAM >= 1812
421 submesh->setCellSubset(uniqueMagicPoints());
422#else
423 submesh->setLargeCellSubset(uniqueMagicPoints());
424#endif
425 submesh->subMesh().fvSchemes::readOpt() = mesh.fvSchemes::readOpt();
426 submesh->subMesh().fvSolution::readOpt() = mesh.fvSolution::readOpt();
427 submesh->subMesh().fvSchemes::read();
428 submesh->subMesh().fvSolution::read();
429 std::cout.clear();
430 S f = submesh->interpolate(field).ref();
431 scalar zerodot25 = 0.25;
432 ITHACAutilities::assignIF(Indici, zerodot25,
433 uniqueMagicPoints().List<label>::clone()());
434 ITHACAutilities::assignONE(Indici, magicPoints());
435
436 if (!secondTime)
437 {
438 localMagicPoints = global2local(magicPoints(), submesh());
439 ITHACAstream::exportSolution(Indici, "1", "./ITHACAoutput/DEIM/" + FunctionName
440 );
441 }
442
443 totalMagicPoints().write();
444 uniqueMagicPoints().write();
445 return f;
446}
447
448template<typename T>
449template<typename S>
450S DEIM<T>::generateSubmeshMatrix(label layers, const fvMesh& mesh, S field,
451 label secondTime)
452{
454 totalMagicPointsA = autoPtr<IOList<labelList>>
455 (
456 new IOList<labelList>
457 (
458 IOobject
459 (
460 "totalMagicPointsA",
461 para->runTime.time().constant(),
462 "../" + FolderM,
463 para->mesh,
464 IOobject::READ_IF_PRESENT,
465 IOobject::NO_WRITE
466 )
467 )
468 );
469 uniqueMagicPointsA = autoPtr<IOList<label>>
470 (
471 new IOList<label>
472 (
473 IOobject
474 (
475 "uniqueMagicPointsA",
476 para->runTime.time().constant(),
477 "../" + FolderM,
478 para->mesh,
479 IOobject::READ_IF_PRESENT,
480 IOobject::NO_WRITE
481 )
482 )
483 );
484 List<label> indices;
485 volScalarField Indici
486 (
487 IOobject
488 (
489 MatrixName + "_A_indices",
490 mesh.time().timeName(),
491 mesh,
492 IOobject::NO_READ,
493 IOobject::NO_WRITE
494 ),
495 mesh,
496 dimensionedScalar(MatrixName + "_A_indices", dimensionSet(0, 0, 0, 0, 0, 0, 0),
497 Foam::scalar(0))
498 );
499 submeshA = autoPtr<fvMeshSubset>(new fvMeshSubset(mesh));
500
501 for (label i = 0; i < magicPointsArow().size(); i++)
502 {
503 indices = ITHACAutilities::getIndices(mesh, magicPointsArow()[i], layers);
504 indices.append(ITHACAutilities::getIndices(mesh, magicPointsAcol()[i],
505 layers));
506 totalMagicPointsA().append(indices);
507 }
508
509 uniqueMagicPointsA() = ITHACAutilities::combineList(totalMagicPointsA());
510#if OPENFOAM >= 1812
511 submeshA->setCellSubset(uniqueMagicPointsA());
512#else
513 submeshA->setLargeCellSubset(uniqueMagicPointsA());
514#endif
515 submeshA->subMesh().fvSchemes::readOpt() = mesh.fvSchemes::readOpt();
516 submeshA->subMesh().fvSolution::readOpt() = mesh.fvSolution::readOpt();
517 submeshA->subMesh().fvSchemes::read();
518 submeshA->subMesh().fvSolution::read();
519 std::cout.clear();
520 S f = submeshA->interpolate(field).ref();
521 scalar zerodot25 = 0.25;
522 ITHACAutilities::assignIF(Indici, zerodot25,
523 uniqueMagicPointsA().List<label>::clone()());
524 ITHACAutilities::assignONE(Indici, magicPointsArow());
525 ITHACAutilities::assignONE(Indici, magicPointsAcol());
526
527 if (!secondTime)
528 {
529 localMagicPointsArow = global2local(magicPointsArow(), submeshA());
530 localMagicPointsAcol = global2local(magicPointsAcol(), submeshA());
531 ITHACAstream::exportSolution(Indici, "1", "./ITHACAoutput/DEIM/" + MatrixName
532 );
533 totalMagicPointsA().write();
534 uniqueMagicPointsA().write();
535 }
536
537 runSubMeshA = true;
538 return f;
539}
540
541template<typename T>
542template<typename S>
543S DEIM<T>::generateSubmeshVector(label layers, const fvMesh& mesh, S field,
544 label secondTime)
545{
547 totalMagicPointsB = autoPtr<IOList<labelList>>
548 (
549 new IOList<labelList>
550 (
551 IOobject
552 (
553 "totalMagicPointsB",
554 para->runTime.time().constant(),
555 "../" + FolderM,
556 para->mesh,
557 IOobject::READ_IF_PRESENT,
558 IOobject::NO_WRITE
559 )
560 )
561 );
562 uniqueMagicPointsB = autoPtr<IOList<label>>
563 (
564 new IOList<label>
565 (
566 IOobject
567 (
568 "uniqueMagicPointsB",
569 para->runTime.time().constant(),
570 "../" + FolderM,
571 para->mesh,
572 IOobject::READ_IF_PRESENT,
573 IOobject::NO_WRITE
574 )
575 )
576 );
577 List<label> indices;
578 volScalarField Indici
579 (
580 IOobject
581 (
582 MatrixName + "_B_indices",
583 mesh.time().timeName(),
584 mesh,
585 IOobject::NO_READ,
586 IOobject::NO_WRITE
587 ),
588 mesh,
589 dimensionedScalar(MatrixName + "_B_indices", dimensionSet(0, 0, 0, 0, 0, 0, 0),
590 Foam::scalar(0))
591 );
592 submeshB = autoPtr<fvMeshSubset>(new fvMeshSubset(mesh));
593
594 for (label i = 0; i < magicPointsB().size(); i++)
595 {
596 indices = ITHACAutilities::getIndices(mesh, magicPointsB()[i], layers);
597 totalMagicPointsB().append(indices);
598
599 if (!secondTime)
600 {
601 ITHACAutilities::assignONE(Indici, indices);
602 }
603 }
604
605 uniqueMagicPointsB() = ITHACAutilities::combineList(totalMagicPointsB());
606 std::cout.setstate(std::ios_base::failbit);
607#if OPENFOAM >= 1812
608 submeshB->setCellSubset(uniqueMagicPointsB());
609#else
610 submeshB->setLargeCellSubset(uniqueMagicPointsB());
611#endif
612 submeshB->subMesh().fvSchemes::readOpt() = mesh.fvSchemes::readOpt();
613 submeshB->subMesh().fvSolution::readOpt() = mesh.fvSolution::readOpt();
614 submeshB->subMesh().fvSchemes::read();
615 submeshB->subMesh().fvSolution::read();
616 std::cout.clear();
617 S f = submeshB->interpolate(field).ref();
618 scalar zerodot25 = 0.25;
619 ITHACAutilities::assignIF(Indici, zerodot25,
620 uniqueMagicPointsB().List<label>::clone()());
621 ITHACAutilities::assignONE(Indici, magicPointsB());
622
623 if (!secondTime)
624 {
625 localMagicPointsB = global2local(magicPointsB(), submeshB());
626 ITHACAstream::exportSolution(Indici, "1", "./ITHACAoutput/DEIM/" + MatrixName
627 );
628 totalMagicPointsB().write();
629 uniqueMagicPointsB().write();
630 }
631
632 runSubMeshB = true;
633 return f;
634}
635
636
637template<typename T>
638List<label> DEIM<T>::global2local(List<label>& points,
639 fvMeshSubset& submesh)
640{
641 List<label> localPoints;
642
643 for (label i = 0; i < points.size(); i++)
644 {
645 for (label j = 0; j < submesh.cellMap().size(); j++)
646 {
647 if (submesh.cellMap()[j] == points[i])
648 {
649 localPoints.append(j);
650 break;
651 }
652 }
653 }
654
655 return localPoints;
656}
657
658template<typename T>
659void DEIM<T>::check3DIndices(label& ind_rowA, label& ind_colA, label& xyz_rowA,
660 label& xyz_colA)
661{
662 if (ind_rowA < Ncells)
663 {
664 xyz_rowA = 0;
665 }
666 else if (ind_rowA < Ncells * 2)
667 {
668 xyz_rowA = 1;
669 ind_rowA = ind_rowA - Ncells;
670 }
671 else
672 {
673 xyz_rowA = 2;
674 ind_rowA = ind_rowA - 2 * Ncells;
675 }
676
677 if (ind_colA < Ncells )
678 {
679 xyz_colA = 0;
680 }
681 else if (ind_colA < Ncells * 2)
682 {
683 xyz_colA = 1;
684 ind_colA = ind_colA - Ncells;
685 }
686 else
687 {
688 xyz_colA = 2;
689 ind_colA = ind_colA - 2 * Ncells;
690 }
691};
692
693template<typename T>
694void DEIM<T>::check3DIndices(label& ind_rowA, label& xyz_rowA)
695{
696 if (ind_rowA < Ncells)
697 {
698 xyz_rowA = 0;
699 }
700 else if (ind_rowA < Ncells * 2)
701 {
702 xyz_rowA = 1;
703 ind_rowA = ind_rowA - Ncells;
704 }
705 else
706 {
707 xyz_rowA = 2;
708 ind_rowA = ind_rowA - 2 * Ncells;
709 }
710};
711
712template<class T>
713template <class F>
715{
716 M_Assert(runSubMesh == true,
717 "You have to compute the magicPoints before calling this function, try to rerun generateSubmesh");
718 F f = submesh().interpolate(field).ref();
719 return f;
720}
721
722template<class T>
723template <class F>
725{
726 M_Assert(runSubMeshA == true,
727 "You have to compute the magicPointsA before calling this function, try to rerun generateSubmeshMatrix");
728 F f = submeshA().interpolate(field).ref();
729 return f;
730}
731
732template<typename T>
733template <class F>
735{
736 M_Assert(runSubMeshB == true,
737 "You have to compute the magicPointsB before calling this function, try to rerun generateSubmeshVector");
738 F f = submeshB().interpolate(field).ref();
739 return f;
740}
741
742
743template<> label DEIM<fvScalarMatrix>::getNcells(label sizeM)
744{
745 label Ncells = sizeM;
746 return Ncells;
747}
748
749template<> label DEIM<fvVectorMatrix>::getNcells(label sizeM)
750{
751 label Ncells = sizeM / 3;
752 return Ncells;
753}
754
755template<typename T>
756void DEIM<T>::setMagicPoints(labelList& newMagicPoints, labelList& newxyz)
757{
758 Folder = "ITHACAoutput/DEIM/" + FunctionName;
760 word command = "rm -r " + Folder + "/magicPoints";
761 system(command);
762 command = "rm -r " + Folder + "/xyz";
763 system(command);
764 command = "rm -r " + Folder + "/totalMagicPoints";
765 system(command);
766 command = "rm -r " + Folder + "/uniqueMagicPoints";
767 system(command);
768 magicPoints.reset
769 (
770 autoPtr<IOList<label>>
771 (
772 new IOList<label>
773 (
774 IOobject
775 (
776 "magicPoints",
777 para->runTime.time().constant(),
778 "../" + Folder,
779 para->mesh,
780 IOobject::READ_IF_PRESENT,
781 IOobject::NO_WRITE
782 )
783 )
784 )
785 );
786
787 for (label i = 0; i < newMagicPoints.size(); ++i)
788 {
789 magicPoints().append(newMagicPoints[i]);
790 }
791
792 xyz.reset
793 (
794 autoPtr<IOList<label>>
795 (
796 new IOList<label>
797 (
798 IOobject
799 (
800 "xyz",
801 para->runTime.time().constant(),
802 "../" + Folder,
803 para->mesh,
804 IOobject::READ_IF_PRESENT,
805 IOobject::NO_WRITE
806 )
807 )
808 )
809 );
810
811 for (label i = 0; i < newMagicPoints.size(); ++i)
812 {
813 xyz().append(newxyz[i]);
814 }
815
816 magicPoints().write();
817 xyz().write();
818}
819
820// Specialization of the constructor
821template DEIM<fvScalarMatrix>::DEIM (PtrList<fvScalarMatrix>& s,
822 label MaxModesA,
823 label MaxModesB, word MatrixName);
824template DEIM<fvVectorMatrix>::DEIM (PtrList<fvVectorMatrix>& s,
825 label MaxModesA,
826 label MaxModesB, word MatrixName);
827template DEIM<volScalarField>::DEIM(PtrList<volScalarField>& s, label MaxModes,
828 word FunctionName, word FieldName);
829template DEIM<volVectorField>::DEIM(PtrList<volVectorField>& s, label MaxModes,
830 word FunctionName, word FieldName);
831
832// Specialization for generateSubField
833template volVectorField DEIM<volScalarField>::generateSubField(
834 volVectorField& field);
835template volVectorField DEIM<volVectorField>::generateSubField(
836 volVectorField& field);
837template volScalarField DEIM<volScalarField>::generateSubField(
838 volScalarField& field);
839template volScalarField DEIM<volVectorField>::generateSubField(
840 volScalarField& field);
841template volTensorField DEIM<volScalarField>::generateSubField(
842 volTensorField& field);
843template volTensorField DEIM<volVectorField>::generateSubField(
844 volTensorField& field);
845template const volVectorField DEIM<volScalarField>::generateSubField(
846 const volVectorField& field);
847template const volVectorField DEIM<volVectorField>::generateSubField(
848 const volVectorField& field);
849template const volScalarField DEIM<volScalarField>::generateSubField(
850 const volScalarField& field);
851template const volScalarField DEIM<volVectorField>::generateSubField(
852 const volScalarField& field);
853template const volTensorField DEIM<volScalarField>::generateSubField(
854 const volTensorField& field);
855template const volTensorField DEIM<volVectorField>::generateSubField(
856 const volTensorField& field);
857
858// Specialization for generateSubFieldMatrix
860 volScalarField& field);
862 volVectorField& field);
863template surfaceScalarField
864DEIM<fvScalarMatrix>::generateSubFieldMatrix(surfaceScalarField& field);
865template surfaceVectorField
866DEIM<fvScalarMatrix>::generateSubFieldMatrix(surfaceVectorField& field);
868 volScalarField& field);
870 volVectorField& field);
871template surfaceScalarField
872DEIM<fvVectorMatrix>::generateSubFieldMatrix(surfaceScalarField& field);
873template surfaceVectorField
874DEIM<fvVectorMatrix>::generateSubFieldMatrix(surfaceVectorField& field);
875
876// Specialization for generateSubFieldVector
878 volScalarField& field);
880 volVectorField& field);
881template surfaceScalarField
882DEIM<fvScalarMatrix>::generateSubFieldVector(surfaceScalarField& field);
883template surfaceVectorField
884DEIM<fvScalarMatrix>::generateSubFieldVector(surfaceVectorField& field);
886 volScalarField& field);
888 volVectorField& field);
889template surfaceScalarField
890DEIM<fvVectorMatrix>::generateSubFieldVector(surfaceScalarField& field);
891template surfaceVectorField
892DEIM<fvVectorMatrix>::generateSubFieldVector(surfaceVectorField& field);
893
894// Specialization for generateSubmesh
895template volScalarField DEIM<volScalarField>::generateSubmesh(
896 label layers, const fvMesh& mesh, volScalarField field,
897 label secondTime);
898template volVectorField DEIM<volVectorField>::generateSubmesh(
899 label layers, const fvMesh& mesh, volVectorField field,
900 label secondTime);
901template volScalarField DEIM<volVectorField>::generateSubmesh(
902 label layers, const fvMesh& mesh, volScalarField field,
903 label secondTime);
904template volVectorField DEIM<volScalarField>::generateSubmesh(
905 label layers, const fvMesh& mesh, volVectorField field,
906 label secondTime);
907
908// Specialization for generateSubmeshVector
909template volScalarField DEIM<fvScalarMatrix>::generateSubmeshVector(
910 label layers, const fvMesh& mesh, volScalarField field, label secondTime);
911template volVectorField DEIM<fvScalarMatrix>::generateSubmeshVector(
912 label layers, const fvMesh& mesh, volVectorField field, label secondTime);
913template surfaceScalarField
914DEIM<fvScalarMatrix>::generateSubmeshVector(label layers, const fvMesh& mesh,
915 surfaceScalarField field, label secondTime);
916template surfaceVectorField
917DEIM<fvScalarMatrix>::generateSubmeshVector(label layers, const fvMesh& mesh,
918 surfaceVectorField field, label secondTime);
919template volScalarField DEIM<fvVectorMatrix>::generateSubmeshVector(
920 label layers, const fvMesh& mesh, volScalarField field, label secondTime);
921template volVectorField DEIM<fvVectorMatrix>::generateSubmeshVector(
922 label layers, const fvMesh& mesh, volVectorField field, label secondTime);
923template surfaceScalarField
924DEIM<fvVectorMatrix>::generateSubmeshVector(label layers, const fvMesh& mesh,
925 surfaceScalarField field, label secondTime);
926template surfaceVectorField
927DEIM<fvVectorMatrix>::generateSubmeshVector(label layers, const fvMesh& mesh,
928 surfaceVectorField field, label secondTime);
929
930// Specialization for generateSubmeshMatrix
931template volScalarField DEIM<fvScalarMatrix>::generateSubmeshMatrix(
932 label layers, const fvMesh& mesh, volScalarField field, label secondTime);
933template volVectorField DEIM<fvScalarMatrix>::generateSubmeshMatrix(
934 label layers, const fvMesh& mesh, volVectorField field, label secondTime);
935template surfaceScalarField
936DEIM<fvScalarMatrix>::generateSubmeshMatrix(label layers, const fvMesh& mesh,
937 surfaceScalarField field, label secondTime);
938template surfaceVectorField
939DEIM<fvScalarMatrix>::generateSubmeshMatrix(label layers, const fvMesh& mesh,
940 surfaceVectorField field, label secondTime);
941template volScalarField DEIM<fvVectorMatrix>::generateSubmeshMatrix(
942 label layers, const fvMesh& mesh, volScalarField field, label secondTime);
943template volVectorField DEIM<fvVectorMatrix>::generateSubmeshMatrix(
944 label layers, const fvMesh& mesh, volVectorField field, label secondTime);
945template surfaceScalarField
946DEIM<fvVectorMatrix>::generateSubmeshMatrix(label layers, const fvMesh& mesh,
947 surfaceScalarField field, label secondTime);
948template surfaceVectorField
949DEIM<fvVectorMatrix>::generateSubmeshMatrix(label layers, const fvMesh& mesh,
950 surfaceVectorField field, label secondTime);
951
952// specialization for setMagicPoints
953template void DEIM<volScalarField>::setMagicPoints(labelList& newMagicPoints,
954 labelList& newxyz);
955template void DEIM<volVectorField>::setMagicPoints(labelList& newMagicPoints,
956 labelList& newxyz);
Foam::fvMesh & mesh
Definition createMesh.H:47
#define M_Assert(Expr, Msg)
autoPtr< IOList< label > > magicPointsAcol
Magic points in the case of the a matrix function (cols indices)
Definition DEIM.H:125
autoPtr< IOList< label > > xyz
Definition of the x, y, z coordinate of the identified element in the matrix or source term 0 for x,...
Definition DEIM.H:151
label MaxModes
The maximum number of modes to be considered.
Definition DEIM.H:77
S generateSubmesh(label layers, const fvMesh &mesh, S field, label secondTime=0)
Function to generate the submesh for the nonlinear function case.
Definition DEIM.C:357
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
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 > > magicPoints
Magic points in the case of the a nonlinear function.
Definition DEIM.H:121
List< Eigen::SparseMatrix< double > > UA
The U matrix of the DEIM method.
Definition DEIM.H:116
autoPtr< IOList< label > > xyz_Acol
Definition of the x, y, z coordinate of the identified element in the matrix or source term 0 for x,...
Definition DEIM.H:153
word MatrixName
The name of the matrix.
Definition DEIM.H:101
Eigen::MatrixXd U
The U matrix of the DEIM method.
Definition DEIM.H:115
DEIM(PtrList< T > &SnapShotsMatrix, label MaxModes, word FunctionName, word FieldName)
Construct DEIM for non-linear function.
Definition DEIM.C:34
F generateSubFieldVector(F &field)
Function to generate a a subfield in the location of the magic points computed for the Matrix (RHS)
Definition DEIM.C:734
label MaxModesA
Definition DEIM.H:78
Eigen::MatrixXd MatrixModes
The matrix containing the modes.
Definition DEIM.H:104
autoPtr< IOList< label > > magicPointsArow
Magic points in the case of the a matrix function (rows indices)
Definition DEIM.H:123
word Folder
Folder for nonlinear functions.
Definition DEIM.H:144
std::tuple< List< Eigen::SparseMatrix< double > >, List< Eigen::VectorXd > > Matrix_Modes
Matrix Modes.
Definition DEIM.H:74
PtrList< T > modes
The POD modes of the DEIM procedure that can be.
Definition DEIM.H:70
word FolderM
Folder in the matrix case.
Definition DEIM.H:146
void check3DIndices(label &ind_rowA, label &ind_colA, label &xyz_rowA, label &xyz_colA)
check3DIndices in case of three dimensional fields
Definition DEIM.C:659
Eigen::SparseMatrix< double > P
The P matrix of the DEIM method.
Definition DEIM.H:167
F generateSubField(F &field)
Function to generate a a subfield in the location of the magic points.
Definition DEIM.C:714
autoPtr< IOList< label > > xyz_Arow
Definition of the x, y, z coordinate of the identified element in the matrix or source term 0 for x,...
Definition DEIM.H:152
Eigen::MatrixXd MatrixOnline
Online Matrices.
Definition DEIM.H:108
List< Eigen::SparseMatrix< double > > MatrixOnlineA
Online Matrices.
Definition DEIM.H:109
label MaxModesB
Definition DEIM.H:79
label getNcells(label sizeM)
get the number of cells from the dimension of a LHS matrix
autoPtr< IOList< label > > xyz_B
Definition of the x, y, z coordinate of the identified element in the matrix or source term 0 for x,...
Definition DEIM.H:154
F generateSubFieldMatrix(F &field)
Function to generate a a subfield in the location of the magic points computed for the Matrix (LHS)
Definition DEIM.C:724
autoPtr< IOList< label > > magicPointsB
Magic points in the case of the a matrix function, right hand side.
Definition DEIM.H:127
Eigen::SparseMatrix< double > PB
The P matrix of the DEIM method.
Definition DEIM.H:169
void setMagicPoints(labelList &newMagicPoints, labelList &newxyz)
Function to set the magic points externally.
Definition DEIM.C:756
PtrList< T > SnapShotsMatrix
The snapshots matrix containing the nonlinear function or operator.
Definition DEIM.H:67
label Ncells
Int Number of Cells;.
Definition DEIM.H:95
Eigen::MatrixXd UB
The U matrix of the DEIM method.
Definition DEIM.H:117
word FunctionName
The name of the non-linear function.
Definition DEIM.H:82
List< Eigen::SparseMatrix< double > > PA
The P matrix of the DEIM method.
Definition DEIM.H:168
List< label > global2local(List< label > &points, fvMeshSubset &submesh)
Get local indices in the submeshe from indices in the global ones.
Definition DEIM.C:638
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:649
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
Time & runTime
runTime defined locally
fvMesh & mesh
Mesh object defined locally.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
volScalarField & A
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.
T max(Eigen::SparseMatrix< T > &mat, label &ind_row, label &ind_col)
Find the maximum of a sparse Matrix (Useful for DEIM)
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > innerProduct(List< Eigen::SparseMatrix< T > > &A, List< Eigen::SparseMatrix< T > > &B)
Perform Frobenius inner Product between two list of sparse matrices A and B.
List< Eigen::SparseMatrix< T > > MMproduct(List< Eigen::SparseMatrix< T > > &A, Eigen::DenseBase< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > > &C)
Sparse Matrix - Dense Matrix product between a list of sparse matrices and a dense matrix.
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 save(const List< Eigen::SparseMatrix< T > > &MatrixList, word folder, word MatrixName)
void load(List< Eigen::SparseMatrix< T > > &MatrixList, word folder, word MatrixName)
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 assignIF(GeometricField< Type, fvPatchField, volMesh > &s, Type value)
Assign internal field.
List< label > getIndices(const fvMesh &mesh, int index, int layers)
Gets the indices of the cells around a certain cell.
void assignONE(volScalarField &s, List< label > &L)
Assign one to volScalarField.
List< T > combineList(List< List< T > > &doubleList)
Combine a list of list into a single list with unique elements.
void save(const Eigen::Matrix< T, -1, dim > &mat, const std::string fname)
Eigen::Matrix< T, -1, dim > load(Eigen::Matrix< T, -1, dim > &mat, const std::string fname, std::string order="rowMajor")
label i
Definition pEqn.H:46
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))