PetscSolverFetiOperators.cc 17.4 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
//
// 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;
  }


37
38
39
40
41
42
43
  int petscMultMatSchurPrimalAugmented(Mat mat, Vec x, Vec y)
  {
    void *ctx;
    MatShellGetContext(mat, &ctx);
    SchurPrimalAugmentedData* data = 
      static_cast<SchurPrimalAugmentedData*>(ctx);

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
    Vec x_primal, x_mu, y_primal, y_mu;
    
    if (data->nestedVec) {
      VecNestGetSubVec(x, 0, &x_primal);
      VecNestGetSubVec(x, 1, &x_mu);
      VecNestGetSubVec(y, 0, &y_primal);
      VecNestGetSubVec(y, 1, &y_mu);
    } else {
      VecDuplicate(data->tmp_vec_primal, &x_primal);
      VecDuplicate(data->tmp_vec_primal, &y_primal);

      PetscInt allLocalSize, allSize;
      VecGetLocalSize(x, &allLocalSize);
      VecGetSize(x, &allSize);

      PetscInt primalLocalSize, primalSize;
      VecGetLocalSize(x_primal, &primalLocalSize);
      VecGetSize(x_primal, &primalSize);

      TEST_EXIT_DBG(allSize > primalSize)("Should not happen!\n");
      TEST_EXIT_DBG(allLocalSize >= primalLocalSize)("Should not happen!\n");

      PetscInt muLocalSize = allLocalSize - primalLocalSize;
      PetscInt muSize = allSize - primalSize;

      VecCreateMPI(PETSC_COMM_WORLD, muLocalSize, muSize, &x_mu);
      VecCreateMPI(PETSC_COMM_WORLD, muLocalSize, muSize, &y_mu);

      PetscScalar *allValue;
      PetscScalar *primalValue;
      PetscScalar *muValue;
      VecGetArray(x, &allValue);
      VecGetArray(x_primal, &primalValue);
      VecGetArray(x_mu, &muValue);

      for (int i = 0; i < primalLocalSize; i++)
	primalValue[i] = allValue[i];
      for (int i = 0; i < muLocalSize; i++)
	muValue[i] = allValue[primalLocalSize + i];

      VecRestoreArray(x, &allValue);
      VecRestoreArray(x_primal, &primalValue);
      VecRestoreArray(x_mu, &muValue);
    }
88
89
90
91
92

    // inv(K_BB) K_BPi x_Pi
    MatMult(data->subSolver->getMatInteriorCoarse(), x_primal, data->tmp_vec_b0);
    data->subSolver->solveGlobal(data->tmp_vec_b0, data->tmp_vec_b0);

93
94
    // inv(K_BB) trans(J) trans(Q) x_mu
    MatMultTranspose(*(data->mat_augmented_lagrange), x_mu, data->tmp_vec_lagrange);
95
96
97
98
99
100
101
102
103
    MatMultTranspose(*(data->mat_lagrange), data->tmp_vec_lagrange, data->tmp_vec_b1);
    data->subSolver->solveGlobal(data->tmp_vec_b1, data->tmp_vec_b1);

    // y_Pi = (K_PiPi - K_PiB inv(K_BB) K_BPi) x_pi
    MatMult(data->subSolver->getMatCoarseInterior(), data->tmp_vec_b0,
	    data->tmp_vec_primal);
    MatMult(data->subSolver->getMatCoarse(), x_primal, y_primal);
    VecAXPY(y_primal, -1.0, data->tmp_vec_primal);

104
    // y_Pi += (-K_PiB inv(K_BB) J^T Q^T) x_mu
105
106
107
108
    MatMult(data->subSolver->getMatCoarseInterior(), data->tmp_vec_b1,
	    data->tmp_vec_primal);
    VecAXPY(y_primal, -1.0, data->tmp_vec_primal);

109
110
    // y_mu = (-Q J inv(K_BB) K_BPi) x_Pi + (-Q J inv(K_BB) J^T Q) x_mu
    //      = -Q J (inv(K_BB) K_BPi x_Pi + inv(K_BB) J^T Q x_mu)
111
112
    VecAXPY(data->tmp_vec_b0, 1.0, data->tmp_vec_b1);
    MatMult(*(data->mat_lagrange), data->tmp_vec_b0, data->tmp_vec_lagrange);
113
    MatMult(*(data->mat_augmented_lagrange), data->tmp_vec_lagrange, y_mu);
114
115
    VecScale(y_mu, -1.0);

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
    if (!data->nestedVec) {
      PetscInt allLocalSize;
      VecGetLocalSize(x, &allLocalSize);
      PetscInt primalLocalSize;
      VecGetLocalSize(x_primal, &primalLocalSize);
      PetscInt muLocalSize = allLocalSize - primalLocalSize;

      PetscScalar *allValue;
      PetscScalar *primalValue;
      PetscScalar *muValue;
      VecGetArray(y, &allValue);
      VecGetArray(y_primal, &primalValue);
      VecGetArray(y_mu, &muValue);

      for (int i = 0; i < primalLocalSize; i++)
	allValue[i] = primalValue[i];

      for (int i = 0; i < muLocalSize; i++)
	allValue[primalLocalSize + i] = muValue[i];

      VecRestoreArray(y, &allValue);
      VecRestoreArray(y_primal, &primalValue);
      VecRestoreArray(y_mu, &muValue);

      VecDestroy(&x_primal);
      VecDestroy(&y_primal);
      VecDestroy(&x_mu);
      VecDestroy(&y_mu);
    }

146
147
148
149
    return 0;
  }


