PetscSolverFetiOperators.cc 9.67 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
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
//
// 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/PetscSolverFetiOperators.h"
#include "parallel/PetscSolverFetiStructs.h"
#include "parallel/PetscSolverFetiTimings.h"

namespace AMDiS {

  int petscMultMatSchurPrimal(Mat mat, Vec x, Vec y)
  {
    // S_PiPi = K_PiPi - K_PiB inv(K_BB) K_BPi

    void *ctx;
    MatShellGetContext(mat, &ctx);
    SchurPrimalData* data = static_cast<SchurPrimalData*>(ctx);

    MatMult(data->subSolver->getMatInteriorCoarse(), x, data->tmp_vec_b);
    data->subSolver->solveGlobal(data->tmp_vec_b, data->tmp_vec_b);
    MatMult(data->subSolver->getMatCoarseInterior(), data->tmp_vec_b, 
	    data->tmp_vec_primal);
    MatMult(data->subSolver->getMatCoarse(), x, y);
    VecAXPBY(y, -1.0, 1.0, data->tmp_vec_primal);

    return 0;
  }


  // y = mat * x
  int petscMultMatFeti(Mat mat, Vec x, Vec y)
  {
    FUNCNAME("petscMultMatFeti()");

    //    F = L inv(K_BB) trans(L) + L inv(K_BB) K_BPi inv(S_PiPi) K_PiB inv(K_BB) trans(L)
    // => F = L [I + inv(K_BB) K_BPi inv(S_PiPi) K_PiB] inv(K_BB) trans(L)

    double wtime = MPI::Wtime();

    void *ctx;
    MatShellGetContext(mat, &ctx);
    FetiData* data = static_cast<FetiData*>(ctx);

Thomas Witkowski's avatar
Thomas Witkowski committed
51
    MatMultTranspose(*(data->mat_lagrange), x, data->tmp_vec_b0);
52
53

    double wtime01 = MPI::Wtime();
Thomas Witkowski's avatar
Thomas Witkowski committed
54
    data->subSolver->solveGlobal(data->tmp_vec_b0, data->tmp_vec_b0);
55
56
57

    FetiTimings::fetiSolve01 += (MPI::Wtime() - wtime01);

Thomas Witkowski's avatar
Thomas Witkowski committed
58
    MatMult(*(data->mat_lagrange), data->tmp_vec_b0, data->tmp_vec_lagrange);
59
60

    MatMult(data->subSolver->getMatCoarseInterior(), 
Thomas Witkowski's avatar
Thomas Witkowski committed
61
	    data->tmp_vec_b0, data->tmp_vec_primal0);
62
63

    wtime01 = MPI::Wtime();
Thomas Witkowski's avatar
Thomas Witkowski committed
64
    KSPSolve(*(data->ksp_schur_primal), data->tmp_vec_primal0, data->tmp_vec_primal0);
65
66
67
    FetiTimings::fetiSolve02 += (MPI::Wtime() - wtime01);

    MatMult(data->subSolver->getMatInteriorCoarse(), 
Thomas Witkowski's avatar
Thomas Witkowski committed
68
	    data->tmp_vec_primal0, data->tmp_vec_b0);
69
70

    wtime01 = MPI::Wtime();
Thomas Witkowski's avatar
Thomas Witkowski committed
71
    data->subSolver->solveGlobal(data->tmp_vec_b0, data->tmp_vec_b0);
72
73
    FetiTimings::fetiSolve01 += (MPI::Wtime() - wtime01);

Thomas Witkowski's avatar
Thomas Witkowski committed
74
    MatMult(*(data->mat_lagrange), data->tmp_vec_b0, y);
75
76
77
78
79
80
81
82
83

    VecAXPBY(y, 1.0, 1.0, data->tmp_vec_lagrange);

    FetiTimings::fetiSolve += (MPI::Wtime() - wtime);

    return 0;
  }


84
85
86
87
88
89
  int petscMultMatFetiInterface(Mat mat, Vec x, Vec y)
  {
    FUNCNAME("petscMultMatFetiInterface()");

    double wtime = MPI::Wtime();

90
91
92
93
94
95
    Vec x_interface, x_lagrange, y_interface, y_lagrange;    
    VecNestGetSubVec(x, 0, &x_interface);
    VecNestGetSubVec(x, 1, &x_lagrange);
    VecNestGetSubVec(y, 0, &y_interface);
    VecNestGetSubVec(y, 1, &y_lagrange);

96
97
    void *ctx;
    MatShellGetContext(mat, &ctx);
98
    FetiData* data = static_cast<FetiData*>(ctx);
99
100


Thomas Witkowski's avatar
Thomas Witkowski committed
101
    // === Calculation of v_{B} ===
102

Thomas Witkowski's avatar
Thomas Witkowski committed
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
    // tmp_vec_b0 = J^{T} \lambda
    MatMultTranspose(*(data->mat_lagrange), x_lagrange, data->tmp_vec_b0);
    // tmp_vec_b1 = A_{B\Gamma} u_{\Gamma}
    MatMult(data->subSolver->getMatInteriorCoarse(1), x_interface, data->tmp_vec_b1);
    // tmp_vec_b0 = A_{B\Gamma} u_{\Gamma} + J^{T} \lambda
    VecAXPY(data->tmp_vec_b0, 1.0, data->tmp_vec_b1);


    // == Calculation of v_{\Pi}
    
    // tmp_vec_primal0 = A_{\Pi\Gamma} u_{\Gamma}
    MatMult(data->subSolver->getMatCoarse(0, 1), x_interface, data->tmp_vec_primal0);
    

    // === Calculate action of FETI-DP operator ===
118

Thomas Witkowski's avatar
Thomas Witkowski committed
119
120
121
122
    double wtime01 = MPI::Wtime();
    // tmp_vec_b0 = A_{BB}^{-1} v_{B}
    data->subSolver->solveGlobal(data->tmp_vec_b0, data->tmp_vec_b0);
    FetiTimings::fetiSolve01 += (MPI::Wtime() - wtime01);
123

Thomas Witkowski's avatar
Thomas Witkowski committed
124
    // tmp_vec_primal1 = A_{\Pi B} A_{BB}^{-1} v_{B}
125
    MatMult(data->subSolver->getMatCoarseInterior(), 
Thomas Witkowski's avatar
Thomas Witkowski committed
126
127
128
129
	    data->tmp_vec_b0, data->tmp_vec_primal1);

    // tmp_vec_primal0 = v_{\Pi} - A_{\Pi B} A_{BB}^{-1} v_{B}
    VecAXPY(data->tmp_vec_primal0, -1.0, data->tmp_vec_primal1);
130
131

    wtime01 = MPI::Wtime();
Thomas Witkowski's avatar
Thomas Witkowski committed
132
133
    // tmp_vec_primal0 = S_{\Pi\Pi}^{-1} (v_{\Pi} - A_{\Pi B} A_{BB}^{-1} v_{B})
    KSPSolve(*(data->ksp_schur_primal), data->tmp_vec_primal0, data->tmp_vec_primal0);
134
135
    FetiTimings::fetiSolve02 += (MPI::Wtime() - wtime01);

Thomas Witkowski's avatar
Thomas Witkowski committed
136
    // tmp_vec_b1 = A_{B\Pi} S_{\Pi\Pi}^{-1} (v_{\Pi} - A_{\Pi B} A_{BB}^{-1} v_{B})
137
    MatMult(data->subSolver->getMatInteriorCoarse(), 
Thomas Witkowski's avatar
Thomas Witkowski committed
138
	    data->tmp_vec_primal0, data->tmp_vec_b1);
139
140

    wtime01 = MPI::Wtime();
Thomas Witkowski's avatar
Thomas Witkowski committed
141
142
    // tmp_vec_b1 = A_{BB}^{-1} A_{B\Pi} S_{\Pi\Pi}^{-1} (v_{\Pi} - A_{\Pi B} A_{BB}^{-1} v_{B})
    data->subSolver->solveGlobal(data->tmp_vec_b1, data->tmp_vec_b1);
143
144
    FetiTimings::fetiSolve01 += (MPI::Wtime() - wtime01);

Thomas Witkowski's avatar
Thomas Witkowski committed
145
146
147
148
149
150
151
152
153
154
155
156
157
158
    // tmp_vec_b0 = A_{BB}^{-1} v_{B} - A_{BB}^{-1} A_{B\Pi} S_{\Pi\Pi}^{-1} (v_{\Pi} - A_{\Pi B} A_{BB}^{-1} v_{B})
    VecAXPY(data->tmp_vec_b0, -1.0, data->tmp_vec_b1);


    // === Calculate projection to interface and constraint variable ===

    // y_interface = A_{\Gamma B} tmp_vec_b0
    MatMult(data->subSolver->getMatCoarseInterior(1), data->tmp_vec_b0, y_interface);

    // tmp_vec_primal0 = S_{\Pi\Pi}^{-1} (v_{\Pi} - A_{\Pi B} A_{BB}^{-1} v_{B})
    // tmp_vec_interface = A_{\Gamma \Pi} tmp_vec_primal0
    MatMult(data->subSolver->getMatCoarse(1, 0), data->tmp_vec_primal0, data->tmp_vec_interface);
    // y_interface = A_{\Gamma B} tmp_vec_b0 + A_{\Gamma \Pi} tmp_vec_primal0
    VecAXPY(y_interface, 1.0, data->tmp_vec_interface);
159

Thomas Witkowski's avatar
Thomas Witkowski committed
160
161
    // y_lagrange = J tmp_vec_b0
    MatMult(*(data->mat_lagrange), data->tmp_vec_b0, y_lagrange);
162
163
164
165
166
167
168

    FetiTimings::fetiSolve += (MPI::Wtime() - wtime);

    return 0;
  }


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
238
239
240
  PetscErrorCode petscApplyFetiDirichletPrecon(PC pc, Vec x, Vec y)
  {
    double wtime = MPI::Wtime();

    // Get data for the preconditioner
    void *ctx;
    PCShellGetContext(pc, &ctx);
    FetiDirichletPreconData* data = static_cast<FetiDirichletPreconData*>(ctx);

    // Multiply with scaled Lagrange constraint matrix.
    MatMultTranspose(*(data->mat_lagrange_scaled), x, data->tmp_vec_b);


    // === Restriction of the B nodes to the boundary nodes. ===

    int nLocalB;
    int nLocalDuals;
    VecGetLocalSize(data->tmp_vec_b, &nLocalB);
    VecGetLocalSize(data->tmp_vec_duals0, &nLocalDuals);

    PetscScalar *local_b, *local_duals;
    VecGetArray(data->tmp_vec_b, &local_b);
    VecGetArray(data->tmp_vec_duals0, &local_duals);

    for (map<int, int>::iterator it = data->localToDualMap.begin();
	 it != data->localToDualMap.end(); ++it)
      local_duals[it->second] = local_b[it->first];

    VecRestoreArray(data->tmp_vec_b, &local_b);
    VecRestoreArray(data->tmp_vec_duals0, &local_duals);


    // === K_DD - K_DI inv(K_II) K_ID ===

    MatMult(*(data->mat_duals_duals), data->tmp_vec_duals0, data->tmp_vec_duals1);

    MatMult(*(data->mat_interior_duals), data->tmp_vec_duals0, data->tmp_vec_interior);
    KSPSolve(*(data->ksp_interior), data->tmp_vec_interior, data->tmp_vec_interior);
    MatMult(*(data->mat_duals_interior), data->tmp_vec_interior, data->tmp_vec_duals0);

    VecAXPBY(data->tmp_vec_duals0, 1.0, -1.0, data->tmp_vec_duals1);


    // === Prolongation from local dual nodes to B nodes.

    VecGetArray(data->tmp_vec_b, &local_b);
    VecGetArray(data->tmp_vec_duals0, &local_duals);

    for (map<int, int>::iterator it = data->localToDualMap.begin();
	 it != data->localToDualMap.end(); ++it)
      local_b[it->first] = local_duals[it->second];

    VecRestoreArray(data->tmp_vec_b, &local_b);
    VecRestoreArray(data->tmp_vec_duals0, &local_duals);


    // Multiply with scaled Lagrange constraint matrix.
    MatMult(*(data->mat_lagrange_scaled), data->tmp_vec_b, y);

    FetiTimings::fetiPreconditioner += (MPI::Wtime() - wtime);

    return 0;
  }


