PetscHelper.cc 7.93 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/PetscHelper.h"
14
#include "parallel/PetscSolver.h"
15
#include "Global.h"
16 17 18 19 20 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 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74

namespace AMDiS {

  namespace petsc_helper {

    using namespace std;
    
    void getMatLocalColumn(Mat mat, PetscMatCol &matCol)
    {
      PetscInt firstIndex, lastIndex;
      MatGetOwnershipRange(mat, &firstIndex, &lastIndex);
    
      PetscInt nCols;
      const PetscInt *cols;
      const PetscScalar *values;
      
      for (int row = firstIndex; row < lastIndex; row++) {
	MatGetRow(mat, row, &nCols, &cols, &values);
	
	for (int i = 0; i < nCols; i++) {
	  if (values[i] != 0.0) {
	    matCol[cols[i]].first.push_back(row - firstIndex);
	    matCol[cols[i]].second.push_back(values[i]);
	  }
	}

	MatRestoreRow(mat, row, &nCols, &cols, &values);
      }
    }


    void setMatLocalColumn(Mat mat, int column, Vec vec)
    {
      PetscInt firstIndex;
      MatGetOwnershipRange(mat, &firstIndex, PETSC_NULL);

      PetscInt vecSize;
      VecGetLocalSize(vec, &vecSize);

      PetscScalar *tmpValues;
      VecGetArray(vec, &tmpValues);
      for (int i  = 0; i < vecSize; i++)
	MatSetValue(mat, 
		    firstIndex + i,
		    column,
		    tmpValues[i],
		    ADD_VALUES);
      VecRestoreArray(vec, &tmpValues);
    }


    void getColumnVec(const SparseCol &matCol, Vec vec)
    {
      VecSet(vec, 0.0);
      VecSetValues(vec, matCol.first.size(), 
		   &(matCol.first[0]), &(matCol.second[0]), INSERT_VALUES);
      VecAssemblyBegin(vec);
      VecAssemblyEnd(vec);
    }
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


