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 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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"
34#include "emptyFvPatchFields.H"
35#include "volFields.H"
36#include "surfaceFields.H"
37
38// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39
40template<class Type>
43(
45) const
46{
47 // Create and map the internal field values
48 Field<Type> mappedField(field, cellAddressing_);
49
50 // Create the field for the processor
51 return
53 (
55 (
56 field.name(),
57 procMesh_.thisDb().time().timeName(),
58 procMesh_.thisDb(),
61 false
62 ),
63 procMesh_,
64 field.dimensions(),
65 std::move(mappedField)
66 );
67}
68
69
70template<class Type>
73(
75 const bool allowUnknownPatchFields
76) const
77{
79
80 // 1. Create the complete field with dummy patch fields
81 PtrList<fvPatchField<Type>> patchFields(boundaryAddressing_.size());
82
83 forAll(boundaryAddressing_, patchi)
84 {
85 patchFields.set
86 (
87 patchi,
89 (
91 procMesh_.boundary()[patchi],
93 )
94 );
95 }
96
97 // Create the field for the processor
99 (
100 new VolFieldType
101 (
103 (
104 field.name(),
105 procMesh_.thisDb().time().timeName(),
106 procMesh_.thisDb(),
109 ),
110 procMesh_,
111 field.dimensions(),
112 Field<Type>(field.primitiveField(), cellAddressing_),
113 patchFields
114 )
115 );
116 VolFieldType& resF = tresF.ref();
117 resF.oriented() = field().oriented();
118
119
120 // 2. Change the fvPatchFields to the correct type using a mapper
121 // constructor (with reference to the now correct internal field)
122
123 auto& bf = resF.boundaryFieldRef();
124
125 forAll(bf, patchi)
126 {
127 if (patchFieldDecomposerPtrs_.set(patchi))
128 {
129 bf.set
130 (
131 patchi,
133 (
134 field.boundaryField()[boundaryAddressing_[patchi]],
135 procMesh_.boundary()[patchi],
136 resF(),
137 patchFieldDecomposerPtrs_[patchi]
138 )
139 );
140 }
141 else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
142 {
143 bf.set
144 (
145 patchi,
147 (
148 procMesh_.boundary()[patchi],
149 resF(),
151 (
152 field.primitiveField(),
153 processorVolPatchFieldDecomposerPtrs_[patchi]
154 )
155 )
156 );
157 }
158 else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
159 {
160 bf.set
161 (
162 patchi,
164 (
165 procMesh_.boundary()[patchi],
166 resF(),
168 (
169 field.primitiveField(),
170 processorVolPatchFieldDecomposerPtrs_[patchi]
171 )
172 )
173 );
174 }
175 else if (allowUnknownPatchFields)
176 {
177 bf.set
178 (
179 patchi,
181 (
182 procMesh_.boundary()[patchi],
183 resF()
184 )
185 );
186 }
187 else
188 {
190 << "Unknown type." << abort(FatalError);
191 }
192 }
193
194 // Create the field for the processor
195 return tresF;
196}
197
198
199template<class Type>
202(
204) const
205{
207
208 labelList mapAddr
209 (
211 (
212 faceAddressing_,
213 procMesh_.nInternalFaces()
214 )
215 );
216 forAll(mapAddr, i)
217 {
218 mapAddr[i] -= 1;
219 }
220
221 // Create and map the internal field values
222 Field<Type> internalField
223 (
224 field.primitiveField(),
225 mapAddr
226 );
227
228 // Problem with addressing when a processor patch picks up both internal
229 // faces and faces from cyclic boundaries. This is a bit of a hack, but
230 // I cannot find a better solution without making the internal storage
231 // mechanism for surfaceFields correspond to the one of faces in polyMesh
232 // (i.e. using slices)
233 Field<Type> allFaceField(field.mesh().nFaces());
234
235 forAll(field.primitiveField(), i)
236 {
237 allFaceField[i] = field.primitiveField()[i];
238 }
239
240 forAll(field.boundaryField(), patchi)
241 {
242 const Field<Type>& p = field.boundaryField()[patchi];
243
244 const label patchStart = field.mesh().boundaryMesh()[patchi].start();
245
246 forAll(p, i)
247 {
248 allFaceField[patchStart + i] = p[i];
249 }
250 }
251
252
253 // 1. Create the complete field with dummy patch fields
254 PtrList<fvsPatchField<Type>> patchFields(boundaryAddressing_.size());
255
256 forAll(boundaryAddressing_, patchi)
257 {
258 patchFields.set
259 (
260 patchi,
262 (
264 procMesh_.boundary()[patchi],
266 )
267 );
268 }
269
271 (
272 new SurfaceFieldType
273 (
275 (
276 field.name(),
277 procMesh_.thisDb().time().timeName(),
278 procMesh_.thisDb(),
281 ),
282 procMesh_,
283 field.dimensions(),
284 Field<Type>(field.primitiveField(), mapAddr),
285 patchFields
286 )
287 );
288 SurfaceFieldType& resF = tresF.ref();
289 resF.oriented() = field().oriented();
290
291 // 2. Change the fvsPatchFields to the correct type using a mapper
292 // constructor (with reference to the now correct internal field)
293
294 auto& bf = resF.boundaryFieldRef();
295
296 forAll(boundaryAddressing_, patchi)
297 {
298 if (patchFieldDecomposerPtrs_.set(patchi))
299 {
300 bf.set
301 (
302 patchi,
304 (
305 field.boundaryField()[boundaryAddressing_[patchi]],
306 procMesh_.boundary()[patchi],
307 resF(),
308 patchFieldDecomposerPtrs_[patchi]
309 )
310 );
311 }
312 else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
313 {
314 bf.set
315 (
316 patchi,
318 (
319 procMesh_.boundary()[patchi],
320 resF(),
322 (
323 allFaceField,
324 processorSurfacePatchFieldDecomposerPtrs_[patchi]
325 )
326 )
327 );
328
329 if (resF.oriented()())
330 {
331 bf[patchi] *= faceSign_[patchi];
332 }
333 }
334 else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
335 {
336 bf.set
337 (
338 patchi,
340 (
341 procMesh_.boundary()[patchi],
342 resF(),
344 (
345 allFaceField,
346 processorSurfacePatchFieldDecomposerPtrs_[patchi]
347 )
348 )
349 );
350
351 if (resF.oriented()())
352 {
353 bf[patchi] *= faceSign_[patchi];
354 }
355 }
356 else
357 {
359 << "Unknown type." << abort(FatalError);
360 }
361 }
362
363 // Create the field for the processor
364 return tresF;
365}
366
367
368template<class GeoField>
370(
372) const
373{
374 for (const auto& fld : fields)
375 {
376 decomposeField(fld)().write();
377 }
378}
379
380
381// ************************************************************************* //
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
A List obtained as a section of another List.
Definition: SubList.H:70
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Foam::calculatedFvsPatchField.
This boundary condition provides an 'empty' condition for reduced dimensions cases,...
tmp< DimensionedField< Type, volMesh > > decomposeField(const DimensionedField< Type, volMesh > &field) const
Decompose internal field.
void decomposeFields(const PtrList< GeoField > &fields) const
Decompose list of fields.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:302
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
const Time & time() const noexcept
Return time registry.
This boundary condition enables processor communication across cyclic patches.
Foam::processorCyclicFvsPatchField.
This boundary condition enables processor communication across patches.
Foam::processorFvsPatchField.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
volScalarField & p
rDeltaTY field()
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
errorManip< error > abort(error &err)
Definition: errorManip.H:144
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
runTime write()
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.