ParticleI.H 10.3 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
6
7
8
9
10
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

11
12
13
14
    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.
15
16
17
18
19
20
21

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    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 "polyMesh.H"
27
#include "wallPolyPatch.H"
28
29
30
31

// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //

template<class ParticleType>
32
inline Foam::scalar Foam::Particle<ParticleType>::lambda
33
34
35
36
37
(
    const vector& from,
    const vector& to,
    const label facei,
    const scalar stepFraction
38
) const
39
40
41
42
43
44
45
46
47
48
{
    const polyMesh& mesh = cloud_.polyMesh_;
    bool movingMesh = mesh.moving();

    if (movingMesh)
    {
        vector Sf = mesh.faceAreas()[facei];
        Sf /= mag(Sf);
        vector Cf = mesh.faceCentres()[facei];

49
        // patch interaction
50
        if (!cloud_.internalFace(facei))
51
        {
52
53
54
55
56
            label patchi = patch(facei);
            const polyPatch& patch = mesh.boundaryMesh()[patchi];

            // move reference point for wall
            if (isA<wallPolyPatch>(patch))
57
            {
58
59
60
61
62
63
64
65
66
67
                const vector& C = mesh.cellCentres()[celli_];
                scalar CCf = mag((C - Cf) & Sf);
                // check if distance between cell centre and face centre
                // is larger than wallImpactDistance
                const ParticleType& p =
                    static_cast<const ParticleType&>(*this);
                if (CCf > p.wallImpactDistance(Sf))
                {
                    Cf -=p.wallImpactDistance(Sf)*Sf;
                }
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
            }
        }

        // for a moving mesh we need to reconstruct the old
        // Sf and Cf from oldPoints (they aren't stored)

        const vectorField& oldPoints = mesh.oldPoints();

        vector Cf00 = mesh.faces()[facei].centre(oldPoints);
        vector Cf0 = Cf00 + stepFraction*(Cf - Cf00);

        vector Sf00 = mesh.faces()[facei].normal(oldPoints);

        // for layer addition Sf00 = vector::zero and we use Sf
        if (mag(Sf00) > SMALL)
        {
            Sf00 /= mag(Sf00);
        }
        else
        {
            Sf00 = Sf;
        }

        scalar magSfDiff = mag(Sf - Sf00);

        // check if the face is rotating
        if (magSfDiff > SMALL)
        {
            vector Sf0 = Sf00 + stepFraction*(Sf - Sf00);

            // find center of rotation
            vector omega = Sf0 ^ Sf;
            scalar magOmega = mag(omega);
101
            omega /= magOmega + SMALL;
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
            vector n0 = omega ^ Sf0;
            scalar lam = ((Cf - Cf0) & Sf)/(n0 & Sf);
            vector r0 = Cf0 + lam*n0;

            // solve '(p - r0) & Sfp = 0', where
            // p = from + lambda*(to - from)
            // Sfp = Sf0 + lambda*(Sf - Sf0)
            // which results in the quadratic eq.
            // a*lambda^2 + b*lambda + c = 0
            vector alpha = from - r0;
            vector beta = to - from;
            scalar a = beta & (Sf - Sf0);
            scalar b = (alpha & (Sf - Sf0)) + (beta & Sf0);
            scalar c = alpha & Sf0;

            if (mag(a) > SMALL)
            {
                // solve the second order polynomial
                scalar ap = b/a;
                scalar bp = c/a;
                scalar cp = ap*ap - 4.0*bp;

                if (cp < 0.0)
                {
                    // imaginary roots only
                    return GREAT;
                }
                else
                {
                    scalar l1 = -0.5*(ap - ::sqrt(cp));
                    scalar l2 = -0.5*(ap + ::sqrt(cp));

                    // one root is around 0-1, while
                    // the other is very large in mag
                    if (mag(l1) < mag(l2))
                    {
                        return l1;
                    }
                    else
                    {
                        return l2;
                    }
                }
            }
            else
            {
                // when a==0, solve the first order polynomial
                return (-c/b);
            }
        }
        else // no rotation
        {
            vector alpha = from - Cf0;
            vector beta = to - from - (Cf - Cf0);
            scalar lambdaNominator = alpha & Sf;
            scalar lambdaDenominator = beta & Sf;

            // check if trajectory is parallel to face
            if (mag(lambdaDenominator) < SMALL)
            {
                if (lambdaDenominator < 0.0)
                {
                    lambdaDenominator = -SMALL;
                }
                else
                {
                    lambdaDenominator = SMALL;
                }
            }

            return (-lambdaNominator/lambdaDenominator);
        }
    }
    else
    {
        // mesh is static and stepFraction is not needed
        return lambda(from, to, facei);
    }
}


