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 -------------------------------------------------------------------------------
10 License
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"
30 #include "fixedValueFvPatchField.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 template<class Type>
35 void 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 auto& i : indices)
47  {
48  lstWithin[j] = lst[i];
49  ++j;
50  }
51  lst.transfer(lstWithin);
52 }
53 
54 
55 template<class MatrixType>
56 void 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 auto& i : indices)
68  {
69  matWithin.subColumn(j) = mat.subColumn(i);
70  ++j;
71  }
72  mat.transfer(matWithin);
73 }
74 
75 
76 template<class Type>
77 bool Foam::DMDModels::STDMD::modes()
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 
95 template<class GeoFieldType>
96 bool 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 i : magsi_)
113  {
114  GeoFieldType modeRe
115  (
116  IOobject
117  (
118  "modeRe_" + std::to_string(modei)
119  + "_" + fieldName_ + "_" + name_,
120  mesh_.time().timeName(),
121  mesh_,
124  false
125  ),
126  mesh_,
127  dimensioned<Type>(dimless, Zero),
128  fixedValueFvPatchField<Type>::typeName
129  );
130 
131  GeoFieldType modeIm
132  (
133  IOobject
134  (
135  "modeIm_" + std::to_string(modei)
136  + "_" + fieldName_ + "_" + name_,
137  mesh_.time().timeName(),
138  mesh_,
141  false
142  ),
143  mesh_,
144  dimensioned<Type>(dimless, Zero),
145  fixedValueFvPatchField<Type>::typeName
146  );
147 
148  if (modeRe.size() != 0 && !empty_)
149  {
150  if (patch_.empty())
151  {
152  auto& inModeRe = modeRe.primitiveFieldRef();
153  auto& inModeIm = modeIm.primitiveFieldRef();
154 
155  calcMode(inModeRe, inModeIm, primitiveMode, i);
156  }
157  else
158  {
159  const label patchi = mesh_.boundaryMesh().findPatchID(patch_);
160 
161  auto& bfModeRe = modeRe.boundaryFieldRef()[patchi];
162  auto& bfModeIm = modeIm.boundaryFieldRef()[patchi];
163 
164  calcMode(bfModeRe, bfModeIm, primitiveMode, i);
165  }
166  }
167 
168  modeRe.write();
169  modeIm.write();
170  ++modei;
171  }
172 
173  return true;
174 }
175 
176 
177 template<class GeoFieldType>
178 typename std::enable_if
179 <
180  std::is_same<Foam::scalar, typename GeoFieldType::value_type>::value,
181  void
182 >::type Foam::DMDModels::STDMD::calcMode
183 (
184  GeoFieldType& modeRe,
185  GeoFieldType& modeIm,
186  const RMatrix& primitiveMode,
187  const label i
188 )
189 {
190  const label fieldSize = modeRe.size();
191 
192  for (label p = 0; p < primitiveMode.m(); ++p)
193  {
194  complex mode(Zero);
195  for (label q = 0; q < evecs_.m(); ++q)
196  {
197  mode += primitiveMode(p, q)*evecs_(q, i);
198  }
199  label p1 = p%fieldSize;
200  modeRe[p1] = mode.real();
201  modeIm[p1] = mode.imag();
202  }
203 }
204 
205 
206 template<class GeoFieldType>
207 typename std::enable_if
208 <
209  !std::is_same<Foam::scalar, typename GeoFieldType::value_type>::value,
210  void
211 >::type Foam::DMDModels::STDMD::calcMode
212 (
213  GeoFieldType& modeRe,
214  GeoFieldType& modeIm,
215  const RMatrix& primitiveMode,
216  const label i
217 )
218 {
219  const label fieldSize = modeRe.size();
220 
221  for (label p = 0; p < primitiveMode.m(); ++p)
222  {
223  complex mode(Zero);
224  for (label q = 0; q < evecs_.m(); ++q)
225  {
226  mode += primitiveMode(p, q)*evecs_(q, i);
227  }
228  label p1 = p%fieldSize;
229  label p2 = p/fieldSize;
230  modeRe[p1][p2] = mode.real();
231  modeIm[p1][p2] = mode.imag();
232  }
233 }
234 
235 
236 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
volFields.H
Foam::roots::complex
Definition: Roots.H:57
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
surfaceFields.H
Foam::surfaceFields.
Foam::labelPair
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:54
Foam::mode
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition: MSwindows.C:564
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
fixedValueFvPatchField.H
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189