PetscSolverCahnHilliard_DD.cc 7.68 KB
Newer Older
1 2 3 4 5 6 7
/******************************************************************************
 *
 * AMDiS - Adaptive multidimensional simulations
 *
 * Copyright (C) 2013 Dresden University of Technology. All Rights Reserved.
 * Web: https://fusionforge.zih.tu-dresden.de/projects/amdis
 *
8
 * Authors:
9 10 11 12 13 14 15 16 17
 * Simon Vey, Thomas Witkowski, Andreas Naumann, Simon Praetorius, et al.
 *
 * This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
 * WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
 *
 *
 * This file is part of AMDiS
 *
 * See also license.opensource.txt in the distribution.
18
 *
19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38
 ******************************************************************************/



#include "PetscSolverCahnHilliard_DD.h"
#include "parallel/PetscHelper.h"
#include "TransformDOF.h"
#include "Expressions.h"

namespace AMDiS { namespace Parallel {


  using namespace std;


  PetscErrorCode pcChShell_dd_b(PC pc, Vec b, Vec x) // solve Px=b
  {FUNCNAME("PCApply()");
    void *ctx;
    PCShellGetContext(pc, &ctx);
    CahnHilliardData2* data = static_cast<CahnHilliardData2*>(ctx);
39

40 41 42
    Vec y1, y2;
    VecDuplicate(b, &y1);
    VecDuplicate(b, &y2);
43

44
    KSPSolve(data->kspMplusK, b, y1);
45
    MatMult(data->matMass, y1, y2);
46
    KSPSolve(data->kspMplusK, y2, x);
47

48 49 50 51 52 53 54 55 56
    PetscFunctionReturn(0);
  }

  /// solve Cahn-Hilliard Preconditioner
  PetscErrorCode pcChSchurShell_dd(PC pc, Vec b, Vec x) // solve Px=b
  {FUNCNAME("PCApply()");
    void *ctx;
    PCShellGetContext(pc, &ctx);
    CahnHilliardData2* data = static_cast<CahnHilliardData2*>(ctx);
57

58 59 60
    /// create S = M + eps^2*delta*K*D^(-1)*K where D=diag(M)
    Mat K, S;
    MatDuplicate(data->matMinusDeltaK, MAT_COPY_VALUES, &K);
61

62 63 64 65 66
    MatGetDiagonal(data->matMass, x);
    VecReciprocal(x);
    MatDiagonalScale(K, x, PETSC_NULL); 						// K' := -delta*D^(-1)*K
    MatMatMult(data->matMinusDeltaK, K, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &S); 		// S := -delta*K*K'
    MatAYPX(S, sqr(*data->eps)/(*data->delta), data->matMass, DIFFERENT_NONZERO_PATTERN); 	// S = eps^2/delta*S + M
67

68 69 70 71 72 73 74 75 76 77 78 79 80 81 82
    /// create new solver for S
    KSP kspS;
    KSPCreate(*data->mpiCommGlobal, &kspS);
#if (PETSC_VERSION_MINOR >= 5)
    KSPSetOperators(kspS, S, S);
#else
    KSPSetOperators(kspS, S, S, SAME_NONZERO_PATTERN);
#endif
    petsc_helper::setSolver(kspS, "S_", KSPFGMRES, PCSHELL, 1e-6, 1e-8, 5);
    {
      PC pc;
      KSPGetPC(kspS, &pc);
      PCShellSetApply(pc, pcChShell_dd_b);
      PCShellSetContext(pc, data);
    }
83

84
    KSPSolve(kspS, b, x); 				// S*x2 = x1
85

86 87 88
    MatDestroy(&S);
    MatDestroy(&K);
    KSPDestroy(&kspS);
89

90 91
    PetscFunctionReturn(0);
  }
92 93


94 95

  PetscSolverCahnHilliard_DD::PetscSolverCahnHilliard_DD(string name)
96 97
    : PetscSolverGlobalMatrix(name),
      useOldInitialGuess(false),
98 99
      laplaceSolutionMode(0),
      massMatrixSolver(NULL),
100
      laplaceMatrixSolver(NULL),
101 102 103 104 105 106 107
      deltaKMatrixSolver(NULL),
      eps(NULL),
      delta(NULL),
      tau(NULL),
      solution(NULL),
      phase(NULL)
  {
108 109
    Parameters::get(initFileStr + "->use old initial guess",
		    useOldInitialGuess);
110 111 112 113 114 115 116 117
  }


  void PetscSolverCahnHilliard_DD::initSolver(KSP &ksp)
  {
    FUNCNAME("PetscSolverCahnHilliard_DD::initSolver()");

    // Create FGMRES based outer solver
118

119
    MSG("CREATE POS 1: %p\n", &ksp);
120
    KSPCreate(domainComm, &ksp);
121 122 123 124 125 126 127 128 129 130
#if (PETSC_VERSION_MINOR >= 5)
    KSPSetOperators(ksp, getMatInterior(), getMatInterior());
#else
    KSPSetOperators(ksp, getMatInterior(), getMatInterior(), SAME_NONZERO_PATTERN);
#endif
    if (getInfo() >= 10)
      KSPMonitorSet(ksp, KSPMonitorDefault, PETSC_NULL, PETSC_NULL);
    else if (getInfo() >= 20)
      KSPMonitorSet(ksp, KSPMonitorTrueResidualNorm, PETSC_NULL, PETSC_NULL);
    petsc_helper::setSolver(ksp, "ch_", KSPFGMRES, PCNONE, getRelative(), getTolerance(), getMaxIterations());
131

132 133
    if (useOldInitialGuess)
      KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
134

135 136 137 138 139 140 141 142 143
  }


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

