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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2022-2023 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 "SlicedGeometricField.H"
30 
31 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
32 
33 template
34 <
35  class Type,
36  template<class> class PatchField,
37  template<class> class SlicedPatchField,
38  class GeoMesh
39 >
40 bool
43 (
44  const Mesh& mesh,
45  const label fieldSize
46 )
47 {
48  label maxAddress(0);
49 
50  if (!mesh.boundary().empty())
51  {
52  const auto& p = mesh.boundary().back();
53  maxAddress = (p.start() + p.size());
54  }
55 
56  // If field size appear to not include internal field
57  return (fieldSize < maxAddress);
58 }
59 
60 
61 template
62 <
63  class Type,
64  template<class> class PatchField,
65  template<class> class SlicedPatchField,
66  class GeoMesh
67 >
71 (
72  const Mesh& mesh,
73  const Field<Type>& completeOrBoundaryField,
74  const bool preserveCouples,
75  const bool preserveProcessorOnly,
76  const bool isBoundaryOnly
77 ) const
78 {
79  typedef typename
80  SlicedPatchField<Type>::processorPatchType
81  processorPatchType;
82 
83  auto tbf = tmp<FieldField<PatchField, Type>>::New(mesh.boundary().size());
84  auto& bf = tbf.ref();
85 
86  forAll(mesh.boundary(), patchi)
87  {
88  const auto& p = mesh.boundary()[patchi];
89 
90  if
91  (
92  preserveCouples && p.coupled()
93  && (!preserveProcessorOnly || isA<processorPatchType>(p))
94  )
95  {
96  // For coupled patched construct the correct patch field type
97  bf.set
98  (
99  patchi,
100  PatchField<Type>::New(p.type(), p, *this)
101  );
102 
103  // Initialize the values on the coupled patch to those of the slice
104  // of the given field.
105  // Note: these will usually be over-ridden by the boundary field
106  // evaluation e.g. in the case of processor and cyclic patches.
107  bf[patchi] = SlicedPatchField<Type>
108  (
109  p,
110  DimensionedField<Type, GeoMesh>::null(),
111  completeOrBoundaryField,
112  isBoundaryOnly
113  );
114  }
115  else
116  {
117  bf.set
118  (
119  patchi,
120  new SlicedPatchField<Type>
121  (
122  p,
123  DimensionedField<Type, GeoMesh>::null(),
124  completeOrBoundaryField,
125  isBoundaryOnly
126  )
127  );
128  }
129  }
130 
131  return tbf;
132 }
133 
134 
135 template
136 <
137  class Type,
138  template<class> class PatchField,
139  template<class> class SlicedPatchField,
140  class GeoMesh
141 >
145 (
146  const Mesh& mesh,
147  const FieldField<PatchField, Type>& bField,
148  const bool preserveCouples
149 ) const
150 {
151  auto tbf = tmp<FieldField<PatchField, Type>>::New(mesh.boundary().size());
152  auto& bf = tbf.ref();
153 
154  forAll(mesh.boundary(), patchi)
155  {
156  const auto& p = mesh.boundary()[patchi];
157 
158  if (preserveCouples && p.coupled())
159  {
160  // For coupled patched construct the correct patch field type
161  bf.set
162  (
163  patchi,
164  PatchField<Type>::New(p.type(), p, *this)
165  );
166 
167  // Assign field
168  bf[patchi] == bField[patchi];
169  }
170  else
171  {
172  // Create unallocated copy of patch field
173  bf.set
174  (
175  patchi,
176  new SlicedPatchField<Type>
177  (
178  p,
179  DimensionedField<Type, GeoMesh>::null()
180  )
181  );
182 
183  bf[patchi].UList<Type>::shallowCopy(bField[patchi]);
184  }
185  }
186 
187  return tbf;
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
192 
193 template
194 <
195  class Type,
196  template<class> class PatchField,
197  template<class> class SlicedPatchField,
198  class GeoMesh
199 >
202 (
203  const IOobject& io,
204  const Mesh& mesh,
205  const dimensionSet& dims,
206  const Field<Type>& completeField,
207  const bool preserveCouples
208 )
209 :
210  GeometricField<Type, PatchField, GeoMesh>
211  (
212  io,
213  mesh,
214  dims,
215  Field<Type>(),
216  // preserveProcessorOnly = false
217  // isBoundaryOnly = false
218  makeBoundary(mesh, completeField, preserveCouples)
219  )
220 {
221  // Set internalField to the slice of the complete field
223  (
224  SubList<Type>(completeField, GeoMesh::size(mesh))
225  );
226 
228 }
229 
230 
231 template
232 <
233  class Type,
234  template<class> class PatchField,
235  template<class> class SlicedPatchField,
236  class GeoMesh
237 >
240 (
241  const IOobject& io,
242  const Mesh& mesh,
243  const dimensionSet& dims,
244  const Field<Type>& completeIField,
245  const Field<Type>& completeBField,
246  const bool preserveCouples,
247  const bool preserveProcessorOnly
248 )
249 :
250  GeometricField<Type, PatchField, GeoMesh>
251  (
252  io,
253  mesh,
254  dims,
255  Field<Type>(),
256  makeBoundary
257  (
258  mesh,
259  completeBField,
260  preserveCouples,
261  preserveProcessorOnly,
262  isBoundaryAddressing(mesh, completeBField.size())
263  )
264  )
265 {
266  // Set internalField to the slice of the complete field
268  (
269  SubList<Type>(completeIField, GeoMesh::size(mesh))
270  );
271 
273 }
274 
275 
276 template
277 <
278  class Type,
279  template<class> class PatchField,
280  template<class> class SlicedPatchField,
281  class GeoMesh
282 >
285 (
286  const IOobject& io,
288  const bool preserveCouples
289 )
290 :
291  GeometricField<Type, PatchField, GeoMesh>
292  (
293  io,
294  gf.mesh(),
295  gf.dimensions(),
296  Field<Type>(),
297  makeBoundary(gf.mesh(), gf.boundaryField(), preserveCouples)
298  )
299 {
300  // Set internalField to the internal field
302 
304 }
305 
306 
307 template
308 <
309  class Type,
310  template<class> class PatchField,
311  template<class> class SlicedPatchField,
312  class GeoMesh
313 >
316 (
318 )
319 :
320  GeometricField<Type, PatchField, GeoMesh>
321  (
322  gf,
323  gf.mesh(),
324  gf.dimensions(),
325  Field<Type>(),
326  // preserveCouples = true
327  makeBoundary(gf.mesh(), gf.boundaryField(), true)
328  )
329 {
330  // Set internalField to the internal field
332 }
333 
334 
335 template
336 <
337  class Type,
338  template<class> class PatchField,
339  template<class> class SlicedPatchField,
340  class GeoMesh
341 >
342 Foam::tmp
343 <
345 >
347 clone() const
348 {
349  return tmp
350  <
352  >::New
353  (
354  *this
355  );
356 }
357 
358 
359 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
360 
361 template
362 <
363  class Type,
364  template<class> class PatchField,
365  template<class> class SlicedPatchField,
366  class GeoMesh
367 >
370 {
371  // Set internalField to nullptr to avoid deletion of underlying field
372  UList<Type>::shallowCopy(nullptr);
373 }
374 
375 
376 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
377 
378 template
379 <
380  class Type,
381  template<class> class PatchField,
382  template<class> class SlicedPatchField,
383  class GeoMesh
384 >
387 {
389 }
390 
391 
392 // ************************************************************************* //
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
SlicedGeometricField(const IOobject &, const Mesh &, const dimensionSet &dims, const Field< Type > &completeField, const bool preserveCouples=true)
Construct from components and field to slice.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Generic GeometricField class.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dynamicFvMesh & mesh
Generic templated field type.
Definition: Field.H:62
GeoMesh::Mesh Mesh
The mesh type for the DimensionedField.
const Mesh & mesh() const noexcept
Return mesh.
tmp< SlicedGeometricField< Type, PatchField, SlicedPatchField, GeoMesh > > clone() const
Clone.
void correctBoundaryConditions()
Correct boundary field.
void shallowCopy(T *__restrict__ ptr, const label len) noexcept
Copy the pointer and size.
Definition: UListI.H:309
void correctBoundaryConditions()
Correct boundary field.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:42
volScalarField & p
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180