Loading...
Searching...
No Matches
ITHACAcoeffsMass.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-------------------------------------------------------------------------------
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/>.
24\*---------------------------------------------------------------------------*/
25#include "ITHACAcoeffsMass.H"
26
27namespace ITHACAutilities
28{
29
30template<class Type, template<class> class PatchField, class GeoMesh>
31PtrList<GeometricField<Type, PatchField, GeoMesh>> reconstructFromCoeff(
32 PtrList<GeometricField<Type, PatchField, GeoMesh>>& modes,
33 Eigen::MatrixXd& coeff_matrix, label Nmodes)
34{
35 PtrList<GeometricField<Type, PatchField, GeoMesh>> rec_field;
36 rec_field.resize(0);
37
38 for (label k = 0; k < coeff_matrix.cols(); k++)
39 {
40 for (label i = 0; i < Nmodes; i++)
41 if ( i == 0)
42 {
43 rec_field.append(modes[i]*coeff_matrix(i, k));
44 }
45 else
46 {
47 rec_field[k] += modes[i] * coeff_matrix(i, k);
48 }
49 }
50
51 return rec_field;
52}
53
54template PtrList<GeometricField<scalar, fvPatchField, volMesh>>
56 PtrList<GeometricField<scalar, fvPatchField, volMesh>>& modes,
57 Eigen::MatrixXd& coeff_matrix, label Nmodes);
58template PtrList<GeometricField<vector, fvPatchField, volMesh>>
60 PtrList<GeometricField<vector, fvPatchField, volMesh>>& modes,
61 Eigen::MatrixXd& coeff_matrix, label Nmodes);
62template PtrList<GeometricField<tensor, fvPatchField, volMesh>>
64 PtrList<GeometricField<tensor, fvPatchField, volMesh>>& modes,
65 Eigen::MatrixXd& coeff_matrix, label Nmodes);
66
67template<class Type, template<class> class PatchField, class GeoMesh>
68Eigen::MatrixXd getMassMatrix(
69 PtrList<GeometricField<Type, PatchField, GeoMesh>>& modes, label Nmodes,
70 bool consider_volumes)
71{
72 label Msize;
73
74 if (Nmodes == 0)
75 {
76 Msize = modes.size();
77 }
78 else
79 {
80 Msize = Nmodes;
81 }
82
83 M_Assert(modes.size() >= Msize,
84 "The Number of requested modes is larger then the available quantity.");
85 Eigen::MatrixXd F = Foam2Eigen::PtrList2Eigen(modes);
86 Eigen::MatrixXd M;
87
88 if (consider_volumes)
89 {
90 Eigen::VectorXd V = getMassMatrixFV(modes[0]);
91 M = F.transpose().topRows(Msize) * V.asDiagonal() * F.leftCols(Msize);
92 }
93 else
94 {
95 M = F.transpose().topRows(Msize) * F.leftCols(Msize);
96 }
97
98 if (Pstream::parRun())
99 {
100 reduce(M, sumOp<Eigen::MatrixXd>());
101 }
102
103 return M;
104}
105
106template Eigen::MatrixXd getMassMatrix(
107 PtrList<GeometricField<scalar, fvPatchField, volMesh>>& modes, label Nmodes,
108 bool consider_volumes);
109
110template Eigen::MatrixXd getMassMatrix(
111 PtrList<GeometricField<scalar, fvsPatchField, surfaceMesh>>& modes,
112 label Nmodes,
113 bool consider_volumes = false);
114
115template Eigen::MatrixXd getMassMatrix(
116 PtrList<GeometricField<vector, fvPatchField, volMesh>>& modes, label Nmodes,
117 bool consider_volumes);
118
119template<class Type, template<class> class PatchField, class GeoMesh>
120Eigen::VectorXd getMassMatrixFV(
121 GeometricField<Type, PatchField, GeoMesh>& snapshot)
122{
123 Eigen::MatrixXd snapEigen = Foam2Eigen::field2Eigen(snapshot);
124 label dim = std::nearbyint(snapEigen.rows() / (snapshot.mesh().V()).size());
125 Eigen::VectorXd volumes = Foam2Eigen::field2Eigen(snapshot.mesh().V());
126 Eigen::VectorXd vol3 = volumes.replicate(dim, 1);
127 return vol3;
128}
129
130template Eigen::VectorXd getMassMatrixFV(
131 GeometricField<scalar, fvPatchField, volMesh>& snapshot);
132template Eigen::VectorXd getMassMatrixFV(
133 GeometricField<vector, fvPatchField, volMesh>& snapshot);
134
135template<class Type, template<class> class PatchField, class GeoMesh>
136Eigen::VectorXd getCoeffs(GeometricField<Type, PatchField, GeoMesh>&
137 snapshot,
138 PtrList<GeometricField<Type, PatchField, GeoMesh>>& modes, label Nmodes,
139 bool consider_volumes)
140{
141 label Msize;
142
143 if (Nmodes == 0)
144 {
145 Msize = modes.size();
146 }
147 else
148 {
149 Msize = Nmodes;
150 }
151
152 M_Assert(modes.size() >= Msize,
153 "The Number of requested modes is larger then the available quantity.");
154 Eigen::MatrixXd F = Foam2Eigen::PtrList2Eigen(modes);
155 Eigen::MatrixXd M_matrix = getMassMatrix(modes, Nmodes, consider_volumes);
156 Eigen::MatrixXd snapEigen = Foam2Eigen::field2Eigen(snapshot);
157 Eigen::VectorXd a(Msize);
158 Eigen::VectorXd b(Msize);
159
160 if (consider_volumes)
161 {
162 Eigen::VectorXd V = getMassMatrixFV(modes[0]);
163 b = F.transpose().topRows(Msize) * V.asDiagonal() * snapEigen;
164 }
165 else
166 {
167 b = F.transpose().topRows(Msize) * snapEigen;
168 }
169
170 if (Pstream::parRun())
171 {
172 reduce(b, sumOp<Eigen::VectorXd>());
173 }
174
175 a = M_matrix.fullPivLu().solve(b);
176 return a;
177}
178
179template Eigen::VectorXd getCoeffs(
180 GeometricField<scalar, fvPatchField, volMesh>&
181 snapshot, PtrList<GeometricField<scalar, fvPatchField, volMesh>>& modes,
182 label Nmodes,
183 bool consider_volumes);
184
185template Eigen::VectorXd getCoeffs(
186 GeometricField<scalar, fvsPatchField, surfaceMesh>&
187 snapshot, PtrList<GeometricField<scalar, fvsPatchField, surfaceMesh>>& modes,
188 label Nmodes,
189 bool consider_volumes = false);
190
191template Eigen::VectorXd getCoeffs(
192 GeometricField<vector, fvPatchField, volMesh>&
193 snapshot, PtrList<GeometricField<vector, fvPatchField, volMesh>>& modes,
194 label Nmodes,
195 bool consider_volumes);
196
197template<class Type, template<class> class PatchField, class GeoMesh>
198Eigen::MatrixXd getCoeffs(PtrList<GeometricField<Type, PatchField, GeoMesh>>&
199 snapshots,
200 PtrList<GeometricField<Type, PatchField, GeoMesh>>& modes, label Nmodes,
201 bool consider_volumes)
202{
203 label Msize;
204
205 if (Nmodes == 0)
206 {
207 Msize = modes.size();
208 }
209 else
210 {
211 Msize = Nmodes;
212 }
213
214 M_Assert(modes.size() >= Msize,
215 "The Number of requested modes is larger then the available quantity.");
216 Eigen::MatrixXd coeff(Msize, snapshots.size());
217
218 for (auto i = 0; i < snapshots.size(); i++)
219 {
220 coeff.col(i) = getCoeffs(snapshots[i], modes, Nmodes, consider_volumes);
221 }
222
223 return coeff;
224}
225
226template Eigen::MatrixXd getCoeffs(
227 PtrList<GeometricField<scalar, fvPatchField, volMesh>>&
228 snapshot, PtrList<GeometricField<scalar, fvPatchField, volMesh>>& modes,
229 label Nmodes,
230 bool consider_volumes);
231
232template Eigen::MatrixXd getCoeffs(
233 PtrList<GeometricField<vector, fvPatchField, volMesh>>&
234 snapshot, PtrList<GeometricField<vector, fvPatchField, volMesh>>& modes,
235 label Nmodes,
236 bool consider_volumes);
237
238template Eigen::MatrixXd getCoeffs(
239 PtrList<GeometricField<scalar, fvsPatchField, surfaceMesh>>&
240 snapshot, PtrList<GeometricField<scalar, fvsPatchField, surfaceMesh>>& modes,
241 label Nmodes,
242 bool consider_volumes);
243
244Eigen::MatrixXd parTimeCombMat(List<Eigen::VectorXd>
245 acquiredSnapshotsTimes,
246 Eigen::MatrixXd parameters)
247{
248 label parsNum = parameters.cols();
249 label parsSamplesNum = parameters.rows();
250 M_Assert(parsSamplesNum == acquiredSnapshotsTimes.size(),
251 "The list of time instants does not have the same number of vectors as the number of parameters samples");
252 Eigen::MatrixXd comb;
253 label totalSnapshotsNum = 0;
254
255 for (label k = 0; k < acquiredSnapshotsTimes.size(); k++)
256 {
257 totalSnapshotsNum += acquiredSnapshotsTimes[k].size();
258 }
259
260 comb.resize(totalSnapshotsNum, parsNum + 1);
261 label i = 0;
262
263 for (label j = 0; j < acquiredSnapshotsTimes.size(); j++)
264 {
265 for (label k = 0; k < acquiredSnapshotsTimes[j].size(); k++)
266 {
267 comb(i, parsNum) = (acquiredSnapshotsTimes[j])(k, 0);
268 comb.block(i, 0, 1, parsNum) = parameters.row(j);
269 i = i + 1;
270 }
271 }
272
273 return comb;
274}
275
276}
#define M_Assert(Expr, Msg)
Header file of the ITHACAcoeffsMass file.
static Eigen::VectorXd field2Eigen(GeometricField< tensor, PatchField, GeoMesh > &field)
Convert a vector OpenFOAM field into an Eigen Vector.
Definition Foam2Eigen.C:40
static Eigen::MatrixXd PtrList2Eigen(PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, label Nfields=-1)
Convert a PtrList of snapshots to Eigen matrix (only internal field)
Definition Foam2Eigen.C:638
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.
label i
Definition pEqn.H:46