EllipticDomain.cpp 18.8 KB
Newer Older
1
//
frey_m's avatar
frey_m committed
2 3
// Class EllipticDomain
//   :FIXME: add brief description
4
//
5
// Copyright (c) 2008,        Yves Ineichen, ETH Zürich,
frey_m's avatar
frey_m committed
6
//               2013 - 2015, Tülin Kaman, Paul Scherrer Institut, Villigen PSI, Switzerland
7
//               2017 - 2020, Paul Scherrer Institut, Villigen PSI, Switzerland
frey_m's avatar
frey_m committed
8
// All rights reserved
9
//
10 11 12 13 14
// Implemented as part of the master thesis
// "A Parallel Multigrid Solver for Beam Dynamics"
// and the paper
// "A fast parallel Poisson solver on irregular domains applied to beam dynamics simulations"
// (https://doi.org/10.1016/j.jcp.2010.02.022)
frey_m's avatar
frey_m committed
15 16 17 18 19 20 21 22 23 24
//
// 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/>.
25
//
kraus's avatar
kraus committed
26 27
#include "Solvers/EllipticDomain.h"

gsell's avatar
gsell committed
28 29 30
#include <map>
#include <cmath>
#include <iostream>
31
#include <cassert>
gsell's avatar
gsell committed
32 33 34 35 36 37 38 39 40 41 42 43 44 45 46

//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);
}

EllipticDomain::EllipticDomain(double semimajor, double semiminor, Vector_t nr, Vector_t hr, std::string interpl) {
    SemiMajor = semimajor;
    SemiMinor = semiminor;
    setNr(nr);
    setHr(hr);

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

55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73
EllipticDomain::EllipticDomain(BoundaryGeometry *bgeom, Vector_t nr, Vector_t hr, std::string interpl) {
    SemiMajor = bgeom->getA();
    SemiMinor = bgeom->getB();
    setMinMaxZ(bgeom->getS(), bgeom->getLength());
    setNr(nr);
    setHr(hr);

    if(interpl == "CONSTANT")
        interpolationMethod = CONSTANT;
    else if(interpl == "LINEAR")
        interpolationMethod = LINEAR;
    else if(interpl == "QUADRATIC")
        interpolationMethod = QUADRATIC;
}

EllipticDomain::EllipticDomain(BoundaryGeometry *bgeom, Vector_t nr, std::string interpl) {
    SemiMajor = bgeom->getA();
    SemiMinor = bgeom->getB();
    setMinMaxZ(bgeom->getS(), bgeom->getLength());
kraus's avatar
kraus committed
74
    Vector_t hr_m;
75 76 77 78 79 80 81 82 83 84 85 86 87
    hr_m[0] = (getXRangeMax()-getXRangeMin())/nr[0];
    hr_m[1] = (getYRangeMax()-getYRangeMin())/nr[1];
    hr_m[2] = (getZRangeMax()-getZRangeMin())/nr[2];
    setHr(hr_m);

    if(interpl == "CONSTANT")
        interpolationMethod = CONSTANT;
    else if(interpl == "LINEAR")
        interpolationMethod = LINEAR;
    else if(interpl == "QUADRATIC")
        interpolationMethod = QUADRATIC;
}

gsell's avatar
gsell committed
88 89 90 91 92 93 94 95
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
96
void EllipticDomain::compute(Vector_t hr){
gsell's avatar
gsell committed
97 98 99 100 101 102
    //there is nothing to be done if the mesh spacings have not changed
    if(hr[0] == getHr()[0] && hr[1] == getHr()[1] && hr[2] == getHr()[2]) {
        hasGeometryChanged_m = false;
        return;
    }
    setHr(hr);
103
    hasGeometryChanged_m = true;
gsell's avatar
gsell committed
104 105 106 107 108 109 110 111 112 113 114
    //reset number of points inside domain
    nxy_m = 0;

    // clear previous coordinate maps
    IdxMap.clear();
    CoordMap.clear();
    //clear previous intersection points
    IntersectYDir.clear();
    IntersectXDir.clear();

    // build a index and coordinate map
115 116
    int idx = 0;
    int x, y;
117

gsell's avatar
gsell committed
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
    for(x = 0; x < nr[0]; x++) {
        for(y = 0; y < nr[1]; y++) {
            if(isInside(x, y, 1)) {
                IdxMap[toCoordIdx(x, y)] = idx;
                CoordMap[idx++] = toCoordIdx(x, y);
                nxy_m++;
            }

        }
    }

    switch(interpolationMethod) {
        case CONSTANT:
            break;
        case LINEAR:
        case QUADRATIC:

            double smajsq = SemiMajor * SemiMajor;
            double sminsq = SemiMinor * SemiMinor;
            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 = 0; x < nr[0]; x++) {
                pos = x * hr[0] - mx;
146 147
                if (pos <= -SemiMajor || pos >= SemiMajor)
                {
kraus's avatar
kraus committed
148 149 150 151 152 153 154
                    IntersectYDir.insert(std::pair<int, double>(x, 0));
                    IntersectYDir.insert(std::pair<int, double>(x, 0));
                }else{
                    yd = std::abs(sqrt(sminsq - sminsq * pos * pos / smajsq)); // + 0.5*nr[1]*hr[1]);
                    IntersectYDir.insert(std::pair<int, double>(x, yd));
                    IntersectYDir.insert(std::pair<int, double>(x, -yd));
                }
155

gsell's avatar
gsell committed
156 157 158 159
            }

            for(y = 0; y < nr[1]; y++) {
                pos = y * hr[1] - my;
kraus's avatar
kraus committed
160
                if (pos <= -SemiMinor || pos >= SemiMinor)
161
                {
kraus's avatar
kraus committed
162 163 164 165 166 167 168
                    IntersectXDir.insert(std::pair<int, double>(y, 0));
                    IntersectXDir.insert(std::pair<int, double>(y, 0));
                }else{
                    xd = std::abs(sqrt(smajsq - smajsq * pos * pos / sminsq)); // + 0.5*nr[0]*hr[0]);
                    IntersectXDir.insert(std::pair<int, double>(y, xd));
                    IntersectXDir.insert(std::pair<int, double>(y, -xd));
                }
gsell's avatar
gsell committed
169 170 171 172
            }
    }
}

