strainRateFunction.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) 2018-2021 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 
26 #include "strainRateFunction.H"
27 #include "volFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace laminarModels
35 {
36 namespace generalisedNewtonianViscosityModels
37 {
38  defineTypeNameAndDebug(strainRateFunction, 0);
39 
41  (
42  generalisedNewtonianViscosityModel,
43  strainRateFunction,
44  dictionary
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
55 (
56  const dictionary& viscosityProperties
57 )
58 :
59  generalisedNewtonianViscosityModel(viscosityProperties),
60  strainRateFunction_
61  (
63  (
64  "function",
65  viscosityProperties.optionalSubDict(typeName + "Coeffs")
66  )
67  )
68 {}
69 
70 
71 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
72 
75 (
76  const dictionary& viscosityProperties
77 )
78 {
79  generalisedNewtonianViscosityModel::read(viscosityProperties);
80 
81  strainRateFunction_.clear();
82  strainRateFunction_ = Function1<scalar>::New
83  (
84  "function",
85  viscosityProperties.optionalSubDict
86  (
87  typeName + "Coeffs"
88  )
89  );
90 
91  return true;
92 }
93 
94 
97 nu
98 (
99  const volScalarField& nu0,
100  const volScalarField& strainRate
101 ) const
102 {
104  (
106  (
107  IOobject::groupName(type() + ":nu", nu0.group()),
108  nu0.mesh(),
110  )
111  );
112 
113  tnu.ref().primitiveFieldRef() = strainRateFunction_->value(strainRate);
114 
115  volScalarField::Boundary& nuBf = tnu.ref().boundaryFieldRef();
116  const volScalarField::Boundary& sigmaBf = strainRate.boundaryField();
117 
118  forAll(nuBf, patchi)
119  {
120  nuBf[patchi] = strainRateFunction_->value(sigmaBf[patchi]);
121  }
122 
123  return tnu;
124 }
125 
126 
127 // ************************************************************************* //
const dimensionSet dimViscosity
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
Run-time selectable general function of one variable.
Definition: Function1.H:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
strainRateFunction(const dictionary &viscosityProperties)
Construct from components.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual tmp< volScalarField > nu(const volScalarField &nu0, const volScalarField &strainRate) const
Return the laminar viscosity.
Macros for easy insertion into run-time selection tables.
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1040
static word groupName(Name name, const word &group)
addToRunTimeSelectionTable(generalisedNewtonianViscosityModel, BirdCarreau, dictionary)
virtual bool read(const dictionary &viscosityProperties)
Read transportProperties dictionary.
const Mesh & mesh() const
Return mesh.
An abstract base class for generalised Newtonian viscosity models.
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
virtual bool read(const dictionary &viscosityProperties)=0
Read transportProperties dictionary.
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32