fixedFluxPressureFvPatchScalarField.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-2018 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 
27 #include "fvPatchFieldMapper.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 (
36  const fvPatch& p,
38 )
39 :
40  fixedGradientFvPatchScalarField(p, iF),
41  curTimeIndex_(-1)
42 {}
43 
44 
46 (
47  const fvPatch& p,
49  const dictionary& dict
50 )
51 :
52  fixedGradientFvPatchScalarField(p, iF),
53  curTimeIndex_(-1)
54 {
55  if (dict.found("value") && dict.found("gradient"))
56  {
58  (
59  scalarField("value", dict, p.size())
60  );
61  gradient() = scalarField("gradient", dict, p.size());
62  }
63  else
64  {
65  fvPatchField<scalar>::operator=(patchInternalField());
66  gradient() = 0.0;
67  }
68 }
69 
70 
72 (
74  const fvPatch& p,
76  const fvPatchFieldMapper& mapper
77 )
78 :
79  fixedGradientFvPatchScalarField(p, iF),
80  curTimeIndex_(-1)
81 {
82  patchType() = ptf.patchType();
83 
84  // Map gradient. Set unmapped values and overwrite with mapped ptf
85  gradient() = 0.0;
86  gradient().map(ptf.gradient(), mapper);
87 
88  // Evaluate the value field from the gradient if the internal field is valid
89  if (notNull(iF) && iF.size())
90  {
91  scalarField::operator=
92  (
93  // patchInternalField() + gradient()/patch().deltaCoeffs()
94  // ***HGW Hack to avoid the construction of mesh.deltaCoeffs
95  // which fails for AMI patches for some mapping operations
96  patchInternalField() + gradient()*(patch().nf() & patch().delta())
97  );
98  }
99  else
100  {
101  // Enforce mapping of values so we have a valid starting value. This
102  // constructor is used when reconstructing fields
103  this->map(ptf, mapper);
104  }
105 }
106 
107 
109 (
111 )
112 :
113  fixedGradientFvPatchScalarField(wbppsf),
114  curTimeIndex_(-1)
115 {}
116 
117 
119 (
122 )
123 :
124  fixedGradientFvPatchScalarField(wbppsf, iF),
125  curTimeIndex_(-1)
126 {}
127 
128 
129 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130 
132 (
133  const scalarField& snGradp
134 )
135 {
136  if (updated())
137  {
138  return;
139  }
140 
141  curTimeIndex_ = this->db().time().timeIndex();
142 
143  gradient() = snGradp;
144  fixedGradientFvPatchScalarField::updateCoeffs();
145 }
146 
147 
149 {
150  if (updated())
151  {
152  return;
153  }
154 
155  if (curTimeIndex_ != this->db().time().timeIndex())
156  {
158  << "updateCoeffs(const scalarField& snGradp) MUST be called before"
159  " updateCoeffs() or evaluate() to set the boundary gradient."
160  << exit(FatalError);
161  }
162 }
163 
164 
166 {
168  writeEntry("value", os);
169 }
170 
171 
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 
174 namespace Foam
175 {
177  (
180  );
181 }
182 
183 
184 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
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
Macros for easy insertion into run-time selection tables.
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
virtual void updateCoeffs()
Update the patch pressure gradient field.
Foam::fvPatchFieldMapper.
This boundary condition sets the pressure gradient to the provided value such that the flux on the bo...
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)
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:46
fixedFluxPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
label timeIndex
Definition: getTimeIndex.H:4
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.