fvFieldDecomposerTemplates.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2021-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU 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 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "fvFieldDecomposer.H"
30 #include "processorFvPatchField.H"
31 #include "processorFvsPatchField.H"
34 #include "emptyFvPatchFields.H"
35 #include "volFields.H"
36 #include "surfaceFields.H"
37 
38 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39 
40 template<class Type>
43 (
45 ) const
46 {
47  // Create the field for the processor
48 
50  (
51  field.name(),
53  procMesh_,
54  field.dimensions(),
55  // Internal field - mapped values
56  Field<Type>(field.field(), cellAddressing_)
57  );
58 }
59 
60 
61 template<class Type>
64 (
66  const bool allowUnknownPatchFields
67 ) const
68 {
69  // Create the field for the processor
70  // - with dummy patch fields
72  (
73  field.name(),
75  procMesh_,
76  field.dimensions(),
77  // Internal field - mapped values
78  Field<Type>(field.primitiveField(), cellAddressing_),
80  );
81  auto& result = tresult.ref();
82  result.oriented() = field.oriented();
83 
84 
85  // 2. Change the fvPatchFields to the correct type using a mapper
86  // constructor (with reference to the now correct internal field)
87 
88  auto& bf = result.boundaryFieldRef();
89 
90  forAll(bf, patchi)
91  {
92  if (patchFieldDecomposerPtrs_.set(patchi))
93  {
94  bf.set
95  (
96  patchi,
98  (
99  field.boundaryField()[boundaryAddressing_[patchi]],
100  procMesh_.boundary()[patchi],
101  result.internalField(),
102  patchFieldDecomposerPtrs_[patchi]
103  )
104  );
105  }
106  else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
107  {
108  bf.set
109  (
110  patchi,
111  new processorCyclicFvPatchField<Type>
112  (
113  procMesh_.boundary()[patchi],
114  result.internalField(),
115  Field<Type>
116  (
117  field.primitiveField(),
118  processorVolPatchFieldDecomposerPtrs_[patchi]
119  )
120  )
121  );
122  }
123  else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
124  {
125  bf.set
126  (
127  patchi,
128  new processorFvPatchField<Type>
129  (
130  procMesh_.boundary()[patchi],
131  result.internalField(),
132  Field<Type>
133  (
134  field.primitiveField(),
135  processorVolPatchFieldDecomposerPtrs_[patchi]
136  )
137  )
138  );
139  }
140  else if (allowUnknownPatchFields)
141  {
142  bf.set
143  (
144  patchi,
145  new emptyFvPatchField<Type>
146  (
147  procMesh_.boundary()[patchi],
148  result.internalField()
149  )
150  );
151  }
152  else
153  {
155  << "Unknown type." << abort(FatalError);
156  }
157  }
158 
159  // Create the field for the processor
160  return tresult;
161 }
162 
163 
164 template<class Type>
167 (
169 ) const
170 {
171  labelList mapAddr
172  (
173  labelList::subList(faceAddressing_, procMesh_.nInternalFaces())
174  );
175  forAll(mapAddr, i)
176  {
177  mapAddr[i] -= 1;
178  }
179 
180 
181  // Problem with addressing when a processor patch picks up both internal
182  // faces and faces from cyclic boundaries. This is a bit of a hack, but
183  // I cannot find a better solution without making the internal storage
184  // mechanism for surfaceFields correspond to the one of faces in polyMesh
185  // (i.e. using slices)
186  Field<Type> allFaceField(field.mesh().nFaces());
187 
188  {
189  SubList<Type>(allFaceField, field.primitiveField().size()) =
190  field.primitiveField();
191 
192  forAll(field.boundaryField(), patchi)
193  {
194  const Field<Type>& pfld = field.boundaryField()[patchi];
195 
196  const label start = field.mesh().boundaryMesh()[patchi].start();
197 
198  SubList<Type>(allFaceField, pfld.size(), start) = pfld;
199  }
200  }
201 
202 
203  // 1. Create the complete field with dummy patch fields
204 
206  (
207  field.name(),
209  procMesh_,
210  field.dimensions(),
211  // Internal field - mapped values
212  Field<Type>(field.primitiveField(), mapAddr),
214  );
215  auto& result = tresult.ref();
216  result.oriented() = field.oriented();
217 
218  // 2. Change the fvsPatchFields to the correct type using a mapper
219  // constructor (with reference to the now correct internal field)
220 
221  auto& bf = result.boundaryFieldRef();
222 
223  forAll(boundaryAddressing_, patchi)
224  {
225  if (patchFieldDecomposerPtrs_.set(patchi))
226  {
227  bf.set
228  (
229  patchi,
231  (
232  field.boundaryField()[boundaryAddressing_[patchi]],
233  procMesh_.boundary()[patchi],
234  result.internalField(),
235  patchFieldDecomposerPtrs_[patchi]
236  )
237  );
238  }
239  else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
240  {
241  bf.set
242  (
243  patchi,
244  new processorCyclicFvsPatchField<Type>
245  (
246  procMesh_.boundary()[patchi],
247  result.internalField(),
248  Field<Type>
249  (
250  allFaceField,
251  processorSurfacePatchFieldDecomposerPtrs_[patchi]
252  )
253  )
254  );
255 
256  if (result.is_oriented())
257  {
258  bf[patchi] *= faceSign_[patchi];
259  }
260  }
261  else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
262  {
263  bf.set
264  (
265  patchi,
266  new processorFvsPatchField<Type>
267  (
268  procMesh_.boundary()[patchi],
269  result.internalField(),
270  Field<Type>
271  (
272  allFaceField,
273  processorSurfacePatchFieldDecomposerPtrs_[patchi]
274  )
275  )
276  );
277 
278  if (result.is_oriented())
279  {
280  bf[patchi] *= faceSign_[patchi];
281  }
282  }
283  else
284  {
286  << "Unknown type." << abort(FatalError);
287  }
288  }
289 
290  // Create the field for the processor
291  return tresult;
292 }
293 
294 
295 template<class GeoField>
297 (
299 ) const
300 {
301  for (const auto& fld : fields)
302  {
303  decomposeField(fld)().write();
304  }
305 }
306 
307 
308 // ************************************************************************* //
Foam::surfaceFields.
rDeltaTY field()
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
orientedType oriented() const noexcept
Return oriented type.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Generic GeometricField class.
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
tmp< DimensionedField< Type, volMesh > > decomposeField(const DimensionedField< Type, volMesh > &field) const
Decompose internal field.
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
A List obtained as a section of another List.
Definition: SubList.H:50
Generic templated field type.
Definition: Field.H:62
errorManip< error > abort(error &err)
Definition: errorManip.H:139
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=PatchField< Type >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const Field< Type > &iField)
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions, copy of internal field...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:204
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void decomposeFields(const PtrList< GeoField > &fields) const
Decompose list of fields.
Do not request registration (bool: false)
static tmp< fvsPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, surfaceMesh > &)
Return a pointer to a new patchField created on freestore given.