template<class ParticleType>
184
inline Foam::scalar Foam::Particle<ParticleType>::lambda
185
186
187
188
(
    const vector& from,
    const vector& to,
    const label facei
189
) const
190
191
192
193
194
195
196
{
    const polyMesh& mesh = cloud_.polyMesh_;

    vector Sf = mesh.faceAreas()[facei];
    Sf /= mag(Sf);
    vector Cf = mesh.faceCentres()[facei];

197
    // patch interaction
198
    if (!cloud_.internalFace(facei))
199
    {
200
        label patchi = patch(facei);
201
202
        const polyPatch& patch = mesh.boundaryMesh()[patchi];

203
        // move reference point for wall
204
        if (isA<wallPolyPatch>(patch))
205
        {
206
207
208
209
            const vector& C = mesh.cellCentres()[celli_];
            scalar CCf = mag((C - Cf) & Sf);
            // check if distance between cell centre and face centre
            // is larger than wallImpactDistance
210
211
            const ParticleType& p = static_cast<const ParticleType&>(*this);
            if (CCf > p.wallImpactDistance(Sf))
212
            {
213
                Cf -=p.wallImpactDistance(Sf)*Sf;
214
            }
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
        }
    }

    scalar lambdaNominator = (Cf - from) & Sf;
    scalar lambdaDenominator = (to - from) & Sf;

    // check if trajectory is parallel to face
    if (mag(lambdaDenominator) < SMALL)
    {
        if (lambdaDenominator < 0.0)
        {
            lambdaDenominator = -SMALL;
        }
        else
        {
            lambdaDenominator = SMALL;
        }
    }

    return lambdaNominator/lambdaDenominator;
}


template<class ParticleType>
239
inline bool Foam::Particle<ParticleType>::inCell() const
240
{
mattijs's avatar
mattijs committed
241
242
    DynamicList<label>& faces = cloud_.labels_;
    findFaces(position_, faces);
243
244
245
246
247
248

    return (!faces.size());
}


template<class ParticleType>
249
inline bool Foam::Particle<ParticleType>::inCell
250
251
252
253
(
    const vector& position,
    const label celli,
    const scalar stepFraction
254
) const
255
{
mattijs's avatar
mattijs committed
256
257
    DynamicList<label>& faces = cloud_.labels_;
    findFaces(position, celli, stepFraction, faces);
258
259
260
261
262
263
264
265

    return (!faces.size());
}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

template<class ParticleType>
266
inline Foam::Particle<ParticleType>::trackData::trackData
267
268
269
270
271
272
273
(
    Cloud<ParticleType>& cloud
)
:
    cloud_(cloud)
{}

274

275
template<class ParticleType>
276
277
inline Foam::Cloud<ParticleType>&
Foam::Particle<ParticleType>::trackData::cloud()
278
279
280
281
282
283
284
285
{
    return cloud_;
}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

template<class ParticleType>
286
287
inline const Foam::Cloud<ParticleType>&
Foam::Particle<ParticleType>::cloud() const
288
289
290
291
292
293
{
    return cloud_;
}


template<class ParticleType>
294
inline const Foam::vector& Foam::Particle<ParticleType>::position() const
295
296
297
298
299
300
{
    return position_;
}


template<class ParticleType>
301
inline Foam::vector& Foam::Particle<ParticleType>::position()
302
303
304
305
306
307
{
    return position_;
}


template<class ParticleType>
308
inline Foam::label Foam::Particle<ParticleType>::cell() const
309
310
311
312
{
    return celli_;
}

313

314
template<class ParticleType>
315
inline Foam::label& Foam::Particle<ParticleType>::cell()
316
317
318
319
{
    return celli_;
}

320
321

template<class ParticleType>
322
inline Foam::label Foam::Particle<ParticleType>::face() const
323
324
325
326
327
328
{
    return facei_;
}


template<class ParticleType>
329
inline bool Foam::Particle<ParticleType>::onBoundary() const
330
331
332
333
334
335
{
    return facei_ != -1 && facei_ >= cloud_.pMesh().nInternalFaces();
}


template<class ParticleType>
336
inline Foam::scalar& Foam::Particle<ParticleType>::stepFraction()
337
338
339
340
341
342
{
    return stepFraction_;
}


template<class ParticleType>
343
inline Foam::scalar Foam::Particle<ParticleType>::stepFraction() const
344
345
346
347
348
349
{
    return stepFraction_;
}


template<class ParticleType>
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
inline Foam::label Foam::Particle<ParticleType>::origProc() const
{
    return origProc_;
}


template<class ParticleType>
inline Foam::label Foam::Particle<ParticleType>::origId() const
{
    return origId_;
}


template<class ParticleType>
inline bool Foam::Particle<ParticleType>::softImpact() const
365
366
367
368
369
{
    return false;
}


370
template<class ParticleType>
371
inline Foam::scalar Foam::Particle<ParticleType>::currentTime() const
372
373
374
{
    return
        cloud_.pMesh().time().value()
graham's avatar
graham committed
375
      + stepFraction_*cloud_.pMesh().time().deltaTValue();
376
377
378
}


379
template<class ParticleType>
380
inline Foam::label Foam::Particle<ParticleType>::patch(const label facei) const
381
382
383
384
385
386
{
    return cloud_.facePatch(facei);
}


template<class ParticleType>
387
inline Foam::label Foam::Particle<ParticleType>::patchFace
388
389
390
391
392
393
394
395
396
397
(
    const label patchi,
    const label facei
) const
{
    return cloud_.patchFace(patchi, facei);
}


template<class ParticleType>
398
399
inline Foam::scalar
Foam::Particle<ParticleType>::wallImpactDistance(const vector&) const
400
401
402
403
404
405
{
    return 0.0;
}


template<class ParticleType>
406
inline Foam::label Foam::Particle<ParticleType>::faceInterpolation() const
407
408
409
410
411
412
{
    return facei_;
}


// ************************************************************************* //