SubDomainSolver.cc 9.19 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
21
22

namespace AMDiS {
  
  using namespace std;


  void SubDomainSolver::fillPetscMatrix(Matrix<DOFMatrix*> *mat)
  {
23
24
25
26
    FUNCNAME("SubDomainSolver::fillPetscMatrix()");

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

27
28
29
30
31
32
33
34
35
36
37
38
39
40
    int nRowsRankInterior = interiorMap->getRankDofs(level);
    int nRowsOverallInterior = interiorMap->getOverallDofs(level);
    int nRowsRankCoarse = coarseSpaceMap->getRankDofs(level);
    int nRowsOverallCoarse = coarseSpaceMap->getOverallDofs(level);

    if (mpiCommInterior.Get_size() == 1) {
      MatCreateSeqAIJ(mpiCommInterior, nRowsRankInterior, nRowsRankInterior,
		      60, PETSC_NULL, &matIntInt);
    } else {
      MatCreateMPIAIJ(mpiCommInterior, 
		      nRowsRankInterior, nRowsRankInterior,
		      nRowsOverallInterior, nRowsOverallInterior,
		      60, PETSC_NULL, 60, PETSC_NULL, &matIntInt);
    }
41
42
43
44
45
46
47
48
49
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
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
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186

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

    MatCreateMPIAIJ(mpiCommCoarseSpace,
		    nRowsRankCoarse, nRowsRankInterior,
		    nRowsOverallCoarse, nRowsOverallInterior,
		    60, PETSC_NULL, 60, PETSC_NULL, &matCoarseInt);

    MatCreateMPIAIJ(mpiCommCoarseSpace,
		    nRowsRankInterior, nRowsRankCoarse,
		    nRowsOverallInterior, nRowsOverallCoarse,
		    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) {
	    int rowIndex = coarseSpaceMap->getMatIndex(i, *cursor);
	    for (unsigned int k = 0; k < cols.size(); k++)
	      cols[k] = coarseSpaceMap->getMatIndex(j, cols[k]);

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

	    if (colsOther.size()) {
	      for (unsigned int k = 0; k < colsOther.size(); k++)
		colsOther[k] = interiorMap->getMatIndex(j, colsOther[k]);
 	      
	      MatSetValues(matCoarseInt, 1, &rowIndex, colsOther.size(),
 			   &(colsOther[0]), &(valuesOther[0]), ADD_VALUES);
	    }
	  } else {
	    int localRowIndex = interiorMap->getLocalMatIndex(i, *cursor);

	    for (unsigned int k = 0; k < cols.size(); k++)
	      cols[k] = interiorMap->getLocalMatIndex(j, cols[k]);
	    
  	    MatSetValues(matIntInt, 1, &localRowIndex, cols.size(),
  			 &(cols[0]), &(values[0]), ADD_VALUES);

	    if (colsOther.size()) {
	      int globalRowIndex = interiorMap->getMatIndex(i, *cursor);
	      for (unsigned int k = 0; k < colsOther.size(); k++)
		colsOther[k] = coarseSpaceMap->getMatIndex(j, colsOther[k]);

  	      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);
    PCFactorSetMatSolverPackage(pcInterior, MATSOLVERUMFPACK);
    KSPSetFromOptions(kspInterior);
187
188
189
190
191
  }


  void SubDomainSolver::fillPetscRhs(SystemVector *vec)
  {
192
193
194
    FUNCNAME("SubDomainSolver::fillPetscRhs()");

    VecCreateMPI(mpiCommCoarseSpace, 
195
196
		 coarseSpaceMap->getRankDofs(level), 
		 coarseSpaceMap->getOverallDofs(level),
197
198
199
		 &rhsCoarseSpace);

    VecCreateMPI(mpiCommCoarseSpace, 
200
201
		 interiorMap->getRankDofs(level), 
		 interiorMap->getOverallDofs(level),
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
		 &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)) {
	  index = coarseSpaceMap->getMatIndex(i, index);
	  VecSetValue(rhsCoarseSpace, index, *dofIt, ADD_VALUES);
	} else {
	  index = interiorMap->getMatIndex(i, index);
	  VecSetValue(rhsInterior, index, *dofIt, INSERT_VALUES);
	}      
      }
    }

    VecAssemblyBegin(rhsCoarseSpace);
    VecAssemblyEnd(rhsCoarseSpace);

    VecAssemblyBegin(rhsInterior);
    VecAssemblyEnd(rhsInterior);
224
225
226
227
228
  }


  void SubDomainSolver::solvePetscMatrix(SystemVector &vec, AdaptInfo *adaptInfo)
  {
229
    FUNCNAME("SubDomainSolver::solvePetscMatrix()");
230
231
232
233
234
235
236
  }


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

237
238
    VecDestroy(&rhsCoarseSpace);
    VecDestroy(&rhsInterior);
239
240
241
242
243
  }


  void SubDomainSolver::destroyMatrixData()
  {
244
245
246
247
248
249
250
251
    FUNCNAME("SubDomainSolver::destroyMatrixData()");

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

    KSPDestroy(&kspInterior);
252
253
254
255
256
  }


  void SubDomainSolver::solve(Vec &rhs, Vec &sol)
  {
257
258
259
260
    KSPSolve(kspInterior, rhs, sol);
  }


261
262
263
264
  void SubDomainSolver::solveGlobal(Vec &rhs, Vec &sol)
  {
    FUNCNAME("SubDomainSolver::solveGlobal()");

265
266
    ERROR_EXIT("BLUB!\n");

267
    Vec tmp;
268
    VecCreateSeq(PETSC_COMM_SELF, interiorMap->getRankDofs(level), &tmp);
269
270
271
272
273

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

274
    for (int i = 0; i < interiorMap->getRankDofs(level); i++)
275
276
277
278
279
280
281
282
283
284
      tmpValues[i] = rhsValues[i];

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

    KSPSolve(kspInterior, tmp, tmp);

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

285
    for (int i = 0; i < interiorMap->getRankDofs(level); i++) 
286
287
288
289
290
291
292
293
294
      rhsValues[i] = tmpValues[i];

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

    VecDestroy(&tmp);
  }


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
  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;
324
  }
325

326
}