pressureFvPatchScalarField.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) 2018-2020 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"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 (
36  const fvPatch& p,
38 )
39 :
40  fixedValueFvPatchScalarField(p, iF),
41  p_(p.size(), 0.0)
42 {}
43 
44 
46 (
47  const fvPatch& p,
49  const dictionary& dict
50 )
51 :
52  fixedValueFvPatchScalarField(p, iF, dict, false),
53  p_("p", dict, p.size())
54 {
55  if (dict.found("value"))
56  {
57  fvPatchScalarField::operator=
58  (
59  scalarField("value", dict, p.size())
60  );
61  }
62  else
63  {
65  }
66 }
67 
68 
70 (
71  const pressureFvPatchScalarField& ptf,
72  const fvPatch& p,
74  const fvPatchFieldMapper& mapper
75 )
76 :
77  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
78  p_(mapper(ptf.p_))
79 {}
80 
81 
83 (
84  const pressureFvPatchScalarField& ptf,
86 )
87 :
88  fixedValueFvPatchScalarField(ptf, iF),
89  p_(ptf.p_)
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 (
97  const fvPatchFieldMapper& m
98 )
99 {
100  fixedValueFvPatchScalarField::autoMap(m);
101  m(p_, p_);
102 }
103 
104 
106 (
107  const fvPatchScalarField& ptf,
108  const labelList& addr
109 )
110 {
111  fixedValueFvPatchScalarField::rmap(ptf, addr);
112 
113  const pressureFvPatchScalarField& tiptf =
114  refCast<const pressureFvPatchScalarField>(ptf);
115 
116  p_.rmap(tiptf.p_, addr);
117 }
118 
119 
121 {
122  if (updated())
123  {
124  return;
125  }
126 
127  operator==(p_);
128 
129  fixedValueFvPatchScalarField::updateCoeffs();
130 }
131 
132 
134 {
136  writeEntry(os, "p", p_);
137  writeEntry(os, "value", *this);
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142 
143 namespace Foam
144 {
146  (
149  );
150 }
151 
152 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
pressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:260
Macros for easy insertion into run-time selection tables.
virtual void write(Ostream &) const
Write.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:156
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
volScalarField scalarField(fieldObject, mesh)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Static pressure boundary condition.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.