Loading...
Searching...
No Matches
strainRateFunction.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) 2018-2020 OpenFOAM Foundation
9 Copyright (C) 2020 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
29#include "strainRateFunction.H"
30#include "volFields.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35namespace Foam
37namespace laminarModels
38{
40{
42
44 (
48 );
50}
51}
52
53
54// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55
58(
59 const dictionary& viscosityProperties
60)
61:
62 generalizedNewtonianViscosityModel(viscosityProperties),
63 strainRateFunctionCoeffs_
64 (
65 viscosityProperties.optionalSubDict(typeName + "Coeffs")
66 ),
67 strainRateFunction_
68 (
69 Function1<scalar>::New("function", strainRateFunctionCoeffs_)
70 )
71{}
72
73
74// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
75
78(
79 const dictionary& viscosityProperties
80)
81{
83
84 strainRateFunctionCoeffs_ =
85 viscosityProperties.optionalSubDict(typeName + "Coeffs");
86
87 strainRateFunction_.reset();
88 strainRateFunction_ = Function1<scalar>::New
89 (
90 "function",
91 strainRateFunctionCoeffs_
92 );
93
94 return true;
95}
96
97
100nu
101(
102 const volScalarField& nu0,
103 const volScalarField& strainRate
104) const
105{
106 auto tnu = volScalarField::New
107 (
109 nu0.mesh(),
111 );
112
113 // Set the database - could not be done in constructor
114 const_cast<Function1<scalar>&>(strainRateFunction_()).resetDb
115 (
116 nu0.mesh().thisDb()
117 );
118
119 tnu.ref().primitiveFieldRef() = strainRateFunction_->value(strainRate);
120
121 volScalarField::Boundary& nuBf = tnu.ref().boundaryFieldRef();
122 const volScalarField::Boundary& sigmaBf = strainRate.boundaryField();
123
124 forAll(nuBf, patchi)
125 {
126 nuBf[patchi] = strainRateFunction_->value(sigmaBf[patchi]);
127 }
128
129 return tnu;
130}
131
132
133// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const Mesh & mesh() const noexcept
Return const reference to mesh.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
static word group(const word &name)
Return group (extension part of name).
Definition IOobject.C:300
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition dictionary.C:560
An abstract base class for generalized Newtonian viscosity models.
generalizedNewtonianViscosityModel(const generalizedNewtonianViscosityModel &)=delete
No copy construct.
static autoPtr< generalizedNewtonianViscosityModel > New(const dictionary &dict)
Select a viscosity model.
const dictionary & viscosityProperties() const
Return the phase transport properties dictionary.
virtual bool read(const dictionary &viscosityProperties)=0
Read transportProperties dictionary.
Run-time selected strain-rate function generalized Newtonian viscosity model.
virtual bool read(const dictionary &viscosityProperties)
Read transportProperties dictionary.
virtual tmp< volScalarField > nu(const volScalarField &nu0, const volScalarField &strainRate) const
Return the laminar viscosity.
strainRateFunction(const dictionary &viscosityProperties)
Construct from dictionary (components).
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
A namespace for various generalized Newtonian viscosity model implementations.
Definition BirdCarreau.C:33
Namespace for OpenFOAM.
const dimensionSet dimViscosity
GeometricField< scalar, fvPatchField, volMesh > volScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299