DSMCParcel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "DSMCParcel.H"
27 #include "meshTools.H"
28 
29 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30 
31 template<class ParcelType>
32 template<class TrackData>
33 bool Foam::DSMCParcel<ParcelType>::move(TrackData& td, const scalar trackTime)
34 {
35  typename TrackData::cloudType::parcelType& p =
36  static_cast<typename TrackData::cloudType::parcelType&>(*this);
37 
38  td.switchProcessor = false;
39  td.keepParticle = true;
40 
41  const polyMesh& mesh = td.cloud().pMesh();
42  const polyBoundaryMesh& pbMesh = mesh.boundaryMesh();
43 
44  scalar tEnd = (1.0 - p.stepFraction())*trackTime;
45  const scalar dtMax = tEnd;
46 
47  // For reduced-D cases, the velocity used to track needs to be
48  // constrained, but the actual U_ of the parcel must not be
49  // altered or used, as it is altered by patch interactions an
50  // needs to retain its 3D value for collision purposes.
51  vector Utracking = U_;
52 
53  while (td.keepParticle && !td.switchProcessor && tEnd > ROOTVSMALL)
54  {
55  // Apply correction to position for reduced-D cases
56  meshTools::constrainToMeshCentre(mesh, p.position());
57 
58  Utracking = U_;
59 
60  // Apply correction to velocity to constrain tracking for
61  // reduced-D cases
62  meshTools::constrainDirection(mesh, mesh.solutionD(), Utracking);
63 
64  // Set the Lagrangian time-step
65  scalar dt = min(dtMax, tEnd);
66 
67  dt *= p.trackToFace(p.position() + dt*Utracking, td);
68 
69  tEnd -= dt;
70 
71  p.stepFraction() = 1.0 - tEnd/trackTime;
72 
73  if (p.onBoundary() && td.keepParticle)
74  {
75  if (isA<processorPolyPatch>(pbMesh[p.patch(p.face())]))
76  {
77  td.switchProcessor = true;
78  }
79  }
80  }
81 
82  return td.keepParticle;
83 }
84 
85 
86 template<class ParcelType>
87 template<class TrackData>
89 (
90  const polyPatch&,
91  TrackData& td,
92  const label,
93  const scalar,
94  const tetIndices&
95 )
96 {
97  return false;
98 }
99 
100 
101 template<class ParcelType>
102 template<class TrackData>
104 (
105  const processorPolyPatch&,
106  TrackData& td
107 )
108 {
109  td.switchProcessor = true;
110 }
111 
112 
113 template<class ParcelType>
114 template<class TrackData>
116 (
117  const wallPolyPatch& wpp,
118  TrackData& td,
119  const tetIndices& tetIs
120 )
121 {
122  label wppIndex = wpp.index();
123 
124  label wppLocalFace = wpp.whichFace(this->face());
125 
126  const scalar fA = mag(wpp.faceAreas()[wppLocalFace]);
127 
128  const scalar deltaT = td.cloud().pMesh().time().deltaTValue();
129 
130  const constantProperties& constProps(td.cloud().constProps(typeId_));
131 
132  scalar m = constProps.mass();
133 
134  vector nw = wpp.faceAreas()[wppLocalFace];
135  nw /= mag(nw);
136 
137  scalar U_dot_nw = U_ & nw;
138 
139  vector Ut = U_ - U_dot_nw*nw;
140 
141  scalar invMagUnfA = 1/max(mag(U_dot_nw)*fA, VSMALL);
142 
143  td.cloud().rhoNBF()[wppIndex][wppLocalFace] += invMagUnfA;
144 
145  td.cloud().rhoMBF()[wppIndex][wppLocalFace] += m*invMagUnfA;
146 
147  td.cloud().linearKEBF()[wppIndex][wppLocalFace] +=
148  0.5*m*(U_ & U_)*invMagUnfA;
149 
150  td.cloud().internalEBF()[wppIndex][wppLocalFace] += Ei_*invMagUnfA;
151 
152  td.cloud().iDofBF()[wppIndex][wppLocalFace] +=
153  constProps.internalDegreesOfFreedom()*invMagUnfA;
154 
155  td.cloud().momentumBF()[wppIndex][wppLocalFace] += m*Ut*invMagUnfA;
156 
157  // pre-interaction energy
158  scalar preIE = 0.5*m*(U_ & U_) + Ei_;
159 
160  // pre-interaction momentum
161  vector preIMom = m*U_;
162 
163  td.cloud().wallInteraction().correct
164  (
165  static_cast<DSMCParcel<ParcelType> &>(*this),
166  wpp
167  );
168 
169  U_dot_nw = U_ & nw;
170 
171  Ut = U_ - U_dot_nw*nw;
172 
173  invMagUnfA = 1/max(mag(U_dot_nw)*fA, VSMALL);
174 
175  td.cloud().rhoNBF()[wppIndex][wppLocalFace] += invMagUnfA;
176 
177  td.cloud().rhoMBF()[wppIndex][wppLocalFace] += m*invMagUnfA;
178 
179  td.cloud().linearKEBF()[wppIndex][wppLocalFace] +=
180  0.5*m*(U_ & U_)*invMagUnfA;
181 
182  td.cloud().internalEBF()[wppIndex][wppLocalFace] += Ei_*invMagUnfA;
183 
184  td.cloud().iDofBF()[wppIndex][wppLocalFace] +=
185  constProps.internalDegreesOfFreedom()*invMagUnfA;
186 
187  td.cloud().momentumBF()[wppIndex][wppLocalFace] += m*Ut*invMagUnfA;
188 
189  // post-interaction energy
190  scalar postIE = 0.5*m*(U_ & U_) + Ei_;
191 
192  // post-interaction momentum
193  vector postIMom = m*U_;
194 
195  scalar deltaQ = td.cloud().nParticle()*(preIE - postIE)/(deltaT*fA);
196 
197  vector deltaFD = td.cloud().nParticle()*(preIMom - postIMom)/(deltaT*fA);
198 
199  td.cloud().qBF()[wppIndex][wppLocalFace] += deltaQ;
200 
201  td.cloud().fDBF()[wppIndex][wppLocalFace] += deltaFD;
202 
203 }
204 
205 
206 template<class ParcelType>
207 template<class TrackData>
209 {
210  td.keepParticle = false;
211 }
212 
213 
214 template<class ParcelType>
216 {
217  ParcelType::transformProperties(T);
218  U_ = transform(T, U_);
219 }
220 
221 
222 template<class ParcelType>
224 (
225  const vector& separation
226 )
227 {
228  ParcelType::transformProperties(separation);
229 }
230 
231 
232 // * * * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * //
233 
234 #include "DSMCParcelIO.C"
235 
236 // ************************************************************************* //
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:290
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:377
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:801
void hitWallPatch(const wallPolyPatch &, TrackData &td, const tetIndices &)
Overridable function to handle the particle hitting a wallPatch.
Definition: DSMCParcel.C:116
bool move(TrackData &td, const scalar trackTime)
Move the parcel.
Definition: DSMCParcel.C:33
DSMC parcel class.
Definition: DSMCParcel.H:51
void hitProcessorPatch(const processorPolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a.
Definition: DSMCParcel.C:104
bool hitPatch(const polyPatch &, TrackData &td, const label patchi, const scalar trackFraction, const tetIndices &tetIs)
Overridable function to handle the particle hitting a patch.
Definition: DSMCParcel.C:89
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Neighbour processor patch.
dynamicFvMesh & mesh
Class to hold DSMC particle constant properties.
Definition: DSMCParcel.H:80
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
Foam::polyBoundaryMesh.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const volScalarField & T
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: DSMCParcel.C:215
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar mass() const
Return const access to the particle mass [kg].
Definition: DSMCParcelI.H:78
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition: meshTools.C:671
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
void constrainToMeshCentre(const polyMesh &mesh, point &pt)
Set the constrained components of position to mesh centre.
Definition: meshTools.C:613
label index() const
Return the index of this patch in the boundaryMesh.
label nw
Definition: createFields.H:25
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465