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