ParallelDomainProblem.cc 17.2 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
#include "ParallelDomainProblem.h"
#include "ProblemScal.h"
#include "ProblemInstat.h"
#include "ParMetisPartitioner.h"
#include "Mesh.h"
#include "Traverse.h"
#include "ElInfo.h"
#include "Element.h"
#include "MacroElement.h"
#include "PartitionElementData.h"
11
12
#include "DOFMatrix.h"
#include "DOFVector.h"
13
14
15
#include "VtkWriter.h"

#include "petscksp.h"
16
17
18
19

namespace AMDiS {

  ParallelDomainProblemBase::ParallelDomainProblemBase(const std::string& name,
20
21
						       ProblemIterationInterface *iIF,
						       ProblemTimeInterface *tIF,
22
						       FiniteElemSpace *fe)
23
24
    : iterationIF(iIF),
      timeIF(tIF),
25
26
      feSpace(fe),
      mesh(fe->getMesh()),
27
      initialPartitionMesh(true),
28
      nRankDOFs(0)
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
  {
    mpiRank = MPI::COMM_WORLD.Get_rank();
    mpiSize = MPI::COMM_WORLD.Get_size();
    mpiComm = MPI::COMM_WORLD;
    partitioner = new ParMetisPartitioner(mesh, &mpiComm);
  }

  void ParallelDomainProblemBase::initParallelization(AdaptInfo *adaptInfo)
  {
    if (mpiSize <= 1)
      return;

    // create an initial partitioning of the mesh
    partitioner->createPartitionData();
    // set the element weights, which are 1 at the very first begin
    setElemWeights(adaptInfo);
    // and now partition the mesh
    partitionMesh(adaptInfo);   


49
    // === Create new global and local DOF numbering. ===
50

51
52
53
    int nRankDOFs = 0;
    int nOverallDOFs = 0;
    createLocalGlobalNumbering(nRankDOFs, nOverallDOFs);
54

Thomas Witkowski's avatar
Thomas Witkowski committed
55
56
57

    // === Create interior boundary information ===

58
    createInteriorBoundaryInfo();
Thomas Witkowski's avatar
Thomas Witkowski committed
59

Thomas Witkowski's avatar
Thomas Witkowski committed
60

61
62
    // === Remove all macro elements that are not part of the rank partition. ===

63
    removeMacroElements();
64
      
65
66
67
68
69
70

    /// === Reset all DOFAdmins of the mesh. ===

    int nAdmins = mesh->getNumberOfDOFAdmin();
    for (int i = 0; i < nAdmins; i++) {
      for (int j = 0; j < mesh->getDOFAdmin(i).getSize(); j++)
71
72
73
 	const_cast<DOFAdmin&>(mesh->getDOFAdmin(i)).setDOFFree(j, true);
      for (int j = 0; j < static_cast<int>(mapLocalGlobalDOFs.size()); j++)
 	const_cast<DOFAdmin&>(mesh->getDOFAdmin(i)).setDOFFree(j, false);
Thomas Witkowski's avatar
Thomas Witkowski committed
74
75
76
    }


77
    /// === Create petsc matrix. ===
78

79
80
    int ierr;
    ierr = MatCreate(PETSC_COMM_WORLD, &petscMatrix);
81
    ierr = MatSetSizes(petscMatrix, nRankDOFs, nRankDOFs, nOverallDOFs, nOverallDOFs);
82
83
84
    ierr = MatSetType(petscMatrix, MATAIJ);

    ierr = VecCreate(PETSC_COMM_WORLD, &petscRhsVec);
85
    ierr = VecSetSizes(petscRhsVec, nRankDOFs, nOverallDOFs);
86
    ierr = VecSetType(petscRhsVec, VECMPI);
87
88

    ierr = VecCreate(PETSC_COMM_WORLD, &petscSolVec);
89
    ierr = VecSetSizes(petscSolVec, nRankDOFs, nOverallDOFs);
90
    ierr = VecSetType(petscSolVec, VECMPI);
91
92
93
  }

