SubDomainSolver.cc 10.9 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
#include "AMDiS.h"
17
18
19
20
21
22

namespace AMDiS {
  
  using namespace std;


23
24
  void SubDomainSolver::setDofMapping(ParallelDofMapping *interiorDofs,
				      ParallelDofMapping *coarseDofs)
25
26
  {
    interiorMap = interiorDofs;
27
    coarseSpaceMap = coarseDofs;
28

29
    if (mpiCommLocal.Get_size() == 1) {
30
31
32
33
      rStartInterior = 0;
      nGlobalOverallInterior = interiorMap->getOverallDofs();
    } else {
      int groupRowsInterior = 0;
34
      if (mpiCommLocal.Get_rank() == 0)
35
36
	groupRowsInterior = interiorMap->getOverallDofs();

37
      mpi::getDofNumbering(mpiCommGlobal, groupRowsInterior,
38
39
40
			   rStartInterior, nGlobalOverallInterior);

      int tmp = 0;
41
      if (mpiCommLocal.Get_rank() == 0)
42
43
	tmp = rStartInterior;

44
      mpiCommLocal.Allreduce(&tmp, &rStartInterior, 1, MPI_INT, MPI_SUM);
45
46
47
48
    }
  }


49
50
  void SubDomainSolver::fillPetscMatrix(Matrix<DOFMatrix*> *mat)
  {
51
52
53
54
    FUNCNAME("SubDomainSolver::fillPetscMatrix()");

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

55
56
    int nRowsRankInterior = interiorMap->getRankDofs();
    int nRowsOverallInterior = interiorMap->getOverallDofs();
57

58
    if (subdomainLevel == 0) {
59
60
      nGlobalOverallInterior = nRowsOverallInterior;

61
      MatCreateSeqAIJ(mpiCommLocal, nRowsRankInterior, nRowsRankInterior,
62
63
		      60, PETSC_NULL, &matIntInt);
    } else {
64
      MatCreateMPIAIJ(mpiCommLocal, 
65
66
67
68
		      nRowsRankInterior, nRowsRankInterior,
		      nRowsOverallInterior, nRowsOverallInterior,
		      60, PETSC_NULL, 60, PETSC_NULL, &matIntInt);
    }
69

70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
    if (coarseSpaceMap) {
      int nRowsRankCoarse = coarseSpaceMap->getRankDofs();
      int nRowsOverallCoarse = coarseSpaceMap->getOverallDofs();

      MatCreateMPIAIJ(mpiCommGlobal,
		      nRowsRankCoarse, nRowsRankCoarse,
		      nRowsOverallCoarse, nRowsOverallCoarse,
		      60, PETSC_NULL, 60, PETSC_NULL, &matCoarseCoarse);
      
      MatCreateMPIAIJ(mpiCommGlobal,
		      nRowsRankCoarse, nRowsRankInterior,
		      nRowsOverallCoarse, nGlobalOverallInterior,
		      60, PETSC_NULL, 60, PETSC_NULL, &matCoarseInt);
      
      MatCreateMPIAIJ(mpiCommGlobal,
		      nRowsRankInterior, nRowsRankCoarse,
		      nGlobalOverallInterior, nRowsOverallCoarse,
		      60, PETSC_NULL, 60, PETSC_NULL, &matIntCoarse);
    }
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

    // === 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) {
158
	    int rowIndex = coarseSpaceMap->getMatIndex(i, *cursor);
159
	    for (unsigned int k = 0; k < cols.size(); k++)
160
	      cols[k] = coarseSpaceMap->getMatIndex(j, cols[k]);
161
162
163
164
165

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

	    if (colsOther.size()) {
166
	      if (subdomainLevel == 0) {
167
168
169
		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
170
		for (unsigned int k = 0; k < colsOther.size(); k++)
171
172
173
		  colsOther[k] = 
		    interiorMap->getMatIndex(j, colsOther[k]) + rStartInterior;
	      }
174
175
176
177
178
 	      
	      MatSetValues(matCoarseInt, 1, &rowIndex, colsOther.size(),
 			   &(colsOther[0]), &(valuesOther[0]), ADD_VALUES);
	    }
	  } else {
179
	    int localRowIndex = 
180
	      (subdomainLevel == 0 ? interiorMap->getLocalMatIndex(i, *cursor) :
181
182
183
	       interiorMap->getMatIndex(i, *cursor));

	    for (unsigned int k = 0; k < cols.size(); k++) {
184
	      if (subdomainLevel == 0)
185
		cols[k] = interiorMap->getLocalMatIndex(j, cols[k]);
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
186
	      else
187
188
		cols[k] = interiorMap->getMatIndex(j, cols[k]);
	    }
189
190
191
192
193
	    
  	    MatSetValues(matIntInt, 1, &localRowIndex, cols.size(),
  			 &(cols[0]), &(values[0]), ADD_VALUES);

	    if (colsOther.size()) {
194
	      int globalRowIndex = interiorMap->getMatIndex(i, *cursor);
195

196
	      if (subdomainLevel != 0)
197
198
		globalRowIndex += rStartInterior;

199
	      for (unsigned int k = 0; k < colsOther.size(); k++)
200
		colsOther[k] = coarseSpaceMap->getMatIndex(j, colsOther[k]);
201
202
203
204
205
206
207
208
209
210
211
212
213
214

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

215
216
217
218
219
220
221
222
223
224
    if (coarseSpaceMap) {
      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);
    }
225
226
227
228


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

229
    KSPCreate(mpiCommLocal, &kspInterior);
230
231
232
233
234
235
    KSPSetOperators(kspInterior, matIntInt, matIntInt, SAME_NONZERO_PATTERN);
    KSPSetOptionsPrefix(kspInterior, "interior_");
    KSPSetType(kspInterior, KSPPREONLY);
    PC pcInterior;
    KSPGetPC(kspInterior, &pcInterior);
    PCSetType(pcInterior, PCLU);
236
    if (subdomainLevel == 0)
237
238
239
      PCFactorSetMatSolverPackage(pcInterior, MATSOLVERUMFPACK);
    else
      PCFactorSetMatSolverPackage(pcInterior, MATSOLVERMUMPS);
240
    KSPSetFromOptions(kspInterior);
241
242
243
244
245
  }


