kaqRWallFunctionFvPatchScalarField.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) 2014-2022 PCOpt/NTUA
9 Copyright (C) 2014-2022 FOSS GP
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
30#include "fvPatchFieldMapper.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
37(
38 const fvPatch& p,
40)
41:
42 kqRWallFunctionFvPatchField<scalar>(p, iF),
44{}
45
46
48(
49 const fvPatch& p,
51 const dictionary& dict
52)
53:
55 adjointScalarBoundaryCondition(p, iF, dict.get<word>("solverName"))
56{}
57
58
60(
62 const fvPatch& p,
64 const fvPatchFieldMapper& mapper
65)
66:
67 kqRWallFunctionFvPatchField<scalar>(ptf, p, iF, mapper),
68 adjointScalarBoundaryCondition(p, iF, ptf.adjointSolverName_)
69{}
70
71
73(
75)
76:
77 kqRWallFunctionFvPatchField<scalar>(tkqrwfpf),
79{}
80
81
83(
86)
87:
88 kqRWallFunctionFvPatchField<scalar>(tkqrwfpf, iF),
90{}
91
92// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93
95(
96 fvMatrix<scalar>& matrix
97)
98{
99 scalarField& source = matrix.source();
100 tmp<fvPatchScalarField> tnutWall(boundaryContrPtr_->turbulentDiffusivity());
101 fvPatchScalarField& nutWall = tnutWall.constCast();
102 if (isA<nutkWallFunctionFvPatchScalarField>(nutWall))
103 {
104 const label patchi(patch().index());
105 const scalarField& magSf = patch().magSf();
106 const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
107 (
109 (
111 internalField().group()
112 )
113 );
114 const scalarField& y = turbModel.y()[patchi];
115 const tmp<scalarField> tnuw = turbModel.nu(patchi);
116 const scalarField& nuw = tnuw();
118 refCast<nutWallFunctionFvPatchScalarField>(nutWall);
119 const wallFunctionCoefficients& wallCoeffs = nutWF.wallCoeffs();
120 const scalar Cmu = wallCoeffs.Cmu();
121 const scalar kappa = wallCoeffs.kappa();
122 const scalar E = wallCoeffs.E();
123 const scalar yPlusLam = wallCoeffs.yPlusLam();
124
125 const scalar Cmu25 = pow025(Cmu);
126
127 const labelList& faceCells = patch().faceCells();
128 const fvPatchVectorField& Uw = boundaryContrPtr_->Ub();
129 const scalarField magGradUw(mag(Uw.snGrad()));
130
131 tmp<scalarField> tdJdnut = boundaryContrPtr_->dJdnut();
132 const scalarField& dJdnut = tdJdnut();
133
134 const tmp<volScalarField> tk = turbModel.k();
135 const volScalarField& k = tk();
136 forAll(dJdnut, facei)
137 {
138 const label celli = faceCells[facei];
139
140 const scalar sqrtkCell(sqrt(k[celli]));
141 const scalar yPlus = Cmu25*y[facei]*sqrtkCell/nuw[facei];
142 if (yPlusLam < yPlus)
143 {
144 const scalar logEyPlus = log(E*yPlus);
145 const scalar dnut_dyPlus =
146 nuw[facei]*kappa*(logEyPlus - 1)/sqr(logEyPlus);
147 const scalar dyPlus_dk =
148 Cmu25*y[facei]/(2*nuw[facei]*sqrtkCell);
149 const scalar dnut_dk = dnut_dyPlus*dyPlus_dk;
150 source[celli] -= dJdnut[facei]*dnut_dk*magSf[facei];
151 }
152 }
153 }
154}
155
156
158{
160 os.writeEntry("solverName", adjointSolverName_);
161}
162
163
164// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165
166namespace Foam
167{
169 (
172 );
173}
174
175
176// ************************************************************************* //
scalar y
label k
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
Base class for solution control classes.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
virtual bool write()
Write the output fields.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Field< Type > & source() noexcept
Definition: fvMatrix.H:458
A FieldMapper for finite-volume patch fields.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:229
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Add source terms to the rhs of the first cell centre.
This boundary condition provides a simple wrapper around the zero-gradient condition,...
The class nutWallFunction is an abstract base class that hosts calculation methods and common functi...
const wallFunctionCoefficients & wallCoeffs() const noexcept
Return wallFunctionCoefficients.
A class for managing temporary objects.
Definition: tmp.H:65
T & constCast() const
Definition: tmpI.H:248
Abstract base class for turbulence models (RAS, LES and laminar).
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.
const nearWallDist & y() const
Return the near wall distances.
Class to host the wall-function coefficients being used in the wall function boundary conditions.
scalar kappa() const noexcept
Return the object: kappa.
scalar E() const noexcept
Return the object: E.
scalar yPlusLam() const noexcept
Return the object: yPlusLam.
scalar Cmu() const noexcept
Return the object: Cmu.
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
scalar yPlus
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar pow025(const dimensionedScalar &ds)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333