fvFieldDecomposerDecomposeFields.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fvFieldDecomposer.H"
27 #include "processorFvPatchField.H"
28 #include "processorFvsPatchField.H"
31 #include "emptyFvPatchFields.H"
32 
33 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
34 
35 template<class Type>
36 Foam::tmp<Foam::Field<Type>> Foam::fvFieldDecomposer::mapField
37 (
38  const Field<Type>& field,
39  const labelUList& mapAndSign,
40  const bool applyFlip
41 )
42 {
43  tmp<Field<Type>> tfld(new Field<Type>(mapAndSign.size()));
44  Field<Type>& fld = tfld.ref();
45 
46  if (applyFlip)
47  {
48  forAll(mapAndSign, i)
49  {
50  if (mapAndSign[i] < 0)
51  {
52  fld[i] = -field[-mapAndSign[i] - 1];
53  }
54  else
55  {
56  fld[i] = field[mapAndSign[i] - 1];
57  }
58  }
59  }
60  else
61  {
62  // Ignore face flipping
63  fld.map(field, mag(mapAndSign) - 1);
64  }
65  return tfld;
66 }
67 
68 
69 template<class Type>
72 (
74  const bool allowUnknownPatchFields
75 ) const
76 {
77  // 1. Create the complete field with dummy patch fields
78  PtrList<fvPatchField<Type>> patchFields(boundaryAddressing_.size());
79 
80  forAll(boundaryAddressing_, patchi)
81  {
82  patchFields.set
83  (
84  patchi,
86  (
88  procMesh_.boundary()[patchi],
90  )
91  );
92  }
93 
94  // Create the field for the processor
96  (
98  (
99  IOobject
100  (
101  field.name(),
102  procMesh_.time().timeName(),
103  procMesh_,
104  IOobject::NO_READ,
105  IOobject::NO_WRITE,
106  false
107  ),
108  procMesh_,
109  field.dimensions(),
110  Field<Type>(field.primitiveField(), cellAddressing_),
111  patchFields
112  )
113  );
115 
116 
117  // 2. Change the fvPatchFields to the correct type using a mapper
118  // constructor (with reference to the now correct internal field)
119 
121  Boundary& bf = resF.boundaryFieldRef();
122 
123  forAll(bf, patchi)
124  {
125  const fvPatch& procPatch = procMesh_.boundary()[patchi];
126 
127  label fromPatchi = boundaryAddressing_[patchi];
128  if (fromPatchi < 0 && isA<processorCyclicFvPatch>(procPatch))
129  {
130  const label referPatchi =
131  refCast<const processorCyclicPolyPatch>
132  (procPatch.patch()).referPatchID();
133  if (field.boundaryField()[referPatchi].overridesConstraint())
134  {
135  fromPatchi = boundaryAddressing_[referPatchi];
136  }
137  }
138 
139  if (fromPatchi >= 0)
140  {
141  bf.set
142  (
143  patchi,
145  (
146  field.boundaryField()[fromPatchi],
147  procPatch,
148  resF(),
149  *patchFieldDecomposerPtrs_[patchi]
150  )
151  );
152  }
153  else if (isA<processorCyclicFvPatch>(procPatch))
154  {
155  bf.set
156  (
157  patchi,
159  (
160  procPatch,
161  resF(),
162  (*processorVolPatchFieldDecomposerPtrs_[patchi])
163  (
164  field.primitiveField()
165  )
166  )
167  );
168  }
169  else if (isA<processorFvPatch>(procPatch))
170  {
171  bf.set
172  (
173  patchi,
175  (
176  procPatch,
177  resF(),
178  (*processorVolPatchFieldDecomposerPtrs_[patchi])
179  (
180  field.primitiveField()
181  )
182  )
183  );
184  }
185  else if (allowUnknownPatchFields)
186  {
187  bf.set
188  (
189  patchi,
191  (
192  procPatch,
193  resF()
194  )
195  );
196  }
197  else
198  {
200  << "Unknown type." << abort(FatalError);
201  }
202  }
203 
204  // Create the field for the processor
205  return tresF;
206 }
207 
208 
209 template<class Type>
212 (
214 ) const
215 {
216  // Apply flipping to surfaceScalarFields only
217  const bool doFlip = (pTraits<Type>::nComponents == 1);
218 
219 
220  // Problem with addressing when a processor patch picks up both internal
221  // faces and faces from cyclic boundaries. This is a bit of a hack, but
222  // I cannot find a better solution without making the internal storage
223  // mechanism for surfaceFields correspond to the one of faces in polyMesh
224  // (i.e. using slices)
225  Field<Type> allFaceField(field.mesh().nFaces());
226 
227  forAll(field.primitiveField(), i)
228  {
229  allFaceField[i] = field.primitiveField()[i];
230  }
231 
232  forAll(field.boundaryField(), patchi)
233  {
234  const Field<Type> & p = field.boundaryField()[patchi];
235 
236  const label patchStart = field.mesh().boundaryMesh()[patchi].start();
237 
238  forAll(p, i)
239  {
240  allFaceField[patchStart + i] = p[i];
241  }
242  }
243 
244 
245  // 1. Create the complete field with dummy patch fields
246  PtrList<fvsPatchField<Type>> patchFields(boundaryAddressing_.size());
247 
248  forAll(boundaryAddressing_, patchi)
249  {
250  patchFields.set
251  (
252  patchi,
254  (
256  procMesh_.boundary()[patchi],
258  )
259  );
260  }
261 
263  (
265  (
266  IOobject
267  (
268  field.name(),
269  procMesh_.time().timeName(),
270  procMesh_,
271  IOobject::NO_READ,
272  IOobject::NO_WRITE,
273  false
274  ),
275  procMesh_,
276  field.dimensions(),
277  mapField
278  (
279  field,
281  (
282  faceAddressing_,
283  procMesh_.nInternalFaces()
284  ),
285  doFlip
286  ),
287  patchFields
288  )
289  );
291 
292 
293  // 2. Change the fvsPatchFields to the correct type using a mapper
294  // constructor (with reference to the now correct internal field)
295 
297  Boundary& bf = resF.boundaryFieldRef();
298 
299  forAll(boundaryAddressing_, patchi)
300  {
301  const fvPatch& procPatch = procMesh_.boundary()[patchi];
302 
303  if (boundaryAddressing_[patchi] >= 0)
304  {
305  bf.set
306  (
307  patchi,
309  (
310  field.boundaryField()[boundaryAddressing_[patchi]],
311  procPatch,
312  resF(),
313  *patchFieldDecomposerPtrs_[patchi]
314  )
315  );
316  }
317  else if (isA<processorCyclicFvPatch>(procPatch))
318  {
319  // Do our own mapping. Avoids a lot of mapping complexity.
320  bf.set
321  (
322  patchi,
324  (
325  procPatch,
326  resF(),
327  mapField
328  (
329  allFaceField,
330  procPatch.patchSlice(faceAddressing_),
331  doFlip
332  )
333  )
334  );
335  }
336  else if (isA<processorFvPatch>(procPatch))
337  {
338  // Do our own mapping. Avoids a lot of mapping complexity.
339  bf.set
340  (
341  patchi,
343  (
344  procPatch,
345  resF(),
346  mapField
347  (
348  allFaceField,
349  procPatch.patchSlice(faceAddressing_),
350  doFlip
351  )
352  )
353  );
354  }
355  else
356  {
358  << "Unknown type." << abort(FatalError);
359  }
360  }
361 
362  // Create the field for the processor
363  return tresF;
364 }
365 
366 
367 template<class GeoField>
369 (
370  const PtrList<GeoField>& fields
371 ) const
372 {
373  forAll(fields, fieldi)
374  {
375  decomposeField(fields[fieldi])().write();
376  }
377 }
378 
379 
380 // ************************************************************************* //
Foam::processorCyclicFvsPatchField.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
Definition: IOobject.H:303
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
Traits class for primitives.
Definition: pTraits.H:50
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Generic GeometricField class.
tmp< GeometricField< Type, fvPatchField, volMesh > > decomposeField(const GeometricField< Type, fvPatchField, volMesh > &field, const bool allowUnknownPatchFields=false) const
Decompose volume field.
const dimensionSet & dimensions() const
Return dimensions.
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::processorFvsPatchField.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Pre-declare SubField and related Field type.
Definition: Field.H:56
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:138
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition: Field.C:265
const Mesh & mesh() const
Return mesh.
const List< T >::subList patchSlice(const List< T > &l) const
Slice list to patch.
Definition: fvPatch.H:190
This boundary condition enables processor communication across cyclic patches.
Internal & ref()
Return a reference to the dimensioned internal field.
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
This boundary condition enables processor communication across patches.
This boundary condition provides an &#39;empty&#39; condition for reduced dimensions cases, i.e. 1- and 2-D geometries. Apply this condition to patches whose normal is aligned to geometric directions that do not constitue solution directions.
Foam::calculatedFvsPatchField.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
void decomposeFields(const PtrList< GeoField > &fields) const
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65