All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
varScModel7.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 #include "varScModel7.H"
30 #include "psiQGDThermo.H"
31 #include "addToRunTimeSelectionTable.H"
32 #include "fvcSmooth.H"
33 #include "fvcGrad.H"
34 #include "fvcSnGrad.H"
35 #include "primitiveMeshTools.H"
36 #include "linear.H"
37 #include "cellSet.H"
38 #include "emptyFvPatch.H"
39 #include "wedgeFvPatch.H"
40 
41 namespace Foam
42 {
43 namespace qgd
44 {
45  defineTypeNameAndDebug(varScModel7,0);
47  (
48  QGDCoeffs,
49  varScModel7,
50  dictionary
51  );
52 }
53 }
54 
57 (
58  const IOobject& io,
59  const fvMesh& mesh,
60  const dictionary& dict
61 )
62 :
63  QGDCoeffs(io, mesh, dict),
64  smoothCoeff_(0.1),
65  rC_(0.5),
66  minSc_(-1.0),
67  maxSc_(-1.0),
68  badQualitySc_(0.05),
69  qgdAspectRatioThreshold_(1.5),
70  constSc_(0.5),
71  cSc1_(1.0),
72  cqSc_(mesh.V().size(), 0.0),
73  constScCellSetPtr_(nullptr)
74 {
75  scalar ScQGD = 0.0, PrQGD = 1.0;
76 
77  dict.lookup("ScQGD") >> ScQGD;
78  dict.lookup("PrQGD") >> PrQGD;
79 
80  ScQGD_.primitiveFieldRef() = ScQGD;
81  PrQGD_.primitiveFieldRef() = PrQGD;
82 
83  ScQGD_.boundaryFieldRef() = ScQGD;
84  PrQGD_.boundaryFieldRef() = PrQGD;
85 
86  if (dict.found("smoothCoeff"))
87  {
88  dict.lookup("smoothCoeff") >> smoothCoeff_;
89  }
90 
91  if (dict.found("rC"))
92  {
93  dict.lookup("rC") >> rC_;
94  }
95 
96  if (dict.found("minSc"))
97  {
98  dict.lookup("minSc") >> minSc_;
99  }
100 
101  if (dict.found("maxSc"))
102  {
103  dict.lookup("maxSc") >> maxSc_;
104  }
105 
106  if (dict.found("badQualitySc"))
107  {
108  dict.lookup("badQualitySc") >> badQualitySc_;
109  }
110 
111  if (dict.found("maxAspectRatio"))
112  {
113  dict.lookup("maxAspectRatio") >> qgdAspectRatioThreshold_;
114  }
115 
116  if (dict.found("cSc1"))
117  {
118  dict.lookup("cSc1") >> cSc1_;
119  }
120 
121  scalarField openness(mesh.V().size(), 0);
122  scalarField aspectRatio(mesh.V().size(), 1);
123 
124  primitiveMeshTools::cellClosedness
125  (
126  mesh_,
127  mesh_.geometricD(),
128  mesh_.faceAreas(),
129  mesh_.cellVolumes(),
130  openness,
131  aspectRatio
132  );
133 
134  forAll(aspectRatio, iCell)
135  {
136  if (aspectRatio[iCell] > qgdAspectRatioThreshold_)
137  {
138  cqSc_[iCell] = badQualitySc_ * aspectRatio[iCell] / qgdAspectRatioThreshold_;
139  //cqSc_[iCell] = badQualitySc_;
140  }
141  }
142 
143  if (dict.found("constScCellSet"))
144  {
145  word constScCellSet (dict.lookup("constScCellSet"));
146  constSc_ = ScQGD;
147 
148  constScCellSetPtr_.reset
149  (
150  new cellSet
151  (
152  mesh,
153  constScCellSet,
154  IOobject::MUST_READ,
155  IOobject::NO_WRITE
156  )
157  );
158  }
159 }
163 {
164 }
165 
166 void Foam::qgd::
167 varScModel7::correct(const Foam::QGDThermo& qgdThermo)
168 {
169  const volScalarField& cSound = qgdThermo.c();
170  const volScalarField& p = qgdThermo.p();
171  const volScalarField rho(qgdThermo.rho()*1.0);
172 
173  this->tauQGDf_= linearInterpolate(this->aQGD_ / cSound) * hQGDf_;
174  this->tauQGD_ = this->aQGD_ * this->hQGD_ / cSound;
175 
176  const surfaceScalarField pf =
177  linearInterpolate(p);
178  const surfaceScalarField dpf =
179  fvc::snGrad(p)/mesh_.deltaCoeffs();
180 
181  forAll(ScQGD_, icell)
182  {
183  const cell& facesOfCell = mesh_.cells()[icell];
184  ScQGD_[icell] = 0.0;
185  scalar sumDpF = 0.0;
186  scalar sumpf = 0.0;
187  label faceid = -1;
188  label patchid = -1;
189  scalar n = 0.0;
190 
191  forAll(facesOfCell, iface)
192  {
193  faceid = facesOfCell[iface];
194  if (mesh_.isInternalFace(faceid))
195  {
196  sumpf += pf.primitiveField()[faceid];
197  if (mesh_.owner()[faceid] == icell)
198  {
199  sumDpF += dpf.primitiveField()[faceid];
200  }
201  else
202  {
203  sumDpF -= dpf.primitiveField()[faceid];
204  }
205  n = n + 1;
206  }
207  else
208  {
209  patchid = mesh_.boundaryMesh().whichPatch(faceid);
210  if (isA<emptyFvPatch>(mesh_.boundary()[patchid]))
211  {
212  continue;
213  }
214  if (isA<wedgeFvPatch>(mesh_.boundary()[patchid]))
215  {
216  continue;
217  }
218  faceid -= mesh_.boundaryMesh()[patchid].start();
219 
220  if (p.boundaryField()[patchid].coupled())
221  {
222  sumDpF +=
223  (p.boundaryField()[patchid].patchNeighbourField()()[faceid] - p[icell]);
224  }
225  else
226  {
227  sumDpF += dpf.boundaryField()[patchid][faceid];
228  }
229  sumpf += pf.boundaryField()[patchid][faceid];
230  n = n + 1;
231  }
232  }
233  sumpf /= n;
234  ScQGD_[icell] = cSc1_*(mag(sumDpF) / sumpf);
235  }
236 
237  if (minSc_ >= 0)
238  {
239  this->ScQGD_ = max(this->ScQGD_, minSc_);
240  }
241  if (maxSc_ >= 0)
242  {
243  this->ScQGD_ = min(this->ScQGD_, maxSc_);
244  }
245 
246  if (constScCellSetPtr_.valid())
247  {
248  const cellSet& constScCells = constScCellSetPtr_();
249 
250  forAllConstIter(cellSet, constScCells, iter)
251  {
252  ScQGD_[iter.key()] = constSc_;
253  }
254  }
255 
256 // if (mesh_.time().outputTime())
257 // {
258 // sumdp.write();
259 // avrp.write();
260 // }
261 
262 // this->ScQGD_.primitiveFieldRef() =
263 // max(this->ScQGD_.primitiveField(), cqSc_);
264 
265 // fvc::smooth(this->ScQGD_, smoothCoeff_);
266 
267  Info<< "max/min ScQGD: "
268  << max(this->ScQGD_).value() << "/"
269  << min(this->ScQGD_).value() << endl;
270 
271  if (runTime_.outputTime())
272  {
273  this->ScQGD_.write();
274  this->hQGD_.write();
275  }
276 
277  forAll(p.primitiveField(), celli)
278  {
279  muQGD_.primitiveFieldRef()[celli] =
280  p.primitiveField()[celli] *
281  ScQGD_.primitiveField()[celli] *
282  tauQGD_.primitiveField()[celli];
283 
284  alphauQGD_.primitiveFieldRef()[celli] = muQGD_.primitiveField()[celli] /
285  PrQGD_.primitiveField()[celli];
286  }
287 
288  forAll(p.boundaryField(), patchi)
289  {
290  forAll(p.boundaryField()[patchi], facei)
291  {
292  muQGD_.boundaryFieldRef()[patchi][facei] =
293  p.boundaryField()[patchi][facei] *
294  ScQGD_.boundaryField()[patchi][facei] *
295  tauQGD_.boundaryField()[patchi][facei];
296 
297  alphauQGD_.boundaryFieldRef()[patchi][facei] =
298  muQGD_.boundaryFieldRef()[patchi][facei] /
299  PrQGD_.boundaryField()[patchi][facei];
300  }
301  }
302 }
303 
304 //
305 //END-OF-FILE
306 //
volScalarField PrQGD_
Definition: QGDCoeffs.H:115
pf
Definition: updateFields.H:21
varScModel7(const IOobject &io, const fvMesh &mesh, const dictionary &dict)
Definition: varScModel7.C:50
virtual const volScalarField & c() const =0
virtual tmp< volScalarField > rho() const =0
scalar qgdAspectRatioThreshold_
Definition: varScModel7.H:69
forAll(Y, i)
Definition: QGDYEqn.H:36
autoPtr< cellSet > constScCellSetPtr_
Definition: varScModel7.H:75
Abstract base class for classes implementing thermophysical properties of gases and fluids governed b...
Definition: QGDThermo.H:53
const fvMesh & mesh_
Definition: QGDCoeffs.H:65
volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho())
void correct(const QGDThermo &)
Definition: varScModel7.C:160
virtual const volScalarField & p() const =0
addToRunTimeSelectionTable(QGDCoeffs, constScPrModel1, dictionary)
defineTypeNameAndDebug(constScPrModel1, 0)
volScalarField & p
Definition: createFields.H:52
volScalarField ScQGD_
Definition: QGDCoeffs.H:120