mappedFixedValueFvPatchField.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) 2011-2016 OpenFOAM Foundation
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 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const fvPatch& p,
38 )
39 :
41  mappedPatchFieldBase<Type>(mappedPatchFieldBase<Type>::mapper(p, iF), *this)
42 {}
43 
44 
45 template<class Type>
47 (
48  const fvPatch& p,
50  const dictionary& dict
51 )
52 :
53  fixedValueFvPatchField<Type>(p, iF, dict),
55  (
56  mappedPatchFieldBase<Type>::mapper(p, iF),
57  *this,
58  dict,
59  *this // initial value for database operation
60  )
61 {}
62 
63 
64 template<class Type>
66 (
68  const fvPatch& p,
70  const fvPatchFieldMapper& mapper
71 )
72 :
73  fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
75  (
76  mappedPatchFieldBase<Type>::mapper(p, iF),
77  *this,
78  ptf
79  )
80 {}
81 
82 
83 template<class Type>
85 (
87 )
88 :
90  mappedPatchFieldBase<Type>(ptf)
91 {}
92 
93 
94 template<class Type>
96 (
99 )
100 :
101  fixedValueFvPatchField<Type>(ptf, iF),
103  (
104  mappedPatchFieldBase<Type>::mapper(this->patch(), iF),
105  *this,
106  ptf
107  )
108 {}
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
113 template<class Type>
115 {
116  if (this->updated())
117  {
118  return;
119  }
120 
121  this->operator==(this->mappedField());
122 
123  if (debug)
124  {
125  Info<< "mapped on field:"
126  << this->internalField().name()
127  << " patch:" << this->patch().name()
128  << " avg:" << gAverage(*this)
129  << " min:" << gMin(*this)
130  << " max:" << gMax(*this)
131  << endl;
132  }
133 
135 }
136 
137 
138 template<class Type>
140 {
144 }
145 
146 
147 // ************************************************************************* //
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
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.
mappedFixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
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 maps the value at a set of cells or patch faces back to *this.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
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)
volScalarField & p
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)