fixedFluxPressureFvPatchScalarField.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  Copyright (C) 2015-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
38 (
39  const fvPatch& p,
41 )
42 :
43  fixedGradientFvPatchScalarField(p, iF),
44  curTimeIndex_(-1)
45 {}
46 
47 
49 (
50  const fvPatch& p,
52  const dictionary& dict
53 )
54 :
55  fixedGradientFvPatchScalarField(p, iF),
56  curTimeIndex_(-1)
57 {
58  patchType() = dict.getOrDefault<word>("patchType", word::null);
59  if (dict.found("value") && dict.found("gradient"))
60  {
62  (
63  scalarField("value", dict, p.size())
64  );
65  gradient() = scalarField("gradient", dict, p.size());
66  }
67  else
68  {
69  fvPatchField<scalar>::operator=(patchInternalField());
70  gradient() = 0.0;
71  }
72 }
73 
74 
76 (
78  const fvPatch& p,
80  const fvPatchFieldMapper& mapper
81 )
82 :
83  fixedGradientFvPatchScalarField(p, iF),
84  curTimeIndex_(-1)
85 {
86  patchType() = ptf.patchType();
87 
88  // Map gradient. Set unmapped values and overwrite with mapped ptf
89  gradient() = 0.0;
90  gradient().map(ptf.gradient(), mapper);
91 
92  // Evaluate the value field from the gradient if the internal field is valid
93  if (notNull(iF))
94  {
95  if (iF.size())
96  {
97  // Note: cannot ask for nf() if zero faces
98 
99  scalarField::operator=
100  (
101  //patchInternalField() + gradient()/patch().deltaCoeffs()
102  // ***HGW Hack to avoid the construction of mesh.deltaCoeffs
103  // which fails for AMI patches for some mapping operations
104  patchInternalField()
105  + gradient()*(patch().nf() & patch().delta())
106  );
107  }
108  }
109  else
110  {
111  // Enforce mapping of values so we have a valid starting value
112  this->map(ptf, mapper);
113  }
114 }
115 
116 
118 (
120 )
121 :
122  fixedGradientFvPatchScalarField(wbppsf),
123  curTimeIndex_(-1)
124 {}
125 
126 
128 (
131 )
132 :
133  fixedGradientFvPatchScalarField(wbppsf, iF),
134  curTimeIndex_(-1)
135 {}
136 
137 
138 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139 
141 (
142  const scalarField& snGradp
143 )
144 {
145  if (updated())
146  {
147  return;
148  }
149 
150  curTimeIndex_ = this->db().time().timeIndex();
151 
152  gradient() = snGradp;
153  fixedGradientFvPatchScalarField::updateCoeffs();
154 }
155 
156 
158 {
159  if (updated())
160  {
161  return;
162  }
163 
164  if (curTimeIndex_ != this->db().time().timeIndex())
165  {
167  << "updateCoeffs(const scalarField& snGradp) MUST be called before"
168  " updateCoeffs() or evaluate() to set the boundary gradient."
169  << exit(FatalError);
170  }
171 }
172 
173 
175 {
177  writeEntry("value", os);
178 }
179 
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 namespace Foam
184 {
186  (
189  );
190 }
191 
192 
193 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
volFields.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::fixedFluxPressureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the patch pressure gradient field.
Definition: fixedFluxPressureFvPatchScalarField.C:157
Foam::Field< scalar >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::notNull
bool notNull(const T *ptr)
True if ptr is not a pointer (of type T) to the nullObject.
Definition: nullObject.H:207
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::fixedFluxPressureFvPatchScalarField::updateSnGrad
virtual void updateSnGrad(const scalarField &snGradp)
Update the patch pressure gradient field from the given snGradp.
Definition: fixedFluxPressureFvPatchScalarField.C:141
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::fixedFluxPressureFvPatchScalarField
This boundary condition sets the pressure gradient to the provided value such that the flux on the bo...
Definition: fixedFluxPressureFvPatchScalarField.H:69
Foam::fixedFluxPressureFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: fixedFluxPressureFvPatchScalarField.C:174
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
timeIndex
label timeIndex
Definition: getTimeIndex.H:30
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
fixedFluxPressureFvPatchScalarField.H
Foam::fixedFluxPressureFvPatchScalarField::fixedFluxPressureFvPatchScalarField
fixedFluxPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: fixedFluxPressureFvPatchScalarField.C:38
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54