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) 2016-2018 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 #include "surfaceFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace viscosityModels
35 {
36  defineTypeNameAndDebug(strainRateFunction, 0);
37 
39  (
40  viscosityModel,
41  strainRateFunction,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const word& name,
53  const dictionary& viscosityProperties,
54  const volVectorField& U,
55  const surfaceScalarField& phi
56 )
57 :
58  viscosityModel(name, viscosityProperties, U, phi),
59  strainRateFunctionCoeffs_
60  (
61  viscosityProperties.optionalSubDict(typeName + "Coeffs")
62  ),
63  strainRateFunction_
64  (
65  Function1<scalar>::New("function", strainRateFunctionCoeffs_)
66  ),
67  nu_
68  (
69  IOobject
70  (
71  name,
72  U_.time().timeName(),
73  U_.db(),
76  ),
77  U_.mesh(),
79  )
80 {
81  correct();
82 }
83 
84 
85 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
86 
89 {
90  return nu_;
91 }
92 
93 
96 {
97  return nu_.boundaryField()[patchi];
98 }
99 
100 
102 {
103  tmp<volScalarField> tsigma = strainRate();
104  const volScalarField& sigma = tsigma();
105 
106  nu_.primitiveFieldRef() = strainRateFunction_->value(sigma());
107 
108  volScalarField::Boundary& nuBf = nu_.boundaryFieldRef();
109  const volScalarField::Boundary& sigmaBf = sigma.boundaryField();
110 
111  forAll(nuBf, patchi)
112  {
113  nuBf[patchi] = strainRateFunction_->value(sigmaBf[patchi]);
114  }
115 }
116 
117 
119 (
120  const dictionary& viscosityProperties
121 )
122 {
123  viscosityModel::read(viscosityProperties);
124 
125  strainRateFunctionCoeffs_ = viscosityProperties.optionalSubDict
126  (
127  typeName + "Coeffs"
128  );
129 
130  strainRateFunction_.clear();
131  strainRateFunction_ = Function1<scalar>::New
132  (
133  "function",
134  strainRateFunctionCoeffs_
135  );
136 
137  return true;
138 }
139 
140 
141 // ************************************************************************* //
Foam::surfaceFields.
virtual bool read(const dictionary &viscosityProperties)
Read transportProperties dictionary.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const dimensionSet dimViscosity
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
An abstract base class for incompressible viscosityModels.
Macros for easy insertion into run-time selection tables.
addToRunTimeSelectionTable(viscosityModel, BirdCarreau, dictionary)
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:759
virtual bool read(const dictionary &viscosityProperties)=0
Read transportProperties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-6/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
defineTypeNameAndDebug(BirdCarreau, 0)
strainRateFunction(const word &name, const dictionary &viscosityProperties, const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:53
void clear()
Clear the dictionary.
Definition: dictionary.C:1135
virtual void correct()
Correct the laminar viscosity.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static autoPtr< Function1< Type > > New(const word &entryName, const dictionary &dict)
Selector.
Definition: Function1New.C:32
Namespace for OpenFOAM.