fixed writing out entries in advective bc
[OpenFOAM-1.6-ext.git] / src / lagrangian / intermediate / submodels / Kinematic / InjectionModel / ManualInjection / ManualInjection.C
blobf231d316c020c11e93fe0981d793aec4d59490d2
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "ManualInjection.H"
28 #include "mathematicalConstants.H"
30 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
32 template<class CloudType>
33 Foam::label Foam::ManualInjection<CloudType>::parcelsToInject
35     const scalar time0,
36     const scalar time1
37 ) const
39     if ((0.0 >= time0) && (0.0 < time1))
40     {
41         return positions_.size();
42     }
43     else
44     {
45         return 0;
46     }
50 template<class CloudType>
51 Foam::scalar Foam::ManualInjection<CloudType>::volumeToInject
53     const scalar time0,
54     const scalar time1
55 ) const
57     // All parcels introduced at SOI
58     if ((0.0 >= time0) && (0.0 < time1))
59     {
60         return this->volumeTotal_;
61     }
62     else
63     {
64         return 0.0;
65     }
69 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
71 template<class CloudType>
72 Foam::ManualInjection<CloudType>::ManualInjection
74     const dictionary& dict,
75     CloudType& owner
78     InjectionModel<CloudType>(dict, owner, typeName),
79     positionsFile_(this->coeffDict().lookup("positionsFile")),
80     positions_
81     (
82         IOobject
83         (
84             positionsFile_,
85             owner.db().time().constant(),
86             owner.mesh(),
87             IOobject::MUST_READ,
88             IOobject::NO_WRITE
89         )
90     ),
91     diameters_(positions_.size()),
92     U0_(this->coeffDict().lookup("U0")),
93     parcelPDF_
94     (
95         pdf::New
96         (
97             this->coeffDict().subDict("parcelPDF"),
98             owner.rndGen()
99         )
100     )
102     // Construct parcel diameters
103     forAll(diameters_, i)
104     {
105         diameters_[i] = parcelPDF_->sample();
106     }
108     // Determine volume of particles to inject
109     this->volumeTotal_ = sum(pow3(diameters_))*mathematicalConstant::pi/6.0;
113 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
115 template<class CloudType>
116 Foam::ManualInjection<CloudType>::~ManualInjection()
120 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
122 template<class CloudType>
123 bool Foam::ManualInjection<CloudType>::active() const
125     return true;
129 template<class CloudType>
130 Foam::scalar Foam::ManualInjection<CloudType>::timeEnd() const
132     // Not used
133     return this->SOI_;
137 template<class CloudType>
138 void Foam::ManualInjection<CloudType>::setPositionAndCell
140     const label parcelI,
141     const label,
142     const scalar,
143     vector& position,
144     label& cellOwner
147     position = positions_[parcelI];
148     this->findCellAtPosition(cellOwner, position);
152 template<class CloudType>
153 void Foam::ManualInjection<CloudType>::setProperties
155     const label parcelI,
156     const label,
157     const scalar,
158     typename CloudType::parcelType& parcel
161     // set particle velocity
162     parcel.U() = U0_;
164     // set particle diameter
165     parcel.d() = diameters_[parcelI];
169 template<class CloudType>
170 bool Foam::ManualInjection<CloudType>::fullyDescribed() const
172     return false;
176 template<class CloudType>
177 bool Foam::ManualInjection<CloudType>::validInjection(const label)
179     return true;
183 // ************************************************************************* //