PetscSolver.cc 6.54 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(string name)
24
    : ParallelCoarseSpaceSolver(name),
25
      dofMap(FESPACE_WISE, true),
26
      dofMapSubDomain(FESPACE_WISE, true),
Thomas Witkowski's avatar
Thomas Witkowski committed
27
      parallelDofMappingsRegistered(false),
28
      kspPrefix(""),
29 30
      removeRhsNullspace(false),    
      hasConstantNullspace(false),
31
      isSymmetric(false),
32
      handleDirichletRows(true)
33
  {
34 35
    setDofMapping(&dofMap);

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

Thomas Witkowski's avatar
Thomas Witkowski committed
41 42
    Parameters::get("parallel->remove rhs null space", removeRhsNullspace);
    Parameters::get("parallel->has constant null space", hasConstantNullspace);
43 44
    Parameters::get("parallel->nullspace->const in comp", 
		    constNullspaceComponent);
45 46 47
  }


Thomas Witkowski's avatar
Thomas Witkowski committed
48 49
  PetscSolver::~PetscSolver()
  {
50
    if (parallelDofMappingsRegistered)
51
      meshDistributor->removeDofMap(dofMap);
Thomas Witkowski's avatar
Thomas Witkowski committed
52 53 54
  }


55
  void PetscSolver::init(vector<const FiniteElemSpace*> &fe0,
Thomas Witkowski's avatar
bo eh  
Thomas Witkowski committed
56 57
			 vector<const FiniteElemSpace*> &fe1,
			 bool createGlobalMapping)
58 59 60 61 62 63 64 65 66 67 68
  {
    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();
69
    int nLevels = levelData.getNumberOfLevels();
70 71
    TEST_EXIT_DBG(nLevels >= 1)("Should not happen!\n");

Thomas Witkowski's avatar
bo eh  
Thomas Witkowski committed
72
    if (createGlobalMapping) {
Thomas Witkowski's avatar
Thomas Witkowski committed
73 74
      parallelDofMappingsRegistered = true;

75
      dofMap.init(componentSpaces, feSpaces);
76 77
      dofMap.setMpiComm(levelData.getMpiComm(meshLevel));
      dofMap.setDofComm(meshDistributor->getDofComm(meshLevel));
Thomas Witkowski's avatar
bo eh  
Thomas Witkowski committed
78
      dofMap.clear();
79 80
      meshDistributor->registerDofMap(dofMap);     

81 82 83 84 85 86 87
      if (meshLevel + 1 < nLevels && 
	  levelData.getMpiComm(meshLevel + 1) != MPI::COMM_SELF) {
	dofMapSubDomain.init(componentSpaces, feSpaces);
	dofMapSubDomain.setMpiComm(levelData.getMpiComm(meshLevel + 1));
	dofMapSubDomain.setDofComm(meshDistributor->getDofComm(meshLevel + 1));
	dofMapSubDomain.clear();
	meshDistributor->registerDofMap(dofMapSubDomain);
Thomas Witkowski's avatar
bo eh  
Thomas Witkowski committed
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
 int PetscSolver::solveSystem(const SolverMatrix<Matrix<DOFMatrix*> >& A,
			   SystemVector& x,
			   SystemVector& b,
			   bool createMatrixData,
			   bool storeMatrixData)
 {
	 MPI::COMM_WORLD.Barrier();
	 double wtime = MPI::Wtime();

	 if (createMatrixData)
		 fillPetscMatrix(const_cast< Matrix< DOFMatrix* >* >(A.getOriginalMat()));

	 fillPetscRhs(&b);

	 INFO(info, 8)("creation of parallel data structures needed %.5f seconds\n", 
			 MPI::Wtime() - wtime);
	 wtime = MPI::Wtime();

	 solvePetscMatrix(x, NULL);

	 INFO(info, 8)("solution of discrete system needed %.5f seconds\n", 
			 MPI::Wtime() - wtime);

	 destroyVectorData();
116 117 118

	 if (!storeMatrixData)
		 destroyMatrixData();
119 120
	 return 0;
 }
121

122 123 124 125 126 127 128 129
  void PetscSolver::fillPetscMatrix(DOFMatrix* mat)
  {
    Matrix<DOFMatrix*> sysMat(1, 1);
    sysMat[0][0] = mat;
    fillPetscMatrix(&sysMat);
  }


130 131 132
  void PetscSolver::solve(Vec &rhs, Vec &sol)
  {
    FUNCNAME("PetscSolver::solve()");
133

134 135 136 137
    PetscErrorCode solverError = KSPSolve(kspInterior, rhs, sol);
    if (solverError != 0) {
      AMDiS::finalize();
      exit(-1);
138 139 140
    }
  }

141

142 143 144 145
  void PetscSolver::solveGlobal(Vec &rhs, Vec &sol)
  {
    FUNCNAME("PetscSolver::solveGlobal()");

146 147 148 149 150 151
    int s, ls;
    VecGetSize(rhs, &s);
    VecGetLocalSize(rhs, &ls);

    MSG("Solve global %d %d\n", ls, s);

152 153 154 155
    ERROR_EXIT("Not implemented!\n");
  }


156 157 158 159 160 161
  void PetscSolver::copyVec(Vec& originVec, Vec& destVec, 
			    vector<int>& originIndex, vector<int>& destIndex)
  {
    FUNCNAME("PetscSolver::copyVec()");

    IS originIs, destIs;
162
    ISCreateGeneral(domainComm, 
163 164 165 166 167
		    originIndex.size(), 
		    &(originIndex[0]),
		    PETSC_USE_POINTER,
		    &originIs);

168
    ISCreateGeneral(domainComm, 
169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185
		    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
186 187 188 189 190 191 192 193

  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
194 195
    bool isSym = true;

Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
196 197 198 199 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 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249
    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
250
    return (isSym && norm1 < 1e-10 && norm2 < 1e-10);
Thomas Witkowski's avatar
Blub  
Thomas Witkowski committed
251
  }
252
}