Loading...
Searching...
No Matches
STDMDTemplates.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2020-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "volFields.H"
29#include "surfaceFields.H"
31
32// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33
34template<class Type>
35void Foam::DMDModels::STDMD::filterIndexed
36(
37 List<Type>& lst,
38 const UList<label>& indices
39)
40{
41 // Elements within [a, b]
42 List<Type> lstWithin(indices.size());
43
44 // Copy if frequency of element is within [a, b]
45 label j = 0;
46 for (const label i : indices)
47 {
48 lstWithin[j] = lst[i];
49 ++j;
50 }
51 lst.transfer(lstWithin);
52}
53
54
55template<class MatrixType>
56void Foam::DMDModels::STDMD::filterIndexed
57(
58 MatrixType& mat,
59 const UList<label>& indices
60)
61{
62 // Elements within [a, b]
63 MatrixType matWithin(labelPair(mat.m(), indices.size()));
64
65 // Copy if frequency of element is within [a, b]
66 label j = 0;
67 for (const label i : indices)
68 {
69 matWithin.subColumn(j) = mat.subColumn(i);
70 ++j;
71 }
72 mat.transfer(matWithin);
73}
74
75
76template<class Type>
78{
79 typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
80 typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfaceFieldType;
81
82 if (mesh_.foundObject<VolFieldType>(fieldName_))
83 {
84 return calcModes<VolFieldType>();
85 }
86 else if (mesh_.foundObject<SurfaceFieldType>(fieldName_))
87 {
88 return calcModes<SurfaceFieldType>();
89 }
90
91 return false;
92}
93
94
95template<class GeoFieldType>
96bool Foam::DMDModels::STDMD::calcModes()
97{
98 typedef typename GeoFieldType::value_type Type;
99
100 // Resize the number of output variables to "nModes" if requested
101 if (magsi_.size() > nModes_)
102 {
103 magsi_.resize(nModes_);
104 }
105
106 // Compute and write modes one by one
107 const RMatrix primitiveMode(Qupper_*RxInv_);
108 Qupper_.clear();
109 RxInv_.clear();
110
111 label modei = 0;
112 for (const label magi : magsi_)
113 {
114 GeoFieldType modeRe
115 (
116 IOobject
117 (
118 "modeRe_" + std::to_string(modei)
119 + "_" + fieldName_ + "_" + name_,
120 mesh_.time().timeName(),
121 mesh_,
125 ),
126 mesh_,
127 dimensioned<Type>(dimless, Zero),
129 );
130
131 GeoFieldType modeIm
132 (
133 IOobject
134 (
135 "modeIm_" + std::to_string(modei)
136 + "_" + fieldName_ + "_" + name_,
137 mesh_.time().timeName(),
138 mesh_,
142 ),
143 mesh_,
144 dimensioned<Type>(dimless, Zero),
146 );
147
148 if (modeRe.size() != 0 && !empty_)
149 {
150 if (patches_.empty())
151 {
152 auto& re = modeRe.primitiveFieldRef();
153 auto& im = modeIm.primitiveFieldRef();
154
155 calcMode(re, im, primitiveMode, magi);
156 }
157 else
158 {
159 label rowi = 0;
160 const labelList patchis
161 (
162 mesh_.boundaryMesh().patchSet(patches_).sortedToc()
163 );
164
165 for (const label patchi : patchis)
166 {
167 auto& re = modeRe.boundaryFieldRef()[patchi];
168 auto& im = modeIm.boundaryFieldRef()[patchi];
169
170 calcMode(re, im, primitiveMode, magi, rowi);
171
172 rowi += re.size()*pTraits<Type>::nComponents;
173 }
174 }
175 }
176
177 modeRe.write();
178 modeIm.write();
179 ++modei;
180 }
181
182 return true;
183}
184
185
186template<class GeoFieldType>
187void Foam::DMDModels::STDMD::calcMode
188(
189 GeoFieldType& modeRe,
190 GeoFieldType& modeIm,
191 const RMatrix& primitiveMode,
192 const label magi,
193 const label rowi
194)
195{
196 const label szfld = modeRe.size();
197 const label szfldcmps =
198 szfld*pTraits<typename GeoFieldType::value_type>::nComponents;
199
200 for (label i = rowi; i < szfldcmps + rowi; ++i)
201 {
203 for (label j = 0; j < evecs_.m(); ++j)
204 {
205 mode += primitiveMode(i, j)*evecs_(j, magi);
206 }
207
208 const label k = (i-rowi)%szfld;
209
210 if constexpr
211 (
212 std::is_same_v<scalar, typename GeoFieldType::value_type>
213 )
214 {
215 modeRe[k] = mode.real();
216 modeIm[k] = mode.imag();
217 }
218 else
219 {
220 const label m = (i-rowi)/szfld;
221 modeRe[k][m] = mode.real();
222 modeIm[k][m] = mode.imag();
223 }
224 }
225}
226
227
228// ************************************************************************* //
label k
virtual bool modes()=0
Compute and write modes.
static const char *const typeName
Typename for Field.
Definition Field.H:93
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
const dimensionedScalar re
Classical electron radius: default SI units: [m].
@ complex
Definition Roots.H:55
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Definition List.H:62
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition POSIX.C:775
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Foam::surfaceFields.