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-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 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;
64 class fvPatchFieldMapper;
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>
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  //- Debug switch to disallow the use of generic fvPatchField
117  static int disallowGenericPatchField;
118 
119  //- Runtime type information
120  TypeName("fvPatchField");
121 
122 
123  // Constructors
124 
125  //- Construct from patch
126  explicit fvPatchFieldBase(const fvPatch& p);
127 
128  //- Construct from patch and patch type
129  explicit fvPatchFieldBase(const fvPatch& p, const word& patchType);
130 
131  //- Construct from patch and dictionary
132  fvPatchFieldBase(const fvPatch& p, const dictionary& dict);
133 
134  //- Copy construct with new patch
136 
137  //- Copy construct
139 
140 
141  //- Destructor
142  virtual ~fvPatchFieldBase() = default;
143 
144 
145  // Static Member Functions
146 
147  //- The type name for \c empty patch fields
148  static const word& emptyType() noexcept
149  {
151  }
152 
153  //- The type name for \c calculated patch fields
154  static const word& calculatedType() noexcept
155  {
157  }
158 
159  //- The type name for \c extrapolatedCalculated patch fields
160  //- combines \c zero-gradient and \c calculated
161  static const word& extrapolatedCalculatedType() noexcept
162  {
164  }
165 
166  //- The type name for \c zeroGradient patch fields
167  static const word& zeroGradientType() noexcept
168  {
170  }
171 
172 
173  // Member Functions
174 
175  // Attributes
176 
177  //- True if the value of the patch field is altered by assignment
178  virtual bool assignable() const
179  {
180  return true;
181  }
182 
183  //- True if the patch field fixes a value.
184  // Needed to check if a level has to be specified while solving
185  // Poissons equations.
186  virtual bool fixesValue() const
187  {
188  return false;
189  }
190 
191  //- True if the patch field is coupled
192  virtual bool coupled() const
193  {
194  return false;
195  }
196 
197 
198  // Access
199 
200  //- The associated objectRegistry
201  const objectRegistry& db() const;
202 
203  //- Return the patch
204  const fvPatch& patch() const noexcept
205  {
206  return patch_;
207  }
208 
209  //- The optional patch type
210  const word& patchType() const noexcept
211  {
212  return patchType_;
213  }
214 
215  //- The optional patch type
217  {
218  return patchType_;
219  }
220 
221  //- True if the type does not correspond to the constraint type
222  virtual bool constraintOverride() const
223  {
224  return !patchType_.empty() && patchType_ != type();
225  }
226 
227 
228  // Solution
229 
230  //- True if the boundary condition has already been updated
231  bool updated() const noexcept
232  {
233  return updated_;
234  }
235 
236  //- Set updated state
237  void setUpdated(bool state) noexcept
238  {
239  updated_ = state;
240  }
241 
242  //- True if the matrix has already been manipulated
243  bool manipulatedMatrix() const noexcept
244  {
245  return manipulatedMatrix_;
246  }
247 
248  //- Set matrix manipulated state
249  void setManipulated(bool state) noexcept
250  {
251  manipulatedMatrix_ = state;
252  }
253 
254  //- Use implicit formulation for coupled patches only
255  bool useImplicit() const noexcept
256  {
257  return useImplicit_;
258  }
259 
260  //- Set useImplicit on/off
261  // \return old value
262  bool useImplicit(bool on) noexcept
263  {
264  bool old(useImplicit_);
265  useImplicit_ = on;
266  return old;
267  }
268 
269 
270  // Check
271 
272  //- Check that patches are identical
273  void checkPatch(const fvPatchFieldBase& rhs) const;
274 };
275 
276 
277 /*---------------------------------------------------------------------------*\
278  Class fvPatchField Declaration
279 \*---------------------------------------------------------------------------*/
280 
281 template<class Type>
282 class fvPatchField
283 :
284  public fvPatchFieldBase,
285  public Field<Type>
286 {
287 public:
288 
289  // Public Data Types
290 
291  //- The patch type for the patch field
292  typedef fvPatch Patch;
293 
294  //- The value_type for the patch field
295  typedef Type value_type;
296 
297  //- The component type for patch field
298  typedef typename pTraits<Type>::cmptType cmptType;
299 
300  //- The internal field type associated with the patch field
302 
303  //- Type for a \em calculated patch
305 
306 
307 private:
308 
309  // Private Data
310 
311  //- Reference to internal field
312  const DimensionedField<Type, volMesh>& internalField_;
313 
314 
315 protected:
316 
317  // Protected Member Functions
318 
319  //- Read the "value" entry into \c *this.
320  // The reading can be optional (default), mandatory etc.
321  // \returns True on success
322  bool readValueEntry
323  (
324  const dictionary& dict,
326  );
327 
328  //- Write \c *this field as a "value" entry
329  void writeValueEntry(Ostream& os) const
330  {
331  Field<Type>::writeEntry("value", os);
332  }
333 
334  //- Assign the patch field from the internal field
335  void extrapolateInternal();
336 
337 
338 public:
339 
340  // Declare run-time constructor selection tables
341 
343  (
344  tmp,
345  fvPatchField,
346  patch,
347  (
348  const fvPatch& p,
350  ),
351  (p, iF)
352  );
353 
355  (
356  tmp,
357  fvPatchField,
358  patchMapper,
359  (
360  const fvPatchField<Type>& ptf,
361  const fvPatch& p,
363  const fvPatchFieldMapper& m
364  ),
365  (dynamic_cast<const fvPatchFieldType&>(ptf), p, iF, m)
366  );
367 
369  (
370  tmp,
372  dictionary,
373  (
374  const fvPatch& p,
377  ),
378  (p, iF, dict)
379  );
380 
381 
382  // Constructors
383 
384  //- Construct from patch and internal field
386  (
387  const fvPatch&,
389  );
390 
391  //- Construct from patch, internal field and patch type
393  (
394  const fvPatch&,
396  const word& patchType
397  );
398 
399  //- Construct from patch, internal field and value
401  (
402  const fvPatch&,
404  const Type& value
405  );
406 
407  //- Construct from patch, internal field and patch field
409  (
410  const fvPatch&,
412  const Field<Type>& pfld
413  );
414 
415  //- Construct from patch, internal field and patch field
417  (
418  const fvPatch&,
420  Field<Type>&& pfld
421  );
422 
423  //- Construct from patch, internal field and dictionary
425  (
426  const fvPatch&,
428  const dictionary& dict,
431  );
432 
433  //- Construct, forwarding to readOption variant
435  (
436  const fvPatch& p,
438  const dictionary& dict,
439  const bool needValue
440  )
441  :
443  (
444  p, iF, dict,
445  (needValue? IOobjectOption::MUST_READ : IOobjectOption::NO_READ)
446  )
447  {}
448 
449  //- Construct by mapping the given fvPatchField onto a new patch
451  (
452  const fvPatchField<Type>&,
453  const fvPatch&,
455  const fvPatchFieldMapper&
456  );
457 
458  //- Construct as copy
460 
461  //- Construct as copy setting internal field reference
463  (
464  const fvPatchField<Type>&,
466  );
467 
468  //- Clone patch field with its own internal field reference
469  virtual tmp<fvPatchField<Type>> clone() const
470  {
471  return tmp<fvPatchField<Type>>
472  (
473  new fvPatchField<Type>(*this, this->internalField_)
474  );
475  }
476 
477  //- Clone patch field with given internal field reference
479  (
481  ) const
482  {
483  return tmp<fvPatchField<Type>>
484  (
485  new fvPatchField<Type>(*this, iF)
486  );
487  }
488 
489 
490  // Factory Methods
491 
492  //- Clone a patch field, optionally with internal field reference etc.
493  template<class DerivedPatchField, class... Args>
494  static tmp<fvPatchField<Type>> Clone
495  (
496  const DerivedPatchField& pf,
497  Args&&... args
498  )
499  {
500  return tmp<fvPatchField<Type>>
501  (
502  new DerivedPatchField(pf, std::forward<Args>(args)...)
503  );
504  }
505 
506  //- Return a pointer to a new patchField created on freestore given
507  // patch and internal field
508  // (does not set the patch field values)
509  static tmp<fvPatchField<Type>> New
510  (
511  const word& patchFieldType,
512  const fvPatch&,
513  const DimensionedField<Type, volMesh>&
514  );
515 
516  //- Return a pointer to a new patchField created on freestore given
517  // patch and internal field
518  // (does not set the patch field values).
519  // Allows override of constraint type
520  static tmp<fvPatchField<Type>> New
521  (
522  const word& patchFieldType,
523  const word& actualPatchType,
524  const fvPatch&,
525  const DimensionedField<Type, volMesh>&
526  );
527 
528  //- Return a pointer to a new patchField created on freestore from
529  // a given fvPatchField mapped onto a new patch
530  static tmp<fvPatchField<Type>> New
531  (
532  const fvPatchField<Type>&,
533  const fvPatch&,
535  const fvPatchFieldMapper&
536  );
537 
538  //- Return a pointer to a new patchField created on freestore
539  // from dictionary
541  (
542  const fvPatch&,
544  const dictionary&
545  );
546 
547  //- Return a pointer to a new calculatedFvPatchField created on
548  // freestore without setting patchField values
550  (
551  const fvPatch& p
552  );
553 
554  //- Return a pointer to a new calculatedFvPatchField created on
555  // freestore without setting patchField values
556  template<class AnyType>
558  (
559  const fvPatchField<AnyType>& pf
560  );
561 
562 
563  //- Destructor
564  virtual ~fvPatchField() = default;
565 
566 
567  // Member Functions
568 
569  // Access
570 
571  //- Return const-reference to the dimensioned internal field
573  {
574  return internalField_;
575  }
576 
577  //- Return const-reference to the internal field values
578  const Field<Type>& primitiveField() const noexcept
579  {
580  return internalField_;
581  }
582 
583 
584  // Mapping Functions
585 
586  //- Map (and resize as needed) from self given a mapping object
587  virtual void autoMap
588  (
589  const fvPatchFieldMapper&
590  );
591 
592  //- Reverse map the given fvPatchField onto this fvPatchField
593  virtual void rmap
594  (
595  const fvPatchField<Type>&,
596  const labelList&
597  );
598 
599 
600  // Evaluation Functions
601 
602  //- Return patch-normal gradient
603  virtual tmp<Field<Type>> snGrad() const;
604 
605  //- Return patch-normal gradient for coupled-patches
606  // using the deltaCoeffs provided
607  virtual tmp<Field<Type>> snGrad
608  (
609  const scalarField& deltaCoeffs
610  ) const
611  {
613  return *this;
614  }
615 
616  //- Update the coefficients associated with the patch field
617  // Sets Updated to true
618  virtual void updateCoeffs();
619 
620  //- Update the coefficients associated with the patch field
621  // with a weight field (0..1). This weight field is usually
622  // provided as the amount of geometric overlap for 'duplicate'
623  // patches. Sets Updated to true
624  virtual void updateWeightedCoeffs(const scalarField& weights);
625 
626  //- Return internal field next to patch
627  virtual tmp<Field<Type>> patchInternalField() const;
628 
629  //- Retrieve internal field next to patch
630  // \param [out] pfld The extracted patch field.
631  virtual void patchInternalField(UList<Type>& pfld) const;
632 
633  //- Return patchField on the opposite patch of a coupled patch
634  virtual tmp<Field<Type>> patchNeighbourField() const
635  {
637  return *this;
638  }
639 
640  //- Initialise the evaluation of the patch field
641  virtual void initEvaluate
642  (
643  const Pstream::commsTypes commsType =
645  )
646  {}
647 
648  //- Evaluate the patch field, sets updated() to false
649  virtual void evaluate
650  (
651  const Pstream::commsTypes commsType =
653  );
654 
655  //- Initialise the evaluation of the patch field after a local
656  // operation
657  virtual void initEvaluateLocal
658  (
659  const Pstream::commsTypes commsType =
661  )
662  {}
663 
664  //- Evaluate the patch field after a local operation (e.g. *=)
665  virtual void evaluateLocal
666  (
667  const Pstream::commsTypes commsType =
669  )
670  {}
671 
672  //- Return the matrix diagonal coefficients corresponding to the
673  // evaluation of the value of this patchField with given weights
674  virtual tmp<Field<Type>> valueInternalCoeffs
675  (
676  const tmp<Field<scalar>>&
677  ) const
678  {
680  return *this;
681  }
682 
683  //- Return the matrix source coefficients corresponding to the
684  // evaluation of the value of this patchField with given weights
685  virtual tmp<Field<Type>> valueBoundaryCoeffs
686  (
688  ) const
689  {
691  return *this;
692  }
693 
694  //- Return the matrix diagonal coefficients corresponding to the
695  // evaluation of the gradient of this patchField
697  {
699  return *this;
700  }
701 
702  //- Return the matrix diagonal coefficients corresponding to the
703  // evaluation of the gradient of this coupled patchField
704  // using the deltaCoeffs provided
706  (
707  const scalarField& deltaCoeffs
708  ) const
709  {
711  return *this;
712  }
713 
714  //- Return the matrix source coefficients corresponding to the
715  // evaluation of the gradient of this patchField
716  virtual tmp<Field<Type>> gradientBoundaryCoeffs() const
717  {
719  return *this;
720  }
721 
722  //- Return the matrix source coefficients corresponding to the
723  // evaluation of the gradient of this coupled patchField
724  // using the deltaCoeffs provided
726  (
727  const scalarField& deltaCoeffs
728  ) const
729  {
731  return *this;
732  }
733 
734 
735  //- Manipulate matrix
736  virtual void manipulateMatrix(fvMatrix<Type>& matrix);
737 
738  //- Manipulate matrix with given weights
739  virtual void manipulateMatrix
740  (
741  fvMatrix<Type>& matrix,
742  const scalarField& weights
743  );
744 
745  //- Manipulate fvMatrix
746  virtual void manipulateMatrix
747  (
748  fvMatrix<Type>& matrix,
749  const label iMatrix,
750  const direction cmp
751  );
752 
753 
754  // Other
755 
756  //- Write
757  virtual void write(Ostream&) const;
758 
759  //- Check against given patch field
760  void check(const fvPatchField<Type>&) const;
761 
762 
763  // Member Operators
764 
765  virtual void operator=(const UList<Type>&);
766 
767  virtual void operator=(const fvPatchField<Type>&);
768  virtual void operator+=(const fvPatchField<Type>&);
769  virtual void operator-=(const fvPatchField<Type>&);
770  virtual void operator*=(const fvPatchField<scalar>&);
771  virtual void operator/=(const fvPatchField<scalar>&);
772 
773  virtual void operator+=(const Field<Type>&);
774  virtual void operator-=(const Field<Type>&);
775 
776  virtual void operator*=(const Field<scalar>&);
777  virtual void operator/=(const Field<scalar>&);
778 
779  virtual void operator=(const Type&);
780  virtual void operator+=(const Type&);
781  virtual void operator-=(const Type&);
782  virtual void operator*=(const scalar);
783  virtual void operator/=(const scalar);
784 
785 
786  // Force an assignment irrespective of form of patch
787 
788  virtual void operator==(const fvPatchField<Type>&);
789  virtual void operator==(const Field<Type>&);
790  virtual void operator==(const Type&);
791 
792  // Prevent automatic comparison rewriting (c++20)
793  bool operator!=(const fvPatchField<Type>&) const = delete;
794  bool operator!=(const Field<Type>&) const = delete;
795  bool operator!=(const Type&) const = delete;
796 
797 
798  // Ostream Operator
799 
800  friend Ostream& operator<< <Type>(Ostream&, const fvPatchField<Type>&);
801 };
802 
803 
804 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
805 
806 } // End namespace Foam
807 
808 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
809 
810 #ifdef NoRepository
811  #include "fvPatchField.C"
812  #include "fvPatchFieldNew.C"
813  #include "calculatedFvPatchField.H"
814  #include "zeroGradientFvPatchField.H"
815 #endif
816 
817 // Runtime selection macros
818 #include "fvPatchFieldMacros.H"
819 
820 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
821 
822 #endif
823 
824 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type.
Definition: Field.H:166
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:754
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
readOption
Enumeration defining read preferences.
@ MUST_READ
Reading required.
@ LAZY_READ
Reading is optional [identical to READ_IF_PRESENT].
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:59
commsTypes
Communications types.
Definition: UPstream.H:78
@ buffered
"buffered" : (MPI_Bsend, MPI_Recv)
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
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
fvPatchFieldBase(const fvPatch &p)
Construct from patch.
virtual ~fvPatchFieldBase()=default
Destructor.
virtual bool fixesValue() const
True if the patch field fixes a value.
Definition: fvPatchField.H:226
void checkPatch(const fvPatchFieldBase &rhs) const
Check that patches are identical.
virtual bool coupled() const
True if the patch field is coupled.
Definition: fvPatchField.H:234
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:250
const objectRegistry & db() const
The associated objectRegistry.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:185
virtual void readDict(const dictionary &dict)
Read dictionary entries.
static const word & emptyType() noexcept
The type name for empty patch fields.
Definition: fvPatchField.H:177
void setManipulated(bool state) noexcept
Set matrix manipulated state.
Definition: fvPatchField.H:309
bool useImplicit() const noexcept
Use implicit formulation for coupled patches only.
Definition: fvPatchField.H:317
void setUpdated(bool state) noexcept
Set updated state.
Definition: fvPatchField.H:293
static int disallowGenericPatchField
Debug switch to disallow the use of generic fvPatchField.
Definition: fvPatchField.H:130
TypeName("fvPatchField")
Runtime type information.
const word & patchType() const noexcept
The optional patch type.
Definition: fvPatchField.H:258
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Definition: fvPatchField.H:194
bool manipulatedMatrix() const noexcept
True if the matrix has already been manipulated.
Definition: fvPatchField.H:301
virtual bool constraintOverride() const
True if the type does not correspond to the constraint type.
Definition: fvPatchField.H:274
virtual bool assignable() const
True if the value of the patch field is altered by assignment.
Definition: fvPatchField.H:215
bool updated() const noexcept
True if the boundary condition has already been updated.
Definition: fvPatchField.H:285
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:202
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 tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
Definition: fvPatchField.H:898
DimensionedField< Type, volMesh > Internal
The internal field type associated with the patch field.
Definition: fvPatchField.H:376
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 void write(Ostream &) const
Write.
Definition: fvPatchField.C:374
Type value_type
The value_type for the patch field.
Definition: fvPatchField.H:366
virtual void initEvaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Initialise the evaluation of the patch field.
Definition: fvPatchField.H:802
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:548
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
Definition: fvPatchField.H:706
declareRunTimeSelectionTable(tmp, fvPatchField, patch,(const fvPatch &p, const DimensionedField< Type, volMesh > &iF),(p, iF))
virtual tmp< fvPatchField< Type > > clone() const
Clone patch field with its own internal field reference.
Definition: fvPatchField.H:577
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
Definition: fvPatchField.C:223
fvPatch Patch
The patch type for the patch field.
Definition: fvPatchField.H:361
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:792
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:413
static tmp< fvPatchField< Type > > NewCalculatedType(const fvPatch &p)
Return a pointer to a new calculatedFvPatchField created on.
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< Field< scalar >> &) const
Return the matrix source coefficients corresponding to the.
Definition: fvPatchField.H:859
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:215
const Field< Type > & primitiveField() const noexcept
Return const-reference to the internal field values.
Definition: fvPatchField.H:714
virtual void initEvaluateLocal(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Initialise the evaluation of the patch field after a local.
Definition: fvPatchField.H:823
virtual void operator-=(const fvPatchField< Type > &)
Definition: fvPatchField.C:425
virtual void evaluateLocal(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field after a local operation (e.g. *=)
Definition: fvPatchField.H:833
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:393
static tmp< fvPatchField< Type > > Clone(const DerivedPatchField &pf, Args &&... args)
Clone a patch field, optionally with internal field reference etc.
Definition: fvPatchField.H:607
pTraits< Type >::cmptType cmptType
The component type for patch field.
Definition: fvPatchField.H:371
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:311
virtual ~fvPatchField()=default
Destructor.
calculatedFvPatchField< Type > Calculated
Type for a calculated patch.
Definition: fvPatchField.H:381
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:344
virtual void operator*=(const fvPatchField< scalar > &)
Definition: fvPatchField.C:436
bool operator!=(const fvPatchField< Type > &) const =delete
virtual void operator+=(const fvPatchField< Type > &)
Definition: fvPatchField.C:414
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
Definition: fvPatchField.H:872
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
virtual void operator/=(const fvPatchField< scalar > &)
Definition: fvPatchField.C:447
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< Field< scalar >> &) const
Return the matrix diagonal coefficients corresponding to the.
Definition: fvPatchField.H:845
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
Registry of regIOobjects.
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:64
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 NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:688
OBJstream os(runTime.globalPath()/outputName)
Macros for creating fvPatchField types.
const word zeroGradientType
A zeroGradient patch field type.
const word extrapolatedCalculatedType
A combined zero-gradient and calculated patch field type.
const word calculatedType
A calculated patch field type.
const word emptyType
An empty patch field type.
Namespace for OpenFOAM.
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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
const direction noexcept
Definition: scalarImpl.H:255
dictionary dict
Foam::argList args(argc, argv)