  PetscErrorCode petscApplyFetiLumpedPrecon(PC pc, Vec x, Vec y)
  {
    // Get data for the preconditioner
    void *ctx;
    PCShellGetContext(pc, &ctx);
    FetiLumpedPreconData* data = static_cast<FetiLumpedPreconData*>(ctx);

241
242
243
244
245
246
247
248
249
250
251
252
    Vec xvec, yvec;
    if (data->enableStokesMode) {
      Vec x_interface, x_lagrange, y_interface, y_lagrange;    
      VecNestGetSubVec(x, 0, &x_interface);
      VecNestGetSubVec(x, 1, &x_lagrange);
      VecNestGetSubVec(y, 0, &y_interface);
      VecNestGetSubVec(y, 1, &y_lagrange);
      
      xvec = x_lagrange;
      yvec = y_lagrange;

      VecCopy(x_interface, y_interface);
Thomas Witkowski's avatar
So    
Thomas Witkowski committed
253
      //      VecScale(y_interface, 0.2519);
254
255
256
257
258
    } else {
      xvec = x;
      yvec = y;
    }

259
    // Multiply with scaled Lagrange constraint matrix.
260
    MatMultTranspose(*(data->mat_lagrange_scaled), xvec, data->tmp_vec_b);
261
262
263
264
265
266
267
268
269
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


    // === Restriction of the B nodes to the boundary nodes. ===

    int nLocalB;
    int nLocalDuals;
    VecGetLocalSize(data->tmp_vec_b, &nLocalB);
    VecGetLocalSize(data->tmp_vec_duals0, &nLocalDuals);

    PetscScalar *local_b, *local_duals;
    VecGetArray(data->tmp_vec_b, &local_b);
    VecGetArray(data->tmp_vec_duals0, &local_duals);

    for (map<int, int>::iterator it = data->localToDualMap.begin();
	 it != data->localToDualMap.end(); ++it)
      local_duals[it->second] = local_b[it->first];

    VecRestoreArray(data->tmp_vec_b, &local_b);
    VecRestoreArray(data->tmp_vec_duals0, &local_duals);


    // === K_DD ===

    MatMult(*(data->mat_duals_duals), data->tmp_vec_duals0, data->tmp_vec_duals1);


    // === Prolongation from local dual nodes to B nodes.

    VecGetArray(data->tmp_vec_b, &local_b);
    VecGetArray(data->tmp_vec_duals1, &local_duals);

    for (map<int, int>::iterator it = data->localToDualMap.begin();
	 it != data->localToDualMap.end(); ++it)
      local_b[it->first] = local_duals[it->second];

    VecRestoreArray(data->tmp_vec_b, &local_b);
    VecRestoreArray(data->tmp_vec_duals0, &local_duals);


    // Multiply with scaled Lagrange constraint matrix.
301
    MatMult(*(data->mat_lagrange_scaled), data->tmp_vec_b, yvec);
302
303
304
305
306

    return 0;
  }

}