SubDomainSolver.cc 11.2 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"
14
#include "SystemVector.h"
15
16
17
18
19
20

namespace AMDiS {
  
  using namespace std;


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
  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);
                 
      MSG("COMM TEST: %d %d %d %d %d\n",
	  mpiCommInterior.Get_size(), 
	  interiorMap->getRankDofs(),
	  interiorMap->getOverallDofs(),
	  nGlobalOverallInterior, rStartInterior);
    }
  }


53
54
  void SubDomainSolver::fillPetscMatrix(Matrix<DOFMatrix*> *mat)
  {
55
56
57
58
    FUNCNAME("SubDomainSolver::fillPetscMatrix()");

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

59
60
61
62
    int nRowsRankInterior = interiorMap->getRankDofs();
    int nRowsOverallInterior = interiorMap->getOverallDofs();
    int nRowsRankCoarse = coarseSpaceMap->getRankDofs();
    int nRowsOverallCoarse = coarseSpaceMap->getOverallDofs();
63

64
    bool multilevel = false;
65
    if (mpiCommInterior.Get_size() == 1) {
66
67
      nGlobalOverallInterior = nRowsOverallInterior;

68
69
70
      MatCreateSeqAIJ(mpiCommInterior, nRowsRankInterior, nRowsRankInterior,
		      60, PETSC_NULL, &matIntInt);
    } else {
71
72
      multilevel = true;

73
74
75
76
77
      MatCreateMPIAIJ(mpiCommInterior, 
		      nRowsRankInterior, nRowsRankInterior,
		      nRowsOverallInterior, nRowsOverallInterior,
		      60, PETSC_NULL, 60, PETSC_NULL, &matIntInt);
    }
78
79
80
81
82
83
84
85

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

    MatCreateMPIAIJ(mpiCommCoarseSpace,
		    nRowsRankCoarse, nRowsRankInterior,
86
		    nRowsOverallCoarse, nGlobalOverallInterior,
87
88
89
90
		    60, PETSC_NULL, 60, PETSC_NULL, &matCoarseInt);

    MatCreateMPIAIJ(mpiCommCoarseSpace,
		    nRowsRankInterior, nRowsRankCoarse,
91
		    nGlobalOverallInterior, nRowsOverallCoarse,
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
		    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) {
162
	    int rowIndex = coarseSpaceMap->getMatIndex(i, *cursor);
163
	    for (unsigned int k = 0; k < cols.size(); k++)
164
	      cols[k] = coarseSpaceMap->getMatIndex(j, cols[k]);
165
166
167
168
169

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

	    if (colsOther.size()) {
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
	      if (multilevel == false) {
		for (unsigned int k = 0; k < colsOther.size(); k++)
		  colsOther[k] = interiorMap->getMatIndex(j, colsOther[k]);
	      } else {
		for (unsigned int k = 0; k < colsOther.size(); k++) {
		  colsOther[k] = 
		    interiorMap->getMatIndex(j, colsOther[k]) + rStartInterior;
		  if (colsOther[k] == 50723) {
		    MSG("RRTEST A: %d %d %d\n",
			nRowsOverallInterior,
			nGlobalOverallInterior,
			rStartInterior);
			
		  }
		}
	      }
186
187
188
189
190
 	      
	      MatSetValues(matCoarseInt, 1, &rowIndex, colsOther.size(),
 			   &(colsOther[0]), &(valuesOther[0]), ADD_VALUES);
	    }
	  } else {
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
	    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]);
	      else  {
		cols[k] = interiorMap->getMatIndex(j, cols[k]);
		  if (colsOther[k] == 50723) {
		    MSG("RRTEST B: %d %d %d\n",
			nRowsOverallInterior,
			nGlobalOverallInterior,
			rStartInterior);
			
		  }
207

208
209
	      }
	    }
210
211
212
213
214
	    
  	    MatSetValues(matIntInt, 1, &localRowIndex, cols.size(),
  			 &(cols[0]), &(values[0]), ADD_VALUES);

	    if (colsOther.size()) {
215
	      int globalRowIndex = interiorMap->getMatIndex(i, *cursor);
216
217
218
219

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

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

  	      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);
255
256
257
258
    if (multilevel == false)
      PCFactorSetMatSolverPackage(pcInterior, MATSOLVERUMFPACK);
    else
      PCFactorSetMatSolverPackage(pcInterior, MATSOLVERMUMPS);
259
    KSPSetFromOptions(kspInterior);
260
261
262
263
264
  }


  void SubDomainSolver::fillPetscRhs(SystemVector *vec)
  {
265
266
267
    FUNCNAME("SubDomainSolver::fillPetscRhs()");

    VecCreateMPI(mpiCommCoarseSpace, 
268
269
		 coarseSpaceMap->getRankDofs(), 
		 coarseSpaceMap->getOverallDofs(),
270
271
272
		 &rhsCoarseSpace);

    VecCreateMPI(mpiCommCoarseSpace, 
273
		 interiorMap->getRankDofs(), 
274
		 nGlobalOverallInterior,
275
276
277
278
279
280
281
282
		 &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)) {
283
	  index = coarseSpaceMap->getMatIndex(i, index);
284
285
	  VecSetValue(rhsCoarseSpace, index, *dofIt, ADD_VALUES);
	} else {
286
	  index = interiorMap->getMatIndex(i, index) + rStartInterior;
287
288
289
290
291
292
293
294
295
296
	  VecSetValue(rhsInterior, index, *dofIt, INSERT_VALUES);
	}      
      }
    }

    VecAssemblyBegin(rhsCoarseSpace);
    VecAssemblyEnd(rhsCoarseSpace);

    VecAssemblyBegin(rhsInterior);
    VecAssemblyEnd(rhsInterior);
297
298
299
300
301
  }


  void SubDomainSolver::solvePetscMatrix(SystemVector &vec, AdaptInfo *adaptInfo)
  {
302
    FUNCNAME("SubDomainSolver::solvePetscMatrix()");
303
304
305
306
307
308
309
  }


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

310
311
    VecDestroy(&rhsCoarseSpace);
    VecDestroy(&rhsInterior);
312
313
314
315
316
  }


  void SubDomainSolver::destroyMatrixData()
  {
317
318
319
320
321
322
323
324
    FUNCNAME("SubDomainSolver::destroyMatrixData()");

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

    KSPDestroy(&kspInterior);
325
326
327
328
329
  }


  void SubDomainSolver::solve(Vec &rhs, Vec &sol)
  {
330
331
332
333
    KSPSolve(kspInterior, rhs, sol);
  }


334
335
336
337
  void SubDomainSolver::solveGlobal(Vec &rhs, Vec &sol)
  {
    FUNCNAME("SubDomainSolver::solveGlobal()");

338
339
340
    int ml = 0;
    Parameters::get("parallel->multi level test", ml);

341
    Vec tmp;
342
343
344
345
346
347
348
    if (ml == 0)
      VecCreateSeq(PETSC_COMM_SELF, interiorMap->getRankDofs(), &tmp);
    else
      VecCreateMPI(mpiCommInterior,
		   interiorMap->getRankDofs(),
		   interiorMap->getOverallDofs(),
		   &tmp);
349
350
351
352
353

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

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

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

    KSPSolve(kspInterior, tmp, tmp);

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

365
    for (int i = 0; i < interiorMap->getRankDofs(); i++) 
366
367
368
369
370
371
372
373
374
      rhsValues[i] = tmpValues[i];

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

    VecDestroy(&tmp);
  }


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
  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;
404
  }
405

406
}