PetscSolverNavierStokes.cc 9.72 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19
//
// 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/PetscSolverNavierStokes.h"

namespace AMDiS {

  using namespace std;

20

Thomas Witkowski's avatar
Thomas Witkowski committed
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
  PetscErrorCode pcNs(PC pc, Vec x, Vec y)
  {
    void *ctx;
    PCShellGetContext(pc, &ctx);
    MatShellNavierStokesSchur* data = static_cast<MatShellNavierStokesSchur*>(ctx);
    
    Vec velocity, pressure;
    VecGetSubVector(x, data->isVelocity, &velocity);
    VecGetSubVector(x, data->isPressure, &pressure);

    Vec velocityTmp;
    VecDuplicate(velocity, &velocityTmp);

    Vec pressureTmp;
    VecDuplicate(pressure, &pressureTmp);

#if 0
    MatMult(data->A01, pressure, velocityTmp);
    VecAXPY(velocity, 1.0, velocityTmp);
    KSPSolve(data->kspVelocity, velocity, velocity);
    VecScale(pressure, -1.0);
#else
    KSPSolve(data->kspVelocity, velocity, velocityTmp);
    MatMult(data->A10, velocityTmp, pressureTmp);
    VecAXPY(pressure, -1.0, pressureTmp);

    KSPSolve(data->kspLaplace, pressure, pressure);
    MatMult(data->matConDif, pressure, pressureTmp);
    KSPSolve(data->kspMass, pressureTmp, pressure);

    MatMult(data->A01, pressure, velocityTmp);
    VecAXPY(velocity, -1.0, velocityTmp);
    KSPSolve(data->kspVelocity, velocity, velocity);
#endif


    VecRestoreSubVector(x, data->isVelocity, &velocity);
    VecRestoreSubVector(x, data->isPressure, &pressure);

    VecCopy(x, y);

    VecDestroy(&velocityTmp);
    VecDestroy(&pressureTmp);
  }

66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94
  int petscMultNavierStokesSchur(Mat mat, Vec x, Vec y)
  {
    void *ctx;
    MatShellGetContext(mat, &ctx);
    MatShellNavierStokesSchur* data = static_cast<MatShellNavierStokesSchur*>(ctx);

    switch (data->solverMode) {
    case 0:
      {
	Vec vec0, vec1;
	VecDuplicate(x, &vec0);
	MatGetVecs(data->A00, &vec1, PETSC_NULL);
	
	MatMult(data->A11, x, y);
	MatMult(data->A01, x, vec1);
	KSPSolve(data->kspVelocity, vec1, vec1);
	MatMult(data->A10, vec1, y);
	VecAYPX(y, -1.0, vec0);
	
	VecDestroy(&vec0);
	VecDestroy(&vec1);
      }
      break;
    case 1:
      {
	Vec vec0, vec1;
	VecDuplicate(x, &vec0);
	VecDuplicate(x, &vec1);

Thomas Witkowski's avatar
Thomas Witkowski committed
95 96 97 98
	//	KSPSolve(data->kspLaplace, x, vec0);
	//	MatMult(data->matConDif, vec0, vec1);
	//KSPSolve(data->kspMass, vec0, y);
	VecSet(y, 0.0);
99 100 101 102 103 104 105 106 107

	VecDestroy(&vec0);
	VecDestroy(&vec1);
      }
      break;
    default:
      ERROR_EXIT("Wrong solver mode %d\n", data->solverMode);
    }
  }
108
  
109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125

  PetscSolverNavierStokes::PetscSolverNavierStokes(string name)
    : PetscSolverGlobalMatrix(name),
      pressureComponent(-1),
      massMatrixSolver(NULL),
      laplaceMatrixSolver(NULL)
  {
    Parameters::get(initFileStr + "->stokes->pressure component", 
		    pressureComponent);

    TEST_EXIT(pressureComponent >= 0)
      ("For using PETSc stokes solver you must define a pressure component!\n");

    TEST_EXIT(pressureComponent == 2)("Fixed for pressure component = 2\n");
  }


126
  void PetscSolverNavierStokes::initSolver(KSP &ksp)
127 128 129 130 131 132 133 134 135 136
  {
    FUNCNAME("PetscSolverNavierStokes::initSolver()");

    KSPCreate(mpiCommGlobal, &ksp);
    KSPSetOperators(ksp, getMatInterior(), getMatInterior(), 
		    SAME_NONZERO_PATTERN); 
    KSPSetTolerances(ksp, 0.0, 1e-8, PETSC_DEFAULT, PETSC_DEFAULT);
    KSPSetType(ksp, KSPGMRES);
    KSPSetOptionsPrefix(ksp, "ns_");
    KSPSetFromOptions(ksp);
137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155
    KSPMonitorSet(ksp, KSPMonitorTrueResidualNorm, PETSC_NULL, PETSC_NULL);

    
    // === Create null space information. ===

    Vec nullSpaceBasis;
    VecDuplicate(getVecSolInterior(), &nullSpaceBasis);

    SystemVector basisVec("tmp", componentSpaces, componentSpaces.size(), true);
    basisVec.set(0.0);
    basisVec.getDOFVector(pressureComponent)->set(1.0);

    setDofVector(nullSpaceBasis, basisVec, true);
    VecAssemblyBegin(nullSpaceBasis);
    VecAssemblyEnd(nullSpaceBasis);
    VecNormalize(nullSpaceBasis, PETSC_NULL);

    MatNullSpace matNullSpace;
    MatNullSpaceCreate(mpiCommGlobal, PETSC_FALSE, 1, &nullSpaceBasis, &matNullSpace);
156 157 158 159 160 161 162
  }


