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 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
30#include "processorFvPatch.H"
31
32// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
33
34template
35<
36 class Type,
37 template<class> class PatchField,
38 template<class> class SlicedPatchField,
39 class GeoMesh
40>
44(
45 const Mesh& mesh,
46 const Field<Type>& completeField,
47 const bool preserveCouples,
48 const bool preserveProcessorOnly
49)
50{
51 auto tbf = tmp<FieldField<PatchField, Type>>::New(mesh.boundary().size());
52 auto& 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,
70 PatchField<Type>::New
71 (
72 mesh.boundary()[patchi].type(),
73 mesh.boundary()[patchi],
74 *this
75 )
76 );
77
78 // Initialize 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
108template
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 auto tbf = tmp<FieldField<PatchField, Type>>::New(mesh.boundary().size());
125 auto& bf = tbf.ref();
126
127 forAll(mesh.boundary(), patchi)
128 {
129 if (preserveCouples && mesh.boundary()[patchi].coupled())
130 {
131 // For coupled patched construct the correct patch field type
132 bf.set
133 (
134 patchi,
135 PatchField<Type>::New
136 (
137 mesh.boundary()[patchi].type(),
138 mesh.boundary()[patchi],
139 *this
140 )
141 );
142
143 // Assign field
144 bf[patchi] == bField[patchi];
145 }
146 else
147 {
148 // Create unallocated copy of patch field
149 bf.set
150 (
151 patchi,
152 new SlicedPatchField<Type>
153 (
154 mesh.boundary()[patchi],
155 DimensionedField<Type, GeoMesh>::null(),
156 bField[patchi]
157 )
158 );
159 }
160 }
161
162 return tbf;
163}
164
165
166// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
167
168template
169<
170 class Type,
171 template<class> class PatchField,
172 template<class> class SlicedPatchField,
173 class GeoMesh
174>
177(
178 const IOobject& io,
179 const Mesh& mesh,
180 const dimensionSet& ds,
181 const Field<Type>& completeField,
182 const bool preserveCouples
183)
184:
185 GeometricField<Type, PatchField, GeoMesh>
186 (
187 io,
188 mesh,
189 ds,
190 Field<Type>(),
191 slicedBoundaryField(mesh, completeField, preserveCouples)
192 )
193{
194 // Set internalField to the slice of the complete field
196 (
197 SubList<Type>(completeField, GeoMesh::size(mesh))
198 );
199
201}
202
203
204template
205<
206 class Type,
207 template<class> class PatchField,
208 template<class> class SlicedPatchField,
209 class GeoMesh
210>
213(
214 const IOobject& io,
215 const Mesh& mesh,
216 const dimensionSet& ds,
217 const Field<Type>& completeIField,
218 const Field<Type>& completeBField,
219 const bool preserveCouples,
220 const bool preserveProcessorOnly
221)
222:
223 GeometricField<Type, PatchField, GeoMesh>
224 (
225 io,
226 mesh,
227 ds,
228 Field<Type>(),
229 slicedBoundaryField
230 (
231 mesh,
232 completeBField,
233 preserveCouples,
234 preserveProcessorOnly
235 )
236 )
237{
238 // Set internalField to the slice of the complete field
240 (
241 SubList<Type>(completeIField, GeoMesh::size(mesh))
242 );
243
245}
246
247
248template
249<
250 class Type,
251 template<class> class PatchField,
252 template<class> class SlicedPatchField,
253 class GeoMesh
254>
257(
258 const IOobject& io,
260 const bool preserveCouples
261)
262:
263 GeometricField<Type, PatchField, GeoMesh>
264 (
265 io,
266 gf.mesh(),
267 gf.dimensions(),
268 Field<Type>(),
269 slicedBoundaryField(gf.mesh(), gf.boundaryField(), preserveCouples)
270 )
271{
272 // Set internalField to the internal field
274
276}
277
278
279template
280<
281 class Type,
282 template<class> class PatchField,
283 template<class> class SlicedPatchField,
284 class GeoMesh
285>
288(
290)
291:
292 GeometricField<Type, PatchField, GeoMesh>
293 (
294 gf,
295 gf.mesh(),
296 gf.dimensions(),
297 Field<Type>(),
298 slicedBoundaryField(gf.mesh(), gf.boundaryField(), true)
299 )
300{
301 // Set internalField to the internal field
303}
304
305
306template
307<
308 class Type,
309 template<class> class PatchField,
310 template<class> class SlicedPatchField,
311 class GeoMesh
312>
314<
316>
318clone() const
319{
320 return tmp
321 <
323 >::New
324 (
325 *this
326 );
327}
328
329
330// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
331
332template
333<
334 class Type,
335 template<class> class PatchField,
336 template<class> class SlicedPatchField,
337 class GeoMesh
338>
341{
342 // Set internalField to nullptr to avoid deletion of underlying field
344}
345
346
347// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
348
349template
350<
351 class Type,
352 template<class> class PatchField,
353 template<class> class SlicedPatchField,
354 class GeoMesh
355>
358{
360}
361
362
363// ************************************************************************* //
const Mesh & mesh() const
Return mesh.
Generic templated field type.
Definition: Field.H:82
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:49
Generic GeometricField class.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
tmp< SlicedGeometricField< Type, PatchField, SlicedPatchField, GeoMesh > > clone() const
Clone.
void correctBoundaryConditions()
Correct boundary field.
A List obtained as a section of another List.
Definition: SubList.H:70
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
void shallowCopy(const UList< T > &list)
Copy the pointer and size held by the given UList.
Definition: UListI.H:272
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
A class for managing temporary objects.
Definition: tmp.H:65
dynamicFvMesh & mesh
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333