fvFieldDecomposerFields.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "fvFieldDecomposer.H"
30 #include "processorFvPatchField.H"
31 #include "processorFvsPatchField.H"
34 #include "emptyFvPatchFields.H"
35 
36 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
37 
38 template<class Type>
41 (
43 ) const
44 {
45  // Create and map the internal field values
46  Field<Type> mappedField(field, cellAddressing_);
47 
48  // Create the field for the processor
49  return
51  (
52  IOobject
53  (
54  field.name(),
55  procMesh_.time().timeName(),
56  procMesh_,
57  IOobject::NO_READ,
58  IOobject::NO_WRITE,
59  false
60  ),
61  procMesh_,
62  field.dimensions(),
63  std::move(mappedField)
64  );
65 }
66 
67 
68 template<class Type>
71 (
73  const bool allowUnknownPatchFields
74 ) const
75 {
77 
78  // 1. Create the complete field with dummy patch fields
79  PtrList<fvPatchField<Type>> patchFields(boundaryAddressing_.size());
80 
81  forAll(boundaryAddressing_, patchi)
82  {
83  patchFields.set
84  (
85  patchi,
87  (
89  procMesh_.boundary()[patchi],
91  )
92  );
93  }
94 
95  // Create the field for the processor
96  tmp<VolFieldType> tresF
97  (
98  new VolFieldType
99  (
100  IOobject
101  (
102  field.name(),
103  procMesh_.time().timeName(),
104  procMesh_,
105  IOobject::NO_READ,
106  IOobject::NO_WRITE
107  ),
108  procMesh_,
109  field.dimensions(),
110  Field<Type>(field.primitiveField(), cellAddressing_),
111  patchFields
112  )
113  );
114  VolFieldType& resF = tresF.ref();
115  resF.oriented() = field().oriented();
116 
117 
118  // 2. Change the fvPatchFields to the correct type using a mapper
119  // constructor (with reference to the now correct internal field)
120 
121  typename VolFieldType::Boundary& bf = resF.boundaryFieldRef();
122 
123  forAll(bf, patchi)
124  {
125  if (patchFieldDecomposerPtrs_.set(patchi))
126  {
127  bf.set
128  (
129  patchi,
131  (
132  field.boundaryField()[boundaryAddressing_[patchi]],
133  procMesh_.boundary()[patchi],
134  resF(),
135  patchFieldDecomposerPtrs_[patchi]
136  )
137  );
138  }
139  else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
140  {
141  bf.set
142  (
143  patchi,
145  (
146  procMesh_.boundary()[patchi],
147  resF(),
149  (
150  field.primitiveField(),
151  processorVolPatchFieldDecomposerPtrs_[patchi]
152  )
153  )
154  );
155  }
156  else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
157  {
158  bf.set
159  (
160  patchi,
162  (
163  procMesh_.boundary()[patchi],
164  resF(),
166  (
167  field.primitiveField(),
168  processorVolPatchFieldDecomposerPtrs_[patchi]
169  )
170  )
171  );
172  }
173  else if (allowUnknownPatchFields)
174  {
175  bf.set
176  (
177  patchi,
179  (
180  procMesh_.boundary()[patchi],
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  typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfaceFieldType;
205 
206  labelList mapAddr
207  (
209  (
210  faceAddressing_,
211  procMesh_.nInternalFaces()
212  )
213  );
214  forAll(mapAddr, i)
215  {
216  mapAddr[i] -= 1;
217  }
218 
219  // Create and map the internal field values
220  Field<Type> internalField
221  (
222  field.primitiveField(),
223  mapAddr
224  );
225 
226  // Problem with addressing when a processor patch picks up both internal
227  // faces and faces from cyclic boundaries. This is a bit of a hack, but
228  // I cannot find a better solution without making the internal storage
229  // mechanism for surfaceFields correspond to the one of faces in polyMesh
230  // (i.e. using slices)
231  Field<Type> allFaceField(field.mesh().nFaces());
232 
233  forAll(field.primitiveField(), i)
234  {
235  allFaceField[i] = field.primitiveField()[i];
236  }
237 
238  forAll(field.boundaryField(), patchi)
239  {
240  const Field<Type>& p = field.boundaryField()[patchi];
241 
242  const label patchStart = field.mesh().boundaryMesh()[patchi].start();
243 
244  forAll(p, i)
245  {
246  allFaceField[patchStart + i] = p[i];
247  }
248  }
249 
250 
251  // 1. Create the complete field with dummy patch fields
252  PtrList<fvsPatchField<Type>> patchFields(boundaryAddressing_.size());
253 
254  forAll(boundaryAddressing_, patchi)
255  {
256  patchFields.set
257  (
258  patchi,
260  (
262  procMesh_.boundary()[patchi],
264  )
265  );
266  }
267 
269  (
270  new SurfaceFieldType
271  (
272  IOobject
273  (
274  field.name(),
275  procMesh_.time().timeName(),
276  procMesh_,
277  IOobject::NO_READ,
278  IOobject::NO_WRITE
279  ),
280  procMesh_,
281  field.dimensions(),
282  Field<Type>(field.primitiveField(), mapAddr),
283  patchFields
284  )
285  );
286  SurfaceFieldType& resF = tresF.ref();
287  resF.oriented() = field().oriented();
288 
289  // 2. Change the fvsPatchFields to the correct type using a mapper
290  // constructor (with reference to the now correct internal field)
291 
292  typename SurfaceFieldType::Boundary& bf = resF.boundaryFieldRef();
293 
294  forAll(boundaryAddressing_, patchi)
295  {
296  if (patchFieldDecomposerPtrs_.set(patchi))
297  {
298  bf.set
299  (
300  patchi,
302  (
303  field.boundaryField()[boundaryAddressing_[patchi]],
304  procMesh_.boundary()[patchi],
305  resF(),
306  patchFieldDecomposerPtrs_[patchi]
307  )
308  );
309  }
310  else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
311  {
312  bf.set
313  (
314  patchi,
316  (
317  procMesh_.boundary()[patchi],
318  resF(),
320  (
321  allFaceField,
322  processorSurfacePatchFieldDecomposerPtrs_[patchi]
323  )
324  )
325  );
326 
327  if (resF.oriented()())
328  {
329  bf[patchi] *= faceSign_[patchi];
330  }
331  }
332  else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
333  {
334  bf.set
335  (
336  patchi,
338  (
339  procMesh_.boundary()[patchi],
340  resF(),
342  (
343  allFaceField,
344  processorSurfacePatchFieldDecomposerPtrs_[patchi]
345  )
346  )
347  );
348 
349  if (resF.oriented()())
350  {
351  bf[patchi] *= faceSign_[patchi];
352  }
353  }
354  else
355  {
357  << "Unknown type." << abort(FatalError);
358  }
359  }
360 
361  // Create the field for the processor
362  return tresF;
363 }
364 
365 
366 template<class GeoField>
368 (
370 ) const
371 {
372  for (const auto& fld : fields)
373  {
374  decomposeField(fld)().write();
375  }
376 }
377 
378 
379 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
processorFvPatchField.H
Foam::processorFvPatchField
This boundary condition enables processor communication across patches.
Definition: processorFvPatchField.H:66
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:54
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
Foam::calculatedFvsPatchField
Foam::calculatedFvsPatchField.
Definition: calculatedFvsPatchField.H:51
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
processorFvsPatchField.H
Foam::processorCyclicFvPatchField
This boundary condition enables processor communication across cyclic patches.
Definition: processorCyclicFvPatchField.H:71
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::fvFieldDecomposer::decomposeFields
void decomposeFields(const PtrList< GeoField > &fields) const
Decompose list of fields.
Definition: fvFieldDecomposerFields.C:368
processorCyclicFvPatchField.H
field
rDeltaTY field()
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
emptyFvPatchFields.H
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::emptyFvPatchField
This boundary condition provides an 'empty' condition for reduced dimensions cases,...
Definition: emptyFvPatchField.H:67
Foam::FatalError
error FatalError
Foam::processorCyclicFvsPatchField
Foam::processorCyclicFvsPatchField.
Definition: processorCyclicFvsPatchField.H:52
Foam::calculatedFvPatchField
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Definition: calculatedFvPatchField.H:66
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
fvFieldDecomposer.H
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::List< label >
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::processorFvsPatchField
Foam::processorFvsPatchField.
Definition: processorFvsPatchField.H:52
Foam::GeometricField< Type, fvPatchField, volMesh >
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::fvFieldDecomposer::decomposeField
tmp< DimensionedField< Type, volMesh > > decomposeField(const DimensionedField< Type, volMesh > &field) const
Decompose internal field.
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
processorCyclicFvsPatchField.H