SubDomainSolver.cc 13.7 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
//
// 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/SubDomainSolver.h"
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
14
#include "parallel/ParallelDebug.h"
15
#include "SystemVector.h"
16
17
18
19
20
21

namespace AMDiS {
  
  using namespace std;


22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
  void SubDomainSolver::setDofMapping(ParallelDofMapping *coarseDofs,
				      ParallelDofMapping *interiorDofs)
  {
    coarseSpaceMap = coarseDofs;
    interiorMap = interiorDofs;

    if (mpiCommInterior.Get_size() == 1) {
      rStartInterior = 0;
      nGlobalOverallInterior = interiorMap->getOverallDofs();
    } else {
      int groupRowsInterior = 0;
      if (mpiCommInterior.Get_rank() == 0)
	groupRowsInterior = interiorMap->getOverallDofs();

      mpi::getDofNumbering(mpiCommCoarseSpace, groupRowsInterior,
			   rStartInterior, nGlobalOverallInterior);

      int tmp = 0;
      if (mpiCommInterior.Get_rank() == 0)
	tmp = rStartInterior;

      mpiCommInterior.Allreduce(&tmp, &rStartInterior, 1, MPI_INT, MPI_SUM);
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
44
45

      /*                 
46
47
48
49
50
      MSG("COMM TEST: %d %d %d %d %d\n",
	  mpiCommInterior.Get_size(), 
	  interiorMap->getRankDofs(),
	  interiorMap->getOverallDofs(),
	  nGlobalOverallInterior, rStartInterior);
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
51
      */
52
53
54
55
    }
  }


56
57
  void SubDomainSolver::fillPetscMatrix(Matrix<DOFMatrix*> *mat)
  {
58
59
60
61
    FUNCNAME("SubDomainSolver::fillPetscMatrix()");

    vector<const FiniteElemSpace*> feSpaces = getFeSpaces(mat);

62
63
64
65
    int nRowsRankInterior = interiorMap->getRankDofs();
    int nRowsOverallInterior = interiorMap->getOverallDofs();
    int nRowsRankCoarse = coarseSpaceMap->getRankDofs();
    int nRowsOverallCoarse = coarseSpaceMap->getOverallDofs();
66

67
    bool multilevel = false;
68
    if (mpiCommInterior.Get_size() == 1) {
69
70
      nGlobalOverallInterior = nRowsOverallInterior;

71
72
73
      MatCreateSeqAIJ(mpiCommInterior, nRowsRankInterior, nRowsRankInterior,
		      60, PETSC_NULL, &matIntInt);
    } else {
74
75
      multilevel = true;

76
77
78
79
80
      MatCreateMPIAIJ(mpiCommInterior, 
		      nRowsRankInterior, nRowsRankInterior,
		      nRowsOverallInterior, nRowsOverallInterior,
		      60, PETSC_NULL, 60, PETSC_NULL, &matIntInt);
    }
81
82
83
84
85
86
87
88

    MatCreateMPIAIJ(mpiCommCoarseSpace,
		    nRowsRankCoarse, nRowsRankCoarse,
		    nRowsOverallCoarse, nRowsOverallCoarse,
		    60, PETSC_NULL, 60, PETSC_NULL, &matCoarseCoarse);

    MatCreateMPIAIJ(mpiCommCoarseSpace,
		    nRowsRankCoarse, nRowsRankInterior,
89
		    nRowsOverallCoarse, nGlobalOverallInterior,
90
91
92
93
		    60, PETSC_NULL, 60, PETSC_NULL, &matCoarseInt);

    MatCreateMPIAIJ(mpiCommCoarseSpace,
		    nRowsRankInterior, nRowsRankCoarse,
94
		    nGlobalOverallInterior, nRowsOverallCoarse,
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
		    60, PETSC_NULL, 60, PETSC_NULL, &matIntCoarse);

    // === Prepare traverse of sequentially created matrices. ===

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

    // === Traverse all sequentially created matrices and add the values to ===
    // === the global PETSc matrices.                                       ===

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

	traits::col<Matrix>::type col((*mat)[i][j]->getBaseMatrix());
	traits::const_value<Matrix>::type value((*mat)[i][j]->getBaseMatrix());
	
	// Traverse all rows.
	for (cursor_type cursor = begin<row>((*mat)[i][j]->getBaseMatrix()), 
	       cend = end<row>((*mat)[i][j]->getBaseMatrix()); cursor != cend; ++cursor) {

	  bool rowPrimal = isCoarseSpace(feSpaces[i], *cursor);
  
	  cols.clear();
	  colsOther.clear();
	  values.clear();	  
	  valuesOther.clear();

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

	    bool colPrimal = isCoarseSpace(feSpaces[j], col(*icursor));

	    if (colPrimal) {
	      if (rowPrimal) {
		cols.push_back(col(*icursor));
		values.push_back(value(*icursor));
	      } else {
		colsOther.push_back(col(*icursor));
		valuesOther.push_back(value(*icursor));
	      }
	    } else {
	      if (rowPrimal) {
		colsOther.push_back(col(*icursor));
		valuesOther.push_back(value(*icursor));
	      } else {
		cols.push_back(col(*icursor));
		values.push_back(value(*icursor));
	      }
	    }
	  }  // for each nnz in row


	  // === Set matrix values. ===

	  if (rowPrimal) {
165
	    int rowIndex = coarseSpaceMap->getMatIndex(i, *cursor);
166
	    for (unsigned int k = 0; k < cols.size(); k++)
167
	      cols[k] = coarseSpaceMap->getMatIndex(j, cols[k]);
168
169
170
171
172

	    MatSetValues(matCoarseCoarse, 1, &rowIndex, cols.size(),
			 &(cols[0]), &(values[0]), ADD_VALUES);

	    if (colsOther.size()) {
173
174
175
176
	      if (multilevel == false) {
		for (unsigned int k = 0; k < colsOther.size(); k++)
		  colsOther[k] = interiorMap->getMatIndex(j, colsOther[k]);
	      } else {
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
177
		for (unsigned int k = 0; k < colsOther.size(); k++)
178
179
180
		  colsOther[k] = 
		    interiorMap->getMatIndex(j, colsOther[k]) + rStartInterior;
	      }
181
182
183
184
185
 	      
	      MatSetValues(matCoarseInt, 1, &rowIndex, colsOther.size(),
 			   &(colsOther[0]), &(valuesOther[0]), ADD_VALUES);
	    }
	  } else {
186
187
188
189
190
191
192
	    int localRowIndex = 
	      (multilevel == false ? interiorMap->getLocalMatIndex(i, *cursor) :
	       interiorMap->getMatIndex(i, *cursor));

	    for (unsigned int k = 0; k < cols.size(); k++) {
	      if (multilevel == false)
		cols[k] = interiorMap->getLocalMatIndex(j, cols[k]);
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
193
	      else
194
195
		cols[k] = interiorMap->getMatIndex(j, cols[k]);
	    }
196
197
198
199
200
	    
  	    MatSetValues(matIntInt, 1, &localRowIndex, cols.size(),
  			 &(cols[0]), &(values[0]), ADD_VALUES);

	    if (colsOther.size()) {
201
	      int globalRowIndex = interiorMap->getMatIndex(i, *cursor);
202
203
204
205

	      if (multilevel == false)
		globalRowIndex += rStartInterior;

206
	      for (unsigned int k = 0; k < colsOther.size(); k++)
207
		colsOther[k] = coarseSpaceMap->getMatIndex(j, colsOther[k]);
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

  	      MatSetValues(matIntCoarse, 1, &globalRowIndex, colsOther.size(),
  			   &(colsOther[0]), &(valuesOther[0]), ADD_VALUES);
	    }
	  }
	} 
      }
    }

    // === Start global assembly procedure. ===

    MatAssemblyBegin(matIntInt, MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(matIntInt, MAT_FINAL_ASSEMBLY);

    MatAssemblyBegin(matCoarseCoarse, MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(matCoarseCoarse, MAT_FINAL_ASSEMBLY);

    MatAssemblyBegin(matIntCoarse, MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(matIntCoarse, MAT_FINAL_ASSEMBLY);

    MatAssemblyBegin(matCoarseInt, MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(matCoarseInt, MAT_FINAL_ASSEMBLY);


    // === Create solver for the non primal (thus local) variables. ===

    KSPCreate(mpiCommInterior, &kspInterior);
    KSPSetOperators(kspInterior, matIntInt, matIntInt, SAME_NONZERO_PATTERN);
    KSPSetOptionsPrefix(kspInterior, "interior_");
    KSPSetType(kspInterior, KSPPREONLY);
    PC pcInterior;
    KSPGetPC(kspInterior, &pcInterior);
    PCSetType(pcInterior, PCLU);
241
242
243
244
    if (multilevel == false)
      PCFactorSetMatSolverPackage(pcInterior, MATSOLVERUMFPACK);
    else
      PCFactorSetMatSolverPackage(pcInterior, MATSOLVERMUMPS);
245
    KSPSetFromOptions(kspInterior);
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
246
247
248
249
250
251
252
253


#if 0
    PetscViewer matview;
    PetscViewerBinaryOpen(mpiCommCoarseSpace, "mat_primal.dat", 
			  FILE_MODE_WRITE, &matview);
    MatView(matCoarseCoarse, matview);
    PetscViewerDestroy(&matview);
Thomas Witkowski's avatar
Thomas Witkowski committed
254
255
256
257

    ParallelDebug::writeDebugFile(feSpaces[0],
				  *coarseSpaceMap,
				  "coarsespace", "dat");
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
258
259
260
#endif

#if 0
Thomas Witkowski's avatar
Thomas Witkowski committed
261
262
263
264
    if (MPI::COMM_WORLD.Get_rank() == 4 ||
	MPI::COMM_WORLD.Get_rank() == 5 ||
	MPI::COMM_WORLD.Get_rank() == 6 ||
	MPI::COMM_WORLD.Get_rank() == 7) {
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
265
266
267
268
269
270
271
272
273
274
275
276
      PetscViewerBinaryOpen(mpiCommInterior, "mat_interior.dat", 
			    FILE_MODE_WRITE, &matview);
      MatView(matIntInt, matview);
      PetscViewerDestroy(&matview);

      ParallelDebug::writeDebugFile(feSpaces[0],
				    *interiorMap,
				    "interiorspace", "dat");
    }
#endif

#if 0
Thomas Witkowski's avatar
Thomas Witkowski committed
277
    if (MPI::COMM_WORLD.Get_rank() == 1) {
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
278
279
280
281
282
283
284
285
286
287
288
289
      PetscViewerBinaryOpen(mpiCommInterior, "mat_interior.dat", 
			    FILE_MODE_WRITE, &matview);
      MatView(matIntInt, matview);
      PetscViewerDestroy(&matview);

      ParallelDebug::writeDebugFile(feSpaces[0],
				    *interiorMap,
				    "interiorspace", "dat");
    }
#endif


290
291
292
293
294
  }


  void SubDomainSolver::fillPetscRhs(SystemVector *vec)
  {
295
296
297
    FUNCNAME("SubDomainSolver::fillPetscRhs()");

    VecCreateMPI(mpiCommCoarseSpace, 
298
299
		 coarseSpaceMap->getRankDofs(), 
		 coarseSpaceMap->getOverallDofs(),
300
301
302
		 &rhsCoarseSpace);

    VecCreateMPI(mpiCommCoarseSpace, 
303
		 interiorMap->getRankDofs(), 
304
		 nGlobalOverallInterior,
305
306
307
308
309
310
311
312
		 &rhsInterior);

    for (int i = 0; i < vec->getSize(); i++) {
      const FiniteElemSpace *feSpace = vec->getDOFVector(i)->getFeSpace();
      DOFVector<double>::Iterator dofIt(vec->getDOFVector(i), USED_DOFS);
      for (dofIt.reset(); !dofIt.end(); ++dofIt) {
	int index = dofIt.getDOFIndex();
	if (isCoarseSpace(feSpace, index)) {
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
313
314
315
316
	  if ((*coarseSpaceMap)[feSpace].isRankDof(index)) {
	    index = coarseSpaceMap->getMatIndex(i, index);
	    VecSetValue(rhsCoarseSpace, index, *dofIt, INSERT_VALUES);
	  }
317
	} else {
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
318
319
320
321
322
323
	  if ((*interiorMap)[feSpace].isRankDof(index)) {
	    int d = index;
	    index = interiorMap->getMatIndex(i, index) + rStartInterior;
	    int t = interiorMap->getLocalMatIndex(i, d);
	    VecSetValue(rhsInterior, index, *dofIt, INSERT_VALUES);
	  }
324
325
326
327
328
329
330
331
332
	}      
      }
    }

    VecAssemblyBegin(rhsCoarseSpace);
    VecAssemblyEnd(rhsCoarseSpace);

    VecAssemblyBegin(rhsInterior);
    VecAssemblyEnd(rhsInterior);
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
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

#if 0
    PetscViewer matview;
    PetscViewerBinaryOpen(mpiCommCoarseSpace, "vec_interior.dat", 
			  FILE_MODE_WRITE, &matview);
    VecView(rhsInterior, matview);
    PetscViewerDestroy(&matview);

#if 1
    {
      DofMap &m = (*interiorMap)[vec->getDOFVector(0)->getFeSpace()].getMap();

      int groupRowsInterior = 0;
      if (mpiCommInterior.Get_rank() == 0)
	groupRowsInterior = 
	  (*interiorMap)[vec->getDOFVector(0)->getFeSpace()].nOverallDofs;

      int rStart, nGlobal;
      mpi::getDofNumbering(mpiCommCoarseSpace, groupRowsInterior, rStart, nGlobal);

      int tmp = 0;
      if (mpiCommInterior.Get_rank() == 0)
	tmp = rStart;

      int a = 0;
      mpiCommInterior.Allreduce(&tmp, &a, 1, MPI_INT, MPI_SUM);

      for (DofMap::iterator it = m.begin(); it != m.end(); it++)
	it->second.global += a;
    }
#else
    {
      int groupRowsInterior = 0;
      groupRowsInterior = 
	(*interiorMap)[vec->getDOFVector(0)->getFeSpace()].nOverallDofs;
      
      int rStart, nGlobal;
      mpi::getDofNumbering(mpiCommCoarseSpace, groupRowsInterior, rStart, nGlobal);

      DofMap &m = (*interiorMap)[vec->getDOFVector(0)->getFeSpace()].getMap();
      for (DofMap::iterator it = m.begin(); it != m.end(); it++)
	it->second.global = it->second.local + 
	  (*interiorMap)[vec->getDOFVector(0)->getFeSpace()].rStartDofs;
    }
#endif

    ParallelDebug::writeDebugFile(vec->getDOFVector(0)->getFeSpace(),
 				  *interiorMap,
 				  "interior", "dat");

    exit(0);
#endif
385
386
387
388
389
  }


  void SubDomainSolver::solvePetscMatrix(SystemVector &vec, AdaptInfo *adaptInfo)
  {
390
    FUNCNAME("SubDomainSolver::solvePetscMatrix()");
391
392
393
394
395
396
397
  }


  void SubDomainSolver::destroyVectorData()
  {
    FUNCNAME("SubDomainSolver::destroyVectorData()");

398
399
    VecDestroy(&rhsCoarseSpace);
    VecDestroy(&rhsInterior);
400
401
402
403
404
  }


  void SubDomainSolver::destroyMatrixData()
  {
405
406
407
408
409
410
411
412
    FUNCNAME("SubDomainSolver::destroyMatrixData()");

    MatDestroy(&matIntInt);
    MatDestroy(&matCoarseCoarse);
    MatDestroy(&matCoarseInt);
    MatDestroy(&matIntCoarse);

    KSPDestroy(&kspInterior);
413
414
415
416
417
  }


  void SubDomainSolver::solve(Vec &rhs, Vec &sol)
  {
418
419
420
421
    KSPSolve(kspInterior, rhs, sol);
  }


422
423
424
425
  void SubDomainSolver::solveGlobal(Vec &rhs, Vec &sol)
  {
    FUNCNAME("SubDomainSolver::solveGlobal()");

426
427
428
    int ml = 0;
    Parameters::get("parallel->multi level test", ml);

429
    Vec tmp;
430
431
432
433
434
435
436
    if (ml == 0)
      VecCreateSeq(PETSC_COMM_SELF, interiorMap->getRankDofs(), &tmp);
    else
      VecCreateMPI(mpiCommInterior,
		   interiorMap->getRankDofs(),
		   interiorMap->getOverallDofs(),
		   &tmp);
437
438
439
440
441

    PetscScalar *tmpValues, *rhsValues;
    VecGetArray(tmp, &tmpValues);
    VecGetArray(rhs, &rhsValues);

442
    for (int i = 0; i < interiorMap->getRankDofs(); i++)
443
444
445
446
447
448
449
450
451
452
      tmpValues[i] = rhsValues[i];

    VecRestoreArray(rhs, &rhsValues);
    VecRestoreArray(tmp, &tmpValues);

    KSPSolve(kspInterior, tmp, tmp);

    VecGetArray(tmp, &tmpValues);
    VecGetArray(sol, &rhsValues);

453
    for (int i = 0; i < interiorMap->getRankDofs(); i++) 
454
455
456
457
458
459
460
461
462
      rhsValues[i] = tmpValues[i];

    VecRestoreArray(sol, &rhsValues);
    VecRestoreArray(tmp, &tmpValues);

    VecDestroy(&tmp);
  }


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
  vector<const FiniteElemSpace*> SubDomainSolver::getFeSpaces(Matrix<DOFMatrix*> *mat)
  {
    FUNCNAME("SubDomainSolver::getFeSpaces()");

    int nComponents = mat->getNumRows();
    vector<const FiniteElemSpace*> result(nComponents);

    for (int i = 0; i < nComponents; i++) 
      for (int j = 0; j < nComponents; j++)
	if ((*mat)[i][j]) {
	  result[i] = (*mat)[i][j]->getRowFeSpace();
	  break;
	}

    return result;
  }


  vector<const FiniteElemSpace*> SubDomainSolver::getFeSpaces(SystemVector *vec)
  {
    FUNCNAME("SubDomainSolver::getFeSpaces()");

    int nComponents = vec->getSize();
    vector<const FiniteElemSpace*> result(nComponents);

    for (int i = 0; i < nComponents; i++) 
      result[i] = vec->getFeSpace(i);

    return result;
492
  }
493

494
}