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

fvcMeshPhi.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 "fvcMeshPhi.H"
00027 #include <finiteVolume/fvMesh.H>
00028 #include <finiteVolume/ddtScheme.H>
00029 
00030 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00031 
00032 Foam::tmp<Foam::surfaceScalarField> Foam::fvc::meshPhi
00033 (
00034     const volVectorField& vf
00035 )
00036 {
00037     return fv::ddtScheme<vector>::New
00038     (
00039         vf.mesh(),
00040         vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
00041     )().meshPhi(vf);
00042 }
00043 
00044 
00045 Foam::tmp<Foam::surfaceScalarField> Foam::fvc::meshPhi
00046 (
00047     const dimensionedScalar& rho,
00048     const volVectorField& vf
00049 )
00050 {
00051     return fv::ddtScheme<vector>::New
00052     (
00053         vf.mesh(),
00054         vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
00055     )().meshPhi(vf);
00056 }
00057 
00058 
00059 Foam::tmp<Foam::surfaceScalarField> Foam::fvc::meshPhi
00060 (
00061     const volScalarField& rho,
00062     const volVectorField& vf
00063 )
00064 {
00065     return fv::ddtScheme<vector>::New
00066     (
00067         vf.mesh(),
00068         vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
00069     )().meshPhi(vf);
00070 }
00071 
00072 
00073 void Foam::fvc::makeRelative
00074 (
00075     surfaceScalarField& phi,
00076     const volVectorField& U
00077 )
00078 {
00079     if (phi.mesh().moving())
00080     {
00081         phi -= fvc::meshPhi(U);
00082     }
00083 }
00084 
00085 void Foam::fvc::makeRelative
00086 (
00087     surfaceScalarField& phi,
00088     const dimensionedScalar& rho,
00089     const volVectorField& U
00090 )
00091 {
00092     if (phi.mesh().moving())
00093     {
00094         phi -= rho*fvc::meshPhi(rho, U);
00095     }
00096 }
00097 
00098 void Foam::fvc::makeRelative
00099 (
00100     surfaceScalarField& phi,
00101     const volScalarField& rho,
00102     const volVectorField& U
00103 )
00104 {
00105     if (phi.mesh().moving())
00106     {
00107         phi -= fvc::interpolate(rho)*fvc::meshPhi(rho, U);
00108     }
00109 }
00110 
00111 
00112 void Foam::fvc::makeAbsolute
00113 (
00114     surfaceScalarField& phi,
00115     const volVectorField& U
00116 )
00117 {
00118     if (phi.mesh().moving())
00119     {
00120         phi += fvc::meshPhi(U);
00121     }
00122 }
00123 
00124 void Foam::fvc::makeAbsolute
00125 (
00126     surfaceScalarField& phi,
00127     const dimensionedScalar& rho,
00128     const volVectorField& U
00129 )
00130 {
00131     if (phi.mesh().moving())
00132     {
00133         phi += rho*fvc::meshPhi(rho, U);
00134     }
00135 }
00136 
00137 void Foam::fvc::makeAbsolute
00138 (
00139     surfaceScalarField& phi,
00140     const volScalarField& rho,
00141     const volVectorField& U
00142 )
00143 {
00144     if (phi.mesh().moving())
00145     {
00146         phi += fvc::interpolate(rho)*fvc::meshPhi(rho, U);
00147     }
00148 }
00149 
00150 
00151 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines