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
41
43{
44 _args = autoPtr<argList>
45 (
46 new argList(argc, argv)
47 );
48
49 if (!_args->checkRootCase())
50 {
51 Foam::FatalError.exit();
52 }
53
54 argList& args = _args();
55#include "createTime.H"
56#include "createMesh.H"
57#include "createFields.H"
60}
61
62
63// * * * * * * * * * * * * * * Full Order Methods * * * * * * * * * * * * * * //
64// Method to performa a truthSolve
65void laplacianProblem::truthSolve(List<scalar> mu_now, word folder)
66{
67 volScalarField& T = _T();
68 volScalarField& S = _S();
69 fvScalarMatrix lhs(theta[0] * operator_list[0]);
70
71 for (label i = 1; i < operator_list.size(); i++)
72 {
73 lhs += theta[i] * operator_list[i];
74 }
75
76 solve(lhs == -S);
78 Tfield.append(T.clone());
79 counter++;
80 writeMu(mu_now);
81 // --- Fill in the mu_samples with parameters (mu) to be used for the PODI sample points
82 mu_samples.conservativeResize(mu_samples.rows() + 1, mu_now.size());
83
84 for (label i = 0; i < mu_now.size(); i++)
85 {
86 mu_samples(mu_samples.rows() - 1, i) = mu_now[i];
87 }
88
89 // Resize to Unitary if not initialized by user (i.e. non-parametric problem)
90 if (mu.cols() == 0)
91 {
92 mu.resize(1, 1);
93 }
94
95 if (mu_samples.rows() == mu.cols())
96 {
97 ITHACAstream::exportMatrix(mu_samples, "mu_samples", "eigen",
98 folder);
99 }
100}
101
102// Perform the projection onto the POD modes
104{
105 NTmodes = Nmodes;
106 A_matrices.resize(operator_list.size());
107 source.resize(Nmodes, 1);
108 volScalarField& S = _S();
109
110 for (label i = 0; i < operator_list.size(); i++)
111 {
112 A_matrices[i].resize(Nmodes, Nmodes);
113
114 for (label j = 0; j < Nmodes; j++)
115 {
116 if (i == 0)
117 {
118 source(j, 0) = fvc::domainIntegrate( Tmodes[j] * S).value();
119 }
120
121 for (label k = 0; k < Nmodes; k++)
122 {
123 A_matrices[i](j, k) = fvc::domainIntegrate( Tmodes[j] * fvc::laplacian(
124 nu_list[i], Tmodes[k])).value();
125 }
126 }
127
128 if (Pstream::parRun())
129 {
130 reduce(A_matrices[i], sumOp<Eigen::MatrixXd>());
131 }
132 }
133
134 if (Pstream::parRun())
135 {
136 reduce(source, sumOp<Eigen::MatrixXd>());
137 }
138
141 "./ITHACAoutput/Matrices/");
143 "./ITHACAoutput/Matrices/");
145 "./ITHACAoutput/Matrices/A_matrices");
147 ITHACAstream::exportMatrix(source, "S", "python", "./ITHACAoutput/Matrices/");
148 ITHACAstream::exportMatrix(source, "S", "matlab", "./ITHACAoutput/Matrices/");
149 ITHACAstream::exportMatrix(source, "S", "eigen", "./ITHACAoutput/Matrices/");
150}
151
152
153
154
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))