EllipticDomain.cpp 16.9 KB
Newer Older
1
//
frey_m's avatar
frey_m committed
2 3
// Class EllipticDomain
//   :FIXME: add brief description
4
//
frey_m's avatar
frey_m committed
5 6 7
// Copyright (c) 2010 - 2013, Yves Ineichen, ETH Zürich,
//               2013 - 2015, Tülin Kaman, Paul Scherrer Institut, Villigen PSI, Switzerland
// All rights reserved
8
//
frey_m's avatar
frey_m committed
9
// Implemented as part of the PhD thesis
frey_m's avatar
frey_m committed
10
// "Toward massively parallel multi-objective optimization with application to
frey_m's avatar
frey_m committed
11 12 13 14 15 16 17 18 19 20 21
// particle accelerators" (https://doi.org/10.3929/ethz-a-009792359)
//
// This file is part of OPAL.
//
// OPAL 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.
//
// You should have received a copy of the GNU General Public License
// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
22
//
kraus's avatar
kraus committed
23 24
#include "Solvers/EllipticDomain.h"

gsell's avatar
gsell committed
25 26 27
#include <map>
#include <cmath>
#include <iostream>
28
#include <cassert>
gsell's avatar
gsell committed
29 30 31 32 33 34 35 36 37

//FIXME: ORDER HOW TO TRAVERSE NODES IS FIXED, THIS SHOULD BE MORE GENERIC! (PLACES MARKED)


EllipticDomain::EllipticDomain(Vector_t nr, Vector_t hr) {
    setNr(nr);
    setHr(hr);
}

38 39 40 41 42
EllipticDomain::EllipticDomain(double semimajor, double semiminor, Vector_t nr,
                               Vector_t hr, std::string interpl)
    : semiMajor_m(semimajor)
    , semiMinor_m(semiminor)
{
gsell's avatar
gsell committed
43 44 45
    setNr(nr);
    setHr(hr);

46
    if(interpl == "CONSTANT")
47
        interpolationMethod_m = CONSTANT;
48
    else if(interpl == "LINEAR")
49
        interpolationMethod_m = LINEAR;
50
    else if(interpl == "QUADRATIC")
51
        interpolationMethod_m = QUADRATIC;
gsell's avatar
gsell committed
52 53
}

54 55 56 57
EllipticDomain::EllipticDomain(BoundaryGeometry *bgeom, Vector_t nr, Vector_t hr,
                               std::string interpl)
    : EllipticDomain(bgeom->getA(), bgeom->getB(), nr, hr, interpl)
{
58 59 60
    setMinMaxZ(bgeom->getS(), bgeom->getLength());
}

61 62 63 64 65 66 67
EllipticDomain::EllipticDomain(BoundaryGeometry *bgeom, Vector_t nr, std::string interpl)
    : EllipticDomain(bgeom, nr,
                     Vector_t(2.0 * bgeom->getA() / nr[0],
                              2.0 * bgeom->getB() / nr[1],
                              (bgeom->getLength() - bgeom->getS()) / nr[2]),
                     interpl)
{ }
68

gsell's avatar
gsell committed
69 70 71 72 73 74 75 76
EllipticDomain::~EllipticDomain() {
    //nothing so far
}


