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-2024 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 public:
286  // Public Data Types
287 
288  //- The patch type for the patch field
289  typedef fvPatch Patch;
290 
291  //- The value_type for the patch field
292  typedef Type value_type;
293 
294  //- The internal field type associated with the patch field
297  //- Type for a \em calculated patch
299 
300 
301 private:
302 
303  // Private Data
305  //- Reference to internal field
306  const DimensionedField<Type, volMesh>& internalField_;
307 
308 
309 protected:
310 
311  // Protected Member Functions
313  //- Read the "value" entry into \c *this.
314  // The reading can be optional (default), mandatory etc.
315  // \returns True on success
316  bool readValueEntry
317  (
318  const dictionary& dict,
320  );
321 
322  //- Write \c *this field as a "value" entry
323  void writeValueEntry(Ostream& os) const
324  {
325  Field<Type>::writeEntry("value", os);
326  }
327 
328  //- Assign the patch field from the internal field
329  void extrapolateInternal();
330 
331 
332 public:
333 
334  // Declare run-time constructor selection tables
335 
337  (
338  tmp,
339  fvPatchField,
340  patch,
341  (
342  const fvPatch& p,
344  ),
345  (p, iF)
346  );
347 
349  (
350  tmp,
351  fvPatchField,
352  patchMapper,
353  (
354  const fvPatchField<Type>& ptf,
355  const fvPatch& p,
357  const fvPatchFieldMapper& m
358  ),
359  (dynamic_cast<const fvPatchFieldType&>(ptf), p, iF, m)
360  );
363  (
364  tmp,
365  fvPatchField,
367  (
368  const fvPatch& p,
370  const dictionary& dict
371  ),
372  (p, iF, dict)
373  );
374 
375 
376  // Constructors
377 
378  //- Construct from patch and internal field
380  (
381  const fvPatch&,
383  );
384 
385  //- Construct from patch, internal field and patch type
387  (
388  const fvPatch&,
390  const word& patchType
391  );
392 
393  //- Construct from patch, internal field and value
395  (
396  const fvPatch&,
398  const Type& value
399  );
400 
401  //- Construct from patch, internal field and patch field
403  (
404  const fvPatch&,
406  const Field<Type>& pfld
407  );
408 
409  //- Construct from patch, internal field and patch field
411  (
412  const fvPatch&,
414  Field<Type>&& pfld
415  );
416 
417  //- Construct from patch, internal field and dictionary
419  (
420  const fvPatch&,
422  const dictionary& dict,
425  );
426 
427  //- Construct, forwarding to readOption variant
429  (
430  const fvPatch& p,
432  const dictionary& dict,
433  const bool needValue
434  )
435  :
437  (
438  p, iF, dict,
439  (needValue? IOobjectOption::MUST_READ : IOobjectOption::NO_READ)
440  )
441  {}
442 
443  //- Construct by mapping the given fvPatchField onto a new patch
445  (
446  const fvPatchField<Type>&,
447  const fvPatch&,
449  const fvPatchFieldMapper&
450  );
451 
452  //- Construct as copy
454 
455  //- Construct as copy setting internal field reference
457  (
458  const fvPatchField<Type>&,
460  );
461 
462  //- Clone patch field with its own internal field reference
463  virtual tmp<fvPatchField<Type>> clone() const
464  {
465  return tmp<fvPatchField<Type>>
466  (
467  new fvPatchField<Type>(*this, this->internalField_)
468  );
469  }
470 
471  //- Clone patch field with given internal field reference
472  virtual tmp<fvPatchField<Type>> clone
473  (
474  const DimensionedField<Type, volMesh>& iF
475  ) const
476  {
477  return tmp<fvPatchField<Type>>
478  (
479  new fvPatchField<Type>(*this, iF)
480  );
481  }
482 
483 
484  // Factory Methods
485 
486  //- Clone a patch field, optionally with internal field reference etc.
487  template<class DerivedPatchField, class... Args>
488  static tmp<fvPatchField<Type>> Clone
489  (
490  const DerivedPatchField& pf,
491  Args&&... args
492  )
493  {
494  return tmp<fvPatchField<Type>>
495  (
496  new DerivedPatchField(pf, std::forward<Args>(args)...)
497  );
498  }
499 
500  //- Return a pointer to a new patchField created on freestore given
501  // patch and internal field
502  // (does not set the patch field values)
503  static tmp<fvPatchField<Type>> New
504  (
505  const word& patchFieldType,
506  const fvPatch&,
507  const DimensionedField<Type, volMesh>&
508  );
509 
510  //- Return a pointer to a new patchField created on freestore given
511  // patch and internal field
512  // (does not set the patch field values).
513  // Allows override of constraint type
514  static tmp<fvPatchField<Type>> New
515  (
516  const word& patchFieldType,
517  const word& actualPatchType,
518  const fvPatch&,
519  const DimensionedField<Type, volMesh>&
520  );
521 
522  //- Return a pointer to a new patchField created on freestore from
523  // a given fvPatchField mapped onto a new patch
524  static tmp<fvPatchField<Type>> New
525  (
526  const fvPatchField<Type>&,
527  const fvPatch&,
529  const fvPatchFieldMapper&
530  );
531 
532  //- Return a pointer to a new patchField created on freestore
533  // from dictionary
535  (
536  const fvPatch&,
538  const dictionary&
539  );
540 
541  //- Return a pointer to a new calculatedFvPatchField created on
542  // freestore without setting patchField values
544  (
545  const fvPatch& p
546  );
547 
548  //- Return a pointer to a new calculatedFvPatchField created on
549  // freestore without setting patchField values
550  template<class AnyType>
552  (
553  const fvPatchField<AnyType>& pf
554  );
555 
556 
557  //- Destructor
558  virtual ~fvPatchField() = default;
559 
560 
561  // Member Functions
562 
563  // Access
564 
565  //- Return const-reference to the dimensioned internal field
567  {
568  return internalField_;
569  }
570 
571  //- Return const-reference to the internal field values
572  const Field<Type>& primitiveField() const noexcept
573  {
574  return internalField_;
575  }
576 
577 
578  // Mapping Functions
580  //- Map (and resize as needed) from self given a mapping object
581  virtual void autoMap
582  (
583  const fvPatchFieldMapper&
584  );
585 
586  //- Reverse map the given fvPatchField onto this fvPatchField
587  virtual void rmap
588  (
589  const fvPatchField<Type>&,
590  const labelList&
591  );
592 
593 
594  // Evaluation Functions
595 
596  //- Return patch-normal gradient
597  virtual tmp<Field<Type>> snGrad() const;
598 
599  //- Return patch-normal gradient for coupled-patches
600  // using the deltaCoeffs provided
601  virtual tmp<Field<Type>> snGrad
602  (
603  const scalarField& deltaCoeffs
604  ) const
605  {
607  return *this;
608  }
609 
610  //- Update the coefficients associated with the patch field
611  // Sets Updated to true
612  virtual void updateCoeffs();
613 
614  //- Update the coefficients associated with the patch field
615  // with a weight field (0..1). This weight field is usually
616  // provided as the amount of geometric overlap for 'duplicate'
617  // patches. Sets Updated to true
618  virtual void updateWeightedCoeffs(const scalarField& weights);
619 
620  //- Return internal field next to patch
621  virtual tmp<Field<Type>> patchInternalField() const;
622 
623  //- Extract internal field next to patch
624  // \param [out] pfld The extracted patch field.
625  virtual void patchInternalField(Field<Type>& pfld) const;
626 
627  //- Return patchField on the opposite patch of a coupled patch
628  virtual tmp<Field<Type>> patchNeighbourField() const
629  {
631  return *this;
632  }
633 
634  //- Initialise the evaluation of the patch field
635  virtual void initEvaluate
636  (
637  const Pstream::commsTypes commsType =
639  )
640  {}
641 
642  //- Evaluate the patch field, sets updated() to false
643  virtual void evaluate
644  (
645  const Pstream::commsTypes commsType =
647  );
648 
649  //- Initialise the evaluation of the patch field after a local
650  // operation
651  virtual void initEvaluateLocal
652  (
653  const Pstream::commsTypes commsType =
655  )
656  {}
657 
658  //- Evaluate the patch field after a local operation (e.g. *=)
659  virtual void evaluateLocal
660  (
661  const Pstream::commsTypes commsType =
663  )
664  {}
665 
666  //- Return the matrix diagonal coefficients corresponding to the
667  // evaluation of the value of this patchField with given weights
668  virtual tmp<Field<Type>> valueInternalCoeffs
669  (
670  const tmp<Field<scalar>>&
671  ) const
672  {
674  return *this;
675  }
676 
677  //- Return the matrix source coefficients corresponding to the
678  // evaluation of the value of this patchField with given weights
679  virtual tmp<Field<Type>> valueBoundaryCoeffs
680  (
681  const tmp<Field<scalar>>&
682  ) const
683  {
685  return *this;
686  }
687 
688  //- Return the matrix diagonal coefficients corresponding to the
689  // evaluation of the gradient of this patchField
690  virtual tmp<Field<Type>> gradientInternalCoeffs() const
691  {
693  return *this;
694  }
695 
696  //- Return the matrix diagonal coefficients corresponding to the
697  // evaluation of the gradient of this coupled patchField
698  // using the deltaCoeffs provided
700  (
701  const scalarField& deltaCoeffs
702  ) const
703  {
705  return *this;
706  }
707 
708  //- Return the matrix source coefficients corresponding to the
709  // evaluation of the gradient of this patchField
710  virtual tmp<Field<Type>> gradientBoundaryCoeffs() const
711  {
713  return *this;
714  }
715 
716  //- Return the matrix source coefficients corresponding to the
717  // evaluation of the gradient of this coupled patchField
718  // using the deltaCoeffs provided
719  virtual tmp<Field<Type>> gradientBoundaryCoeffs
720  (
721  const scalarField& deltaCoeffs
722  ) const
723  {
725  return *this;
726  }
727 
728 
729  //- Manipulate matrix
730  virtual void manipulateMatrix(fvMatrix<Type>& matrix);
731 
732  //- Manipulate matrix with given weights
733  virtual void manipulateMatrix
734  (
735  fvMatrix<Type>& matrix,
736  const scalarField& weights
737  );
738 
739  //- Manipulate fvMatrix
740  virtual void manipulateMatrix
741  (
742  fvMatrix<Type>& matrix,
743  const label iMatrix,
744  const direction cmp
745  );
746 
747 
748  // Other
749 
750  //- Write
751  virtual void write(Ostream&) const;
752 
753  //- Check against given patch field
754  void check(const fvPatchField<Type>&) const;
755 
756 
757  // Member Operators
758 
759  virtual void operator=(const UList<Type>&);
760 
761  virtual void operator=(const fvPatchField<Type>&);
762  virtual void operator+=(const fvPatchField<Type>&);
763  virtual void operator-=(const fvPatchField<Type>&);
764  virtual void operator*=(const fvPatchField<scalar>&);
765  virtual void operator/=(const fvPatchField<scalar>&);
766 
767  virtual void operator+=(const Field<Type>&);
768  virtual void operator-=(const Field<Type>&);
769 
770  virtual void operator*=(const Field<scalar>&);
771  virtual void operator/=(const Field<scalar>&);
772 
773  virtual void operator=(const Type&);
774  virtual void operator+=(const Type&);
775  virtual void operator-=(const Type&);
776  virtual void operator*=(const scalar);
777  virtual void operator/=(const scalar);
778 
779 
780  // Force an assignment irrespective of form of patch
781 
782  virtual void operator==(const fvPatchField<Type>&);
783  virtual void operator==(const Field<Type>&);
784  virtual void operator==(const Type&);
785 
786  // Prevent automatic comparison rewriting (c++20)
787  bool operator!=(const fvPatchField<Type>&) const = delete;
788  bool operator!=(const Field<Type>&) const = delete;
789  bool operator!=(const Type&) const = delete;
790 
791 
792  // Ostream Operator
793 
794  friend Ostream& operator<< <Type>(Ostream&, const fvPatchField<Type>&);
795 };
796 
797 
798 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
799 
800 } // End namespace Foam
801 
802 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
803 
804 #ifdef NoRepository
805  #include "fvPatchField.C"
806  #include "fvPatchFieldNew.C"
807  #include "calculatedFvPatchField.H"
808  #include "zeroGradientFvPatchField.H"
809 #endif
810 
811 // Runtime selection macros
812 #include "fvPatchFieldMacros.H"
814 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
815 
816 #endif
817 
818 // ************************************************************************* //
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
virtual void initEvaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Initialise the evaluation of the patch field.
Definition: fvPatchField.H:792
dictionary dict
virtual bool coupled() const
True if the patch field is coupled.
Definition: fvPatchField.H:253
const objectRegistry & db() const
The associated objectRegistry.
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:77
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 evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field, sets updated() to false.
Definition: fvPatchField.C:329
virtual void operator*=(const fvPatchField< scalar > &)
Definition: fvPatchField.C:434
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 tmp< Field< Type > > valueBoundaryCoeffs(const tmp< Field< scalar >> &) const
Return the matrix source coefficients corresponding to the.
Definition: fvPatchField.H:849
Template invariant parts for fvPatchField.
Definition: fvPatchField.H:77
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:403
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
virtual void readDict(const dictionary &dict)
Read dictionary entries.
static tmp< fvPatchField< Type > > NewCalculatedType(const fvPatch &p)
Return a pointer to a new calculatedFvPatchField created on.
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.
static tmp< fvPatchField< Type > > Clone(const DerivedPatchField &pf, Args &&... args)
Clone a patch field, optionally with internal field reference etc.
Definition: fvPatchField.H:597
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:782
bool operator!=(const fvPatchField< Type > &) const =delete
fvPatch Patch
The patch type for the patch field.
Definition: fvPatchField.H:356
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:862
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 assumed 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
virtual void initEvaluateLocal(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Initialise the evaluation of the patch field after a local.
Definition: fvPatchField.H:813
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:835
const Field< Type > & primitiveField() const noexcept
Return const-reference to the internal field values.
Definition: fvPatchField.H:704
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 ~fvPatchFieldBase()=default
Destructor.
const word & patchType() const noexcept
The optional patch type.
Definition: fvPatchField.H:277
calculatedFvPatchField< Type > Calculated
Type for a calculated patch.
Definition: fvPatchField.H:371
Type value_type
The value_type for the patch field.
Definition: fvPatchField.H:361
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:888
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...
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:696
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.
"buffered" : (MPI_Bsend, MPI_Recv)
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Definition: fvPatchField.H:213
Foam::argList args(argc, argv)
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:696
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
virtual void evaluateLocal(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field after a local operation (e.g. *=)
Definition: fvPatchField.H:823
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
Clone patch field with its own internal field reference.
Definition: fvPatchField.H:567
DimensionedField< Type, volMesh > Internal
The internal field type associated with the patch field.
Definition: fvPatchField.H:366
readOption
Enumeration defining read preferences.