SlicedGeometricField.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 "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  FieldField<PatchField, Type>& bf = tbf.ref();
53 
54  forAll(mesh.boundary(), patchi)
55  {
56  if
57  (
58  preserveCouples
59  && mesh.boundary()[patchi].coupled()
60  && (
61  !preserveProcessorOnly
62  || isA<processorFvPatch>(mesh.boundary()[patchi])
63  )
64  )
65  {
66  // For coupled patched construct the correct patch field type
67  bf.set
68  (
69  patchi,
71  (
72  mesh.boundary()[patchi].type(),
73  mesh.boundary()[patchi],
74  *this
75  )
76  );
77 
78  // Initialise the values on the coupled patch to those of the slice
79  // of the given field.
80  // Note: these will usually be over-ridden by the boundary field
81  // evaluation e.g. in the case of processor and cyclic patches.
82  bf[patchi] = SlicedPatchField<Type>
83  (
84  mesh.boundary()[patchi],
85  DimensionedField<Type, GeoMesh>::null(),
86  completeField
87  );
88  }
89  else
90  {
91  bf.set
92  (
93  patchi,
94  new SlicedPatchField<Type>
95  (
96  mesh.boundary()[patchi],
97  DimensionedField<Type, GeoMesh>::null(),
98  completeField
99  )
100  );
101  }
102  }
103 
104  return tbf;
105 }
106 
107 
108 template
109 <
110  class Type,
111  template<class> class PatchField,
112  template<class> class SlicedPatchField,
113  class GeoMesh
114 >
118 (
119  const Mesh& mesh,
120  const FieldField<PatchField, Type>& bField,
121  const bool preserveCouples
122 )
123 {
124  tmp<FieldField<PatchField, Type>> tbf
125  (
126  new FieldField<PatchField, Type>(mesh.boundary().size())
127  );
128  FieldField<PatchField, Type>& bf = tbf.ref();
129 
130  forAll(mesh.boundary(), patchi)
131  {
132  if (preserveCouples && mesh.boundary()[patchi].coupled())
133  {
134  // For coupled patched construct the correct patch field type
135  bf.set
136  (
137  patchi,
139  (
140  mesh.boundary()[patchi].type(),
141  mesh.boundary()[patchi],
142  *this
143  )
144  );
145 
146  // Assign field
147  bf[patchi] == bField[patchi];
148  }
149  else
150  {
151  // Create unallocated copy of patch field
152  bf.set
153  (
154  patchi,
155  new SlicedPatchField<Type>
156  (
157  mesh.boundary()[patchi],
158  DimensionedField<Type, GeoMesh>::null(),
159  bField[patchi]
160  )
161  );
162  }
163  }
164 
165  return tbf;
166 }
167 
168 
169 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
170 
171 template
172 <
173  class Type,
174  template<class> class PatchField,
175  template<class> class SlicedPatchField,
176  class GeoMesh
177 >
180 (
181  const IOobject& io,
182  const Mesh& mesh,
183  const dimensionSet& ds,
184  const Field<Type>& completeField,
185  const bool preserveCouples
186 )
187 :
188  GeometricField<Type, PatchField, GeoMesh>
189  (
190  io,
191  mesh,
192  ds,
193  Field<Type>(),
194  slicedBoundaryField(mesh, completeField, preserveCouples)
195  )
196 {
197  // Set the internalField to the slice of the complete field
199  (
200  typename Field<Type>::subField(completeField, GeoMesh::size(mesh))
201  );
202 
204 }
205 
206 
207 template
208 <
209  class Type,
210  template<class> class PatchField,
211  template<class> class SlicedPatchField,
212  class GeoMesh
213 >
216 (
217  const IOobject& io,
218  const Mesh& mesh,
219  const dimensionSet& ds,
220  const Field<Type>& completeIField,
221  const Field<Type>& completeBField,
222  const bool preserveCouples,
223  const bool preserveProcessorOnly
224 )
225 :
226  GeometricField<Type, PatchField, GeoMesh>
227  (
228  io,
229  mesh,
230  ds,
231  Field<Type>(),
232  slicedBoundaryField
233  (
234  mesh,
235  completeBField,
236  preserveCouples,
237  preserveProcessorOnly
238  )
239  )
240 {
241  // Set the internalField to the slice of the complete field
243  (
244  typename Field<Type>::subField(completeIField, GeoMesh::size(mesh))
245  );
246 
248 }
249 
250 
251 template
252 <
253  class Type,
254  template<class> class PatchField,
255  template<class> class SlicedPatchField,
256  class GeoMesh
257 >
260 (
261  const IOobject& io,
263  const bool preserveCouples
264 )
265 :
266  GeometricField<Type, PatchField, GeoMesh>
267  (
268  io,
269  gf.mesh(),
270  gf.dimensions(),
271  Field<Type>(),
272  slicedBoundaryField(gf.mesh(), gf.boundaryField(), preserveCouples)
273  )
274 {
275  // Set the internalField to the supplied internal field
277 
279 }
280 
281 
282 template
283 <
284  class Type,
285  template<class> class PatchField,
286  template<class> class SlicedPatchField,
287  class GeoMesh
288 >
291 (
293 )
294 :
295  GeometricField<Type, PatchField, GeoMesh>
296  (
297  gf,
298  gf.mesh(),
299  gf.dimensions(),
300  Field<Type>(),
301  slicedBoundaryField(gf.mesh(), gf.boundaryField(), true)
302  )
303 {
304  // Set the internalField to the supplied internal field
306 }
307 
308 
309 template
310 <
311  class Type,
312  template<class> class PatchField,
313  template<class> class SlicedPatchField,
314  class GeoMesh
315 >
316 Foam::tmp
317 <
319 >
321 clone() const
322 {
323  return tmp
324  <
326  >
327  (
329  (
330  *this
331  )
332  );
333 }
334 
335 
336 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
337 
338 template
339 <
340  class Type,
341  template<class> class PatchField,
342  template<class> class SlicedPatchField,
343  class GeoMesh
344 >
347 {
348  // Set the internalField storage pointer to nullptr before its destruction
349  // to protect the field it a slice of.
351 }
352 
353 
354 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
355 
356 template
357 <
358  class Type,
359  template<class> class PatchField,
360  template<class> class SlicedPatchField,
361  class GeoMesh
362 >
365 splice() const
366 {
367  const Mesh& mesh = this->mesh();
368 
369  label completeSize = GeoMesh::size(mesh);
370 
371  forAll(mesh.boundaryMesh(), patchi)
372  {
373  completeSize += mesh.boundaryMesh()[patchi].size();
374  }
375 
376  tmp<Field<Type>> tCompleteField(new Field<Type>(completeSize));
377  Field<Type>& completeField(tCompleteField.ref());
378 
379  typename Field<Type>::subField(completeField, GeoMesh::size(mesh))
380  = this->primitiveField();
381 
382  label start = GeoMesh::size(mesh);
383 
384  forAll(mesh.boundaryMesh(), patchi)
385  {
386  if
387  (
388  mesh.boundary()[patchi].size()
389  == mesh.boundaryMesh()[patchi].size()
390  )
391  {
392  typename Field<Type>::subField
393  (
394  completeField,
395  mesh.boundary()[patchi].size(),
396  start
397  ) = this->boundaryField()[patchi];
398  }
399  else
400  {
401  typename Field<Type>::subField
402  (
403  completeField,
404  mesh.boundaryMesh()[patchi].size(),
405  start
406  ) = Zero;
407  }
408 
409  start += mesh.boundaryMesh()[patchi].size();
410  }
411 
412  return tCompleteField;
413 }
414 
415 
416 template
417 <
418  class Type,
419  template<class> class PatchField,
420  template<class> class SlicedPatchField,
421  class GeoMesh
422 >
425 {
427 }
428 
429 
430 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
GeoMesh::Mesh Mesh
Type of mesh on which this DimensionedField is instantiated.
const Mesh & mesh() const
Return mesh.
Pre-declare SubField and related Field type.
Definition: Field.H:82
SubField< Type > subField
Declare type of subField.
Definition: Field.H:100
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
Generic GeometricField class.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
void correctBoundaryConditions()
Correct boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Specialisation of GeometricField which holds slices of given complete fields in a form that they act ...
tmp< SlicedGeometricField< Type, PatchField, SlicedPatchField, GeoMesh > > clone() const
Clone.
SlicedGeometricField(const IOobject &, const Mesh &, const dimensionSet &, const Field< Type > &completeField, const bool preserveCouples=true)
Construct from components and field to slice.
tmp< Field< Type > > splice() const
Splice the sliced field and return the complete field.
void correctBoundaryConditions()
Correct boundary field.
Pre-declare related SubField type.
Definition: SubField.H:63
void shallowCopy(const UList< T > &)
Copy the pointer held by the given UList.
Definition: UListI.H:156
Dimension set for the base types.
Definition: dimensionSet.H:122
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
label patchi
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static const zero Zero
Definition: zero.H:97
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