GlobalMatrixSolver.cc 18.2 KB
Newer Older
1
2
3
4
5
#include "GlobalMatrixSolver.h"
#include "DOFVector.h"
#include "Debug.h"
#include "SystemVector.h"
#include "parallel/StdMpi.h"
6
#include "parallel/ParallelDomainDbg.h"
7
8
9
10
11
12
13

#include "petscksp.h"

namespace AMDiS {

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

    return 0;
  }
 
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37

  void GlobalMatrixSolver::addToMeshDistributor(MeshDistributor& m)
  {
    meshDistributor = &m;
    m.addProblemStat(this);
  }


  void GlobalMatrixSolver::buildAfterCoarsen(AdaptInfo *adaptInfo, Flag flag,
					     bool assembleMatrix,
					     bool assembleVector)
  {
    meshDistributor->checkMeshChange();
    ProblemVec::buildAfterCoarsen(adaptInfo, flag, assembleMatrix, assembleVector);
  }


  void GlobalMatrixSolver::solve(AdaptInfo *adaptInfo, bool fixedMatrix)
38
39
40
  {
    FUNCNAME("GlobalMatrixSolver::solve()");

41
42
    TEST_EXIT(meshDistributor)("Should not happen!\n");

43
44
45
46
47
#ifdef _OPENMP
    double wtime = omp_get_wtime();
#endif
    clock_t first = clock();

48
49
    fillPetscMatrix(systemMatrix, rhs);
    solvePetscMatrix(*solution);
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
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93

#ifdef _OPENMP
    INFO(info, 8)("solution of discrete system needed %.5f seconds system time / %.5f seconds wallclock time\n",
		   TIME_USED(first, clock()),
		   omp_get_wtime() - wtime);
#else
    INFO(info, 8)("solution of discrete system needed %.5f seconds\n",
		   TIME_USED(first, clock()));
#endif    
  }


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

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

    using mtl::tag::row; using mtl::tag::nz; using mtl::begin; using mtl::end;
    namespace traits= mtl::traits;
    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;

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

    // === 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) {

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

      // Global index of the current row dof.
94
      DegreeOfFreedom globalRowDof = meshDistributor->mapLocalToGlobal(*cursor);
95
      // Test if the current row dof is a periodic dof.
96
      bool periodicRow = (meshDistributor->getPeriodicDofMap().count(globalRowDof) > 0);
