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-2020 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"
30 #include "interpolationCell.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class Type>
36 (
37  const fvPatch& p,
39 )
40 :
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 
52 template<class Type>
54 (
55  const fvPatch& p,
57  const dictionary& dict
58 )
59 :
61  mappedPatchBase(p.patch(), dict),
62  mappedPatchFieldBase<Type>(*this, *this, dict),
63  weightFieldName_(dict.getOrDefault<word>("weightField", word::null))
64 {
66  (
67  Field<Type>("value", dict, p.size())
68  );
69 
70  if (dict.found("refValue"))
71  {
72  // Full restart
73  this->refValue() = Field<Type>("refValue", dict, p.size());
74  this->refGrad() = Field<Type>("refGradient", dict, p.size());
75  this->valueFraction() = scalarField("valueFraction", dict, p.size());
76  }
77  else
78  {
79  // Start from user entered data. Assume fixedValue.
80  this->refValue() = *this;
81  this->refGrad() = Zero;
82  this->valueFraction() = 1.0;
83  }
84 
85  // Store patch value as initial guess when running in database mode
87  (
88  this->internalField().name(),
89  *this
90  );
92  (
93  this->internalField().name() + "_weights",
94  this->patch().deltaCoeffs()
95  );
96 }
97 
98 
99 template<class Type>
101 (
103  const fvPatch& p,
105  const fvPatchFieldMapper& mapper
106 )
107 :
108  mixedFvPatchField<Type>(ptf, p, iF, mapper),
109  mappedPatchBase(p.patch(), ptf),
110  mappedPatchFieldBase<Type>(*this, *this, ptf),
111  weightFieldName_(ptf.weightFieldName_)
112 {}
113 
114 
115 template<class Type>
117 (
119 )
120 :
122  mappedPatchBase(ptf.patch().patch(), ptf),
124  weightFieldName_(ptf.weightFieldName_)
125 {}
126 
127 
128 template<class Type>
130 (
133 )
134 :
135  mixedFvPatchField<Type>(ptf, iF),
136  mappedPatchBase(ptf.patch().patch(), ptf),
137  mappedPatchFieldBase<Type>(*this, *this, ptf),
138  weightFieldName_(ptf.weightFieldName_)
139 {}
140 
141 
142 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
143 
144 template<class Type>
146 (
147  const fvPatchFieldMapper& m
148 )
149 {
151  mappedPatchBase::clearOut();
152 }
153 
154 
155 template<class Type>
157 (
158  const fvPatchField<Type>& ptf,
159  const labelList& addr
160 )
161 {
163  mappedPatchBase::clearOut();
164 }
165 
166 
167 template<class Type>
169 {
170  if (this->updated())
171  {
172  return;
173  }
174 
175  const tmp<Field<Type>> nbrIntFld(this->mappedInternalField());
176 
177  //- Unweighted
178  //const tmp<scalarField> nbrKDelta(this->mappedWeightField());
179 
180  //- Weighted
181  tmp<scalarField> myKDelta;
182  tmp<scalarField> nbrKDelta;
183  this->mappedWeightField(weightFieldName_, myKDelta, nbrKDelta);
184 
185  // Both sides agree on
186  // - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
187  // - gradient : (temperature-fld)*delta
188  // We've got a degree of freedom in how to implement this in a mixed bc.
189  // (what gradient, what fixedValue and mixing coefficient)
190  // Two reasonable choices:
191  // 1. specify above temperature on one side (preferentially the high side)
192  // and above gradient on the other. So this will switch between pure
193  // fixedvalue and pure fixedgradient
194  // 2. specify gradient and temperature such that the equations are the
195  // same on both sides. This leads to the choice of
196  // - refGradient = zero gradient
197  // - refValue = neighbour value
198  // - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
199 
200  this->refValue() = nbrIntFld;
201  this->refGrad() = Zero;
202  this->valueFraction() = nbrKDelta()/(nbrKDelta() + myKDelta());
203 
205 
206  if (debug)
207  {
208  Info<< this->patch().boundaryMesh().mesh().name() << ':'
209  << this->patch().name() << ':'
210  << this->internalField().name() << " <- "
211  << this->mapper_.sampleRegion() << ':'
212  << this->mapper_.samplePatch() << ':'
213  << this->fieldName_ << " :"
214  << " value "
215  << " min:" << gMin(*this)
216  << " max:" << gMax(*this)
217  << " avg:" << gAverage(*this)
218  << endl;
219  }
220 }
221 
222 
223 template<class Type>
225 {
228  os.writeEntryIfDifferent<word>("weightField", word::null, weightFieldName_);
230 }
231 
232 
233 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
volFields.H
Foam::mappedMixedFieldFvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: mappedMixedFieldFvPatchField.C:146
Foam::Ostream::writeEntryIfDifferent
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:244
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
interpolationCell.H
Foam::mappedMixedFieldFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: mappedMixedFieldFvPatchField.C:224
Foam::mappedMixedFieldFvPatchField::mappedMixedFieldFvPatchField
mappedMixedFieldFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: mappedMixedFieldFvPatchField.C:36
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:112
mappedMixedFieldFvPatchField.H
Foam::mappedMixedFieldFvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: mappedMixedFieldFvPatchField.C:157
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::mappedPatchFieldBase
Functionality for sampling fields using mappedPatchBase. Every call to mappedField() returns a sample...
Definition: mappedPatchFieldBase.H:104
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::mixedFvPatchField
This boundary condition provides a base class for 'mixed' type boundary conditions,...
Definition: mixedFvPatchField.H:123
Foam::mappedMixedFieldFvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: mappedMixedFieldFvPatchField.C:168
Foam::mappedMixedFieldFvPatchField
This boundary condition provides a self-contained version of e.g. mapped boundary conditions.
Definition: mappedMixedFieldFvPatchField.H:129
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:343
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54