mixedFvPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "mixedFvPatchField.H"
27 
28 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
30 template<class Type>
32 (
33  const fvPatch& p,
35 )
36 :
37  fvPatchField<Type>(p, iF),
38  refValue_(p.size()),
39  refGrad_(p.size()),
40  valueFraction_(p.size())
41 {}
42 
43 
44 template<class Type>
46 (
47  const fvPatch& p,
49  const dictionary& dict
50 )
51 :
52  fvPatchField<Type>(p, iF, dict, false),
53  refValue_("refValue", dict, p.size()),
54  refGrad_("refGradient", dict, p.size()),
55  valueFraction_("valueFraction", dict, p.size())
56 {
57  evaluate();
58 }
59 
60 
61 template<class Type>
63 (
64  const mixedFvPatchField<Type>& ptf,
65  const fvPatch& p,
67  const fvPatchFieldMapper& mapper,
68  const bool mappingRequired
69 )
70 :
71  fvPatchField<Type>(ptf, p, iF, mapper, mappingRequired),
72  refValue_(mapper(ptf.refValue_)),
73  refGrad_(mapper(ptf.refGrad_)),
74  valueFraction_(mapper(ptf.valueFraction_))
75 {}
76 
77 
78 template<class Type>
80 (
81  const mixedFvPatchField<Type>& ptf,
83 )
84 :
85  fvPatchField<Type>(ptf, iF),
86  refValue_(ptf.refValue_),
87  refGrad_(ptf.refGrad_),
88  valueFraction_(ptf.valueFraction_)
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
94 template<class Type>
96 (
97  const fvPatchFieldMapper& m
98 )
99 {
101  m(refValue_, refValue_);
102  m(refGrad_, refGrad_);
103  m(valueFraction_, valueFraction_);
104 }
105 
106 
107 template<class Type>
109 (
110  const fvPatchField<Type>& ptf,
111  const labelList& addr
112 )
113 {
114  fvPatchField<Type>::rmap(ptf, addr);
115 
116  const mixedFvPatchField<Type>& mptf =
117  refCast<const mixedFvPatchField<Type>>(ptf);
118 
119  refValue_.rmap(mptf.refValue_, addr);
120  refGrad_.rmap(mptf.refGrad_, addr);
121  valueFraction_.rmap(mptf.valueFraction_, addr);
122 }
123 
124 
125 template<class Type>
127 {
128  if (!this->updated())
129  {
130  this->updateCoeffs();
131  }
132 
134  (
135  valueFraction_*refValue_
136  +
137  (1.0 - valueFraction_)*
138  (
139  this->patchInternalField()
140  + refGrad_/this->patch().deltaCoeffs()
141  )
142  );
143 
145 }
146 
147 
148 template<class Type>
151 {
152  return
153  valueFraction_
154  *(refValue_ - this->patchInternalField())
155  *this->patch().deltaCoeffs()
156  + (1.0 - valueFraction_)*refGrad_;
157 }
158 
159 
160 template<class Type>
163 (
164  const tmp<scalarField>&
165 ) const
166 {
167  return Type(pTraits<Type>::one)*(1.0 - valueFraction_);
168 
169 }
170 
171 
172 template<class Type>
175 (
176  const tmp<scalarField>&
177 ) const
178 {
179  return
180  valueFraction_*refValue_
181  + (1.0 - valueFraction_)*refGrad_/this->patch().deltaCoeffs();
182 }
183 
184 
185 template<class Type>
188 {
189  return -Type(pTraits<Type>::one)*valueFraction_*this->patch().deltaCoeffs();
190 }
191 
192 
193 template<class Type>
196 {
197  return
198  valueFraction_*this->patch().deltaCoeffs()*refValue_
199  + (1.0 - valueFraction_)*refGrad_;
200 }
201 
202 
203 template<class Type>
205 {
207  writeEntry(os, "refValue", refValue_);
208  writeEntry(os, "refGradient", refGrad_);
209  writeEntry(os, "valueFraction", valueFraction_);
210  writeEntry(os, "value", *this);
211 }
212 
213 
214 // ************************************************************************* //
dictionary dict
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
commsTypes
Types of communications.
Definition: UPstream.H:64
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
Traits class for primitives.
Definition: pTraits.H:50
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
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.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
Pre-declare SubField and related Field type.
Definition: Field.H:56
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
Foam::fvPatchFieldMapper.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
mixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void write(Ostream &) const
Write.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.