97
98
99
100
101
102
103
104
105
106
107
108


      // === Traverse all non zero entries of the row and produce vector cols ===
      // === with the column indices of all row entries and vector values     ===
      // === with the corresponding values.                                   ===

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

	// Set only non null values.
	if (value(*icursor) != 0.0) {
	  // Global index of the current column index.
109
	  int globalColDof = meshDistributor->mapLocalToGlobal(col(*icursor));
110
111
112
113
114
	  // Calculate the exact position of the column index in the petsc matrix.
	  int colIndex = globalColDof * dispMult + dispAddCol;

	  // If the current row is not periodic, but the current dof index is periodic,
	  // we have to duplicate the value to the other corresponding periodic columns.
115
116
 	  if (periodicRow == false && 
	      meshDistributor->getPeriodicDofMap().count(globalColDof) > 0) {
117
118
	    // The value is assign to n matrix entries, therefore, every entry 
	    // has only 1/n value of the original entry.
119
120
	    double scalFactor = 
	      1.0 / (meshDistributor->getPeriodicDof(globalColDof).size() + 1.0);
121
122
123
124
125
126
127

	    // Insert original entry.
 	    cols.push_back(colIndex);
 	    values.push_back(value(*icursor) * scalFactor);

	    // Insert the periodic entries.
 	    for (std::set<DegreeOfFreedom>::iterator it = 
128
129
		   meshDistributor->getPeriodicDof(globalColDof).begin();
 		 it != meshDistributor->getPeriodicDof(globalColDof).end(); ++it) {
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
 	      cols.push_back(*it * dispMult + dispAddCol);
 	      values.push_back(value(*icursor) * scalFactor);
	    }
 	  } else {
	    // Neigher row nor column dof index is periodic, simple add entry.
	    cols.push_back(colIndex);
	    values.push_back(value(*icursor));
	  }
	}
      }


      // === Up to now we have assembled on row. Now, the row must be send to the ===
      // === corresponding rows to the petsc matrix.                              ===

      // Calculate petsc row index.
      int rowIndex = globalRowDof * dispMult + dispAddRow;
      
      if (periodicRow) {
	// The row dof is periodic, so send dof to all the corresponding rows.

151
152
	double scalFactor = 
	  1.0 / (meshDistributor->getPeriodicDof(globalRowDof).size() + 1.0);
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
	
	int diagIndex = -1;
	for (int i = 0; i < static_cast<int>(values.size()); i++) {
	  // Change only the non diagonal values in the col. For the diagonal test
	  // we compare the global dof indices of the dof matrix (not of the petsc
	  // matrix!).
	  if ((cols[i] - dispAddCol) / dispMult != globalRowDof)
	    values[i] *= scalFactor;
	  else
	    diagIndex = i;
	}
	
	// Send the main row to the petsc matrix.
	MatSetValues(petscMatrix, 1, &rowIndex, cols.size(), 
		     &(cols[0]), &(values[0]), ADD_VALUES);	
 
	// Set diagonal element to zero, i.e., the diagonal element of the current
	// row is not send to the periodic row indices.
	if (diagIndex != -1)
	  values[diagIndex] = 0.0;

	// Send the row to all periodic row indices.
175
176
	for (std::set<DegreeOfFreedom>::iterator it = meshDistributor->getPeriodicDof(globalRowDof).begin();
	     it != meshDistributor->getPeriodicDof(globalRowDof).end(); ++it) {
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
	  int perRowIndex = *it * dispMult + dispAddRow;
	  MatSetValues(petscMatrix, 1, &perRowIndex, cols.size(), 
		       &(cols[0]), &(values[0]), ADD_VALUES);
	}

      } else {
	// The row dof is not periodic, simply send the row to the petsc matrix.
	MatSetValues(petscMatrix, 1, &rowIndex, cols.size(), 
		     &(cols[0]), &(values[0]), ADD_VALUES);
      }    
    }
  }


  void GlobalMatrixSolver::setDofVector(Vec& petscVec, DOFVector<double>* vec, 
					int dispMult, int dispAdd)
  {
194
195
    FUNCNAME("GlobalMatrixSolver::setDofVector()");

196
197
198
199
    // Traverse all used dofs in the dof vector.
    DOFVector<double>::Iterator dofIt(vec, USED_DOFS);
    for (dofIt.reset(); !dofIt.end(); ++dofIt) {
      // Calculate global row index of the dof.
200
201
      DegreeOfFreedom globalRow = 
	meshDistributor->mapLocalToGlobal(dofIt.getDOFIndex());
202
203
204
      // Calculate petsc index of the row dof.
      int index = globalRow * dispMult + dispAdd;

205
      if (meshDistributor->getPeriodicDofMap().count(globalRow) > 0) {
206
207
	// The dof index is periodic, so devide the value to all dof entries.

208
	double value = *dofIt / (meshDistributor->getPeriodicDof(globalRow).size() + 1.0);
209
210
	VecSetValues(petscVec, 1, &index, &value, ADD_VALUES);

211
212
	for (std::set<DegreeOfFreedom>::iterator it = meshDistributor->getPeriodicDof(globalRow).begin();
	     it != meshDistributor->getPeriodicDof(globalRow).end(); ++it) {
213
214
215
216
217
218
219
220
221
	  index = *it * dispMult + dispAdd;
	  VecSetValues(petscVec, 1, &index, &value, ADD_VALUES);
	}

      } else {
	// The dof index is not periodic.
	double value = *dofIt;
	VecSetValues(petscVec, 1, &index, &value, ADD_VALUES);
      }
222
    }
223
224
225
226
227
228
229
230
231
232
  }


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

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