  void PetscSolverNavierStokes::initPreconditioner(PC pc)
  {
    FUNCNAME("PetscSolverNavierStokes::initPreconditioner()");

163

Thomas Witkowski's avatar
Thomas Witkowski committed
164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182
    PCSetType(pc, PCSHELL);
    PCShellSetApply(pc, pcNs);
    PCShellSetContext(pc, &matShellContext);

    interiorMap->createIndexSet(matShellContext.isVelocity, 0, 2);
    interiorMap->createIndexSet(matShellContext.isPressure, 2, 1);
    MatGetSubMatrix(getMatInterior(), matShellContext.isVelocity, 
		    matShellContext.isVelocity, MAT_INITIAL_MATRIX, &matShellContext.A00);
    MatGetSubMatrix(getMatInterior(), matShellContext.isVelocity, 
		    matShellContext.isPressure, MAT_INITIAL_MATRIX, &matShellContext.A01);
    MatGetSubMatrix(getMatInterior(), matShellContext.isPressure, 
		    matShellContext.isVelocity, MAT_INITIAL_MATRIX, &matShellContext.A10);
    MatGetSubMatrix(getMatInterior(), matShellContext.isPressure, 
		    matShellContext.isPressure, MAT_INITIAL_MATRIX, &matShellContext.A11);

    KSPCreate(mpiCommGlobal, &(matShellContext.kspVelocity));
    KSPSetOperators(matShellContext.kspVelocity, matShellContext.A00, matShellContext.A00, SAME_NONZERO_PATTERN);
    //  KSPSetType(matShellContext.kspVelocity, KSPPREONLY);
    KSPSetInitialGuessNonzero(matShellContext.kspVelocity, PETSC_TRUE);
183
    PC pcSub;
Thomas Witkowski's avatar
Thomas Witkowski committed
184
    KSPGetPC(matShellContext.kspVelocity, &pcSub);
185 186
    PCSetType(pcSub, PCLU);
    PCFactorSetMatSolverPackage(pcSub, MATSOLVERMUMPS);
Thomas Witkowski's avatar
Thomas Witkowski committed
187
    
188

189
    // === Mass matrix solver ===
190

191 192 193 194 195 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
    MSG("Initialize mass matrix solver ...\n");

    const FiniteElemSpace *pressureFeSpace = componentSpaces[pressureComponent];
    DOFMatrix massMatrix(pressureFeSpace, pressureFeSpace);
    Operator massOp(pressureFeSpace, pressureFeSpace);
    Simple_ZOT zot;
    massOp.addTerm(&zot);
    massMatrix.assembleOperator(massOp);

    massMatrixSolver = createSubSolver(pressureComponent, "mass_");
    massMatrixSolver->fillPetscMatrix(&massMatrix);

    MSG("... OK!\n");


    // === Laplace matrix solver ===

    MSG("Initialize laplace matrix solver ...\n");

    DOFMatrix laplaceMatrix(pressureFeSpace, pressureFeSpace);
    Operator laplaceOp(pressureFeSpace, pressureFeSpace);
    Simple_SOT sot;
    laplaceOp.addTerm(&sot);
    laplaceMatrix.assembleOperator(laplaceOp);

    laplaceMatrixSolver = createSubSolver(pressureComponent, "laplace_");
    laplaceMatrixSolver->fillPetscMatrix(&laplaceMatrix);

    MSG("... OK!\n");


    // === Create convection-diffusion operator ===

    MSG("Initialize pressure convection-diffusion operator ...\n");

    double timestep = 1.0;
    Parameters::get("navierstokes->adapt->timestep", timestep);
    timestep = 1.0 / timestep;

    double nu = 0.01;
    Parameters::get("nu", nu);

    DOFMatrix conDifMatrix(pressureFeSpace, pressureFeSpace);
    Operator conDifOp(pressureFeSpace, pressureFeSpace);
    Simple_ZOT conDif0(timestep);
    conDifOp.addTerm(&conDif0);
    Simple_SOT conDif1(nu);
    conDifOp.addTerm(&conDif1);
    Vector_FOT conDif2(0);
Thomas Witkowski's avatar
Thomas Witkowski committed
240
    conDifOp.addTerm(&conDif2, GRD_PHI);
241
    Vector_FOT conDif3(1);
Thomas Witkowski's avatar
Thomas Witkowski committed
242
    conDifOp.addTerm(&conDif3, GRD_PHI);
243 244 245 246 247 248 249 250 251 252 253 254

    conDifMatrix.assembleOperator(conDifOp);

    conDifMatrixSolver = createSubSolver(pressureComponent, "condif_");
    conDifMatrixSolver->fillPetscMatrix(&conDifMatrix);

    MSG("... OK!\n");

    matShellContext.kspMass = massMatrixSolver->getSolver();
    matShellContext.kspLaplace = laplaceMatrixSolver->getSolver();
    matShellContext.matConDif = conDifMatrixSolver->getMatInterior();    

Thomas Witkowski's avatar
Thomas Witkowski committed
255

256 257 258
    KSPSetType(matShellContext.kspMass, KSPRICHARDSON);
    KSPSetTolerances(matShellContext.kspMass, 0, 1e-13, 1e+3, 1);
    KSPGetPC(matShellContext.kspMass, &pcSub);
259 260
    PCSetType(pcSub, PCLU);
    PCFactorSetMatSolverPackage(pcSub, MATSOLVERMUMPS);
261 262 263
    KSPSetFromOptions(matShellContext.kspMass);

    KSPSetType(matShellContext.kspLaplace, KSPCG);
Thomas Witkowski's avatar
Thomas Witkowski committed
264
    KSPSetInitialGuessNonzero(matShellContext.kspLaplace, PETSC_TRUE);
265 266 267 268 269
    KSPSetTolerances(matShellContext.kspLaplace, 0, 1e-13, 1e+3, 10000);
    KSPGetPC(matShellContext.kspLaplace, &pcSub);
    PCSetType(pcSub, PCJACOBI);
    KSPSetFromOptions(matShellContext.kspLaplace);

Thomas Witkowski's avatar
Thomas Witkowski committed
270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323



#if 0
    vector<int> velocityComponents;
    velocityComponents.push_back(0);
    velocityComponents.push_back(1);

    PCSetType(pc, PCFIELDSPLIT);
    PCFieldSplitSetType(pc, PC_COMPOSITE_SCHUR);
    PCFieldSplitSetSchurFactType(pc, PC_FIELDSPLIT_SCHUR_FACT_FULL);
    createFieldSplit(pc, "velocity", velocityComponents);
    createFieldSplit(pc, "pressure", pressureComponent);

    KSPSetUp(kspInterior);

    KSP *subKsp;
    int nSubKsp;
    PCFieldSplitGetSubKSP(pc, &nSubKsp, &subKsp);

    TEST_EXIT(nSubKsp == 2)
      ("Wrong numer of KSPs inside of the fieldsplit preconditioner!\n");

    KSP kspVelocity = subKsp[0];
    KSP kspSchur = subKsp[1];
    PetscFree(subKsp);

    Mat A00, A01, A10, A11;
    PCFieldSplitGetSchurBlocks(pc, &A00, &A01, &A10, &A11);
    matShellContext.A00 = A00;
    matShellContext.A01 = A01;
    matShellContext.A10 = A10;
    matShellContext.A11 = A11;
    matShellContext.kspVelocity = kspVelocity;

    Mat matShell;
    MatDuplicate(A11, MAT_DO_NOT_COPY_VALUES, &matShell);
    MatSetType(matShell, MATSHELL);
    MatShellSetContext(matShell, &matShellContext);
    MatShellSetOperation(matShell, MATOP_MULT, (void(*)(void))petscMultNavierStokesSchur);

    KSPSetType(kspVelocity, KSPPREONLY);
    PC pcSub;
    KSPGetPC(kspVelocity, &pcSub);
    PCSetType(pcSub, PCLU);
    PCFactorSetMatSolverPackage(pcSub, MATSOLVERMUMPS);


    KSPSetType(kspSchur, KSPGMRES);
    KSPSetTolerances(kspSchur, 0, 1e-14, 1e+3, 1000);
    KSPSetOperators(kspSchur, matShell, matShell, SAME_NONZERO_PATTERN);
    KSPGetPC(kspSchur, &pcSub);
    PCSetType(pcSub, PCNONE);

324
    matShellContext.solverMode = 1;
325 326


327 328
    KSPSetFromOptions(kspVelocity);
    KSPSetFromOptions(kspSchur);
Thomas Witkowski's avatar
Thomas Witkowski committed
329
#endif
330 331 332
  }

}