epsilonmWallFunctionFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2022-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "fvMatrix.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 Foam::scalar Foam::epsilonmWallFunctionFvPatchScalarField::tolerance_ = 1e-1;
33 
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42  const dictionary& dict
43 )
44 :
45  fixedValueFvPatchField<scalar>(p, iF, dict)
46 {
47  // Apply zero-gradient condition on start-up
49 }
50 
51 
54 (
56  const fvPatch& p,
58  const fvPatchFieldMapper& mapper
59 )
60 :
61  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper)
62 {}
63 
64 
67 (
70 )
71 :
72  fixedValueFvPatchField<scalar>(ewfpsf, iF)
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
79 (
80  fvMatrix<scalar>& matrix
81 )
82 {
83  if (manipulatedMatrix())
84  {
85  return;
86  }
87 
88  const DimensionedField<scalar, volMesh>& epsilon = internalField();
89 
90  scalarField weights(patch().magSf()/patch().patch().magFaceAreas());
91  forAll(weights, facei)
92  {
93  scalar& w = weights[facei];
94  w = w <= tolerance_ ? 0 : (w - tolerance_)/(1 - tolerance_);
95  }
96 
97  matrix.setValues
98  (
99  patch().faceCells(),
100  UIndirectList<scalar>(epsilon, patch().faceCells()),
101  weights
102  );
103 
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
109 
110 namespace Foam
111 {
113  (
116  );
117 }
118 
119 
120 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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...
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
This boundary condition provides a turbulence dissipation wall constraint for the Foam::mixtureKEpsil...
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
epsilonmWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void setValues(const labelUList &cells, const ListType< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:502
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:417
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:172
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:224
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
const scalar epsilon
const dimensionedScalar e
Elementary charge.
Namespace for OpenFOAM.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p