fixed writing out entries in advective bc
[OpenFOAM-1.6-ext.git] / src / lagrangian / intermediate / clouds / Templates / KinematicCloud / KinematicCloud.C
blob5f25f60fdb135216a5dba9d5b75e6f488bf76993
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 // * * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * //
39 template<class ParcelType>
40 void Foam::KinematicCloud<ParcelType>::preEvolve()
42     this->dispersion().cacheFields(true);
43     forces_.cacheFields(true);
47 template<class ParcelType>
48 void Foam::KinematicCloud<ParcelType>::evolveCloud()
50     autoPtr<interpolation<scalar> > rhoInterpolator =
51         interpolation<scalar>::New
52         (
53             interpolationSchemes_,
54             rho_
55         );
57     autoPtr<interpolation<vector> > UInterpolator =
58         interpolation<vector>::New
59         (
60             interpolationSchemes_,
61             U_
62         );
64     autoPtr<interpolation<scalar> > muInterpolator =
65         interpolation<scalar>::New
66         (
67             interpolationSchemes_,
68             mu_
69         );
71     typename ParcelType::trackData td
72     (
73         *this,
74         constProps_,
75         rhoInterpolator(),
76         UInterpolator(),
77         muInterpolator(),
78         g_.value()
79     );
81     this->injection().inject(td);
83     if (coupled_)
84     {
85         resetSourceTerms();
86     }
88     Cloud<ParcelType>::move(td);
92 template<class ParcelType>
93 void Foam::KinematicCloud<ParcelType>::postEvolve()
95     if (debug)
96     {
97         this->writePositions();
98     }
100     this->dispersion().cacheFields(false);
101     forces_.cacheFields(false);
103     this->postProcessing().post();
107 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
109 template<class ParcelType>
110 Foam::KinematicCloud<ParcelType>::KinematicCloud
112     const word& cloudName,
113     const volScalarField& rho,
114     const volVectorField& U,
115     const volScalarField& mu,
116     const dimensionedVector& g,
117     bool readFields
120     Cloud<ParcelType>(rho.mesh(), cloudName, false),
121     kinematicCloud(),
122     mesh_(rho.mesh()),
123     particleProperties_
124     (
125         IOobject
126         (
127             cloudName + "Properties",
128             rho.mesh().time().constant(),
129             rho.mesh(),
130             IOobject::MUST_READ,
131             IOobject::NO_WRITE
132         )
133     ),
134     constProps_(particleProperties_),
135     active_(particleProperties_.lookup("active")),
136     parcelTypeId_(readLabel(particleProperties_.lookup("parcelTypeId"))),
137     coupled_(particleProperties_.lookup("coupled")),
138     cellValueSourceCorrection_
139     (
140         particleProperties_.lookup("cellValueSourceCorrection")
141     ),
142     rndGen_(label(0)),
143     rho_(rho),
144     U_(U),
145     mu_(mu),
146     g_(g),
147     forces_(mesh_, particleProperties_, g_.value()),
148     interpolationSchemes_(particleProperties_.subDict("interpolationSchemes")),
149     dispersionModel_
150     (
151         DispersionModel<KinematicCloud<ParcelType> >::New
152         (
153             particleProperties_,
154             *this
155         )
156     ),
157     dragModel_
158     (
159         DragModel<KinematicCloud<ParcelType> >::New
160         (
161             particleProperties_,
162             *this
163         )
164     ),
165     injectionModel_
166     (
167         InjectionModel<KinematicCloud<ParcelType> >::New
168         (
169             particleProperties_,
170             *this
171         )
172     ),
173     patchInteractionModel_
174     (
175         PatchInteractionModel<KinematicCloud<ParcelType> >::New
176         (
177             particleProperties_,
178             *this
179         )
180     ),
181     postProcessingModel_
182     (
183         PostProcessingModel<KinematicCloud<ParcelType> >::New
184         (
185             this->particleProperties_,
186             *this
187         )
188     ),
189     UIntegrator_
190     (
191         vectorIntegrationScheme::New
192         (
193             "U",
194             particleProperties_.subDict("integrationSchemes")
195         )
196     ),
197     UTrans_
198     (
199         IOobject
200         (
201             this->name() + "UTrans",
202             this->db().time().timeName(),
203             this->db(),
204             IOobject::NO_READ,
205             IOobject::NO_WRITE,
206             false
207         ),
208         mesh_,
209         dimensionedVector("zero", dimMass*dimVelocity, vector::zero)
210     )
212     if (readFields)
213     {
214         ParcelType::readFields(*this);
215     }
219 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
221 template<class ParcelType>
222 Foam::KinematicCloud<ParcelType>::~KinematicCloud()
226 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
228 template<class ParcelType>
229 void Foam::KinematicCloud<ParcelType>::checkParcelProperties
231     ParcelType& parcel,
232     const scalar lagrangianDt,
233     const bool fullyDescribed
236     if (!fullyDescribed)
237     {
238         parcel.rho() = constProps_.rho0();
239     }
241     scalar carrierDt = this->db().time().deltaTValue();
242     parcel.stepFraction() = (carrierDt - lagrangianDt)/carrierDt;
246 template<class ParcelType>
247 void Foam::KinematicCloud<ParcelType>::resetSourceTerms()
249     UTrans_.field() = vector::zero;
253 template<class ParcelType>
254 void Foam::KinematicCloud<ParcelType>::evolve()
256     if (active_)
257     {
258         preEvolve();
260         evolveCloud();
262         postEvolve();
264         info();
265         Info<< endl;
266     }
270 template<class ParcelType>
271 void Foam::KinematicCloud<ParcelType>::info() const
273     Info<< "Cloud: " << this->name() << nl
274         << "    Total number of parcels added   = "
275         << this->injection().parcelsAddedTotal() << nl
276         << "    Total mass introduced           = "
277         << this->injection().massInjected() << nl
278         << "    Current number of parcels       = "
279         << returnReduce(this->size(), sumOp<label>()) << nl
280         << "    Current mass in system          = "
281         << returnReduce(massInSystem(), sumOp<scalar>()) << nl;
285 // ************************************************************************* //