150
151
152
153
154
  // y = mat * x
  int petscMultMatFeti(Mat mat, Vec x, Vec y)
  {
    FUNCNAME("petscMultMatFeti()");

155
156
    //    F = J inv(K_BB) trans(J) + J inv(K_BB) K_BPi inv(S_PiPi) K_PiB inv(K_BB) trans(J)
    // => F = J [I + inv(K_BB) K_BPi inv(S_PiPi) K_PiB] inv(K_BB) trans(J)
157
158
159
160
161
162
163

    double wtime = MPI::Wtime();

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

Thomas Witkowski's avatar
Thomas Witkowski committed
164
    MatMultTranspose(*(data->mat_lagrange), x, data->tmp_vec_b0);
165
166

    double wtime01 = MPI::Wtime();
Thomas Witkowski's avatar
Thomas Witkowski committed
167
    data->subSolver->solveGlobal(data->tmp_vec_b0, data->tmp_vec_b0);
168
169
170

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

Thomas Witkowski's avatar
Thomas Witkowski committed
171
    MatMult(*(data->mat_lagrange), data->tmp_vec_b0, data->tmp_vec_lagrange);
172
173

    MatMult(data->subSolver->getMatCoarseInterior(), 
Thomas Witkowski's avatar
Thomas Witkowski committed
174
	    data->tmp_vec_b0, data->tmp_vec_primal0);
175
176

    wtime01 = MPI::Wtime();
Thomas Witkowski's avatar
Thomas Witkowski committed
177
    KSPSolve(*(data->ksp_schur_primal), data->tmp_vec_primal0, data->tmp_vec_primal0);
178
179
180
    FetiTimings::fetiSolve02 += (MPI::Wtime() - wtime01);

    MatMult(data->subSolver->getMatInteriorCoarse(), 
Thomas Witkowski's avatar
Thomas Witkowski committed
181
	    data->tmp_vec_primal0, data->tmp_vec_b0);
182
183

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

Thomas Witkowski's avatar
Thomas Witkowski committed
187
    MatMult(*(data->mat_lagrange), data->tmp_vec_b0, y);
188
189
190
191
192
193
194
195
196

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

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

    return 0;
  }


197
198
199
200
201
202
203
204
205
206
  // y = mat * x
  int petscMultMatFetiAugmented(Mat mat, Vec x, Vec y)
  {
    FUNCNAME("petscMultMatFetiAugmented()");

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

    Vec vec_mu0, vec_mu1;
207
    MatGetVecs(*(data->mat_augmented_lagrange), PETSC_NULL, &vec_mu0);
208
209
210
211
212
213
214
    VecDuplicate(vec_mu0, &vec_mu1);

    MatMultTranspose(*(data->mat_lagrange), x, data->tmp_vec_b0);
    data->subSolver->solveGlobal(data->tmp_vec_b0, data->tmp_vec_b0);

    MatMult(data->subSolver->getMatCoarseInterior(), data->tmp_vec_b0, data->tmp_vec_primal0);
    MatMult(*(data->mat_lagrange), data->tmp_vec_b0, data->tmp_vec_lagrange);
215
    MatMult(*(data->mat_augmented_lagrange), data->tmp_vec_lagrange, vec_mu0);
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234

    Vec vec_array0[2] = {data->tmp_vec_primal0, vec_mu0};
    Vec vec_array1[2] = {data->tmp_vec_primal1, vec_mu1};
    Vec vec_nest0, vec_nest1;
    VecCreateNest(PETSC_COMM_WORLD, 2, PETSC_NULL, vec_array0, &vec_nest0);
    VecCreateNest(PETSC_COMM_WORLD, 2, PETSC_NULL, vec_array1, &vec_nest1);

    KSPSolve(*(data->ksp_schur_primal), vec_nest0, vec_nest1);

    // Step 1
    MatMult(*(data->mat_lagrange), data->tmp_vec_b0, y);

    // Step 2
    MatMult(data->subSolver->getMatInteriorCoarse(), data->tmp_vec_primal1, data->tmp_vec_b0);
    data->subSolver->solveGlobal(data->tmp_vec_b0, data->tmp_vec_b0);
    MatMult(*(data->mat_lagrange), data->tmp_vec_b0, data->tmp_vec_lagrange);
    VecAXPY(y, 1.0, data->tmp_vec_lagrange);

    // Step 3
235
    MatMultTranspose(*(data->mat_augmented_lagrange), vec_mu1, data->tmp_vec_lagrange);
236
237
238
239
240
241
242
243
244
245
246
247
248
    MatMultTranspose(*(data->mat_lagrange), data->tmp_vec_lagrange, data->tmp_vec_b0);
    data->subSolver->solveGlobal(data->tmp_vec_b0, data->tmp_vec_b0);
    MatMult(*(data->mat_lagrange), data->tmp_vec_b0, data->tmp_vec_lagrange);
    VecAXPY(y, 1.0, data->tmp_vec_lagrange);

    VecDestroy(&vec_mu0);
    VecDestroy(&vec_mu1);
    VecDestroy(&vec_nest0);
    VecDestroy(&vec_nest1);
    return 0;
  }


