SlicedGeometricField.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 | Copyright (C) 2011-2012 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 "SlicedGeometricField.H"
27 #include "processorFvPatch.H"
28 
29 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
30 
31 template
32 <
33  class Type,
34  template<class> class PatchField,
35  template<class> class SlicedPatchField,
36  class GeoMesh
37 >
41 (
42  const Mesh& mesh,
43  const Field<Type>& completeField,
44  const bool preserveCouples,
45  const bool preserveProcessorOnly
46 )
47 {
48  tmp<FieldField<PatchField, Type> > tbf
49  (
50  new FieldField<PatchField, Type>(mesh.boundary().size())
51  );
52 
53  FieldField<PatchField, Type>& bf = tbf();
54 
55  forAll(mesh.boundary(), patchi)
56  {
57  if
58  (
59  preserveCouples
60  && mesh.boundary()[patchi].coupled()
61  && (
62  !preserveProcessorOnly
63  || isA<processorFvPatch>(mesh.boundary()[patchi])
64  )
65  )
66  {
67  // For coupled patched construct the correct patch field type
68  bf.set
69  (
70  patchi,
72  (
73  mesh.boundary()[patchi].type(),
74  mesh.boundary()[patchi],
75  *this
76  )
77  );
78 
79  // Initialize the values on the coupled patch to those of the slice
80  // of the given field.
81  // Note: these will usually be over-ridden by the boundary field
82  // evaluation e.g. in the case of processor and cyclic patches.
83  bf[patchi] = SlicedPatchField<Type>
84  (
85  mesh.boundary()[patchi],
86  DimensionedField<Type, GeoMesh>::null(),
87  completeField
88  );
89  }
90  else
91  {
92  bf.set
93  (
94  patchi,
95  new SlicedPatchField<Type>
96  (
97  mesh.boundary()[patchi],
98  DimensionedField<Type, GeoMesh>::null(),
99  completeField
100  )
101  );
102  }
103  }
104 
105  return tbf;
106 }
107 
108 
109 template
110 <
111  class Type,
112  template<class> class PatchField,
113  template<class> class SlicedPatchField,
114  class GeoMesh
115 >
119 (
120  const Mesh& mesh,
121  const FieldField<PatchField, Type>& bField,
122  const bool preserveCouples
123 )
124 {
126  (
127  new FieldField<PatchField, Type>(mesh.boundary().size())
128  );
129 
130  FieldField<PatchField, Type>& bf = tbf();
131 
132  forAll(mesh.boundary(), patchi)
133  {
134  if (preserveCouples && mesh.boundary()[patchi].coupled())
135  {
136  // For coupled patched construct the correct patch field type
137  bf.set
138  (
139  patchi,
141  (
142  mesh.boundary()[patchi].type(),
143  mesh.boundary()[patchi],
144  *this
145  )
146  );
147 
148  // Assign field
149  bf[patchi] == bField[patchi];
150  }
151  else
152  {
153  // Create unallocated copy of patch field
154  bf.set
155  (
156  patchi,
157  new SlicedPatchField<Type>
158  (
159  mesh.boundary()[patchi],
161  )
162  );
163  bf[patchi].UList<Type>::operator=(bField[patchi]);
164  }
165  }
166 
167  return tbf;
168 }
169 
170 
171 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
172 
173 template
174 <
175  class Type,
176  template<class> class PatchField,
177  template<class> class SlicedPatchField,
178  class GeoMesh
179 >
182 (
183  const IOobject& io,
184  const Mesh& mesh,
185  const dimensionSet& ds,
186  const Field<Type>& iField
187 )
188 :
190  (
191  io,
192  mesh,
193  ds,
194  Field<Type>()
195  )
196 {
197  // Set the internalField to the slice of the complete field
199  (
200  typename Field<Type>::subField(iField, GeoMesh::size(mesh))
201  );
202 }
203 
204 
205 template
206 <
207  class Type,
208  template<class> class PatchField,
209  template<class> class SlicedPatchField,
210  class GeoMesh
211 >
214 (
215  const IOobject& io,
216  const Mesh& mesh,
217  const dimensionSet& ds,
218  const Field<Type>& completeField,
219  const bool preserveCouples
220 )
221 :
222  GeometricField<Type, PatchField, GeoMesh>
223  (
224  io,
225  mesh,
226  ds,
227  Field<Type>(),
228  slicedBoundaryField(mesh, completeField, preserveCouples)
229  )
230 {
231  // Set the internalField to the slice of the complete field
233  (
234  typename Field<Type>::subField(completeField, GeoMesh::size(mesh))
235  );
236 
238 }
239 
240 
241 template
242 <
243  class Type,
244  template<class> class PatchField,
245  template<class> class SlicedPatchField,
246  class GeoMesh
247 >
250 (
251  const IOobject& io,
252  const Mesh& mesh,
253  const dimensionSet& ds,
254  const Field<Type>& completeIField,
255  const Field<Type>& completeBField,
256  const bool preserveCouples,
257  const bool preserveProcessorOnly
258 )
259 :
260  GeometricField<Type, PatchField, GeoMesh>
261  (
262  io,
263  mesh,
264  ds,
265  Field<Type>(),
267  (
268  mesh,
269  completeBField,
270  preserveCouples,
271  preserveProcessorOnly
272  )
273  )
274 {
275  // Set the internalField to the slice of the complete field
277  (
278  typename Field<Type>::subField(completeIField, GeoMesh::size(mesh))
279  );
280 
282 }
283 
284 
285 template
286 <
287  class Type,
288  template<class> class PatchField,
289  template<class> class SlicedPatchField,
290  class GeoMesh
291 >
294 (
295  const IOobject& io,
296  const GeometricField<Type, PatchField, GeoMesh>& gf,
297  const bool preserveCouples
298 )
299 :
300  GeometricField<Type, PatchField, GeoMesh>
301  (
302  io,
303  gf.mesh(),
304  gf.dimensions(),
305  Field<Type>(),
306  slicedBoundaryField(gf.mesh(), gf.boundaryField(), preserveCouples)
307  )
308 {
309  // Set the internalField to the supplied internal field
310  UList<Type>::operator=(gf.internalField());
311 
313 }
314 
315 
316 template
317 <
318  class Type,
319  template<class> class PatchField,
320  template<class> class SlicedPatchField,
321  class GeoMesh
322 >
325 (
326  const SlicedGeometricField<Type, PatchField, SlicedPatchField, GeoMesh>& gf
327 )
328 :
329  GeometricField<Type, PatchField, GeoMesh>
330  (
331  gf,
332  gf.mesh(),
333  gf.dimensions(),
334  Field<Type>(),
335  slicedBoundaryField(gf.mesh(), gf.boundaryField(), true)
336  )
337 {
338  // Set the internalField to the supplied internal field
339  UList<Type>::operator=(gf.internalField());
340 }
341 
342 
343 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
344 
345 template
346 <
347  class Type,
348  template<class> class PatchField,
349  template<class> class SlicedPatchField,
350  class GeoMesh
351 >
354 {
355  // Set the internalField storage pointer to NULL before its destruction
356  // to protect the field it a slice of.
358 }
359 
360 
361 template
362 <
363  class Type,
364  template<class> class PatchField,
365  template<class> class SlicedPatchField,
366  class GeoMesh
367 >
370 {
371  // Set the internalField storage pointer to NULL before its destruction
372  // to protect the field it a slice of.
374 }
375 
376 
377 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
378 
379 template
380 <
381  class Type,
382  template<class> class PatchField,
383  template<class> class SlicedPatchField,
384  class GeoMesh
385 >
388 {
390 }
391 
392 
393 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
processorFvPatch.H
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::FieldField< PatchField, Type >
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
SlicedGeometricField.H
Foam::SlicedGeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Foam::SlicedGeometricField::slicedBoundaryField
tmp< FieldField< PatchField, Type > > slicedBoundaryField(const Mesh &mesh, const Field< Type > &completeField, const bool preserveCouples, const bool preserveProcessorOnly=false)
Slice the given field and a create a PtrList of SlicedPatchField.
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:116
Foam::SlicedGeometricField::SlicedGeometricField
SlicedGeometricField(const IOobject &, const Mesh &, const dimensionSet &, const Field< Type > &completeField, const bool preserveCouples=true)
Construct from components and field to slice.
Foam::SlicedGeometricField::DimensionedInternalField::DimensionedInternalField
DimensionedInternalField(const IOobject &, const Mesh &, const dimensionSet &, const Field< Type > &iField)
Construct from components and field to slice.
Definition: SlicedGeometricField.C:182
Foam::DimensionedField::Mesh
GeoMesh::Mesh Mesh
Definition: DimensionedField.H:81
Foam::Field< Type >
Foam::DimensionedField::mesh
const Mesh & mesh() const
Return mesh.
Definition: DimensionedFieldI.H:38
Foam::SlicedGeometricField::DimensionedInternalField::~DimensionedInternalField
~DimensionedInternalField()
Destructor.
Definition: SlicedGeometricField.C:369
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
boundaryField
cellIbMask *cellIbMaskExt *faceIbMask *cellIbMask boundaryField().evaluateCoupled()
Foam::GeoMesh
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
Foam::SlicedGeometricField::Mesh
GeoMesh::Mesh Mesh
Definition: SlicedGeometricField.H:69
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
Foam::UList
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
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::SlicedGeometricField::~SlicedGeometricField
~SlicedGeometricField()
Destructor.
Definition: SlicedGeometricField.C:353
Foam::Field::subField
SubField< Type > subField
Declare type of subField.
Definition: Field.H:89
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::UList::operator=
void operator=(const T &)
Assignment of all entries to the given value.
Definition: UList.C:70
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51