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