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 "KinematicCloud_.H"
00027 #include <lagrangianIntermediate/IntegrationScheme.H>
00028 #include <finiteVolume/interpolation.H>
00029
00030 #include <lagrangianIntermediate/DispersionModel.H>
00031 #include <lagrangianIntermediate/DragModel.H>
00032 #include <lagrangianIntermediate/InjectionModel.H>
00033 #include <lagrangianIntermediate/PatchInteractionModel.H>
00034 #include <lagrangianIntermediate/PostProcessingModel.H>
00035
00036
00037
00038 template<class ParcelType>
00039 void Foam::KinematicCloud<ParcelType>::preEvolve()
00040 {
00041 this->dispersion().cacheFields(true);
00042 forces_.cacheFields(true);
00043 }
00044
00045
00046 template<class ParcelType>
00047 void Foam::KinematicCloud<ParcelType>::evolveCloud()
00048 {
00049 autoPtr<interpolation<scalar> > rhoInterpolator =
00050 interpolation<scalar>::New
00051 (
00052 interpolationSchemes_,
00053 rho_
00054 );
00055
00056 autoPtr<interpolation<vector> > UInterpolator =
00057 interpolation<vector>::New
00058 (
00059 interpolationSchemes_,
00060 U_
00061 );
00062
00063 autoPtr<interpolation<scalar> > muInterpolator =
00064 interpolation<scalar>::New
00065 (
00066 interpolationSchemes_,
00067 mu_
00068 );
00069
00070 typename ParcelType::trackData td
00071 (
00072 *this,
00073 constProps_,
00074 rhoInterpolator(),
00075 UInterpolator(),
00076 muInterpolator(),
00077 g_.value()
00078 );
00079
00080 this->injection().inject(td);
00081
00082 if (coupled_)
00083 {
00084 resetSourceTerms();
00085 }
00086
00087 Cloud<ParcelType>::move(td);
00088 }
00089
00090
00091 template<class ParcelType>
00092 void Foam::KinematicCloud<ParcelType>::postEvolve()
00093 {
00094 if (debug)
00095 {
00096 this->writePositions();
00097 }
00098
00099 this->dispersion().cacheFields(false);
00100 forces_.cacheFields(false);
00101
00102 this->postProcessing().post();
00103 }
00104
00105
00106
00107
00108 template<class ParcelType>
00109 Foam::KinematicCloud<ParcelType>::KinematicCloud
00110 (
00111 const word& cloudName,
00112 const volScalarField& rho,
00113 const volVectorField& U,
00114 const volScalarField& mu,
00115 const dimensionedVector& g,
00116 bool readFields
00117 )
00118 :
00119 Cloud<ParcelType>(rho.mesh(), cloudName, false),
00120 kinematicCloud(),
00121 mesh_(rho.mesh()),
00122 particleProperties_
00123 (
00124 IOobject
00125 (
00126 cloudName + "Properties",
00127 rho.mesh().time().constant(),
00128 rho.mesh(),
00129 IOobject::MUST_READ,
00130 IOobject::NO_WRITE
00131 )
00132 ),
00133 constProps_(particleProperties_),
00134 active_(particleProperties_.lookup("active")),
00135 parcelTypeId_(readLabel(particleProperties_.lookup("parcelTypeId"))),
00136 coupled_(particleProperties_.lookup("coupled")),
00137 cellValueSourceCorrection_
00138 (
00139 particleProperties_.lookup("cellValueSourceCorrection")
00140 ),
00141 rndGen_(label(0)),
00142 rho_(rho),
00143 U_(U),
00144 mu_(mu),
00145 g_(g),
00146 forces_(mesh_, particleProperties_, g_.value()),
00147 interpolationSchemes_(particleProperties_.subDict("interpolationSchemes")),
00148 dispersionModel_
00149 (
00150 DispersionModel<KinematicCloud<ParcelType> >::New
00151 (
00152 particleProperties_,
00153 *this
00154 )
00155 ),
00156 dragModel_
00157 (
00158 DragModel<KinematicCloud<ParcelType> >::New
00159 (
00160 particleProperties_,
00161 *this
00162 )
00163 ),
00164 injectionModel_
00165 (
00166 InjectionModel<KinematicCloud<ParcelType> >::New
00167 (
00168 particleProperties_,
00169 *this
00170 )
00171 ),
00172 patchInteractionModel_
00173 (
00174 PatchInteractionModel<KinematicCloud<ParcelType> >::New
00175 (
00176 particleProperties_,
00177 *this
00178 )
00179 ),
00180 postProcessingModel_
00181 (
00182 PostProcessingModel<KinematicCloud<ParcelType> >::New
00183 (
00184 this->particleProperties_,
00185 *this
00186 )
00187 ),
00188 UIntegrator_
00189 (
00190 vectorIntegrationScheme::New
00191 (
00192 "U",
00193 particleProperties_.subDict("integrationSchemes")
00194 )
00195 ),
00196 UTrans_
00197 (
00198 IOobject
00199 (
00200 this->name() + "UTrans",
00201 this->db().time().timeName(),
00202 this->db(),
00203 IOobject::NO_READ,
00204 IOobject::NO_WRITE,
00205 false
00206 ),
00207 mesh_,
00208 dimensionedVector("zero", dimMass*dimVelocity, vector::zero)
00209 )
00210 {
00211 if (readFields)
00212 {
00213 ParcelType::readFields(*this);
00214 }
00215 }
00216
00217
00218
00219
00220 template<class ParcelType>
00221 Foam::KinematicCloud<ParcelType>::~KinematicCloud()
00222 {}
00223
00224
00225
00226
00227 template<class ParcelType>
00228 void Foam::KinematicCloud<ParcelType>::checkParcelProperties
00229 (
00230 ParcelType& parcel,
00231 const scalar lagrangianDt,
00232 const bool fullyDescribed
00233 )
00234 {
00235 if (!fullyDescribed)
00236 {
00237 parcel.rho() = constProps_.rho0();
00238 }
00239
00240 scalar carrierDt = this->db().time().deltaTValue();
00241 parcel.stepFraction() = (carrierDt - lagrangianDt)/carrierDt;
00242 }
00243
00244
00245 template<class ParcelType>
00246 void Foam::KinematicCloud<ParcelType>::resetSourceTerms()
00247 {
00248 UTrans_.field() = vector::zero;
00249 }
00250
00251
00252 template<class ParcelType>
00253 void Foam::KinematicCloud<ParcelType>::evolve()
00254 {
00255 if (active_)
00256 {
00257 preEvolve();
00258
00259 evolveCloud();
00260
00261 postEvolve();
00262
00263 info();
00264 Info<< endl;
00265 }
00266 }
00267
00268
00269 template<class ParcelType>
00270 void Foam::KinematicCloud<ParcelType>::info() const
00271 {
00272 Info<< "Cloud: " << this->name() << nl
00273 << " Total number of parcels added = "
00274 << this->injection().parcelsAddedTotal() << nl
00275 << " Total mass introduced = "
00276 << this->injection().massInjected() << nl
00277 << " Current number of parcels = "
00278 << returnReduce(this->size(), sumOp<label>()) << nl
00279 << " Current mass in system = "
00280 << returnReduce(massInSystem(), sumOp<scalar>()) << nl;
00281 }
00282
00283
00284