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 -------------------------------------------------------------------------------
10 License
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"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<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_()
49 {
50  this->refValue() = Zero;
51  this->refGrad() = Zero;
52  this->valueFraction() = 0.0;
53 }
54 
55 
56 template<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 //
104 }
105 
106 
107 template<class Type>
109 (
110  const mappedMixedFvPatchField<Type>& ptf,
111  const fvPatch& p,
112  const DimensionedField<Type, volMesh>& iF,
113  const fvPatchFieldMapper& mapper
114 )
115 :
116  mixedFvPatchField<Type>(ptf, p, iF, mapper),
117  mappedPatchFieldBase<Type>
118  (
119  mappedFixedValueFvPatchField<Type>::mapper(p, iF),
120  *this,
121  ptf
122  ),
123  weightFieldName_(ptf.weightFieldName_)
124 {}
125 
126 
127 template<class Type>
129 (
131 )
132 :
133  mixedFvPatchField<Type>(ptf),
135  weightFieldName_(ptf.weightFieldName_)
136 {}
137 
138 
139 template<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 
159 template<class Type>
161 {
163 }
164 
165 
166 template<class Type>
168 (
169  const fvPatchField<Type>& ptf,
170  const labelList& addr
171 )
172 {
174 }
175 
176 
177 template<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  Info<< this->patch().boundaryMesh().mesh().name() << ':'
219  << this->patch().name() << ':'
220  << this->internalField().name() << " <- "
221  << this->mapper_.sampleRegion() << ':'
222  << this->mapper_.samplePatch() << ':'
223  << this->fieldName_ << " :"
224  << " value "
225  << " min:" << gMin(*this)
226  << " max:" << gMax(*this)
227  << " avg:" << gAverage(*this)
228  << endl;
229  }
230 }
231 
233 template<class Type>
235 {
237  os.writeEntryIfDifferent<word>("weightField", word::null, weightFieldName_);
239 }
240 
241 
242 // ************************************************************************* //
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
Definition: fvPatchField.C:32
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Type gMin(const FieldField< Field, Type > &f)
virtual void write(Ostream &os) const
Write.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
This boundary condition provides a base class for &#39;mixed&#39; type boundary conditions, i.e. conditions that mix fixed value and patch-normal gradient conditions.
This boundary condition maps the value at a set of cells or patch faces back to *this.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void readDict(const dictionary &dict)
Read dictionary entries.
virtual Field< Type > & refValue()
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
bool readMixedEntries(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the &#39;refValue&#39;, &#39;refGradient&#39; and &#39;valueFraction&#39; entries into their respective places...
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual scalarField & valueFraction()
A FieldMapper for finite-volume patch fields.
mappedMixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
static const word null
An empty word.
Definition: word.H:84
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
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:336
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
OBJstream os(runTime.globalPath()/outputName)
Functionality for sampling fields using mappedPatchBase. Every call to mappedField() returns a sample...
Type gAverage(const FieldField< Field, Type > &f)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual void write(Ostream &) const
Write.
volScalarField & p
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
This boundary condition maps the value at a set of cells or patch faces back to *this.
virtual Field< Type > & refGrad()
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual void write(Ostream &) const
Write.