initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / lagrangian / dsmc / parcels / Templates / DsmcParcel / DsmcParcel.C
blobc52d92d20fef43858cfdd7075883003e6a043b76
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2009-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 "DsmcParcel.H"
28 #include "dimensionedConstants.H"
30 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
32 template<class ParcelType>
33 template<class TrackData>
34 bool Foam::DsmcParcel<ParcelType>::move
36     TrackData& td
39     ParcelType& p = static_cast<ParcelType&>(*this);
41     td.switchProcessor = false;
42     td.keepParticle = true;
44     const polyMesh& mesh = td.cloud().pMesh();
45     const polyBoundaryMesh& pbMesh = mesh.boundaryMesh();
47     const scalar deltaT = td.cloud().cachedDeltaT();
48     scalar tEnd = (1.0 - p.stepFraction())*deltaT;
49     const scalar dtMax = tEnd;
51     while (td.keepParticle && !td.switchProcessor && tEnd > ROOTVSMALL)
52     {
53         // Set the Lagrangian time-step
54         scalar dt = min(dtMax, tEnd);
56         dt *= p.trackToFace(p.position() + dt*U_, td);
58         tEnd -= dt;
60         p.stepFraction() = 1.0 - tEnd/deltaT;
62         if (p.onBoundary() && td.keepParticle)
63         {
64             if (isType<processorPolyPatch>(pbMesh[p.patch(p.face())]))
65             {
66                 td.switchProcessor = true;
67             }
68         }
69     }
71     return td.keepParticle;
75 template<class ParcelType>
76 template<class TrackData>
77 bool Foam::DsmcParcel<ParcelType>::hitPatch
79     const polyPatch&,
80     TrackData& td,
81     const label patchI
84     return false;
88 template<class ParcelType>
89 template<class TrackData>
90 void Foam::DsmcParcel<ParcelType>::hitProcessorPatch
92     const processorPolyPatch&,
93     TrackData& td
96     td.switchProcessor = true;
100 template<class ParcelType>
101 void Foam::DsmcParcel<ParcelType>::hitProcessorPatch
103     const processorPolyPatch&,
104     int&
109 template<class ParcelType>
110 template<class TrackData>
111 void Foam::DsmcParcel<ParcelType>::hitWallPatch
113     const wallPolyPatch& wpp,
114     TrackData& td
117     const constantProperties& constProps(td.cloud().constProps(typeId_));
119     scalar m = constProps.mass();
121     // pre-interaction energy
122     scalar preIE = 0.5*m*(U_ & U_) + Ei_;
124     // pre-interaction momentum
125     vector preIMom = m*U_;
127     td.cloud().wallInteraction().correct
128     (
129         wpp,
130         this->face(),
131         U_,
132         Ei_,
133         typeId_
134     );
136     // post-interaction energy
137     scalar postIE = 0.5*m*(U_ & U_) + Ei_;
139     // post-interaction momentum
140     vector postIMom = m*U_;
142     label wppIndex = wpp.index();
144     label wppLocalFace = wpp.whichFace(this->face());
146     const scalar fA = mag(wpp.faceAreas()[wppLocalFace]);
148     const scalar deltaT = td.cloud().cachedDeltaT();
150     scalar deltaQ = td.cloud().nParticle()*(preIE - postIE)/(deltaT*fA);
152     vector deltaFD = td.cloud().nParticle()*(preIMom - postIMom)/(deltaT*fA);
154     td.cloud().q().boundaryField()[wppIndex][wppLocalFace] += deltaQ;
156     td.cloud().fD().boundaryField()[wppIndex][wppLocalFace] += deltaFD;
160 template<class ParcelType>
161 void Foam::DsmcParcel<ParcelType>::hitWallPatch
163     const wallPolyPatch&,
164     int&
169 template<class ParcelType>
170 template<class TrackData>
171 void Foam::DsmcParcel<ParcelType>::hitPatch
173     const polyPatch&,
174     TrackData& td
177     td.keepParticle = false;
181 template<class ParcelType>
182 void Foam::DsmcParcel<ParcelType>::hitPatch
184     const polyPatch&,
185     int&
190 template<class ParcelType>
191 void Foam::DsmcParcel<ParcelType>::transformProperties
193     const tensor& T
196     Particle<ParcelType>::transformProperties(T);
197     U_ = transform(T, U_);
201 template<class ParcelType>
202 void Foam::DsmcParcel<ParcelType>::transformProperties
204     const vector& separation
207     Particle<ParcelType>::transformProperties(separation);
211 // * * * * * * * * * * * * * * * *  IOStream operators * * * * * * * * * * * //
213 #include "DsmcParcelIO.C"
215 // ************************************************************************* //