Liebe Gitlab-Nutzer, lieber Gitlab-Nutzer, es ist nun möglich sich mittels des ZIH-Logins/LDAP an unserem Dienst anzumelden. Ein Anmelden über dieses erzeugt ein neues Konto. Das alte Konto ist über den Reiter "Standard" erreichbar. Die Administratoren

Dear Gitlab user, it is now possible to log in to our service using the ZIH login/LDAP. Logging in via this will create a new account. The old account can be accessed via the "Standard" tab. The administrators

PetscSolverFeti.cc 29.6 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22
//
// Software License for AMDiS
//
// Copyright (c) 2010 Dresden University of Technology 
// All rights reserved.
// Authors: Simon Vey, Thomas Witkowski et al.
//
// This file is part of AMDiS
//
// See also license.opensource.txt in the distribution.


#include "parallel/PetscSolverFeti.h"
#include "parallel/StdMpi.h"
#include "parallel/MpiHelper.h"

namespace AMDiS {

  using namespace std;


#ifdef HAVE_PETSC_DEV 
23 24 25 26 27 28 29 30 31 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 65 66 67 68 69 70 71 72
  // y = mat * x
  int petscMultMatSchurPrimal(Mat mat, Vec x, Vec y)
  {
    // S_PiPi = K_PiPi - K_PiB inv(K_BB) K_BPi

    void *ctx;
    MatShellGetContext(mat, &ctx);
    PetscSchurPrimalData* data = static_cast<PetscSchurPrimalData*>(ctx);

    MatMult(*(data->mat_b_primal), x, data->tmp_vec_b);
    KSPSolve(*(data->ksp_b), data->tmp_vec_b, data->tmp_vec_b);

    MatMult(*(data->mat_primal_b), data->tmp_vec_b, data->tmp_vec_primal);
    MatMult(*(data->mat_primal_primal), x, y);
    VecAXPBY(y, -1.0, 1.0, data->tmp_vec_primal);

    return 0;
  }


  // y = mat * x
  int petscMultMatFeti(Mat mat, Vec x, Vec y)
  {
    // F = L inv(K_BB) trans(L) + L inv(K_BB) K_BPi inv(S_PiPi) K_PiB inv(K_BB) trans(L)

    void *ctx;
    MatShellGetContext(mat, &ctx);
    PetscFetiData* data = static_cast<PetscFetiData*>(ctx);

    // y = L inv(K_BB) trans(L) x
    MatMultTranspose(*(data->mat_lagrange), x, data->tmp_vec_b);
    KSPSolve(*(data->ksp_b), data->tmp_vec_b, data->tmp_vec_b);
    MatMult(*(data->mat_lagrange), data->tmp_vec_b, y);

    // tmp_vec_primal = inv(S_PiPi) K_PiB inv(K_BB) trans(L)
    MatMult(*(data->mat_primal_b), data->tmp_vec_b, data->tmp_vec_primal);
    KSPSolve(*(data->ksp_schur_primal), data->tmp_vec_primal, data->tmp_vec_primal);

    // tmp_vec_lagrange = L inv(K_BB) K_BPi tmp_vec_primal
    //                  = L inv(K_BB) K_BPi inv(S_PiPi) K_PiB inv(K_BB) trans(L)
    MatMult(*(data->mat_b_primal), data->tmp_vec_primal, data->tmp_vec_b);
    KSPSolve(*(data->ksp_b), data->tmp_vec_b, data->tmp_vec_b);
    MatMult(*(data->mat_lagrange), data->tmp_vec_b, data->tmp_vec_lagrange);

    VecAXPBY(y, 1.0, 1.0, data->tmp_vec_lagrange);

    return 0;
  }


73
  void PetscSolverFeti::updateDofData()
74 75
  {
    FUNCNAME("PetscSolverFeti::updateDofData()");
76 77 78 79 80 81

    TEST_EXIT(meshDistributor->getMesh()->getDim() == 2)
      ("Works for 2D problems only!");

    TEST_EXIT(meshDistributor->getFeSpace()->getBasisFcts()->getDegree() == 1)
      ("Works for linear basis functions only!\n");
82
   
83 84 85 86 87 88 89
    createPrimals();

    createDuals();

    createLagrange();

    createIndexB();
90 91 92
  }


93
  void PetscSolverFeti::createPrimals()
94
  {
95
    FUNCNAME("PetscSolverFeti::createPrimals()");  
96

97 98 99 100 101 102 103
    primals.clear();
    DofContainerSet& vertices = 
      meshDistributor->getBoundaryDofInfo().geoDofs[VERTEX];
    TEST_EXIT_DBG(vertices.size())("No primal vertices on this rank!\n");
    for (DofContainerSet::iterator it = vertices.begin(); 
	 it != vertices.end(); ++it)
      primals.insert(**it);
104 105

    globalPrimalIndex.clear();
106 107 108 109
    nRankPrimals = 0;
    for (DofIndexSet::iterator it = primals.begin(); it != primals.end(); ++it)
      if (meshDistributor->getIsRankDof(*it)) {
	globalPrimalIndex[*it] = nRankPrimals;
110 111 112
	nRankPrimals++;
      }

113 114

    nOverallPrimals = 0;
115
    rStartPrimals = 0;
116 117 118 119 120 121 122 123 124
    mpi::getDofNumbering(meshDistributor->getMpiComm(),
			 nRankPrimals, rStartPrimals, nOverallPrimals);

    for (DofMapping::iterator it = globalPrimalIndex.begin();
	 it != globalPrimalIndex.end(); ++it)
      it->second += rStartPrimals;

    MSG_DBG("nRankPrimals = %d   nOverallPrimals = %d\n",
	    nRankPrimals, nOverallPrimals);
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 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 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 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

    StdMpi<vector<int> > stdMpi(meshDistributor->getMpiComm());
    RankToDofContainer& sendDofs = meshDistributor->getSendDofs();
    for (RankToDofContainer::iterator it = sendDofs.begin();
	 it != sendDofs.end(); ++it)
      for (DofContainer::iterator dofIt = it->second.begin();
	   dofIt != it->second.end(); ++dofIt)
	if (globalPrimalIndex.count(**dofIt))
	  stdMpi.getSendData(it->first).push_back(globalPrimalIndex[**dofIt]);
    stdMpi.updateSendDataSize();

    RankToDofContainer& recvDofs = meshDistributor->getRecvDofs();
    for (RankToDofContainer::iterator it = recvDofs.begin();
	 it != recvDofs.end(); ++it) {
      bool recvFromRank = false;
      for (DofContainer::iterator dofIt = it->second.begin();
	   dofIt != it->second.end(); ++dofIt)
	if (primals.count(**dofIt) && 
	    meshDistributor->getIsRankDof(**dofIt) == false) {
	  recvFromRank = true;
	  break;
	}

      if (recvFromRank) 
	stdMpi.recv(it->first);
    }
    stdMpi.startCommunication();

    for (RankToDofContainer::iterator it = recvDofs.begin();
	 it != recvDofs.end(); ++it) {
      int i = 0;
      for (DofContainer::iterator dofIt = it->second.begin();
	   dofIt != it->second.end(); ++dofIt) {
	if (primals.count(**dofIt) && 
	    meshDistributor->getIsRankDof(**dofIt) == false)
	  globalPrimalIndex[**dofIt] = stdMpi.getRecvData(it->first)[i++];
      }
    }

    TEST_EXIT_DBG(primals.size() == globalPrimalIndex.size())
      ("Number of primals %d, but number of global primals on this rank is %d!\n",
       primals.size(), globalPrimalIndex.size());


    TEST_EXIT_DBG(nOverallPrimals > 0)
      ("There are zero primal nodes in domain!\n");
  }


  void PetscSolverFeti::createDuals()
  {
    FUNCNAME("PetscSolverFeti::createDuals()");
    
    // === Create for each dual node the set of ranks that contain this ===
    // === node (denoted by W(x_j)).                                    ===

    boundaryDofRanks.clear();

    RankToDofContainer& sendDofs = meshDistributor->getSendDofs();
    for (RankToDofContainer::iterator it = sendDofs.begin();
	 it != sendDofs.end(); ++it) {
      for (DofContainer::iterator dofIt = it->second.begin();
	   dofIt != it->second.end(); ++dofIt) {
	// If DOF is not primal, i.e., its a dual node
	if (primals.count(**dofIt) == 0) {
	  boundaryDofRanks[**dofIt].insert(mpiRank);
	  boundaryDofRanks[**dofIt].insert(it->first);
	}
      }
    }

    StdMpi<vector<std::set<int> > > stdMpi(meshDistributor->getMpiComm());
    for (RankToDofContainer::iterator it = sendDofs.begin();
	 it != sendDofs.end(); ++it)
      for (DofContainer::iterator dofIt = it->second.begin();
	   dofIt != it->second.end(); ++dofIt)
	if (primals.count(**dofIt) == 0)
	  stdMpi.getSendData(it->first).push_back(boundaryDofRanks[**dofIt]);

    stdMpi.updateSendDataSize();

    RankToDofContainer& recvDofs = meshDistributor->getRecvDofs();
    for (RankToDofContainer::iterator it = recvDofs.begin();
	 it != recvDofs.end(); ++it) {
      bool recvFromRank = false;
      for (DofContainer::iterator dofIt = it->second.begin();
	   dofIt != it->second.end(); ++dofIt)
	if (primals.count(**dofIt) == 0) {
	  recvFromRank = true;
	  break;
	}

      if (recvFromRank)
	stdMpi.recv(it->first);
    }
    stdMpi.startCommunication();

    for (RankToDofContainer::iterator it = recvDofs.begin();
	 it != recvDofs.end(); ++it) {
      int i = 0;
      for (DofContainer::iterator dofIt = it->second.begin();
	   dofIt != it->second.end(); ++dofIt)	
	if (primals.count(**dofIt) == 0)
	  boundaryDofRanks[**dofIt] = stdMpi.getRecvData(it->first)[i++];	      
    }


    // === Create global index of the dual nodes on each rank. ===

    duals.clear();
    globalDualIndex.clear();

    int nRankAllDofs = meshDistributor->getFeSpace()->getAdmin()->getUsedDofs();
    nRankB = nRankAllDofs - primals.size();
    nOverallB = 0;
    rStartB = 0;
    mpi::getDofNumbering(meshDistributor->getMpiComm(),
			 nRankB, rStartB, nOverallB);
    DofContainer allBoundaryDofs;
    meshDistributor->getAllBoundaryDofs(allBoundaryDofs);
    int nRankInteriorDofs = nRankAllDofs - allBoundaryDofs.size();

    int nRankDuals = 0;
    for (DofContainer::iterator it = allBoundaryDofs.begin();
	 it != allBoundaryDofs.end(); ++it) {
      if (primals.count(**it) == 0) {
	duals.insert(**it);
	globalDualIndex[**it] = rStartB + nRankInteriorDofs + nRankDuals;
	nRankDuals++;
      }
    }

    int nOverallDuals = nRankDuals;
    mpi::globalAdd(nOverallDuals);

    MSG_DBG("nRankDuals = %d   nOverallDuals = %d\n",
	    nRankDuals, nOverallDuals);
  }

  
  void PetscSolverFeti::createLagrange()
  {
    FUNCNAME("PetscSolverFeti::createLagrange()");

    nRankLagrange = 0;
    for (DofIndexSet::iterator it = duals.begin(); it != duals.end(); ++it) {
      if (meshDistributor->getIsRankDof(*it)) {
	dofFirstLagrange[*it] = nRankLagrange;
	int degree = boundaryDofRanks[*it].size();
	nRankLagrange += (degree * (degree - 1)) / 2;
      }
    }

    nOverallLagrange = 0;
279
    rStartLagrange = 0;
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 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359
    mpi::getDofNumbering(meshDistributor->getMpiComm(),
			 nRankLagrange, rStartLagrange, nOverallLagrange);

    for (DofIndexSet::iterator it = duals.begin(); it != duals.end(); ++it)
      if (meshDistributor->getIsRankDof(*it))
	dofFirstLagrange[*it] += rStartLagrange;

    MSG_DBG("nRankLagrange = %d  nOverallLagrange = %d\n",
	    nRankLagrange, nOverallLagrange);


    // === ===

    StdMpi<vector<int> > stdMpi(meshDistributor->getMpiComm());
    RankToDofContainer& sendDofs = meshDistributor->getSendDofs();
    for (RankToDofContainer::iterator it = sendDofs.begin();
	 it != sendDofs.end(); ++it)
      for (DofContainer::iterator dofIt = it->second.begin();
	   dofIt != it->second.end(); ++dofIt) {
	if (primals.count(**dofIt) == 0) {
	  TEST_EXIT_DBG(dofFirstLagrange.count(**dofIt))("Should not happen!\n");
	  stdMpi.getSendData(it->first).push_back(dofFirstLagrange[**dofIt]);
	}
      }
    stdMpi.updateSendDataSize();

    RankToDofContainer& recvDofs = meshDistributor->getRecvDofs();
    for (RankToDofContainer::iterator it = recvDofs.begin();
	 it != recvDofs.end(); ++it) {
      bool recvData = false;
      for (DofContainer::iterator dofIt = it->second.begin();
	   dofIt != it->second.end(); ++dofIt)
	if (primals.count(**dofIt) == 0) {
	  recvData = true;
	  break;
	}
	  
      if (recvData)
	stdMpi.recv(it->first);
    }

    stdMpi.startCommunication();

    for (RankToDofContainer::iterator it = recvDofs.begin();
	 it != recvDofs.end(); ++it) {
      int counter = 0;
      for (unsigned int i = 0; i < it->second.size(); i++)
	if (primals.count(*(it->second[i])) == 0)
	  dofFirstLagrange[*(it->second[i])] = stdMpi.getRecvData(it->first)[counter++];
    }
     
  }


  void PetscSolverFeti::createIndexB()
  {
    FUNCNAME("PetscSolverFeti::createIndeB()");

    globalIndexB.clear();
    DOFAdmin* admin = meshDistributor->getFeSpace()->getAdmin();
    
    for (int i = 0; i < admin->getUsedSize(); i++)
      if (admin->isDofFree(i) == false && primals.count(i) == 0)
	if (duals.count(i) == 0 && primals.count(i) == 0)
	  globalIndexB[i] = -1;

    int nInterior = 0;
    for (DofMapping::iterator it = globalIndexB.begin(); 
	 it != globalIndexB.end(); ++it) {
      it->second = nInterior + rStartB;
      nInterior++;
    }

    TEST_EXIT_DBG(nInterior + primals.size() + duals.size() == 
		  static_cast<unsigned int>(admin->getUsedDofs()))
      ("Should not happen!\n");

    for (DofIndexSet::iterator it = duals.begin();
	 it != duals.end(); ++it)
      globalIndexB[*it] = globalDualIndex[*it];
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
  void PetscSolverFeti::createMatLagrange(int nComponents)
  {
    FUNCNAME("PetscSolverFeti::createMatLagrange()");

    MatCreateMPIAIJ(PETSC_COMM_WORLD,
		    nRankLagrange * nComponents, nRankB * nComponents,
		    nOverallLagrange * nComponents, nOverallB * nComponents,
		    2, PETSC_NULL, 2, PETSC_NULL,
		    &mat_lagrange);

    for (DofIndexSet::iterator it = duals.begin(); it != duals.end(); ++it) {
      TEST_EXIT_DBG(dofFirstLagrange.count(*it))("Should not happen!\n");
      TEST_EXIT_DBG(boundaryDofRanks.count(*it))("Should not happen!\n");

      int index = dofFirstLagrange[*it];
      vector<int> W(boundaryDofRanks[*it].begin(), boundaryDofRanks[*it].end());
      int degree = W.size();

      TEST_EXIT_DBG(globalDualIndex.count(*it))("Should not happen!\n");
      int dualCol = globalDualIndex[*it];

      for (int i = 0; i < degree; i++) {
	for (int j = i + 1; j < degree; j++) {
	  if (W[i] == mpiRank || W[j] == mpiRank) {	      
	    for (int k = 0; k < nComponents; k++) {
	      int rowIndex = index * nComponents + k;
	      int colIndex = dualCol * nComponents + k;
	      double value = (W[i] == mpiRank ? 1.0 : -1.0);
	      MatSetValue(mat_lagrange, rowIndex, colIndex, value, 
			  INSERT_VALUES);
	    }
	  }

	  index++;
	}
      }
    }

    MatAssemblyBegin(mat_lagrange, MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(mat_lagrange, MAT_FINAL_ASSEMBLY);
  }


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 444 445 446 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
  void PetscSolverFeti::createSchurPrimalKsp(int nComponents)
  {
    FUNCNAME("PetscSolverFeti::createSchurPrimal()");

    petscSchurPrimalData.mat_primal_primal = &mat_primal_primal;
    petscSchurPrimalData.mat_primal_b = &mat_primal_b;
    petscSchurPrimalData.mat_b_primal = &mat_b_primal;
    petscSchurPrimalData.ksp_b = &ksp_b;

    MatGetVecs(mat_b_b, 
	       PETSC_NULL, &(petscSchurPrimalData.tmp_vec_b));
    MatGetVecs(mat_primal_primal, 
	       PETSC_NULL, &(petscSchurPrimalData.tmp_vec_primal));

    MatCreateShell(PETSC_COMM_WORLD,
		   nRankPrimals * nComponents, nRankPrimals * nComponents,
		   nOverallPrimals * nComponents, nOverallPrimals * nComponents,
		   &petscSchurPrimalData, 
		   &mat_schur_primal);
    MatShellSetOperation(mat_schur_primal, MATOP_MULT, 
			 (void(*)(void))petscMultMatSchurPrimal);

    KSPCreate(PETSC_COMM_WORLD, &ksp_schur_primal);
    KSPSetOperators(ksp_schur_primal, mat_schur_primal, mat_schur_primal, SAME_NONZERO_PATTERN);
    KSPSetOptionsPrefix(ksp_schur_primal, "solver_sp_");
    KSPSetFromOptions(ksp_schur_primal);
  }


  void PetscSolverFeti::destroySchurPrimalKsp()
  {
    FUNCNAME("PetscSolverFeti::destroySchurPrimal()");

    petscSchurPrimalData.mat_primal_primal = PETSC_NULL;
    petscSchurPrimalData.mat_primal_b = PETSC_NULL;
    petscSchurPrimalData.mat_b_primal = PETSC_NULL;
    petscSchurPrimalData.ksp_b = PETSC_NULL;

    VecDestroy(petscSchurPrimalData.tmp_vec_b);
    VecDestroy(petscSchurPrimalData.tmp_vec_primal);

    MatDestroy(mat_schur_primal);
    KSPDestroy(ksp_schur_primal);
  }


  void PetscSolverFeti::createFetiKsp()
  {
    FUNCNAME("PetscSolverFeti::createFetiKsp()");

    petscFetiData.mat_primal_primal = &mat_primal_primal;
    petscFetiData.mat_primal_b = &mat_primal_b;
    petscFetiData.mat_b_primal = &mat_b_primal;
    petscFetiData.mat_lagrange = &mat_lagrange;
    petscFetiData.ksp_b = &ksp_b;
    petscFetiData.ksp_schur_primal = &ksp_schur_primal;


    MatGetVecs(mat_b_b, PETSC_NULL, &(petscFetiData.tmp_vec_b));
    MatGetVecs(mat_primal_primal, PETSC_NULL, &(petscFetiData.tmp_vec_primal));
    MatGetVecs(mat_lagrange, PETSC_NULL, &(petscFetiData.tmp_vec_lagrange));


    MatCreateShell(PETSC_COMM_WORLD,
		   nRankLagrange, nRankLagrange,
		   nOverallLagrange, nOverallLagrange,
		   &petscFetiData, &mat_feti);
    MatShellSetOperation(mat_feti, MATOP_MULT, (void(*)(void))petscMultMatFeti);


    KSPCreate(PETSC_COMM_WORLD, &ksp_feti);
    KSPSetOperators(ksp_feti, mat_feti, mat_feti, SAME_NONZERO_PATTERN);
    KSPSetOptionsPrefix(ksp_feti, "solver_feti_");
    KSPSetFromOptions(ksp_feti);
  }
  

  void PetscSolverFeti::destroyFetiKsp()
  {
    FUNCNAME("PetscSolverFeti::destroyFetiKsp()");

    petscFetiData.mat_primal_primal = PETSC_NULL;
    petscFetiData.mat_primal_b = PETSC_NULL;
    petscFetiData.mat_b_primal = PETSC_NULL;
    petscFetiData.mat_lagrange = PETSC_NULL;
    petscFetiData.ksp_b = PETSC_NULL;
    petscFetiData.ksp_schur_primal = PETSC_NULL;


    VecDestroy(petscFetiData.tmp_vec_b);
    VecDestroy(petscFetiData.tmp_vec_primal);
    VecDestroy(petscFetiData.tmp_vec_lagrange);


    MatDestroy(mat_feti);
    KSPDestroy(ksp_feti);
  }


  void PetscSolverFeti::recoverSolution(Vec &vec_sol_b,
					Vec &vec_sol_primal,
					SystemVector &vec)
  {
    FUNCNAME("PetscSolverFeti::recoverSolution()");

    int nComponents = vec.getSize();

    // === ===

    PetscScalar *localSolB;
    VecGetArray(vec_sol_b, &localSolB);


    // === ===
    
    vector<PetscInt> globalIsIndex, localIsIndex;
    globalIsIndex.reserve(globalPrimalIndex.size() * nComponents);
    localIsIndex.reserve(globalPrimalIndex.size() * nComponents);

    {
      int counter = 0;
      for (DofMapping::iterator it = globalPrimalIndex.begin();
	   it != globalPrimalIndex.end(); ++it) {
	for (int i = 0; i < nComponents; i++) {
	  globalIsIndex.push_back(it->second * nComponents + i);
	  localIsIndex.push_back(counter++);
	}
      }
    }
    
    IS globalIs, localIs;
    ISCreateGeneral(PETSC_COMM_SELF, 
		    globalIsIndex.size(), 
		    &(globalIsIndex[0]),
		    PETSC_USE_POINTER,
		    &globalIs);

    ISCreateGeneral(PETSC_COMM_SELF, 
		    localIsIndex.size(), 
		    &(localIsIndex[0]),
		    PETSC_USE_POINTER,
		    &localIs);

    Vec local_sol_primal;
    VecCreateSeq(PETSC_COMM_SELF, localIsIndex.size(), &local_sol_primal);


    VecScatter primalScatter;
    VecScatterCreate(vec_sol_primal, globalIs, local_sol_primal, localIs, &primalScatter);
    VecScatterBegin(primalScatter, vec_sol_primal, local_sol_primal, 
		    INSERT_VALUES, SCATTER_FORWARD);
    VecScatterEnd(primalScatter, vec_sol_primal, local_sol_primal, 
		  INSERT_VALUES, SCATTER_FORWARD);

    ISDestroy(globalIs);
    ISDestroy(localIs);    
    VecScatterDestroy(primalScatter);    


    PetscScalar *localSolPrimal;
    VecGetArray(local_sol_primal, &localSolPrimal);




    for (int i = 0; i < nComponents; i++) {
      DOFVector<double>& dofVec = *(vec.getDOFVector(i));

      for (DofMapping::iterator it = globalIndexB.begin();
	   it != globalIndexB.end(); ++it) {
	int petscIndex = (it->second - rStartB) * nComponents + i;
	dofVec[it->first] = localSolB[petscIndex];
      }

      int counter = 0;
      for (DofMapping::iterator it = globalPrimalIndex.begin();
	   it != globalPrimalIndex.end(); ++it) {
	dofVec[it->first] = localSolPrimal[counter * nComponents + i];
	counter++;
      }
    }



    VecRestoreArray(vec_sol_b, &localSolB);

    VecRestoreArray(local_sol_primal, &localSolPrimal);
    VecDestroy(local_sol_primal);
  }


597 598
  void PetscSolverFeti::fillPetscMatrix(Matrix<DOFMatrix*> *mat, 
					SystemVector *vec)
599 600
  {
    FUNCNAME("PetscSolverFeti::fillPetscMatrix()");   
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 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732

    updateDofData();

    int nComponents = vec->getSize();

    MatCreateMPIAIJ(PETSC_COMM_WORLD,
		    nRankB * nComponents, nRankB * nComponents,
		    nOverallB * nComponents, nOverallB * nComponents,
		    100, PETSC_NULL, 100, PETSC_NULL, 
		    &mat_b_b);

    MatCreateMPIAIJ(PETSC_COMM_WORLD,
		    nRankPrimals * nComponents, nRankPrimals * nComponents,
		    nOverallPrimals * nComponents, nOverallPrimals * nComponents,
		    10, PETSC_NULL, 10, PETSC_NULL, 
		    &mat_primal_primal);

    MatCreateMPIAIJ(PETSC_COMM_WORLD,
		    nRankB * nComponents, nRankPrimals * nComponents,
		    nOverallB * nComponents, nOverallPrimals * nComponents,
		    100, PETSC_NULL, 100, PETSC_NULL, 
		    &mat_b_primal);

    MatCreateMPIAIJ(PETSC_COMM_WORLD,
		    nRankPrimals * nComponents, nRankB * nComponents,
		    nOverallPrimals * nComponents, nOverallB * nComponents,
		    100, PETSC_NULL, 100, PETSC_NULL, 
		    &mat_primal_b);

    using mtl::tag::row; using mtl::tag::nz; using mtl::begin; using mtl::end;
    namespace traits = mtl::traits;
    typedef DOFMatrix::base_matrix_type Matrix;

    typedef traits::range_generator<row, Matrix>::type cursor_type;
    typedef traits::range_generator<nz, cursor_type>::type icursor_type;

    vector<int> cols, colsOther;
    vector<double> values, valuesOther;
    cols.reserve(300);
    colsOther.reserve(300);
    values.reserve(300);
    valuesOther.reserve(300);

    for (int i = 0; i < nComponents; i++)
      for (int j = 0; j < nComponents; j++)
	if ((*mat)[i][j]) {
	  traits::col<Matrix>::type col((*mat)[i][j]->getBaseMatrix());
	  traits::const_value<Matrix>::type value((*mat)[i][j]->getBaseMatrix());

	  for (cursor_type cursor = begin<row>((*mat)[i][j]->getBaseMatrix()), 
		 cend = end<row>((*mat)[i][j]->getBaseMatrix()); cursor != cend; ++cursor) {
	    bool rowPrimal = primals.count(*cursor) != 0;

	    cols.clear();
	    values.clear();

	    colsOther.clear();
	    valuesOther.clear();

	    for (icursor_type icursor = begin<nz>(cursor), icend = end<nz>(cursor); 
		 icursor != icend; ++icursor) {
	      if (primals.count(col(*icursor)) != 0) {
		TEST_EXIT_DBG(globalPrimalIndex.count(col(*icursor)))
		  ("No global primal index for DOF %d!\n", col(*icursor));

		int colIndex = globalPrimalIndex[col(*icursor)] * nComponents + j;

		if (rowPrimal) {
		  cols.push_back(colIndex);
		  values.push_back(value(*icursor));
		} else {
		  colsOther.push_back(colIndex);
		  valuesOther.push_back(value(*icursor));
		}
	      } else {
		TEST_EXIT_DBG(globalIndexB.count(col(*icursor)))
		  ("No global B index for DOF %d!\n", col(*icursor));
	      
		int colIndex = globalIndexB[col(*icursor)] * nComponents + j;

		if (rowPrimal) {
		  colsOther.push_back(colIndex);
		  valuesOther.push_back(value(*icursor));
		} else {
		  cols.push_back(colIndex);
		  values.push_back(value(*icursor));
		}
	      }
	    }

	    if (rowPrimal) {
	      TEST_EXIT_DBG(globalPrimalIndex.count(*cursor))
		("Should not happen!\n");

 	      int rowIndex = globalPrimalIndex[*cursor] * nComponents + i;
 	      MatSetValues(mat_primal_primal, 1, &rowIndex, cols.size(),
 			   &(cols[0]), &(values[0]), ADD_VALUES);

   	      if (colsOther.size())
  		MatSetValues(mat_primal_b, 1, &rowIndex, colsOther.size(),
  			     &(colsOther[0]), &(valuesOther[0]), ADD_VALUES);
	    } else {
	      TEST_EXIT_DBG(globalIndexB.count(*cursor))
		("Should not happen!\n");

	      int rowIndex = globalIndexB[*cursor] * nComponents + i;
	      MatSetValues(mat_b_b, 1, &rowIndex, cols.size(),
			   &(cols[0]), &(values[0]), ADD_VALUES);

 	      if (colsOther.size())
 		MatSetValues(mat_b_primal, 1, &rowIndex, colsOther.size(),
 			     &(colsOther[0]), &(valuesOther[0]), ADD_VALUES);
	    }
	  } 
	}


    MatAssemblyBegin(mat_b_b, MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(mat_b_b, MAT_FINAL_ASSEMBLY);

    MatAssemblyBegin(mat_primal_primal, MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(mat_primal_primal, MAT_FINAL_ASSEMBLY);

    MatAssemblyBegin(mat_b_primal, MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(mat_b_primal, MAT_FINAL_ASSEMBLY);

    MatAssemblyBegin(mat_primal_b, MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(mat_primal_b, MAT_FINAL_ASSEMBLY);
	  

    // === ===

733 734 735
    VecCreate(PETSC_COMM_WORLD, &f_b);
    VecSetSizes(f_b, nRankB * nComponents, nOverallB * nComponents);
    VecSetType(f_b, VECMPI);
736

737 738
    VecCreate(PETSC_COMM_WORLD, &f_primal);
    VecSetSizes(f_primal, nRankPrimals * nComponents, 
739
		nOverallPrimals * nComponents);
740
    VecSetType(f_primal, VECMPI);
741 742 743 744 745 746 747 748 749 750 751
    
    for (int i = 0; i < nComponents; i++) {
      DOFVector<double>::Iterator dofIt(vec->getDOFVector(i), USED_DOFS);
      for (dofIt.reset(); !dofIt.end(); ++dofIt) {
	int index = dofIt.getDOFIndex();
	if (primals.count(index)) {
	  TEST_EXIT_DBG(globalPrimalIndex.count(index))
	    ("Should not happen!\n");

	  index = globalPrimalIndex[index] * nComponents + i;
	  double value = *dofIt;
752
	  VecSetValues(f_primal, 1, &index, &value, ADD_VALUES);
753 754 755 756 757 758
	} else {
	  TEST_EXIT_DBG(globalIndexB.count(index))
	    ("Should not happen!\n");

	  index = globalIndexB[index] * nComponents + i;
	  double value = *dofIt;
759
	  VecSetValues(f_b, 1, &index, &value, ADD_VALUES);
760 761 762 763
	}      
      }
    }

764 765
    VecAssemblyBegin(f_b);
    VecAssemblyEnd(f_b);
766

767 768
    VecAssemblyBegin(f_primal);
    VecAssemblyEnd(f_primal);
769 770 771 772


    createMatLagrange(nComponents);

773 774 775 776 777

    
    createSchurPrimalKsp(nComponents);

    createFetiKsp();
778 779 780 781 782 783
  }


  void PetscSolverFeti::solvePetscMatrix(SystemVector &vec, AdaptInfo *adaptInfo)
  {
    FUNCNAME("PetscSolverFeti::solvePetscMatrix()");
784

785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025
    int debug = 0;
    Parameters::get("parallel->debug feti", debug);

    if (debug) {
      int nComponents = vec.getSize();

      Mat mat_lagrange_transpose;
      MatTranspose(mat_lagrange, MAT_INITIAL_MATRIX, &mat_lagrange_transpose);

      Mat A;
      Mat nestedA[3][3];
      nestedA[0][0] = mat_b_b;
      nestedA[0][1] = mat_b_primal;
      nestedA[0][2] = mat_lagrange_transpose;
      nestedA[1][0] = mat_primal_b;
      nestedA[1][1] = mat_primal_primal;
      nestedA[1][2] = PETSC_NULL;
      nestedA[2][0] = mat_lagrange;
      nestedA[2][1] = PETSC_NULL;
      nestedA[2][2] = PETSC_NULL;

      MatCreateNest(PETSC_COMM_WORLD, 3, PETSC_NULL, 3, PETSC_NULL, &(nestedA[0][0]), &A);

      MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
      MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
  


      int nRankNest = (nRankB + nRankPrimals) * nComponents + nRankLagrange;    
      int nOverallNest = (nOverallB + nOverallPrimals) * nComponents + nOverallLagrange;
      int rStartNest = (rStartB + rStartPrimals) * nComponents + rStartLagrange;

      {
	int matRow, matCol;
	MatGetLocalSize(A, &matRow, &matCol);
	TEST_EXIT_DBG(matRow == nRankNest)("Should not happen!\n");
	mpi::globalAdd(matRow);
	TEST_EXIT_DBG(matRow == nOverallNest)("Should not happen!\n");

	MatGetOwnershipRange(A, &matRow, &matCol);
	TEST_EXIT_DBG(matRow == rStartNest)("Should not happen!\n");
      }

      Vec f;
      VecCreate(PETSC_COMM_WORLD, &f);
      VecSetSizes(f, nRankNest, nOverallNest);
      VecSetType(f, VECMPI);

      Vec b;
      VecDuplicate(f, &b);

      PetscScalar *local_f_b;
      VecGetArray(f_b, &local_f_b);

      {
	int size;
	VecGetLocalSize(f_b, &size);
	TEST_EXIT_DBG(size == nRankB * nComponents)("Should not happen!\n");
      }

      PetscScalar *local_f_primal;
      VecGetArray(f_primal, &local_f_primal);

      {
	int size;
	VecGetLocalSize(f_primal, &size);
	TEST_EXIT_DBG(size == nRankPrimals * nComponents)("Should not happen!\n");
      }

      PetscScalar *local_f;
      VecGetArray(f, &local_f);

      int index = 0;
      for (int i = 0; i < nRankB * nComponents; i++)
	local_f[index++] = local_f_b[i];
      for (int i = 0; i < nRankPrimals * nComponents; i++)
	local_f[index++] = local_f_primal[i];

      VecRestoreArray(f, &local_f);  
      VecRestoreArray(f_b, &local_f_b);
      VecRestoreArray(f_primal, &local_f_primal);


      KSP ksp;
      KSPCreate(PETSC_COMM_WORLD, &ksp);
      KSPSetOperators(ksp, A, A, SAME_NONZERO_PATTERN);
      KSPSetFromOptions(ksp);



      KSPSolve(ksp, f, b);




      Vec u_b, u_primal;
      VecDuplicate(f_b, &u_b);
      VecDuplicate(f_primal, &u_primal);
    

      PetscScalar *local_b;
      VecGetArray(b, &local_b);

      PetscScalar *local_u_b;
      VecGetArray(u_b, &local_u_b);

      PetscScalar *local_u_primal;
      VecGetArray(u_primal, &local_u_primal);

      index = 0;
      for (int i = 0; i < nRankB * nComponents; i++)
	local_u_b[i] = local_b[index++];
      for (int i = 0; i < nRankPrimals * nComponents; i++)
	local_u_primal[i] = local_b[index++];
      for (int i = 0; i < nRankLagrange; i++) {
	MSG("MY %d-ith Lagrange: %f\n", i, local_b[index++]);
      }

      VecRestoreArray(f, &local_b);
      VecRestoreArray(u_b, &local_u_b);
      VecRestoreArray(u_primal, &local_u_primal);


      recoverSolution(u_b, u_primal, vec);

      VecDestroy(u_b);
      VecDestroy(u_primal);
      VecDestroy(b);
      VecDestroy(f);

      KSPDestroy(ksp);

    } else {

      KSPCreate(PETSC_COMM_WORLD, &ksp_b);
      KSPSetOperators(ksp_b, mat_b_b, mat_b_b, SAME_NONZERO_PATTERN);
      KSPSetOptionsPrefix(ksp_b, "solver_b_");
      KSPSetFromOptions(ksp_b);

    
      Vec vec_rhs;

      Vec tmp_b0, tmp_b1, tmp_lagrange0, tmp_primal0, tmp_primal1;
      MatGetVecs(mat_lagrange, PETSC_NULL, &tmp_lagrange0);
      MatGetVecs(mat_lagrange, PETSC_NULL, &vec_rhs);
      MatGetVecs(mat_b_b, PETSC_NULL, &tmp_b0);
      MatGetVecs(mat_b_b, PETSC_NULL, &tmp_b1);
      MatGetVecs(mat_primal_primal, PETSC_NULL, &tmp_primal0);
      MatGetVecs(mat_primal_primal, PETSC_NULL, &tmp_primal1);


      // === Create new rhs ===

      // vec_rhs = L * inv(K_BB) * f_b
      KSPSolve(ksp_b, f_b, tmp_b0);
      MatMult(mat_lagrange, tmp_b0, vec_rhs);

      // tmp_primal0 = M_PiB * inv(K_BB) * f_b
      MatMult(mat_primal_b, tmp_b0, tmp_primal0);

      // tmp_primal0 = f_Pi - M_PiB * inv(K_BB) * f_b
      //      VecAXPBY(tmp_primal0, 1.0, -1.0, f_primal);
      VecAXPBY(tmp_primal0, -1.0, 1.0, f_primal);

      // tmp_primal0 = inv(S_PiPi) (f_Pi - M_PiB * inv(K_BB) * f_b)
      KSPSolve(ksp_schur_primal, tmp_primal0, tmp_primal0);

      //
      MatMult(mat_b_primal, tmp_primal0, tmp_b0);
      KSPSolve(ksp_b, tmp_b0, tmp_b0);
      MatMult(mat_lagrange, tmp_b0, tmp_lagrange0);

      //
      VecAXPBY(vec_rhs, 1.0, 1.0, tmp_lagrange0);


      // === Solve with FETI-DP operator. ===

      KSPSolve(ksp_feti, vec_rhs, vec_rhs);

   
      // === Solve for u_primals. ===

      VecCopy(f_primal, tmp_primal0);

      KSPSolve(ksp_b, f_b, tmp_b0);
      MatMult(mat_primal_b, tmp_b0, tmp_primal1);

      VecAXPBY(tmp_primal0, -1.0, 1.0, tmp_primal1);

      MatMultTranspose(mat_lagrange, vec_rhs, tmp_b0);
      KSPSolve(ksp_b, tmp_b0, tmp_b0);
      MatMult(mat_primal_b, tmp_b0, tmp_primal1);

      VecAXPBY(tmp_primal0, 1.0, 1.0, tmp_primal1);
      KSPSolve(ksp_schur_primal, tmp_primal0, tmp_primal0);

    
      // === Solve for u_b. ===

      VecCopy(f_b, tmp_b0);
      MatMultTranspose(mat_lagrange, vec_rhs, tmp_b1);
      VecAXPBY(tmp_b0, -1.0, 1.0, tmp_b1);

      MatMult(mat_b_primal, tmp_primal0, tmp_b1);
      VecAXPBY(tmp_b0, -1.0, 1.0, tmp_b1);

      KSPSolve(ksp_b, tmp_b0, tmp_b0);


      // === And recover AMDiS solution vectors. ===
    
      recoverSolution(tmp_b0, tmp_primal0, vec);


      // === Destroy all data structures. ===
    

      VecDestroy(vec_rhs);
      VecDestroy(tmp_b0);
      VecDestroy(tmp_b1);
      VecDestroy(tmp_lagrange0);
      VecDestroy(tmp_primal0);
      VecDestroy(tmp_primal1);
	    

      KSPDestroy(ksp_b);

      MatDestroy(mat_b_b);
      MatDestroy(mat_primal_primal);
      MatDestroy(mat_b_primal);
      MatDestroy(mat_primal_b);
      MatDestroy(mat_lagrange);

      VecDestroy(f_b);
      VecDestroy(f_primal);

      destroySchurPrimalKsp();

      destroyFetiKsp();   
    }
1026

1027 1028 1029 1030
  }
#endif

}