fvPatchField.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 "IOobject.H"
27 #include "dictionary.H"
28 #include "fvMesh.H"
29 #include "fieldMapper.H"
30 #include "volMesh.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class Type>
36 (
37  const fvPatch& p,
39 )
40 :
41  Field<Type>(p.size()),
42  patch_(p),
43  internalField_(iF),
44  updated_(false),
45  manipulatedMatrix_(false)
46 {}
47 
48 
49 template<class Type>
51 (
52  const fvPatch& p,
54  const Field<Type>& f
55 )
56 :
57  Field<Type>(f),
58  patch_(p),
59  internalField_(iF),
60  updated_(false),
61  manipulatedMatrix_(false)
62 {}
63 
64 
65 template<class Type>
67 (
68  const fvPatch& p,
70  const dictionary& dict,
71  const bool valueRequired
72 )
73 :
74  Field<Type>(p.size()),
75  libs_
76  (
77  dict.lookupOrDefault
78  (
79  "libs",
80  fileNameList::null(),
81  dictionary::writeOptionalEntries > 1
82  )
83  ),
84  patch_(p),
85  internalField_(iF),
86  updated_(false),
87  manipulatedMatrix_(false)
88 {
89  if (valueRequired)
90  {
91  if (dict.found("value"))
92  {
94  (
95  Field<Type>("value", iF.dimensions(), dict, p.size())
96  );
97  }
98  else
99  {
101  << "Essential entry 'value' missing"
102  << exit(FatalIOError);
103  }
104  }
105 }
106 
107 
108 template<class Type>
110 (
111  const fvPatchField<Type>& ptf,
112  const fvPatch& p,
114  const fieldMapper& mapper,
115  const bool mappingRequired
116 )
117 :
118  Field<Type>(p.size()),
119  libs_(ptf.libs_),
120  patch_(p),
121  internalField_(iF),
122  updated_(false),
123  manipulatedMatrix_(false)
124 {
125  if (mappingRequired)
126  {
127  mapper(*this, ptf);
128  }
129 }
130 
131 
132 template<class Type>
134 (
135  const fvPatchField<Type>& ptf,
137 )
138 :
139  Field<Type>(ptf),
140  libs_(ptf.libs_),
141  patch_(ptf.patch_),
142  internalField_(iF),
143  updated_(false),
144  manipulatedMatrix_(false)
145 {}
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
150 template<class Type>
152 {
153  return patch_.boundaryMesh().mesh();
154 }
155 
156 
157 template<class Type>
159 {
160  if (&patch_ != &(ptf.patch_))
161  {
163  << "different patches for fvPatchField<Type>s"
164  << abort(FatalError);
165  }
166 }
167 
168 
169 template<class Type>
171 {
172  return patch_.deltaCoeffs()*(*this - patchInternalField());
173 }
174 
175 
176 template<class Type>
179 {
180  return patch_.patchInternalField(internalField_);
181 }
182 
183 
184 template<class Type>
186 {
187  patch_.patchInternalField(internalField_, pif);
188 }
189 
190 
191 template<class Type>
193 (
194  const fvPatchField<Type>& ptf,
195  const fieldMapper& mapper
196 )
197 {
198  mapper(*this, ptf);
199 }
200 
201 
202 template<class Type>
204 {
205  Field<Type>::reset(ptf);
206 }
207 
208 
209 template<class Type>
211 {
212  updated_ = true;
213 }
214 
215 
216 template<class Type>
218 {
219  if (!updated_)
220  {
221  updateCoeffs();
222  }
223 
224  updated_ = false;
225  manipulatedMatrix_ = false;
226 }
227 
228 
229 template<class Type>
231 {
232  manipulatedMatrix_ = true;
233 }
234 
235 
236 template<class Type>
238 {
239  writeEntry(os, "type", type());
240 
241  if (overridesConstraint())
242  {
243  writeEntry(os, "patchType", patch().type());
244  }
245 
246  if (libs_.size())
247  {
248  writeEntry(os, "libs", libs_);
249  }
250 }
251 
252 
253 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
254 
255 template<class Type>
257 (
258  const UList<Type>& ul
259 )
260 {
262 }
263 
264 
265 template<class Type>
267 (
268  const fvPatchField<Type>& ptf
269 )
270 {
271  check(ptf);
273 }
274 
275 
276 template<class Type>
278 (
279  const fvPatchField<Type>& ptf
280 )
281 {
282  check(ptf);
284 }
285 
286 
287 template<class Type>
289 (
290  const fvPatchField<Type>& ptf
291 )
292 {
293  check(ptf);
295 }
296 
297 
298 template<class Type>
300 (
301  const fvPatchField<scalar>& ptf
302 )
303 {
304  if (&patch_ != &ptf.patch())
305  {
307  << "incompatible patches for patch fields"
308  << abort(FatalError);
309  }
310 
312 }
313 
314 
315 template<class Type>
317 (
318  const fvPatchField<scalar>& ptf
319 )
320 {
321  if (&patch_ != &ptf.patch())
322  {
324  << abort(FatalError);
325  }
326 
328 }
329 
330 
331 template<class Type>
333 (
334  const Field<Type>& tf
335 )
336 {
338 }
339 
340 
341 template<class Type>
343 (
344  const Field<Type>& tf
345 )
346 {
348 }
349 
350 
351 template<class Type>
353 (
354  const scalarField& tf
355 )
356 {
358 }
359 
360 
361 template<class Type>
363 (
364  const scalarField& tf
365 )
366 {
368 }
369 
370 
371 template<class Type>
373 (
374  const Type& t
375 )
376 {
378 }
379 
380 
381 template<class Type>
383 (
384  const Type& t
385 )
386 {
388 }
389 
390 
391 template<class Type>
393 (
394  const Type& t
395 )
396 {
398 }
399 
400 
401 template<class Type>
403 (
404  const scalar s
405 )
406 {
408 }
409 
410 
411 template<class Type>
413 (
414  const scalar s
415 )
416 {
418 }
419 
420 
421 template<class Type>
423 (
424  const fvPatchField<Type>& ptf
425 )
426 {
428 }
429 
430 
431 template<class Type>
433 (
434  const Field<Type>& tf
435 )
436 {
438 }
439 
440 
441 template<class Type>
443 (
444  const Type& t
445 )
446 {
448 }
449 
450 
451 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
452 
453 template<class Type>
455 {
456  ptf.write(os);
457 
458  os.check("Ostream& operator<<(Ostream&, const fvPatchField<Type>&");
459 
460  return os;
461 }
462 
463 
464 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
465 
466 #include "fvPatchFieldNew.C"
467 
468 // ************************************************************************* //
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:557
void operator+=(const UList< Type > &)
Definition: Field.C:655
void operator-=(const UList< Type > &)
Definition: Field.C:656
void operator*=(const UList< scalar > &)
Definition: Field.C:657
void operator/=(const UList< scalar > &)
Definition: Field.C:658
void reset(const Field< Type > &)
Reset the field values to the given field.
Definition: Field.C:456
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:217
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:237
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:178
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:170
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:151
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:210
virtual void reset(const fvPatchField< Type > &)
Reset the fvPatchField to the given fvPatchField.
Definition: fvPatchField.C:203
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:230
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:368
virtual void map(const fvPatchField< Type > &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:193
void check(const fvPatchField< Type > &) const
Check fvPatchField<Type> against given fvPatchField<Type>
Definition: fvPatchField.C:158
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: fvPatchField.C:36
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Registry of regIOobjects.
A class for managing temporary objects.
Definition: tmp.H:55
#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
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
error FatalError
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList f(nPoints)
dictionary dict
volScalarField & p