Loading...
Searching...
No Matches
fixedMeanOutletInletFvPatchField.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2018 OpenFOAM Foundation
9 Copyright (C) 2020-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
27\*---------------------------------------------------------------------------*/
28
30#include "volFields.H"
31#include "surfaceFields.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
35template<class Type>
37(
38 const fvPatch& p,
40)
42 outletInletFvPatchField<Type>(p, iF),
43 meanValue_()
44{}
45
46
47template<class Type>
49(
50 const fvPatch& p,
52 const dictionary& dict
53)
54:
55 outletInletFvPatchField<Type>(p, iF),
56 meanValue_(Function1<Type>::New("meanValue", dict, &this->db()))
57{
58 this->phiName_ = dict.getOrDefault<word>("phi", "phi");
59
61
62 this->refValue() = *this;
63 this->refGrad() = Zero;
64 this->valueFraction() = 0.0;
65}
66
67
68template<class Type>
70(
72 const fvPatch& p,
74 const fvPatchFieldMapper& mapper
75)
77 outletInletFvPatchField<Type>(ptf, p, iF, mapper),
78 meanValue_(ptf.meanValue_.clone())
79{}
80
81
82template<class Type>
84(
86)
88 outletInletFvPatchField<Type>(ptf),
89 meanValue_(ptf.meanValue_.clone())
90{}
91
92
93template<class Type>
95(
98)
99:
100 outletInletFvPatchField<Type>(ptf, iF),
101 meanValue_(ptf.meanValue_.clone())
102{}
103
104
105// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106
107template<class Type>
109{
110 if (this->updated())
111 {
112 return;
113 }
114
115 const scalar t = this->db().time().timeOutputValue();
116 Type meanValue = meanValue_->value(t);
117
118 Field<Type> newValues(this->patchInternalField());
119
120 Type meanValuePsi = gWeightedAverage(this->patch().magSf(), newValues);
121
122 if (mag(meanValue) > SMALL && mag(meanValuePsi) > 0.5*mag(meanValue))
123 {
124 newValues *= mag(meanValue)/mag(meanValuePsi);
125 }
126 else
127 {
128 newValues += (meanValue - meanValuePsi);
129 }
130
131 this->refValue() = newValues;
132
134}
135
136
137template<class Type>
139{
141 meanValue_->writeData(os);
143}
144
145
146// ************************************************************************* //
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
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
This boundary condition extrapolates field to the patch using the near-cell values and adjusts the di...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
fixedMeanOutletInletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual tmp< fvPatchField< Type > > clone() const
Return a clone.
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual Field< Type > & refGrad()
virtual Field< Type > & refValue()
virtual scalarField & valueFraction()
This boundary condition provides a generic inflow condition, with specified outflow for the case of r...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
word phiName_
Name of flux field (default: phi).
outletInletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dictionary dict
Foam::surfaceFields.