plane.C 11.3 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
6
     \\/     M anipulation  | Copyright (C) 2016-2017 OpenCFD Ltd.
7
8
9
10
-------------------------------------------------------------------------------
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
27
28
29
30
31
32
33
34

\*---------------------------------------------------------------------------*/

#include "plane.H"
#include "tensor.H"

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

void Foam::plane::calcPntAndVec(const scalarList& C)
{
    if (mag(C[0]) > VSMALL)
    {
35
        point_ = vector((-C[3]/C[0]), 0, 0);
36
37
38
39
40
    }
    else
    {
        if (mag(C[1]) > VSMALL)
        {
41
            point_ = vector(0, (-C[3]/C[1]), 0);
42
43
44
45
46
        }
        else
        {
            if (mag(C[2]) > VSMALL)
            {
47
                point_ = vector(0, 0, (-C[3]/C[2]));
48
49
50
            }
            else
            {
51
                FatalErrorInFunction
52
53
54
55
56
57
                    << "At least one plane coefficient must have a value"
                    << abort(FatalError);
            }
        }
    }

58
59
    normal_ = vector(C[0], C[1], C[2]);
    scalar magUnitVector(mag(normal_));
60
61
62

    if (magUnitVector < VSMALL)
    {
63
        FatalErrorInFunction
64
65
66
67
            << "Plane normal defined with zero length"
            << abort(FatalError);
    }

68
    normal_ /= magUnitVector;
69
70
71
72
73
74
75
76
77
78
}


