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

PatchInjection.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) 2009-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 "PatchInjection.H"
00027 #include <lagrangianIntermediate/DataEntry.H>
00028 #include <pdf/pdf.H>
00029 
00030 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
00031 
00032 template<class CloudType>
00033 Foam::label Foam::PatchInjection<CloudType>::parcelsToInject
00034 (
00035     const scalar time0,
00036     const scalar time1
00037 ) const
00038 {
00039     if ((time0 >= 0.0) && (time0 < duration_))
00040     {
00041         return round(fraction_*(time1 - time0)*parcelsPerSecond_);
00042     }
00043     else
00044     {
00045         return 0;
00046     }
00047 }
00048 
00049 
00050 template<class CloudType>
00051 Foam::scalar Foam::PatchInjection<CloudType>::volumeToInject
00052 (
00053     const scalar time0,
00054     const scalar time1
00055 ) const
00056 {
00057     if ((time0 >= 0.0) && (time0 < duration_))
00058     {
00059         return fraction_*volumeFlowRate_().integrate(time0, time1);
00060     }
00061     else
00062     {
00063         return 0.0;
00064     }
00065 }
00066 
00067 
00068 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00069 
00070 template<class CloudType>
00071 Foam::PatchInjection<CloudType>::PatchInjection
00072 (
00073     const dictionary& dict,
00074     CloudType& owner
00075 )
00076 :
00077     InjectionModel<CloudType>(dict, owner, typeName),
00078     patchName_(this->coeffDict().lookup("patchName")),
00079     duration_(readScalar(this->coeffDict().lookup("duration"))),
00080     parcelsPerSecond_
00081     (
00082         readScalar(this->coeffDict().lookup("parcelsPerSecond"))
00083     ),
00084     U0_(this->coeffDict().lookup("U0")),
00085     volumeFlowRate_
00086     (
00087         DataEntry<scalar>::New
00088         (
00089             "volumeFlowRate",
00090             this->coeffDict()
00091         )
00092     ),
00093     parcelPDF_
00094     (
00095         pdfs::pdf::New
00096         (
00097             this->coeffDict().subDict("parcelPDF"),
00098             owner.rndGen()
00099         )
00100     ),
00101     cellOwners_(),
00102     fraction_(1.0)
00103 {
00104     label patchId = owner.mesh().boundaryMesh().findPatchID(patchName_);
00105 
00106     if (patchId < 0)
00107     {
00108         FatalErrorIn
00109         (
00110             "PatchInjection<CloudType>::PatchInjection"
00111             "("
00112                 "const dictionary&, "
00113                 "CloudType&"
00114             ")"
00115         )   << "Requested patch " << patchName_ << " not found" << nl
00116             << "Available patches are: " << owner.mesh().boundaryMesh().names()
00117             << nl << exit(FatalError);
00118     }
00119 
00120     const polyPatch& patch = owner.mesh().boundaryMesh()[patchId];
00121 
00122     cellOwners_ = patch.faceCells();
00123 
00124     label patchSize = cellOwners_.size();
00125     label totalPatchSize = patchSize;
00126     reduce(totalPatchSize, sumOp<label>());
00127     fraction_ = scalar(patchSize)/totalPatchSize;
00128 
00129     // Set total volume/mass to inject
00130     this->volumeTotal_ = fraction_*volumeFlowRate_().integrate(0.0, duration_);
00131     this->massTotal_ *= fraction_;
00132 }
00133 
00134 
00135 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00136 
00137 template<class CloudType>
00138 Foam::PatchInjection<CloudType>::~PatchInjection()
00139 {}
00140 
00141 
00142 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00143 
00144 template<class CloudType>
00145 bool Foam::PatchInjection<CloudType>::active() const
00146 {
00147     return true;
00148 }
00149 
00150 
00151 template<class CloudType>
00152 Foam::scalar Foam::PatchInjection<CloudType>::timeEnd() const
00153 {
00154     return this->SOI_ + duration_;
00155 }
00156 
00157 
00158 template<class CloudType>
00159 void Foam::PatchInjection<CloudType>::setPositionAndCell
00160 (
00161     const label,
00162     const label,
00163     const scalar,
00164     vector& position,
00165     label& cellOwner
00166 )
00167 {
00168     if (cellOwners_.size() > 0)
00169     {
00170         label cellI = this->owner().rndGen().integer(0, cellOwners_.size() - 1);
00171 
00172         cellOwner = cellOwners_[cellI];
00173         position = this->owner().mesh().C()[cellOwner];
00174     }
00175     else
00176     {
00177         cellOwner = -1;
00178         // dummy position
00179         position = pTraits<vector>::max;
00180     }
00181 }
00182 
00183 
00184 template<class CloudType>
00185 void Foam::PatchInjection<CloudType>::setProperties
00186 (
00187     const label,
00188     const label,
00189     const scalar,
00190     typename CloudType::parcelType& parcel
00191 )
00192 {
00193     // set particle velocity
00194     parcel.U() = U0_;
00195 
00196     // set particle diameter
00197     parcel.d() = parcelPDF_->sample();
00198 }
00199 
00200 
00201 template<class CloudType>
00202 bool Foam::PatchInjection<CloudType>::fullyDescribed() const
00203 {
00204     return false;
00205 }
00206 
00207 
00208 template<class CloudType>
00209 bool Foam::PatchInjection<CloudType>::validInjection(const label)
00210 {
00211     return true;
00212 }
00213 
00214 
00215 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines