Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::fvPatchField
29
30Description
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
41SourceFiles
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
58namespace Foam
59{
60
61// Forward Declarations
62class dictionary;
63class objectRegistry;
65class volMesh;
66
67template<class Type> class fvPatchField;
68template<class Type> class calculatedFvPatchField;
69template<class Type> class fvMatrix;
70
71template<class Type>
73
74
75/*---------------------------------------------------------------------------*\
76 Class fvPatchFieldBase Declaration
77\*---------------------------------------------------------------------------*/
78
79//- Template invariant parts for fvPatchField
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
105protected:
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
114public:
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);
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 //- The type name for \c zeroValue patch fields
173 static const word& zeroValueType() noexcept
174 {
176 }
178
179 // Member Functions
180
181 // Attributes
182
183 //- True if the value of the patch field is altered by assignment
184 virtual bool assignable() const
186 return true;
187 }
188
189 //- True if the patch field fixes a value.
190 // Needed to check if a level has to be specified while solving
191 // Poissons equations.
192 virtual bool fixesValue() const
193 {
194 return false;
195 }
196
197 //- True if the patch field is coupled
198 virtual bool coupled() const
199 {
200 return false;
201 }
203
204 // Access
205
206 //- The associated objectRegistry
207 const objectRegistry& db() const;
208
209 //- Return the patch
210 const fvPatch& patch() const noexcept
211 {
212 return patch_;
213 }
214
215 //- The optional patch type
216 const word& patchType() const noexcept
217 {
218 return patchType_;
219 }
220
221 //- The optional patch type
222 word& patchType() noexcept
224 return patchType_;
225 }
226
227 //- True if the type does not correspond to the constraint type
228 virtual bool constraintOverride() const
229 {
230 return !patchType_.empty() && patchType_ != type();
231 }
232
233
234 // Solution
235
236 //- True if the boundary condition has already been updated
237 bool updated() const noexcept
238 {
239 return updated_;
240 }
241
242 //- Set updated state
243 void setUpdated(bool state) noexcept
244 {
245 updated_ = state;
246 }
247
248 //- True if the matrix has already been manipulated
249 bool manipulatedMatrix() const noexcept
250 {
251 return manipulatedMatrix_;
252 }
253
254 //- Set matrix manipulated state
255 void setManipulated(bool state) noexcept
256 {
257 manipulatedMatrix_ = state;
259
260 //- Use implicit formulation for coupled patches only
261 bool useImplicit() const noexcept
262 {
263 return useImplicit_;
264 }
265
266 //- Set useImplicit on/off
267 // \return old value
268 bool useImplicit(bool on) noexcept
269 {
270 bool old(useImplicit_);
271 useImplicit_ = on;
272 return old;
273 }
275
276 // Check
277
278 //- Check that patches are identical
279 void checkPatch(const fvPatchFieldBase& rhs) const;
280};
281
283/*---------------------------------------------------------------------------*\
284 Class fvPatchField Declaration
285\*---------------------------------------------------------------------------*/
286
287template<class Type>
288class fvPatchField
289:
290 public fvPatchFieldBase,
291 public Field<Type>
292{
293public:
294
295 // Public Data Types
296
297 //- The patch type for the patch field
298 typedef fvPatch Patch;
299
300 //- The value_type for the patch field
301 typedef Type value_type;
302
303 //- The component type for patch field
304 typedef typename pTraits<Type>::cmptType cmptType;
305
306 //- The internal field type associated with the patch field
308
309 //- Type for a \em calculated patch
311
312
313private:
314
315 // Private Data
316
317 //- Reference to internal field
318 const DimensionedField<Type, volMesh>& internalField_;
319
320
321protected:
322
323 // Protected Member Functions
324
325 //- Read the "value" entry into \c *this.
326 // The reading can be optional (default), mandatory etc.
327 // \returns True on success
328 bool readValueEntry
329 (
330 const dictionary& dict,
332 );
333
334 //- Write \c *this field as a "value" entry
336 {
337 Field<Type>::writeEntry("value", os);
338 }
339
340 //- Assign the patch field from the internal field
341 void extrapolateInternal();
342
343
344public:
345
346 // Declare run-time constructor selection tables
347
349 (
350 tmp,
352 patch,
353 (
354 const fvPatch& p,
356 ),
357 (p, iF)
358 );
359
361 (
362 tmp,
364 patchMapper,
365 (
366 const fvPatchField<Type>& ptf,
367 const fvPatch& p,
370 ),
371 (dynamic_cast<const fvPatchFieldType&>(ptf), p, iF, m)
372 );
373
375 (
376 tmp,
380 const fvPatch& p,
382 const dictionary& dict
383 ),
384 (p, iF, dict)
385 );
386
387
388 // Constructors
390 //- Construct from patch and internal field
392 (
393 const fvPatch&,
395 );
396
397 //- Construct from patch, internal field and patch type
399 (
400 const fvPatch&,
402 const word& patchType
403 );
404
405 //- Construct from patch, internal field and value
407 (
408 const fvPatch&,
410 const Type& value
411 );
413 //- Construct from patch, internal field and patch field
415 (
416 const fvPatch&,
418 const Field<Type>& pfld
419 );
420
421 //- Construct from patch, internal field and patch field
423 (
424 const fvPatch&,
426 Field<Type>&& pfld
427 );
428
429 //- Construct from patch, internal field and dictionary
431 (
432 const fvPatch&,
434 const dictionary& dict,
437 );
438
439 //- Construct, forwarding to readOption variant
441 (
442 const fvPatch& p,
444 const dictionary& dict,
445 const bool needValue
446 )
447 :
449 (
450 p, iF, dict,
451 (needValue? IOobjectOption::MUST_READ : IOobjectOption::NO_READ)
452 )
453 {}
454
455 //- Construct by mapping the given fvPatchField onto a new patch
457 (
458 const fvPatchField<Type>&,
459 const fvPatch&,
461 const fvPatchFieldMapper&
462 );
463
464 //- Copy construct onto a new patch with internal field reference
465 //- and specified value
467 (
468 const fvPatchField<Type>& pfld,
469 const fvPatch& p,
471 const Type& value
472 );
473
474 //- Copy construct with internal field reference
476 (
477 const fvPatchField<Type>& pfld,
479 );
480
481 //- Copy construct
483 :
484 fvPatchField<Type>(pfld, pfld.internalField())
485 {}
486
487
488 //- Clone patch field with its own internal field reference
489 virtual tmp<fvPatchField<Type>> clone() const
492 (
493 new fvPatchField<Type>(*this, this->internalField_)
494 );
495 }
496
497 //- Clone patch field with given internal field reference
499 (
501 ) const
502 {
504 (
505 new fvPatchField<Type>(*this, iF)
506 );
507 }
508
509
510 // Factory Methods
511
512 //- Clone a patch field, optionally with internal field reference etc.
513 template<class DerivedPatchField, class... Args>
515 (
516 const DerivedPatchField& pf,
517 Args&&... args
518 )
519 {
521 (
522 new DerivedPatchField(pf, std::forward<Args>(args)...)
523 );
524 }
525
526 //- Return a pointer to a new patchField created on freestore given
527 // patch and internal field
528 // (does not set the patch field values)
531 const word& patchFieldType,
532 const fvPatch&,
534 );
535
536 //- Return a pointer to a new patchField created on freestore given
537 // patch and internal field
538 // (does not set the patch field values).
539 // Allows override of constraint type
541 (
542 const word& patchFieldType,
543 const word& actualPatchType,
544 const fvPatch&,
546 );
547
548 //- Return a pointer to a new patchField created on freestore from
549 // a given fvPatchField mapped onto a new patch
551 (
552 const fvPatchField<Type>&,
553 const fvPatch&,
555 const fvPatchFieldMapper&
556 );
557
558 //- Return a pointer to a new patchField created on freestore
559 // from dictionary
561 (
562 const fvPatch&,
564 const dictionary&
565 );
566
567 //- Return a pointer to a new calculatedFvPatchField created on
568 // freestore without setting patchField values
570 (
571 const fvPatch& p
572 );
573
574 //- Return a pointer to a new calculatedFvPatchField created on
575 // freestore without setting patchField values
576 template<class AnyType>
578 (
579 const fvPatchField<AnyType>& pf
580 );
581
582
583 //- Destructor
584 virtual ~fvPatchField() = default;
585
586
587 // Member Functions
588
589 // Access
590
591 //- Return const-reference to the dimensioned internal field
593 {
594 return internalField_;
595 }
596
597 //- Return const-reference to the internal field values
598 const Field<Type>& primitiveField() const noexcept
599 {
600 return internalField_;
602
603
604 // Mapping Functions
605
606 //- Map (and resize as needed) from self given a mapping object
607 virtual void autoMap
608 (
609 const fvPatchFieldMapper&
610 );
611
612 //- Reverse map the given fvPatchField onto this fvPatchField
613 virtual void rmap
614 (
615 const fvPatchField<Type>&,
616 const labelList&
617 );
618
619
620 // Evaluation Functions
621
622 //- Return patch-normal gradient
623 virtual tmp<Field<Type>> snGrad() const;
624
625 //- Retrieve patch-normal gradient [contiguous storage]
626 virtual void snGrad(UList<Type>& result) const;
627
628 //- Return patch-normal gradient for coupled-patches
629 //- using the deltaCoeffs provided
631 (
632 const scalarField& deltaCoeffs
633 ) const
634 {
636 return *this;
637 }
638
639 //- Retrieve patch-normal gradient for coupled-patches
640 //- using the deltaCoeffs provided [contiguous storage]
641 virtual void snGrad
642 (
643 const scalarField& deltaCoeffs,
645 ) const
646 {
649
650 //- Return internal field next to patch
651 virtual tmp<Field<Type>> patchInternalField() const;
652
653 //- Retrieve internal field next to patch
654 // \param [out] pfld The extracted patch field.
655 virtual void patchInternalField(UList<Type>& pfld) const;
656
657 //- Return patchField on the opposite patch of a coupled patch
658 virtual tmp<Field<Type>> patchNeighbourField() const
659 {
661 return *this;
663
664 //- Retrieve patchField on the opposite patch of a coupled patch
665 virtual void patchNeighbourField(UList<Type>&) const
666 {
668 }
669
670 //- Update the coefficients associated with the patch field
671 // Sets Updated to true
672 virtual void updateCoeffs();
673
674 //- Update the coefficients associated with the patch field
675 // with a weight field (0..1). This weight field is usually
676 // provided as the amount of geometric overlap for 'duplicate'
677 // patches. Sets Updated to true
678 virtual void updateWeightedCoeffs(const scalarField& weights);
679
680 //- Initialise the evaluation of the patch field
681 virtual void initEvaluate
682 (
683 const Pstream::commsTypes commsType =
685 )
686 {}
687
688 //- Evaluate the patch field, sets updated() to false
689 virtual void evaluate
690 (
691 const Pstream::commsTypes commsType =
693 );
694
695 //- Initialise the evaluation of the patch field after a local
696 // operation
697 virtual void initEvaluateLocal
698 (
699 const Pstream::commsTypes commsType =
701 )
702 {}
703
704 //- Evaluate the patch field after a local operation (e.g. *=)
705 virtual void evaluateLocal
706 (
707 const Pstream::commsTypes commsType =
709 )
710 {}
712 //- Return the matrix diagonal coefficients corresponding to the
713 //- evaluation of the value of this patchField with given weights
715 (
716 const tmp<scalarField>&
717 ) const
718 {
720 return *this;
721 }
722
723 //- Return the matrix source coefficients corresponding to the
724 //- evaluation of the value of this patchField with given weights
726 (
727 const tmp<scalarField>&
728 ) const
729 {
731 return *this;
732 }
733
734 //- Return the matrix diagonal coefficients corresponding to the
735 //- evaluation of the gradient of this patchField
737 {
739 return *this;
740 }
741
742 //- Return the matrix diagonal coefficients corresponding to the
743 //- evaluation of the gradient of this coupled patchField
744 //- using the deltaCoeffs provided
746 (
747 const scalarField& deltaCoeffs
748 ) const
751 return *this;
752 }
753
754 //- Return the matrix source coefficients corresponding to the
755 //- evaluation of the gradient of this patchField
759 return *this;
760 }
761
762 //- Return the matrix source coefficients corresponding to the
763 //- evaluation of the gradient of this coupled patchField
764 //- using the deltaCoeffs provided
766 (
767 const scalarField& deltaCoeffs
768 ) const
771 return *this;
772 }
773
774 //- Manipulate matrix
775 virtual void manipulateMatrix(fvMatrix<Type>& matrix);
776
777 //- Manipulate matrix with given weights
778 virtual void manipulateMatrix
779 (
781 const scalarField& weights
782 );
783
784 //- Manipulate fvMatrix
785 virtual void manipulateMatrix
786 (
787 fvMatrix<Type>& matrix,
788 const label iMatrix,
789 const direction cmp
790 );
791
792
793 // Contiguous storage
794
795 //- Retrieve the matrix diagonal coefficients corresponding to the
796 //- evaluation of the value of this patchField with given weights
797 virtual void valueInternalCoeffs
798 (
799 const tmp<Field<scalar>>&,
801 ) const
802 {
804 }
806 //- Retrieve the matrix source coefficients corresponding to the
807 //- evaluation of the value of this patchField with given weights
808 virtual void valueBoundaryCoeffs
809 (
810 const tmp<Field<scalar>>&,
812 ) const
813 {
815 }
816
817 //- Retrieve the matrix diagonal coefficients corresponding to the
818 //- evaluation of the gradient of this patchField
819 virtual void gradientInternalCoeffs(UList<Type>&) const
820 {
822 }
823
824 //- Retrieve the matrix diagonal coefficients corresponding to the
825 //- evaluation of the gradient of this coupled patchField
826 //- using the deltaCoeffs provided
827 virtual void gradientInternalCoeffs
828 (
829 const scalarField& deltaCoeffs,
831 ) const
832 {
834 }
835
836 //- Retrieve the matrix source coefficients corresponding to the
837 //- evaluation of the gradient of this patchField
838 virtual void gradientBoundaryCoeffs(UList<Type>&) const
839 {
841 }
842
843 //- Retrieve the matrix source coefficients corresponding to the
844 //- evaluation of the gradient of this coupled patchField
845 //- using the deltaCoeffs provided
846 virtual void gradientBoundaryCoeffs
847 (
848 const scalarField& deltaCoeffs,
850 ) const
851 {
853 }
854
855
856 // Other
857
858 //- Write
859 virtual void write(Ostream&) const;
861 //- Check against given patch field
862 void check(const fvPatchField<Type>&) const;
863
864
865 // Member Operators
866
867 virtual void operator=(const UList<Type>&);
868
869 virtual void operator=(const fvPatchField<Type>&);
870 virtual void operator+=(const fvPatchField<Type>&);
871 virtual void operator-=(const fvPatchField<Type>&);
872 virtual void operator*=(const fvPatchField<scalar>&);
873 virtual void operator/=(const fvPatchField<scalar>&);
874
875 virtual void operator+=(const Field<Type>&);
876 virtual void operator-=(const Field<Type>&);
877
878 virtual void operator*=(const Field<scalar>&);
879 virtual void operator/=(const Field<scalar>&);
880
881 virtual void operator=(const Type&);
882 virtual void operator+=(const Type&);
883 virtual void operator-=(const Type&);
884 virtual void operator*=(const scalar);
885 virtual void operator/=(const scalar);
886
887
888 // Force an assignment irrespective of form of patch
889
890 virtual void operator==(const fvPatchField<Type>&);
891 virtual void operator==(const Field<Type>&);
892 virtual void operator==(const Type&);
893
894 // Prevent automatic comparison rewriting (c++20)
895 bool operator!=(const fvPatchField<Type>&) const = delete;
896 bool operator!=(const Field<Type>&) const = delete;
897 bool operator!=(const Type&) const = delete;
898
899
900 // Ostream Operator
901
902 friend Ostream& operator<< <Type>(Ostream&, const fvPatchField<Type>&);
903};
904
906// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
907
908} // End namespace Foam
909
910// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
911
912#ifdef NoRepository
913 #include "fvPatchField.C"
914 #include "fvPatchFieldNew.C"
915 #include "calculatedFvPatchField.H"
917#endif
919// Runtime selection macros
920#include "fvPatchFieldMacros.H"
921
922// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
923
924#endif
925
926// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition Field.C:754
constexpr Field() noexcept
Default construct.
Definition FieldI.H:24
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
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
commsTypes
Communications types.
Definition UPstream.H:81
@ buffered
"buffered" : (MPI_Bsend, MPI_Recv)
Definition UPstream.H:82
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:133
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
Template invariant parts for fvPatchField.
fvPatchFieldBase(const fvPatch &p)
Construct from patch.
virtual ~fvPatchFieldBase()=default
Destructor.
virtual bool fixesValue() const
True if the patch field fixes a value.
void checkPatch(const fvPatchFieldBase &rhs) const
Check that patches are identical.
virtual bool coupled() const
True if the patch field is coupled.
static const word & zeroValueType() noexcept
The type name for zeroValue patch fields.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
const objectRegistry & db() const
The associated objectRegistry.
virtual void readDict(const dictionary &dict)
Read dictionary entries.
const fvPatch & patch() const noexcept
Return the patch.
void setManipulated(bool state) noexcept
Set matrix manipulated state.
bool useImplicit() const noexcept
Use implicit formulation for coupled patches only.
void setUpdated(bool state) noexcept
Set updated state.
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
static int disallowGenericPatchField
Debug switch to disallow the use of generic fvPatchField.
TypeName("fvPatchField")
Runtime type information.
bool useImplicit(bool on) noexcept
Set useImplicit on/off.
const word & patchType() const noexcept
The optional patch type.
bool manipulatedMatrix() const noexcept
True if the matrix has already been manipulated.
virtual bool constraintOverride() const
True if the type does not correspond to the constraint type.
virtual bool assignable() const
True if the value of the patch field is altered by assignment.
bool updated() const noexcept
True if the boundary condition has already been updated.
word & patchType() noexcept
The optional patch type.
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
static const word & emptyType() noexcept
The type name for empty patch fields.
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...
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the evaluation of the value of this patchField...
fvPatchField(const fvPatch &p, const DimensionedField< Type, volMesh > &iF, const dictionary &dict, const bool needValue)
Construct, forwarding to readOption variant.
virtual void valueBoundaryCoeffs(const tmp< Field< scalar > > &, UList< Type > &) const
Retrieve the matrix source coefficients corresponding to the evaluation of the value of this patchFie...
virtual void operator/=(const Field< scalar > &)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
static tmp< fvPatchField< Type > > NewCalculatedType(const fvPatchField< AnyType > &pf)
Return a pointer to a new calculatedFvPatchField created on.
DimensionedField< Type, volMesh > Internal
The internal field type associated with the patch field.
static tmp< fvPatchField< Type > > New(const fvPatchField< Type > &, const fvPatch &, const DimensionedField< Type, volMesh > &, const fvPatchFieldMapper &)
Return a pointer to a new patchField created on freestore from.
virtual tmp< Field< Type > > snGrad(const scalarField &deltaCoeffs) const
Return patch-normal gradient for coupled-patches using the deltaCoeffs provided.
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
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 operator=(const Type &)
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the evaluation of the value of this patchFie...
virtual void write(Ostream &) const
Write.
virtual void operator-=(const Type &)
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &dict, IOobjectOption::readOption requireValue=IOobjectOption::MUST_READ)
Construct from patch, internal field and dictionary.
declareRunTimeSelectionTable(tmp, fvPatchField, dictionary,(const fvPatch &p, const DimensionedField< Type, volMesh > &iF, const dictionary &dict),(p, iF, dict))
Type value_type
The value_type for the patch field.
virtual void gradientBoundaryCoeffs(const scalarField &deltaCoeffs, UList< Type > &) const
Retrieve the matrix source coefficients corresponding to the evaluation of the gradient of this coupl...
virtual void initEvaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Initialise the evaluation of the patch field.
virtual void operator==(const fvPatchField< Type > &)
declareRunTimeSelectionTable(tmp, fvPatchField, patch,(const fvPatch &p, const DimensionedField< Type, volMesh > &iF),(p, iF))
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
fvPatch Patch
The patch type for the patch field.
virtual void patchNeighbourField(UList< Type > &) const
Retrieve patchField on the opposite patch of a coupled patch.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
static tmp< fvPatchField< Type > > New(const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &)
Return a pointer to a new patchField created on freestore.
static tmp< fvPatchField< Type > > NewCalculatedType(const fvPatch &p)
Return a pointer to a new calculatedFvPatchField created on.
virtual void operator*=(const Field< scalar > &)
virtual void operator==(const Field< Type > &)
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const word &patchType)
Construct from patch, internal field and patch type.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
virtual void operator+=(const Field< Type > &)
virtual void initEvaluateLocal(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Initialise the evaluation of the patch field after a local.
virtual void operator-=(const fvPatchField< Type > &)
const Field< Type > & primitiveField() const noexcept
Return const-reference to the internal field values.
virtual void operator+=(const Type &)
const DimensionedField< scalar, volMesh > & internalField() const noexcept
virtual void evaluateLocal(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field after a local operation (e.g. *=).
virtual void gradientBoundaryCoeffs(UList< Type > &) const
Retrieve the matrix source coefficients corresponding to the evaluation of the gradient of this patch...
virtual tmp< Field< Type > > gradientBoundaryCoeffs(const scalarField &deltaCoeffs) const
Return the matrix source coefficients corresponding to the evaluation of the gradient of this coupled...
virtual void operator-=(const Field< Type > &)
virtual void operator=(const UList< Type > &)
pTraits< Type >::cmptType cmptType
The component type for patch field.
fvPatchField(const fvPatchField< Type > &, const fvPatch &, const DimensionedField< Type, volMesh > &, const fvPatchFieldMapper &)
Construct by mapping the given fvPatchField onto a new patch.
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the evaluation of the gradient of this patch...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
fvPatchField(const fvPatchField< Type > &pfld)
Copy construct.
virtual void gradientInternalCoeffs(const scalarField &deltaCoeffs, UList< Type > &) const
Retrieve the matrix diagonal coefficients corresponding to the evaluation of the gradient of this cou...
virtual ~fvPatchField()=default
Destructor.
calculatedFvPatchField< Type > Calculated
Type for a calculated patch.
virtual void patchInternalField(UList< Type > &pfld) const
Retrieve internal field next to patch.
virtual void snGrad(UList< Type > &result) const
Retrieve patch-normal gradient [contiguous storage].
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
virtual void operator/=(const scalar)
virtual void operator*=(const fvPatchField< scalar > &)
virtual void operator=(const fvPatchField< Type > &)
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const Type &value)
Construct from patch, internal field and value.
virtual void snGrad(const scalarField &deltaCoeffs, UList< Type > &) const
Retrieve patch-normal gradient for coupled-patches using the deltaCoeffs provided [contiguous storage...
bool operator!=(const fvPatchField< Type > &) const =delete
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, Field< Type > &&pfld)
Construct from patch, internal field and patch field.
bool operator!=(const Field< Type > &) const =delete
virtual void operator+=(const fvPatchField< Type > &)
declareRunTimeSelectionTable(tmp, fvPatchField, patchMapper,(const fvPatchField< Type > &ptf, const fvPatch &p, const DimensionedField< Type, volMesh > &iF, const fvPatchFieldMapper &m),(dynamic_cast< const fvPatchFieldType & >(ptf), p, iF, m))
fvPatchField(const fvPatchField< Type > &pfld, const DimensionedField< Type, volMesh > &iF)
Copy construct with internal field reference.
virtual void valueInternalCoeffs(const tmp< Field< scalar > > &, UList< Type > &) const
Retrieve the matrix diagonal coefficients corresponding to the evaluation of the value of this patchF...
bool operator!=(const Type &) const =delete
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
static tmp< fvPatchField< Type > > Clone(const DerivedPatchField &pf, Args &&... args)
Clone a patch field, optionally with internal field reference etc.
virtual void gradientInternalCoeffs(UList< Type > &) const
Retrieve the matrix diagonal coefficients corresponding to the evaluation of the gradient of this pat...
virtual tmp< fvPatchField< Type > > clone() const
Clone patch field with its own internal field reference.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
void extrapolateInternal()
Assign the patch field from the internal field.
virtual void operator/=(const fvPatchField< scalar > &)
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const word &actualPatchType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
void check(const fvPatchField< Type > &) const
Check against given patch field.
virtual void operator*=(const scalar)
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field, sets updated() to false.
virtual void manipulateMatrix(fvMatrix< Type > &matrix, const label iMatrix, const direction cmp)
Manipulate fvMatrix.
virtual void manipulateMatrix(fvMatrix< Type > &matrix, const scalarField &weights)
Manipulate matrix with given weights.
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual void operator==(const Type &)
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const Field< Type > &pfld)
Construct from patch, internal field and patch field.
virtual tmp< fvPatchField< Type > > clone(const DimensionedField< Type, volMesh > &iF) const
Clone patch field with given internal field reference.
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the evaluation of the gradient of this patchFi...
fvPatchField(const fvPatchField< Type > &pfld, const fvPatch &p, const DimensionedField< Type, volMesh > &iF, const Type &value)
Copy construct onto a new patch with internal field reference and specified value.
virtual tmp< Field< Type > > gradientInternalCoeffs(const scalarField &deltaCoeffs) const
Return the matrix diagonal coefficients corresponding to the evaluation of the gradient of this coupl...
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:75
Mesh data needed to do the Finite Volume discretisation.
Definition volMesh.H:47
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.
const word zeroValueType
A zeroValue patch field type.
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
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).
uint8_t direction
Definition direction.H:49
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
const direction noexcept
Definition scalarImpl.H:265
runTime write()
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
dictionary dict
Foam::argList args(argc, argv)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68