PetscSolverGlobalBlockMatrix.cc 8.27 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23
//
// 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/PetscSolverGlobalBlockMatrix.h"
#include "parallel/StdMpi.h"
#include "parallel/MpiHelper.h"

namespace AMDiS {

  void PetscSolverGlobalBlockMatrix::fillPetscMatrix(Matrix<DOFMatrix*> *mat)
  {
    FUNCNAME("PetscSolverGlobalBlockMatrix::fillPetscMatrix()");

    TEST_EXIT_DBG(meshDistributor)("No mesh distributor object defined!\n");
24
    TEST_EXIT_DBG(interiorMap)("No parallel mapping object defined!\n");
25 26 27
    TEST_EXIT_DBG(mat)("No DOF matrix defined!\n");

    double wtime = MPI::Wtime();
28
    const FiniteElemSpace *feSpace = meshDistributor->getFeSpace(0);
29
    nComponents = mat->getNumRows();
30 31
    int nRankRows = (*interiorMap)[feSpace].nRankDofs;
    int nOverallRows = (*interiorMap)[feSpace].nOverallDofs;
32 33 34 35 36

#if (DEBUG != 0)
    MSG("Fill petsc matrix 1 needed %.5f seconds\n", MPI::Wtime() - wtime);
#endif

37 38 39 40 41 42 43 44 45 46 47 48 49 50 51
    if (nBlocks == -1) {
      nBlocks = nComponents;
      for (int i = 0; i < nBlocks; i++)
	componentInBlock[i] = i;
    }

    vector<int> compNthInBlock(nComponents, 0);
    vector<int> blockSize(nBlocks, 0);

    for (int i = 0; i < nComponents; i++) {
      compNthInBlock[i] = blockSize[componentInBlock[i]];
      blockSize[componentInBlock[i]]++;
    }

    nestMat.resize(nBlocks * nBlocks);
52 53 54

    // === Transfer values from DOF matrices to the PETSc matrix. === 

55 56
    for (int i = 0; i < nBlocks; i++)
      for (int j = 0; j < nBlocks; j++)
57 58 59 60 61 62
	MatCreateAIJ(mpiCommGlobal,
		     nRankRows * blockSize[i], nRankRows * blockSize[j],
		     nOverallRows * blockSize[i], nOverallRows * blockSize[j],
		     30 * blockSize[i], PETSC_NULL, 
		     30 * blockSize[j], PETSC_NULL,
		     &(nestMat[i * nBlocks + j]));
63
			
64 65 66
    for (int i = 0; i < nComponents; i++)
      for (int j = 0; j < nComponents; j++)
	if ((*mat)[i][j]) {
67 68 69
	  int idx = componentInBlock[i] * nBlocks + componentInBlock[j];
	  setDofMatrix(nestMat[idx], (*mat)[i][j], 
		       compNthInBlock[i], compNthInBlock[j]);
70 71
	}

72 73 74 75 76 77 78 79 80 81 82
    for (int i = 0; i < nBlocks; i++) {
      for (int j = 0; j < nBlocks; j++) {
	int idx = i * nBlocks + j;
	if (nestMat[idx]) {
	  MatAssemblyBegin(nestMat[idx], MAT_FINAL_ASSEMBLY);
	  MatAssemblyEnd(nestMat[idx], MAT_FINAL_ASSEMBLY);
	}
      }
    }	  
	

83
    MatCreateNest(mpiCommGlobal, 
84
		  nBlocks, PETSC_NULL, nBlocks, PETSC_NULL, 
85
		  &(nestMat[0]), &matIntInt);
86 87 88 89 90

#if (DEBUG != 0)
    MSG("Fill petsc matrix 2 needed %.5f seconds\n", MPI::Wtime() - wtime);
#endif

91 92
    MatAssemblyBegin(matIntInt, MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(matIntInt, MAT_FINAL_ASSEMBLY);
93 94

    // === Init PETSc solver. ===
95 96 97
    KSPCreate(mpiCommGlobal, &kspInterior);
    KSPSetOperators(kspInterior, matIntInt, matIntInt, SAME_NONZERO_PATTERN); 
    KSPSetFromOptions(kspInterior);
98

99 100 101 102 103 104 105 106 107 108 109
    MSG("Fill petsc matrix needed %.5f seconds\n", MPI::Wtime() - wtime);
  }


  void PetscSolverGlobalBlockMatrix::fillPetscRhs(SystemVector *vec)
  {
    FUNCNAME("PetscSolverGlobalBlockMatrix::fillPetscRhs()");

    TEST_EXIT_DBG(vec)("NO DOF vector defined!\n");

    nComponents = vec->getSize();
110
    const FiniteElemSpace *feSpace = meshDistributor->getFeSpace(0);
111 112
    int nRankRows = (*interiorMap)[feSpace].nRankDofs;
    int nOverallRows = (*interiorMap)[feSpace].nOverallDofs;
113

114
    nestVec.resize(nComponents);
115

116
    for (int i = 0; i < nComponents; i++) {
117
      VecCreateMPI(mpiCommGlobal, nRankRows, nOverallRows, &(nestVec[i]));
118

119 120 121 122 123 124
      setDofVector(nestVec[i], vec->getDOFVector(i));
      
      VecAssemblyBegin(nestVec[i]);
      VecAssemblyEnd(nestVec[i]);
    }

125 126
    VecCreateNest(mpiCommGlobal, nComponents, PETSC_NULL, 
		  &(nestVec[0]), &rhsInterior);
127

128 129
    VecAssemblyBegin(rhsInterior);
    VecAssemblyEnd(rhsInterior);
130 131 132 133 134 135 136 137
  }


  void PetscSolverGlobalBlockMatrix::solvePetscMatrix(SystemVector &vec, 
						      AdaptInfo *adaptInfo)
  {
    FUNCNAME("PetscSolverGlobalBlockMatrix::solvePetscMatrix()");

138 139
    KSPGetPC(kspInterior, &pcInterior);
    setBlockPreconditioner(pcInterior);
140

141
    const FiniteElemSpace *feSpace = meshDistributor->getFeSpace(0);
142
    VecDuplicate(rhsInterior, &petscSolVec);
Thomas Witkowski's avatar
Thomas Witkowski committed
143

144
    // PETSc.
145
    solve(rhsInterior, petscSolVec);
146 147 148 149

    // === Transfere values from PETSc's solution vectors to the DOF vectors. ===
    for (int i = 0; i < nComponents; i++) {
      DOFVector<double> &dofvec = *(vec.getDOFVector(i));
150 151

      Vec tmp;
Thomas Witkowski's avatar
Thomas Witkowski committed
152
      VecNestGetSubVec(petscSolVec, i, &tmp);
153

154
      int nRankDofs = (*interiorMap)[feSpace].nRankDofs;
155 156 157
      PetscScalar *vecPointer;
      VecGetArray(tmp, &vecPointer);

158
      DofMap& d = (*interiorMap)[feSpace].getMap();
159 160 161
      for (DofMap::iterator it = d.begin(); it != d.end(); ++it)
	if (it->second.local != -1)
	  dofvec[it->first] = vecPointer[it->second.local];
162

163 164
      VecRestoreArray(tmp, &vecPointer);
    }
165 166 167 168 169 170 171 172 173 174 175 176


    // === Synchronize DOFs at common DOFs, i.e., DOFs that correspond to ===
    // === more than one partition.                                       ===
    meshDistributor->synchVector(vec);
  }


  void PetscSolverGlobalBlockMatrix::destroyMatrixData()
  {
    FUNCNAME("PetscSolverGlobalBlockMatrix::destroyMatrixData()");

177 178 179 180
    for (unsigned int i = 0; i < nestMat.size(); i++)
      if (nestMat[i] != PETSC_NULL)
	MatDestroy(&(nestMat[i]));

181 182
    MatDestroy(&matIntInt);
    KSPDestroy(&kspInterior);
183 184 185
  }


186 187 188 189
  void PetscSolverGlobalBlockMatrix::destroyVectorData()
  {
    FUNCNAME("PetscSolverGlobalBlockMatrix::destroyVectorData()");

190
    VecDestroy(&rhsInterior);
191 192 193 194 195 196 197
    for (int i = 0; i < nComponents; i++)
      VecDestroy(&(nestVec[i]));

    VecDestroy(&petscSolVec);
  }


198 199 200 201
  void PetscSolverGlobalBlockMatrix::setDofMatrix(Mat& petscMat, 
						  DOFMatrix* mat,
						  int dispRowBlock, 
						  int dispColBlock)
202 203 204 205 206 207
  {
    FUNCNAME("PetscSolverGlobalBlockMatrix::setDofMatrix()");

    TEST_EXIT(mat)("No DOFMatrix!\n");
    TEST_EXIT(petscMat)("No PETSc matrix!\n");

208 209
    const FiniteElemSpace *feSpace = meshDistributor->getFeSpace(0);

210 211 212 213 214 215 216 217 218 219
    using mtl::tag::row; using mtl::tag::nz; using mtl::begin; using mtl::end;
    namespace traits = mtl::traits;
    typedef DOFMatrix::base_matrix_type Matrix;

    traits::col<Matrix>::type col(mat->getBaseMatrix());
    traits::const_value<Matrix>::type value(mat->getBaseMatrix());

    typedef traits::range_generator<row, Matrix>::type cursor_type;
    typedef traits::range_generator<nz, cursor_type>::type icursor_type;

220 221
    int dispRowIndex = (*interiorMap)[feSpace].nRankDofs * dispRowBlock;
    int dispColIndex = (*interiorMap)[feSpace].nRankDofs * dispColBlock;
222

223 224 225 226 227 228 229 230 231 232 233 234
    vector<int> cols;
    vector<double> values;
    cols.reserve(300);
    values.reserve(300);
    
    // === Traverse all rows of the dof matrix and insert row wise the values ===
    // === to the PETSc matrix.                                               ===

    for (cursor_type cursor = begin<row>(mat->getBaseMatrix()), 
	   cend = end<row>(mat->getBaseMatrix()); cursor != cend; ++cursor) {

      // Global index of the current row DOF.
235
      int rowIndex = (*interiorMap)[feSpace][*cursor].global + dispRowIndex;
236 237 238 239 240 241 242

      cols.clear();
      values.clear();
      
      for (icursor_type icursor = begin<nz>(cursor), icend = end<nz>(cursor); 
	   icursor != icend; ++icursor) {	
	// Global index of the current column index.
243
	int colIndex = (*interiorMap)[feSpace][col(*icursor)].global + dispColIndex;
244 245 246 247 248 249 250 251 252
	
	// Ignore all zero entries, expect it is a diagonal entry.
	if (value(*icursor) == 0.0 && rowIndex != colIndex)
	  continue;
	
	// Calculate the exact position of the column index in the PETSc matrix.
	cols.push_back(colIndex);
	values.push_back(value(*icursor));
      }
253 254 255

      MatSetValues(petscMat, 1, &rowIndex, cols.size(), 
 		   &(cols[0]), &(values[0]), ADD_VALUES);	
256 257 258 259 260
    }
  }
  
  
  void PetscSolverGlobalBlockMatrix::setDofVector(Vec& petscVec, 
261
						  DOFVector<double>* vec)
262 263 264
  {
    FUNCNAME("PetscSolverGlobalBlockMatrix::setDofVector()");

265 266
    const FiniteElemSpace *feSpace = meshDistributor->getFeSpace(0);

267 268 269
    // Traverse all used DOFs in the dof vector.
    DOFVector<double>::Iterator dofIt(vec, USED_DOFS);
    for (dofIt.reset(); !dofIt.end(); ++dofIt) {
270
      int index = (*interiorMap)[feSpace][dofIt.getDOFIndex()].global;
271 272 273 274 275 276 277
      double value = *dofIt;

      VecSetValues(petscVec, 1, &index, &value, ADD_VALUES);
    }
  }

}