249
250
251
252
253
254
  int petscMultMatFetiInterface(Mat mat, Vec x, Vec y)
  {
    FUNCNAME("petscMultMatFetiInterface()");

    double wtime = MPI::Wtime();

255
256
257
258
259
260
    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);

261
262
    void *ctx;
    MatShellGetContext(mat, &ctx);
263
    FetiData* data = static_cast<FetiData*>(ctx);
264
265


Thomas Witkowski's avatar
Thomas Witkowski committed
266
    // === Calculation of v_{B} ===
267

Thomas Witkowski's avatar
Thomas Witkowski committed
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
    // 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 ===
283

Thomas Witkowski's avatar
Thomas Witkowski committed
284
285
286
287
    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);
288

Thomas Witkowski's avatar
Thomas Witkowski committed
289
    // tmp_vec_primal1 = A_{\Pi B} A_{BB}^{-1} v_{B}
290
    MatMult(data->subSolver->getMatCoarseInterior(), 
Thomas Witkowski's avatar
Thomas Witkowski committed
291
292
293
294
	    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);
295
296

    wtime01 = MPI::Wtime();
Thomas Witkowski's avatar
Thomas Witkowski committed
297
298
    // 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);
299
300
    FetiTimings::fetiSolve02 += (MPI::Wtime() - wtime01);

Thomas Witkowski's avatar
Thomas Witkowski committed
301
    // tmp_vec_b1 = A_{B\Pi} S_{\Pi\Pi}^{-1} (v_{\Pi} - A_{\Pi B} A_{BB}^{-1} v_{B})
302
    MatMult(data->subSolver->getMatInteriorCoarse(), 
Thomas Witkowski's avatar
Thomas Witkowski committed
303
	    data->tmp_vec_primal0, data->tmp_vec_b1);
304
305

    wtime01 = MPI::Wtime();
Thomas Witkowski's avatar
Thomas Witkowski committed
306
307
    // 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);
308
309
    FetiTimings::fetiSolve01 += (MPI::Wtime() - wtime01);

Thomas Witkowski's avatar
Thomas Witkowski committed
310
311
312
313
314
315
316
317
318
319
320
321
322
323
    // 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);
324

Thomas Witkowski's avatar
Thomas Witkowski committed
325
326
    // y_lagrange = J tmp_vec_b0
    MatMult(*(data->mat_lagrange), data->tmp_vec_b0, y_lagrange);
327
328
329
330
331
332
333

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

    return 0;
  }


334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
  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;
  }


Thomas Witkowski's avatar
Thomas Witkowski committed
399
  PetscErrorCode petscApplyFetiLumpedPrecon(PC pc, Vec xvec, Vec yvec)
