Loading...
Searching...
No Matches
Modes.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 EigenFunctions
26Description
27 Container that contains a list of modes with additional operations
28SourceFiles
29 Modes.C
30\*---------------------------------------------------------------------------*/
31
36
37
38#ifndef Modes_H
39#define Modes_H
40#pragma GCC diagnostic push
41#pragma GCC diagnostic ignored "-Wold-style-cast"
42#include <Eigen/Eigen>
43#include <unsupported/Eigen/SparseExtra>
44#include <unsupported/Eigen/CXX11/Tensor>
45#pragma GCC diagnostic pop
46#include "fvCFD.H"
47#include "Foam2Eigen.H"
48#include "ITHACAutilities.H"
49#include "ITHACAstream.H"
50
51
52// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53
54/*---------------------------------------------------------------------------*\
55 Class Modes Declaration
56\*---------------------------------------------------------------------------*/
57
58//--------------------------------------------------------------------------
67template<class Type, template<class> class PatchField, class GeoMesh>
68class Modes : public PtrList<GeometricField<Type, PatchField, GeoMesh>>
69{
70 public:
71
73 List<Eigen::MatrixXd> EigenModes;
74
76 label NBC;
77
79 List<Eigen::MatrixXd> toEigen();
80
81 //--------------------------------------------------------------------------
86 PtrList<GeometricField<Type, PatchField, GeoMesh>>& toPtrList()
87 {
88 return static_cast<PtrList<GeometricField<Type, PatchField, GeoMesh>>&>(*this);
89 }
90
91 //----------------------------------------------------------------------
113 List<Eigen::MatrixXd> project(fvMatrix<Type>& Af, label numberOfModes = 0,
114 word projType = "G");
115
116 //----------------------------------------------------------------------
133 Eigen::MatrixXd project(GeometricField<Type, PatchField, GeoMesh>& field,
134 label numberOfModes = 0, word projType = "G", fvMatrix<Type>* Af = NULL);
135
136 //----------------------------------------------------------------------
153 Eigen::MatrixXd project(PtrList<GeometricField<Type, PatchField, GeoMesh>>&
154 fields,
155 label numberOfModes = 0, word projType = "G", fvMatrix<Type>* Af = NULL);
156
157 //----------------------------------------------------------------------
167 GeometricField<Type, PatchField, GeoMesh> projectSnapshot(
168 GeometricField<Type, PatchField, GeoMesh>& field,
169 label numberOfModes = 0, word projType = "G", fvMatrix<Type>* Af = NULL);
170
171 //----------------------------------------------------------------------
183 GeometricField<Type, PatchField, GeoMesh> reconstruct(
184 GeometricField<Type, PatchField, GeoMesh>& inputField, Eigen::MatrixXd Coeff,
185 word Name);
186
187 //----------------------------------------------------------------------
199 PtrList<GeometricField<Type, PatchField, GeoMesh>> reconstruct(
200 GeometricField<Type, PatchField, GeoMesh>& inputField,
201 List < Eigen::MatrixXd> Coeff,
202 word Name);
203
204 //----------------------------------------------------------------------
211 void projectSnapshots(PtrList<GeometricField<Type, PatchField, GeoMesh>>
212 snapshots, PtrList<GeometricField<Type, PatchField, GeoMesh>>& projSnapshots,
213 label numberOfModes = 0, word innerProduct = "L2");
214
215 //----------------------------------------------------------------------
221 void projectSnapshots(PtrList<GeometricField<Type, PatchField, GeoMesh>>
222 snapshots, PtrList<GeometricField<Type, PatchField, GeoMesh>>& projSnapshots,
223 word innerProduct = "L2");
224
225 //----------------------------------------------------------------------
233 void projectSnapshots(PtrList<GeometricField<Type, PatchField, GeoMesh>>
234 snapshots, PtrList<GeometricField<Type, PatchField, GeoMesh>>& projSnapshots,
235 PtrList<volScalarField> Volumes,
236 label numberOfModes = 0, word innerProduct = "L2");
237
238 //----------------------------------------------------------------------
245 void projectSnapshots(PtrList<GeometricField<Type, PatchField, GeoMesh>>
246 snapshots, PtrList<GeometricField<Type, PatchField, GeoMesh>>& projSnapshots,
247 PtrList<volScalarField> Volumes,
248 word innerProduct);
249
250 void operator=(const PtrList<GeometricField<Type, PatchField, GeoMesh>>& modes);
251};
252
256
257#endif
258
259
Header file of the Foam2Eigen class.
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the ITHACAutilities namespace.
Modes< vector, fvPatchField, volMesh > volVectorModes
Definition Modes.H:254
Modes< scalar, fvPatchField, volMesh > volScalarModes
Definition Modes.H:253
Modes< scalar, fvsPatchField, surfaceMesh > surfaceScalarModes
Definition Modes.H:255
Implementation of a container class derived from PtrList.
Definition Modes.H:69
GeometricField< Type, PatchField, GeoMesh > reconstruct(GeometricField< Type, PatchField, GeoMesh > &inputField, Eigen::MatrixXd Coeff, word Name)
Function to reconstruct the solution starting from the coefficients, in this case the field is passed...
Definition Modes.C:275
List< Eigen::MatrixXd > EigenModes
List of Matrices that contains the internalField and the additional matrices for the boundary patches...
Definition Modes.H:73
List< Eigen::MatrixXd > toEigen()
Method that convert a PtrList of modes into Eigen matrices filling the EigenModes object.
Definition Modes.C:38
label NBC
Number of patches.
Definition Modes.H:76
List< Eigen::MatrixXd > project(fvMatrix< Type > &Af, label numberOfModes=0, word projType="G")
A function that project an FvMatrix (OpenFoam linear System) on the modes.
Definition Modes.C:63
void projectSnapshots(PtrList< GeometricField< Type, PatchField, GeoMesh > > snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &projSnapshots, label numberOfModes=0, word innerProduct="L2")
Function to project a list of fields into the modes manifold.
Definition Modes.C:408
PtrList< GeometricField< Type, PatchField, GeoMesh > > & toPtrList()
Function that returns the Modes object as a standard PtrList.
Definition Modes.H:86
GeometricField< Type, PatchField, GeoMesh > projectSnapshot(GeometricField< Type, PatchField, GeoMesh > &field, label numberOfModes=0, word projType="G", fvMatrix< Type > *Af=NULL)
A function that project and reconstruct a snapshot on the modes.
Definition Modes.C:195
void operator=(const PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes)
Definition Modes.C:477