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-2023 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 #include "stringOps.H"
33 
34 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 
36 template<class Type>
37 Foam::tmp<Foam::Field<Type>> Foam::fvFieldDecomposer::mapCellToFace
38 (
39  const labelUList& owner,
40  const labelUList& neighbour,
41  const Field<Type>& field,
42  const labelUList& addressing
43 )
44 {
45  tmp<Field<Type>> tfld(new Field<Type>(addressing.size()));
46  Field<Type>& fld = tfld.ref();
47 
48  forAll(addressing, i)
49  {
50  fld[i] =
51  field
52  [
53  addressing[i] > 0
54  ? neighbour[addressing[i] - 1]
55  : owner[- addressing[i] - 1]
56  ];
57  }
58 
59  return tfld;
60 }
61 
62 
63 template<class Type>
64 Foam::tmp<Foam::Field<Type>> Foam::fvFieldDecomposer::mapFaceToFace
65 (
66  const Field<Type>& field,
67  const labelUList& addressing,
68  const bool isFlux
69 )
70 {
71  tmp<Field<Type>> tfld(new Field<Type>(addressing.size()));
72  Field<Type>& fld = tfld.ref();
73 
74  forAll(addressing, i)
75  {
76  fld[i] =
77  (isFlux && addressing[i] < 0 ? -1 : +1)
78  *field[mag(addressing[i]) - 1];
79  }
80 
81  return tfld;
82 }
83 
84 
85 template<class Type>
88 (
89  const VolField<Type>& field
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<VolField<Type>> tresF
110  (
111  new VolField<Type>
112  (
113  IOobject
114  (
115  field.name(),
116  procMesh_.time().name(),
117  procMesh_,
120  false
121  ),
122  procMesh_,
123  field.dimensions(),
124  Field<Type>(field.primitiveField(), cellAddressing_),
125  patchFields
126  )
127  );
128  VolField<Type>& 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 VolField<Type>::
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  {
158  OStringStream str;
159  str << "\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.";
169  << stringOps::breakIntoIndentedLines(str.str()).c_str()
170  << exit(FatalError);
171  }
172 
173  const label nbrCompletePatchi =
174  refCast<const processorCyclicFvPatch>(procPatch)
175  .referPatch().nbrPatchID();
176 
177  bf.set
178  (
179  procPatchi,
180  new processorCyclicFvPatchField<Type>
181  (
182  procPatch,
183  resF(),
184  mapCellToFace
185  (
186  labelUList(),
187  completeMesh_.lduAddr().patchAddr(nbrCompletePatchi),
188  field.primitiveField(),
189  faceAddressingBf_[procPatchi]
190  )
191  )
192  );
193  }
194  else if (isA<processorFvPatch>(procPatch))
195  {
196  bf.set
197  (
198  procPatchi,
199  new processorFvPatchField<Type>
200  (
201  procPatch,
202  resF(),
203  mapCellToFace
204  (
205  completeMesh_.owner(),
206  completeMesh_.neighbour(),
207  field.primitiveField(),
208  faceAddressingBf_[procPatchi]
209  )
210  )
211  );
212  }
213  else
214  {
216  << "Unknown type." << abort(FatalError);
217  }
218  }
219 
220  return tresF;
221 }
222 
223 
224 template<class Type>
227 (
228  const SurfaceField<Type>& field
229 ) const
230 {
231  const SubList<label> faceAddressingIf
232  (
233  faceAddressing_,
234  procMesh_.nInternalFaces()
235  );
236 
237  // Create dummy patch fields
238  PtrList<fvsPatchField<Type>> patchFields(procMesh_.boundary().size());
239  forAll(procMesh_.boundary(), procPatchi)
240  {
241  patchFields.set
242  (
243  procPatchi,
245  (
247  procMesh_.boundary()[procPatchi],
249  )
250  );
251  }
252 
253  // Create the processor field with the dummy patch fields
254  tmp<SurfaceField<Type>> tresF
255  (
256  new SurfaceField<Type>
257  (
258  IOobject
259  (
260  field.name(),
261  procMesh_.time().name(),
262  procMesh_,
265  false
266  ),
267  procMesh_,
268  field.dimensions(),
269  mapFaceToFace
270  (
271  field,
272  faceAddressingIf,
273  isFlux(field)
274  ),
275  patchFields
276  )
277  );
278  SurfaceField<Type>& resF = tresF.ref();
279 
280  // Change the patch fields to the correct type using a mapper constructor
281  // (with reference to the now correct internal field)
282  typename SurfaceField<Type>::
283  Boundary& bf = resF.boundaryFieldRef();
284  forAll(procMesh_.boundary(), procPatchi)
285  {
286  const fvPatch& procPatch = procMesh_.boundary()[procPatchi];
287 
288  const label completePatchi = completePatchID(procPatchi);
289 
290  if (completePatchi == procPatchi)
291  {
292  bf.set
293  (
294  procPatchi,
296  (
297  field.boundaryField()[procPatchi],
298  procPatch,
299  resF(),
300  patchFieldDecomposers_[procPatchi]
301  )
302  );
303  }
304  else if (isA<processorCyclicFvPatch>(procPatch))
305  {
306  bf.set
307  (
308  procPatchi,
309  new processorCyclicFvsPatchField<Type>
310  (
311  procPatch,
312  resF(),
313  mapFaceToFace
314  (
315  field.boundaryField()[completePatchi],
316  faceAddressingBf_[procPatchi],
317  isFlux(field)
318  )
319  )
320  );
321  }
322  else if (isA<processorFvPatch>(procPatch))
323  {
324  bf.set
325  (
326  procPatchi,
327  new processorFvsPatchField<Type>
328  (
329  procPatch,
330  resF(),
331  mapFaceToFace
332  (
333  field.primitiveField(),
334  faceAddressingBf_[procPatchi],
335  isFlux(field)
336  )
337  )
338  );
339  }
340  else
341  {
343  << "Unknown type." << abort(FatalError);
344  }
345  }
346 
347  return tresF;
348 }
349 
350 
351 template<class GeoField>
353 (
354  const PtrList<GeoField>& fields
355 ) const
356 {
357  forAll(fields, fieldi)
358  {
359  decomposeField(fields[fieldi])().write();
360  }
361 }
362 
363 
364 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
static const char *const typeName
Definition: Field.H:105
void decomposeFields(const PtrList< GeoField > &fields) const
Decompose a list of fields.
tmp< VolField< Type > > decomposeField(const VolField< Type > &field) const
Decompose volume field.
static tmp< fvPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
static tmp< fvsPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, surfaceMesh > &)
Return a pointer to a new patchField created on freestore given.
A class for managing temporary objects.
Definition: tmp.H:55
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
string breakIntoIndentedLines(const string &str, const string::size_type nLength=80, const string::size_type nIndent=0)
Break a string up into indented lines.
Definition: stringOps.C:956
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
errorManip< error > abort(error &err)
Definition: errorManip.H:131
bool isFlux(const DimensionedField< Type, surfaceMesh > &df)
Check if surfaceField is a flux.
Definition: surfaceFields.H:55
dimensioned< scalar > mag(const dimensioned< Type > &)
error FatalError
UList< label > labelUList
Definition: UList.H:65