Loading...
Searching...
No Matches
ReducedLaplacian.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 "ReducedLaplacian.H"
35
36
37// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
38
39// Constructor
43
45 :
46 problem(&problem)
47{
48}
49
50void reducedLaplacian::solveOnline(Eigen::MatrixXd mu)
51{
52 if (mu.cols() != problem->A_matrices.size())
53 {
54 Info << "wrong dimension of online parameters" << endl;
55 exit(0);
56 }
57
58 Eigen::MatrixXd A;
59 A.setZero(problem->NTmodes, problem->NTmodes);
60
61 for (int i = 0; i < problem->A_matrices.size() ; i++)
62 {
63 A += problem->A_matrices[i] * mu(0, i);
64 }
65
66 Eigen::MatrixXd x;
67 x = A.fullPivLu().solve(-problem->source);
68 online_solution.conservativeResize(count_online_solve, problem->NTmodes + 1);
71 x.transpose();
73}
74
75void reducedLaplacian::reconstruct(fileName folder, int printevery)
76{
77 mkDir(folder);
79 int counter = 0;
80 int nextwrite = 0;
81
82 for (int i = 0; i < online_solution.rows(); i++)
83 {
84 if (counter == nextwrite)
85 {
86 volScalarField T_rec("T_rec", problem->Tmodes[0] * 0);
87
88 for (int j = 0; j < problem->NTmodes; j++)
89 {
90 T_rec += problem->Tmodes[j] * online_solution(i, j + 1);
91 }
92
93 Trec.append((T_rec).clone());
94 ITHACAstream::exportSolution(T_rec, name(online_solution(i, 0)), folder);
95 nextwrite += printevery;
96 }
97
98 counter++;
99 }
100}
101
102
103
104
105// ************************************************************************* //
106
Header file of the reducedLaplacian class.
Class to implement a full order laplacian parametrized problem.
List< Eigen::MatrixXd > A_matrices
A matrices.
volScalarModes Tmodes
List of POD modes.
Eigen::MatrixXd source
Source vector.
label NTmodes
Number of modes reduced problem.
void reconstruct(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Function to recover the solution given the online solution.
PtrList< volScalarField > Trec
Reduced reconstructed Solution.
laplacianProblem * problem
Problem object.
reducedLaplacian()
Construct Null.
int count_online_solve
Counter for online sol.
Eigen::MatrixXd online_solution
Online solution.
virtual void solveOnline()
Virtual Method to perform and online Solve.
volScalarField & A
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 createSymLink(word folder)
Creates symbolic links to 0, system and constant.
label i
Definition pEqn.H:46