Loading...
Searching...
No Matches
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),
43 mappedPatchFieldBase<Type>(*this, *this),
44 weightFieldName_()
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),
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),
78 mappedPatchFieldBase<Type>(*this, *this, ptf),
79 weightFieldName_(ptf.weightFieldName_)
80{}
81
82
83template<class Type>
85(
87)
88:
89 mixedFvPatchField<Type>(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)
126 const fvPatchField<Type>& ptf,
127 const labelList& addr
128)
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 auto limits = gMinMax(*this);
177 auto avg = gAverage(*this);
178
179 Info<< this->patch().boundaryMesh().mesh().name() << ':'
180 << this->patch().name() << ':'
181 << this->internalField().name() << " <- "
182 << this->mapper_.sampleRegion() << ':'
183 << this->mapper_.samplePatch() << ':'
184 << this->fieldName_ << " :"
185 << " value "
186 << " min:" << limits.min()
187 << " max:" << limits.max()
188 << " avg:" << avg
189 << endl;
190 }
191}
192
194template<class Type>
196{
199 os.writeEntryIfDifferent<word>("weightField", word::null, weightFieldName_);
201}
202
203
204// ************************************************************************* //
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
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
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.
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 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.
mappedMixedFieldFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
virtual void write(Ostream &os) const
Write as a dictionary.
mappedPatchBase(const polyPatch &)
Construct from patch.
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.
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