valuePointPatchField.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 "valuePointPatchField.H"
27 #include "fieldMapper.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 template<class Type>
33 {
34  if (this->size() != this->patch().size())
35  {
37  << "field does not correspond to patch. " << endl
38  << "Field size: " << size() << " patch size: "
39  << this->patch().size()
40  << abort(FatalError);
41  }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
46 
47 template<class Type>
49 (
50  const pointPatch& p,
52 )
53 :
54  pointPatchField<Type>(p, iF),
55  Field<Type>(p.size())
56 {}
57 
58 
59 template<class Type>
61 (
62  const pointPatch& p,
64  const dictionary& dict,
65  const bool valueRequired
66 )
67 :
68  pointPatchField<Type>(p, iF, dict),
69  Field<Type>(p.size())
70 {
71  if (dict.found("value"))
72  {
74  (
75  Field<Type>("value", iF.dimensions(), dict, p.size())
76  );
77  }
78  else if (!valueRequired)
79  {
81  }
82  else
83  {
85  << "Essential entry 'value' missing"
86  << exit(FatalIOError);
87  }
88 }
89 
90 
91 template<class Type>
93 (
94  const valuePointPatchField<Type>& ptf,
95  const pointPatch& p,
97  const fieldMapper& mapper
98 )
99 :
100  pointPatchField<Type>(ptf, p, iF, mapper),
101  Field<Type>(mapper(ptf))
102 {}
103 
104 
105 template<class Type>
107 (
108  const valuePointPatchField<Type>& ptf,
110 )
111 :
112  pointPatchField<Type>(ptf, iF),
113  Field<Type>(ptf)
114 {}
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
119 template<class Type>
121 (
122  const pointPatchField<Type>& ptf,
123  const fieldMapper& mapper
124 )
125 {
126  const valuePointPatchField<Type>& vptf =
127  refCast<const valuePointPatchField<Type>>(ptf);
128 
129  mapper(*this, vptf);
130 }
131 
132 
133 template<class Type>
135 (
136  const pointPatchField<Type>& ptf
137 )
138 {
140 }
141 
142 
143 template<class Type>
145 {
146  if (this->updated())
147  {
148  return;
149  }
150 
151  // Get internal field to insert values into
152  Field<Type>& iF = const_cast<Field<Type>&>(this->primitiveField());
153 
154  this->setInternalField(iF, *this);
155 
157 }
158 
159 
160 template<class Type>
162 {
163  // Get internal field to insert values into
164  Field<Type>& iF = const_cast<Field<Type>&>(this->primitiveField());
165 
166  this->setInternalField(iF, *this);
167 
169 }
170 
171 
172 template<class Type>
174 {
176  writeEntry(os, "value", static_cast<const Field<Type>&>(*this));
177 }
178 
179 
180 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
181 
182 template<class Type>
184 (
185  const valuePointPatchField<Type>& ptf
186 )
187 {
189 }
190 
191 
192 template<class Type>
194 (
195  const pointPatchField<Type>& ptf
196 )
197 {
198  Field<Type>::operator=(this->patchInternalField());
199 }
200 
201 
202 template<class Type>
204 (
205  const Field<Type>& tf
206 )
207 {
209 }
210 
211 
212 template<class Type>
214 (
215  const Type& t
216 )
217 {
219 }
220 
221 
222 template<class Type>
224 (
225  const valuePointPatchField<Type>& ptf
226 )
227 {
229 }
230 
231 
232 template<class Type>
234 (
235  const pointPatchField<Type>& ptf
236 )
237 {
238  Field<Type>::operator=(this->patchInternalField());
239 }
240 
241 
242 template<class Type>
244 (
245  const Field<Type>& tf
246 )
247 {
249 }
250 
251 
252 template<class Type>
254 (
255  const Type& t
256 )
257 {
259 }
260 
261 
262 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Pre-declare SubField and related Field type.
Definition: Field.H:83
void operator=(const Field< Type > &)
Definition: Field.C:550
void reset(const Field< Type > &)
Reset the field values to the given field.
Definition: Field.C:456
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
commsTypes
Types of communications.
Definition: UPstream.H:65
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class for point-mesh patch fields.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:57
Foam::valuePointPatchField.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
virtual void write(Ostream &) const
Write.
valuePointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void reset(const pointPatchField< Type > &)
Reset the pointPatchField to the given pointPatchField.
virtual void map(const pointPatchField< Type > &, const fieldMapper &)
Map the given pointPatchField onto this pointPatchField.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const tensorField & tf
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:129
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
IOerror FatalIOError
error FatalError
dictionary dict
volScalarField & p