phaseCompressibleMomentumTransportModel.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) 2013-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 
27 #include "surfaceInterpolate.H"
28 
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  defineTypeNameAndDebug(phaseCompressibleMomentumTransportModel, 0);
35 }
36 
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
42 (
43  const word& type,
44  const volScalarField& alpha,
45  const volScalarField& rho,
46  const volVectorField& U,
47  const surfaceScalarField& alphaRhoPhi,
48  const surfaceScalarField& phi,
49  const viscosity& viscosity
50 )
51 :
53  (
54  type,
56  rho,
57  U,
58  alphaRhoPhi,
59  phi,
60  viscosity
61  ),
62  alpha_(alpha)
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
67 
70 (
71  const volScalarField& alpha,
72  const volScalarField& rho,
73  const volVectorField& U,
74  const surfaceScalarField& alphaRhoPhi,
75  const surfaceScalarField& phi,
76  const viscosity& viscosity
77 )
78 {
80  <
82  >
83  (
84  alpha,
85  rho,
86  U,
87  alphaRhoPhi,
88  phi,
89  viscosity
90  );
91 }
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
98 {
99  return volScalarField::New
100  (
101  IOobject::groupName("pPrime", this->alphaRhoPhi_.group()),
102  this->mesh_,
104  );
105 }
106 
107 
110 {
112  (
113  IOobject::groupName("pPrimef", this->alphaRhoPhi_.group()),
114  this->mesh_,
116  );
117 }
118 
119 
120 // ************************************************************************* //
Templated abstract base class for multiphase compressible turbulence models.
U
Definition: pEqn.H:72
const dimensionSet dimPressure
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,.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure&#39;.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A class for handling words, derived from string.
Definition: word.H:59
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))
phaseCompressibleMomentumTransportModel(const word &type, const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct.
phi
Definition: correctPhi.H:3
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:53
static autoPtr< MomentumTransportModel > New(const typename MomentumTransportModel::alphaField &alpha, const typename MomentumTransportModel::rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
static autoPtr< phaseCompressibleMomentumTransportModel > New(const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Return a reference to the selected turbulence model.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
Base class for single-phase compressible turbulence models.
Namespace for OpenFOAM.