Loading...
Searching...
No Matches
laplacianProblem.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
13 License
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
35#include "laplacianProblem.H"
36
37// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
38
39// Constructors
42{
43 _args = autoPtr<argList>
44 (
45 new argList(argc, argv)
46 );
47
48 if (!_args->checkRootCase())
49 {
50 Foam::FatalError.exit();
51 }
52
53 argList& args = _args();
54#include "createTime.H"
55#include "createMesh.H"
56#include "createFields.H"
59}
60
61
62// * * * * * * * * * * * * * * Full Order Methods * * * * * * * * * * * * * * //
63// Method to performa a truthSolve
64void laplacianProblem::truthSolve(List<scalar> mu_now, word folder)
65{
66 volScalarField& T = _T();
67 volScalarField& S = _S();
68 fvScalarMatrix lhs(theta[0]*operator_list[0]);
69
70 for (label i = 1; i < operator_list.size(); i++)
71 {
72 lhs += theta[i] * operator_list[i];
73 }
74
75 solve(lhs == -S);
77 Tfield.append(T.clone());
78 counter++;
79 writeMu(mu_now);
80 // --- Fill in the mu_samples with parameters (mu) to be used for the PODI sample points
81 mu_samples.conservativeResize(mu_samples.rows() + 1, mu_now.size());
82
83 for (label i = 0; i < mu_now.size(); i++)
84 {
85 mu_samples(mu_samples.rows() - 1, i) = mu_now[i];
86 }
87
88 // Resize to Unitary if not initialized by user (i.e. non-parametric problem)
89 if (mu.cols() == 0)
90 {
91 mu.resize(1, 1);
92 }
93
94 if (mu_samples.rows() == mu.cols())
95 {
96 ITHACAstream::exportMatrix(mu_samples, "mu_samples", "eigen",
97 folder);
98 }
99}
100// Perform the projection onto the POD modes
102{
103 NTmodes = Nmodes;
104 A_matrices.resize(operator_list.size());
105 source.resize(Nmodes, 1);
106 volScalarField& S = _S();
107
108 for (label i = 0; i < operator_list.size(); i++)
109 {
110 A_matrices[i].resize(Nmodes, Nmodes);
111
112 for (label j = 0; j < Nmodes; j++)
113 {
114 if (i == 0)
115 {
116 source(j, 0) = fvc::domainIntegrate( Tmodes[j] * S).value();
117 }
118
119 for (label k = 0; k < Nmodes; k++)
120 {
121 A_matrices[i](j, k) = fvc::domainIntegrate( Tmodes[j] * fvc::laplacian(
122 nu_list[i], Tmodes[k])).value();
123 }
124 }
125
126 if (Pstream::parRun())
127 {
128 reduce(A_matrices[i], sumOp<Eigen::MatrixXd>());
129 }
130 }
131
132 if (Pstream::parRun())
133 {
134 reduce(source, sumOp<Eigen::MatrixXd>());
135 }
136
139 "./ITHACAoutput/Matrices/");
141 "./ITHACAoutput/Matrices/");
143 "./ITHACAoutput/Matrices/A_matrices");
145 ITHACAstream::exportMatrix(source, "S", "python", "./ITHACAoutput/Matrices/");
146 ITHACAstream::exportMatrix(source, "S", "matlab", "./ITHACAoutput/Matrices/");
147 ITHACAstream::exportMatrix(source, "S", "eigen", "./ITHACAoutput/Matrices/");
148}
149
150
151
152
TEqn solve()
PtrList< volScalarField > nu_list
Nu (diffusivity)
List< scalar > theta
Theta (coefficients of the affine expansion)
PtrList< fvScalarMatrix > operator_list
List of operators.
autoPtr< volScalarField > _T
Temperature field.
List< Eigen::MatrixXd > A_matrices
A matrices.
volScalarModes Tmodes
List of POD modes.
Eigen::MatrixXd source
Source vector.
autoPtr< volScalarField > _S
Source Term.
PtrList< volScalarField > Tfield
List of snapshots for the solution.
label NTmodes
Number of modes reduced problem.
void project()
General projection operation.
void writeMu(List< scalar > mu_now)
Write out a list of scalar corresponding to the parameters used in the truthSolve.
Eigen::MatrixXd mu_samples
Matrix of parameters to be used for PODI, where each row corresponds to a sample point....
label counter
Counter used for the output of the full order solutions.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
Eigen::MatrixXd mu
Row matrix of parameters.
autoPtr< argList > _args
argList
bool podex
Boolean variable, it is 1 if the POD has already been computed, else 0.
void truthSolve()
Perform a TruthSolve.
volScalarField & T
Definition createT.H:46
Header file of the laplacianProblem class.
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 ...
bool check_pod()
Check if the POD data folder "./ITHACAoutput/POD" exists.
bool check_off()
Check if the offline data folder "./ITHACAoutput/Offline" exists.
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))