void Foam::plane::calcPntAndVec
(
    const point& point1,
    const point& point2,
    const point& point3
)
{
79
    point_ = (point1 + point2 + point3)/3;
80
81
82
83
84
85
86
87
88
89
    vector line12 = point1 - point2;
    vector line23 = point2 - point3;

    if
    (
        mag(line12) < VSMALL
     || mag(line23) < VSMALL
     || mag(point3-point1) < VSMALL
    )
    {
90
91
        FatalErrorInFunction
            << "Bad points:" << point1 << ' ' << point2 << ' ' << point3
mattijs's avatar
mattijs committed
92
            << abort(FatalError);
93
94
    }

95
96
    normal_ = line12 ^ line23;
    scalar magUnitVector(mag(normal_));
97
98
99

    if (magUnitVector < VSMALL)
    {
100
101
        FatalErrorInFunction
            << "Plane normal defined with zero length" << nl
mattijs's avatar
mattijs committed
102
            << "Bad points:" << point1 << ' ' << point2 << ' ' << point3
103
104
105
            << abort(FatalError);
    }

106
    normal_ /= magUnitVector;
107
108
109
110
111
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

112
113
114
115
Foam::plane::plane()
{}


116
117
Foam::plane::plane(const vector& normalVector)
:
118
119
    normal_(normalVector),
    point_(Zero)
120
{
121
    scalar magUnitVector(mag(normal_));
122
123
124

    if (magUnitVector > VSMALL)
    {
125
        normal_ /= magUnitVector;
126
127
128
    }
    else
    {
129
        FatalErrorInFunction
130
            << "plane normal has zero length. basePoint:" << point_
131
132
133
            << abort(FatalError);
    }
}
134
135


136
137
138
139
140
141
Foam::plane::plane
(
    const point& basePoint,
    const vector& normalVector,
    const bool normalise
)
142
:
143
144
    normal_(normalVector),
    point_(basePoint)
145
{
146
    scalar magSqrNormalVector(magSqr(normal_));
147

148
    if (magSqrNormalVector < VSMALL)
149
    {
150
        FatalErrorInFunction
151
            << "plane normal has zero length. basePoint:" << point_
152
            << abort(FatalError);
153
    }
154
155
156

    if (normalise)
    {
157
        normal_ /= Foam::sqrt(magSqrNormalVector);
158
    }
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
}


Foam::plane::plane(const scalarList& C)
{
    calcPntAndVec(C);
}


Foam::plane::plane
(
    const point& a,
    const point& b,
    const point& c
)
{
    calcPntAndVec(a, b, c);
}


Foam::plane::plane(const dictionary& dict)
:
181
182
    normal_(Zero),
    point_(Zero)
183
{
184
    const word planeType(dict.lookup("planeType"));
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200

    if (planeType == "planeEquation")
    {
        const dictionary& subDict = dict.subDict("planeEquationDict");
        scalarList C(4);

        C[0] = readScalar(subDict.lookup("a"));
        C[1] = readScalar(subDict.lookup("b"));
        C[2] = readScalar(subDict.lookup("c"));
        C[3] = readScalar(subDict.lookup("d"));

        calcPntAndVec(C);

    }
    else if (planeType == "embeddedPoints")
    {
201
        const dictionary& subDict = dict.subDict("embeddedPointsDict");
202
203
204
205
206
207
208
209
210
211
212

        point point1(subDict.lookup("point1"));
        point point2(subDict.lookup("point2"));
        point point3(subDict.lookup("point3"));

        calcPntAndVec(point1, point2, point3);
    }
    else if (planeType == "pointAndNormal")
    {
        const dictionary& subDict = dict.subDict("pointAndNormalDict");

213
        point_ =
214
215
216
            subDict.found("basePoint")
          ? subDict.lookup("basePoint")
          : subDict.lookup("point");
217
218

        normal_ =
219
220
221
            subDict.found("normalVector")
          ? subDict.lookup("normalVector")
          : subDict.lookup("normal");
222
223

        normal_ /= mag(normal_);
224
225
226
    }
    else
    {
227
        FatalIOErrorInFunction(dict)
228
229
230
            << "Invalid plane type: " << planeType << nl
            << "Valid options include: planeEquation, embeddedPoints and "
            << "pointAndNormal"
231
            << abort(FatalIOError);
232
233
234
235
236
237
    }
}


Foam::plane::plane(Istream& is)
:
238
239
    normal_(is),
    point_(is)
240
{
241
    scalar magUnitVector(mag(normal_));
242
243
244

    if (magUnitVector > VSMALL)
    {
245
        normal_ /= magUnitVector;
246
247
248
    }
    else
    {
249
        FatalErrorInFunction
250
            << "plane normal has zero length. basePoint:" << point_
251
252
253
254
255
256
257
258
259
            << abort(FatalError);
    }
}


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

const Foam::vector& Foam::plane::normal() const
{
260
    return normal_;
261
262
263
264
265
}


const Foam::point& Foam::plane::refPoint() const
{
266
    return point_;
267
268
269
}


270
Foam::FixedList<Foam::scalar, 4> Foam::plane::planeCoeffs() const
271
{
272
    FixedList<scalar, 4> C(4);
273

274
275
276
    scalar magX = mag(normal_.x());
    scalar magY = mag(normal_.y());
    scalar magZ = mag(normal_.z());
277
278
279
280
281
282

    if (magX > magY)
    {
        if (magX > magZ)
        {
            C[0] = 1;
283
284
            C[1] = normal_.y()/normal_.x();
            C[2] = normal_.z()/normal_.x();
285
286
287
        }
        else
        {
288
289
            C[0] = normal_.x()/normal_.z();
            C[1] = normal_.y()/normal_.z();
290
291
292
293
294
295
296
            C[2] = 1;
        }
    }
    else
    {
        if (magY > magZ)
        {
297
            C[0] = normal_.x()/normal_.y();
298
            C[1] = 1;
299
            C[2] = normal_.z()/normal_.y();
300
301
302
        }
        else
        {
303
304
            C[0] = normal_.x()/normal_.z();
            C[1] = normal_.y()/normal_.z();
305
306
307
308
            C[2] = 1;
        }
    }

309
310
311
    C[3] = - C[0] * point_.x()
           - C[1] * point_.y()
           - C[2] * point_.z();
312
313
314
315
316
317
318

    return C;
}


Foam::point Foam::plane::nearestPoint(const point& p) const
{
319
    return p - normal_*((p - point_) & normal_);
320
321
322
323
324
}


Foam::scalar Foam::plane::distance(const point& p) const
{
325
    return mag((p - point_) & normal_);
326
327
328
329
330
331
332
333
334
}


Foam::scalar Foam::plane::normalIntersect
(
    const point& pnt0,
    const vector& dir
) const
{
335
    scalar denom = stabilise((dir & normal_), VSMALL);
336

337
    return ((point_ - pnt0) & normal_)/denom;
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
}


Foam::plane::ray Foam::plane::planeIntersect(const plane& plane2) const
{
    // Mathworld plane-plane intersection. Assume there is a point on the
    // intersection line with z=0 and solve the two plane equations
    // for that (now 2x2 equation in x and y)
    // Better: use either z=0 or x=0 or y=0.

    const vector& n1 = normal();
    const vector& n2 = plane2.normal();

    const point& p1 = refPoint();
    const point& p2 = plane2.refPoint();

    scalar n1p1 = n1&p1;
    scalar n2p2 = n2&p2;

    vector dir = n1 ^ n2;

    // Determine zeroed out direction (can be x,y or z) by looking at which
    // has the largest component in dir.
    scalar magX = mag(dir.x());
    scalar magY = mag(dir.y());
    scalar magZ = mag(dir.z());

    direction iZero, i1, i2;

    if (magX > magY)
    {
        if (magX > magZ)
        {
            iZero = 0;
            i1 = 1;
            i2 = 2;
        }
        else
        {
            iZero = 2;
            i1 = 0;
            i2 = 1;
        }
    }
    else
    {
        if (magY > magZ)
        {
            iZero = 1;
            i1 = 2;
            i2 = 0;
        }
        else
        {
            iZero = 2;
            i1 = 0;
            i2 = 1;
        }
    }

    vector pt;

    pt[iZero] = 0;
    pt[i1] = (n2[i2]*n1p1 - n1[i2]*n2p2) / (n1[i1]*n2[i2] - n2[i1]*n1[i2]);
    pt[i2] = (n2[i1]*n1p1 - n1[i1]*n2p2) / (n1[i2]*n2[i1] - n1[i1]*n2[i2]);

    return ray(pt, dir);
}


Foam::point Foam::plane::planePlaneIntersect
(
    const plane& plane2,
    const plane& plane3
) const
{
414
415
416
    FixedList<scalar, 4> coeffs1(planeCoeffs());
    FixedList<scalar, 4> coeffs2(plane2.planeCoeffs());
    FixedList<scalar, 4> coeffs3(plane3.planeCoeffs());
417
418
419

    tensor a
    (
420
421
422
        coeffs1[0],coeffs1[1],coeffs1[2],
        coeffs2[0],coeffs2[1],coeffs2[2],
        coeffs3[0],coeffs3[1],coeffs3[2]
423
424
    );

425
    vector b(coeffs1[3],coeffs2[3],coeffs3[3]);
426
427
428
429
430

    return (inv(a) & (-b));
}


431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
Foam::point Foam::plane::somePointInPlane(const scalar dist) const
{
    // ax + by + cz + d = 0
    const FixedList<scalar, 4> coeff(planeCoeffs());

    // Perturb the base-point
    point p = refPoint() + point::uniform(dist);

    if (Foam::mag(coeff[2]) < SMALL)
    {
        if (Foam::mag(coeff[1]) < SMALL)
        {
            p[0] = -1.0*(coeff[1]*p[1] + coeff[2]*p[2] + coeff[3])/coeff[0];
        }
        else
        {
            p[1] = -1.0*(coeff[0]*p[0] + coeff[2]*p[2] + coeff[3])/coeff[1];
        }
    }
    else
    {
        p[2] = -1.0*(coeff[0]*p[0] + coeff[1]*p[1] + coeff[3])/coeff[2];
    }

    return p;
}


459
460
Foam::plane::side Foam::plane::sideOfPlane(const point& p) const
{
461
    const scalar angle((p - point_) & normal_);
462
463
464
465
466

    return (angle < 0 ? FLIP : NORMAL);
}


467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
Foam::point Foam::plane::mirror(const point& p) const
{
    const vector mirroredPtDir = p - nearestPoint(p);

    if ((normal() & mirroredPtDir) > 0)
    {
        return p - 2.0*distance(p)*normal();
    }
    else
    {
        return p + 2.0*distance(p)*normal();
    }
}


482
483
void Foam::plane::writeDict(Ostream& os) const
{
484
    os.writeEntry("planeType", "pointAndNormal");
485

486
    os.beginBlock("pointAndNormalDict");
487

488
489
    os.writeEntry("point",  point_);
    os.writeEntry("normal", normal_);
490

491
    os.endBlock() << flush;
492
493
}

494
// * * * * * * * * * * * * * * * Global Operators  * * * * * * * * * * * * * //
495
496
497

bool Foam::operator==(const plane& a, const plane& b)
{
498
    return (a.refPoint() == b.refPoint() && a.normal() == b.normal());
499
500
}

501

502
503
504
505
506
507
508
509
bool Foam::operator!=(const plane& a, const plane& b)
{
    return !(a == b);
}


// * * * * * * * * * * * * * * * Friend Functions  * * * * * * * * * * * * * //

510
Foam::Ostream& Foam::operator<<(Ostream& os, const plane& pln)
511
{
512
    os  << pln.normal_ << token::SPACE << pln.point_;
513
514
515
516
517
518

    return os;
}


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