    void blockMatMatSolve(KSP ksp, Mat mat0, Mat &mat1)
    {
      // === We have to  calculate mat1 = ksp mat0:                       ===
      // ===    - get all local column vectors from mat0                  ===
      // ===    - solve with ksp for this column vector as the rhs vector ===
      // ===    - set the result to the corresponding column vector of    ===
      // ===      matrix mat1                                             ===
      
      PetscInt localRows, localCols, globalRows, globalCols;
      MatGetLocalSize(mat0, &localRows, &localCols);
      MatGetSize(mat0, &globalRows, &globalCols);

      MatCreateAIJ(PETSC_COMM_WORLD,
		   localRows, localCols, globalRows, globalCols,
		   150, PETSC_NULL, 150, PETSC_NULL, &mat1);
      MatSetOption(mat1, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
      
      // Transform matrix mat0 into a sparse column format.
      PetscMatCol mat0_cols;
      getMatLocalColumn(mat0, mat0_cols);
      
      Vec tmpVec;
      VecCreateSeq(PETSC_COMM_SELF, localRows, &tmpVec);
      
      // Solve for all column vectors of mat A_BPi and create matrix matK
      for (PetscMatCol::iterator it = mat0_cols.begin(); 
	   it != mat0_cols.end(); ++it) {
	getColumnVec(it->second, tmpVec);
	KSPSolve(ksp, tmpVec, tmpVec);
	setMatLocalColumn(mat1, it->first, tmpVec);
      }
      
      VecDestroy(&tmpVec);
      
      MatAssemblyBegin(mat1, MAT_FINAL_ASSEMBLY);
      MatAssemblyEnd(mat1, MAT_FINAL_ASSEMBLY);
    }


    void matNestConvert(Mat matNest, Mat &mat)  
    {
      FUNCNAME("matNestConvert()");

      PetscInt nestRows, nestCols;
      MatNestGetSize(matNest, &nestRows, &nestCols);

      TEST_EXIT(nestRows == 2 && nestCols == 2)
	("This funciton is only implemented for 2x2 nested matrices!\n");

      Mat mat00, mat01, mat10, mat11;
      MatNestGetSubMat(matNest, 0, 0, &mat00);
      MatNestGetSubMat(matNest, 0, 1, &mat01);
      MatNestGetSubMat(matNest, 1, 0, &mat10);
      MatNestGetSubMat(matNest, 1, 1, &mat11);

      int nRankRows0, nOverallRows0;
      MatGetLocalSize(mat00, &nRankRows0, PETSC_NULL);
      MatGetSize(mat00, &nOverallRows0, PETSC_NULL);

      int nRankRows1, nOverallRows1;
      MatGetLocalSize(mat11, &nRankRows1, PETSC_NULL);
      MatGetSize(mat11, &nOverallRows1, PETSC_NULL);

      int firstRow0;
      MatGetOwnershipRange(mat00, &firstRow0, PETSC_NULL);

      int firstRow1;
      MatGetOwnershipRange(mat11, &firstRow1, PETSC_NULL);

      int nRankRows = nRankRows0 + nRankRows1;
      int nOverallRows = nOverallRows0 + nOverallRows1;
      int firstRow = firstRow0 + firstRow1;

150 151 152 153 154 155 156 157 158 159 160 161 162 163 164
      int mpiSize = MPI::COMM_WORLD.Get_size();
      int mpiRank = MPI::COMM_WORLD.Get_rank();
      vector<int> allFirstRow0(mpiSize + 1, 0);
      vector<int> allFirstRow1(mpiSize + 1, 0);
      MPI::COMM_WORLD.Allgather(&nRankRows0, 1, MPI_INT, &(allFirstRow0[1]), 1, MPI_INT);
      MPI::COMM_WORLD.Allgather(&nRankRows1, 1, MPI_INT, &(allFirstRow1[1]), 1, MPI_INT);

      for (int i = 1; i < mpiSize + 1; i++) {
	allFirstRow0[i] += allFirstRow0[i - 1];
	allFirstRow1[i] += allFirstRow1[i - 1];
      }

      TEST_EXIT_DBG(allFirstRow0[mpiRank] == firstRow0)("Should not happen!\n");
      TEST_EXIT_DBG(allFirstRow1[mpiRank] == firstRow1)("Should not happen!\n");

165 166 167
      MatCreateAIJ(PETSC_COMM_WORLD, 
		   nRankRows, nRankRows,
		   nOverallRows, nOverallRows,
168 169 170 171 172 173
		   25, PETSC_NULL, 25, PETSC_NULL, &mat);
      MatSetOption(mat, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);

      PetscInt nCols;
      const PetscInt *cols;
      const PetscScalar *vals;
174 175

      for (int i = 0; i < nRankRows0; i++) {
176
	PetscInt newRowIndex = firstRow + i;
177

178 179 180 181 182 183 184 185 186 187 188 189 190 191 192
	MatGetRow(mat00, firstRow0 + i, &nCols, &cols, &vals);
	for (int j = 0; j < nCols; j++) {
	  int rank = -1;
	  for (int k = 0; k < mpiSize; k++) {
	    if (cols[j] >= allFirstRow0[k] && cols[j] < allFirstRow0[k + 1]) {
	      rank = k;
	      break;
	    }
	  }
	  TEST_EXIT_DBG(rank != -1)("Should not happen!\n");

	  PetscInt newColIndex = cols[j] + allFirstRow1[rank];
	  MatSetValue(mat, newRowIndex, newColIndex, vals[j], INSERT_VALUES);
	}
	MatRestoreRow(mat00, firstRow0 + i, &nCols, &cols, &vals);
193

194
	MatGetRow(mat01, firstRow0 + i, &nCols, &cols, &vals);
195
	for (int j = 0; j < nCols; j++) {
196 197 198 199 200 201 202 203 204 205 206
	  int rank = -1;
	  for (int k = 0; k < mpiSize; k++) {
	    if (cols[j] >= allFirstRow1[k] && cols[j] < allFirstRow1[k + 1]) {
	      rank = k;
	      break;
	    }
	  }
	  TEST_EXIT_DBG(rank != -1)("Should not happen!\n");

	  PetscInt newColIndex = cols[j] + allFirstRow0[rank + 1];
	  MatSetValue(mat, newRowIndex, newColIndex, vals[j], INSERT_VALUES);
207
	}
208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223
	MatRestoreRow(mat01, firstRow0 + i, &nCols, &cols, &vals);	
      }

      for (int i = 0; i < nRankRows1; i++) {
	PetscInt newRowIndex = firstRow + nRankRows0 + i;

	MatGetRow(mat10, firstRow1 + i, &nCols, &cols, &vals);
	for (int j = 0; j < nCols; j++) {
	  int rank = -1;
	  for (int k = 0; k < mpiSize; k++) {
	    if (cols[j] >= allFirstRow0[k] && cols[j] < allFirstRow0[k + 1]) {
	      rank = k;
	      break;
	    }
	  }
	  TEST_EXIT_DBG(rank != -1)("Should not happen!\n");
224

225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244
	  PetscInt newColIndex = cols[j] + allFirstRow1[rank];
	  MatSetValue(mat, newRowIndex, newColIndex, vals[j], INSERT_VALUES);
	}
	MatRestoreRow(mat10, firstRow1 + i, &nCols, &cols, &vals);

	MatGetRow(mat11, firstRow1 + i, &nCols, &cols, &vals);
	for (int j = 0; j < nCols; j++) {
	  int rank = -1;
	  for (int k = 0; k < mpiSize; k++) {
	    if (cols[j] >= allFirstRow1[k] && cols[j] < allFirstRow1[k + 1]) {
	      rank = k;
	      break;
	    }
	  }
	  TEST_EXIT_DBG(rank != -1)("Should not happen!\n");

	  PetscInt newColIndex = cols[j] + allFirstRow0[rank + 1];
	  MatSetValue(mat, newRowIndex, newColIndex, vals[j], INSERT_VALUES);
	}
	MatRestoreRow(mat11, firstRow1 + i, &nCols, &cols, &vals);
245 246 247 248 249 250
      }

      MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);
      MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);
    }

251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269

    void setSolver(KSP ksp, 
		   string kspPrefix,
		   KSPType kspType, 
		   PCType pcType, 
		   PetscReal rtol,
		   PetscReal atol,
		   PetscInt maxIt)
    {
      KSPSetType(ksp, kspType);
      KSPSetTolerances(ksp, rtol, atol, PETSC_DEFAULT, maxIt);
      if (kspPrefix != "")
	KSPSetOptionsPrefix(ksp, kspPrefix.c_str());
      KSPSetFromOptions(ksp);

      PC pc;
      KSPGetPC(ksp, &pc);
      PCSetType(pc, pcType);
    }
270
  }
271

272
}