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-2024 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"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace laminarModels
34 {
35 namespace generalisedNewtonianViscosityModels
36 {
38 
40  (
44  );
45 }
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
54 (
55  const dictionary& viscosityProperties,
57  const volVectorField& U
58 )
59 :
60  strainRateViscosityModel(viscosityProperties, viscosity, U),
61  strainRateFunction_
62  (
63  Function1<scalar>::New
64  (
65  "function",
68  viscosityProperties.optionalSubDict(typeName + "Coeffs")
69  )
70  )
71 {
72  correct();
73 }
74 
75 
76 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
77 
80 (
81  const dictionary& viscosityProperties
82 )
83 {
84  strainRateViscosityModel::read(viscosityProperties);
85 
86  strainRateFunction_.clear();
87  strainRateFunction_ = Function1<scalar>::New
88  (
89  "function",
92  viscosityProperties.optionalSubDict
93  (
94  typeName + "Coeffs"
95  )
96  );
97 
98  return true;
99 }
100 
101 
104 nu
105 (
106  const volScalarField& nu0,
107  const volScalarField& strainRate
108 ) const
109 {
111  (
113  (
114  IOobject::groupName(typedName("nu"), nu0.group()),
115  nu0.mesh(),
117  )
118  );
119 
120  tnu.ref().primitiveFieldRef() = strainRateFunction_->value(strainRate);
121 
122  volScalarField::Boundary& nuBf = tnu.ref().boundaryFieldRef();
123  const volScalarField::Boundary& sigmaBf = strainRate.boundaryField();
124 
125  forAll(nuBf, patchi)
126  {
127  nuBf[patchi] = strainRateFunction_->value(sigmaBf[patchi]);
128  }
129 
130  return tnu;
131 }
132 
133 
134 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
const Mesh & mesh() const
Return mesh.
Run-time selectable general function of one variable.
Definition: Function1.H:125
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
static word groupName(Name name, const word &group)
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:921
An abstract base class for generalised Newtonian viscosity models.
virtual bool read(const dictionary &viscosityProperties)=0
Read transportProperties dictionary.
Run-time selected strain-rate function generalised Newtonian viscosity model.
virtual bool read(const dictionary &viscosityProperties)
Read transportProperties dictionary.
strainRateFunction(const dictionary &viscosityProperties, const Foam::viscosity &viscosity, const volVectorField &U)
Construct from components.
An abstract base class for strain-rate dependent generalised Newtonian viscosity models.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
label patchi
U
Definition: pEqn.H:72
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
addToRunTimeSelectionTable(generalisedNewtonianViscosityModel, Newtonian, dictionary)
Namespace for OpenFOAM.
const dimensionSet dimKinematicViscosity
const dimensionSet dimless
const dimensionSet dimTime
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176