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-2024 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<class Type, class GeoMesh>
34 (
35  const Mesh& mesh,
36  const Field<Type>& completeField,
37  const bool preserveCouples,
38  const bool preserveProcessorOnly
39 )
40 {
41  tmp<FieldField<GeoMesh::template PatchField, Type>> tbf
42  (
43  new FieldField<GeoMesh::template PatchField, Type>
44  (
45  mesh.boundary().size()
46  )
47  );
48  FieldField<GeoMesh::template PatchField, Type>& bf = tbf.ref();
49 
51  {
52  if
53  (
54  preserveCouples
55  && mesh.boundary()[patchi].coupled()
56  && (
57  !preserveProcessorOnly
58  || isA<processorFvPatch>(mesh.boundary()[patchi])
59  )
60  )
61  {
62  // For coupled patched construct the correct patch field type
63  bf.set
64  (
65  patchi,
67  (
68  mesh.boundary()[patchi].type(),
69  mesh.boundary()[patchi],
70  *this
71  )
72  );
73 
74  // Initialise the values on the coupled patch to those of the slice
75  // of the given field.
76  // Note: these will usually be over-ridden by the boundary field
77  // evaluation e.g. in the case of processor and cyclic patches.
78  bf[patchi] = SlicedPatch
79  (
80  mesh.boundary()[patchi],
82  completeField
83  );
84  }
85  else
86  {
87  bf.set
88  (
89  patchi,
90  new SlicedPatch
91  (
92  mesh.boundary()[patchi],
94  completeField
95  )
96  );
97  }
98  }
99 
100  return tbf;
101 }
102 
103 
104 template<class Type, class GeoMesh>
107 (
108  const Mesh& mesh,
109  const FieldField<GeoMesh::template PatchField, Type>& bField,
110  const bool preserveCouples
111 )
112 {
113  tmp<FieldField<GeoMesh::template PatchField, Type>> tbf
114  (
115  new FieldField<GeoMesh::template PatchField, Type>
116  (
117  mesh.boundary().size()
118  )
119  );
120  FieldField<GeoMesh::template PatchField, Type>& bf = tbf.ref();
121 
123  {
124  if (preserveCouples && mesh.boundary()[patchi].coupled())
125  {
126  // For coupled patched construct the correct patch field type
127  bf.set
128  (
129  patchi,
130  Patch::New
131  (
132  mesh.boundary()[patchi].type(),
133  mesh.boundary()[patchi],
134  *this
135  )
136  );
137 
138  // Assign field
139  bf[patchi] == bField[patchi];
140  }
141  else
142  {
143  // Create unallocated copy of patch field
144  bf.set
145  (
146  patchi,
147  new SlicedPatch
148  (
149  mesh.boundary()[patchi],
151  bField[patchi]
152  )
153  );
154  }
155  }
156 
157  return tbf;
158 }
159 
160 
161 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
162 
163 template<class Type, class GeoMesh>
165 (
166  const IOobject& io,
167  const Mesh& mesh,
168  const dimensionSet& ds,
169  const Field<Type>& completeField,
170  const bool preserveCouples
171 )
172 :
173  GeometricField<Type, GeoMesh>
174  (
175  io,
176  mesh,
177  ds,
178  Field<Type>(),
179  slicedBoundaryField(mesh, completeField, preserveCouples)
180  )
181 {
182  // Set the internalField to the slice of the complete field
184  (
185  typename Field<Type>::subField(completeField, GeoMesh::size(mesh))
186  );
187 
189 }
190 
191 
192 template<class Type, class GeoMesh>
194 (
195  const IOobject& io,
196  const Mesh& mesh,
197  const dimensionSet& ds,
198  const Field<Type>& completeIField,
199  const Field<Type>& completeBField,
200  const bool preserveCouples,
201  const bool preserveProcessorOnly
202 )
203 :
204  GeometricField<Type, GeoMesh>
205  (
206  io,
207  mesh,
208  ds,
209  Field<Type>(),
210  slicedBoundaryField
211  (
212  mesh,
213  completeBField,
214  preserveCouples,
215  preserveProcessorOnly
216  )
217  )
218 {
219  // Set the internalField to the slice of the complete field
221  (
222  typename Field<Type>::subField(completeIField, GeoMesh::size(mesh))
223  );
224 
226 }
227 
228 
229 template<class Type, class GeoMesh>
231 (
232  const IOobject& io,
234  const bool preserveCouples
235 )
236 :
237  GeometricField<Type, GeoMesh>
238  (
239  io,
240  gf.mesh(),
241  gf.dimensions(),
242  Field<Type>(),
243  slicedBoundaryField(gf.mesh(), gf.boundaryField(), preserveCouples)
244  )
245 {
246  // Set the internalField to the supplied internal field
248 
250 }
251 
252 
253 template<class Type, class GeoMesh>
255 (
257 )
258 :
259  GeometricField<Type, GeoMesh>
260  (
261  gf,
262  gf.mesh(),
263  gf.dimensions(),
264  Field<Type>(),
265  slicedBoundaryField(gf.mesh(), gf.boundaryField(), true)
266  )
267 {
268  // Set the internalField to the supplied internal field
270 }
271 
272 
273 template<class Type, class GeoMesh>
276 {
278  (
280  );
281 }
282 
283 
284 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
285 
286 template<class Type, class GeoMesh>
288 {
289  // Set the internalField storage pointer to nullptr before its destruction
290  // to protect the field it a slice of.
292 }
293 
294 
295 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
296 
297 template<class Type, class GeoMesh>
300 {
301  const Mesh& mesh = this->mesh();
302 
303  label completeSize = GeoMesh::size(mesh);
304 
306  {
307  completeSize += mesh.boundaryMesh()[patchi].size();
308  }
309 
310  tmp<Field<Type>> tCompleteField(new Field<Type>(completeSize));
311  Field<Type>& completeField(tCompleteField.ref());
312 
313  typename Field<Type>::subField(completeField, GeoMesh::size(mesh))
314  = this->primitiveField();
315 
316  label start = GeoMesh::size(mesh);
317 
319  {
320  if
321  (
322  mesh.boundary()[patchi].size()
323  == mesh.boundaryMesh()[patchi].size()
324  )
325  {
326  typename Field<Type>::subField
327  (
328  completeField,
329  mesh.boundary()[patchi].size(),
330  start
331  ) = this->boundaryField()[patchi];
332  }
333  else
334  {
335  typename Field<Type>::subField
336  (
337  completeField,
339  start
340  ) = Zero;
341  }
342 
343  start += mesh.boundaryMesh()[patchi].size();
344  }
345 
346  return tCompleteField;
347 }
348 
349 
350 template<class Type, class GeoMesh>
352 {
354 }
355 
356 
357 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
static const DimensionedField< Type, GeoMesh, PrimitiveField > & null()
Return a null DimensionedField.
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:83
SubField< Type > subField
Declare type of subField.
Definition: Field.H:101
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 primitive field.
void correctBoundaryConditions()
Correct boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Specialisation of GeometricField which holds slices of given complete fields in a form that they act ...
tmp< SlicedGeometricField< Type, 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
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
Dimension set for the base types.
Definition: dimensionSet.H:125
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:962
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
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:197
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96