All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
generalisedNewtonian.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 "generalisedNewtonian.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "fvcGrad.H"
30 #include "fvcDiv.H"
31 #include "fvmLaplacian.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace laminarModels
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class BasicMomentumTransportModel>
44 (
45  const alphaField& alpha,
46  const rhoField& rho,
47  const volVectorField& U,
48  const surfaceScalarField& alphaRhoPhi,
49  const surfaceScalarField& phi,
50  const viscosity& viscosity
51 )
52 :
54  (
55  typeName,
56  alpha,
57  rho,
58  U,
59  alphaRhoPhi,
60  phi,
61  viscosity
62  ),
63 
64  viscosityModel_
65  (
67  (
68  this->coeffDict_,
69  viscosity,
70  U
71  )
72  )
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
78 template<class BasicMomentumTransportModel>
80 {
81  viscosityModel_->read(this->coeffDict_);
82 
83  return true;
84 }
85 
86 
87 template<class BasicMomentumTransportModel>
90 {
91  return volScalarField::New
92  (
93  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
94  viscosityModel_->nu()
95  );
96 }
97 
98 
99 template<class BasicMomentumTransportModel>
102 (
103  const label patchi
104 ) const
105 {
106  return viscosityModel_->nu(patchi);
107 }
108 
109 
110 template<class BasicMomentumTransportModel>
112 {
113  viscosityModel_->correct();
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
119 
120 } // End namespace laminarModels
121 } // End namespace Foam
122 
123 // ************************************************************************* //
Foam::surfaceFields.
U
Definition: pEqn.H:72
virtual bool read()
Read momentumTransport dictionary.
static autoPtr< generalisedNewtonianViscosityModel > New(const dictionary &viscosityProperties, const viscosity &viscosity, const volVectorField &U)
Return a reference to the selected viscosity model.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
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,.
Calculate the matrix for the laplacian of the field.
Linear viscous stress turbulence model base class.
BasicMomentumTransportModel::rhoField rhoField
virtual void correct()
Correct the generalisedNewtonian viscosity.
Calculate the gradient of the given field.
static word groupName(Name name, const word &group)
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))
BasicMomentumTransportModel::alphaField alphaField
phi
Definition: correctPhi.H:3
Calculate the divergence of the given field.
generalisedNewtonian(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
label patchi
virtual void correct()
Correct the laminar viscosity.
Definition: laminarModel.C:260
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.