initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / lagrangian / intermediate / clouds / Templates / KinematicCloud / KinematicCloud.C
blob650cb7ddcbdc9f47bb7ca2782d67d99e0d743a84
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
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 "KinematicCloud.H"
28 #include "IntegrationScheme.H"
29 #include "interpolation.H"
31 #include "DispersionModel.H"
32 #include "DragModel.H"
33 #include "InjectionModel.H"
34 #include "PatchInteractionModel.H"
35 #include "PostProcessingModel.H"
37 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
39 template<class ParcelType>
40 Foam::KinematicCloud<ParcelType>::KinematicCloud
42     const word& cloudName,
43     const volScalarField& rho,
44     const volVectorField& U,
45     const volScalarField& mu,
46     const dimensionedVector& g
49     Cloud<ParcelType>(rho.mesh(), cloudName, false),
50     kinematicCloud(),
51     mesh_(rho.mesh()),
52     particleProperties_
53     (
54         IOobject
55         (
56             cloudName + "Properties",
57             rho.mesh().time().constant(),
58             rho.mesh(),
59             IOobject::MUST_READ,
60             IOobject::NO_WRITE
61         )
62     ),
63     constProps_(particleProperties_),
64     parcelTypeId_(readLabel(particleProperties_.lookup("parcelTypeId"))),
65     coupled_(particleProperties_.lookup("coupled")),
66     cellValueSourceCorrection_
67     (
68         particleProperties_.lookup("cellValueSourceCorrection")
69     ),
70     rndGen_(label(0)),
71     rho_(rho),
72     U_(U),
73     mu_(mu),
74     g_(g),
75     forces_(mesh_, particleProperties_, g_.value()),
76     interpolationSchemes_(particleProperties_.subDict("interpolationSchemes")),
77     dispersionModel_
78     (
79         DispersionModel<KinematicCloud<ParcelType> >::New
80         (
81             particleProperties_,
82             *this
83         )
84     ),
85     dragModel_
86     (
87         DragModel<KinematicCloud<ParcelType> >::New
88         (
89             particleProperties_,
90             *this
91         )
92     ),
93     injectionModel_
94     (
95         InjectionModel<KinematicCloud<ParcelType> >::New
96         (
97             particleProperties_,
98             *this
99         )
100     ),
101     patchInteractionModel_
102     (
103         PatchInteractionModel<KinematicCloud<ParcelType> >::New
104         (
105             particleProperties_,
106             *this
107         )
108     ),
109     postProcessingModel_
110     (
111         PostProcessingModel<KinematicCloud<ParcelType> >::New
112         (
113             this->particleProperties_,
114             *this
115         )
116     ),
117     UIntegrator_
118     (
119         vectorIntegrationScheme::New
120         (
121             "U",
122             particleProperties_.subDict("integrationSchemes")
123         )
124     ),
125     UTrans_
126     (
127         IOobject
128         (
129             this->name() + "UTrans",
130             this->db().time().timeName(),
131             this->db(),
132             IOobject::NO_READ,
133             IOobject::NO_WRITE,
134             false
135         ),
136         mesh_,
137         dimensionedVector("zero", dimMass*dimVelocity, vector::zero)
138     )
142 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
144 template<class ParcelType>
145 Foam::KinematicCloud<ParcelType>::~KinematicCloud()
149 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
151 template<class ParcelType>
152 void Foam::KinematicCloud<ParcelType>::checkParcelProperties
154     ParcelType& parcel,
155     const scalar lagrangianDt,
156     const bool fullyDescribed
159     if (!fullyDescribed)
160     {
161         parcel.rho() = constProps_.rho0();
162     }
164     scalar carrierDt = this->db().time().deltaT().value();
165     parcel.stepFraction() = (carrierDt - lagrangianDt)/carrierDt;
169 template<class ParcelType>
170 void Foam::KinematicCloud<ParcelType>::resetSourceTerms()
172     UTrans_.field() = vector::zero;
176 template<class ParcelType>
177 void Foam::KinematicCloud<ParcelType>::preEvolve()
179     this->dispersion().cacheFields(true);
180     forces_.cacheFields(true);
184 template<class ParcelType>
185 void Foam::KinematicCloud<ParcelType>::postEvolve()
187     if (debug)
188     {
189         this->writePositions();
190     }
192     this->dispersion().cacheFields(false);
193     forces_.cacheFields(false);
195     this->postProcessing().post();
199 template<class ParcelType>
200 void Foam::KinematicCloud<ParcelType>::evolve()
202     preEvolve();
204     autoPtr<interpolation<scalar> > rhoInterpolator =
205         interpolation<scalar>::New
206         (
207             interpolationSchemes_,
208             rho_
209         );
211     autoPtr<interpolation<vector> > UInterpolator =
212         interpolation<vector>::New
213         (
214             interpolationSchemes_,
215             U_
216         );
218     autoPtr<interpolation<scalar> > muInterpolator =
219         interpolation<scalar>::New
220         (
221             interpolationSchemes_,
222             mu_
223         );
225     typename ParcelType::trackData td
226     (
227         *this,
228         constProps_,
229         rhoInterpolator(),
230         UInterpolator(),
231         muInterpolator(),
232         g_.value()
233     );
235     this->injection().inject(td);
237     if (coupled_)
238     {
239         resetSourceTerms();
240     }
242     Cloud<ParcelType>::move(td);
244     postEvolve();
248 template<class ParcelType>
249 void Foam::KinematicCloud<ParcelType>::info() const
251     Info<< "Cloud: " << this->name() << nl
252         << "    Total number of parcels added   = "
253         << returnReduce(this->injection().parcelsAddedTotal(), sumOp<label>())
254             << nl
255         << "    Total mass introduced           = "
256         << returnReduce(this->injection().massInjected(), sumOp<scalar>())
257             << nl
258         << "    Current number of parcels       = "
259         << returnReduce(this->size(), sumOp<label>()) << nl
260         << "    Current mass in system          = "
261         << returnReduce(massInSystem(), sumOp<scalar>()) << nl;
265 // ************************************************************************* //