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 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"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const fvPatch& p,
40 )
41 :
44  (
46  *this
47  ),
48  weightFieldName_(word::null)
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 :
66  (
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  // Store patch value as initial guess when running in database mode
95  (
96  this->internalField().name(),
97  *this
98  );
100  (
101  this->internalField().name() + "_weights",
102  this->patch().deltaCoeffs()
103  );
104 }
105 
106 
107 template<class Type>
109 (
111  const fvPatch& p,
113  const fvPatchFieldMapper& mapper
114 )
115 :
116  mixedFvPatchField<Type>(ptf, p, iF, mapper),
118  (
120  *this,
121  ptf
122  ),
123  weightFieldName_(ptf.weightFieldName_)
124 {}
125 
126 
127 template<class Type>
129 (
131 )
132 :
135  weightFieldName_(ptf.weightFieldName_)
136 {}
137 
138 
139 template<class Type>
141 (
144 )
145 :
146  mixedFvPatchField<Type>(ptf, iF),
148  (
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 
232 
233 template<class Type>
235 {
237  os.writeEntryIfDifferent<word>("weightField", word::null, weightFieldName_);
239 }
240 
241 
242 // ************************************************************************* //
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::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::mappedMixedFvPatchField::mappedMixedFvPatchField
mappedMixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: mappedMixedFvPatchField.C:37
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::mappedFixedValueFvPatchField
This boundary condition maps the value at a set of cells or patch faces back to *this.
Definition: mappedFixedValueFvPatchField.H:125
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::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::mappedMixedFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: mappedMixedFvPatchField.C:234
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::mappedMixedFvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: mappedMixedFvPatchField.C:178
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
mappedFixedValueFvPatchField.H
Foam::mappedMixedFvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: mappedMixedFvPatchField.C:168
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::mappedMixedFvPatchField
This boundary condition maps the value at a set of cells or patch faces back to *this.
Definition: mappedMixedFvPatchField.H:128
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::mappedMixedFvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: mappedMixedFvPatchField.C:160
Foam::mixedFvPatchField
This boundary condition provides a base class for 'mixed' type boundary conditions,...
Definition: mixedFvPatchField.H:123
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
mappedMixedFvPatchField.H
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