fvPatchField.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) 2015-2025 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
29 #include "IOobject.H"
30 #include "dictionary.H"
31 #include "fvMesh.H"
32 #include "fvPatchFieldMapper.H"
33 #include "volMesh.H"
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
37 template<class Type>
39 (
40  const dictionary& dict,
42 )
43 {
44  if (!IOobjectOption::isAnyRead(readOpt)) return false;
45  const auto& p = fvPatchFieldBase::patch();
46 
47 
48  const auto* eptr = dict.findEntry("value", keyType::LITERAL);
49 
50  if (eptr)
51  {
52  Field<Type>::assign(*eptr, p.size());
53  return true;
54  }
55 
56  if (IOobjectOption::isReadRequired(readOpt))
57  {
59  << "Required entry 'value' : missing for patch " << p.name()
60  << " in dictionary " << dict.relativeName() << nl
61  << exit(FatalIOError);
62  }
63 
64  return false;
65 }
66 
67 
68 template<class Type>
70 {
71  const auto& p = fvPatchFieldBase::patch();
72  this->resize_nocopy(p.size()); // In general this is a no-op
73  p.patchInternalField(internalField_, *this);
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78 
79 template<class Type>
81 (
82  const fvPatch& p,
84 )
85 :
87  Field<Type>(p.size()),
88  internalField_(iF)
89 {}
90 
91 
92 template<class Type>
94 (
95  const fvPatch& p,
97  const word& patchType
98 )
99 :
100  fvPatchFieldBase(p, patchType),
101  Field<Type>(p.size()),
102  internalField_(iF)
103 {}
104 
105 
106 template<class Type>
108 (
109  const fvPatch& p,
111  const Type& value
112 )
113 :
115  Field<Type>(p.size(), value),
116  internalField_(iF)
117 {}
118 
119 
120 template<class Type>
122 (
123  const fvPatch& p,
125  const Field<Type>& pfld
126 )
127 :
129  Field<Type>(pfld),
130  internalField_(iF)
131 {}
132 
133 
134 template<class Type>
136 (
137  const fvPatch& p,
139  Field<Type>&& pfld
140 )
141 :
143  Field<Type>(std::move(pfld)),
144  internalField_(iF)
145 {}
146 
147 
148 template<class Type>
150 (
151  const fvPatch& p,
153  const dictionary& dict,
154  IOobjectOption::readOption requireValue
155 )
156 :
158  Field<Type>(p.size()),
159  internalField_(iF)
160 {
161  readValueEntry(dict, requireValue);
162 }
163 
164 
165 template<class Type>
167 (
168  const fvPatchField<Type>& ptf,
169  const fvPatch& p,
171  const fvPatchFieldMapper& mapper
172 )
173 :
174  fvPatchFieldBase(ptf, p),
175  Field<Type>(p.size()),
176  internalField_(iF)
177 {
178  // For unmapped faces set to internal field value (zero-gradient)
179  if (notNull(iF) && mapper.hasUnmapped())
180  {
182  }
183  this->map(ptf, mapper);
184 }
185 
186 
187 template<class Type>
189 (
190  const fvPatchField<Type>& ptf
191 )
192 :
194  Field<Type>(ptf),
195  internalField_(ptf.internalField_)
196 {}
197 
198 
199 template<class Type>
201 (
202  const fvPatchField<Type>& ptf,
204 )
205 :
206  fvPatchFieldBase(ptf),
207  Field<Type>(ptf),
208  internalField_(iF)
209 {}
210 
211 
212 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
213 
214 template<class Type>
216 {
218 }
219 
220 
221 template<class Type>
223 {
224  return patch().deltaCoeffs()*(*this - patchInternalField());
225 }
226 
227 
228 template<class Type>
231 {
232  return patch().patchInternalField(internalField_);
233 }
234 
235 
236 template<class Type>
238 {
239  patch().patchInternalField(internalField_, pfld);
240 }
241 
242 
243 template<class Type>
245 (
246  const fvPatchFieldMapper& mapper
247 )
248 {
249  Field<Type>& f = *this;
250 
251  if (!this->size() && !mapper.distributed())
252  {
253  f.resize_nocopy(mapper.size());
254  if (f.size())
255  {
256  f = this->patchInternalField();
257  }
258  }
259  else
260  {
261  // Map all faces provided with mapping data
262  Field<Type>::autoMap(mapper);
263 
264 
265  // For unmapped faces set to internal field value (zero-gradient)
266  if (mapper.hasUnmapped())
267  {
268  Field<Type> pif(this->patchInternalField());
269 
270  if
271  (
272  mapper.direct()
273  && notNull(mapper.directAddressing())
274  && mapper.directAddressing().size()
275  )
276  {
277  const labelUList& mapAddressing = mapper.directAddressing();
278 
279  forAll(mapAddressing, i)
280  {
281  if (mapAddressing[i] < 0)
282  {
283  f[i] = pif[i];
284  }
285  }
286  }
287  else if (!mapper.direct() && mapper.addressing().size())
288  {
289  const labelListList& mapAddressing = mapper.addressing();
290 
291  forAll(mapAddressing, i)
292  {
293  const labelList& localAddrs = mapAddressing[i];
294 
295  if (!localAddrs.size())
296  {
297  f[i] = pif[i];
298  }
299  }
300  }
301  }
302  }
303 }
304 
305 
306 template<class Type>
308 (
309  const fvPatchField<Type>& ptf,
310  const labelList& addr
311 )
312 {
313  Field<Type>::rmap(ptf, addr);
314 }
315 
316 
317 template<class Type>
319 {
321 }
322 
323 
324 template<class Type>
326 {
327  // Default behaviour ignores the weights
328  if (!updated())
329  {
330  updateCoeffs();
331 
333  }
334 }
335 
336 
337 template<class Type>
339 {
340  if (!updated())
341  {
342  updateCoeffs();
343  }
344 
347 }
348 
349 
350 template<class Type>
352 {
354 }
355 
356 
357 template<class Type>
359 (
360  fvMatrix<Type>& matrix,
361  const scalarField& weights
362 )
363 {
365 }
366 
367 
368 template<class Type>
370 (
371  fvMatrix<Type>& matrix,
372  const label iMatrix,
373  const direction cmp
374 )
375 {
377 }
378 
379 
380 template<class Type>
382 {
383  os.writeEntry("type", type());
384 
385  if (!patchType().empty())
386  {
387  os.writeEntry("patchType", patchType());
388  }
389  if (useImplicit())
390  {
391  os.writeEntry("useImplicit", "true");
392  }
393 }
394 
395 
396 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
397 
398 template<class Type>
400 (
401  const UList<Type>& ul
402 )
403 {
405 }
406 
407 
408 template<class Type>
410 (
411  const fvPatchField<Type>& ptf
412 )
413 {
416 }
417 
418 
419 template<class Type>
421 (
422  const fvPatchField<Type>& ptf
423 )
424 {
427 }
428 
429 
430 template<class Type>
432 (
433  const fvPatchField<Type>& ptf
434 )
435 {
438 }
439 
440 
441 template<class Type>
443 (
444  const fvPatchField<scalar>& ptf
445 )
446 {
449 }
450 
451 
452 template<class Type>
454 (
455  const fvPatchField<scalar>& ptf
456 )
457 {
460 }
461 
462 
463 template<class Type>
465 (
466  const Field<Type>& tf
467 )
468 {
470 }
471 
472 
473 template<class Type>
475 (
476  const Field<Type>& tf
477 )
478 {
480 }
481 
482 
483 template<class Type>
485 (
486  const scalarField& tf
487 )
488 {
490 }
491 
492 
493 template<class Type>
495 (
496  const scalarField& tf
497 )
498 {
500 }
501 
502 
503 template<class Type>
505 (
506  const Type& t
507 )
508 {
510 }
511 
512 
513 template<class Type>
515 (
516  const Type& t
517 )
518 {
520 }
521 
522 
523 template<class Type>
525 (
526  const Type& t
527 )
528 {
530 }
531 
532 
533 template<class Type>
535 (
536  const scalar s
537 )
538 {
540 }
541 
542 
543 template<class Type>
545 (
546  const scalar s
547 )
548 {
550 }
551 
552 
553 template<class Type>
555 (
556  const fvPatchField<Type>& ptf
557 )
558 {
560 }
561 
562 
563 template<class Type>
565 (
566  const Field<Type>& tf
567 )
568 {
570 }
571 
572 
573 template<class Type>
575 (
576  const Type& t
577 )
578 {
580 }
581 
582 
583 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
584 
585 template<class Type>
586 Foam::Ostream& Foam::operator<<(Ostream& os, const fvPatchField<Type>& ptf)
587 {
588  ptf.write(os);
589 
591 
592  return os;
593 }
594 
595 
596 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual bool direct() const =0
Is it a direct (non-interpolating) mapper?
virtual const labelListList & addressing() const
Return the interpolation addressing.
Definition: FieldMapper.H:115
virtual const labelUList & directAddressing() const
Return the direct addressing values.
Definition: FieldMapper.H:91
virtual label size() const =0
The size of the mapper.
virtual bool distributed() const
Does the mapper have remote contributions?
Definition: FieldMapper.H:76
virtual bool hasUnmapped() const =0
Any unmapped values?
Generic templated field type.
Definition: Field.H:166
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
void operator+=(const UList< Type > &)
Definition: Field.C:833
void operator-=(const UList< Type > &)
Definition: Field.C:834
void operator*=(const UList< scalar > &)
Definition: Field.C:835
void operator/=(const UList< scalar > &)
Definition: Field.C:836
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition: Field.C:302
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:528
readOption
Enumeration defining read preferences.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition: ListI.H:168
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:59
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
commsTypes
Communications types.
Definition: UPstream.H:78
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:132
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:116
Template invariant parts for fvPatchField.
Definition: fvPatchField.H:78
void checkPatch(const fvPatchFieldBase &rhs) const
Check that patches are identical.
void setManipulated(bool state) noexcept
Set matrix manipulated state.
Definition: fvPatchField.H:309
void setUpdated(bool state) noexcept
Set updated state.
Definition: fvPatchField.H:293
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:353
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:238
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:374
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
Definition: fvPatchField.C:223
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:215
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:311
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:344
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:318
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:301
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
Definition: fvPatchField.C:32
void extrapolateInternal()
Assign the patch field from the internal field.
Definition: fvPatchField.C:62
void check(const fvPatchField< Type > &) const
Check against given patch field.
Definition: fvPatchField.C:208
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field, sets updated() to false.
Definition: fvPatchField.C:331
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: fvPatchField.C:74
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
A class for managing temporary objects.
Definition: tmp.H:66
A class for handling words, derived from Foam::string.
Definition: word.H:66
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:629
OBJstream os(runTime.globalPath()/outputName)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define FUNCTION_NAME
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:28
const std::string patch
OpenFOAM patch number as a std::string.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
bool notNull(const T *ptr) noexcept
True if ptr is not a pointer (of type T) to the nullObject.
Definition: nullObject.H:267
List< label > labelList
A List of labels.
Definition: List.H:61
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:801
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
uint8_t direction
Definition: direction.H:46
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
UList< label > labelUList
A UList of labels.
Definition: UList.H:76
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:50
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:286