Loading...
Searching...
No Matches
valuePointPatchField.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-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2025 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\*---------------------------------------------------------------------------*/
31
32// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
33
34template<class Type>
36(
37 const dictionary& dict,
39)
40{
41 if (!IOobjectOption::isAnyRead(readOpt)) return false;
42 const auto& p = pointPatchFieldBase::patch();
43
44
45 const auto* eptr = dict.findEntry("value", keyType::LITERAL);
46
47 if (eptr)
48 {
49 Field<Type>::assign(*eptr, p.size());
50 return true;
51 }
52
54 {
56 << "Required entry 'value' : missing for patch " << p.name()
57 << " in dictionary " << dict.relativeName() << nl
59 }
60
61 return false;
62}
64
65template<class Type>
67{
68 const labelUList& meshPoints = pointPatchFieldBase::patch().meshPoints();
69
70 const Field<Type>& iF = this->primitiveField();
71 Field<Type>& pfld = *this;
72
73 pfld.resize_nocopy(meshPoints.size()); // In general this is a no-op
74
75 forAll(meshPoints, i)
76 {
77 pfld[i] = iF[meshPoints[i]];
78 }
79}
81
82// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
83
84template<class Type>
86(
87 const pointPatch& p,
89)
91 pointPatchField<Type>(p, iF),
92 Field<Type>(p.size())
93{}
94
95
96template<class Type>
98(
99 const pointPatch& p,
101 const Type& value
102)
104 pointPatchField<Type>(p, iF),
105 Field<Type>(p.size(), value)
106{}
107
108
109template<class Type>
111(
112 const pointPatch& p,
114 const dictionary& dict,
115 IOobjectOption::readOption requireValue
116)
117:
118 pointPatchField<Type>(p, iF, dict),
119 Field<Type>(p.size())
120{
121 if (!readValueEntry(dict, requireValue))
122 {
123 // Not read (eg, optional and missing): define zero
125 }
126}
127
128
129template<class Type>
131(
132 const valuePointPatchField<Type>& pfld,
133 const pointPatch& p,
134 const DimensionedField<Type, pointMesh>& iF,
135 const pointPatchFieldMapper& mapper
136)
138 pointPatchField<Type>(pfld, p, iF, mapper),
139 Field<Type>(pfld, mapper)
140{}
141
142
143template<class Type>
145(
146 const valuePointPatchField<Type>& pfld,
147 const pointPatch& p,
149 const Type& value
150)
152 pointPatchField<Type>(pfld, p, iF, value),
153 Field<Type>(p.size(), value)
154{}
155
156
157template<class Type>
159(
160 const valuePointPatchField<Type>& pfld,
162)
163:
164 pointPatchField<Type>(pfld, iF),
165 Field<Type>(pfld)
166{}
167
168
169// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
170
171template<class Type>
173(
174 const pointPatchFieldMapper& m
176{
178}
179
180
181template<class Type>
183(
184 const pointPatchField<Type>& ptf,
185 const labelList& addr
186)
187{
189 (
191 (
192 ptf
193 ),
194 addr
195 );
196}
197
198
199template<class Type>
201{
202 if (this->updated())
203 {
204 return;
205 }
206
207 // Get internal field to insert values into
208 Field<Type>& iF = const_cast<Field<Type>&>(this->primitiveField());
209
210 this->setInInternalField(iF, *this);
211
213}
214
215
216template<class Type>
218{
219 // Get internal field to insert values into
220 Field<Type>& iF = const_cast<Field<Type>&>(this->primitiveField());
221
222 this->setInInternalField(iF, *this);
223
225}
226
227
228template<class Type>
230{
233}
234
235
236// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
237
238template<class Type>
240(
243{
245}
246
247
248template<class Type>
250(
251 const pointPatchField<Type>& ptf
252)
254 // pointPatchField has no values to copy, assign internal values
255 this->extrapolateInternal();
256}
257
258
259template<class Type>
261(
262 const Field<Type>& tf
264{
266}
267
268
269template<class Type>
271(
272 const Type& t
274{
277
278
279template<class Type>
281(
284{
286}
287
288
289template<class Type>
291(
292 const pointPatchField<Type>& ptf
293)
295 // pointPatchField has no values to copy, assign internal values
296 this->extrapolateInternal();
297}
298
299
300template<class Type>
302(
303 const Field<Type>& tf
305{
307}
308
309
310template<class Type>
312(
313 const Type& t
314)
315{
317}
318
319
320// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void operator=(const Field< Type > &)
Copy assignment.
Definition Field.C:781
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition Field.C:465
constexpr Field() noexcept
void assign(const entry &e, const label len)
Assign from a primitive dictionary entry with the following behaviour:
Definition Field.C:206
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition Field.C:528
static bool isReadRequired(readOption opt) noexcept
True if (MUST_READ | READ_MODIFIED) bits are set.
readOption
Enumeration defining read preferences.
static bool isAnyRead(readOption opt) noexcept
True if any reading may be required (ie, != NO_READ).
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
commsTypes
Communications types.
Definition UPstream.H:81
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
@ LITERAL
String literal.
Definition keyType.H:82
const pointPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
void setInInternalField(Field< Type1 > &iF, const Field< Type1 > &pF, const labelUList &meshPoints) const
Given the internal field and a patch field, set the patch field in the internal field.
const Field< vector > & primitiveField() const noexcept
virtual void write(Ostream &os) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field, sets updated() to false.
pointPatchField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Basic pointPatch represents a set of points from the mesh.
Definition pointPatch.H:67
Foam::valuePointPatchField.
valuePointPatchField(const valuePointPatchField &)=default
Copy construct.
virtual void write(Ostream &) const
Write.
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
void extrapolateInternal()
Assign the patch field from the internal field.
virtual void rmap(const pointPatchField< Type > &, const labelList &)
Reverse map the given PointPatchField onto.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field.
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
OBJstream os(runTime.globalPath()/outputName)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< label > labelList
A List of labels.
Definition List.H:62
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299