lambdaThixotropic.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) 2020-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 "lambdaThixotropic.H"
27 #include "fvModels.H"
28 #include "fvConstraints.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace laminarModels
36 {
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 template<class BasicMomentumTransportModel>
42 (
43  const alphaField& alpha,
44  const rhoField& rho,
45  const volVectorField& U,
46  const surfaceScalarField& alphaRhoPhi,
47  const surfaceScalarField& phi,
48  const viscosity& viscosity
49 )
50 :
52  (
53  typeName,
54  alpha,
55  rho,
56  U,
57  alphaRhoPhi,
58  phi,
59  viscosity
60  ),
61 
62  a_("a", dimless/dimTime, this->coeffDict_),
63  b_("b", dimless, this->coeffDict_),
64  d_("d", dimless, this->coeffDict_),
65  c_("c", pow(dimTime, d_.value() - scalar(1)), this->coeffDict_),
66  nu0_("nu0", dimViscosity, this->coeffDict_),
67  nuInf_("nuInf", dimViscosity, this->coeffDict_),
68  K_(1 - sqrt(nuInf_/nu0_)),
69 
70  lambda_
71  (
72  IOobject
73  (
75  (
76  IOobject::modelName("lambda", typeName),
77  alphaRhoPhi.group()
78  ),
79  this->runTime_.timeName(),
80  this->mesh_,
83  ),
84  this->mesh_
85  ),
86 
87  nu_
88  (
89  IOobject
90  (
92  (
93  IOobject::modelName("nu", typeName),
94  alphaRhoPhi.group()
95  ),
96  this->runTime_.timeName(),
97  this->mesh_,
100  ),
101  calcNu()
102  )
103 {}
104 
105 
106 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
107 
108 template<class BasicMomentumTransportModel>
111 {
112  return nuInf_/(sqr(1 - K_*lambda_) + rootVSmall);
113 }
114 
115 
116 template<class BasicMomentumTransportModel>
119 {
120  return sqrt(2.0)*mag(symm(fvc::grad(this->U())()()));
121 }
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
126 template<class BasicMomentumTransportModel>
128 {
130  {
131  a_.read(this->coeffDict());
132  b_.read(this->coeffDict());
133  d_.read(this->coeffDict());
134 
135  c_ = dimensionedScalar
136  (
137  "c",
138  pow(dimTime, d_.value() - scalar(1)),
139  this->coeffDict_
140  );
141 
142  nu0_.read(this->coeffDict());
143  nuInf_.read(this->coeffDict());
144 
145  K_ = (1 - sqrt(nuInf_/nu0_));
146 
147  return true;
148  }
149  else
150  {
151  return false;
152  }
153 }
154 
155 
156 template<class BasicMomentumTransportModel>
159 {
160  return volScalarField::New
161  (
162  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
163  nu_
164  );
165 }
166 
167 
168 template<class BasicMomentumTransportModel>
171 (
172  const label patchi
173 ) const
174 {
175  return nu_.boundaryField()[patchi];
176 }
177 
178 
179 template<class BasicMomentumTransportModel>
181 {
182  // Local references
183  const surfaceScalarField& phi = this->phi_;
184  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
186  (
187  Foam::fvConstraints::New(this->mesh_)
188  );
189 
190  tmp<fvScalarMatrix> lambdaEqn
191  (
192  fvm::ddt(lambda_) + fvm::div(phi, lambda_)
193  - fvm::Sp(fvc::div(phi), lambda_)
194  ==
195  a_*pow(1 - lambda_(), b_)
196  - fvm::Sp(c_*pow(strainRate(), d_), lambda_)
197  + fvModels.source(lambda_)
198  );
199 
200  lambdaEqn.ref().relax();
201  fvConstraints.constrain(lambdaEqn.ref());
202  solve(lambdaEqn);
203  fvConstraints.constrain(lambda_);
204 
205  lambda_.maxMin(scalar(0), scalar(1));
206 
207  nu_ = calcNu();
208 
210 }
211 
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 } // End namespace laminarModels
216 } // End namespace Foam
217 
218 // ************************************************************************* //
const dimensionSet dimViscosity
virtual void correct()
Correct the lambdaThixotropic viscosity.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Templated abstract base class for laminar transport models.
Definition: laminarModel.H:49
U
Definition: pEqn.H:72
BasicMomentumTransportModel::alphaField alphaField
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Linear viscous stress turbulence model base class.
const dimensionSet dimless
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
const dimensionSet dimTime
Foam::fvConstraints & fvConstraints
static word groupName(Name name, const word &group)
Finite volume models.
Definition: fvModels.H:60
Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, simple.dict());mesh.schemes().setFluxRequired(p.name());Info<< "Reading field pa\"<< endl;volScalarField pa(IOobject("pa", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Ua\"<< endl;volVectorField Ua(IOobject("Ua", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);# 65 "/home/ubuntu/OpenFOAM-10/applications/solvers/incompressible/adjointShapeOptimisationFoam/createFields.H" 2label paRefCell=0;scalar paRefValue=0.0;setRefCell(pa, simple.dict(), paRefCell, paRefValue);mesh.schemes().setFluxRequired(pa.name());autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
phi
Definition: correctPhi.H:3
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
Finite volume constraints.
Definition: fvConstraints.H:57
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
virtual void correct()
Correct the laminar viscosity.
Definition: laminarModel.C:260
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read()
Read momentumTransport dictionary.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dimensioned< scalar > mag(const dimensioned< Type > &)
Thixotropic viscosity momentum transport model based on the evolution of the structural parameter : ...
static word modelName(Name name, const word &model)
Return the name of the object within the given model.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
BasicMomentumTransportModel::rhoField rhoField
lambdaThixotropic(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
Namespace for OpenFOAM.