mappedMixedFvPatchField.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) 2020-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#include "volFields.H"
30#include "interpolationCell.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
35template<class Type>
37(
38 const fvPatch& p,
40)
41:
42 mixedFvPatchField<Type>(p, iF),
44 (
45 mappedFixedValueFvPatchField<Type>::mapper(p, iF),
46 *this
47 ),
48 weightFieldName_(word::null)
49{
50 this->refValue() = Zero;
51 this->refGrad() = Zero;
52 this->valueFraction() = 0.0;
53}
54
55
56template<class Type>
58(
59 const fvPatch& p,
61 const dictionary& dict
62)
63:
64 mixedFvPatchField<Type>(p, iF, dict),
66 (
67 mappedFixedValueFvPatchField<Type>::mapper(p, iF),
68 *this,
69 dict
70 ),
71 weightFieldName_(dict.getOrDefault<word>("weightField", word::null))
72{
74 (
75 Field<Type>("value", dict, p.size())
76 );
77
78 if (dict.found("refValue"))
79 {
80 // Full restart
81 this->refValue() = Field<Type>("refValue", dict, p.size());
82 this->refGrad() = Field<Type>("refGradient", dict, p.size());
83 this->valueFraction() = scalarField("valueFraction", dict, p.size());
84 }
85 else
86 {
87 // Start from user entered data. Assume fixedValue.
88 this->refValue() = *this;
89 this->refGrad() = Zero;
90 this->valueFraction() = 1.0;
91 }
92
93// This blocks (crashes) with more than two worlds!
94//
106}
107
108
109template<class Type>
111(
113 const fvPatch& p,
115 const fvPatchFieldMapper& mapper
116)
117:
118 mixedFvPatchField<Type>(ptf, p, iF, mapper),
120 (
121 mappedFixedValueFvPatchField<Type>::mapper(p, iF),
122 *this,
123 ptf
124 ),
125 weightFieldName_(ptf.weightFieldName_)
126{}
127
128
129template<class Type>
131(
133)
134:
135 mixedFvPatchField<Type>(ptf),
136 mappedPatchFieldBase<Type>(ptf),
137 weightFieldName_(ptf.weightFieldName_)
138{}
139
140
141template<class Type>
143(
146)
147:
148 mixedFvPatchField<Type>(ptf, iF),
150 (
151 mappedFixedValueFvPatchField<Type>::mapper(ptf.patch(), iF),
152 *this,
153 ptf
154 ),
155 weightFieldName_(ptf.weightFieldName_)
156{}
157
158
159// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
160
161template<class Type>
163{
165}
166
167
168template<class Type>
170(
171 const fvPatchField<Type>& ptf,
172 const labelList& addr
173)
174{
176}
177
178
179template<class Type>
181{
182 if (this->updated())
183 {
184 return;
185 }
186
187 const tmp<Field<Type>> nbrIntFld(this->mappedInternalField());
188
189 //- Unweighted
190 //const tmp<scalarField> nbrKDelta(this->mappedWeightField());
191
192 //- Weighted
193 tmp<scalarField> myKDelta;
194 tmp<scalarField> nbrKDelta;
195 this->mappedWeightField(weightFieldName_, myKDelta, nbrKDelta);
196
197 // Both sides agree on
198 // - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
199 // - gradient : (temperature-fld)*delta
200 // We've got a degree of freedom in how to implement this in a mixed bc.
201 // (what gradient, what fixedValue and mixing coefficient)
202 // Two reasonable choices:
203 // 1. specify above temperature on one side (preferentially the high side)
204 // and above gradient on the other. So this will switch between pure
205 // fixedvalue and pure fixedgradient
206 // 2. specify gradient and temperature such that the equations are the
207 // same on both sides. This leads to the choice of
208 // - refGradient = zero gradient
209 // - refValue = neighbour value
210 // - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
211
212 this->refValue() = nbrIntFld;
213 this->refGrad() = Zero;
214 this->valueFraction() = nbrKDelta()/(nbrKDelta() + myKDelta());
215
217
218 if (debug)
219 {
220 Info<< this->patch().boundaryMesh().mesh().name() << ':'
221 << this->patch().name() << ':'
222 << this->internalField().name() << " <- "
223 << this->mapper_.sampleRegion() << ':'
224 << this->mapper_.samplePatch() << ':'
225 << this->fieldName_ << " :"
226 << " value "
227 << " min:" << gMin(*this)
228 << " max:" << gMax(*this)
229 << " avg:" << gAverage(*this)
230 << endl;
231 }
232}
233
234
235template<class Type>
237{
239 os.writeEntryIfDifferent<word>("weightField", word::null, weightFieldName_);
241}
242
243
244// ************************************************************************* //
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:82
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:251
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
virtual bool write()
Write the output 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...
Definition: fvPatchField.H:82
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
This boundary condition maps the value at a set of cells or patch faces back to *this.
This boundary condition maps the value at a set of cells or patch faces back to *this.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Functionality for sampling fields using mappedPatchBase. Every call to mappedField() returns a sample...
This boundary condition provides a base class for 'mixed' type boundary conditions,...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual Field< Type > & refGrad()
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual Field< Type > & refValue()
virtual scalarField & valueFraction()
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
static const word null
An empty word.
Definition: word.H:80
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
dictionary dict