Loading...
Searching...
No Matches
TurbDiffusionHyperreduction.H
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-------------------------------------------------------------------------------
12License
13 This file is part of ITHACA-FV
14 ITHACA-FV is free software: you can redistribute it and/or modify
15 it under the terms of the GNU Lesser General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18 ITHACA-FV is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU Lesser General Public License for more details.
22 You should have received a copy of the GNU Lesser General Public License
23 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
24Class
25 TurbDiffusionHyperreduction
26Description
27 Class to compute and hyperreduce the non-polynomial diffusion for a NON-Stationary turbulent NS problem
28SourceFiles
29 TurbDiffusionHyperreduction.C
30
31\*---------------------------------------------------------------------------*/
32
37
38#ifndef TurbDiffusionHyperreduction_H
39#define TurbDiffusionHyperreduction_H
40
41#include "UnsteadyNSTurb.H"
42#include "PODTemplate.H"
43#include "PODTemplateH1.H"
44#include "hyperReduction.templates.H"
45
46
47template<typename volF, typename T, typename S>
49{
50 protected :
53
55 autoPtr<UnsteadyNSTurb> ithacaUnsteadyNSTurb;
56
58 autoPtr<HyperReduction<PtrList<S>>> ithacaHyperreduction;
59
62
65
68
71 label POD_nSnapshot = -1;
72 label POD_startTime = -1;
73 label POD_endTime = -1;
74
77
80
82 PtrList<volF> f_spatialModesU; // TO DO : need to be duplicated (local => global) if //
83
85 volF* f_meanU; // TO DO : need to be duplicated (local => global) if //
86
87 // Spatial POD-HR modes
88 PtrList<S> f_spatialModesHR;
89
90 // Temporal modes of HR field
91 Eigen::MatrixXd m_temporalModesHR;
92
93 // Test temporal modes of HR field
94 Eigen::MatrixXd m_temporalModesHRSimulation;
95
96 // Covariance matrix the resolved modes of HR field
97 Eigen::MatrixXd covMatrixHR;
98
101
103 Foam::Time runTime2;
104
108
109 // Name of the HR method: DEIM or ECP
110 word HRMethod;
111
112 // Name of the interpolated field: fullStressFunction or nut or their reduced version
113 word HRInterpField;
114
115 // Name of the snapshots' field if ECP
116 word HRSnapshotsField;
117
118 // Name of the algorithm choice if ECP: Global or EachMode
119 word ECPAlgo;
120
121 // Centered field to perform the hyperreduction
122 bool interpFieldCenteredOrNot;
123
124
125
126 public :
127
128 //--------------------------------------------------------------------------
140 volF*& meanU, PtrList<volF>& spatialModesU);
141
142 //--------------------------------------------------------------------------
145
146
147 //--------------------------------------------------------------------------
153 void precomputeTurbDiffusionFunctions(word& fieldToCompute);
154
155 //--------------------------------------------------------------------------
159 template<typename volGradF>
160 void computeDefTensor(volGradF tensorMeanU);
161
162 //--------------------------------------------------------------------------
167
168 //--------------------------------------------------------------------------
173
174};
175
176
177#endif
Class that contains all parameters of the stochastic POD.
void computeTurbDiffusionHyperreduction()
Compute the hyperreduction by callng other methods.
PtrList< volF > f_spatialModesU
Pointer list of POD Modes of U.
label l_nSnapshotSimulation
Number of Simulation snapshot, also number of time steps.
S template_HRSnapshotsField
Snapshots field template.
autoPtr< UnsteadyNSTurb > ithacaUnsteadyNSTurb
Parameter object to evaluate interpolated functions.
Foam::Time runTime2
Dummy Foam::Time object for mesh and submesh management.
autoPtr< HyperReduction< PtrList< S > > > ithacaHyperreduction
Hyperreduction objects to choose magic points.
word snapshotsFolder
Path to snapshots, number, start and end time indices for each HR POD. If -1, set automatically by th...
label l_startTime
start time of the POD decomposition
~TurbDiffusionHyperreduction()
Class TurbDiffusionHyperreduction Destructor.
bool mgPoints_0_Neighborhoods_1
Boolean to include or not magic points' neighbourhood in the submesh 0 if only magicPoints and 1 if i...
word folder_HR
Export folder for HR operators.
label l_nSnapshot
Number of snapshot, also number of time step.
StoredParameters * m_parameters
Parameter object containing main information.
volF * f_meanU
Pointer of mean of U.
void precomputeTurbDiffusionFunctions(word &fieldToCompute)
Precompute and save non polynomial field.
void common_MeanHRInterpField()
Compute the interpolated field's mean on the magic points submesh.
T template_HRInterpField
Interpolated field template.
TurbDiffusionHyperreduction(Parameters *myParameters, bool mgPoints_0_Neighborhoods_1, T &template_HRInterpField, S &template_HRSnapshotsField, volF *&meanU, PtrList< volF > &spatialModesU)
Class TurbDiffusionHyperreduction Constructor.
void computeDefTensor(volGradF tensorMeanU)
Compute the symmetric tensor of the velocity modes on the magic points submesh.