  void ParallelDomainProblemBase::exitParallelization(AdaptInfo *adaptInfo)
94
  {}
95

96

97
98
99
  void ParallelDomainProblemBase::fillPetscMatrix(DOFMatrix *mat, 
						  DOFVector<double> *vec)
  {
100
101
102
103
104
105
106
107
108
109
110
111
112
113
    using mtl::tag::major; using mtl::tag::nz; using mtl::begin; using mtl::end;
    namespace traits= mtl::traits;
    typedef DOFMatrix::base_matrix_type Matrix;

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

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

    for (cursor_type cursor = begin<major>(mat->getBaseMatrix()), cend = end<major>(mat->getBaseMatrix()); cursor != cend; ++cursor)
      for (icursor_type icursor = begin<nz>(cursor), icend = end<nz>(cursor); icursor != icend; ++icursor)
	if (value(*icursor) != 0.0) {
114
115
	  int r = mapLocalGlobalDOFs[row(*icursor)];
	  int c = mapLocalGlobalDOFs[col(*icursor)];
116
	  double v = value(*icursor);
117

118
	  MatSetValues(petscMatrix, 1, &r, 1, &c, &v, ADD_VALUES);
119
	}
120

121
122
123
124
125
126

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

    DOFVector<double>::Iterator dofIt(vec, USED_DOFS);
    for (dofIt.reset(); !dofIt.end(); ++dofIt) {
127
      int index = mapLocalGlobalDOFs[dofIt.getDOFIndex()];
128
129
130
      double value = *dofIt;

      VecSetValues(petscRhsVec, 1, &index, &value, ADD_VALUES);
131
132
133
    }
  }

134

135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
  void ParallelDomainProblemBase::solvePetscMatrix(DOFVector<double> *vec)
  {
    KSP ksp;
    PC pc;

    KSPCreate(PETSC_COMM_WORLD, &ksp);
    KSPSetOperators(ksp, petscMatrix, petscMatrix, DIFFERENT_NONZERO_PATTERN);
    KSPGetPC(ksp, &pc);
    PCSetType(pc, PCJACOBI);
    KSPSetTolerances(ksp, 1.e-7, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
    KSPSetType(ksp, KSPBCGS);
    KSPMonitorSet(ksp, KSPMonitorDefault, PETSC_NULL, 0);
    KSPSolve(ksp, petscRhsVec, petscSolVec);

    PetscScalar *vecPointer;
    VecGetArray(petscSolVec, &vecPointer);

    DOFVector<double>::Iterator dofIt(vec, USED_DOFS);
    int counter = 0;
    for (dofIt.reset(); !dofIt.end(); ++dofIt)
155
      *dofIt = vecPointer[counter++];
156
157
158

    VecRestoreArray(petscSolVec, &vecPointer);

159
160
161
162
163
164
165
166
167
    std::vector<double*> sendBuffers(sendDofs.size());
    std::vector<double*> recvBuffers(recvDofs.size());
    
    int i = 0;
    for (std::map<int, std::vector<DegreeOfFreedom> >::iterator sendIt = sendDofs.begin();
	 sendIt != sendDofs.end();
	 ++sendIt, i++) {
      sendBuffers[i] = new double[sendIt->second.size()];

168
      for (int j = 0; j < static_cast<int>(sendIt->second.size()); j++)
169
170
171
172
173
174
175
176
177
	sendBuffers[i][j] = (*vec)[(sendIt->second)[j]];

      mpiComm.Isend(sendBuffers[i], sendIt->second.size(), MPI_DOUBLE, sendIt->first, 0);
    }

    i = 0;
    for (std::map<int, std::vector<DegreeOfFreedom> >::iterator recvIt = recvDofs.begin();
	 recvIt != recvDofs.end();
	 ++recvIt, i++) {
178
      recvBuffers[i] = new double[recvIt->second.size()];      
179
180
181
182
183
184
185
186
187
188
189

      mpiComm.Irecv(recvBuffers[i], recvIt->second.size(), MPI_DOUBLE, recvIt->first, 0);
    }

    
    mpiComm.Barrier();
    
    i = 0;
    for (std::map<int, std::vector<DegreeOfFreedom> >::iterator recvIt = recvDofs.begin();
	 recvIt != recvDofs.end();
	 ++recvIt, i++) {
190
      for (int j = 0; j < static_cast<int>(recvIt->second.size()); j++)
191
192
193
194
195
	(*vec)[(recvIt->second)[j]] = recvBuffers[i][j];

      delete [] recvBuffers[i];
    }
    
196
    for (int i = 0; i < static_cast<int>(sendBuffers.size()); i++)
197
      delete [] sendBuffers[i];    
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
  double ParallelDomainProblemBase::setElemWeights(AdaptInfo *adaptInfo) 
  {
    double localWeightSum = 0.0;
    int elNum = -1;

    elemWeights.clear();

    TraverseStack stack;
    ElInfo *elInfo = stack.traverseFirst(mesh, -1,
					 Mesh::CALL_EVERY_EL_PREORDER);
    while (elInfo) {
      Element *element = elInfo->getElement();

      // get partition data
      PartitionElementData *partitionData = dynamic_cast<PartitionElementData*>
	(element->getElementData(PARTITION_ED));

      if (partitionData && partitionData->getPartitionStatus() == IN) {
	if (partitionData->getLevel() == 0) {
	  elNum = element->getIndex();
	}
	TEST_EXIT(elNum != -1)("invalid element number\n");
	if (element->isLeaf()) {
	  elemWeights[elNum] += 1.0;
	  localWeightSum += 1.0;
	}
      }

      elInfo = stack.traverseNext(elInfo);
    }

    return localWeightSum;
  }

235

236
237
238
239
240
241
242
243
244
245
246
247
248
249
  void ParallelDomainProblemBase::partitionMesh(AdaptInfo *adaptInfo)
  {
    if (initialPartitionMesh) {
      initialPartitionMesh = false;
      partitioner->fillCoarsePartitionVec(&oldPartitionVec);
      partitioner->partition(&elemWeights, INITIAL);
    } else {
      oldPartitionVec = partitionVec;
      partitioner->partition(&elemWeights, ADAPTIVE_REPART, 100.0 /*0.000001*/);
    }    

    partitioner->fillCoarsePartitionVec(&partitionVec);
  }

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

  void ParallelDomainProblemBase::createInteriorBoundaryInfo()
  {
    TraverseStack stack;
    ElInfo *elInfo = stack.traverseFirst(mesh, -1, Mesh::CALL_LEAF_EL | Mesh::FILL_NEIGH);
    while (elInfo) {
      Element *element = elInfo->getElement();

      // Hidde elements which are not part of ranks partition.
      PartitionElementData *partitionData = 
	dynamic_cast<PartitionElementData*>(element->getElementData(PARTITION_ED));   
      if (partitionData->getPartitionStatus() == IN) {
	for (int i = 0; i < 3; i++) {
	  if (!elInfo->getNeighbour(i))
	    continue;

	  PartitionElementData *neighbourPartitionData =
	    dynamic_cast<PartitionElementData*>(elInfo->getNeighbour(i)->getElementData(PARTITION_ED));
 	  if (neighbourPartitionData->getPartitionStatus() == OUT) {
 	    AtomicBoundary& bound = interiorBoundary.
	      getNewAtomicBoundary(partitionVec[elInfo->getNeighbour(i)->getIndex()]);
 	    bound.rankObject.el = element;
 	    bound.rankObject.subObjAtBoundary = EDGE;
 	    bound.rankObject.ithObjAtBoundary = i;
 	    bound.neighbourObject.el = elInfo->getNeighbour(i);
 	    bound.neighbourObject.subObjAtBoundary = EDGE;
 	    bound.neighbourObject.ithObjAtBoundary = -1;
 	  }
	}
      }

      elInfo = stack.traverseNext(elInfo);
    }
  }


  void ParallelDomainProblemBase::removeMacroElements()
  {
    std::vector<MacroElement*> macrosToRemove;
    for (std::deque<MacroElement*>::iterator it = mesh->firstMacroElement();
	 it != mesh->endOfMacroElements();
	 ++it) {
      PartitionElementData *partitionData = 
	dynamic_cast<PartitionElementData*>
	((*it)->getElement()->getElementData(PARTITION_ED));
      if (partitionData->getPartitionStatus() != IN)
	macrosToRemove.push_back(*it);
    }

    mesh->removeMacroElements(macrosToRemove);
  }


  void ParallelDomainProblemBase::createLocalGlobalNumbering(int& nRankDOFs, 
							     int& nOverallDOFs)
  {
    /// === Get rank information about DOFs. ===

    // Stores to each DOF pointer the set of ranks the DOF is part of.
    std::map<const DegreeOfFreedom*, std::set<int> > partitionDOFs;
    // Set of all DOFs of the rank.
    std::vector<const DegreeOfFreedom*> rankDOFs;
    // Set of all interior boundary DOFs in ranks partition which are owned by 
    // another rank.
    std::map<const DegreeOfFreedom*, int> boundaryDOFs;

    createDOFMemberInfo(partitionDOFs, rankDOFs, boundaryDOFs);

    nRankDOFs = rankDOFs.size();
    nOverallDOFs = partitionDOFs.size();

    // === Get starting position for global rank dof ordering ====

    int rstart = 0;
    MPI_Scan(&nRankDOFs, &rstart, 1, MPI_INT, MPI_SUM, PETSC_COMM_WORLD);
    rstart -= nRankDOFs;

   
    /// === Create information which dof indices must be send and which must be received. ===

    std::map<int, std::map<DegreeOfFreedom, DegreeOfFreedom> > sendNewDofs;
    std::map<int, std::vector<DegreeOfFreedom> > recvNewDofs;

    for (std::map<const DegreeOfFreedom*, int>::iterator it = boundaryDOFs.begin();
	 it != boundaryDOFs.end();
	 ++it) {

      if (it->second == mpiRank) {
	// If the boundary dof is a rank dof, it must be send to other ranks.

	// old global index
	int oldDofIndex = (it->first)[0];
	// search for new dof index in ranks partition for this boundary dof
	int newDofIndex = 0;
	for (int i = 0; i < static_cast<int>(rankDOFs.size()); i++) {
	  if (rankDOFs[i] == it->first) {
	    newDofIndex = rstart + i;
	    break;
	  }
	}

	// Search for all ranks that have this dof too.
	for (std::set<int>::iterator itRanks = partitionDOFs[it->first].begin();
	     itRanks != partitionDOFs[it->first].end();
	     ++itRanks) {
	  if (*itRanks != mpiRank)
	    sendNewDofs[*itRanks][oldDofIndex] = newDofIndex;
	}
      } else {
	// If the boundary dof is not a rank dof, its new dof index, and later
	// also the dof values, must be received from another rank.
	recvNewDofs[it->second].push_back((it->first)[0]);
      }
    }


    /// === Send and receive the dof indices at boundary. ===

    std::vector<int*> sendBuffers(sendNewDofs.size());
    std::vector<int*> recvBuffers(recvNewDofs.size());
    
    int i = 0;
    for (std::map<int, std::map<DegreeOfFreedom, DegreeOfFreedom> >::iterator sendIt = sendNewDofs.begin();
	 sendIt != sendNewDofs.end();
	 ++sendIt, i++) {
      sendBuffers[i] = new int[sendIt->second.size() * 2];
      int c = 0;
      for (std::map<DegreeOfFreedom, DegreeOfFreedom>::iterator dofIt = sendIt->second.begin();
	   dofIt != sendIt->second.end();
	   ++dofIt, c += 2) {
	sendBuffers[i][c] = dofIt->first;
	sendBuffers[i][c + 1] = dofIt->second;

	sendDofs[sendIt->first].push_back(dofIt->second);
      }

      mpiComm.Isend(sendBuffers[i], sendIt->second.size() * 2, MPI_INT, sendIt->first, 0);
    }

    i = 0;
    for (std::map<int, std::vector<DegreeOfFreedom> >::iterator recvIt = recvNewDofs.begin();
	 recvIt != recvNewDofs.end();
	 ++recvIt, i++) {
      recvBuffers[i] = new int[recvIt->second.size() * 2];

      mpiComm.Irecv(recvBuffers[i], recvIt->second.size() * 2, MPI_INT, recvIt->first, 0);
    }


    mpiComm.Barrier();

    
    /// === Delete send buffers. ===

    i = 0;
    for (std::map<int, std::map<DegreeOfFreedom, DegreeOfFreedom> >::iterator sendIt = sendNewDofs.begin();
	 sendIt != sendNewDofs.end();
	 ++sendIt, i++) 
      delete [] sendBuffers[i];


    /// === Change dof indices for rank partition. ===

    for (int i = 0; i < static_cast<int>(rankDOFs.size()); i++) {
      const_cast<DegreeOfFreedom*>(rankDOFs[i])[0] = i; 
      mapLocalGlobalDOFs[i] = rstart + i;
      mapGlobalLocalDOFs[rstart + i] = i;
      isRankDOF[i] = true;
    }


    /// === Change dof indices at boundary from other ranks. ===

    i = 0;
    for (std::map<int, std::vector<DegreeOfFreedom> >::iterator recvIt = recvNewDofs.begin();
	 recvIt != recvNewDofs.end();
	 ++recvIt, i++) {

      for (int j = 0; j < static_cast<int>(recvIt->second.size()); j++) {

	int oldDof = recvBuffers[i][j * 2];
	int newDof = recvBuffers[i][j * 2 + 1];
	int newLocalDof = mapLocalGlobalDOFs.size();

	recvDofs[recvIt->first].push_back(newDof);

	for (std::map<const DegreeOfFreedom*, int>::iterator dofIt = boundaryDOFs.begin();
	     dofIt != boundaryDOFs.end();
	     ++dofIt) {
	  if ((dofIt->first)[0] == oldDof) {
	    const_cast<DegreeOfFreedom*>(dofIt->first)[0] = newLocalDof;
	    mapLocalGlobalDOFs[newLocalDof] = newDof;
	    mapGlobalLocalDOFs[newDof] = newLocalDof;
	    isRankDOF[newLocalDof] = false;
	    break;
	  }
	}
      }

      delete [] recvBuffers[i];
    }


    /// === Create local information from sendDofs and recvDofs

    for (std::map<int, std::vector<DegreeOfFreedom> >::iterator it = sendDofs.begin();
	 it != sendDofs.end();
	 ++it)
      for (std::vector<DegreeOfFreedom>::iterator dofIt = it->second.begin();
	   dofIt != it->second.end();
	   dofIt++)
	*dofIt = mapGlobalLocalDOFs[*dofIt];

    for (std::map<int, std::vector<DegreeOfFreedom> >::iterator it = recvDofs.begin();
	 it != recvDofs.end();
	 ++it)
      for (std::vector<DegreeOfFreedom>::iterator dofIt = it->second.begin();
	   dofIt != it->second.end();
	   dofIt++)
	*dofIt = mapGlobalLocalDOFs[*dofIt];  
  }



474
  void ParallelDomainProblemBase::createDOFMemberInfo(
475
		       std::map<const DegreeOfFreedom*, std::set<int> >& partitionDOFs,
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
		       std::vector<const DegreeOfFreedom*>& rankDOFs,
		       std::map<const DegreeOfFreedom*, int>& boundaryDOFs)
  {
    /// === Determine to each dof the set of partitions the dof belongs to. ===

    TraverseStack stack;
    ElInfo *elInfo = stack.traverseFirst(mesh, -1, Mesh::CALL_LEAF_EL);
    while (elInfo) {
      Element *element = elInfo->getElement();

      // Determine to each dof the partition(s) it corresponds to.
      for (int i = 0; i < 3; i++) 
	partitionDOFs[element->getDOF(i)].insert(partitionVec[element->getIndex()]);
          
      elInfo = stack.traverseNext(elInfo);
    }

    /// === Determine the set of ranks dofs and the dofs ownership at the boundary. ===

    for (std::map<const DegreeOfFreedom*, std::set<int> >::iterator it = partitionDOFs.begin();
	 it != partitionDOFs.end();
	 ++it) {
      for (std::set<int>::iterator itpart1 = it->second.begin();
	   itpart1 != it->second.end();
	   ++itpart1) {
	if (*itpart1 == mpiRank) {
	  if (it->second.size() == 1) {
	    rankDOFs.push_back(it->first);
	  } else {	    
	    // This dof is at the ranks boundary. It is owned by the rank only if
	    // the rank number is the highest of all ranks containing this dof.

	    bool insert = true;
	    int highestRank = mpiRank;
	    for (std::set<int>::iterator itpart2 = it->second.begin();
		 itpart2 != it->second.end();
		 ++itpart2) {
	      if (*itpart2 > mpiRank)
		insert = false;

	      if (*itpart2 > highestRank)
		highestRank = *itpart2;
	    }

	    if (insert)
	      rankDOFs.push_back(it->first);

	    boundaryDOFs[it->first] = highestRank;
	  }

	  break;
	}
      }
    }
  }


533
534
535
  ParallelDomainProblemScal::ParallelDomainProblemScal(const std::string& name,
						       ProblemScal *problem,
						       ProblemInstatScal *problemInstat)
536
537
    : ParallelDomainProblemBase(name, problem, problemInstat, problem->getFESpace()),
      probScal(problem)
538
539
540
  {
  }

541
542
543
544
545
546
547
  void ParallelDomainProblemScal::initParallelization(AdaptInfo *adaptInfo)
  {
    ParallelDomainProblemBase::initParallelization(adaptInfo);

    probScal->getSystemMatrix()->setIsRankDOF(isRankDOF);
  }

548
549
550
551
552
553
554
555
  Flag ParallelDomainProblemScal::oneIteration(AdaptInfo *adaptInfo, Flag toDo)
  {
    //    return iterationIF->oneIteration(adaptInfo, toDo);

    Flag flag = dynamic_cast<StandardProblemIteration*>(iterationIF)->buildAndAdapt(adaptInfo, toDo);

    fillPetscMatrix(probScal->getSystemMatrix(), probScal->getRHS());

556
557
    solvePetscMatrix(probScal->getSolution());

558
559
560
561
562
563
564
565
566
567
568
569
//     if (toDo.isSet(SOLVE))
//       iterationIF->getProblem()->solve(adaptInfo, false);

//     if (toDo.isSet(SOLVE_RHS))
//       iterationIF->getProblem()->solve(adaptInfo, true);

//     if (toDo.isSet(ESTIMATE)) 
//       iterationIF->getProblem()->estimate(adaptInfo);

    return flag;

  }
570
571

}