Loading...
Searching...
No Matches
omegaWallFunctionFvPatchScalarField.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, 2019 OpenFOAM Foundation
9 Copyright (C) 2017-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
31#include "turbulenceModel.H"
32#include "fvMatrix.H"
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
38
39// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40
42{
43 if (master_ != -1)
44 {
45 return;
46 }
47
48 const auto& omega =
49 static_cast<const volScalarField&>(this->internalField());
50
51 const volScalarField::Boundary& bf = omega.boundaryField();
52
53 label master = -1;
54 forAll(bf, patchi)
55 {
57 {
59
60 if (master == -1)
61 {
62 master = patchi;
63 }
65 opf.master() = master;
66 }
67 }
68}
69
70
72{
73 const auto& omega =
74 static_cast<const volScalarField&>(this->internalField());
75
76 const volScalarField::Boundary& bf = omega.boundaryField();
77
78 const fvMesh& mesh = omega.mesh();
79
80 if (initialised_ && !mesh.changing())
81 {
82 return;
83 }
84
85 volScalarField weights
86 (
88 (
89 "weights",
90 mesh.time().timeName(),
91 mesh,
95 ),
96 mesh,
98 );
99
100 DynamicList<label> omegaPatches(bf.size());
101 forAll(bf, patchi)
102 {
104 {
105 omegaPatches.append(patchi);
106
107 const labelUList& faceCells = bf[patchi].patch().faceCells();
108 for (const auto& celli : faceCells)
109 {
110 ++weights[celli];
111 }
112 }
113 }
114
115 cornerWeights_.setSize(bf.size());
116 for (const auto& patchi : omegaPatches)
117 {
118 const fvPatchScalarField& wf = weights.boundaryField()[patchi];
119 cornerWeights_[patchi] = 1.0/wf.patchInternalField();
120 }
121
122 G_.setSize(internalField().size(), 0.0);
123 omega_.setSize(internalField().size(), 0.0);
124
125 initialised_ = true;
126}
127
128
131(
132 const label patchi
133)
134{
135 const auto& omega =
136 static_cast<const volScalarField&>(this->internalField());
137
138 const volScalarField::Boundary& bf = omega.boundaryField();
139
140 const auto& opf =
142
143 return const_cast<omegaWallFunctionFvPatchScalarField&>(opf);
144}
145
146
148(
149 const turbulenceModel& turbModel,
150 scalarField& G0,
151 scalarField& omega0
152)
153{
154 // accumulate all of the G and omega contributions
155 forAll(cornerWeights_, patchi)
156 {
157 if (!cornerWeights_[patchi].empty())
158 {
159 omegaWallFunctionFvPatchScalarField& opf = omegaPatch(patchi);
160
161 const List<scalar>& w = cornerWeights_[patchi];
162
163 opf.calculate(turbModel, w, opf.patch(), G0, omega0);
164 }
165 }
166
167 // apply zero-gradient condition for omega
168 forAll(cornerWeights_, patchi)
169 {
170 if (!cornerWeights_[patchi].empty())
171 {
172 omegaWallFunctionFvPatchScalarField& opf = omegaPatch(patchi);
174 opf == scalarField(omega0, opf.patch().faceCells());
175 }
176 }
177}
178
179
181(
182 const turbulenceModel& turbModel,
183 const List<scalar>& cornerWeights,
184 const fvPatch& patch,
185 scalarField& G0,
186 scalarField& omega0
187)
188{
189 const label patchi = patch.index();
190
191 const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
192 const scalar kappa = wallCoeffs_.kappa();
193 const scalar yPlusLam = wallCoeffs_.yPlusLam();
194
195 const scalarField& y = turbModel.y()[patchi];
196
197 const labelUList& faceCells = patch.faceCells();
198
199 const tmp<scalarField> tnutw = turbModel.nut(patchi);
200 const scalarField& nutw = tnutw();
201
202 const tmp<scalarField> tnuw = turbModel.nu(patchi);
203 const scalarField& nuw = tnuw();
204
205 const tmp<volScalarField> tk = turbModel.k();
206 const volScalarField& k = tk();
207
208 // Calculate y-plus
209 const auto yPlus = [&](const label facei) -> scalar
210 {
211 return
212 (
213 Cmu25*y[facei]*sqrt(k[faceCells[facei]])/nuw[facei]
214 );
215 };
216
217 // Contribution from the viscous sublayer
218 const auto omegaVis = [&](const label facei) -> scalar
219 {
220 return
221 (
222 6.0*nuw[facei]/(beta1_*sqr(y[facei]))
223 );
224 };
225
226 // Contribution from the inertial sublayer
227 const auto omegaLog = [&](const label facei) -> scalar
228 {
229 return
230 (
231 sqrt(k[faceCells[facei]])
232 / (Cmu25*kappa*y[facei])
233 );
234 };
235
236 switch (blender_)
237 {
238 case blenderType::STEPWISE:
239 {
240 forAll(faceCells, facei)
241 {
242 if (yPlus(facei) > yPlusLam)
243 {
244 omega0[faceCells[facei]] +=
245 cornerWeights[facei]
246 * omegaLog(facei);
247 }
248 else
249 {
250 omega0[faceCells[facei]] +=
251 cornerWeights[facei]
252 * omegaVis(facei);
253 }
254 }
255 break;
256 }
257
258 case blenderType::BINOMIAL:
259 {
260 forAll(faceCells, facei)
261 {
262 omega0[faceCells[facei]] +=
263 cornerWeights[facei]
264 * pow
265 (
266 pow(omegaVis(facei), n_) + pow(omegaLog(facei), n_),
267 scalar(1)/n_
268 );
269 }
270 break;
271 }
272
273 case blenderType::MAX:
274 {
275 forAll(faceCells, facei)
276 {
277 // (PH:Eq. 27)
278 omega0[faceCells[facei]] +=
279 cornerWeights[facei]
280 * max(omegaVis(facei), omegaLog(facei));
281 }
282 break;
283 }
284
285 case blenderType::EXPONENTIAL:
286 {
287 forAll(faceCells, facei)
288 {
289 // (PH:Eq. 31)
290 const scalar yPlusFace = yPlus(facei);
291 const scalar Gamma = 0.01*pow4(yPlusFace)/(1 + 5*yPlusFace);
292 const scalar invGamma = scalar(1)/(Gamma + ROOTVSMALL);
293
294 omega0[faceCells[facei]] +=
295 cornerWeights[facei]
296 * (
297 omegaVis(facei)*exp(-Gamma)
298 + omegaLog(facei)*exp(-invGamma)
299 );
300 }
301 break;
302 }
303
304 case blenderType::TANH:
305 {
306 forAll(faceCells, facei)
307 {
308 // (KAS:Eqs. 33-34)
309 const scalar omegaVisFace = omegaVis(facei);
310 const scalar omegaLogFace = omegaLog(facei);
311 const scalar b1 = omegaVisFace + omegaLogFace;
312 const scalar b2 =
313 pow
314 (
315 pow(omegaVisFace, 1.2) + pow(omegaLogFace, 1.2),
316 1.0/1.2
317 );
318 const scalar phiTanh = tanh(pow4(0.1*yPlus(facei)));
319
320 omega0[faceCells[facei]] +=
321 cornerWeights[facei]
322 *(phiTanh*b1 + (1 - phiTanh)*b2);
323 }
324 break;
325 }
326 }
327
328 const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
329 const scalarField magGradUw(mag(Uw.snGrad()));
330
331 forAll(faceCells, facei)
332 {
333 if (!(blender_ == blenderType::STEPWISE) || yPlus(facei) > yPlusLam)
334 {
335 G0[faceCells[facei]] +=
336 cornerWeights[facei]
337 *(nutw[facei] + nuw[facei])
338 *magGradUw[facei]
339 *Cmu25*sqrt(k[faceCells[facei]])
340 /(kappa*y[facei]);
341 }
342 }
343}
344
345
347(
348 Ostream& os
349) const
350{
352 os.writeEntryIfDifferent<scalar>("beta1", 0.075, beta1_);
353 wallCoeffs_.writeEntries(os);
354}
355
356
357// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
358
360(
361 const fvPatch& p,
363)
364:
365 fixedValueFvPatchField<scalar>(p, iF),
367 initialised_(false),
368 master_(-1),
369 beta1_(0.075),
371 G_(),
372 omega_(),
374{}
375
376
378(
380 const fvPatch& p,
382 const fvPatchFieldMapper& mapper
383)
384:
385 fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
387 initialised_(false),
388 master_(-1),
389 beta1_(ptf.beta1_),
391 G_(),
392 omega_(),
394{}
395
396
398(
399 const fvPatch& p,
401 const dictionary& dict
402)
403:
404 fixedValueFvPatchField<scalar>(p, iF, dict),
405 wallFunctionBlenders(dict, blenderType::BINOMIAL, scalar(2)),
406 initialised_(false),
407 master_(-1),
408 beta1_(dict.getOrDefault<scalar>("beta1", 0.075)),
409 wallCoeffs_(dict),
410 G_(),
411 omega_(),
413{
414 // Apply zero-gradient condition on start-up
415 this->extrapolateInternal();
416}
417
418
420(
422)
423:
424 fixedValueFvPatchField<scalar>(owfpsf),
425 wallFunctionBlenders(owfpsf),
426 initialised_(false),
427 master_(-1),
428 beta1_(owfpsf.beta1_),
445 master_(-1),
446 beta1_(owfpsf.beta1_),
447 wallCoeffs_(owfpsf.wallCoeffs_),
448 G_(),
451{}
452
453
454// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
455
457(
458 bool init
459)
460{
461 if (patch().index() == master_)
462 {
463 if (init)
464 {
465 G_ = 0.0;
466 }
467
468 return G_;
469 }
470
471 return omegaPatch(master_).G();
472}
473
474
476(
477 bool init
478)
479{
480 if (patch().index() == master_)
481 {
482 if (init)
483 {
484 omega_ = 0.0;
485 }
486
487 return omega_;
488 }
489
490 return omegaPatch(master_).omega(init);
491}
492
493
495{
496 if (updated())
497 {
498 return;
499 }
500
501 const auto& turbModel = db().lookupObject<turbulenceModel>
502 (
504 (
506 internalField().group()
507 )
508 );
509
510 setMaster();
511
512 if (patch().index() == master_)
513 {
514 createAveragingWeights();
515 calculateTurbulenceFields(turbModel, G(true), omega(true));
516 }
517
518 const scalarField& G0 = this->G();
519 const scalarField& omega0 = this->omega();
520
521 typedef DimensionedField<scalar, volMesh> FieldType;
522
523 FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
524
525 FieldType& omega = const_cast<FieldType&>(internalField());
526
527 forAll(*this, facei)
528 {
529 const label celli = patch().faceCells()[facei];
530
531 G[celli] = G0[celli];
532 omega[celli] = omega0[celli];
533 }
534
536}
537
538
540(
541 const scalarField& weights
542)
543{
544 if (updated())
545 {
546 return;
547 }
548
549 const auto& turbModel = db().lookupObject<turbulenceModel>
550 (
552 (
554 internalField().group()
555 )
556 );
557
558 setMaster();
559
560 if (patch().index() == master_)
561 {
562 createAveragingWeights();
563 calculateTurbulenceFields(turbModel, G(true), omega(true));
564 }
565
566 const scalarField& G0 = this->G();
567 const scalarField& omega0 = this->omega();
568
569 typedef DimensionedField<scalar, volMesh> FieldType;
570
571 FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
572
573 FieldType& omega = const_cast<FieldType&>(internalField());
574
575 scalarField& omegaf = *this;
576
577 // only set the values if the weights are > tolerance
578 forAll(weights, facei)
579 {
580 const scalar w = weights[facei];
581
582 if (w > tolerance_)
583 {
584 const label celli = patch().faceCells()[facei];
585
586 G[celli] = (1.0 - w)*G[celli] + w*G0[celli];
587 omega[celli] = (1.0 - w)*omega[celli] + w*omega0[celli];
588 omegaf[facei] = omega[celli];
590 }
591
593}
594
595
597(
598 fvMatrix<scalar>& matrix
599)
600{
601 if (manipulatedMatrix())
602 {
603 return;
604 }
605
606 matrix.setValues(patch().faceCells(), patchInternalField());
607
609}
610
611
613(
614 fvMatrix<scalar>& matrix,
615 const Field<scalar>& weights
616)
617{
618 if (manipulatedMatrix())
619 {
620 return;
621 }
622
623 DynamicList<label> constraintCells(weights.size());
624 DynamicList<scalar> constraintValues(weights.size());
625 const labelUList& faceCells = patch().faceCells();
626
627 const DimensionedField<scalar, volMesh>& fld = internalField();
628
629 forAll(weights, facei)
630 {
631 // only set the values if the weights are > tolerance
632 if (tolerance_ < weights[facei])
633 {
634 const label celli = faceCells[facei];
635
636 constraintCells.append(celli);
637 constraintValues.append(fld[celli]);
638 }
639 }
640
641 if (debug)
642 {
643 Pout<< "Patch: " << patch().name()
644 << ": number of constrained cells = " << constraintCells.size()
645 << " out of " << patch().size()
646 << endl;
647 }
649 matrix.setValues(constraintCells, constraintValues);
650
652}
653
654
656(
657 Ostream& os
658) const
659{
661 writeLocalEntries(os);
664
665
666// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
667
668namespace Foam
669{
671 (
673 omegaWallFunctionFvPatchScalarField
674 );
675}
676
677
678// ************************************************************************* //
scalar y
label k
Macros for easy insertion into run-time selection tables.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
bool empty() const noexcept
Definition UList.H:701
void size(const label n)
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
fixedValueFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
void setValues(const labelUList &cellLabels, const Type &value)
Set solution in given cells to the specified value and eliminate the corresponding equations from the...
Definition fvMatrix.C:971
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool manipulatedMatrix() const noexcept
True if the matrix has already been manipulated.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
const DimensionedField< scalar, volMesh > & internalField() const noexcept
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual const labelUList & faceCells() const
Return faceCells.
Definition fvPatch.C:107
This boundary condition provides a wall function for the specific dissipation rate (i....
scalarField G_
Local copy of turbulence G field.
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
Calculate the omega and G.
List< List< scalar > > cornerWeights_
List of averaging corner weights.
omegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
static scalar tolerance_
Tolerance used in weighted calculations.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by multiple wall function faces.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
virtual label & master()
Return non-const access to the master patch ID.
scalarField & G(bool init=false)
Return non-const access to the master's G field.
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
scalarField & omega(bool init=false)
Return non-const access to the master's omega field.
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0)
Main driver to calculate the turbulence fields.
scalarField omega_
Local copy of turbulence omega field.
virtual void setMaster()
Set the master patch - master is responsible for updating all wall function patches.
virtual omegaWallFunctionFvPatchScalarField & omegaPatch(const label patchi)
Helper function to return non-const access to an omega patch.
bool changing() const noexcept
Is mesh changing (topology changing and/or moving).
Definition polyMesh.H:768
A class for managing temporary objects.
Definition tmp.H:75
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
const nearWallDist & y() const
Return the near wall distances.
The class wallFunctionBlenders is a base class that hosts common entries for various derived wall-fun...
void writeEntries(Ostream &) const
Write wall-function blending data as dictionary entries.
wallFunctionBlenders()
Default construct with default coefficients.
blenderType
Options for the blending treatment of viscous and inertial sublayers.
@ BINOMIAL
"Binomial blending (smooth)"
enum blenderType blender_
Blending treatment.
scalar n_
Binomial blending exponent being used when blenderType is blenderType::BINOMIAL.
volScalarField & p
dynamicFvMesh & mesh
scalar yPlus
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
constexpr const char *const group
Group name for atomic constants.
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar G
Newtonian constant of gravitation.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField
dimensionedScalar pow025(const dimensionedScalar &ds)
static bool initialised_(false)
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299