PetscSolver.cc 5.33 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12
//
// 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.


13
#include "AMDiS.h"
Thomas Witkowski's avatar
Thomas Witkowski committed
14
#include "parallel/PetscSolver.h"
15
#include "parallel/StdMpi.h"
16
#include "parallel/MpiHelper.h"
Thomas Witkowski's avatar
Thomas Witkowski committed
17
#include "parallel/ParallelDofMapping.h"
18

19 20
namespace AMDiS {

Thomas Witkowski's avatar
Thomas Witkowski committed
21 22
  using namespace std;

23
  PetscSolver::PetscSolver()
24 25 26
    : ParallelCoarseSpaceMatVec(),
      dofMap(FESPACE_WISE),
      dofMapSd(FESPACE_WISE),
27
      kspPrefix(""),
28 29
      removeRhsNullspace(false),    
      hasConstantNullspace(false),
30
      isSymmetric(false),
31
      handleDirichletRows(true)
32
  {
33 34
    setDofMapping(&dofMap);

35 36 37 38
    string kspStr = "";
    Parameters::get("parallel->solver->petsc->ksp", kspStr);
    if (kspStr != "")
      PetscOptionsInsertString(kspStr.c_str());
39

Thomas Witkowski's avatar
Thomas Witkowski committed
40 41
    Parameters::get("parallel->remove rhs null space", removeRhsNullspace);
    Parameters::get("parallel->has constant null space", hasConstantNullspace);
42 43
    Parameters::get("parallel->nullspace->const in comp", 
		    constNullspaceComponent);
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
  void PetscSolver::init(vector<const FiniteElemSpace*> &fe0,
			 vector<const FiniteElemSpace*> &fe1)
  {
    FUNCNAME("PetscSolver::init()");

    TEST_EXIT(meshDistributor)("No mesh distributor object defined!\n");
    TEST_EXIT(fe0.size())("No component spaces in PETSc solver object!\n");
    TEST_EXIT(fe1.size())("No FE spaces in PETSc solver object!\n");

    componentSpaces = fe0;
    feSpaces = fe1;

    MeshLevelData& levelData = meshDistributor->getMeshLevelData();
    int nLevels = levelData.getLevelNumber();
    TEST_EXIT_DBG(nLevels >= 1)("Should not happen!\n");

    dofMap.init(levelData, componentSpaces, feSpaces);
    dofMap.setMpiComm(levelData.getMpiComm(0), 0);
    dofMap.setDofComm(meshDistributor->getDofComm());
    dofMap.clear();
    meshDistributor->registerDofMap(dofMap);

    if (nLevels > 1) {
      dofMapSd.init(levelData, componentSpaces, feSpaces);
      dofMapSd.setMpiComm(levelData.getMpiComm(1), 1);
      dofMapSd.setDofComm(meshDistributor->getDofCommSd());
      dofMapSd.clear();
      meshDistributor->registerDofMap(dofMapSd);
    }
  }


79 80 81 82 83 84 85 86
  void PetscSolver::fillPetscMatrix(DOFMatrix* mat)
  {
    Matrix<DOFMatrix*> sysMat(1, 1);
    sysMat[0][0] = mat;
    fillPetscMatrix(&sysMat);
  }


87 88 89
  void PetscSolver::solve(Vec &rhs, Vec &sol)
  {
    FUNCNAME("PetscSolver::solve()");
90

91 92 93 94
    PetscErrorCode solverError = KSPSolve(kspInterior, rhs, sol);
    if (solverError != 0) {
      AMDiS::finalize();
      exit(-1);
95 96 97
    }
  }

98

99 100 101 102 103 104 105 106
  void PetscSolver::solveGlobal(Vec &rhs, Vec &sol)
  {
    FUNCNAME("PetscSolver::solveGlobal()");

    ERROR_EXIT("Not implemented!\n");
  }


107 108 109 110 111 112
  void PetscSolver::copyVec(Vec& originVec, Vec& destVec, 
			    vector<int>& originIndex, vector<int>& destIndex)
  {
    FUNCNAME("PetscSolver::copyVec()");

    IS originIs, destIs;
113
    ISCreateGeneral(mpiCommGlobal, 
114 115 116 117 118
		    originIndex.size(), 
		    &(originIndex[0]),
		    PETSC_USE_POINTER,
		    &originIs);

119
    ISCreateGeneral(mpiCommGlobal, 
120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136
		    destIndex.size(), 
		    &(destIndex[0]),
		    PETSC_USE_POINTER,
		    &destIs);

    VecScatter scatter;
    VecScatterCreate(originVec, originIs, destVec, destIs, &scatter);
    VecScatterBegin(scatter, originVec, destVec,
		    INSERT_VALUES, SCATTER_FORWARD);
    VecScatterEnd(scatter, originVec, destVec,
		  INSERT_VALUES, SCATTER_FORWARD);

    ISDestroy(&originIs);
    ISDestroy(&destIs);    
    VecScatterDestroy(&scatter);
  }

Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
137 138 139 140 141 142 143 144

  bool PetscSolver::testMatrixSymmetric(Mat mat, bool advancedTest)
  {
    FUNCNAME("PetscSolver::testMatrixSymmetric()");

    Mat matTrans;
    MatTranspose(mat, MAT_INITIAL_MATRIX, &matTrans);

Thomas Witkowski's avatar
Thomas Witkowski committed
145 146
    bool isSym = true;

Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
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 187 188 189 190 191 192 193 194 195 196 197 198 199 200
    if (advancedTest) {
      int rowStart, rowEnd;
      MatGetOwnershipRange(mat, &rowStart, &rowEnd);

      PetscInt nCols, nColsTrans;
      const PetscInt *cols, *colsTrans;
      const PetscScalar *vals, *valsTrans;

      for (int i = rowStart; i < rowEnd; i++) {
	MatGetRow(mat, i, &nCols, &cols, &vals);
	MatGetRow(matTrans, i, &nColsTrans, &colsTrans, &valsTrans);

	if (nCols != nColsTrans) {
	  MSG("Symmetric test fails: mat row %d: %d nnz   mat col %d: %d nnz\n",
	      i, nCols, i, nColsTrans);
	  isSym = false;
	} else {
	  for (int j = 0; j < nCols; j++) {
	    if (cols[j] != colsTrans[j]) {
	      if (cols[j] > colsTrans[j]) {
		MSG("mat[%d][%d] does not exists: mat[%d][%d] = %e\n", 
		    i, colsTrans[j],
		    colsTrans[j], i, valsTrans[j]);
	      } else {
		MSG("mat[%d][%d] does not exists: mat[%d][%d] = %e\n", 
		    cols[j], i,
		    i, cols[j], vals[j]);
	      }
	      isSym = false;
	    } else {
	      double n = fabs(vals[j] - valsTrans[j]);
	      if (n > 1e-10) {
		MSG("value diff:  mat[%d][%d] = %e   mat[%d][%d] = %e\n",
		    i, cols[j], vals[j], colsTrans[j], i, valsTrans[j]);

		isSym = false;
	      }
	    }
	  }
	}

	MatRestoreRow(mat, i, &nCols, &cols, &vals);
	MatRestoreRow(matTrans, i, &nColsTrans, &colsTrans, &valsTrans);
      }
    } 
      
    MatAXPY(matTrans, -1.0, mat, DIFFERENT_NONZERO_PATTERN);
    double norm1, norm2;
    MatNorm(matTrans, NORM_FROBENIUS, &norm1);
    MatNorm(matTrans, NORM_INFINITY, &norm2);
    MatDestroy(&matTrans);
    
    MSG("Matrix norm test: %e %e\n", norm1, norm2);
    
Thomas Witkowski's avatar
Thomas Witkowski committed
201
    return (isSym && norm1 < 1e-10 && norm2 < 1e-10);
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
202
  }
203
}