    MPI::COMM_WORLD.Barrier();
    double wtime = MPI::Wtime();
144 145


146 147 148 149
    if (tau) {
      delta = new double;
      *delta = std::sqrt(*tau);
    }
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
    vector<int> chPotentialComponent;
    chPotentialComponent.push_back(0);
    vector<int> chSchurComponent;
    chSchurComponent.push_back(1);

    PCSetType(pc, PCFIELDSPLIT);
    PCFieldSplitSetType(pc, PC_COMPOSITE_SCHUR);
    PCFieldSplitSetSchurFactType(pc, PC_FIELDSPLIT_SCHUR_FACT_FULL);

    createFieldSplit(pc, "ch_potential", chPotentialComponent);
    createFieldSplit(pc, "ch_schur", chSchurComponent);
    PCSetFromOptions(pc);

    KSPSetUp(kspInterior);

    KSP *subKspCH;
    int nSubKspCH;
    PCFieldSplitGetSubKSP(pc, &nSubKspCH, &subKspCH);


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

    KSP kspChPotential = subKspCH[0];
    KSP kspChSchur = subKspCH[1];
    PetscFree(subKspCH);

    KSPSetType(kspChSchur, KSPPREONLY);
    PC pcSub;
    KSPGetPC(kspChSchur, &pcSub);
    PCSetType(pcSub, PCSHELL);
    PCShellSetApply(pcSub, pcChSchurShell_dd);
    PCShellSetContext(pcSub, &matShellContext);

    const FiniteElemSpace *feSpace = componentSpaces[0];
186

187 188 189 190 191 192 193
    // === Mass matrix solver ===
    DOFMatrix massMatrix(feSpace, feSpace);
    Operator massOp(feSpace, feSpace);
    addZOT(massOp, valueOf(phase));
    massMatrix.assembleOperator(massOp);
    massMatrixSolver = createSubSolver(0, "mass_");
    massMatrixSolver->fillPetscMatrix(&massMatrix);
194

195 196 197 198 199 200 201 202
    // === Laplace matrix solver ===
    DOFMatrix laplaceMatrix(feSpace, feSpace);
    Operator laplaceOp(feSpace, feSpace);
    addSOT(laplaceOp, ((*eps)*std::sqrt(*delta)) * valueOf(phase));
    addSOT(laplaceOp, valueOf(phase));
    laplaceMatrix.assembleOperator(laplaceOp);
    laplaceMatrixSolver = createSubSolver(0, "MpK_");
    laplaceMatrixSolver->fillPetscMatrix(&laplaceMatrix);
203

204 205 206 207 208 209 210 211
    // === matrix (-delta*K) ===
    DOFMatrix deltaKMatrix(feSpace, feSpace);
    Operator laplaceOp2(feSpace, feSpace);
    addSOT(laplaceOp2, -(*delta) * valueOf(phase));
    deltaKMatrix.assembleOperator(laplaceOp2);
    deltaKMatrixSolver = createSubSolver(0, "laplace_");
    deltaKMatrixSolver->fillPetscMatrix(&deltaKMatrix);

212

213 214
    // === Setup solver ===
    matShellContext.kspMplusK = laplaceMatrixSolver->getSolver();
215
    matShellContext.matMinusDeltaK = deltaKMatrixSolver->getMatInterior();
216 217 218
    matShellContext.eps = eps;
    matShellContext.delta = delta;
    matShellContext.kspMass = massMatrixSolver->getSolver();
219
    matShellContext.matMass = massMatrixSolver->getMatInterior();
220 221 222 223 224
    matShellContext.mpiCommGlobal = &(meshDistributor->getMpiComm(0));

    petsc_helper::setSolver(matShellContext.kspMass, "mass_", KSPCG, PCJACOBI, 0.0, 1e-14, 5);
    petsc_helper::setSolver(matShellContext.kspMplusK, "MpK_", KSPRICHARDSON, PCHYPRE, 0.0, 1e-14, 1);
    petsc_helper::setSolver(kspChPotential, "chPotential",  KSPRICHARDSON, PCHYPRE, 0.0, 1e-14, 1);
225 226

    MSG("Setup of Cahn-Hilliard preconditioner needed %.5f seconds\n",
227 228 229 230 231 232 233 234 235 236 237
	MPI::Wtime() - wtime);
  }


  void PetscSolverCahnHilliard_DD::exitPreconditioner(PC pc)
  {
    FUNCNAME("PetscSolverCahnHilliard_DD::exitPreconditioner()");

    massMatrixSolver->destroyMatrixData();
    massMatrixSolver->destroyVectorData();
    laplaceMatrixSolver->destroyMatrixData();
238
    laplaceMatrixSolver->destroyVectorData();
239 240
    deltaKMatrixSolver->destroyMatrixData();
    deltaKMatrixSolver->destroyVectorData();
241

242 243 244 245

    massMatrixSolver->destroyVectorData();
    laplaceMatrixSolver->destroyVectorData();
    deltaKMatrixSolver->destroyVectorData();
246

247 248 249 250 251
    delete massMatrixSolver;
    massMatrixSolver = NULL;

    delete laplaceMatrixSolver;
    laplaceMatrixSolver = NULL;
252

253 254 255
    delete deltaKMatrixSolver;
    deltaKMatrixSolver = NULL;

256

257 258 259 260 261
    if (tau) {
      delete delta;
    }
  }
} }