  void SubDomainSolver::fillPetscRhs(SystemVector *vec)
  {
246
247
    FUNCNAME("SubDomainSolver::fillPetscRhs()");

248
    VecCreateMPI(mpiCommGlobal, 
249
		 interiorMap->getRankDofs(), 
250
		 nGlobalOverallInterior,
251
252
		 &rhsInterior);

253
254
255
256
257
258
259
    if (coarseSpaceMap) 
      VecCreateMPI(mpiCommGlobal, 
		   coarseSpaceMap->getRankDofs(), 
		   coarseSpaceMap->getOverallDofs(),
		   &rhsCoarseSpace);


260
261
262
263
264
    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();
Thomas Witkowski's avatar
Thomas Witkowski committed
265
266
267
	if (isCoarseSpace(feSpace, index)) {	  
	  index = coarseSpaceMap->getMatIndex(i, index);
	  VecSetValue(rhsCoarseSpace, index, *dofIt, ADD_VALUES);
268
	} else {
Thomas Witkowski's avatar
Thomas Witkowski committed
269
270
	  index = interiorMap->getMatIndex(i, index) + rStartInterior;
	  VecSetValue(rhsInterior, index, *dofIt, ADD_VALUES);
271
272
273
274
275
276
	}      
      }
    }

    VecAssemblyBegin(rhsInterior);
    VecAssemblyEnd(rhsInterior);
277
278
279
280
281

    if (coarseSpaceMap) {
      VecAssemblyBegin(rhsCoarseSpace);
      VecAssemblyEnd(rhsCoarseSpace);
    }
282
283
284
285
286
  }


  void SubDomainSolver::solvePetscMatrix(SystemVector &vec, AdaptInfo *adaptInfo)
  {
287
    FUNCNAME("SubDomainSolver::solvePetscMatrix()");
288
289
290
291
292
293
294
  }


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

295
296
    VecDestroy(&rhsCoarseSpace);
    VecDestroy(&rhsInterior);
297
298
299
300
301
  }


  void SubDomainSolver::destroyMatrixData()
  {
302
303
304
305
306
307
308
309
    FUNCNAME("SubDomainSolver::destroyMatrixData()");

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

    KSPDestroy(&kspInterior);
310
311
312
313
314
  }


  void SubDomainSolver::solve(Vec &rhs, Vec &sol)
  {
315
316
317
318
319
320
321
    FUNCNAME("SubDomainSolver::solve()");

    PetscErrorCode solverError = KSPSolve(kspInterior, rhs, sol);
    if (solverError != 0) {
      AMDiS::finalize();
      exit(-1);
    }
322
323
324
  }


325
326
327
328
329
  void SubDomainSolver::solveGlobal(Vec &rhs, Vec &sol)
  {
    FUNCNAME("SubDomainSolver::solveGlobal()");

    Vec tmp;
330
331
    if (mpiCommLocal.Get_size() == 1)
      VecCreateSeq(mpiCommLocal, interiorMap->getRankDofs(), &tmp);
332
    else
333
      VecCreateMPI(mpiCommLocal,
334
335
336
		   interiorMap->getRankDofs(),
		   interiorMap->getOverallDofs(),
		   &tmp);
337
338
339
340
341

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

342
    for (int i = 0; i < interiorMap->getRankDofs(); i++)
343
344
345
346
347
348
349
350
351
352
      tmpValues[i] = rhsValues[i];

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

    KSPSolve(kspInterior, tmp, tmp);

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

353
    for (int i = 0; i < interiorMap->getRankDofs(); i++) 
354
355
356
357
358
359
360
361
362
      rhsValues[i] = tmpValues[i];

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

    VecDestroy(&tmp);
  }


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
  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;
392
  }
393

394
}