FieldActivatedInjection.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) 2011-2012 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 
27 #include "volFields.H"
28 #include "mathematicalConstants.H"
29 
30 using namespace Foam::constant::mathematical;
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class CloudType>
36 (
37  const dictionary& dict,
38  CloudType& owner,
39  const word& modelName
40 )
41 :
42  InjectionModel<CloudType>(dict, owner, modelName, typeName),
43  factor_(readScalar(this->coeffDict().lookup("factor"))),
44  referenceField_
45  (
46  owner.db().objectRegistry::template lookupObject<volScalarField>
47  (
48  this->coeffDict().lookup("referenceField")
49  )
50  ),
51  thresholdField_
52  (
53  owner.db().objectRegistry::template lookupObject<volScalarField>
54  (
55  this->coeffDict().lookup("thresholdField")
56  )
57  ),
58  positionsFile_(this->coeffDict().lookup("positionsFile")),
59  positions_
60  (
61  IOobject
62  (
63  positionsFile_,
64  owner.db().time().constant(),
65  owner.mesh(),
68  )
69  ),
70  injectorCells_(positions_.size()),
71  injectorTetFaces_(positions_.size()),
72  injectorTetPts_(positions_.size()),
73  nParcelsPerInjector_
74  (
75  readLabel(this->coeffDict().lookup("parcelsPerInjector"))
76  ),
77  nParcelsInjected_(positions_.size(), 0),
78  U0_(this->coeffDict().lookup("U0")),
79  diameters_(positions_.size()),
80  sizeDistribution_
81  (
83  (
84  this->coeffDict().subDict("sizeDistribution"),
85  owner.rndGen()
86  )
87  )
88 {
89  // Construct parcel diameters - one per injector cell
90  forAll(diameters_, i)
91  {
92  diameters_[i] = sizeDistribution_->sample();
93  }
94 
95  // Determine total volume of particles to inject
96  this->volumeTotal_ =
97  nParcelsPerInjector_*sum(pow3(diameters_))*pi/6.0;
98 
99  updateMesh();
100 }
101 
102 
103 template<class CloudType>
105 (
107 )
108 :
110  factor_(im.factor_),
111  referenceField_(im.referenceField_),
112  thresholdField_(im.thresholdField_),
113  positionsFile_(im.positionsFile_),
114  positions_(im.positions_),
115  injectorCells_(im.injectorCells_),
116  injectorTetFaces_(im.injectorTetFaces_),
117  injectorTetPts_(im.injectorTetPts_),
118  nParcelsPerInjector_(im.nParcelsPerInjector_),
119  nParcelsInjected_(im.nParcelsInjected_),
120  U0_(im.U0_),
121  diameters_(im.diameters_),
122  sizeDistribution_(im.sizeDistribution_, false)
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
127 
128 template<class CloudType>
130 {}
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
135 template<class CloudType>
137 {
138  // Set/cache the injector cells
139  forAll(positions_, i)
140  {
141  this->findCellAtPosition
142  (
143  injectorCells_[i],
144  injectorTetFaces_[i],
145  injectorTetPts_[i],
146  positions_[i]
147  );
148  }
149 }
150 
151 
152 template<class CloudType>
154 {
155  return GREAT;
156 }
157 
158 
159 template<class CloudType>
161 (
162  const scalar time0,
163  const scalar time1
164 )
165 {
166  if (sum(nParcelsInjected_) < nParcelsPerInjector_*positions_.size())
167  {
168  return positions_.size();
169  }
170  else
171  {
172  return 0;
173  }
174 }
175 
176 
177 template<class CloudType>
179 (
180  const scalar time0,
181  const scalar time1
182 )
183 {
184  if (sum(nParcelsInjected_) < nParcelsPerInjector_*positions_.size())
185  {
186  return this->volumeTotal_/nParcelsPerInjector_;
187  }
188  else
189  {
190  return 0;
191  }
192 }
193 
194 
195 template<class CloudType>
197 (
198  const label parcelI,
199  const label,
200  const scalar,
201  vector& position,
202  label& cellOwner,
203  label& tetFaceI,
204  label& tetPtI
205 )
206 {
207  position = positions_[parcelI];
208  cellOwner = injectorCells_[parcelI];
209  tetFaceI = injectorTetFaces_[parcelI];
210  tetPtI = injectorTetPts_[parcelI];
211 }
212 
213 
214 template<class CloudType>
216 (
217  const label parcelI,
218  const label,
219  const scalar,
220  typename CloudType::parcelType& parcel
221 )
222 {
223  // set particle velocity
224  parcel.U() = U0_;
225 
226  // set particle diameter
227  parcel.d() = diameters_[parcelI];
228 }
229 
230 
231 template<class CloudType>
233 {
234  return false;
235 }
236 
237 
238 template<class CloudType>
240 (
241  const label parcelI
242 )
243 {
244  const label cellI = injectorCells_[parcelI];
245 
246  if
247  (
248  nParcelsInjected_[parcelI] < nParcelsPerInjector_
249  && factor_*referenceField_[cellI] > thresholdField_[cellI]
250  )
251  {
252  nParcelsInjected_[parcelI]++;
253  return true;
254  }
255 
256  return false;
257 }
258 
259 
260 // ************************************************************************* //
Foam::FieldActivatedInjection::sizeDistribution_
const autoPtr< distributionModels::distributionModel > sizeDistribution_
Parcel size distribution model.
Definition: FieldActivatedInjection.H:115
volFields.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
mathematicalConstants.H
Foam::FieldActivatedInjection::referenceField_
const volScalarField & referenceField_
Reference field.
Definition: FieldActivatedInjection.H:75
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::FieldActivatedInjection::positionsFile_
const word positionsFile_
Name of file containing positions data.
Definition: FieldActivatedInjection.H:84
Foam::FieldActivatedInjection::U0_
const vector U0_
Initial parcel velocity.
Definition: FieldActivatedInjection.H:108
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::FieldActivatedInjection::injectorTetPts_
labelList injectorTetPts_
List of tetPt labels corresponding to injector positions.
Definition: FieldActivatedInjection.H:96
FieldActivatedInjection.H
Foam::FieldActivatedInjection::~FieldActivatedInjection
virtual ~FieldActivatedInjection()
Destructor.
Definition: FieldActivatedInjection.C:129
Foam::FieldActivatedInjection::thresholdField_
const volScalarField & thresholdField_
Threshold field.
Definition: FieldActivatedInjection.H:78
Foam::FieldActivatedInjection::FieldActivatedInjection
FieldActivatedInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: FieldActivatedInjection.C:36
Foam::InjectionModel
Templated injection model class.
Definition: InjectionModel.H:67
Foam::distributionModels::distributionModel::New
static autoPtr< distributionModel > New(const dictionary &dict, cachedRandom &rndGen)
Selector.
Definition: distributionModelNew.C:32
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::FieldActivatedInjection::fullyDescribed
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Definition: FieldActivatedInjection.C:232
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::DSMCCloud::rndGen
Random & rndGen()
Return refernce to the random object.
Definition: DSMCCloudI.H:121
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
Foam::FieldActivatedInjection::injectorTetFaces_
labelList injectorTetFaces_
List of tetFace labels corresponding to injector positions.
Definition: FieldActivatedInjection.H:93
Foam::FieldActivatedInjection::updateMesh
virtual void updateMesh()
Set injector locations when mesh is updated.
Definition: FieldActivatedInjection.C:136
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
Foam::FieldActivatedInjection::factor_
const scalar factor_
Factor to apply to reference field.
Definition: FieldActivatedInjection.H:72
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::FieldActivatedInjection::parcelsToInject
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
Definition: FieldActivatedInjection.C:161
Foam::FieldActivatedInjection::volumeToInject
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
Definition: FieldActivatedInjection.C:179
Foam::FieldActivatedInjection::timeEnd
scalar timeEnd() const
Return the end-of-injection time.
Definition: FieldActivatedInjection.C:153
Foam::FieldActivatedInjection
Conditional injection at specified positions.
Definition: FieldActivatedInjection.H:63
Foam::constant::mathematical
mathematical constants.
Foam::FieldActivatedInjection::setPositionAndCell
virtual void setPositionAndCell(const label parcelI, const label nParcels, const scalar time, vector &position, label &cellOwner, label &tetFaceI, label &tetPtI)
Set the injection position and owner cell, tetFace and tetPt.
Definition: FieldActivatedInjection.C:197
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::FieldActivatedInjection::positions_
vectorIOField positions_
Field of injector (x,y,z) positions.
Definition: FieldActivatedInjection.H:87
Foam::FieldActivatedInjection::diameters_
scalarList diameters_
List of parcel diameters.
Definition: FieldActivatedInjection.H:111
Foam::Vector< scalar >
Foam::FieldActivatedInjection::setProperties
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
Definition: FieldActivatedInjection.C:216
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
Foam::readLabel
label readLabel(Istream &is)
Definition: label.H:64
Foam::constant::mathematical::pi
const scalar pi(M_PI)
Foam::FieldActivatedInjection::injectorCells_
labelList injectorCells_
List of cell labels corresponding to injector positions.
Definition: FieldActivatedInjection.H:90
Foam::FieldActivatedInjection::validInjection
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
Definition: FieldActivatedInjection.C:240
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
Foam::FieldActivatedInjection::nParcelsInjected_
labelList nParcelsInjected_
List of number of parcels injected for each injector.
Definition: FieldActivatedInjection.H:102
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
Foam::FieldActivatedInjection::nParcelsPerInjector_
const label nParcelsPerInjector_
Number of parcels per injector.
Definition: FieldActivatedInjection.H:99