400
401
402
403
404
405
  {
    // Get data for the preconditioner
    void *ctx;
    PCShellGetContext(pc, &ctx);
    FetiLumpedPreconData* data = static_cast<FetiLumpedPreconData*>(ctx);

Thomas Witkowski's avatar
Thomas Witkowski committed
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
    // Multiply with scaled Lagrange constraint matrix.
    MatMultTranspose(*(data->mat_lagrange_scaled), xvec, data->tmp_vec_b0);


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

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

    PetscScalar *local_b, *local_duals;
    VecGetArray(data->tmp_vec_b0, &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_b0, &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_b0, &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_b0, &local_b);
    VecRestoreArray(data->tmp_vec_duals0, &local_duals);


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

    return 0;
  }


  PetscErrorCode petscApplyFetiInterfaceLumpedPrecon(PC pc, Vec xvec, Vec yvec)
  {
Thomas Witkowski's avatar
Blub    
Thomas Witkowski committed
456
    FUNCNAME("precon");
Thomas Witkowski's avatar
Thomas Witkowski committed
457
458
459
460
461
462
463
464
465
466
467
468
    // Get data for the preconditioner
    void *ctx;
    PCShellGetContext(pc, &ctx);
    FetiInterfaceLumpedPreconData* data = 
      static_cast<FetiInterfaceLumpedPreconData*>(ctx);
    
    Vec x_interface, x_lagrange, y_interface, y_lagrange;    
    VecNestGetSubVec(xvec, 0, &x_interface);
    VecNestGetSubVec(xvec, 1, &x_lagrange);
    VecNestGetSubVec(yvec, 0, &y_interface);
    VecNestGetSubVec(yvec, 1, &y_lagrange);

469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
    bool useMassMatrix = false;
    Parameters::get("lp->mass matrix", useMassMatrix);
    if (useMassMatrix) {
      KSPSolve(data->ksp_mass, x_interface, y_interface);
    } else {
      VecCopy(x_interface, y_interface);
      double scaleFactor = 1.0;
      Parameters::get("lp->scal", scaleFactor);
      bool mult = false;
      Parameters::get("lp->mult", mult);
      if (mult)
	VecPointwiseMult(y_interface, x_interface, data->h2vec);
      if (scaleFactor != 0.0)
	VecScale(y_interface, scaleFactor);
    }
484

Thomas Witkowski's avatar
Thomas Witkowski committed
485
    MatMultTranspose(*(data->mat_lagrange_scaled), x_lagrange, data->tmp_vec_b0);   
486
487

    // === Restriction of the B nodes to the boundary nodes. ===
Thomas Witkowski's avatar
Thomas Witkowski committed
488
    
489
490
    int nLocalB;
    int nLocalDuals;
Thomas Witkowski's avatar
Thomas Witkowski committed
491
    VecGetLocalSize(data->tmp_vec_b0, &nLocalB);
492
    VecGetLocalSize(data->tmp_vec_duals0, &nLocalDuals);
Thomas Witkowski's avatar
Thomas Witkowski committed
493
    
494
    PetscScalar *local_b, *local_duals;
Thomas Witkowski's avatar
Thomas Witkowski committed
495
    VecGetArray(data->tmp_vec_b0, &local_b);
496
    VecGetArray(data->tmp_vec_duals0, &local_duals);
Thomas Witkowski's avatar
Thomas Witkowski committed
497
    
498
499
500
    for (map<int, int>::iterator it = data->localToDualMap.begin();
	 it != data->localToDualMap.end(); ++it)
      local_duals[it->second] = local_b[it->first];
Thomas Witkowski's avatar
Thomas Witkowski committed
501
    
Thomas Witkowski's avatar
Thomas Witkowski committed
502
    VecRestoreArray(data->tmp_vec_b0, &local_b);
503
    VecRestoreArray(data->tmp_vec_duals0, &local_duals);
Thomas Witkowski's avatar
Thomas Witkowski committed
504
505
    
    
506
    // === K_DD ===
Thomas Witkowski's avatar
Thomas Witkowski committed
507
    
508
    MatMult(*(data->mat_duals_duals), data->tmp_vec_duals0, data->tmp_vec_duals1);
Thomas Witkowski's avatar
Thomas Witkowski committed
509
510
    
    
511
    // === Prolongation from local dual nodes to B nodes.
Thomas Witkowski's avatar
Thomas Witkowski committed
512
    
Thomas Witkowski's avatar
Thomas Witkowski committed
513
    VecGetArray(data->tmp_vec_b0, &local_b);
514
    VecGetArray(data->tmp_vec_duals1, &local_duals);
Thomas Witkowski's avatar
Thomas Witkowski committed
515
    
516
517
518
    for (map<int, int>::iterator it = data->localToDualMap.begin();
	 it != data->localToDualMap.end(); ++it)
      local_b[it->first] = local_duals[it->second];
Thomas Witkowski's avatar
Thomas Witkowski committed
519
    
Thomas Witkowski's avatar
Thomas Witkowski committed
520
    VecRestoreArray(data->tmp_vec_b0, &local_b);
521
    VecRestoreArray(data->tmp_vec_duals0, &local_duals);
Thomas Witkowski's avatar
Thomas Witkowski committed
522
        
523
    // Multiply with scaled Lagrange constraint matrix.
Thomas Witkowski's avatar
Thomas Witkowski committed
524
    MatMult(*(data->mat_lagrange_scaled), data->tmp_vec_b0, y_lagrange);
Thomas Witkowski's avatar
Thomas Witkowski committed
525
  
526
527
528
    return 0;
  }

Thomas Witkowski's avatar
Thomas Witkowski committed
529

530
}