SubDomainSolver.cc 11.5 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 44 45 46 47
  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);
    }
  }


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

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

54 55 56 57
    int nRowsRankInterior = interiorMap->getRankDofs();
    int nRowsOverallInterior = interiorMap->getOverallDofs();
    int nRowsRankCoarse = coarseSpaceMap->getRankDofs();
    int nRowsOverallCoarse = coarseSpaceMap->getOverallDofs();
58

59
    bool multilevel = false;
60
    if (mpiCommInterior.Get_size() == 1) {
61 62
      nGlobalOverallInterior = nRowsOverallInterior;

63 64 65
      MatCreateSeqAIJ(mpiCommInterior, nRowsRankInterior, nRowsRankInterior,
		      60, PETSC_NULL, &matIntInt);
    } else {
66 67
      multilevel = true;

68 69 70 71 72
      MatCreateMPIAIJ(mpiCommInterior, 
		      nRowsRankInterior, nRowsRankInterior,
		      nRowsOverallInterior, nRowsOverallInterior,
		      60, PETSC_NULL, 60, PETSC_NULL, &matIntInt);
    }
73 74 75 76 77 78 79 80

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

    MatCreateMPIAIJ(mpiCommCoarseSpace,
		    nRowsRankCoarse, nRowsRankInterior,
81
		    nRowsOverallCoarse, nGlobalOverallInterior,
82 83 84 85
		    60, PETSC_NULL, 60, PETSC_NULL, &matCoarseInt);

    MatCreateMPIAIJ(mpiCommCoarseSpace,
		    nRowsRankInterior, nRowsRankCoarse,
86
		    nGlobalOverallInterior, nRowsOverallCoarse,
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
		    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) {
157
	    int rowIndex = coarseSpaceMap->getMatIndex(i, *cursor);
158
	    for (unsigned int k = 0; k < cols.size(); k++)
159
	      cols[k] = coarseSpaceMap->getMatIndex(j, cols[k]);
160 161 162 163 164

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

	    if (colsOther.size()) {
165 166 167 168
	      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
169
		for (unsigned int k = 0; k < colsOther.size(); k++)
170 171 172
		  colsOther[k] = 
		    interiorMap->getMatIndex(j, colsOther[k]) + rStartInterior;
	      }
173 174 175 176 177
 	      
	      MatSetValues(matCoarseInt, 1, &rowIndex, colsOther.size(),
 			   &(colsOther[0]), &(valuesOther[0]), ADD_VALUES);
	    }
	  } else {
178 179 180 181 182 183 184
	    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
185
	      else
186 187
		cols[k] = interiorMap->getMatIndex(j, cols[k]);
	    }
188 189 190 191 192
	    
  	    MatSetValues(matIntInt, 1, &localRowIndex, cols.size(),
  			 &(cols[0]), &(values[0]), ADD_VALUES);

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

Thomas Witkowski's avatar
Thomas Witkowski committed
195
	      if (multilevel)
196 197
		globalRowIndex += rStartInterior;

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

  	      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);
233 234 235 236
    if (multilevel == false)
      PCFactorSetMatSolverPackage(pcInterior, MATSOLVERUMFPACK);
    else
      PCFactorSetMatSolverPackage(pcInterior, MATSOLVERMUMPS);
237
    KSPSetFromOptions(kspInterior);
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
238 239


Thomas Witkowski's avatar
Thomas Witkowski committed
240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262
    if (multilevel == false) {
      PetscViewer matview;
      PetscViewerBinaryOpen(mpiCommCoarseSpace, "mat_coarse_int.dat", 
			    FILE_MODE_WRITE, &matview);
      MatView(matCoarseInt, matview);
      PetscViewerDestroy(&matview);

      PetscViewerBinaryOpen(mpiCommCoarseSpace, "mat_int_coarse.dat", 
			    FILE_MODE_WRITE, &matview);
      MatView(matIntCoarse, matview);
      PetscViewerDestroy(&matview);      
    } else {
      PetscViewer matview;
      PetscViewerBinaryOpen(mpiCommCoarseSpace, "mat_coarse_int.dat", 
			    FILE_MODE_WRITE, &matview);
      MatView(matCoarseInt, matview);
      PetscViewerDestroy(&matview);

      PetscViewerBinaryOpen(mpiCommCoarseSpace, "mat_int_coarse.dat", 
			    FILE_MODE_WRITE, &matview);
      MatView(matIntCoarse, matview);
      PetscViewerDestroy(&matview);
    }
263 264 265 266 267
  }


  void SubDomainSolver::fillPetscRhs(SystemVector *vec)
  {
268 269 270
    FUNCNAME("SubDomainSolver::fillPetscRhs()");

    VecCreateMPI(mpiCommCoarseSpace, 
271 272
		 coarseSpaceMap->getRankDofs(), 
		 coarseSpaceMap->getOverallDofs(),
273 274 275
		 &rhsCoarseSpace);

    VecCreateMPI(mpiCommCoarseSpace, 
276
		 interiorMap->getRankDofs(), 
277
		 nGlobalOverallInterior,
278 279 280 281 282 283 284
		 &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();
Thomas Witkowski's avatar
Thomas Witkowski committed
285 286 287
	if (isCoarseSpace(feSpace, index)) {	  
	  index = coarseSpaceMap->getMatIndex(i, index);
	  VecSetValue(rhsCoarseSpace, index, *dofIt, ADD_VALUES);
288
	} else {
Thomas Witkowski's avatar
Thomas Witkowski committed
289 290
	  index = interiorMap->getMatIndex(i, index) + rStartInterior;
	  VecSetValue(rhsInterior, index, *dofIt, ADD_VALUES);
291 292 293 294 295 296 297 298 299
	}      
      }
    }

    VecAssemblyBegin(rhsCoarseSpace);
    VecAssemblyEnd(rhsCoarseSpace);

    VecAssemblyBegin(rhsInterior);
    VecAssemblyEnd(rhsInterior);
300 301 302 303 304
  }


  void SubDomainSolver::solvePetscMatrix(SystemVector &vec, AdaptInfo *adaptInfo)
  {
305
    FUNCNAME("SubDomainSolver::solvePetscMatrix()");
306 307 308 309 310 311 312
  }


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

313 314
    VecDestroy(&rhsCoarseSpace);
    VecDestroy(&rhsInterior);
315 316 317 318 319
  }


  void SubDomainSolver::destroyMatrixData()
  {
320 321 322 323 324 325 326 327
    FUNCNAME("SubDomainSolver::destroyMatrixData()");

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

    KSPDestroy(&kspInterior);
328 329 330 331 332
  }


  void SubDomainSolver::solve(Vec &rhs, Vec &sol)
  {
333 334 335 336
    KSPSolve(kspInterior, rhs, sol);
  }


337 338 339 340
  void SubDomainSolver::solveGlobal(Vec &rhs, Vec &sol)
  {
    FUNCNAME("SubDomainSolver::solveGlobal()");

341 342 343
    int ml = 0;
    Parameters::get("parallel->multi level test", ml);

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

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

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

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

    KSPSolve(kspInterior, tmp, tmp);

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

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

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

    VecDestroy(&tmp);
  }


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
  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;
407
  }
408

409
}