All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
mQhdFluxFvPatchScalarField.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 \*---------------------------------------------------------------------------*/
29 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "surfaceFields.H"
34 #include "addToRunTimeSelectionTable.H"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const fvPatch& p,
41  const DimensionedField<scalar, volMesh>& iF
42 )
43 :
44  fixedGradientFvPatchScalarField(p, iF)
45 {}
46 
47 
49 (
50  const fvPatch& p,
51  const DimensionedField<scalar, volMesh>& iF,
52  const dictionary& dict
53 )
54 :
55  fixedGradientFvPatchScalarField(p, iF)
56 {
57  patchType() = dict.lookupOrDefault<word>("patchType", word::null);
58 
59  if (dict.found("value") && dict.found("gradient"))
60  {
61  fvPatchField<scalar>::operator=
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 (
77  const mQhdFluxFvPatchScalarField& ptf,
78  const fvPatch& p,
79  const DimensionedField<scalar, volMesh>& iF,
80  const fvPatchFieldMapper& mapper
81 )
82 :
83  fixedGradientFvPatchScalarField(p, iF)
84 {
85  patchType() = ptf.patchType();
86 
87  // Map gradient. Set unmapped values and overwrite with mapped ptf
88  gradient() = 0.0;
89  gradient().map(ptf.gradient(), mapper);
90 
91  // Evaluate the value field from the gradient if the internal field is valid
92  if (notNull(iF))
93  {
94  if (iF.size())
95  {
96  // Note: cannot ask for nf() if zero faces
97 
98  scalarField::operator=
99  (
100  //patchInternalField() + gradient()/patch().deltaCoeffs()
101  // ***HGW Hack to avoid the construction of mesh.deltaCoeffs
102  // which fails for AMI patches for some mapping operations
103  patchInternalField()
104  + gradient()*(patch().nf() & patch().delta())
105  );
106  }
107  else
108  {
109  // Enforce mapping of values so we have a valid starting value. This
110  // constructor is used when reconstructing fields
111  this->map(ptf, mapper);
112  }
113  }
114  else
115  {
116  // Enforce mapping of values so we have a valid starting value
117  this->map(ptf, mapper);
118  }
119 }
120 
121 
123 (
124  const mQhdFluxFvPatchScalarField& wbppsf
125 )
126 :
127  fixedGradientFvPatchScalarField(wbppsf)
128 {}
129 
130 
132 (
133  const mQhdFluxFvPatchScalarField& wbppsf,
134  const DimensionedField<scalar, volMesh>& iF
135 )
136 :
137  fixedGradientFvPatchScalarField(wbppsf, iF)
138 {}
139 
140 
141 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142 
143 //void Foam::mQhdFluxFvPatchScalarField::updateSnGrad
144 //(
145 // const scalarField& snGradp
146 //)
147 //{
148 // if (updated())
149 // {
150 // return;
151 // }
152 //
153 // curTimeIndex_ = this->db().time().timeIndex();
154 //
155 // gradient() = snGradp;
156 // fixedGradientFvPatchScalarField::updateCoeffs();
157 // }
158 
159 
161 {
162  if (updated())
163  {
164  return;
165  }
166 
167  if (this->patch().boundaryMesh().mesh().thisDb().
168  foundObject<surfaceScalarField>("phiwm")
169  )
170  {
171  const surfaceScalarField& phiwm =
172  this->patch().boundaryMesh().mesh().
173  thisDb().lookupObject<surfaceScalarField>
174  ("phiwm");
175 
176  if (this->patch().boundaryMesh().mesh().thisDb().
177  foundObject<surfaceScalarField>("coeffp")
178  )
179  {
180  const surfaceScalarField& coeffp =
181  this->patch().boundaryMesh().mesh().
182  thisDb().lookupObject<surfaceScalarField>
183  ("coeffp");
184 
185  scalarField fluxSnGrad
186  (
187  phiwm.boundaryField()[patch().index()]
188  /
189  coeffp.boundaryField()[patch().index()]
190  /
191  patch().magSf()
192  );
193  this->gradient() = fluxSnGrad;
194  }
195  }
196 
197  fixedGradientFvPatchScalarField::updateCoeffs();
198 }
199 
200 
201 void Foam::mQhdFluxFvPatchScalarField::write(Ostream& os) const
202 {
203  fixedGradientFvPatchScalarField::write(os);
204  writeEntry("value", os);
205 }
206 
207 
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 
210 namespace Foam
211 {
213  (
214  fvPatchScalarField,
215  mQhdFluxFvPatchScalarField
216  );
217 }
218 
219 
220 // ************************************************************************* //
makePatchTypeField(fvPatchVectorField, cosVelocityFvPatchVectorField)
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the patch pressure gradient field.
mQhdFluxFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
surfaceScalarField phiwm("phiwm", phiwo1 *0.0)
This boundary condition sets the pressure gradient to the provided value such that the flux on the bo...
surfaceScalarField coeffp("coeffp",)
volScalarField & p
Definition: createFields.H:52