233
    int nRankRows = meshDistributor->getNumberRankDofs() * nComponents;
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
    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 std::vector<std::pair<int, int> > MatrixNnzEntry;

    // Stores to each rank a list of nnz entries (i.e. pairs of row and column index)
    // that this rank will send to. This 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.
    std::map<int, MatrixNnzEntry> sendMatrixEntry;

    for (int i = 0; i < nComponents; i++) {
      for (int j = 0; j < nComponents; j++) {
 	if ((*mat)[i][j]) {
	  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) {

	    // Map the local row number to the global DOF index and create from it
	    // the global PETSc row index of this DOF.
268
269
	    int petscRowIdx = 
	      meshDistributor->mapLocalToGlobal(*cursor) * nComponents + i;
270

271
	    if (meshDistributor->getIsRankDof(*cursor)) {
272
273
274
275
276

	      // === 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.
277
278
	      int localPetscRowIdx = 
		petscRowIdx - meshDistributor->getRstart() * nComponents;
279
280
281

#if (DEBUG != 0)    
	      if (localPetscRowIdx < 0 || localPetscRowIdx >= nRankRows) {
282
		std::cout << "ERROR in rank: " << meshDistributor->getMpiRank() << std::endl;
283
		std::cout << "  Wrong r = " << localPetscRowIdx << " " << *cursor 
284
			  << " " << meshDistributor->mapLocalToGlobal(*cursor) << " " 
285
286
287
288
289
290
291
292
293
			  << nRankRows << std::endl;
		ERROR_EXIT("Should not happen!\n");
	      }
#endif
	      
	      // Traverse all non zero entries in this row.
	      for (icursor_type icursor = begin<nz>(cursor), 
		     icend = end<nz>(cursor); icursor != icend; ++icursor) {
		if (value(*icursor) != 0.0) {
294
295
		  int petscColIdx = 
		    meshDistributor->mapLocalToGlobal(col(*icursor)) * nComponents + j;
296
297
298

		  // 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.
299
300
		  if (petscColIdx >= meshDistributor->getRstart() * nComponents && 
		      petscColIdx < meshDistributor->getRstart() * nComponents + nRankRows)
301
302
303
304
305
306
		    d_nnz[localPetscRowIdx]++;
		  else
		    o_nnz[localPetscRowIdx]++;
		}    
	      }
	    } else {
307
308
	      typedef std::map<int, DofContainer> RankToDofContainer;

309
310
311
312
313
314
315
	      // === 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.                          ===

	      // Find out who is the member of this DOF.
	      int sendToRank = -1;
316
317
	      for (RankToDofContainer::iterator it = meshDistributor->getRecvDofs().begin();
		   it != meshDistributor->getRecvDofs().end(); ++it) {
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
		for (DofContainer::iterator dofIt = it->second.begin();
		     dofIt != it->second.end(); ++dofIt) {
		  if (**dofIt == *cursor) {
		    sendToRank = it->first;
		    break;
		  }
		}

		if (sendToRank != -1)
		  break;
	      }

	      TEST_EXIT_DBG(sendToRank != -1)("Should not happen!\n");

	      // Send all non zero entries to the member of the row DOF.
	      for (icursor_type icursor = begin<nz>(cursor), 
		     icend = end<nz>(cursor); icursor != icend; ++icursor) {
		if (value(*icursor) != 0.0) {
336
337
		  int petscColIdx = 
		    meshDistributor->mapLocalToGlobal(col(*icursor)) * nComponents + j;
338
339
340
341
342
343
344
345
346
347
348
349
350
351
		  
		  sendMatrixEntry[sendToRank].
		    push_back(std::make_pair(petscRowIdx, petscColIdx));
		}
	      }

	    } // if (isRankDof[*cursor]) ... else ...
	  } // for each row in mat[i][j]
	} // if mat[i][j]
      } 
    }

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

352
    StdMpi<MatrixNnzEntry> stdMpi(meshDistributor->getMpiComm(), true);
353
    stdMpi.send(sendMatrixEntry);
354
    stdMpi.recv(meshDistributor->getSendDofs());
355
356
357
358
359
360
361
362
363
364
365
366
    stdMpi.startCommunication<int>(MPI_INT);

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

    for (std::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;

367
	  int localRowIdx = r - meshDistributor->getRstart() * nComponents;
368
369
370
371
372

	  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);
	  
373
374
	  if (c < meshDistributor->getRstart() * nComponents || 
	      c >= meshDistributor->getRstart() * nComponents + nRankRows)