// for this geometry we only have to calculate the intersection on
// one x-y-plane
// for the moment we center the ellipse around the center of the grid
77
void EllipticDomain::compute(Vector_t hr, NDIndex<3> localId){
78
    //there is nothing to be done if the mesh spacings have not changed
79 80 81 82
    if (hr[0] == getHr()[0] &&
        hr[1] == getHr()[1] &&
        hr[2] == getHr()[2])
    {
83 84 85
        hasGeometryChanged_m = false;
        return;
    }
86

87 88 89
    setHr(hr);
    hasGeometryChanged_m = true;
    //reset number of points inside domain
90
    nxy_m = 0; //nr[0] * nr[1];
91 92

    // clear previous coordinate maps
93 94
    idxMap_m.clear();
    coordMap_m.clear();
95
    //clear previous intersection points
96 97
    intersectYDir_m.clear();
    intersectXDir_m.clear();
98 99

    // build a index and coordinate map
100 101
    int idx = 0;
    int x, y;
102

103 104 105 106 107 108 109 110 111 112
    /* we need to scan the full x-y-plane on all cores
     * in order to figure out the number of valid
     * grid points per plane --> otherwise we might
     * get not unique global indices in the Tpetra::CrsMatrix
     */
    for (x = 0; x < nr[0]; ++x) {
        for (y = 0; y < nr[1]; ++y) {
            if (isInside(x, y, 1)) {
                idxMap_m[toCoordIdx(x, y)] = idx;
                coordMap_m[idx++] = toCoordIdx(x, y);
113 114 115 116 117
                nxy_m++;
            }
        }
    }

118
    switch (interpolationMethod_m) {
119 120 121 122 123
        case CONSTANT:
            break;
        case LINEAR:
        case QUADRATIC:

124 125
            double smajsq = semiMajor_m * semiMajor_m;
            double sminsq = semiMinor_m * semiMinor_m;
126 127 128 129 130 131 132 133 134
            double yd = 0.0;
            double xd = 0.0;
            double pos = 0.0;
            double mx = (nr[0] - 1) * hr[0] / 2.0;
            double my = (nr[1] - 1) * hr[1] / 2.0;

            //calculate intersection with the ellipse
            for(x = localId[0].first(); x <= localId[0].last(); x++) {
                pos = x * hr[0] - mx;
135
                if (std::abs(pos) >= semiMajor_m)
136
                {
137 138
                    intersectYDir_m.insert(std::pair<int, double>(x, 0));
                    intersectYDir_m.insert(std::pair<int, double>(x, 0));
kraus's avatar
kraus committed
139
                }else{
140 141 142
                    yd = std::abs(std::sqrt(sminsq - sminsq * pos * pos / smajsq)); // + 0.5*nr[1]*hr[1]);
                    intersectYDir_m.insert(std::pair<int, double>(x, yd));
                    intersectYDir_m.insert(std::pair<int, double>(x, -yd));
kraus's avatar
kraus committed
143
                }
144 145 146 147 148

            }

            for(y = localId[0].first(); y < localId[1].last(); y++) {
                pos = y * hr[1] - my;
149
                if (std::abs(pos) >= semiMinor_m)
150
                {
151 152
                    intersectXDir_m.insert(std::pair<int, double>(y, 0));
                    intersectXDir_m.insert(std::pair<int, double>(y, 0));
kraus's avatar
kraus committed
153
                }else{
154 155 156
                    xd = std::abs(std::sqrt(smajsq - smajsq * pos * pos / sminsq)); // + 0.5*nr[0]*hr[0]);
                    intersectXDir_m.insert(std::pair<int, double>(y, xd));
                    intersectXDir_m.insert(std::pair<int, double>(y, -xd));
kraus's avatar
kraus committed
157
                }
158 159
            }
    }
160 161
}

162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179
void EllipticDomain::resizeMesh(Vector_t& origin, Vector_t& hr, PartBunchBase<double, 3>* bunch) {
    Vector_t rmin  = bunch->getLowerBound();
    Vector_t rmax  = bunch->getUpperBound();
    Vector_t mymax = Vector_t(0.0, 0.0, 0.0);

    double dh = bunch->getMeshEnlargement();

    // apply bounding box increment, i.e., "BBOXINCR" input argument
    double zmin = std::signbit(rmin[2]) ? rmin[2] * (1.0 + dh) : rmin[2] * (1.0 - dh);
    double zmax = std::signbit(rmax[2]) ? rmax[2] * (1.0 - dh) : rmax[2] * (1.0 + dh);

    origin = Vector_t(-semiMajor_m, -semiMinor_m, zmin);
    mymax  = Vector_t( semiMajor_m,  semiMinor_m, zmax);

    for (int i = 0; i < 3; ++i)
        hr[i]   = (mymax[i] - origin[i]) / nr[i];
}

