All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
qgdFluxFvPatchScalarField.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) 2019 OpenCFD Ltd.
10  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11 -------------------------------------------------------------------------------
12 License
13  This file is part of QGDsolver library, based on OpenFOAM+.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 \*---------------------------------------------------------------------------*/
28 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "addToRunTimeSelectionTable.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
38 (
39  const fvPatch& p,
40  const DimensionedField<scalar, volMesh>& iF
41 )
42 :
43  fixedGradientFvPatchScalarField(p, iF)
44 {}
45 
46 
48 (
49  const fvPatch& p,
50  const DimensionedField<scalar, volMesh>& iF,
51  const dictionary& dict
52 )
53 :
54  fixedGradientFvPatchScalarField(p, iF)
55 {
56  patchType() = dict.lookupOrDefault<word>("patchType", word::null);
57 
58  if (dict.found("value") && dict.found("gradient"))
59  {
60  fvPatchField<scalar>::operator=
61  (
62  scalarField("value", dict, p.size())
63  );
64  gradient() = scalarField("gradient", dict, p.size());
65  }
66  else
67  {
68  fvPatchField<scalar>::operator=(patchInternalField());
69  gradient() = 0.0;
70  }
71 }
72 
73 
75 (
76  const qgdFluxFvPatchScalarField& ptf,
77  const fvPatch& p,
78  const DimensionedField<scalar, volMesh>& iF,
79  const fvPatchFieldMapper& mapper
80 )
81 :
82  fixedGradientFvPatchScalarField(p, iF)
83 {
84  patchType() = ptf.patchType();
85 
86  // Map gradient. Set unmapped values and overwrite with mapped ptf
87  gradient() = 0.0;
88  gradient().map(ptf.gradient(), mapper);
89 
90  // Evaluate the value field from the gradient if the internal field is valid
91  if (notNull(iF))
92  {
93  if (iF.size())
94  {
95  // Note: cannot ask for nf() if zero faces
96 
97  scalarField::operator=
98  (
99  //patchInternalField() + gradient()/patch().deltaCoeffs()
100  // ***HGW Hack to avoid the construction of mesh.deltaCoeffs
101  // which fails for AMI patches for some mapping operations
102  patchInternalField()
103  + gradient()*(patch().nf() & patch().delta())
104  );
105  }
106  else
107  {
108  // Enforce mapping of values so we have a valid starting value. This
109  // constructor is used when reconstructing fields
110  this->map(ptf, mapper);
111  }
112  }
113  else
114  {
115  // Enforce mapping of values so we have a valid starting value
116  this->map(ptf, mapper);
117  }
118 }
119 
120 
122 (
123  const qgdFluxFvPatchScalarField& wbppsf
124 )
125 :
126  fixedGradientFvPatchScalarField(wbppsf)
127 {}
128 
129 
131 (
132  const qgdFluxFvPatchScalarField& wbppsf,
133  const DimensionedField<scalar, volMesh>& iF
134 )
135 :
136  fixedGradientFvPatchScalarField(wbppsf, iF)
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
142 //void Foam::qgdFluxFvPatchScalarField::updateSnGrad
143 //(
144 // const scalarField& snGradp
145 //)
146 //{
147 // if (updated())
148 // {
149 // return;
150 // }
151 //
152 // curTimeIndex_ = this->db().time().timeIndex();
153 //
154 // gradient() = snGradp;
155 // fixedGradientFvPatchScalarField::updateCoeffs();
156 // }
157 
158 
160 {
161  if (updated())
162  {
163  return;
164  }
165 
166  if (this->patch().boundaryMesh().mesh().thisDb().
167  foundObject<surfaceScalarField>("phiwStar")
168  )
169  {
170  const surfaceScalarField& phiws =
171  this->patch().boundaryMesh().mesh().
172  thisDb().lookupObject<surfaceScalarField>
173  ("phiwStar");
174 
175  if (this->patch().boundaryMesh().mesh().thisDb().
176  foundObject<surfaceScalarField>("tauQGDf")
177  )
178  {
179  const surfaceScalarField& tauQGDf =
180  this->patch().boundaryMesh().mesh().
181  thisDb().lookupObject<surfaceScalarField>
182  ("tauQGDf");
183 
184  scalarField fluxSnGrad
185  (
186  phiws.boundaryField()[patch().index()]
187  /
188  tauQGDf.boundaryField()[patch().index()]
189  /
190  patch().magSf()
191  );
192  this->gradient() = -fluxSnGrad;
193 
194  }
195  }
196 
197  fixedGradientFvPatchScalarField::updateCoeffs();
198 
199  // if (curTimeIndex_ != this->db().time().timeIndex())
200  // {
201  // FatalErrorInFunction
202  // << "updateCoeffs(const scalarField& snGradp) MUST be called before"
203  // " updateCoeffs() or evaluate() to set the boundary gradient."
204  // << exit(FatalError);
205  // }
206 
207 
208 }
209 
210 
211 void Foam::qgdFluxFvPatchScalarField::write(Ostream& os) const
212 {
213  fixedGradientFvPatchScalarField::write(os);
214  writeEntry("value", os);
215 }
216 
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 namespace Foam
221 {
223  (
224  fvPatchScalarField,
225  qgdFluxFvPatchScalarField
226  );
227 }
228 
229 
230 // ************************************************************************* //
makePatchTypeField(fvPatchVectorField, cosVelocityFvPatchVectorField)
const surfaceScalarField & tauQGDf
Definition: createFields.H:55
virtual void write(Ostream &) const
Write.
qgdFluxFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
This boundary condition sets the pressure gradient to the provided value such that the flux on the bo...
virtual void updateCoeffs()
Update the patch pressure gradient field.
volScalarField & p
Definition: createFields.H:52