Loading...
Searching...
No Matches
ITHACAcoeffsMass.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/>.
24Namespace
25 ITHACAutilites
26Description
27 Utilities to compute projection coefficients and mass matrices
28SourceFiles
29 ITHACAcoeffsMass.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef ITHACAcoeffsMass_H
38#define ITHACAcoeffsMass_H
39
40#include "fvCFD.H"
41#include "IOmanip.H"
42#include "freestreamFvPatchField.H"
43#include <sys/stat.h>
44#include <unistd.h>
45#pragma GCC diagnostic push
46#pragma GCC diagnostic ignored "-Wold-style-cast"
47#include <Eigen/Eigen>
48#pragma GCC diagnostic pop
49#include <functional>
50#include "./colormod.H"
51#include "polyMeshTools.H"
52#include <chrono>
53#include "mixedFvPatchFields.H"
54#include "fvMeshSubset.H"
55using namespace std::placeholders;
56#include "Foam2Eigen.H"
57#include "ITHACAstream.H"
58
60namespace ITHACAutilities
61{
62
63//--------------------------------------------------------------------------
77template<class Type, template<class> class PatchField, class GeoMesh>
78PtrList<GeometricField<Type, PatchField, GeoMesh>> reconstructFromCoeff(
79 PtrList<GeometricField<Type, PatchField, GeoMesh>>& modes,
80 Eigen::MatrixXd& coeff_matrix, label Nmodes);
81
82//-----------------------------------------------------------------------------
94template<class Type, template<class> class PatchField, class GeoMesh>
95Eigen::MatrixXd getMassMatrix(PtrList<GeometricField<Type, PatchField, GeoMesh>>&
96 fields, label nModes = 0,
97 bool consider_volumes = true);
98
99//-----------------------------------------------------------------------------
111template<class Type, template<class> class PatchField, class GeoMesh>
112Eigen::MatrixXd getMassMatrix(PtrList<GeometricField<Type, PatchField, GeoMesh>>&
113 fields,
114 PtrList<GeometricField<Type, PatchField, GeoMesh>>&
115 fields2, label nModes = 0,
116 bool consider_volumes = true);
117
118//-----------------------------------------------------------------------------
132template<class Type, template<class> class PatchField, class GeoMesh>
133Eigen::MatrixXd getMassMatrix(PtrList<GeometricField<Type, PatchField, GeoMesh>>& fields,
134 PtrList<GeometricField<Type, PatchField, GeoMesh>>& fields2,
135 Eigen::VectorXd weights,
136 label Nmodes = 0,
137 bool consider_volumes = true);
138
139
140//--------------------------------------------------------------------------
151template<class Type, template<class> class PatchField, class GeoMesh>
152Eigen::VectorXd getMassMatrixFV(
153 GeometricField<Type, PatchField, GeoMesh>& snapshot);
154
155//------------------------------------------------------------------------------
170template<class Type, template<class> class PatchField, class GeoMesh>
171Eigen::VectorXd getCoeffs(GeometricField<Type, PatchField, GeoMesh>&
172 snapshot,
173 PtrList<GeometricField<Type, PatchField, GeoMesh>>& modes, label Nmodes = 0,
174 bool consider_volumes = true);
175
176//------------------------------------------------------------------------------
191template<class Type, template<class> class PatchField, class GeoMesh>
192Eigen::MatrixXd getCoeffs(PtrList<GeometricField<Type, PatchField, GeoMesh>>&
193 snapshot,
194 PtrList<GeometricField<Type, PatchField, GeoMesh>>& modes, label Nmodes = 0,
195 bool consider_volumes = true);
196
197//--------------------------------------------------------------------------
209Eigen::MatrixXd parTimeCombMat(List<Eigen::VectorXd>
210 acquiredSnapshotsTimes, Eigen::MatrixXd parameters);
211
212}
213
214#endif
Header file of the Foam2Eigen class.
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Simple header and source file of the Color::Modifier class to change color to the output stream.
Namespace to implement some useful assign operation of OF fields.
Eigen::VectorXd getCoeffs(GeometricField< Type, PatchField, GeoMesh > &snapshot, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes, bool consider_volumes)
Projects a snapshot on a basis function and gets the coefficients of the projection.
PtrList< GeometricField< Type, PatchField, GeoMesh > > reconstructFromCoeff(PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, Eigen::MatrixXd &coeff_matrix, label Nmodes)
Exact reconstruction using a certain number of modes for a list of fields and the projection coeffici...
Eigen::VectorXd getMassMatrixFV(GeometricField< Type, PatchField, GeoMesh > &snapshot)
Gets a vector containing the volumes of each cell of the mesh.
Eigen::MatrixXd parTimeCombMat(List< Eigen::VectorXd > acquiredSnapshotsTimes, Eigen::MatrixXd parameters)
A method to compute the time-parameter combined matrix whose any single element corresponds to a uniq...
Eigen::MatrixXd getMassMatrix(PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes, bool consider_volumes)
Gets the mass matrix using the eigen routine.