180 181 182 183 184
void EllipticDomain::getBoundaryStencil(int x, int y, int z, double &W,
                                        double &E, double &S, double &N,
                                        double &F, double &B, double &C,
                                        double &scaleFactor)
{
185
    scaleFactor = 1.0;
gsell's avatar
gsell committed
186

187
        // determine which interpolation method we use for points near the boundary
188
    switch (interpolationMethod_m) {
gsell's avatar
gsell committed
189
        case CONSTANT:
190
            constantInterpolation(x, y, z, W, E, S, N, F, B, C, scaleFactor);
gsell's avatar
gsell committed
191 192
            break;
        case LINEAR:
193
            linearInterpolation(x, y, z, W, E, S, N, F, B, C, scaleFactor);
gsell's avatar
gsell committed
194 195
            break;
        case QUADRATIC:
196
            quadraticInterpolation(x, y, z, W, E, S, N, F, B, C, scaleFactor);
gsell's avatar
gsell committed
197 198 199 200 201 202 203
            break;
    }

    // stencil center value has to be positive!
    assert(C > 0);
}

204 205 206 207
void EllipticDomain::getBoundaryStencil(int idx, double &W, double &E,
                                        double &S, double &N, double &F,
                                        double &B, double &C, double &scaleFactor)
{
gsell's avatar
gsell committed
208 209 210 211 212 213 214 215 216 217 218 219 220 221
    int x = 0, y = 0, z = 0;
    getCoord(idx, x, y, z);
    getBoundaryStencil(x, y, z, W, E, S, N, F, B, C, scaleFactor);
}


void EllipticDomain::getNeighbours(int idx, int &W, int &E, int &S, int &N, int &F, int &B) {

    int x = 0, y = 0, z = 0;
    getCoord(idx, x, y, z);
    getNeighbours(x, y, z, W, E, S, N, F, B);

}

222 223 224 225
void EllipticDomain::getNeighbours(int x, int y, int z, int &W, int &E,
                                   int &S, int &N, int &F, int &B)
{
    if (x > 0)
gsell's avatar
gsell committed
226 227 228
        W = getIdx(x - 1, y, z);
    else
        W = -1;
229 230

    if (x < nr[0] - 1)
gsell's avatar
gsell committed
231 232 233 234
        E = getIdx(x + 1, y, z);
    else
        E = -1;

235
    if (y < nr[1] - 1)
gsell's avatar
gsell committed
236 237 238
        N = getIdx(x, y + 1, z);
    else
        N = -1;
239 240

    if (y > 0)
gsell's avatar
gsell committed
241 242 243 244
        S = getIdx(x, y - 1, z);
    else
        S = -1;

245
    if (z > 0)
gsell's avatar
gsell committed
246 247 248
        F = getIdx(x, y, z - 1);
    else
        F = -1;
249 250

    if (z < nr[2] - 1)
gsell's avatar
gsell committed
251 252 253 254 255 256
        B = getIdx(x, y, z + 1);
    else
        B = -1;

}

257 258 259 260 261
void EllipticDomain::constantInterpolation(int x, int y, int z, double &WV,
                                           double &EV, double &SV, double &NV,
                                           double &FV, double &BV, double &CV,
                                           double &scaleFactor)
{
gsell's avatar
gsell committed
262 263
    scaleFactor = 1.0;

264 265 266 267 268 269 270 271 272 273
    WV = -1.0 / (hr[0] * hr[0]);
    EV = -1.0 / (hr[0] * hr[0]);
    NV = -1.0 / (hr[1] * hr[1]);
    SV = -1.0 / (hr[1] * hr[1]);
    FV = -1.0 / (hr[2] * hr[2]);
    BV = -1.0 / (hr[2] * hr[2]);

    CV =  2.0 / (hr[0] * hr[0])
       +  2.0 / (hr[1] * hr[1])
       +  2.0 / (hr[2] * hr[2]);
gsell's avatar
gsell committed
274 275

    // we are a right boundary point
276 277
    if (!isInside(x + 1, y, z))
        EV = 0.0;
gsell's avatar
gsell committed
278 279

    // we are a left boundary point
280
    if (!isInside(x - 1, y, z))
281
        WV = 0.0;
gsell's avatar
gsell committed
282 283

    // we are a upper boundary point
284
    if (!isInside(x, y + 1, z))
285
        NV = 0.0;
gsell's avatar
gsell committed
286 287

    // we are a lower boundary point
288
    if (!isInside(x, y - 1, z))
289
        SV = 0.0;
gsell's avatar
gsell committed
290

291
    if (z == 0 || z == nr[2] - 1) {
gsell's avatar
gsell committed
292 293 294 295 296

        // case where we are on the Robin BC in Z-direction
        // where we distinguish two cases
        // IFF: this values should not matter because they
        // never make it into the discretization matrix
297
        if (z == 0)
298
            FV = 0.0;
gsell's avatar
gsell committed
299
        else
300
            BV = 0.0;
gsell's avatar
gsell committed
301 302 303

        // add contribution of Robin discretization to center point
        // d the distance between the center of the bunch and the boundary
304 305
        double d = hr[2] * (nr[2] - 1) / 2.0;
        CV += 2.0 / (d * hr[2]);
gsell's avatar
gsell committed
306 307 308 309
    }

}

310 311 312 313 314
void EllipticDomain::linearInterpolation(int x, int y, int z, double &W,
                                         double &E, double &S, double &N,
                                         double &F, double &B, double &C,
                                         double &scaleFactor)
{
gsell's avatar
gsell committed
315 316 317 318 319 320
    scaleFactor = 1.0;

    double cx = x * hr[0] - (nr[0] - 1) * hr[0] / 2.0;
    double cy = y * hr[1] - (nr[1] - 1) * hr[1] / 2.0;

    double dx = 0.0;
321
    std::multimap<int, double>::iterator it = intersectXDir_m.find(y);
gsell's avatar
gsell committed
322
    if(cx < 0)
323
        ++it;
gsell's avatar
gsell committed
324 325 326
    dx = it->second;

    double dy = 0.0;
327 328
    it = intersectYDir_m.find(x);
    if (cy < 0)
329
        ++it;
gsell's avatar
gsell committed
330 331 332 333 334 335 336 337 338 339
    dy = it->second;


    double dw = hr[0];
    double de = hr[0];
    double dn = hr[1];
    double ds = hr[1];
    C = 0.0;

    //we are a right boundary point
340
    if (!isInside(x + 1, y, z)) {
gsell's avatar
gsell committed
341 342 343 344 345 346 347 348
        C += 1 / ((dx - cx) * de);
        E = 0.0;
    } else {
        C += 1 / (de * de);
        E = -1 / (de * de);
    }

    //we are a left boundary point
349
    if (!isInside(x - 1, y, z)) {
gsell's avatar
gsell committed
350 351 352 353 354 355 356 357
        C += 1 / ((std::abs(dx) - std::abs(cx)) * dw);
        W = 0.0;
    } else {
        C += 1 / (dw * dw);
        W = -1 / (dw * dw);
    }

    //we are a upper boundary point
358
    if (!isInside(x, y + 1, z)) {
gsell's avatar
gsell committed
359 360 361 362 363 364 365 366
        C += 1 / ((dy - cy) * dn);
        N = 0.0;
    } else {
        C += 1 / (dn * dn);
        N = -1 / (dn * dn);
    }

    //we are a lower boundary point
367
    if (!isInside(x, y - 1, z)) {
gsell's avatar
gsell committed
368 369 370 371 372 373 374 375 376 377 378 379
        C += 1 / ((std::abs(dy) - std::abs(cy)) * ds);
        S = 0.0;
    } else {
        C += 1 / (ds * ds);
        S = -1 / (ds * ds);
    }

    F = -1 / (hr[2] * hr[2]);
    B = -1 / (hr[2] * hr[2]);
    C += 2 / (hr[2] * hr[2]);

    // handle boundary condition in z direction
380
/*
gsell's avatar
gsell committed
381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401
    if(z == 0 || z == nr[2] - 1) {

        // case where we are on the NEUMAN BC in Z-direction
        // where we distinguish two cases
        if(z == 0)
            F = 0.0;
        else
            B = 0.0;

        //hr[2]*(nr2[2]-1)/2 = radius
        double d = hr[2] * (nr[2] - 1) / 2;
        C += 2 / (d * hr[2]);

        W /= 2.0;
        E /= 2.0;
        N /= 2.0;
        S /= 2.0;
        C /= 2.0;
        scaleFactor *= 0.5;

    }
402
*/
gsell's avatar
gsell committed
403 404
}