173
void EllipticDomain::compute(Vector_t hr, NDIndex<3> localId){
174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191
    //there is nothing to be done if the mesh spacings have not changed
    if(hr[0] == getHr()[0] && hr[1] == getHr()[1] && hr[2] == getHr()[2]) {
        hasGeometryChanged_m = false;
        return;
    }
    setHr(hr);
    hasGeometryChanged_m = true;
    //reset number of points inside domain
    nxy_m = 0;

    // clear previous coordinate maps
    IdxMap.clear();
    CoordMap.clear();
    //clear previous intersection points
    IntersectYDir.clear();
    IntersectXDir.clear();

    // build a index and coordinate map
192 193
    int idx = 0;
    int x, y;
194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223

    for(x = localId[0].first();  x<= localId[0].last(); x++) {
        for(y = localId[1].first(); y <= localId[1].last(); y++) {
            if(isInside(x, y, 1)) {
                IdxMap[toCoordIdx(x, y)] = idx;
                CoordMap[idx++] = toCoordIdx(x, y);
                nxy_m++;
            }
        }
    }

    switch(interpolationMethod) {
        case CONSTANT:
            break;
        case LINEAR:
        case QUADRATIC:

            double smajsq = SemiMajor * SemiMajor;
            double sminsq = SemiMinor * SemiMinor;
            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;
                if (pos <= -SemiMajor || pos >= SemiMajor)
                {
kraus's avatar
kraus committed
224 225 226 227 228 229 230
                    IntersectYDir.insert(std::pair<int, double>(x, 0));
                    IntersectYDir.insert(std::pair<int, double>(x, 0));
                }else{
                    yd = std::abs(sqrt(sminsq - sminsq * pos * pos / smajsq)); // + 0.5*nr[1]*hr[1]);
                    IntersectYDir.insert(std::pair<int, double>(x, yd));
                    IntersectYDir.insert(std::pair<int, double>(x, -yd));
                }
231 232 233 234 235

            }

            for(y = localId[0].first(); y < localId[1].last(); y++) {
                pos = y * hr[1] - my;
kraus's avatar
kraus committed
236
                if (pos <= -SemiMinor || pos >= SemiMinor)
237
                {
kraus's avatar
kraus committed
238 239 240 241 242 243 244
                    IntersectXDir.insert(std::pair<int, double>(y, 0));
                    IntersectXDir.insert(std::pair<int, double>(y, 0));
                }else{
                    xd = std::abs(sqrt(smajsq - smajsq * pos * pos / sminsq)); // + 0.5*nr[0]*hr[0]);
                    IntersectXDir.insert(std::pair<int, double>(y, xd));
                    IntersectXDir.insert(std::pair<int, double>(y, -xd));
                }
245 246
            }
    }
247 248
}

gsell's avatar
gsell committed
249
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) {
250
    scaleFactor = 1.0;
gsell's avatar
gsell committed
251

252
        // determine which interpolation method we use for points near the boundary
gsell's avatar
gsell committed
253 254
    switch(interpolationMethod) {
        case CONSTANT:
255
            constantInterpolation(x, y, z, W, E, S, N, F, B, C, scaleFactor);
gsell's avatar
gsell committed
256 257
            break;
        case LINEAR:
258
            linearInterpolation(x, y, z, W, E, S, N, F, B, C, scaleFactor);
gsell's avatar
gsell committed
259 260
            break;
        case QUADRATIC:
261
            quadraticInterpolation(x, y, z, W, E, S, N, F, B, C, scaleFactor);
gsell's avatar
gsell committed
262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316
            break;
    }

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

void EllipticDomain::getBoundaryStencil(int idx, double &W, double &E, double &S, double &N, double &F, double &B, double &C, double &scaleFactor) {


    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);

}

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

    if(x > 0)
        W = getIdx(x - 1, y, z);
    else
        W = -1;
    if(x < nr[0] - 1)
        E = getIdx(x + 1, y, z);
    else
        E = -1;

    if(y < nr[1] - 1)
        N = getIdx(x, y + 1, z);
    else
        N = -1;
    if(y > 0)
        S = getIdx(x, y - 1, z);
    else
        S = -1;

    if(z > 0)
        F = getIdx(x, y, z - 1);
    else
        F = -1;
    if(z < nr[2] - 1)
        B = getIdx(x, y, z + 1);
    else
        B = -1;

}

