FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

phaseModel.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "phaseModel.H"
00027 #include <finiteVolume/fixedValueFvPatchFields.H>
00028 #include <finiteVolume/surfaceInterpolate.H>
00029 
00030 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00031 
00032 Foam::phaseModel::phaseModel
00033 (
00034     const fvMesh& mesh,
00035     const dictionary& transportProperties,
00036     const word& phaseName
00037 )
00038 :
00039     dict_
00040     (
00041         transportProperties.subDict("phase" + phaseName)
00042     ),
00043     name_(phaseName),
00044     d_
00045     (
00046         dict_.lookup("d")
00047     ),
00048     nu_
00049     (
00050         dict_.lookup("nu")
00051     ),
00052     rho_
00053     (
00054         dict_.lookup("rho")
00055     ),
00056     U_
00057     (
00058         IOobject
00059         (
00060             "U" + phaseName,
00061             mesh.time().timeName(),
00062             mesh,
00063             IOobject::MUST_READ,
00064             IOobject::AUTO_WRITE
00065         ),
00066         mesh
00067     )
00068 {
00069     const word phiName = "phi" + phaseName;
00070 
00071     IOobject phiHeader
00072     (
00073         phiName,
00074         mesh.time().timeName(),
00075         mesh,
00076         IOobject::NO_READ
00077     );
00078 
00079     if (phiHeader.headerOk())
00080     {
00081         Info<< "Reading face flux field " << phiName << endl;
00082 
00083         phiPtr_.reset
00084         (
00085             new surfaceScalarField
00086             (
00087                 IOobject
00088                 (
00089                     phiName,
00090                     mesh.time().timeName(),
00091                     mesh,
00092                     IOobject::MUST_READ,
00093                     IOobject::AUTO_WRITE
00094                 ),
00095                 mesh
00096             )
00097         );
00098     }
00099     else
00100     {
00101         Info<< "Calculating face flux field " << phiName << endl;
00102 
00103         wordList phiTypes
00104         (
00105             U_.boundaryField().size(),
00106             calculatedFvPatchScalarField::typeName
00107         );
00108 
00109         for (label i=0; i<U_.boundaryField().size(); i++)
00110         {
00111             if (isA<fixedValueFvPatchVectorField>(U_.boundaryField()[i]))
00112             {
00113                 phiTypes[i] = fixedValueFvPatchScalarField::typeName;
00114             }
00115         }
00116 
00117         phiPtr_.reset
00118         (
00119             new surfaceScalarField
00120             (
00121                 IOobject
00122                 (
00123                     phiName,
00124                     mesh.time().timeName(),
00125                     mesh,
00126                     IOobject::NO_READ,
00127                     IOobject::AUTO_WRITE
00128                 ),
00129                 fvc::interpolate(U_) & mesh.Sf(),
00130                 phiTypes
00131             )
00132         );
00133     }
00134 }
00135 
00136 
00137 Foam::autoPtr<Foam::phaseModel> Foam::phaseModel::New
00138 (
00139     const fvMesh& mesh,
00140     const dictionary& transportProperties,
00141     const word& phaseName
00142 )
00143 {
00144     return autoPtr<phaseModel>
00145     (
00146         new phaseModel(mesh, transportProperties, phaseName)
00147     );
00148 }
00149 
00150 
00151 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00152 
00153 Foam::phaseModel::~phaseModel()
00154 {}
00155 
00156 
00157 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines