PetscSolverGlobalMatrix.cc 21.7 KB
Newer Older
Thomas Witkowski's avatar
Thomas Witkowski committed
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 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 73 74
//
// 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/PetscSolverGlobalMatrix.h"
#include "parallel/StdMpi.h"
#include "parallel/MpiHelper.h"

namespace AMDiS {

  PetscErrorCode myKSPMonitor(KSP ksp, PetscInt iter, PetscReal rnorm, void *)
  {    
    if (iter % 100 == 0 && MPI::COMM_WORLD.Get_rank() == 0)
      cout << "[0]  Petsc-Iteration " << iter << ": " << rnorm << endl;

    return 0;
  }


  void PetscSolverGlobalMatrix::fillPetscMatrix(Matrix<DOFMatrix*> *mat, SystemVector *vec)
  {
    FUNCNAME("PetscSolverGlobalMatrix::fillPetscMatrix()");

    TEST_EXIT_DBG(meshDistributor)("No mesh distributor object defined!\n");
    TEST_EXIT_DBG(mat)("No DOF matrix defined!\n");
    TEST_EXIT_DBG(vec)("NO DOF vector defined!\n");

    double wtime = MPI::Wtime();
    int nComponents = mat->getNumRows();
    int nRankRows = meshDistributor->getNumberRankDofs() * nComponents;
    int nOverallRows = meshDistributor->getNumberOverallDofs() * nComponents;

    // === Create PETSc vector (rhs, solution and a temporary vector). ===

    VecCreate(PETSC_COMM_WORLD, &petscRhsVec);
    VecSetSizes(petscRhsVec, nRankRows, nOverallRows);
    VecSetType(petscRhsVec, VECMPI);

    VecCreate(PETSC_COMM_WORLD, &petscSolVec);
    VecSetSizes(petscSolVec, nRankRows, nOverallRows);
    VecSetType(petscSolVec, VECMPI);

    VecCreate(PETSC_COMM_WORLD, &petscTmpVec);
    VecSetSizes(petscTmpVec, nRankRows, nOverallRows);
    VecSetType(petscTmpVec, VECMPI);

    int recvAllValues = 0;
    int sendValue = static_cast<int>(meshDistributor->getLastMeshChangeIndex() != lastMeshNnz);
    meshDistributor->getMpiComm().Allreduce(&sendValue, &recvAllValues, 1, MPI_INT, MPI_SUM);

    if (!d_nnz || recvAllValues != 0) {
      if (d_nnz) {
	delete [] d_nnz;
	d_nnz = NULL;
	delete [] o_nnz;
	o_nnz = NULL;
      }

      createPetscNnzStructure(mat);
      lastMeshNnz = meshDistributor->getLastMeshChangeIndex();
    }


    // === Create PETSc matrix with the computed nnz data structure. ===

    MatCreateMPIAIJ(PETSC_COMM_WORLD, nRankRows, nRankRows, nOverallRows, nOverallRows,
75 76
		     0, d_nnz, 0, o_nnz, &petscMatrix);

Thomas Witkowski's avatar
Thomas Witkowski committed
77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 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 123 124 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
#if (DEBUG != 0)
    MSG("Fill petsc matrix 1 needed %.5f seconds\n", MPI::Wtime() - wtime);
#endif

#if (DEBUG != 0)
    int a, b;
    MatGetOwnershipRange(petscMatrix, &a, &b);
    TEST_EXIT(a == meshDistributor->getRstart() * nComponents)
      ("Wrong matrix ownership range!\n");
    TEST_EXIT(b == meshDistributor->getRstart() * nComponents + nRankRows)
      ("Wrong matrix ownership range!\n");
#endif


    // === Transfer values from DOF matrices to the PETSc matrix. === 

    for (int i = 0; i < nComponents; i++)
      for (int j = 0; j < nComponents; j++)
	if ((*mat)[i][j])
	  setDofMatrix((*mat)[i][j], nComponents, i, j);

#if (DEBUG != 0)
    MSG("Fill petsc matrix 2 needed %.5f seconds\n", MPI::Wtime() - wtime);
#endif

    MatAssemblyBegin(petscMatrix, MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(petscMatrix, MAT_FINAL_ASSEMBLY);

#if (DEBUG != 0)
    MSG("Fill petsc matrix 3 needed %.5f seconds\n", 
	TIME_USED(MPI::Wtime(), wtime));
#endif

    // === Transfer values from DOF vector to the PETSc vector. === 

    for (int i = 0; i < nComponents; i++)
      setDofVector(petscRhsVec, vec->getDOFVector(i), nComponents, i);

    VecAssemblyBegin(petscRhsVec);
    VecAssemblyEnd(petscRhsVec);

    MSG("Fill petsc matrix needed %.5f seconds\n", MPI::Wtime() - wtime);
  }


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

    int nComponents = vec.getSize();

    // === Set old solution to be initiual guess for PETSc solver. ===
    if (!zeroStartVector) {
      VecSet(petscSolVec, 0.0);
      
      for (int i = 0; i < nComponents; i++)
	setDofVector(petscSolVec, vec.getDOFVector(i), nComponents, i, true);
      
      VecAssemblyBegin(petscSolVec);
      VecAssemblyEnd(petscSolVec);
    }

    // === Init PETSc solver. ===
    KSPCreate(PETSC_COMM_WORLD, &solver);
    KSPGetPC(solver, &pc);
    KSPSetOperators(solver, petscMatrix, petscMatrix, SAME_NONZERO_PATTERN); 
    KSPSetTolerances(solver, 0.0, 1e-8, PETSC_DEFAULT, PETSC_DEFAULT);
    KSPSetType(solver, KSPBCGS);
    KSPMonitorSet(solver, myKSPMonitor, PETSC_NULL, 0);
    KSPSetFromOptions(solver);
    PCSetFromOptions(pc);

    // Do not delete the solution vector, use it for the initial guess.
    if (!zeroStartVector)
      KSPSetInitialGuessNonzero(solver, PETSC_TRUE);

    // PETSc.
    KSPSolve(solver, petscRhsVec, petscSolVec);

    // === Transfere values from PETSc's solution vectors to the DOF vectors. ===
    int nRankDofs = meshDistributor->getNumberRankDofs();
    PetscScalar *vecPointer;
    VecGetArray(petscSolVec, &vecPointer);

    for (int i = 0; i < nComponents; i++) {
      DOFVector<double> &dofvec = *(vec.getDOFVector(i));
      for (int j = 0; j < nRankDofs; j++)
	dofvec[meshDistributor->mapLocalToDofIndex(j)] = 
	  vecPointer[j * nComponents + i]; 
    }

    VecRestoreArray(petscSolVec, &vecPointer);


    // === Synchronize DOFs at common DOFs, i.e., DOFs that correspond to ===
    // === more than one partition.                                       ===
    meshDistributor->synchVector(vec);


    // Print iteration counter and residual norm of the solution.
    printSolutionInfo(adaptInfo);


    // === Destroy PETSc's variables. ===

182 183 184 185 186 187 188
#ifdef HAVE_PETSC_DEV 
    MatDestroy(&petscMatrix);
    VecDestroy(&petscRhsVec);
    VecDestroy(&petscSolVec);
    VecDestroy(&petscTmpVec);
    KSPDestroy(&solver);
#else
189 190 191 192 193
    MatDestroy(&petscMatrix);
    VecDestroy(&petscRhsVec);
    VecDestroy(&petscSolVec);
    VecDestroy(&petscTmpVec);
    KSPDestroy(&solver);
194
#endif
Thomas Witkowski's avatar
Thomas Witkowski committed
195 196 197 198 199 200 201 202 203 204 205
  }


  void PetscSolverGlobalMatrix::setDofMatrix(DOFMatrix* mat, int dispMult, 
				      int dispAddRow, int dispAddCol)
  {
    FUNCNAME("PetscSolverGlobalMatrix::setDofMatrix()");

    TEST_EXIT(mat)("No DOFMatrix!\n");

    using mtl::tag::row; using mtl::tag::nz; using mtl::begin; using mtl::end;
206
    namespace traits = mtl::traits;
Thomas Witkowski's avatar
Thomas Witkowski committed
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 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 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 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
    typedef DOFMatrix::base_matrix_type Matrix;

    traits::col<Matrix>::type col(mat->getBaseMatrix());
    traits::const_value<Matrix>::type value(mat->getBaseMatrix());

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

    vector<int> cols;
    vector<double> values;
    cols.reserve(300);
    values.reserve(300);
    
    vector<int> globalCols;


    // === Traverse all rows of the dof matrix and insert row wise the values ===
    // === to the PETSc matrix.                                               ===

    for (cursor_type cursor = begin<row>(mat->getBaseMatrix()), 
	   cend = end<row>(mat->getBaseMatrix()); cursor != cend; ++cursor) {

      // Global index of the current row DOF.
      int globalRowDof = meshDistributor->mapLocalToGlobal(*cursor);
      // Test if the current row DOF is a periodic DOF.
      bool periodicRow = meshDistributor->isPeriodicDof(globalRowDof);

      if (!periodicRow) {
	// === Row DOF index is not periodic. ===

	// Calculate PETSc row index.

	int rowIndex = globalRowDof * dispMult + dispAddRow;

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

	for (icursor_type icursor = begin<nz>(cursor), icend = end<nz>(cursor); 
	     icursor != icend; ++icursor) {

	  // Global index of the current column index.
	  int globalColDof = meshDistributor->mapLocalToGlobal(col(*icursor));
	  // Test if the current col dof is a periodic dof.
	  bool periodicCol = meshDistributor->isPeriodicDof(globalColDof);
	  // Calculate PETSc col index.
	  int colIndex = globalColDof * dispMult + dispAddCol;

	  // Ignore all zero entries, expect it is a diagonal entry.
 	  if (value(*icursor) == 0.0 && rowIndex != colIndex)
 	    continue;

	  if (!periodicCol) {
	    // Calculate the exact position of the column index in the PETSc matrix.
	    cols.push_back(colIndex);
	    values.push_back(value(*icursor));
	  } else {
	    // === Row index is not periodic, but column index is. ===

	    // Create set of all periodic associations of the column index.
	    std::set<int> perAsc;
	    std::set<int>& perColAsc = 
	      meshDistributor->getPerDofAssociations(globalColDof);
	    for (std::set<int>::iterator it = perColAsc.begin(); 
		 it != perColAsc.end(); ++it)
	      if (*it >= -3)
		perAsc.insert(*it);
    
	    // Scale value to the number of periodic associations of the column index.
	    double scaledValue = 
	      value(*icursor) * pow(0.5, static_cast<double>(perAsc.size()));

	    
	    // === Create set of all matrix column indices due to the periodic ===
	    // === associations of the column DOF index.                       ===

	    vector<int> newCols;

	    // First, add the original matrix index.
	    newCols.push_back(globalColDof);
	    
	    // And add all periodic matrix indices.
	    for (std::set<int>::iterator it = perAsc.begin(); 
		 it != perAsc.end(); ++it) {
	      int nCols = static_cast<int>(newCols.size());

	      for (int i = 0; i < nCols; i++) {
 		TEST_EXIT_DBG(meshDistributor->isPeriodicDof(newCols[i], *it))
 		  ("Wrong periodic DOF associations at boundary %d with DOF %d!\n",
		   *it, newCols[i]);

		newCols.push_back(meshDistributor->getPeriodicMapping(newCols[i], *it));
	      }
	    }

	    for (unsigned int i = 0; i < newCols.size(); i++) {
	      cols.push_back(newCols[i] * dispMult + dispAddCol);
	      values.push_back(scaledValue);	      
	    }
	  }
	}

	MatSetValues(petscMatrix, 1, &rowIndex, cols.size(), 
		     &(cols[0]), &(values[0]), ADD_VALUES);	
      } else {
	// === Row DOF index is periodic. ===

	// Because this row is periodic, we will have to add the entries of this 
	// matrix row to multiple rows. The following maps store to each row an
	// array of column indices and values of the entries that must be added to
	// the PETSc matrix.
	map<int, vector<int> > colsMap;
	map<int, vector<double> > valsMap;

	// Traverse all column entries.
	for (icursor_type icursor = begin<nz>(cursor), icend = end<nz>(cursor); 
	     icursor != icend; ++icursor) {

	  // Global index of the current column index.
	  int globalColDof = meshDistributor->mapLocalToGlobal(col(*icursor));

	  // Ignore all zero entries, expect it is a diagonal entry.
 	  if (value(*icursor) == 0.0 && globalRowDof != globalColDof)
 	    continue;


	  // === Add all periodic associations of both, the row and the column ===
	  // === indices to the set perAsc.                                    ===

	  std::set<int> perAsc;

	  if (meshDistributor->isPeriodicDof(globalColDof)) {
	    std::set<int>& perColAsc = 
	      meshDistributor->getPerDofAssociations(globalColDof);
	    for (std::set<int>::iterator it = perColAsc.begin(); 
		 it != perColAsc.end(); ++it)
	      if (*it >= -3)
		perAsc.insert(*it);
	  }

	  std::set<int>& perRowAsc = 
	    meshDistributor->getPerDofAssociations(globalRowDof);
	  for (std::set<int>::iterator it = perRowAsc.begin(); 
	       it != perRowAsc.end(); ++it)
	    if (*it >= -3)
	      perAsc.insert(*it);

	  // Scale the value with respect to the number of periodic associations.
	  double scaledValue = 
	    value(*icursor) * pow(0.5, static_cast<double>(perAsc.size()));


	  // === Create all matrix entries with respect to the periodic  ===
	  // === associations of the row and column indices.             ===

	  vector<pair<int, int> > entry;
	  
	  // First, add the original entry.
	  entry.push_back(make_pair(globalRowDof, globalColDof));

	  // Then, traverse the periodic associations of the row and column indices
	  // and create the corresponding entries.
	  for (std::set<int>::iterator it = perAsc.begin(); it != perAsc.end(); ++it) {
	    int nEntry = static_cast<int>(entry.size());
	    for (int i = 0; i < nEntry; i++) {
	      int perRowDof = 0;
	      if (meshDistributor->getPeriodicMapping()[*it].count(entry[i].first))
		perRowDof = meshDistributor->getPeriodicMapping(entry[i].first, *it);
	      else
		perRowDof = entry[i].first;

	      int perColDof;
	      if (meshDistributor->getPeriodicMapping()[*it].count(entry[i].second))
		perColDof = meshDistributor->getPeriodicMapping(entry[i].second, *it);
	      else
		perColDof = entry[i].second;	      	      
	      

	      entry.push_back(make_pair(perRowDof, perColDof));
	    }
	  }


	  // === Translate the matrix entries to PETSc's matrix.

	  for (vector<pair<int, int> >::iterator eIt = entry.begin(); 
	       eIt != entry.end(); ++eIt) {
	    // Calculate row and column indices of the PETSc matrix.
	    int rowIndex = eIt->first * dispMult + dispAddRow;
	    int colIndex = eIt->second * dispMult + dispAddCol;

	    colsMap[rowIndex].push_back(colIndex);
	    valsMap[rowIndex].push_back(scaledValue);
	  }
	}


	// === Finally, add all periodic rows to the PETSc matrix. ===

	for (map<int, vector<int> >::iterator rowIt = colsMap.begin();
	     rowIt != colsMap.end(); ++rowIt) {
	  TEST_EXIT_DBG(rowIt->second.size() == valsMap[rowIt->first].size())
	    ("Should not happen!\n");

	  int rowIndex = rowIt->first;
	  MatSetValues(petscMatrix, 1, &rowIndex, rowIt->second.size(),
		       &(rowIt->second[0]), &(valsMap[rowIt->first][0]), ADD_VALUES);
	}
      }
    }
  }


419 420 421 422 423
  void PetscSolverGlobalMatrix::setDofVector(Vec& petscVec, 
					     DOFVector<double>* vec, 
					     int dispMult, 
					     int dispAdd, 
					     bool rankOnly)
Thomas Witkowski's avatar
Thomas Witkowski committed
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 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
  {
    FUNCNAME("PetscSolverGlobalMatrix::setDofVector()");

    // Traverse all used DOFs in the dof vector.
    DOFVector<double>::Iterator dofIt(vec, USED_DOFS);
    for (dofIt.reset(); !dofIt.end(); ++dofIt) {
      if (rankOnly && !meshDistributor->getIsRankDof(dofIt.getDOFIndex()))
	continue;

      // Calculate global row index of the DOF.
      DegreeOfFreedom globalRowDof = 
	meshDistributor->mapLocalToGlobal(dofIt.getDOFIndex());
      // Calculate PETSc index of the row DOF.
      int index = globalRowDof * dispMult + dispAdd;

      if (meshDistributor->isPeriodicDof(globalRowDof)) {
	std::set<int>& perAsc = meshDistributor->getPerDofAssociations(globalRowDof);
	double value = *dofIt / (perAsc.size() + 1.0);
	VecSetValues(petscVec, 1, &index, &value, ADD_VALUES);

	for (std::set<int>::iterator perIt = perAsc.begin(); perIt != perAsc.end(); ++perIt) {
	  int mappedDof = meshDistributor->getPeriodicMapping(globalRowDof, *perIt);
	  int mappedIndex = mappedDof * dispMult + dispAdd;
	  VecSetValues(petscVec, 1, &mappedIndex, &value, ADD_VALUES);
	}
      } else {
	// The DOF index is not periodic.
	double value = *dofIt;
	VecSetValues(petscVec, 1, &index, &value, ADD_VALUES);
      }
    }
  }


  void PetscSolverGlobalMatrix::createPetscNnzStructure(Matrix<DOFMatrix*> *mat)
  {
    FUNCNAME("PetscSolverGlobalMatrix::createPetscNnzStructure()");

    TEST_EXIT_DBG(!d_nnz)("There is something wrong!\n");
    TEST_EXIT_DBG(!o_nnz)("There is something wrong!\n");

    int nComponents = mat->getNumRows();
    int nRankRows = meshDistributor->getNumberRankDofs() * nComponents;
    d_nnz = new int[nRankRows];
    o_nnz = new int[nRankRows];
    for (int i = 0; i < nRankRows; i++) {
      d_nnz[i] = 0;
      o_nnz[i] = 0;
    }

    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 vector<pair<int, int> > MatrixNnzEntry;
    typedef map<int, DofContainer> RankToDofContainer;

    // Stores to each rank a list of nnz entries (i.e. pairs of row and column index)
    // that this rank will send to. These nnz entries will be assembled on this rank,
    // but because the row DOFs are not DOFs of this rank they will be send to the
    // owner of the row DOFs.
    map<int, MatrixNnzEntry> sendMatrixEntry;

    // Maps to each DOF that must be send to another rank the rank number of the
    // receiving rank.
    map<DegreeOfFreedom, int> sendDofToRank;


    // First, create for all ranks we send data to MatrixNnzEntry object with 0 entries.
    for (RankToDofContainer::iterator it = meshDistributor->getRecvDofs().begin();
	 it != meshDistributor->getRecvDofs().end(); ++it) {
      sendMatrixEntry[it->first].resize(0);

      for (DofContainer::iterator dofIt = it->second.begin();
	   dofIt != it->second.end(); ++dofIt)
	sendDofToRank[**dofIt] = it->first;
    }


    std::set<int> recvFromRank;
    for (RankToDofContainer::iterator it = meshDistributor->getSendDofs().begin();
	 it != meshDistributor->getSendDofs().end(); ++it)
      recvFromRank.insert(it->first);


    for (int i = 0; i < nComponents; i++) {
      for (int j = 0; j < nComponents; j++) {
 	if (!(*mat)[i][j])
	  continue;

	Matrix bmat = (*mat)[i][j]->getBaseMatrix();

	traits::col<Matrix>::type col(bmat);
	traits::const_value<Matrix>::type value(bmat);
	  
	typedef traits::range_generator<row, Matrix>::type cursor_type;
	typedef traits::range_generator<nz, cursor_type>::type icursor_type;
	
	for (cursor_type cursor = begin<row>(bmat), 
	       cend = end<row>(bmat); cursor != cend; ++cursor) {
	  
	  int globalRowDof = meshDistributor->mapLocalToGlobal(*cursor);

	  vector<int> rows;
	  rows.push_back(globalRowDof);
	  vector<int> rowRank;
	  if (meshDistributor->getIsRankDof(*cursor)) {
	    rowRank.push_back(meshDistributor->getMpiRank());
	  } else {
	    // Find out who is the member of this DOF.
	    TEST_EXIT_DBG(sendDofToRank.count(*cursor))("DOF %d has no receiver!\n", *cursor);
	    
	    rowRank.push_back(sendDofToRank[*cursor]);
	  }

	  // Map the local row number to the global DOF index and create from it
	  // the global PETSc row index of this DOF.
	  
	  int petscRowIdx = globalRowDof * nComponents + i;
	  
	  if (meshDistributor->getIsRankDof(*cursor)) {
	    	    
	    // === The current row DOF is a rank dof, so create the corresponding ===
	    // === nnz values directly on rank's nnz data.                        ===
	    
	    // This is the local row index of the local PETSc matrix.
	    int localPetscRowIdx = 
	      petscRowIdx - meshDistributor->getRstart() * nComponents;
	    
	    TEST_EXIT_DBG(localPetscRowIdx >= 0 && localPetscRowIdx < nRankRows)
	      ("Should not happen! \n Debug info: localRowIdx = %d   globalRowIndx = %d   petscRowIdx = %d   localPetscRowIdx = %d   rStart = %d   nCompontens = %d   nRankRows = %d\n",
	       *cursor, 
	       meshDistributor->mapLocalToGlobal(*cursor), 
	       petscRowIdx, 
	       localPetscRowIdx, 
	       meshDistributor->getRstart(), 
	       nComponents, 
	       nRankRows);
	    
	    
	    // Traverse all non zero entries in this row.
	    for (icursor_type icursor = begin<nz>(cursor), 
		   icend = end<nz>(cursor); icursor != icend; ++icursor) {
	      int petscColIdx = 
		meshDistributor->mapLocalToGlobal(col(*icursor)) * nComponents + j;
	      
	      if (value(*icursor) != 0.0 || petscRowIdx == petscColIdx) {
		// The row DOF is a rank DOF, if also the column is a rank DOF, 
		// increment the d_nnz values for this row, otherwise the o_nnz value.
		if (petscColIdx >= meshDistributor->getRstart() * nComponents && 
		    petscColIdx < meshDistributor->getRstart() * nComponents + nRankRows)
		  d_nnz[localPetscRowIdx]++;
		else
		  o_nnz[localPetscRowIdx]++;
	      }    
	    }
	  } else {
	    // === The current row DOF is not a rank dof, i.e., it will be created ===
	    // === on this rank, but after this it will be send to another rank    ===
	    // === matrix. So we need to send also the corresponding nnz structure ===
	    // === of this row to the corresponding rank.                          ===
	    
	    // Send all non zero entries to the member of the row DOF.
	    int sendToRank = sendDofToRank[*cursor];
	    
	    for (icursor_type icursor = begin<nz>(cursor), 
		   icend = end<nz>(cursor); icursor != icend; ++icursor) {
	      if (value(*icursor) != 0.0) {
		int petscColIdx = 
		  meshDistributor->mapLocalToGlobal(col(*icursor)) * nComponents + j;
		
		sendMatrixEntry[sendToRank].
		  push_back(make_pair(petscRowIdx, petscColIdx));
	      }
	    }
	    
	  } // if (isRankDof[*cursor]) ... else ...
	} // for each row in mat[i][j]
      } 
    }

    // === Send and recv the nnz row structure to/from other ranks. ===

    StdMpi<MatrixNnzEntry> stdMpi(meshDistributor->getMpiComm(), true);
    stdMpi.send(sendMatrixEntry);
    for (std::set<int>::iterator it = recvFromRank.begin(); 
	 it != recvFromRank.end(); ++it)
      stdMpi.recv(*it);
    stdMpi.startCommunication();


    // === Evaluate the nnz structure this rank got from other ranks and add it to ===
    // === the PETSc nnz data structure.                                           ===

    for (map<int, MatrixNnzEntry>::iterator it = stdMpi.getRecvData().begin();
	 it != stdMpi.getRecvData().end(); ++it) {
      if (it->second.size() > 0) {
	for (unsigned int i = 0; i < it->second.size(); i++) {
	  int r = it->second[i].first;
	  int c = it->second[i].second;

	  int localRowIdx = r - meshDistributor->getRstart() * nComponents;

	  TEST_EXIT_DBG(localRowIdx >= 0 && localRowIdx < nRankRows)
	    ("Got row index %d/%d (nRankRows = %d) from rank %d. Should not happen!\n",
	     r, localRowIdx, nRankRows, it->first);
	  
	  if (c < meshDistributor->getRstart() * nComponents || 
	      c >= meshDistributor->getRstart() * nComponents + nRankRows)
	    o_nnz[localRowIdx]++;
	  else
	    d_nnz[localRowIdx]++;
	}
      }
    }

    // The above algorithm for calculating the number of nnz per row over-
    // approximates the value, i.e., the number is always equal or larger to 
    // the real number of nnz values in the global parallel matrix. For small
    // matrices, the problem may arise, that the result is larger than the
    // number of elements in a row. This is fixed in the following.

    if (nRankRows < 100) 
      for (int i = 0; i < nRankRows; i++)
	d_nnz[i] = std::min(d_nnz[i], nRankRows);
  }

}