405 406 407 408 409
void EllipticDomain::quadraticInterpolation(int x, int y, int z,
                                            double &W, double &E, double &S,
                                            double &N, double &F, double &B, double &C,
                                            double &scaleFactor)
{
gsell's avatar
gsell committed
410 411 412 413 414 415
    double cx = (x - (nr[0] - 1) / 2.0) * hr[0];
    double cy = (y - (nr[1] - 1) / 2.0) * hr[1];

    // since every vector for elliptic domains has ALWAYS size 2 we
    // can catch all cases manually
    double dx = 0.0;
416 417
    std::multimap<int, double>::iterator it = intersectXDir_m.find(y);
    if (cx < 0)
418
        ++it;
gsell's avatar
gsell committed
419 420 421
    dx = it->second;

    double dy = 0.0;
422 423
    it = intersectYDir_m.find(x);
    if (cy < 0)
424
        ++it;
gsell's avatar
gsell committed
425 426 427 428 429 430 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 459 460 461 462 463 464 465 466
    dy = it->second;

    double dw = hr[0];
    double de = hr[0];
    double dn = hr[1];
    double ds = hr[1];
    W = 1.0;
    E = 1.0;
    N = 1.0;
    S = 1.0;
    F = 1.0;
    B = 1.0;
    C = 0.0;

    //TODO: = cx+hr[0] > dx && cx > 0
    //if((x-nr[0]/2.0+1)*hr[0] > dx && cx > 0) {
    ////we are a right boundary point
    ////if(!isInside(x+1,y,z)) {
    //de = dx-cx;
    //}

    //if((x-nr[0]/2.0-1)*hr[0] < dx && cx < 0) {
    ////we are a left boundary point
    ////if(!isInside(x-1,y,z)) {
    //dw = std::abs(dx)-std::abs(cx);
    //}

    //if((y-nr[1]/2.0+1)*hr[1] > dy && cy > 0) {
    ////we are a upper boundary point
    ////if(!isInside(x,y+1,z)) {
    //dn = dy-cy;
    //}

    //if((y-nr[1]/2.0-1)*hr[1] < dy && cy < 0) {
    ////we are a lower boundary point
    ////if(!isInside(x,y-1,z)) {
    //ds = std::abs(dy)-std::abs(cy);
    //}

    //TODO: = cx+hr[0] > dx && cx > 0
    //if((x-nr[0]/2.0+1)*hr[0] > dx && cx > 0) {
    //we are a right boundary point
467
    if (!isInside(x + 1, y, z)) {
gsell's avatar
gsell committed
468 469 470 471 472 473
        de = dx - cx;
        E = 0.0;
    }

    //if((x-nr[0]/2.0-1)*hr[0] < dx && cx < 0) {
    //we are a left boundary point
474
    if (!isInside(x - 1, y, z)) {
gsell's avatar
gsell committed
475 476 477 478 479 480
        dw = std::abs(dx) - std::abs(cx);
        W = 0.0;
    }

    //if((y-nr[1]/2.0+1)*hr[1] > dy && cy > 0) {
    //we are a upper boundary point
481
    if (!isInside(x, y + 1, z)) {
gsell's avatar
gsell committed
482 483 484 485 486 487
        dn = dy - cy;
        N = 0.0;
    }

    //if((y-nr[1]/2.0-1)*hr[1] < dy && cy < 0) {
    //we are a lower boundary point
488
    if (!isInside(x, y - 1, z)) {
gsell's avatar
gsell committed
489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600
        ds = std::abs(dy) - std::abs(cy);
        S = 0.0;
    }

    //2/dw*(dw_de)
    W *= -1.0 / (dw * (dw + de));
    E *= -1.0 / (de * (dw + de));
    N *= -1.0 / (dn * (dn + ds));
    S *= -1.0 / (ds * (dn + ds));
    F = -1 / (hr[2] * (hr[2] + hr[2]));
    B = -1 / (hr[2] * (hr[2] + hr[2]));

    //TODO: problem when de,dw,dn,ds == 0
    //is NOT a regular BOUND PT
    C += 1 / de * 1 / dw;
    C += 1 / dn * 1 / ds;
    C += 1 / hr[2] * 1 / hr[2];
    scaleFactor = 0.5;


    //for regular gridpoints no problem with symmetry, just boundary
    //z direction is right
    //implement isLastInside(dir)
    //we have LOCAL x,y coordinates!

    /*
       if(dw != 0 && !wIsB)
       W = -1/dw * (dn+ds) * 2*hr[2];
       else
       W = 0;
       if(de != 0 && !eIsB)
       E = -1/de * (dn+ds) * 2*hr[2];
       else
       E = 0;
       if(dn != 0 && !nIsB)
       N = -1/dn * (dw+de) * 2*hr[2];
       else
       N = 0;
       if(ds != 0 && !sIsB)
       S = -1/ds * (dw+de) * 2*hr[2];
       else
       S = 0;
       F = -(dw+de)*(dn+ds)/hr[2];
       B = -(dw+de)*(dn+ds)/hr[2];
       */

    //if(dw != 0)
    //W = -2*hr[2]*(dn+ds)/dw;
    //else
    //W = 0;
    //if(de != 0)
    //E = -2*hr[2]*(dn+ds)/de;
    //else
    //E = 0;
    //if(dn != 0)
    //N = -2*hr[2]*(dw+de)/dn;
    //else
    //N = 0;
    //if(ds != 0)
    //S = -2*hr[2]*(dw+de)/ds;
    //else
    //S = 0;
    //F = -(dw+de)*(dn+ds)/hr[2];
    //B = -(dw+de)*(dn+ds)/hr[2];

    //// RHS scaleFactor for current 3D index
    //// Factor 0.5 results from the SW/quadratic extrapolation
    //scaleFactor = 0.5*(dw+de)*(dn+ds)*(2*hr[2]);

    // catch the case where a point lies on the boundary
    //FIXME: do this more elegant!
    //double m1 = dw*de;
    //double m2 = dn*ds;
    //if(de == 0)
    //m1 = dw;
    //if(dw == 0)
    //m1 = de;
    //if(dn == 0)
    //m2 = ds;
    //if(ds == 0)
    //m2 = dn;
    ////XXX: dn+ds || dw+de can be 0
    ////C = 2*(dn+ds)*(dw+de)/hr[2];
    //C = 2/hr[2];
    //if(dw != 0 || de != 0)
    //C *= (dw+de);
    //if(dn != 0 || ds != 0)
    //C *= (dn+ds);
    //if(dw != 0 || de != 0)
    //C += (2*hr[2])*(dn+ds)*(dw+de)/m1;
    //if(dn != 0 || ds != 0)
    //C += (2*hr[2])*(dw+de)*(dn+ds)/m2;

    //handle Neumann case
    //if(z == 0 || z == nr[2]-1) {

    //if(z == 0)
    //F = 0.0;
    //else
    //B = 0.0;

    ////neumann stuff
    //W = W/2.0;
    //E = E/2.0;
    //N = N/2.0;
    //S = S/2.0;
    //C /= 2.0;

    //scaleFactor /= 2.0;
    //}

    // handle boundary condition in z direction
601
    if (z == 0 || z == nr[2] - 1) {
gsell's avatar
gsell committed
602 603 604

        // case where we are on the NEUMAN BC in Z-direction
        // where we distinguish two cases
605
        if (z == 0)
gsell's avatar
gsell committed
606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624
            F = 0.0;
        else
            B = 0.0;

        //C += 2/((hr[2]*(nr[2]-1)/2.0) * hr[2]);
        //hr[2]*(nr2[2]-1)/2 = radius
        double d = hr[2] * (nr[2] - 1) / 2;
        C += 2 / (d * hr[2]);

        W /= 2.0;
        E /= 2.0;
        N /= 2.0;
        S /= 2.0;
        C /= 2.0;
        scaleFactor /= 2.0;

    }
}

625 626 627 628 629 630 631
// vi: set et ts=4 sw=4 sts=4:
// Local Variables:
// mode:c
// c-basic-offset: 4
// indent-tabs-mode: nil
// require-final-newline: nil
// End: