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-2017 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  // For reduced-D cases, the velocity used to track needs to be
45  // constrained, but the actual U_ of the parcel must not be
46  // altered or used, as it is altered by patch interactions an
47  // needs to retain its 3D value for collision purposes.
48  vector Utracking = U_;
49 
50  while (td.keepParticle && !td.switchProcessor && p.stepFraction() < 1)
51  {
52  // Apply correction to position for reduced-D cases
53  p.constrainToMeshCentre();
54 
55  Utracking = U_;
56 
57  // Apply correction to velocity to constrain tracking for
58  // reduced-D cases
59  meshTools::constrainDirection(mesh, mesh.solutionD(), Utracking);
60 
61  const scalar f = 1 - p.stepFraction();
62  p.trackToFace(f*trackTime*Utracking, f, td);
63 
64  if (p.onBoundaryFace() && td.keepParticle)
65  {
66  if (isA<processorPolyPatch>(pbMesh[p.patch()]))
67  {
68  td.switchProcessor = true;
69  }
70  }
71  }
72 
73  return td.keepParticle;
74 }
75 
76 
77 template<class ParcelType>
78 template<class TrackData>
80 (
81  const polyPatch&,
82  TrackData& td,
83  const label,
84  const scalar,
85  const tetIndices&
86 )
87 {
88  return false;
89 }
90 
91 
92 template<class ParcelType>
93 template<class TrackData>
95 (
96  const processorPolyPatch&,
97  TrackData& td
98 )
99 {
100  td.switchProcessor = true;
101 }
102 
103 
104 template<class ParcelType>
105 template<class TrackData>
107 (
108  const wallPolyPatch& wpp,
109  TrackData& td,
110  const tetIndices& tetIs
111 )
112 {
113  label wppIndex = wpp.index();
114 
115  label wppLocalFace = wpp.whichFace(this->face());
116 
117  const scalar fA = mag(wpp.faceAreas()[wppLocalFace]);
118 
119  const scalar deltaT = td.cloud().pMesh().time().deltaTValue();
120 
121  const constantProperties& constProps(td.cloud().constProps(typeId_));
122 
123  scalar m = constProps.mass();
124 
125  vector nw = wpp.faceAreas()[wppLocalFace];
126  nw /= mag(nw);
127 
128  scalar U_dot_nw = U_ & nw;
129 
130  vector Ut = U_ - U_dot_nw*nw;
131 
132  scalar invMagUnfA = 1/max(mag(U_dot_nw)*fA, VSMALL);
133 
134  td.cloud().rhoNBF()[wppIndex][wppLocalFace] += invMagUnfA;
135 
136  td.cloud().rhoMBF()[wppIndex][wppLocalFace] += m*invMagUnfA;
137 
138  td.cloud().linearKEBF()[wppIndex][wppLocalFace] +=
139  0.5*m*(U_ & U_)*invMagUnfA;
140 
141  td.cloud().internalEBF()[wppIndex][wppLocalFace] += Ei_*invMagUnfA;
142 
143  td.cloud().iDofBF()[wppIndex][wppLocalFace] +=
144  constProps.internalDegreesOfFreedom()*invMagUnfA;
145 
146  td.cloud().momentumBF()[wppIndex][wppLocalFace] += m*Ut*invMagUnfA;
147 
148  // pre-interaction energy
149  scalar preIE = 0.5*m*(U_ & U_) + Ei_;
150 
151  // pre-interaction momentum
152  vector preIMom = m*U_;
153 
154  td.cloud().wallInteraction().correct
155  (
156  static_cast<DSMCParcel<ParcelType> &>(*this),
157  wpp
158  );
159 
160  U_dot_nw = U_ & nw;
161 
162  Ut = U_ - U_dot_nw*nw;
163 
164  invMagUnfA = 1/max(mag(U_dot_nw)*fA, VSMALL);
165 
166  td.cloud().rhoNBF()[wppIndex][wppLocalFace] += invMagUnfA;
167 
168  td.cloud().rhoMBF()[wppIndex][wppLocalFace] += m*invMagUnfA;
169 
170  td.cloud().linearKEBF()[wppIndex][wppLocalFace] +=
171  0.5*m*(U_ & U_)*invMagUnfA;
172 
173  td.cloud().internalEBF()[wppIndex][wppLocalFace] += Ei_*invMagUnfA;
174 
175  td.cloud().iDofBF()[wppIndex][wppLocalFace] +=
176  constProps.internalDegreesOfFreedom()*invMagUnfA;
177 
178  td.cloud().momentumBF()[wppIndex][wppLocalFace] += m*Ut*invMagUnfA;
179 
180  // post-interaction energy
181  scalar postIE = 0.5*m*(U_ & U_) + Ei_;
182 
183  // post-interaction momentum
184  vector postIMom = m*U_;
185 
186  scalar deltaQ = td.cloud().nParticle()*(preIE - postIE)/(deltaT*fA);
187 
188  vector deltaFD = td.cloud().nParticle()*(preIMom - postIMom)/(deltaT*fA);
189 
190  td.cloud().qBF()[wppIndex][wppLocalFace] += deltaQ;
191 
192  td.cloud().fDBF()[wppIndex][wppLocalFace] += deltaFD;
193 
194 }
195 
196 
197 template<class ParcelType>
198 template<class TrackData>
200 {
201  td.keepParticle = false;
202 }
203 
204 
205 template<class ParcelType>
207 {
208  ParcelType::transformProperties(T);
209  U_ = transform(T, U_);
210 }
211 
212 
213 template<class ParcelType>
215 (
216  const vector& separation
217 )
218 {
219  ParcelType::transformProperties(separation);
220 }
221 
222 
223 // * * * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * //
224 
225 #include "DSMCParcelIO.C"
226 
227 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
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
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void hitWallPatch(const wallPolyPatch &, TrackData &td, const tetIndices &)
Overridable function to handle the particle hitting a wallPatch.
Definition: DSMCParcel.C:107
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:824
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:95
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:80
Neighbour processor patch.
dynamicFvMesh & mesh
Class to hold DSMC particle constant properties.
Definition: DSMCParcel.H:80
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
scalar mass() const
Return const access to the particle mass [kg].
Definition: DSMCParcelI.H:96
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.
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:290
labelList f(nPoints)
const volScalarField & T
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: DSMCParcel.C:206
label index() const
Return the index of this patch in the boundaryMesh.
dimensioned< scalar > mag(const dimensioned< Type > &)
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
label nw
Definition: createFields.H:25
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:377