phase.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) 2011-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 "phase.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 (
32  const word& phaseName,
33  const dictionary& phaseDict,
34  const volVectorField& U,
35  const surfaceScalarField& phi
36 )
37 :
39  (
40  IOobject
41  (
42  IOobject::groupName("alpha", phaseName),
43  U.mesh().time().timeName(),
44  U.mesh(),
45  IOobject::MUST_READ,
46  IOobject::AUTO_WRITE
47  ),
48  U.mesh()
49  ),
50  name_(phaseName),
51  phaseDict_(phaseDict),
52  nuModel_
53  (
54  viscosityModel::New
55  (
56  IOobject::groupName("nu", phaseName),
57  phaseDict_,
58  U,
59  phi
60  )
61  ),
62  rho_("rho", dimDensity, phaseDict_)
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
69 {
71  return autoPtr<phase>(nullptr);
72 }
73 
74 
76 {
77  nuModel_->correct();
78 }
79 
80 
81 bool Foam::phase::read(const dictionary& phaseDict)
82 {
83  phaseDict_ = phaseDict;
84 
85  if (nuModel_->read(phaseDict_))
86  {
87  phaseDict_.lookup("rho") >> rho_;
88 
89  return true;
90  }
91  else
92  {
93  return false;
94  }
95 }
96 
97 
98 // ************************************************************************* //
virtual bool read()
Read object.
void correct()
Correct the phase properties.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
phi
Definition: pEqn.H:104
dictionary()
Construct top-level dictionary null.
Definition: dictionary.C:464
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
autoPtr< phase > clone() const
Return clone.
phase(const word &name, const dictionary &phaseDict, const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimDensity
U
Definition: pEqn.H:72
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366