375
376
377
378
379
380
381
382
383
384
385
386
387
388
	    o_nnz[localRowIdx]++;
	  else
	    d_nnz[localRowIdx]++;
	}
      }
    }  
  }


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

    clock_t first = clock();
389
390
    int nRankRows = meshDistributor->getNumberRankDofs() * nComponents;
    int nOverallRows = meshDistributor->getNumberOverallDofs() * nComponents;
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405

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

406
407
408
409
410
    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) {
411
412
      if (d_nnz) {
	delete [] d_nnz;
413
	d_nnz = NULL;
414
	delete [] o_nnz;
415
	o_nnz = NULL;
416
417
      }

418
      createPetscNnzStructure(mat);
419
      lastMeshNnz = meshDistributor->getLastMeshChangeIndex();
420
    }
421
422
423
424
425
426
427
428
429
430
431

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

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

    INFO(info, 8)("Fill petsc matrix 1 needed %.5f seconds\n", TIME_USED(first, clock()));

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

438

439
440
441
442
443
    // === 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])
444
	  setDofMatrix((*mat)[i][j], nComponents, i, j);	
445
446
447
448
449
450

    INFO(info, 8)("Fill petsc matrix 2 needed %.5f seconds\n", TIME_USED(first, clock()));

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

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

    INFO(info, 8)("Fill petsc matrix needed %.5f seconds\n", TIME_USED(first, clock()));
  }


  void GlobalMatrixSolver::solvePetscMatrix(SystemVector &vec)
  {
    FUNCNAME("GlobalMatrixSolver::solvePetscMatrix()");

#if 0
    // Set old solution to be initiual guess for petsc solver.
    for (int i = 0; i < nComponents; i++)
      setDofVector(petscSolVec, vec->getDOFVector(i), nComponents, i);

    VecAssemblyBegin(petscSolVec);
    VecAssemblyEnd(petscSolVec);
#endif

    // === Init Petsc solver. ===

    KSP solver;
    KSPCreate(PETSC_COMM_WORLD, &solver);
    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);
    // Do not delete the solution vector, use it for the initial guess.
    //    KSPSetInitialGuessNonzero(solver, PETSC_TRUE);

489
490
491
492
#if (DEBUG != 0)
    //    ParallelDomainDbg::writeCoordsFile(*meshDistributor, "mpi-coords", "dat");
#endif

493
494
495
496
497

    // === Run Petsc. ===

    KSPSolve(solver, petscRhsVec, petscSolVec);

498

499
    // === Transfere values from Petsc's solution vectors to the dof vectors.
500

501
502
503
    PetscScalar *vecPointer;
    VecGetArray(petscSolVec, &vecPointer);

504
    int nRankDofs = meshDistributor->getNumberRankDofs();
505
    for (int i = 0; i < nComponents; i++) {
506
      DOFVector<double> &dofvec = *(vec.getDOFVector(i));
507
      for (int j = 0; j < nRankDofs; j++)
508
	dofvec[meshDistributor->mapLocalToDofIndex(j)] = 
509
	  vecPointer[j * nComponents + i]; 
510
511
512
513
514
515
516
    }

    VecRestoreArray(petscSolVec, &vecPointer);


    // === Synchronize dofs at common dofs, i.e., dofs that correspond to more ===
    // === than one partition.                                                 ===
517
    meshDistributor->synchVector(vec);
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


    // === Print information about solution process. ===

    int iterations = 0;
    KSPGetIterationNumber(solver, &iterations);
    MSG("  Number of iterations: %d\n", iterations);
    
    double norm = 0.0;
    MatMult(petscMatrix, petscSolVec, petscTmpVec);
    VecAXPY(petscTmpVec, -1.0, petscRhsVec);
    VecNorm(petscTmpVec, NORM_2, &norm);
    MSG("  Residual norm: %e\n", norm);


    // === Destroy Petsc's variables. ===

    MatDestroy(petscMatrix);
    VecDestroy(petscRhsVec);
    VecDestroy(petscSolVec);
    VecDestroy(petscTmpVec);
    KSPDestroy(solver);
  }

}