Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "FieldActivatedInjection.H"
00027 #include <finiteVolume/volFields.H>
00028 #include <OpenFOAM/mathematicalConstants.H>
00029
00030
00031
00032 template<class CloudType>
00033 Foam::label Foam::FieldActivatedInjection<CloudType>::parcelsToInject
00034 (
00035 const scalar time0,
00036 const scalar time1
00037 ) const
00038 {
00039 if (sum(nParcelsInjected_) < nParcelsPerInjector_*positions_.size())
00040 {
00041 return positions_.size();
00042 }
00043 else
00044 {
00045 return 0;
00046 }
00047 }
00048
00049
00050 template<class CloudType>
00051 Foam::scalar Foam::FieldActivatedInjection<CloudType>::volumeToInject
00052 (
00053 const scalar time0,
00054 const scalar time1
00055 ) const
00056 {
00057 if (sum(nParcelsInjected_) < nParcelsPerInjector_*positions_.size())
00058 {
00059 return this->volumeTotal_/nParcelsPerInjector_;
00060 }
00061 else
00062 {
00063 return 0;
00064 }
00065 }
00066
00067
00068
00069
00070 template<class CloudType>
00071 Foam::FieldActivatedInjection<CloudType>::FieldActivatedInjection
00072 (
00073 const dictionary& dict,
00074 CloudType& owner
00075 )
00076 :
00077 InjectionModel<CloudType>(dict, owner, typeName),
00078 factor_(readScalar(this->coeffDict().lookup("factor"))),
00079 referenceField_
00080 (
00081 owner.db().objectRegistry::lookupObject<volScalarField>
00082 (
00083 this->coeffDict().lookup("referenceField")
00084 )
00085 ),
00086 thresholdField_
00087 (
00088 owner.db().objectRegistry::lookupObject<volScalarField>
00089 (
00090 this->coeffDict().lookup("thresholdField")
00091 )
00092 ),
00093 positionsFile_(this->coeffDict().lookup("positionsFile")),
00094 positions_
00095 (
00096 IOobject
00097 (
00098 positionsFile_,
00099 owner.db().time().constant(),
00100 owner.mesh(),
00101 IOobject::MUST_READ,
00102 IOobject::NO_WRITE
00103 )
00104 ),
00105 injectorCells_(positions_.size()),
00106 nParcelsPerInjector_
00107 (
00108 readLabel(this->coeffDict().lookup("parcelsPerInjector"))
00109 ),
00110 nParcelsInjected_(positions_.size(), 0),
00111 U0_(this->coeffDict().lookup("U0")),
00112 diameters_(positions_.size()),
00113 parcelPDF_
00114 (
00115 pdfs::pdf::New
00116 (
00117 this->coeffDict().subDict("parcelPDF"),
00118 owner.rndGen()
00119 )
00120 )
00121 {
00122
00123 forAll(diameters_, i)
00124 {
00125 diameters_[i] = parcelPDF_->sample();
00126 }
00127
00128
00129 this->volumeTotal_ =
00130 nParcelsPerInjector_*sum(pow3(diameters_))*mathematicalConstant::pi/6.0;
00131
00132
00133 forAll(positions_, i)
00134 {
00135 this->findCellAtPosition
00136 (
00137 injectorCells_[i],
00138 positions_[i]
00139 );
00140 }
00141 }
00142
00143
00144
00145
00146 template<class CloudType>
00147 Foam::FieldActivatedInjection<CloudType>::~FieldActivatedInjection()
00148 {}
00149
00150
00151
00152
00153 template<class CloudType>
00154 bool Foam::FieldActivatedInjection<CloudType>::active() const
00155 {
00156 return true;
00157 }
00158
00159
00160 template<class CloudType>
00161 Foam::scalar Foam::FieldActivatedInjection<CloudType>::timeEnd() const
00162 {
00163 return GREAT;
00164 }
00165
00166
00167 template<class CloudType>
00168 void Foam::FieldActivatedInjection<CloudType>::setPositionAndCell
00169 (
00170 const label parcelI,
00171 const label,
00172 const scalar,
00173 vector& position,
00174 label& cellOwner
00175 )
00176 {
00177 position = positions_[parcelI];
00178 cellOwner = injectorCells_[parcelI];
00179 }
00180
00181
00182 template<class CloudType>
00183 void Foam::FieldActivatedInjection<CloudType>::setProperties
00184 (
00185 const label parcelI,
00186 const label,
00187 const scalar,
00188 typename CloudType::parcelType& parcel
00189 )
00190 {
00191
00192 parcel.U() = U0_;
00193
00194
00195 parcel.d() = diameters_[parcelI];
00196 }
00197
00198
00199 template<class CloudType>
00200 bool Foam::FieldActivatedInjection<CloudType>::fullyDescribed() const
00201 {
00202 return false;
00203 }
00204
00205
00206 template<class CloudType>
00207 bool Foam::FieldActivatedInjection<CloudType>::validInjection
00208 (
00209 const label parcelI
00210 )
00211 {
00212 const label cellI = injectorCells_[parcelI];
00213
00214 if
00215 (
00216 nParcelsInjected_[parcelI] < nParcelsPerInjector_
00217 && factor_*referenceField_[cellI] > thresholdField_[cellI]
00218 )
00219 {
00220 nParcelsInjected_[parcelI]++;
00221 return true;
00222 }
00223
00224 return false;
00225 }
00226
00227
00228