GeometricFieldOps.H
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) 2019 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26InNamespace
27 Foam::FieldOps
28
29Description
30 Various utility functions to work on geometric fields
31
32SourceFiles
33 GeometricFieldOps.H
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef GeometricFieldOps_H
38#define GeometricFieldOps_H
39
40#include "FieldOps.H"
41#include "GeometricField.H"
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45namespace Foam
46{
47namespace FieldOps
48{
49
50/*---------------------------------------------------------------------------*\
51 Namespace FieldOps Declarations
52\*---------------------------------------------------------------------------*/
53
54//- Populate a geometric field as the result of a unary operation on an input.
55// It is permissible for inputs/outputs to refer to the same field(s),
56// but partially overlapping regions are ill-defined.
57template
58<
59 class Tout, class T1, class UnaryOp,
60 template<class> class PatchField, class GeoMesh
61>
63(
66 const UnaryOp& op
67)
68{
70 (
71 result.primitiveFieldRef(),
73 op
74 );
75
76 auto& bfld = result.boundaryFieldRef();
77
78 const label len = bfld.size();
79
80 for (label i = 0; i < len; ++i)
81 {
83 (
84 bfld[i],
85 a.boundaryField()[i],
86 op
87 );
88 }
89}
90
91
92//- Populate a geometric field from the binary operation on two inputs.
93// It is permissible for inputs/outputs to refer to the same field(s),
94// but partially overlapping regions are ill-defined.
95template
96<
97 class Tout, class T1, class T2, class BinaryOp,
98 template<class> class PatchField, class GeoMesh
99>
101(
105 const BinaryOp& bop
106)
107{
109 (
110 result.primitiveFieldRef(),
111 a.primitiveField(),
112 b.primitiveField(),
113 bop
114 );
115
116 auto& bfld = result.boundaryFieldRef();
117
118 const label len = bfld.size();
119
120 for (label i = 0; i < len; ++i)
121 {
123 (
124 bfld[i],
125 a.boundaryField()[i],
126 b.boundaryField()[i],
127 bop
128 );
129 }
130}
131
132
133//- Emulate a ternary operation, selecting values from a or b
134//- depending on the binary predicate.
135template
136<
137 class T, class BinaryOp,
138 template<class> class PatchField, class GeoMesh
139>
141(
145 const BinaryOp& bop
146)
147{
149 (
150 result.primitiveFieldRef(),
151 a.primitiveField(),
152 b.primitiveField(),
153 bop
154 );
155
156 auto& bfld = result.boundaryFieldRef();
157
158 const label len = bfld.size();
159
160 for (label i = 0; i < len; ++i)
161 {
163 (
164 bfld[i],
165 a.boundaryField()[i],
166 b.boundaryField()[i],
167 bop
168 );
169 }
170}
171
172
173//- Emulate a ternary operation, selecting field values from a or b
174//- depending on the conditional.
175//
176// Since boolean fields are not normally used, a flip operation is
177// a general requirement.
178template
179<
180 class T, class BoolType, class FlipOp,
181 template<class> class PatchField, class GeoMesh
182>
184(
189 const FlipOp& flip
190)
191{
193 (
194 result.primitiveFieldRef(),
195 cond.primitiveField(),
196 a.primitiveField(),
197 b.primitiveField(),
198 flip
199 );
200
201 auto& bfld = result.boundaryFieldRef();
202
203 const label len = bfld.size();
204
205 for (label i = 0; i < len; ++i)
206 {
208 (
209 bfld[i],
210 cond.boundaryField()[i],
211 a.boundaryField()[i],
212 b.boundaryField()[i],
213 flip
214 );
215 }
216}
217
218
219// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220
221} // End namespace FieldOps
222} // End namespace Foam
223
224// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225
226#endif
227
228// ************************************************************************* //
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.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
const volScalarField & T
void assign(Field< Tout > &result, const Field< T1 > &a, const UnaryOp &op)
Populate a field as the result of a unary operation on an input.
Definition: FieldOps.C:36
void ternarySelect(Field< T > &result, const BoolListType &cond, const Field< T > &a, const Field< T > &b, const FlipOp &flip)
Definition: FieldOps.C:107
void ternary(Field< T > &result, const Field< T > &a, const Field< T > &b, const BinaryOp &bop)
Definition: FieldOps.C:81
Namespace for OpenFOAM.
volScalarField & b
Definition: createFields.H:27