PetscSolverPfc.cc 7.11 KB
Newer Older
Praetorius, Simon's avatar
Praetorius, Simon committed
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
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
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
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
186
187
188
189
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
/******************************************************************************
 *
 * Extension of AMDiS - Adaptive multidimensional simulations
 *
 * Copyright (C) 2013 Dresden University of Technology. All Rights Reserved.
 * Web: https://fusionforge.zih.tu-dresden.de/projects/amdis
 *
 * Authors: 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.
 *
 *
 * See also license.opensource.txt in the distribution.
 * 
 ******************************************************************************/

#include "PetscSolverPfc.h"
#include "parallel/PetscHelper.h"
#include "parallel/PetscSolverGlobalMatrix.h"

namespace AMDiS { namespace Parallel {

  using namespace std;
  PetscErrorCode pcPfcMatMultSchur(Mat mat, Vec b, Vec x)
  {
    void *ctx;
    MatShellGetContext(mat, &ctx);
    PfcData* data = static_cast<PfcData*>(ctx);

    Vec y1, y2, y3;
    VecDuplicate(b, &y1);
    VecDuplicate(b, &y2);
    VecDuplicate(b, &y3);
    
    MatMult(data->matK, b, y1);		// y1 := K*b
    KSPSolve(data->kspM, y1, y2);	// M*y2 = y1
    MatMult(data->matK, y2, y3);	// y3 := K*y2
//     KSPSolve(data->kspM, y3, y2);	// M*y2 = y3
//     MatMult(data->matK, y2, y3);	// y3 := K*y2
    
    VecScale(y3, data->delta);		// y3 = delta*y3
    MatMultAdd(data->matM, b, y3, x);	// x = M*b + y3
    VecAXPY(x, -2.0*(data->delta), y1); // x := M*b - 2*delta*K*b + delta*B*b
    
    VecDestroy(&y1);
    VecDestroy(&y2);
    VecDestroy(&y3);
    
    PetscFunctionReturn(0);
  }
  
  /// solve Pfc Preconditioner
  PetscErrorCode pcPfcShell(PC pc, Vec b, Vec x) // solve Px=b
  { FUNCNAME("pcPfcShell()");
  
    void *ctx;
    PCShellGetContext(pc, &ctx);
    PfcData* data = static_cast<PfcData*>(ctx);
    
    Vec b1, b2, b3, x1, x2, x3;
    VecNestGetSubVec(b, 0, &b1);
    VecNestGetSubVec(b, 1, &b2);
    VecNestGetSubVec(b, 2, &b3);

    VecNestGetSubVec(x, 0, &x1);
    VecNestGetSubVec(x, 1, &x2);
    VecNestGetSubVec(x, 2, &x3);


    // Hilfsvariablen
    Vec y1, y2, tmp;
    VecDuplicate(b1, &y1);
    VecDuplicate(b1, &y2);
    VecDuplicate(b1, &tmp);

    KSPSolve(data->kspM, b1, y1); 			// M*y1 = b1    
    MatMult(data->matK, y1, tmp); 			// tmp := K*y1
    VecAYPX(tmp, -(data->tau), b2);			// tmp := b2 - tau*tmp
    
    KSPSolve(data->kspMpK, tmp, y2); 			// (M+sqrt(tau)K)*y2 = tmp
    MatMult(data->matM, y2, tmp);			// tmp := M*y2

    KSPSolve(data->kspSchur, tmp, x2); 			// S*x2 = tmp
    VecCopy(x2, x1); 					// x1 := x2
    VecAXPBYPCZ(x1, 1.0, 1.0/(data->delta), -1.0/(data->delta), y1, y2);	// x1 = 1*y1 + factor*y2 - factor*x1

    MatMult(data->matK, x2, tmp); 			// tmp := K*x2
    VecAYPX(tmp, -1.0, b3);				// tmp := b3 - tmp
    KSPSolve(data->kspM2, tmp, x3); 			// M*x3 = tmp    

    VecDestroy(&y1);
    VecDestroy(&y2);
    VecDestroy(&tmp);
    
    PetscFunctionReturn(0);
  }
  
  void PetscSolverPfc::initSolver(KSP &ksp)
  {
    // Create FGMRES based outer solver
    KSPCreate(meshDistributor->getMpiComm(0), &ksp);
    KSPSetOperators(ksp, getMatInterior(), getMatInterior(), SAME_NONZERO_PATTERN);
    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, "pfc_", KSPFGMRES, PCNONE, getRelative(), getTolerance(), getMaxIterations());
    KSPSetFromOptions(ksp);
    
    
    if (useOldInitialGuess)
      KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
  }
  

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

    TEST_EXIT(tau)("tau pointer not set!\n");
    
