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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 }
100 
101 
103 (
105 )
106 :
107  fixedGradientFvPatchScalarField(wbppsf),
108  curTimeIndex_(-1)
109 {}
110 
111 
113 (
116 )
117 :
118  fixedGradientFvPatchScalarField(wbppsf, iF),
119  curTimeIndex_(-1)
120 {}
121 
122 
123 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124 
126 (
127  const scalarField& snGradp
128 )
129 {
130  if (updated())
131  {
132  return;
133  }
134 
135  curTimeIndex_ = this->db().time().timeIndex();
136 
137  gradient() = snGradp;
138  fixedGradientFvPatchScalarField::updateCoeffs();
139 }
140 
141 
143 {
144  if (updated())
145  {
146  return;
147  }
148 
149  if (curTimeIndex_ != this->db().time().timeIndex())
150  {
152  << "updateCoeffs(const scalarField& snGradp) MUST be called before"
153  " updateCoeffs() or evaluate() to set the boundary gradient."
154  << exit(FatalError);
155  }
156 }
157 
158 
160 {
162  writeEntry("value", os);
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 namespace Foam
169 {
171  (
174  );
175 }
176 
177 
178 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
p
p
Definition: pEqn.H:62
Foam::notNull
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:46
Foam::fixedFluxPressureFvPatchScalarField::curTimeIndex_
label curTimeIndex_
Current time index (used for updating)
Definition: fixedFluxPressureFvPatchScalarField.H:72
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::fixedFluxPressureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the patch pressure gradient field.
Definition: fixedFluxPressureFvPatchScalarField.C:142
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::fixedFluxPressureFvPatchScalarField
This boundary condition sets the pressure gradient to the provided value such that the flux on the bo...
Definition: fixedFluxPressureFvPatchScalarField.H:65
Foam::fixedFluxPressureFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: fixedFluxPressureFvPatchScalarField.C:159
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
fixedFluxPressureFvPatchScalarField.H
Foam::fixedFluxPressureFvPatchScalarField::fixedFluxPressureFvPatchScalarField
fixedFluxPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: fixedFluxPressureFvPatchScalarField.C:35
write
Tcoeff write()
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51