EllipticDomain.cpp 17.5 KB
Newer Older
1
#ifdef HAVE_SAAMG_SOLVER
2

kraus's avatar
kraus committed
3 4
#include "Solvers/EllipticDomain.h"

gsell's avatar
gsell committed
5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23
#include <map>
#include <cmath>
#include <iostream>
#include <assert.h>

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

24
    if(interpl == "CONSTANT")
gsell's avatar
gsell committed
25
        interpolationMethod = CONSTANT;
26
    else if(interpl == "LINEAR")
gsell's avatar
gsell committed
27
        interpolationMethod = LINEAR;
28
    else if(interpl == "QUADRATIC")
gsell's avatar
gsell committed
29 30 31
        interpolationMethod = QUADRATIC;
}

32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64
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());
    Vector_t hr_m; 
    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
65 66 67 68 69 70 71 72
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
73
void EllipticDomain::compute(Vector_t hr){
gsell's avatar
gsell committed
74 75 76 77 78 79
    //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);
80
    hasGeometryChanged_m = true;
gsell's avatar
gsell committed
81 82 83 84 85 86 87 88 89 90 91
    //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
92 93
    int idx = 0;
    int x, y;
94

gsell's avatar
gsell committed
95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122
    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;
123 124 125 126 127 128 129 130 131 132
                if (pos <= -SemiMajor || pos >= SemiMajor)
                {
                	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));
		}

gsell's avatar
gsell committed
133 134 135 136
            }

            for(y = 0; y < nr[1]; y++) {
                pos = y * hr[1] - my;
137 138 139 140 141 142 143 144 145
		if (pos <= -SemiMinor || pos >= SemiMinor)
                {
                	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
146 147 148 149
            }
    }
}

150
void EllipticDomain::compute(Vector_t hr, NDIndex<3> localId){
151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168
    //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
169 170
    int idx = 0;
    int x, y;
171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 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)
                {
                	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));
		}

            }

            for(y = localId[0].first(); y < localId[1].last(); y++) {
                pos = y * hr[1] - my;
		if (pos <= -SemiMinor || pos >= SemiMinor)
                {
                	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));
		}
            }
    }
224 225
}

gsell's avatar
gsell committed
226
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) {
227
    scaleFactor = 1.0;
gsell's avatar
gsell committed
228

229
        // determine which interpolation method we use for points near the boundary
gsell's avatar
gsell committed
230 231
    switch(interpolationMethod) {
        case CONSTANT:
232
            constantInterpolation(x, y, z, W, E, S, N, F, B, C, scaleFactor);
gsell's avatar
gsell committed
233 234
            break;
        case LINEAR:
235
            linearInterpolation(x, y, z, W, E, S, N, F, B, C, scaleFactor);
gsell's avatar
gsell committed
236 237
            break;
        case QUADRATIC:
238
            quadraticInterpolation(x, y, z, W, E, S, N, F, B, C, scaleFactor);
gsell's avatar
gsell committed
239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 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
            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;

}

294
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
295 296 297

    scaleFactor = 1.0;

298 299 300 301 302 303 304
    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
305 306 307

    // we are a right boundary point
    if(!isInside(x + 1, y, z))
308
        EV= 0.0;
gsell's avatar
gsell committed
309 310 311

    // we are a left boundary point
    if(!isInside(x - 1, y, z))
312
        WV = 0.0;
gsell's avatar
gsell committed
313 314 315

    // we are a upper boundary point
    if(!isInside(x, y + 1, z))
316
        NV = 0.0;
gsell's avatar
gsell committed
317 318 319

    // we are a lower boundary point
    if(!isInside(x, y - 1, z))
320
        SV = 0.0;
gsell's avatar
gsell committed
321 322 323 324 325 326 327 328

    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)
329
            FV = 0.0;
gsell's avatar
gsell committed
330
        else
331
            BV = 0.0;
gsell's avatar
gsell committed
332 333 334 335 336 337 338 339

        // 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;
340
        CV += 2 / (d * hr[2]);
gsell's avatar
gsell committed
341 342 343
        //C += 2/((hr[2]*(nr[2]-1)/2.0) * hr[2]);

        // scale all stencil-points in z-plane with 0.5 (Robin discretization)
344 345 346 347 348
        WV /= 2.0;
        EV /= 2.0;
        NV /= 2.0;
        SV /= 2.0;
        CV /= 2.0;
gsell's avatar
gsell committed
349 350 351 352 353
        scaleFactor *= 0.5;
    }

}

354
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
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 414 415 416 417 418 419 420

    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)
        it++;
    dx = it->second;

    double dy = 0.0;
    it = IntersectYDir.find(x);
    if(cy < 0)
        it++;
    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
421
/*
gsell's avatar
gsell committed
422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442
    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;

    }
443
*/
gsell's avatar
gsell committed
444 445
}

446
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
447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 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

    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)
        it++;
    dx = it->second;

    double dy = 0.0;
    it = IntersectYDir.find(x);
    if(cy < 0)
        it++;
    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;

    }
}


664
#endif //#ifdef HAVE_SAAMG_SOLVER