    PCSetType(getPc(), PCSHELL);
    PCShellSetApply(getPc(), pcPfcShell);
    PCShellSetContext(getPc(), &data);
    
    const FiniteElemSpace *feSpace = componentSpaces[0];
    
    // create mass-matrix
    DOFMatrix matrixM(feSpace, feSpace);
    Operator massOp(feSpace, feSpace);
    Simple_ZOT zot;
    massOp.addTerm(&zot);
    matrixM.assembleOperator(massOp);
    
    solverM = createSubSolver(0, "M_");
    solverM->fillPetscMatrix(&matrixM);
    data.matM = solverM->getMatInterior();
    data.kspM = solverM->getSolver();
    
    solverM2 = createSubSolver(0, "M2_");
    solverM2->fillPetscMatrix(&matrixM);
    data.kspM2 = solverM2->getSolver();
    
    // create laplace-matrix
    DOFMatrix matrixMpK(feSpace, feSpace);
    Operator laplaceOp2(feSpace, feSpace);
    Simple_SOT sot2(sqrt(*tau));
    laplaceOp2.addTerm(&zot);
    laplaceOp2.addTerm(&sot2);
    matrixMpK.assembleOperator(laplaceOp2);
    
    solverMpK = createSubSolver(0, "MpK_");
    solverMpK->fillPetscMatrix(&matrixMpK);
    matMpK = solverMpK->getMatInterior();
    data.kspMpK = solverMpK->getSolver();
    
    // create laplace-matrix
    DOFMatrix matrixK(feSpace, feSpace);
    Operator laplaceOp(feSpace, feSpace);
    Simple_SOT sot;
    laplaceOp.addTerm(&sot);
    matrixK.assembleOperator(laplaceOp);
    
    solverK = createSubSolver(0, "K_");
    solverK->fillPetscMatrix(&matrixK);
    data.matK = solverK->getMatInterior();
    
    // === Setup preconditioner data ===    
    data.delta = sqrt(*tau);
    data.tau = (*tau);

    petsc_helper::setSolver(data.kspM, "M_", KSPCG, PCJACOBI, 0.0, 1e-14, 2);
    petsc_helper::setSolver(data.kspM2, "M2_", KSPCG, PCBJACOBI, 0.0, 1e-14, 5);
    petsc_helper::setSolver(data.kspMpK, "MpK_", KSPRICHARDSON, PCHYPRE, 0.0, 1e-14, 1);
    
    PetscInt mLocal, mGlobal;
    MatGetLocalSize(data.matM, &mLocal, PETSC_NULL);
    MatGetSize(data.matM, &mGlobal, PETSC_NULL);
    MatCreateShell(meshDistributor->getMpiComm(0),
		    mLocal, mLocal, 
		    mGlobal, mGlobal,
		    &data, 
		    &matSchur);
    MatShellSetOperation(matSchur, MATOP_MULT, (void(*)(void))pcPfcMatMultSchur);

    KSPCreate(meshDistributor->getMpiComm(0), &data.kspSchur);
    KSPSetOperators(data.kspSchur, matSchur, matMpK, SAME_NONZERO_PATTERN);
    petsc_helper::setSolver(data.kspSchur, "schur_", KSPFGMRES, PCJACOBI, 1e-6, 1e-8, 5);
    KSPSetFromOptions(data.kspSchur);
  }


  void PetscSolverPfc::exitPreconditioner(PC pc)
  {
    FUNCNAME("PetscSolverPfc::exit()");
    
    MatDestroy(&matMpK);
    MatDestroy(&matSchur);
    
    solverM->destroyMatrixData();
    solverK->destroyMatrixData();
    solverMpK->destroyMatrixData();
    
    solverM->destroyVectorData();
    solverK->destroyVectorData();
    solverMpK->destroyVectorData();
    
    delete solverM;
    solverM = NULL;
    delete solverK;
    solverK = NULL;
    delete solverMpK;
    solverMpK = NULL;
  }
  
  
  PetscSolver* PetscSolverPfc::createSubSolver(int component,
					      string kspPrefix)
  {
    FUNCNAME("PetscSolverPfc::createSubSolver()");
    
    vector<const FiniteElemSpace*> fe;
    fe.push_back(componentSpaces[component]);
    
    PetscSolver* subSolver = new PetscSolverGlobalMatrix("");
    subSolver->setKspPrefix(kspPrefix.c_str());
    subSolver->setMeshDistributor(meshDistributor, 0);
    subSolver->init(fe, fe);
    
    ParallelDofMapping &subDofMap = subSolver->getDofMap();
    subDofMap[0] = dofMap[component];
    subDofMap.update();
    
    return subSolver;
  }
} }