fvPatchField.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2023 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 Class
28  Foam::fvPatchField
29 
30 Description
31  Abstract base class with a fat-interface to all derived classes
32  covering all possible ways in which they might be used.
33 
34  The first level of derivation is to basic patchFields which cover
35  zero-gradient, fixed-gradient, fixed-value and mixed conditions.
36 
37  The next level of derivation covers all the specialised types with
38  specific evaluation procedures, particularly with respect to specific
39  fields.
40 
41 SourceFiles
42  fvPatchField.C
43  fvPatchFieldBase.C
44  fvPatchFieldNew.C
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #ifndef Foam_fvPatchField_H
49 #define Foam_fvPatchField_H
50 
51 #include "fvPatch.H"
52 #include "DimensionedField.H"
53 #include "fieldTypes.H"
54 #include "scalarField.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 
61 // Forward Declarations
62 class dictionary;
63 class objectRegistry;
65 class volMesh;
66 
67 template<class Type> class fvPatchField;
68 template<class Type> class calculatedFvPatchField;
69 template<class Type> class fvMatrix;
70 
71 template<class Type>
72 Ostream& operator<<(Ostream&, const fvPatchField<Type>&);
73 
74 
75 /*---------------------------------------------------------------------------*\
76  Class fvPatchFieldBase Declaration
77 \*---------------------------------------------------------------------------*/
78 
79 //- Template invariant parts for fvPatchField
80 class fvPatchFieldBase
81 {
82  // Private Data
83 
84  //- Reference to patch
85  const fvPatch& patch_;
86 
87  //- Update index used so that updateCoeffs is called only once during
88  //- the construction of the matrix
89  bool updated_;
90 
91  //- Update index used so that manipulateMatrix is called only once
92  //- during the construction of the matrix
93  bool manipulatedMatrix_;
94 
95  //- Use implicit formulation
96  bool useImplicit_;
97 
98  //- Optional patch type
99  // Used to allow specified boundary conditions to be applied
100  // to constraint patches by providing the constraint
101  // patch type as 'patchType'
102  word patchType_;
103 
104 
105 protected:
106 
107  // Protected Member Functions
108 
109  //- Read dictionary entries.
110  // Useful when initially constructed without a dictionary
111  virtual void readDict(const dictionary& dict);
112 
113 
114 public:
115 
116  //- Set updated state
117  void setUpdated(bool state) noexcept
118  {
119  updated_ = state;
120  }
121 
122  //- Set matrix manipulated state
123  void setManipulated(bool state) noexcept
124  {
125  manipulatedMatrix_ = state;
126  }
127 
128 
129 public:
131  //- Debug switch to disallow the use of generic fvPatchField
132  static int disallowGenericPatchField;
133 
134  //- Runtime type information
135  TypeName("fvPatchField");
136 
137 
138  // Constructors
139 
140  //- Construct from patch
141  explicit fvPatchFieldBase(const fvPatch& p);
142 
143  //- Construct from patch and patch type
144  explicit fvPatchFieldBase(const fvPatch& p, const word& patchType);
145 
146  //- Construct from patch and dictionary
147  fvPatchFieldBase(const fvPatch& p, const dictionary& dict);
148 
149  //- Copy construct with new patch
150  fvPatchFieldBase(const fvPatchFieldBase& rhs, const fvPatch& p);
151 
152  //- Copy construct
154 
155 
156  //- Destructor
157  virtual ~fvPatchFieldBase() = default;
158 
159 
160  // Static Member Functions
161 
162  //- The type name for \c empty patch fields
163  static const word& emptyType() noexcept
164  {
166  }
167 
168  //- The type name for \c calculated patch fields
169  static const word& calculatedType() noexcept
170  {
172  }
173 
174  //- The type name for \c extrapolatedCalculated patch fields
175  //- combines \c zero-gradient and \c calculated
176  static const word& extrapolatedCalculatedType() noexcept
177  {
179  }
180 
181  //- The type name for \c zeroGradient patch fields
182  static const word& zeroGradientType() noexcept
183  {
185  }
186 
187 
188  // Member Functions
189 
190  // Attributes
191 
192  //- True if the value of the patch field is altered by assignment
193  virtual bool assignable() const
194  {
195  return true;
196  }
197 
198  //- True if the patch field fixes a value.
199  // Needed to check if a level has to be specified while solving
200  // Poissons equations.
201  virtual bool fixesValue() const
202  {
203  return false;
204  }
205 
206  //- True if the patch field is coupled
207  virtual bool coupled() const
208  {
209  return false;
210  }
211 
212 
213  // Access
214 
215  //- The associated objectRegistry
216  const objectRegistry& db() const;
217 
218  //- Return the patch
219  const fvPatch& patch() const noexcept
220  {
221  return patch_;
222  }
223 
224  //- The optional patch type
225  const word& patchType() const noexcept
226  {
227  return patchType_;
228  }
229 
230  //- The optional patch type
231  word& patchType() noexcept
232  {
233  return patchType_;
234  }
235 
236 
237  // Solution
238 
239  //- True if the boundary condition has already been updated
240  bool updated() const noexcept
241  {
242  return updated_;
243  }
244 
245  //- True if the matrix has already been manipulated
246  bool manipulatedMatrix() const noexcept
247  {
248  return manipulatedMatrix_;
249  }
250 
251  //- Use implicit formulation for coupled patches only
252  bool useImplicit() const noexcept
253  {
254  return useImplicit_;
255  }
256 
257  //- Set useImplicit on/off
258  // \return old value
259  bool useImplicit(bool on) noexcept
260  {
261  bool old(useImplicit_);
262  useImplicit_ = on;
263  return old;
264  }
265 
266 
267  // Check
268 
269  //- Check that patches are identical
270  void checkPatch(const fvPatchFieldBase& rhs) const;
271 };
272 
273 
274 /*---------------------------------------------------------------------------*\
275  Class fvPatchField Declaration
276 \*---------------------------------------------------------------------------*/
278 template<class Type>
279 class fvPatchField
280 :
281  public fvPatchFieldBase,
282  public Field<Type>
283 {
284  // Private Data
286  //- Reference to internal field
287  const DimensionedField<Type, volMesh>& internalField_;
288 
289 protected:
290 
291  // Protected Member Functions
292 
293  //- Read the "value" entry into \c *this.
294  // The reading can be optional (default), mandatory etc.
295  // \returns True on success
297  (
298  const dictionary& dict,
300  );
301 
302  //- Write \c *this field as a "value" entry
303  void writeValueEntry(Ostream& os) const
304  {
305  Field<Type>::writeEntry("value", os);
306  }
307 
308  //- Assign the patch field from the internal field
309  void extrapolateInternal();
310 
311 
312 public:
313 
314  //- The internal field type associated with the patch field
316 
317  //- The patch type for the patch field
318  typedef fvPatch Patch;
319 
320  //- Type for a \em calculated patch
323 
324  // Declare run-time constructor selection tables
325 
327  (
328  tmp,
329  fvPatchField,
330  patch,
331  (
332  const fvPatch& p,
334  ),
335  (p, iF)
336  );
337 
339  (
340  tmp,
341  fvPatchField,
342  patchMapper,
343  (
344  const fvPatchField<Type>& ptf,
345  const fvPatch& p,
347  const fvPatchFieldMapper& m
348  ),
349  (dynamic_cast<const fvPatchFieldType&>(ptf), p, iF, m)
350  );
351 
353  (
354  tmp,
355  fvPatchField,
356  dictionary,
357  (
358  const fvPatch& p,
360  const dictionary& dict
361  ),
362  (p, iF, dict)
363  );
364 
365 
366  // Constructors
367 
368  //- Construct from patch and internal field
370  (
371  const fvPatch&,
373  );
374 
375  //- Construct from patch, internal field and patch type
377  (
378  const fvPatch&,
380  const word& patchType
381  );
382 
383  //- Construct from patch, internal field and value
385  (
386  const fvPatch&,
388  const Type& value
389  );
390 
391  //- Construct from patch, internal field and patch field
393  (
394  const fvPatch&,
396  const Field<Type>& pfld
397  );
398 
399  //- Construct from patch, internal field and patch field
401  (
402  const fvPatch&,
404  Field<Type>&& pfld
405  );
406 
407  //- Construct from patch, internal field and dictionary
409  (
410  const fvPatch&,
412  const dictionary& dict,
415  );
416 
417  //- Construct, forwarding to readOption variant
419  (
420  const fvPatch& p,
422  const dictionary& dict,
423  const bool valueReqd
424  )
425  :
427  (
428  p, iF, dict,
429  (valueReqd? IOobjectOption::MUST_READ : IOobjectOption::NO_READ)
430  )
431  {}
432 
433  //- Construct by mapping the given fvPatchField onto a new patch
435  (
436  const fvPatchField<Type>&,
437  const fvPatch&,
439  const fvPatchFieldMapper&
440  );
441 
442  //- Construct as copy
444 
445  //- Construct as copy setting internal field reference
447  (
448  const fvPatchField<Type>&,
450  );
451 
452  //- Construct and return a clone
453  virtual tmp<fvPatchField<Type>> clone() const
454  {
455  return tmp<fvPatchField<Type>>::New(*this);
456  }
457 
458  //- Construct and return a clone setting internal field reference
459  virtual tmp<fvPatchField<Type>> clone
460  (
461  const DimensionedField<Type, volMesh>& iF
462  ) const
463  {
464  return tmp<fvPatchField<Type>>::New(*this, iF);
465  }
466 
467 
468  // Selectors
469 
470  //- Return a pointer to a new patchField created on freestore given
471  // patch and internal field
472  // (does not set the patch field values)
473  static tmp<fvPatchField<Type>> New
474  (
475  const word& patchFieldType,
476  const fvPatch&,
477  const DimensionedField<Type, volMesh>&
478  );
479 
480  //- Return a pointer to a new patchField created on freestore given
481  // patch and internal field
482  // (does not set the patch field values).
483  // Allows override of constraint type
484  static tmp<fvPatchField<Type>> New
485  (
486  const word& patchFieldType,
487  const word& actualPatchType,
488  const fvPatch&,
489  const DimensionedField<Type, volMesh>&
490  );
491 
492  //- Return a pointer to a new patchField created on freestore from
493  // a given fvPatchField mapped onto a new patch
494  static tmp<fvPatchField<Type>> New
495  (
496  const fvPatchField<Type>&,
497  const fvPatch&,
498  const DimensionedField<Type, volMesh>&,
499  const fvPatchFieldMapper&
500  );
501 
502  //- Return a pointer to a new patchField created on freestore
503  // from dictionary
504  static tmp<fvPatchField<Type>> New
505  (
506  const fvPatch&,
507  const DimensionedField<Type, volMesh>&,
508  const dictionary&
509  );
510 
511  //- Return a pointer to a new calculatedFvPatchField created on
512  // freestore without setting patchField values
514  (
515  const fvPatch&
516  );
517 
518  //- Return a pointer to a new calculatedFvPatchField created on
519  // freestore without setting patchField values
520  template<class Type2>
522  (
523  const fvPatchField<Type2>&
524  );
525 
526 
527  //- Destructor
528  virtual ~fvPatchField() = default;
529 
530 
531  // Member Functions
532 
533  // Access
534 
535  //- Return const-reference to the dimensioned internal field
537  {
538  return internalField_;
539  }
540 
541  //- Return const-reference to the internal field values
542  const Field<Type>& primitiveField() const noexcept
543  {
544  return internalField_;
545  }
546 
547 
548  // Mapping Functions
549 
550  //- Map (and resize as needed) from self given a mapping object
551  virtual void autoMap
552  (
553  const fvPatchFieldMapper&
554  );
556  //- Reverse map the given fvPatchField onto this fvPatchField
557  virtual void rmap
558  (
559  const fvPatchField<Type>&,
560  const labelList&
561  );
562 
563 
564  // Evaluation Functions
565 
566  //- Return patch-normal gradient
567  virtual tmp<Field<Type>> snGrad() const;
568 
569  //- Return patch-normal gradient for coupled-patches
570  // using the deltaCoeffs provided
571  virtual tmp<Field<Type>> snGrad
572  (
573  const scalarField& deltaCoeffs
574  ) const
575  {
577  return *this;
578  }
579 
580  //- Update the coefficients associated with the patch field
581  // Sets Updated to true
582  virtual void updateCoeffs();
583 
584  //- Update the coefficients associated with the patch field
585  // with a weight field (0..1). This weight field is usually
586  // provided as the amount of geometric overlap for 'duplicate'
587  // patches. Sets Updated to true
588  virtual void updateWeightedCoeffs(const scalarField& weights);
589 
590  //- Return internal field next to patch
591  virtual tmp<Field<Type>> patchInternalField() const;
592 
593  //- Extract internal field next to patch
594  // \param [out] pfld The extracted patch field.
595  virtual void patchInternalField(Field<Type>& pfld) const;
596 
597  //- Return patchField on the opposite patch of a coupled patch
598  virtual tmp<Field<Type>> patchNeighbourField() const
599  {
601  return *this;
602  }
603 
604  //- Initialise the evaluation of the patch field
605  virtual void initEvaluate
606  (
607  const Pstream::commsTypes commsType =
609  )
610  {}
611 
612  //- Evaluate the patch field, sets updated() to false
613  virtual void evaluate
614  (
615  const Pstream::commsTypes commsType =
617  );
618 
619  //- Initialise the evaluation of the patch field after a local
620  // operation
621  virtual void initEvaluateLocal
622  (
623  const Pstream::commsTypes commsType =
625  )
626  {}
627 
628  //- Evaluate the patch field after a local operation (e.g. *=)
629  virtual void evaluateLocal
630  (
631  const Pstream::commsTypes commsType =
633  )
634  {}
635 
636  //- Return the matrix diagonal coefficients corresponding to the
637  // evaluation of the value of this patchField with given weights
638  virtual tmp<Field<Type>> valueInternalCoeffs
639  (
640  const tmp<Field<scalar>>&
641  ) const
642  {
644  return *this;
645  }
646 
647  //- Return the matrix source coefficients corresponding to the
648  // evaluation of the value of this patchField with given weights
649  virtual tmp<Field<Type>> valueBoundaryCoeffs
650  (
651  const tmp<Field<scalar>>&
652  ) const
653  {
655  return *this;
656  }
657 
658  //- Return the matrix diagonal coefficients corresponding to the
659  // evaluation of the gradient of this patchField
660  virtual tmp<Field<Type>> gradientInternalCoeffs() const
661  {
663  return *this;
664  }
665 
666  //- Return the matrix diagonal coefficients corresponding to the
667  // evaluation of the gradient of this coupled patchField
668  // using the deltaCoeffs provided
670  (
671  const scalarField& deltaCoeffs
672  ) const
673  {
675  return *this;
676  }
677 
678  //- Return the matrix source coefficients corresponding to the
679  // evaluation of the gradient of this patchField
680  virtual tmp<Field<Type>> gradientBoundaryCoeffs() const
681  {
683  return *this;
684  }
685 
686  //- Return the matrix source coefficients corresponding to the
687  // evaluation of the gradient of this coupled patchField
688  // using the deltaCoeffs provided
689  virtual tmp<Field<Type>> gradientBoundaryCoeffs
690  (
691  const scalarField& deltaCoeffs
692  ) const
693  {
695  return *this;
696  }
697 
698 
699  //- Manipulate matrix
700  virtual void manipulateMatrix(fvMatrix<Type>& matrix);
701 
702  //- Manipulate matrix with given weights
703  virtual void manipulateMatrix
704  (
705  fvMatrix<Type>& matrix,
706  const scalarField& weights
707  );
708 
709  //- Manipulate fvMatrix
710  virtual void manipulateMatrix
711  (
712  fvMatrix<Type>& matrix,
713  const label iMatrix,
714  const direction cmp
715  );
716 
717 
718  // Other
719 
720  //- Write
721  virtual void write(Ostream&) const;
722 
723  //- Check against given patch field
724  void check(const fvPatchField<Type>&) const;
725 
726 
727  // Member Operators
728 
729  virtual void operator=(const UList<Type>&);
730 
731  virtual void operator=(const fvPatchField<Type>&);
732  virtual void operator+=(const fvPatchField<Type>&);
733  virtual void operator-=(const fvPatchField<Type>&);
734  virtual void operator*=(const fvPatchField<scalar>&);
735  virtual void operator/=(const fvPatchField<scalar>&);
736 
737  virtual void operator+=(const Field<Type>&);
738  virtual void operator-=(const Field<Type>&);
739 
740  virtual void operator*=(const Field<scalar>&);
741  virtual void operator/=(const Field<scalar>&);
742 
743  virtual void operator=(const Type&);
744  virtual void operator+=(const Type&);
745  virtual void operator-=(const Type&);
746  virtual void operator*=(const scalar);
747  virtual void operator/=(const scalar);
749 
750  // Force an assignment irrespective of form of patch
751 
752  virtual void operator==(const fvPatchField<Type>&);
753  virtual void operator==(const Field<Type>&);
754  virtual void operator==(const Type&);
755 
756  // Prevent automatic comparison rewriting (c++20)
757  bool operator!=(const fvPatchField<Type>&) const = delete;
758  bool operator!=(const Field<Type>&) const = delete;
759  bool operator!=(const Type&) const = delete;
760 
761 
762  // Ostream Operator
763 
764  friend Ostream& operator<< <Type>(Ostream&, const fvPatchField<Type>&);
765 };
766 
767 
768 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
769 
770 } // End namespace Foam
771 
772 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
773 
774 #ifdef NoRepository
775  #include "fvPatchField.C"
776  #include "fvPatchFieldNew.C"
777  #include "calculatedFvPatchField.H"
778  #include "zeroGradientFvPatchField.H"
779 #endif
780 
781 // Runtime selection macros
782 #include "fvPatchFieldMacros.H"
783 
784 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
785 
786 #endif
787 
788 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
Definition: fvPatchField.C:32
dictionary dict
virtual bool coupled() const
True if the patch field is coupled.
Definition: fvPatchField.H:253
const objectRegistry & db() const
The associated objectRegistry.
"blocking" : (MPI_Bsend, MPI_Recv)
uint8_t direction
Definition: direction.H:46
const word zeroGradientType
A zeroGradient patch field type.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:236
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:213
declareRunTimeSelectionTable(tmp, fvPatchField, patch,(const fvPatch &p, const DimensionedField< Type, volMesh > &iF),(p, iF))
commsTypes
Communications types.
Definition: UPstream.H:72
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
static const word & emptyType() noexcept
The type name for empty patch fields.
Definition: fvPatchField.H:196
virtual ~fvPatchField()=default
Destructor.
virtual void operator*=(const fvPatchField< scalar > &)
Definition: fvPatchField.C:434
virtual void initEvaluateLocal(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Initialise the evaluation of the patch field after a local.
Definition: fvPatchField.H:779
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
void extrapolateInternal()
Assign the patch field from the internal field.
Definition: fvPatchField.C:62
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:720
virtual void initEvaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Initialise the evaluation of the patch field.
Definition: fvPatchField.H:758
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< Field< scalar >> &) const
Return the matrix source coefficients corresponding to the.
Definition: fvPatchField.H:815
Template invariant parts for fvPatchField.
Definition: fvPatchField.H:77
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
virtual void readDict(const dictionary &dict)
Read dictionary entries.
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:546
bool useImplicit() const noexcept
Use implicit formulation for coupled patches only.
Definition: fvPatchField.H:312
const word calculatedType
A calculated patch field type.
bool manipulatedMatrix() const noexcept
True if the matrix has already been manipulated.
Definition: fvPatchField.H:304
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:45
Macros for creating fvPatchField types.
bool updated() const noexcept
True if the boundary condition has already been updated.
Definition: fvPatchField.H:296
Generic templated field type.
Definition: Field.H:62
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:748
bool operator!=(const fvPatchField< Type > &) const =delete
fvPatch Patch
The patch type for the patch field.
Definition: fvPatchField.H:396
virtual void operator/=(const fvPatchField< scalar > &)
Definition: fvPatchField.C:445
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
Definition: fvPatchField.H:828
A FieldMapper for finite-volume patch fields.
void setUpdated(bool state) noexcept
Set updated state.
Definition: fvPatchField.H:130
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:299
virtual bool assignable() const
True if the value of the patch field is altered by assignment.
Definition: fvPatchField.H:234
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:342
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< Field< scalar >> &) const
Return the matrix diagonal coefficients corresponding to the.
Definition: fvPatchField.H:801
const Field< Type > & primitiveField() const noexcept
Return const-reference to the internal field values.
Definition: fvPatchField.H:670
void check(const fvPatchField< Type > &) const
Check against given patch field.
Definition: fvPatchField.C:206
fvPatchFieldBase(const fvPatch &p)
Construct from patch.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
virtual void operator-=(const fvPatchField< Type > &)
Definition: fvPatchField.C:423
static int disallowGenericPatchField
Debug switch to disallow the use of generic fvPatchField.
Definition: fvPatchField.H:149
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
Definition: fvPatchField.C:221
const word emptyType
An empty patch field type.
OBJstream os(runTime.globalPath()/outputName)
virtual void operator+=(const fvPatchField< Type > &)
Definition: fvPatchField.C:412
virtual void evaluateLocal(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field after a local operation (e.g. *=)
Definition: fvPatchField.H:789
virtual ~fvPatchFieldBase()=default
Destructor.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets updated() to false.
Definition: fvPatchField.C:329
const word & patchType() const noexcept
The optional patch type.
Definition: fvPatchField.H:277
calculatedFvPatchField< Type > Calculated
Type for a calculated patch.
Definition: fvPatchField.H:401
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
Definition: fvPatchField.H:854
virtual bool fixesValue() const
True if the patch field fixes a value.
Definition: fvPatchField.H:245
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:391
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
static tmp< fvPatchField< Type > > NewCalculatedType(const fvPatch &)
Return a pointer to a new calculatedFvPatchField created on.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:204
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
Definition: fvPatchField.H:662
Reading is optional [identical to READ_IF_PRESENT].
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:316
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void checkPatch(const fvPatchFieldBase &rhs) const
Check that patches are identical.
Registry of regIOobjects.
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Definition: fvPatchField.H:213
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
TypeName("fvPatchField")
Runtime type information.
const word extrapolatedCalculatedType
A combined zero-gradient and calculated patch field type.
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: fvPatchField.C:72
void setManipulated(bool state) noexcept
Set matrix manipulated state.
Definition: fvPatchField.H:138
Namespace for OpenFOAM.
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
virtual tmp< fvPatchField< Type > > clone() const
Construct and return a clone.
Definition: fvPatchField.H:555
DimensionedField< Type, volMesh > Internal
The internal field type associated with the patch field.
Definition: fvPatchField.H:391
readOption
Enumeration defining read preferences.