Loading...
Searching...
No Matches
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 (
46 *this
47 ),
48 weightFieldName_()
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 // Bypass dictionary constructor (all reading handled later)
65 // but cannot use NO_READ since will still trigger an evaluate()
66 mixedFvPatchField<Type>(p, iF),
68 (
69 mappedFixedValueFvPatchField<Type>::mapper(p, iF),
70 *this,
71 dict
72 ),
73 weightFieldName_(dict.getOrDefault<word>("weightField", word::null))
74{
75 fvPatchFieldBase::readDict(dict); // Consistent with a dict constructor
76
78
79 if (this->readMixedEntries(dict))
80 {
81 // Full restart
82 }
83 else
84 {
85 // Start from user entered data. Assume fixedValue.
86 this->refValue() = *this;
87 this->refGrad() = Zero;
88 this->valueFraction() = 1.0;
89 }
90
91// This blocks (crashes) with more than two worlds!
92//
98
101/// this->internalField().name() + "_weights",
102/// this->patch().deltaCoeffs()
103/// );
104}
105
106
107template<class Type>
109(
111 const fvPatch& p,
113 const fvPatchFieldMapper& mapper
114)
115:
116 mixedFvPatchField<Type>(ptf, p, iF, mapper),
118 (
119 mappedFixedValueFvPatchField<Type>::mapper(p, iF),
120 *this,
121 ptf
122 ),
123 weightFieldName_(ptf.weightFieldName_)
124{}
125
126
127template<class Type>
129(
131)
132:
134 mappedPatchFieldBase<Type>(ptf),
135 weightFieldName_(ptf.weightFieldName_)
136{}
137
138
139template<class Type>
141(
144)
145:
146 mixedFvPatchField<Type>(ptf, iF),
148 (
149 mappedFixedValueFvPatchField<Type>::mapper(ptf.patch(), iF),
150 *this,
151 ptf
152 ),
153 weightFieldName_(ptf.weightFieldName_)
154{}
155
156
157// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158
159template<class Type>
168(
169 const fvPatchField<Type>& ptf,
170 const labelList& addr
172{
174}
175
176
177template<class Type>
179{
180 if (this->updated())
181 {
182 return;
183 }
184
185 const tmp<Field<Type>> nbrIntFld(this->mappedInternalField());
186
187 //- Unweighted
188 //const tmp<scalarField> nbrKDelta(this->mappedWeightField());
189
190 //- Weighted
191 tmp<scalarField> myKDelta;
192 tmp<scalarField> nbrKDelta;
193 this->mappedWeightField(weightFieldName_, myKDelta, nbrKDelta);
194
195 // Both sides agree on
196 // - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
197 // - gradient : (temperature-fld)*delta
198 // We've got a degree of freedom in how to implement this in a mixed bc.
199 // (what gradient, what fixedValue and mixing coefficient)
200 // Two reasonable choices:
201 // 1. specify above temperature on one side (preferentially the high side)
202 // and above gradient on the other. So this will switch between pure
203 // fixedvalue and pure fixedgradient
204 // 2. specify gradient and temperature such that the equations are the
205 // same on both sides. This leads to the choice of
206 // - refGradient = zero gradient
207 // - refValue = neighbour value
208 // - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
209
210 this->refValue() = nbrIntFld;
211 this->refGrad() = Zero;
212 this->valueFraction() = nbrKDelta()/(nbrKDelta() + myKDelta());
213
215
216 if (debug)
217 {
218 auto limits = gMinMax(*this);
219 auto avg = gAverage(*this);
220
221 Info<< this->patch().boundaryMesh().mesh().name() << ':'
222 << this->patch().name() << ':'
223 << this->internalField().name() << " <- "
224 << this->mapper_.sampleRegion() << ':'
225 << this->mapper_.samplePatch() << ':'
226 << this->fieldName_ << " :"
227 << " value "
228 << " min:" << limits.min()
229 << " max:" << limits.max()
230 << " avg:" << avg
231 << endl;
232 }
233}
234
236template<class Type>
238{
240 os.writeEntryIfDifferent<word>("weightField", word::null, weightFieldName_);
242}
243
244
245// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static const Field< Type > & null() noexcept
Return a null Field (reference to a nullObject). Behaves like an empty Field.
Definition Field.H:192
@ 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
virtual void readDict(const dictionary &dict)
Read dictionary entries.
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.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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
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 write(Ostream &) const
Write.
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.
mappedMixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Functionality for sampling fields using mappedPatchBase. Every call to mappedField() returns a sample...
const mappedPatchBase & mapper_
Mapping engine.
virtual void write(Ostream &os) const
Write.
static const mappedPatchBase & mapper(const fvPatch &p, const DimensionedField< Type, volMesh > &iF)
Check that patch is of correct type.
word fieldName_
Name of field to sample.
virtual tmp< scalarField > mappedWeightField() const
Map optional weightField onto *this patch. Default is deltaCoeffs.
virtual tmp< Field< Type > > mappedInternalField() const
Map internal of sampleField onto *this patch.
mappedPatchFieldBase(const mappedPatchBase &mapper, const fvPatchField< Type > &patchField, const word &fieldName, const bool setAverage, const Type average, const word &interpolationScheme)
Construct from components.
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.
mixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &) const
Write.
virtual Field< Type > & refGrad()
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
bool readMixedEntries(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the 'refValue', 'refGradient' and 'valueFraction' entries into their respective places.
virtual Field< Type > & refValue()
virtual scalarField & valueFraction()
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
volScalarField & p
auto limits
Definition setRDeltaT.H:186
OBJstream os(runTime.globalPath()/outputName)
Namespace for handling debugging switches.
Definition debug.C:45
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dictionary dict