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-2022 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::mapCellToFace
37 (
38  const labelUList& owner,
39  const labelUList& neighbour,
40  const Field<Type>& field,
41  const labelUList& addressing
42 )
43 {
44  tmp<Field<Type>> tfld(new Field<Type>(addressing.size()));
45  Field<Type>& fld = tfld.ref();
46 
47  forAll(addressing, i)
48  {
49  fld[i] =
50  field
51  [
52  addressing[i] > 0
53  ? neighbour[addressing[i] - 1]
54  : owner[- addressing[i] + 1]
55  ];
56  }
57 
58  return tfld;
59 }
60 
61 
62 template<class Type>
63 Foam::tmp<Foam::Field<Type>> Foam::fvFieldDecomposer::mapFaceToFace
64 (
65  const Field<Type>& field,
66  const labelUList& addressing,
67  const bool isFlux
68 )
69 {
70  tmp<Field<Type>> tfld(new Field<Type>(addressing.size()));
71  Field<Type>& fld = tfld.ref();
72 
73  forAll(addressing, i)
74  {
75  fld[i] =
76  (isFlux && addressing[i] < 0 ? -1 : +1)
77  *field[mag(addressing[i]) - 1];
78  }
79 
80  return tfld;
81 }
82 
83 
84 template<class Type>
87 (
88  const GeometricField<Type, fvPatchField, volMesh>& field,
89  const bool allowUnknownPatchFields
90 ) const
91 {
92  // Create dummy patch fields
93  PtrList<fvPatchField<Type>> patchFields(procMesh_.boundary().size());
94  forAll(procMesh_.boundary(), procPatchi)
95  {
96  patchFields.set
97  (
98  procPatchi,
100  (
102  procMesh_.boundary()[procPatchi],
104  )
105  );
106  }
107 
108  // Create the processor field with the dummy patch fields
109  tmp<GeometricField<Type, fvPatchField, volMesh>> tresF
110  (
111  new GeometricField<Type, fvPatchField, volMesh>
112  (
113  IOobject
114  (
115  field.name(),
116  procMesh_.time().timeName(),
117  procMesh_,
120  false
121  ),
122  procMesh_,
123  field.dimensions(),
124  Field<Type>(field.primitiveField(), cellAddressing_),
125  patchFields
126  )
127  );
128  GeometricField<Type, fvPatchField, volMesh>& resF = tresF.ref();
129 
130  // Change the patch fields to the correct type using a mapper constructor
131  // (with reference to the now correct internal field)
132  typename GeometricField<Type, fvPatchField, volMesh>::
133  Boundary& bf = resF.boundaryFieldRef();
134  forAll(bf, procPatchi)
135  {
136  const fvPatch& procPatch = procMesh_.boundary()[procPatchi];
137 
138  const label completePatchi = completePatchID(procPatchi);
139 
140  if (completePatchi == procPatchi)
141  {
142  bf.set
143  (
144  procPatchi,
146  (
147  field.boundaryField()[completePatchi],
148  procPatch,
149  resF(),
150  patchFieldDecomposers_[procPatchi]
151  )
152  );
153  }
154  else if (isA<processorCyclicFvPatch>(procPatch))
155  {
156  if (field.boundaryField()[completePatchi].overridesConstraint())
157  {
159  << "\nThe field \"" << field.name()
160  << "\" on cyclic patch \""
161  << field.boundaryField()[completePatchi].patch().name()
162  << "\" cannot be decomposed as it is not a cyclic "
163  << "patch field. A \"patchType cyclic;\" setting has "
164  << "been used to override the cyclic patch type.\n\n"
165  << "Cyclic patches like this with non-cyclic boundary "
166  << "conditions should be confined to a single "
167  << "processor using decomposition constraints."
168  << exit(FatalError);
169  }
170 
171  const label nbrCompletePatchi =
172  refCast<const processorCyclicFvPatch>(procPatch)
173  .referPatch().nbrPatchID();
174 
175  bf.set
176  (
177  procPatchi,
178  new processorCyclicFvPatchField<Type>
179  (
180  procPatch,
181  resF(),
182  mapCellToFace
183  (
184  labelUList(),
185  completeMesh_.lduAddr().patchAddr(nbrCompletePatchi),
186  field.primitiveField(),
187  faceAddressingBf_[procPatchi]
188  )
189  )
190  );
191  }
192  else if (isA<processorFvPatch>(procPatch))
193  {
194  bf.set
195  (
196  procPatchi,
197  new processorFvPatchField<Type>
198  (
199  procPatch,
200  resF(),
201  mapCellToFace
202  (
203  completeMesh_.owner(),
204  completeMesh_.neighbour(),
205  field.primitiveField(),
206  faceAddressingBf_[procPatchi]
207  )
208  )
209  );
210  }
211  else if (allowUnknownPatchFields)
212  {
213  bf.set
214  (
215  procPatchi,
216  new emptyFvPatchField<Type>
217  (
218  procPatch,
219  resF()
220  )
221  );
222  }
223  else
224  {
226  << "Unknown type." << abort(FatalError);
227  }
228  }
229 
230  return tresF;
231 }
232 
233 
234 template<class Type>
237 (
238  const GeometricField<Type, fvsPatchField, surfaceMesh>& field
239 ) const
240 {
241  const SubList<label> faceAddressingIf
242  (
243  faceAddressing_,
244  procMesh_.nInternalFaces()
245  );
246 
247  // Create dummy patch fields
248  PtrList<fvsPatchField<Type>> patchFields(procMesh_.boundary().size());
249  forAll(procMesh_.boundary(), procPatchi)
250  {
251  patchFields.set
252  (
253  procPatchi,
255  (
257  procMesh_.boundary()[procPatchi],
259  )
260  );
261  }
262 
263  // Create the processor field with the dummy patch fields
264  tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tresF
265  (
266  new GeometricField<Type, fvsPatchField, surfaceMesh>
267  (
268  IOobject
269  (
270  field.name(),
271  procMesh_.time().timeName(),
272  procMesh_,
275  false
276  ),
277  procMesh_,
278  field.dimensions(),
279  mapFaceToFace
280  (
281  field,
282  faceAddressingIf,
283  isFlux(field)
284  ),
285  patchFields
286  )
287  );
288  GeometricField<Type, fvsPatchField, surfaceMesh>& resF = tresF.ref();
289 
290  // Change the patch fields to the correct type using a mapper constructor
291  // (with reference to the now correct internal field)
292  typename GeometricField<Type, fvsPatchField, surfaceMesh>::
293  Boundary& bf = resF.boundaryFieldRef();
294  forAll(procMesh_.boundary(), procPatchi)
295  {
296  const fvPatch& procPatch = procMesh_.boundary()[procPatchi];
297 
298  const label completePatchi = completePatchID(procPatchi);
299 
300  if (completePatchi == procPatchi)
301  {
302  bf.set
303  (
304  procPatchi,
306  (
307  field.boundaryField()[procPatchi],
308  procPatch,
309  resF(),
310  patchFieldDecomposers_[procPatchi]
311  )
312  );
313  }
314  else if (isA<processorCyclicFvPatch>(procPatch))
315  {
316  bf.set
317  (
318  procPatchi,
319  new processorCyclicFvsPatchField<Type>
320  (
321  procPatch,
322  resF(),
323  mapFaceToFace
324  (
325  field.boundaryField()[completePatchi],
326  faceAddressingBf_[procPatchi],
327  isFlux(field)
328  )
329  )
330  );
331  }
332  else if (isA<processorFvPatch>(procPatch))
333  {
334  bf.set
335  (
336  procPatchi,
337  new processorFvsPatchField<Type>
338  (
339  procPatch,
340  resF(),
341  mapFaceToFace
342  (
343  field.primitiveField(),
344  faceAddressingBf_[procPatchi],
345  isFlux(field)
346  )
347  )
348  );
349  }
350  else
351  {
353  << "Unknown type." << abort(FatalError);
354  }
355  }
356 
357  return tresF;
358 }
359 
360 
361 template<class GeoField>
363 (
364  const PtrList<GeoField>& fields
365 ) const
366 {
367  forAll(fields, fieldi)
368  {
369  decomposeField(fields[fieldi])().write();
370  }
371 }
372 
373 
374 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label nInternalFaces() const
static const char *const typeName
Definition: Field.H:105
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
UList< label > labelUList
Definition: UList.H:65
tmp< GeometricField< Type, fvPatchField, volMesh > > decomposeField(const GeometricField< Type, fvPatchField, volMesh > &field, const bool allowUnknownPatchFields=false) const
Decompose volume field.
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:423
static tmp< fvsPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, surfaceMesh > &)
Return a pointer to a new patchField created on freestore given.
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:806
virtual const labelUList & patchAddr(const label patchNo) const =0
Return patch to internal addressing given patch number.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
static tmp< fvPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:417
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
bool isFlux(const DimensionedField< Type, surfaceMesh > &df)
Check if surfaceField is a flux.
Definition: surfaceFields.H:55
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
void decomposeFields(const PtrList< GeoField > &fields) const
Decompose a list of fields.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:800