317
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
318 319 320

    scaleFactor = 1.0;

321 322 323 324 325 326 327
    WV = -1/(hr[0]*hr[0]);
    EV = -1/(hr[0]*hr[0]);
    NV = -1/(hr[1]*hr[1]);
    SV = -1/(hr[1]*hr[1]);
    FV = -1/(hr[2]*hr[2]);
    BV = -1/(hr[2]*hr[2]);
    CV = 2/(hr[0]*hr[0]) + 2/(hr[1]*hr[1]) + 2/(hr[2]*hr[2]);
gsell's avatar
gsell committed
328 329 330

    // we are a right boundary point
    if(!isInside(x + 1, y, z))
331
        EV= 0.0;
gsell's avatar
gsell committed
332 333 334

    // we are a left boundary point
    if(!isInside(x - 1, y, z))
335
        WV = 0.0;
gsell's avatar
gsell committed
336 337 338

    // we are a upper boundary point
    if(!isInside(x, y + 1, z))
339
        NV = 0.0;
gsell's avatar
gsell committed
340 341 342

    // we are a lower boundary point
    if(!isInside(x, y - 1, z))
343
        SV = 0.0;
gsell's avatar
gsell committed
344 345 346 347 348 349 350 351

    if(z == 1 || z == nr[2] - 2) {

        // 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
        if(z == 1)
352
            FV = 0.0;
gsell's avatar
gsell committed
353
        else
354
            BV = 0.0;
gsell's avatar
gsell committed
355 356 357 358 359 360 361 362

        // add contribution of Robin discretization to center point
        // d the distance between the center of the bunch and the boundary
        //double cx = (x-(nr[0]-1)/2)*hr[0];
        //double cy = (y-(nr[1]-1)/2)*hr[1];
        //double cz = hr[2]*(nr[2]-1);
        //double d = sqrt(cx*cx+cy*cy+cz*cz);
        double d = hr[2] * (nr[2] - 1) / 2;
363
        CV += 2 / (d * hr[2]);
gsell's avatar
gsell committed
364 365 366
        //C += 2/((hr[2]*(nr[2]-1)/2.0) * hr[2]);

        // scale all stencil-points in z-plane with 0.5 (Robin discretization)
367 368 369 370 371
        WV /= 2.0;
        EV /= 2.0;
        NV /= 2.0;
        SV /= 2.0;
        CV /= 2.0;
gsell's avatar
gsell committed
372 373 374 375 376
        scaleFactor *= 0.5;
    }

}

377
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
378 379 380 381 382 383 384 385 386

    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;
    std::multimap<int, double>::iterator it = IntersectXDir.find(y);
    if(cx < 0)
387
        ++it;
gsell's avatar
gsell committed
388 389 390 391 392
    dx = it->second;

    double dy = 0.0;
    it = IntersectYDir.find(x);
    if(cy < 0)
393
        ++it;
gsell's avatar
gsell committed
394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443
    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
    if(!isInside(x + 1, y, z)) {
        C += 1 / ((dx - cx) * de);
        E = 0.0;
    } else {
        C += 1 / (de * de);
        E = -1 / (de * de);
    }

    //we are a left boundary point
    if(!isInside(x - 1, y, z)) {
        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
    if(!isInside(x, y + 1, z)) {
        C += 1 / ((dy - cy) * dn);
        N = 0.0;
    } else {
        C += 1 / (dn * dn);
        N = -1 / (dn * dn);
    }

    //we are a lower boundary point
    if(!isInside(x, y - 1, z)) {
        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
444
/*
gsell's avatar
gsell committed
445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465
    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;

    }
466
*/
gsell's avatar
gsell committed
467 468
}

469
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
470 471 472 473 474 475 476 477 478

    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;
    std::multimap<int, double>::iterator it = IntersectXDir.find(y);
    if(cx < 0)
479
        ++it;
gsell's avatar
gsell committed
480 481 482 483 484
    dx = it->second;

    double dy = 0.0;
    it = IntersectYDir.find(x);
    if(cy < 0)
485
        ++it;
gsell's avatar
gsell committed
486 487 488 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 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685
    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
    if(!isInside(x + 1, y, z)) {
        de = dx - cx;
        E = 0.0;
    }

    //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);
        W = 0.0;
    }

    //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;
        N = 0.0;
    }

    //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);
        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
    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;

        //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;

    }
}

686 687 688 689 690 691 692
// 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: