mappedMixedFieldFvPatchField.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) 2019-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"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
34template<class Type>
36(
37 const fvPatch& p,
39)
40:
41 mixedFvPatchField<Type>(p, iF),
42 mappedPatchBase(p.patch()),
43 mappedPatchFieldBase<Type>(*this, *this),
44 weightFieldName_(word::null)
45{
46 this->refValue() = Zero;
47 this->refGrad() = Zero;
48 this->valueFraction() = 0.0;
49}
50
51
52template<class Type>
54(
55 const fvPatch& p,
57 const dictionary& dict
58)
59:
60 mixedFvPatchField<Type>(p, iF, dict),
61 mappedPatchBase(p.patch(), dict),
62 mappedPatchFieldBase<Type>(*this, *this, dict),
63 weightFieldName_(dict.getOrDefault<word>("weightField", word::null))
64{}
65
66
67template<class Type>
69(
71 const fvPatch& p,
73 const fvPatchFieldMapper& mapper
74)
75:
76 mixedFvPatchField<Type>(ptf, p, iF, mapper),
77 mappedPatchBase(p.patch(), ptf),
78 mappedPatchFieldBase<Type>(*this, *this, ptf),
79 weightFieldName_(ptf.weightFieldName_)
80{}
81
82
83template<class Type>
85(
87)
88:
89 mixedFvPatchField<Type>(ptf),
90 mappedPatchBase(ptf.patch().patch(), ptf),
91 mappedPatchFieldBase<Type>(ptf),
92 weightFieldName_(ptf.weightFieldName_)
93{}
94
95
96template<class Type>
98(
101)
102:
103 mixedFvPatchField<Type>(ptf, iF),
104 mappedPatchBase(ptf.patch().patch(), ptf),
105 mappedPatchFieldBase<Type>(*this, *this, ptf),
106 weightFieldName_(ptf.weightFieldName_)
107{}
108
109
110// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111
112template<class Type>
114(
115 const fvPatchFieldMapper& m
116)
117{
120}
121
122
123template<class Type>
125(
126 const fvPatchField<Type>& ptf,
127 const labelList& addr
128)
129{
132}
133
134
135template<class Type>
137{
138 if (this->updated())
139 {
140 return;
141 }
142
143 const tmp<Field<Type>> nbrIntFld(this->mappedInternalField());
144
145 //- Unweighted
146 //const tmp<scalarField> nbrKDelta(this->mappedWeightField());
147
148 //- Weighted
149 tmp<scalarField> myKDelta;
150 tmp<scalarField> nbrKDelta;
151 this->mappedWeightField(weightFieldName_, myKDelta, nbrKDelta);
152
153 // Both sides agree on
154 // - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
155 // - gradient : (temperature-fld)*delta
156 // We've got a degree of freedom in how to implement this in a mixed bc.
157 // (what gradient, what fixedValue and mixing coefficient)
158 // Two reasonable choices:
159 // 1. specify above temperature on one side (preferentially the high side)
160 // and above gradient on the other. So this will switch between pure
161 // fixedvalue and pure fixedgradient
162 // 2. specify gradient and temperature such that the equations are the
163 // same on both sides. This leads to the choice of
164 // - refGradient = zero gradient
165 // - refValue = neighbour value
166 // - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
167
168 this->refValue() = nbrIntFld;
169 this->refGrad() = Zero;
170 this->valueFraction() = nbrKDelta()/(nbrKDelta() + myKDelta());
171
173
174 if (debug)
175 {
176 Info<< this->patch().boundaryMesh().mesh().name() << ':'
177 << this->patch().name() << ':'
178 << this->internalField().name() << " <- "
179 << this->mapper_.sampleRegion() << ':'
180 << this->mapper_.samplePatch() << ':'
181 << this->fieldName_ << " :"
182 << " value "
183 << " min:" << gMin(*this)
184 << " max:" << gMax(*this)
185 << " avg:" << gAverage(*this)
186 << endl;
187 }
188}
189
190
191template<class Type>
193{
196 os.writeEntryIfDifferent<word>("weightField", word::null, weightFieldName_);
198}
199
200
201// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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
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
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 provides a self-contained version of e.g. mapped boundary conditions.
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.
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
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)
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