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-2019 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  if (patchFieldDecomposerPtrs_[patchi])
128  {
129  bf.set
130  (
131  patchi,
133  (
134  field.boundaryField()[boundaryAddressing_[patchi]],
135  procPatch,
136  resF(),
137  *patchFieldDecomposerPtrs_[patchi]
138  )
139  );
140  }
141  else if (isA<processorCyclicFvPatch>(procPatch))
142  {
143  bf.set
144  (
145  patchi,
147  (
148  procPatch,
149  resF(),
150  (*processorVolPatchFieldDecomposerPtrs_[patchi])
151  (
152  field.primitiveField()
153  )
154  )
155  );
156  }
157  else if (isA<processorFvPatch>(procPatch))
158  {
159  bf.set
160  (
161  patchi,
163  (
164  procPatch,
165  resF(),
166  (*processorVolPatchFieldDecomposerPtrs_[patchi])
167  (
168  field.primitiveField()
169  )
170  )
171  );
172  }
173  else if (allowUnknownPatchFields)
174  {
175  bf.set
176  (
177  patchi,
179  (
180  procPatch,
181  resF()
182  )
183  );
184  }
185  else
186  {
188  << "Unknown type." << abort(FatalError);
189  }
190  }
191 
192  // Create the field for the processor
193  return tresF;
194 }
195 
196 
197 template<class Type>
200 (
202 ) const
203 {
204  // Apply flipping to surfaceScalarFields only
205  const bool doFlip = (pTraits<Type>::nComponents == 1);
206 
207 
208  // Problem with addressing when a processor patch picks up both internal
209  // faces and faces from cyclic boundaries. This is a bit of a hack, but
210  // I cannot find a better solution without making the internal storage
211  // mechanism for surfaceFields correspond to the one of faces in polyMesh
212  // (i.e. using slices)
213  Field<Type> allFaceField(field.mesh().nFaces());
214 
215  forAll(field.primitiveField(), i)
216  {
217  allFaceField[i] = field.primitiveField()[i];
218  }
219 
220  forAll(field.boundaryField(), patchi)
221  {
222  const Field<Type> & p = field.boundaryField()[patchi];
223 
224  const label patchStart = field.mesh().boundaryMesh()[patchi].start();
225 
226  forAll(p, i)
227  {
228  allFaceField[patchStart + i] = p[i];
229  }
230  }
231 
232 
233  // 1. Create the complete field with dummy patch fields
234  PtrList<fvsPatchField<Type>> patchFields(boundaryAddressing_.size());
235 
236  forAll(boundaryAddressing_, patchi)
237  {
238  patchFields.set
239  (
240  patchi,
242  (
244  procMesh_.boundary()[patchi],
246  )
247  );
248  }
249 
251  (
253  (
254  IOobject
255  (
256  field.name(),
257  procMesh_.time().timeName(),
258  procMesh_,
259  IOobject::NO_READ,
260  IOobject::NO_WRITE,
261  false
262  ),
263  procMesh_,
264  field.dimensions(),
265  mapField
266  (
267  field,
269  (
270  faceAddressing_,
271  procMesh_.nInternalFaces()
272  ),
273  doFlip
274  ),
275  patchFields
276  )
277  );
279 
280 
281  // 2. Change the fvsPatchFields to the correct type using a mapper
282  // constructor (with reference to the now correct internal field)
283 
285  Boundary& bf = resF.boundaryFieldRef();
286 
287  forAll(boundaryAddressing_, patchi)
288  {
289  const fvPatch& procPatch = procMesh_.boundary()[patchi];
290 
291  if (patchFieldDecomposerPtrs_[patchi])
292  {
293  bf.set
294  (
295  patchi,
297  (
298  field.boundaryField()[boundaryAddressing_[patchi]],
299  procPatch,
300  resF(),
301  *patchFieldDecomposerPtrs_[patchi]
302  )
303  );
304  }
305  else
306  {
307  // Do our own mapping - avoids a lot of mapping complexity
308 
309  if (isA<processorCyclicFvPatch>(procPatch))
310  {
311  bf.set
312  (
313  patchi,
315  (
316  procPatch,
317  resF(),
318  mapField
319  (
320  allFaceField,
321  procPatch.patchSlice
322  (
323  faceAddressing_
324  ),
325  doFlip
326  )
327  )
328  );
329  }
330  else if (isA<processorFvPatch>(procPatch))
331  {
332  bf.set
333  (
334  patchi,
336  (
337  procPatch,
338  resF(),
339  mapField
340  (
341  allFaceField,
342  procPatch.patchSlice
343  (
344  faceAddressing_
345  ),
346  doFlip
347  )
348  )
349  );
350  }
351  else
352  {
354  << "Unknown type." << abort(FatalError);
355  }
356  }
357  }
358 
359  // Create the field for the processor
360  return tresF;
361 }
362 
363 
364 template<class GeoField>
366 (
367  const PtrList<GeoField>& fields
368 ) const
369 {
370  forAll(fields, fieldi)
371  {
372  decomposeField(fields[fieldi])().write();
373  }
374 }
375 
376 
377 // ************************************************************************* //
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:319
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:174
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
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
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:284
const Mesh & mesh() const
Return mesh.
const List< T >::subList patchSlice(const List< T > &l) const
Slice list to patch.
Definition: fvPatch.H:186
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:299
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