alphaFixedPressureFvPatchScalarField.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 | Copyright (C) 2011-2016 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 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40 )
41 :
42  fixedValueFvPatchScalarField(p, iF),
43  p_(p.size(), 0.0)
44 {}
45 
46 
49 (
51  const fvPatch& p,
53  const fvPatchFieldMapper& mapper
54 )
55 :
56  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
57  p_(ptf.p_, mapper)
58 {}
59 
60 
63 (
64  const fvPatch& p,
66  const dictionary& dict
67 )
68 :
69  fixedValueFvPatchScalarField(p, iF, dict, false),
70  p_("p", dict, p.size())
71 {
72  if (dict.found("value"))
73  {
75  (
76  scalarField("value", dict, p.size())
77  );
78  }
79  else
80  {
82  }
83 }
84 
85 
88 (
90 )
91 :
92  fixedValueFvPatchScalarField(tppsf),
93  p_(tppsf.p_)
94 {}
95 
96 
99 (
102 )
103 :
104  fixedValueFvPatchScalarField(tppsf, iF),
105  p_(tppsf.p_)
106 {}
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
112 (
113  const fvPatchFieldMapper& m
114 )
115 {
116  scalarField::autoMap(m);
117  p_.autoMap(m);
118 }
119 
120 
122 (
123  const fvPatchScalarField& ptf,
124  const labelList& addr
125 )
126 {
127  fixedValueFvPatchScalarField::rmap(ptf, addr);
128 
130  refCast<const alphaFixedPressureFvPatchScalarField>(ptf);
131 
132  p_.rmap(tiptf.p_, addr);
133 }
134 
135 
137 {
138  if (updated())
139  {
140  return;
141  }
142 
144  db().lookupObject<uniformDimensionedVectorField>("g");
145 
146  const fvPatchField<scalar>& rho =
147  patch().lookupPatchField<volScalarField, scalar>("rho");
148 
149  operator==(p_ - rho*(g.value() & patch().Cf()));
150 
151  fixedValueFvPatchScalarField::updateCoeffs();
152 }
153 
154 
156 (
157  Ostream& os
158 ) const
159 {
161  p_.writeEntry("p", os);
162  writeEntry("value", os);
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 namespace Foam
169 {
171  (
174  );
175 }
176 
177 // ************************************************************************* //
Foam::surfaceFields.
alphaFixedPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:726
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::fvPatchFieldMapper.
const Type & value() const
Return const reference to value.
virtual label size() const
Return size.
Definition: fvPatch.H:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
volScalarField scalarField(fieldObject, mesh)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
const dimensionedVector & g
A fixed-pressure alphaContactAngle boundary.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Namespace for OpenFOAM.