fvFieldDecomposerDecomposeFields.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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "fvFieldDecomposer.H"
29 #include "processorFvPatchField.H"
30 #include "processorFvsPatchField.H"
33 #include "emptyFvPatchFields.H"
34 
35 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
36 
37 template<class Type>
40 (
42  const bool allowUnknownPatchFields
43 ) const
44 {
46 
47  // 1. Create the complete field with dummy patch fields
48  PtrList<fvPatchField<Type>> patchFields(boundaryAddressing_.size());
49 
50  forAll(boundaryAddressing_, patchi)
51  {
52  patchFields.set
53  (
54  patchi,
56  (
58  procMesh_.boundary()[patchi],
60  )
61  );
62  }
63 
64  // Create the field for the processor
65  tmp<VolFieldType> tresF
66  (
67  new VolFieldType
68  (
69  IOobject
70  (
71  field.name(),
72  procMesh_.time().timeName(),
73  procMesh_,
74  IOobject::NO_READ,
75  IOobject::NO_WRITE
76  ),
77  procMesh_,
78  field.dimensions(),
79  Field<Type>(field.primitiveField(), cellAddressing_),
80  patchFields
81  )
82  );
83  VolFieldType& resF = tresF.ref();
84  resF.oriented() = field().oriented();
85 
86 
87  // 2. Change the fvPatchFields to the correct type using a mapper
88  // constructor (with reference to the now correct internal field)
89 
90  typename VolFieldType::Boundary& bf = resF.boundaryFieldRef();
91 
92  forAll(bf, patchi)
93  {
94  if (patchFieldDecomposerPtrs_.set(patchi))
95  {
96  bf.set
97  (
98  patchi,
100  (
101  field.boundaryField()[boundaryAddressing_[patchi]],
102  procMesh_.boundary()[patchi],
103  resF(),
104  patchFieldDecomposerPtrs_[patchi]
105  )
106  );
107  }
108  else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
109  {
110  bf.set
111  (
112  patchi,
114  (
115  procMesh_.boundary()[patchi],
116  resF(),
118  (
119  field.primitiveField(),
120  processorVolPatchFieldDecomposerPtrs_[patchi]
121  )
122  )
123  );
124  }
125  else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
126  {
127  bf.set
128  (
129  patchi,
131  (
132  procMesh_.boundary()[patchi],
133  resF(),
135  (
136  field.primitiveField(),
137  processorVolPatchFieldDecomposerPtrs_[patchi]
138  )
139  )
140  );
141  }
142  else if (allowUnknownPatchFields)
143  {
144  bf.set
145  (
146  patchi,
148  (
149  procMesh_.boundary()[patchi],
150  resF()
151  )
152  );
153  }
154  else
155  {
157  << "Unknown type." << abort(FatalError);
158  }
159  }
160 
161  // Create the field for the processor
162  return tresF;
163 }
164 
165 
166 template<class Type>
169 (
171 ) const
172 {
173  typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfaceFieldType;
174 
175  labelList mapAddr
176  (
178  (
179  faceAddressing_,
180  procMesh_.nInternalFaces()
181  )
182  );
183  forAll(mapAddr, i)
184  {
185  mapAddr[i] -= 1;
186  }
187 
188  // Create and map the internal field values
189  Field<Type> internalField
190  (
191  field.primitiveField(),
192  mapAddr
193  );
194 
195  // Problem with addressing when a processor patch picks up both internal
196  // faces and faces from cyclic boundaries. This is a bit of a hack, but
197  // I cannot find a better solution without making the internal storage
198  // mechanism for surfaceFields correspond to the one of faces in polyMesh
199  // (i.e. using slices)
200  Field<Type> allFaceField(field.mesh().nFaces());
201 
202  forAll(field.primitiveField(), i)
203  {
204  allFaceField[i] = field.primitiveField()[i];
205  }
206 
207  forAll(field.boundaryField(), patchi)
208  {
209  const Field<Type>& p = field.boundaryField()[patchi];
210 
211  const label patchStart = field.mesh().boundaryMesh()[patchi].start();
212 
213  forAll(p, i)
214  {
215  allFaceField[patchStart + i] = p[i];
216  }
217  }
218 
219 
220  // 1. Create the complete field with dummy patch fields
221  PtrList<fvsPatchField<Type>> patchFields(boundaryAddressing_.size());
222 
223  forAll(boundaryAddressing_, patchi)
224  {
225  patchFields.set
226  (
227  patchi,
229  (
231  procMesh_.boundary()[patchi],
233  )
234  );
235  }
236 
238  (
239  new SurfaceFieldType
240  (
241  IOobject
242  (
243  field.name(),
244  procMesh_.time().timeName(),
245  procMesh_,
246  IOobject::NO_READ,
247  IOobject::NO_WRITE
248  ),
249  procMesh_,
250  field.dimensions(),
251  Field<Type>(field.primitiveField(), mapAddr),
252  patchFields
253  )
254  );
255  SurfaceFieldType& resF = tresF.ref();
256  resF.oriented() = field().oriented();
257 
258  // 2. Change the fvsPatchFields to the correct type using a mapper
259  // constructor (with reference to the now correct internal field)
260 
261  typename SurfaceFieldType::Boundary& bf = resF.boundaryFieldRef();
262 
263  forAll(boundaryAddressing_, patchi)
264  {
265  if (patchFieldDecomposerPtrs_.set(patchi))
266  {
267  bf.set
268  (
269  patchi,
271  (
272  field.boundaryField()[boundaryAddressing_[patchi]],
273  procMesh_.boundary()[patchi],
274  resF(),
275  patchFieldDecomposerPtrs_[patchi]
276  )
277  );
278  }
279  else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
280  {
281  bf.set
282  (
283  patchi,
285  (
286  procMesh_.boundary()[patchi],
287  resF(),
289  (
290  allFaceField,
291  processorSurfacePatchFieldDecomposerPtrs_[patchi]
292  )
293  )
294  );
295 
296  if (resF.oriented()())
297  {
298  bf[patchi] *= faceSign_[patchi];
299  }
300  }
301  else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
302  {
303  bf.set
304  (
305  patchi,
307  (
308  procMesh_.boundary()[patchi],
309  resF(),
311  (
312  allFaceField,
313  processorSurfacePatchFieldDecomposerPtrs_[patchi]
314  )
315  )
316  );
317 
318  if (resF.oriented()())
319  {
320  bf[patchi] *= faceSign_[patchi];
321  }
322  }
323  else
324  {
326  << "Unknown type." << abort(FatalError);
327  }
328  }
329 
330  // Create the field for the processor
331  return tresF;
332 }
333 
334 
335 template<class GeoField>
337 (
339 ) const
340 {
341  forAll(fields, fieldi)
342  {
343  decomposeField(fields[fieldi])().write();
344  }
345 }
346 
347 
348 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
Foam::fvFieldDecomposer::decomposeField
tmp< GeometricField< Type, fvPatchField, volMesh > > decomposeField(const GeometricField< Type, fvPatchField, volMesh > &field, const bool allowUnknownPatchFields=false) const
Decompose volume field.
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
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:53
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:228
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::fvFieldDecomposer::decomposeFields
void decomposeFields(const PtrList< GeoField > &fields) const
Definition: fvFieldDecomposerDecomposeFields.C:337
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:62
emptyFvPatchFields.H
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
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
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:35
Foam::processorFvsPatchField
Foam::processorFvsPatchField.
Definition: processorFvsPatchField.H:52
Foam::GeometricField< Type, fvPatchField, volMesh >
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
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