FreeFOAM The Cross-Platform CFD Toolkit
ManualInjection.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
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 
28 
29 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
30 
31 template<class CloudType>
33 (
34  const scalar time0,
35  const scalar time1
36 ) const
37 {
38  if ((0.0 >= time0) && (0.0 < time1))
39  {
40  return positions_.size();
41  }
42  else
43  {
44  return 0;
45  }
46 }
47 
48 
49 template<class CloudType>
51 (
52  const scalar time0,
53  const scalar time1
54 ) const
55 {
56  // All parcels introduced at SOI
57  if ((0.0 >= time0) && (0.0 < time1))
58  {
59  return this->volumeTotal_;
60  }
61  else
62  {
63  return 0.0;
64  }
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
69 
70 template<class CloudType>
72 (
73  const dictionary& dict,
74  CloudType& owner
75 )
76 :
77  InjectionModel<CloudType>(dict, owner, typeName),
78  positionsFile_(this->coeffDict().lookup("positionsFile")),
79  positions_
80  (
81  IOobject
82  (
83  positionsFile_,
84  owner.db().time().constant(),
85  owner.mesh(),
86  IOobject::MUST_READ,
87  IOobject::NO_WRITE
88  )
89  ),
90  diameters_(positions_.size()),
91  U0_(this->coeffDict().lookup("U0")),
92  parcelPDF_
93  (
94  pdfs::pdf::New
95  (
96  this->coeffDict().subDict("parcelPDF"),
97  owner.rndGen()
98  )
99  )
100 {
101  // Construct parcel diameters
102  forAll(diameters_, i)
103  {
104  diameters_[i] = parcelPDF_->sample();
105  }
106 
107  // Determine volume of particles to inject
108  this->volumeTotal_ = sum(pow3(diameters_))*mathematicalConstant::pi/6.0;
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
113 
114 template<class CloudType>
116 {}
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
121 template<class CloudType>
123 {
124  return true;
125 }
126 
127 
128 template<class CloudType>
130 {
131  // Not used
132  return this->SOI_;
133 }
134 
135 
136 template<class CloudType>
138 (
139  const label parcelI,
140  const label,
141  const scalar,
142  vector& position,
143  label& cellOwner
144 )
145 {
146  position = positions_[parcelI];
147  this->findCellAtPosition(cellOwner, position);
148 }
149 
150 
151 template<class CloudType>
153 (
154  const label parcelI,
155  const label,
156  const scalar,
157  typename CloudType::parcelType& parcel
158 )
159 {
160  // set particle velocity
161  parcel.U() = U0_;
162 
163  // set particle diameter
164  parcel.d() = diameters_[parcelI];
165 }
166 
167 
168 template<class CloudType>
170 {
171  return false;
172 }
173 
174 
175 template<class CloudType>
177 {
178  return true;
179 }
180 
181 
182 // ************************ vim: set sw=4 sts=4 et: ************************ //