Actual source code: qeplin.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*

  3:    Various types of linearization for quadratic eigenvalue problem.

  5:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  6:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  7:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  9:    This file is part of SLEPc.

 11:    SLEPc is free software: you can redistribute it and/or modify it under  the
 12:    terms of version 3 of the GNU Lesser General Public License as published by
 13:    the Free Software Foundation.

 15:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 16:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 17:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 18:    more details.

 20:    You  should have received a copy of the GNU Lesser General  Public  License
 21:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 22:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 23: */

 25: #include <slepc-private/pepimpl.h>
 26:  #include linearp.h

 28: /*
 29:     Given the quadratic problem (l^2*M + l*C + K)*x = 0 the linearization is
 30:     A*z = l*B*z for z = [  x  ] and A,B defined as follows:
 31:                         [ l*x ]

 33:             N1:
 34:                      A = [  0   I ]     B = [ I  0 ]
 35:                          [ -K  -C ]         [ 0  M ]

 37:             N2:
 38:                      A = [ -K   0 ]     B = [ C  M ]
 39:                          [  0   I ]         [ I  0 ]

 41:             S1:
 42:                      A = [  0  -K ]     B = [-K  0 ]
 43:                          [ -K  -C ]         [ 0  M ]

 45:             S2:
 46:                      A = [ -K   0 ]     B = [ C  M ]
 47:                          [  0   M ]         [ M  0 ]

 49:             H1:
 50:                      A = [  K   0 ]     B = [ 0  K ]
 51:                          [  C   K ]         [-M  0 ]

 53:             H2:
 54:                      A = [  0  -K ]     B = [ M  C ]
 55:                          [  M   0 ]         [ 0  M ]
 56:  */

 58: /* - - - N1 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 62: PetscErrorCode MatMult_Linear_N1A(Mat A,Vec x,Vec y)
 63: {
 64:   PetscErrorCode    ierr;
 65:   PEP_LINEAR        *ctx;
 66:   const PetscScalar *px;
 67:   PetscScalar       *py;
 68:   PetscInt          m;

 71:   MatShellGetContext(A,(void**)&ctx);
 72:   MatGetLocalSize(ctx->M,&m,NULL);
 73:   VecGetArrayRead(x,&px);
 74:   VecGetArray(y,&py);
 75:   VecPlaceArray(ctx->x1,px);
 76:   VecPlaceArray(ctx->x2,px+m);
 77:   VecPlaceArray(ctx->y1,py);
 78:   VecPlaceArray(ctx->y2,py+m);
 79:   /* y2 = -(K*x1 + C*x2) */
 80:   MatMult(ctx->K,ctx->x1,ctx->y2);
 81:   MatMult(ctx->C,ctx->x2,ctx->y1);
 82:   VecAXPY(ctx->y2,ctx->sfactor,ctx->y1);
 83:   VecScale(ctx->y2,-1.0);
 84:   /* y1 = x2 */
 85:   VecCopy(ctx->x2,ctx->y1);
 86:   VecResetArray(ctx->x1);
 87:   VecResetArray(ctx->x2);
 88:   VecResetArray(ctx->y1);
 89:   VecResetArray(ctx->y2);
 90:   VecRestoreArrayRead(x,&px);
 91:   VecRestoreArray(y,&py);
 92:   return(0);
 93: }

 97: PetscErrorCode MatMult_Linear_N1B(Mat B,Vec x,Vec y)
 98: {
 99:   PetscErrorCode    ierr;
100:   PEP_LINEAR        *ctx;
101:   const PetscScalar *px;
102:   PetscScalar       *py;
103:   PetscInt          m;

106:   MatShellGetContext(B,(void**)&ctx);
107:   MatGetLocalSize(ctx->M,&m,NULL);
108:   VecGetArrayRead(x,&px);
109:   VecGetArray(y,&py);
110:   VecPlaceArray(ctx->x1,px);
111:   VecPlaceArray(ctx->x2,px+m);
112:   VecPlaceArray(ctx->y1,py);
113:   VecPlaceArray(ctx->y2,py+m);
114:   /* y1 = x1 */
115:   VecCopy(ctx->x1,ctx->y1);
116:   /* y2 = M*x2 */
117:   MatMult(ctx->M,ctx->x2,ctx->y2);
118:   VecScale(ctx->y2,ctx->sfactor*ctx->sfactor);
119:   VecResetArray(ctx->x1);
120:   VecResetArray(ctx->x2);
121:   VecResetArray(ctx->y1);
122:   VecResetArray(ctx->y2);
123:   VecRestoreArrayRead(x,&px);
124:   VecRestoreArray(y,&py);
125:   return(0);
126: }

130: PetscErrorCode MatGetDiagonal_Linear_N1A(Mat A,Vec diag)
131: {
133:   PEP_LINEAR     *ctx;
134:   PetscScalar    *pd;
135:   PetscInt       m;

138:   MatShellGetContext(A,(void**)&ctx);
139:   MatGetLocalSize(ctx->M,&m,NULL);
140:   VecGetArray(diag,&pd);
141:   VecPlaceArray(ctx->x1,pd);
142:   VecPlaceArray(ctx->x2,pd+m);
143:   VecSet(ctx->x1,0.0);
144:   MatGetDiagonal(ctx->C,ctx->x2);
145:   VecScale(ctx->x2,-ctx->sfactor);
146:   VecResetArray(ctx->x1);
147:   VecResetArray(ctx->x2);
148:   VecRestoreArray(diag,&pd);
149:   return(0);
150: }

154: PetscErrorCode MatGetDiagonal_Linear_N1B(Mat B,Vec diag)
155: {
157:   PEP_LINEAR     *ctx;
158:   PetscScalar    *pd;
159:   PetscInt       m;

162:   MatShellGetContext(B,(void**)&ctx);
163:   MatGetLocalSize(ctx->M,&m,NULL);
164:   VecGetArray(diag,&pd);
165:   VecPlaceArray(ctx->x1,pd);
166:   VecPlaceArray(ctx->x2,pd+m);
167:   VecSet(ctx->x1,1.0);
168:   MatGetDiagonal(ctx->M,ctx->x2);
169:   VecScale(ctx->x2,ctx->sfactor*ctx->sfactor);
170:   VecResetArray(ctx->x1);
171:   VecResetArray(ctx->x2);
172:   VecRestoreArray(diag,&pd);
173:   return(0);
174: }

178: PetscErrorCode MatCreateExplicit_Linear_N1A(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
179: {
181:   PetscInt       M,N,m,n;
182:   Mat            Id;

185:   MatGetSize(ctx->M,&M,&N);
186:   MatGetLocalSize(ctx->M,&m,&n);
187:   MatCreate(PetscObjectComm((PetscObject)ctx->M),&Id);
188:   MatSetSizes(Id,m,n,M,N);
189:   MatSetFromOptions(Id);
190:   MatSetUp(Id);
191:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
192:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
193:   MatShift(Id,1.0);
194:   SlepcMatTile(0.0,Id,1.0,Id,-1.0,ctx->K,-ctx->sfactor,ctx->C,A);
195:   MatDestroy(&Id);
196:   return(0);
197: }

201: PetscErrorCode MatCreateExplicit_Linear_N1B(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
202: {
204:   PetscInt       M,N,m,n;
205:   Mat            Id;

208:   MatGetSize(ctx->M,&M,&N);
209:   MatGetLocalSize(ctx->M,&m,&n);
210:   MatCreate(PetscObjectComm((PetscObject)ctx->M),&Id);
211:   MatSetSizes(Id,m,n,M,N);
212:   MatSetFromOptions(Id);
213:   MatSetUp(Id);
214:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
215:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
216:   MatShift(Id,1.0);
217:   SlepcMatTile(1.0,Id,0.0,Id,0.0,Id,ctx->sfactor*ctx->sfactor,ctx->M,B);
218:   MatDestroy(&Id);
219:   return(0);
220: }

222: /* - - - N2 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

226: PetscErrorCode MatMult_Linear_N2A(Mat A,Vec x,Vec y)
227: {
228:   PetscErrorCode    ierr;
229:   PEP_LINEAR        *ctx;
230:   const PetscScalar *px;
231:   PetscScalar       *py;
232:   PetscInt          m;

235:   MatShellGetContext(A,(void**)&ctx);
236:   MatGetLocalSize(ctx->M,&m,NULL);
237:   VecGetArrayRead(x,&px);
238:   VecGetArray(y,&py);
239:   VecPlaceArray(ctx->x1,px);
240:   VecPlaceArray(ctx->x2,px+m);
241:   VecPlaceArray(ctx->y1,py);
242:   VecPlaceArray(ctx->y2,py+m);
243:   /* y1 = -K*x1 */
244:   MatMult(ctx->K,ctx->x1,ctx->y1);
245:   VecScale(ctx->y1,-1.0);
246:   /* y2 = x2 */
247:   VecCopy(ctx->x2,ctx->y2);
248:   VecResetArray(ctx->x1);
249:   VecResetArray(ctx->x2);
250:   VecResetArray(ctx->y1);
251:   VecResetArray(ctx->y2);
252:   VecRestoreArrayRead(x,&px);
253:   VecRestoreArray(y,&py);
254:   return(0);
255: }

259: PetscErrorCode MatMult_Linear_N2B(Mat B,Vec x,Vec y)
260: {
261:   PetscErrorCode    ierr;
262:   PEP_LINEAR        *ctx;
263:   const PetscScalar *px;
264:   PetscScalar       *py;
265:   PetscInt          m;

268:   MatShellGetContext(B,(void**)&ctx);
269:   MatGetLocalSize(ctx->M,&m,NULL);
270:   VecGetArrayRead(x,&px);
271:   VecGetArray(y,&py);
272:   VecPlaceArray(ctx->x1,px);
273:   VecPlaceArray(ctx->x2,px+m);
274:   VecPlaceArray(ctx->y1,py);
275:   VecPlaceArray(ctx->y2,py+m);
276:   /* y1 = C*x1 + M*x2 */
277:   MatMult(ctx->C,ctx->x1,ctx->y1);
278:   VecScale(ctx->y1,ctx->sfactor);
279:   MatMult(ctx->M,ctx->x2,ctx->y2);
280:   VecAXPY(ctx->y1,ctx->sfactor*ctx->sfactor,ctx->y2);
281:   /* y2 = x1 */
282:   VecCopy(ctx->x1,ctx->y2);
283:   VecResetArray(ctx->x1);
284:   VecResetArray(ctx->x2);
285:   VecResetArray(ctx->y1);
286:   VecResetArray(ctx->y2);
287:   VecRestoreArrayRead(x,&px);
288:   VecRestoreArray(y,&py);
289:   return(0);
290: }

294: PetscErrorCode MatGetDiagonal_Linear_N2A(Mat A,Vec diag)
295: {
297:   PEP_LINEAR     *ctx;
298:   PetscScalar    *pd;
299:   PetscInt       m;

302:   MatShellGetContext(A,(void**)&ctx);
303:   MatGetLocalSize(ctx->M,&m,NULL);
304:   VecGetArray(diag,&pd);
305:   VecPlaceArray(ctx->x1,pd);
306:   VecPlaceArray(ctx->x2,pd+m);
307:   MatGetDiagonal(ctx->K,ctx->x1);
308:   VecScale(ctx->x1,-1.0);
309:   VecSet(ctx->x2,1.0);
310:   VecResetArray(ctx->x1);
311:   VecResetArray(ctx->x2);
312:   VecRestoreArray(diag,&pd);
313:   return(0);
314: }

318: PetscErrorCode MatGetDiagonal_Linear_N2B(Mat B,Vec diag)
319: {
321:   PEP_LINEAR     *ctx;
322:   PetscScalar    *pd;
323:   PetscInt       m;

326:   MatShellGetContext(B,(void**)&ctx);
327:   MatGetLocalSize(ctx->M,&m,NULL);
328:   VecGetArray(diag,&pd);
329:   VecPlaceArray(ctx->x1,pd);
330:   VecPlaceArray(ctx->x2,pd+m);
331:   MatGetDiagonal(ctx->C,ctx->x1);
332:   VecScale(ctx->x1,ctx->sfactor);
333:   VecSet(ctx->x2,0.0);
334:   VecResetArray(ctx->x1);
335:   VecResetArray(ctx->x2);
336:   VecRestoreArray(diag,&pd);
337:   return(0);
338: }

342: PetscErrorCode MatCreateExplicit_Linear_N2A(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
343: {
345:   PetscInt       M,N,m,n;
346:   Mat            Id;

349:   MatGetSize(ctx->M,&M,&N);
350:   MatGetLocalSize(ctx->M,&m,&n);
351:   MatCreate(PetscObjectComm((PetscObject)ctx->M),&Id);
352:   MatSetSizes(Id,m,n,M,N);
353:   MatSetFromOptions(Id);
354:   MatSetUp(Id);
355:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
356:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
357:   MatShift(Id,1.0);
358:   SlepcMatTile(-1.0,ctx->K,0.0,Id,0.0,Id,1.0,Id,A);
359:   MatDestroy(&Id);
360:   return(0);
361: }

365: PetscErrorCode MatCreateExplicit_Linear_N2B(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
366: {
368:   PetscInt       M,N,m,n;
369:   Mat            Id;

372:   MatGetSize(ctx->M,&M,&N);
373:   MatGetLocalSize(ctx->M,&m,&n);
374:   MatCreate(PetscObjectComm((PetscObject)ctx->M),&Id);
375:   MatSetSizes(Id,m,n,M,N);
376:   MatSetFromOptions(Id);
377:   MatSetUp(Id);
378:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
379:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
380:   MatShift(Id,1.0);
381:   SlepcMatTile(ctx->sfactor,ctx->C,ctx->sfactor*ctx->sfactor,ctx->M,1.0,Id,0.0,Id,B);
382:   MatDestroy(&Id);
383:   return(0);
384: }

386: /* - - - S1 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

390: PetscErrorCode MatMult_Linear_S1A(Mat A,Vec x,Vec y)
391: {
392:   PetscErrorCode    ierr;
393:   PEP_LINEAR        *ctx;
394:   const PetscScalar *px;
395:   PetscScalar       *py;
396:   PetscInt          m;

399:   MatShellGetContext(A,(void**)&ctx);
400:   MatGetLocalSize(ctx->M,&m,NULL);
401:   VecGetArrayRead(x,&px);
402:   VecGetArray(y,&py);
403:   VecPlaceArray(ctx->x1,px);
404:   VecPlaceArray(ctx->x2,px+m);
405:   VecPlaceArray(ctx->y1,py);
406:   VecPlaceArray(ctx->y2,py+m);
407:   /* y2 = -(K*x1 + C*x2) */
408:   MatMult(ctx->K,ctx->x1,ctx->y2);
409:   VecScale(ctx->y2,-1.0);
410:   MatMult(ctx->C,ctx->x2,ctx->y1);
411:   VecAXPY(ctx->y2,-ctx->sfactor,ctx->y1);
412:   /* y1 = -K*x2 */
413:   MatMult(ctx->K,ctx->x2,ctx->y1);
414:   VecScale(ctx->y1,-1.0);
415:   VecResetArray(ctx->x1);
416:   VecResetArray(ctx->x2);
417:   VecResetArray(ctx->y1);
418:   VecResetArray(ctx->y2);
419:   VecRestoreArrayRead(x,&px);
420:   VecRestoreArray(y,&py);
421:   return(0);
422: }

426: PetscErrorCode MatMult_Linear_S1B(Mat B,Vec x,Vec y)
427: {
428:   PetscErrorCode    ierr;
429:   PEP_LINEAR        *ctx;
430:   const PetscScalar *px;
431:   PetscScalar       *py;
432:   PetscInt          m;

435:   MatShellGetContext(B,(void**)&ctx);
436:   MatGetLocalSize(ctx->M,&m,NULL);
437:   VecGetArrayRead(x,&px);
438:   VecGetArray(y,&py);
439:   VecPlaceArray(ctx->x1,px);
440:   VecPlaceArray(ctx->x2,px+m);
441:   VecPlaceArray(ctx->y1,py);
442:   VecPlaceArray(ctx->y2,py+m);
443:   /* y1 = -K*x1 */
444:   MatMult(ctx->K,ctx->x1,ctx->y1);
445:   VecScale(ctx->y1,-1.0);
446:   /* y2 = M*x2 */
447:   MatMult(ctx->M,ctx->x2,ctx->y2);
448:   VecScale(ctx->y2,ctx->sfactor*ctx->sfactor);
449:   VecResetArray(ctx->x1);
450:   VecResetArray(ctx->x2);
451:   VecResetArray(ctx->y1);
452:   VecResetArray(ctx->y2);
453:   VecRestoreArrayRead(x,&px);
454:   VecRestoreArray(y,&py);
455:   return(0);
456: }

460: PetscErrorCode MatGetDiagonal_Linear_S1A(Mat A,Vec diag)
461: {
463:   PEP_LINEAR     *ctx;
464:   PetscScalar    *pd;
465:   PetscInt       m;

468:   MatShellGetContext(A,(void**)&ctx);
469:   MatGetLocalSize(ctx->M,&m,NULL);
470:   VecGetArray(diag,&pd);
471:   VecPlaceArray(ctx->x1,pd);
472:   VecPlaceArray(ctx->x2,pd+m);
473:   VecSet(ctx->x1,0.0);
474:   MatGetDiagonal(ctx->C,ctx->x2);
475:   VecScale(ctx->x2,-ctx->sfactor);
476:   VecResetArray(ctx->x1);
477:   VecResetArray(ctx->x2);
478:   VecRestoreArray(diag,&pd);
479:   return(0);
480: }

484: PetscErrorCode MatGetDiagonal_Linear_S1B(Mat B,Vec diag)
485: {
487:   PEP_LINEAR     *ctx;
488:   PetscScalar    *pd;
489:   PetscInt       m;

492:   MatShellGetContext(B,(void**)&ctx);
493:   MatGetLocalSize(ctx->M,&m,NULL);
494:   VecGetArray(diag,&pd);
495:   VecPlaceArray(ctx->x1,pd);
496:   VecPlaceArray(ctx->x2,pd+m);
497:   MatGetDiagonal(ctx->K,ctx->x1);
498:   VecScale(ctx->x1,-1.0);
499:   MatGetDiagonal(ctx->M,ctx->x2);
500:   VecScale(ctx->x2,ctx->sfactor*ctx->sfactor);
501:   VecResetArray(ctx->x1);
502:   VecResetArray(ctx->x2);
503:   VecRestoreArray(diag,&pd);
504:   return(0);
505: }

509: PetscErrorCode MatCreateExplicit_Linear_S1A(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
510: {

514:   SlepcMatTile(0.0,ctx->K,-1.0,ctx->K,-1.0,ctx->K,-ctx->sfactor,ctx->C,A);
515:   return(0);
516: }

520: PetscErrorCode MatCreateExplicit_Linear_S1B(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
521: {

525:   SlepcMatTile(-1.0,ctx->K,0.0,ctx->M,0.0,ctx->M,ctx->sfactor*ctx->sfactor,ctx->M,B);
526:   return(0);
527: }

529: /* - - - S2 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

533: PetscErrorCode MatMult_Linear_S2A(Mat A,Vec x,Vec y)
534: {
535:   PetscErrorCode    ierr;
536:   PEP_LINEAR        *ctx;
537:   const PetscScalar *px;
538:   PetscScalar       *py;
539:   PetscInt          m;

542:   MatShellGetContext(A,(void**)&ctx);
543:   MatGetLocalSize(ctx->M,&m,NULL);
544:   VecGetArrayRead(x,&px);
545:   VecGetArray(y,&py);
546:   VecPlaceArray(ctx->x1,px);
547:   VecPlaceArray(ctx->x2,px+m);
548:   VecPlaceArray(ctx->y1,py);
549:   VecPlaceArray(ctx->y2,py+m);
550:   /* y1 = -K*x1 */
551:   MatMult(ctx->K,ctx->x1,ctx->y1);
552:   VecScale(ctx->y1,-1.0);
553:   /* y2 = M*x2 */
554:   MatMult(ctx->M,ctx->x2,ctx->y2);
555:   VecScale(ctx->y2,ctx->sfactor*ctx->sfactor);
556:   VecResetArray(ctx->x1);
557:   VecResetArray(ctx->x2);
558:   VecResetArray(ctx->y1);
559:   VecResetArray(ctx->y2);
560:   VecRestoreArrayRead(x,&px);
561:   VecRestoreArray(y,&py);
562:   return(0);
563: }

567: PetscErrorCode MatMult_Linear_S2B(Mat B,Vec x,Vec y)
568: {
569:   PetscErrorCode    ierr;
570:   PEP_LINEAR        *ctx;
571:   const PetscScalar *px;
572:   PetscScalar       *py;
573:   PetscInt          m;

576:   MatShellGetContext(B,(void**)&ctx);
577:   MatGetLocalSize(ctx->M,&m,NULL);
578:   VecGetArrayRead(x,&px);
579:   VecGetArray(y,&py);
580:   VecPlaceArray(ctx->x1,px);
581:   VecPlaceArray(ctx->x2,px+m);
582:   VecPlaceArray(ctx->y1,py);
583:   VecPlaceArray(ctx->y2,py+m);
584:   /* y1 = C*x1 + M*x2 */
585:   MatMult(ctx->C,ctx->x1,ctx->y1);
586:   VecScale(ctx->y1,ctx->sfactor);
587:   MatMult(ctx->M,ctx->x2,ctx->y2);
588:   VecAXPY(ctx->y1,ctx->sfactor*ctx->sfactor,ctx->y2);
589:   /* y2 = M*x1 */
590:   MatMult(ctx->M,ctx->x1,ctx->y2);
591:   VecScale(ctx->y2,ctx->sfactor*ctx->sfactor);
592:   VecResetArray(ctx->x1);
593:   VecResetArray(ctx->x2);
594:   VecResetArray(ctx->y1);
595:   VecResetArray(ctx->y2);
596:   VecRestoreArrayRead(x,&px);
597:   VecRestoreArray(y,&py);
598:   return(0);
599: }

603: PetscErrorCode MatGetDiagonal_Linear_S2A(Mat A,Vec diag)
604: {
606:   PEP_LINEAR     *ctx;
607:   PetscScalar    *pd;
608:   PetscInt       m;

611:   MatShellGetContext(A,(void**)&ctx);
612:   MatGetLocalSize(ctx->M,&m,NULL);
613:   VecGetArray(diag,&pd);
614:   VecPlaceArray(ctx->x1,pd);
615:   VecPlaceArray(ctx->x2,pd+m);
616:   MatGetDiagonal(ctx->K,ctx->x1);
617:   VecScale(ctx->x1,-1.0);
618:   MatGetDiagonal(ctx->M,ctx->x2);
619:   VecScale(ctx->x2,ctx->sfactor*ctx->sfactor);
620:   VecResetArray(ctx->x1);
621:   VecResetArray(ctx->x2);
622:   VecRestoreArray(diag,&pd);
623:   return(0);
624: }

628: PetscErrorCode MatGetDiagonal_Linear_S2B(Mat B,Vec diag)
629: {
631:   PEP_LINEAR     *ctx;
632:   PetscScalar    *pd;
633:   PetscInt       m;

636:   MatShellGetContext(B,(void**)&ctx);
637:   MatGetLocalSize(ctx->M,&m,NULL);
638:   VecGetArray(diag,&pd);
639:   VecPlaceArray(ctx->x1,pd);
640:   VecPlaceArray(ctx->x2,pd+m);
641:   MatGetDiagonal(ctx->C,ctx->x1);
642:   VecScale(ctx->x1,ctx->sfactor);
643:   VecSet(ctx->x2,0.0);
644:   VecResetArray(ctx->x1);
645:   VecResetArray(ctx->x2);
646:   VecRestoreArray(diag,&pd);
647:   return(0);
648: }

652: PetscErrorCode MatCreateExplicit_Linear_S2A(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
653: {

657:   SlepcMatTile(-1.0,ctx->K,0.0,ctx->M,0.0,ctx->M,ctx->sfactor*ctx->sfactor,ctx->M,A);
658:   return(0);
659: }

663: PetscErrorCode MatCreateExplicit_Linear_S2B(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
664: {

668:   SlepcMatTile(ctx->sfactor,ctx->C,ctx->sfactor*ctx->sfactor,ctx->M,ctx->sfactor*ctx->sfactor,ctx->M,0.0,ctx->M,B);
669:   return(0);
670: }

672: /* - - - H1 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

676: PetscErrorCode MatMult_Linear_H1A(Mat A,Vec x,Vec y)
677: {
678:   PetscErrorCode    ierr;
679:   PEP_LINEAR        *ctx;
680:   const PetscScalar *px;
681:   PetscScalar       *py;
682:   PetscInt          m;

685:   MatShellGetContext(A,(void**)&ctx);
686:   MatGetLocalSize(ctx->M,&m,NULL);
687:   VecGetArrayRead(x,&px);
688:   VecGetArray(y,&py);
689:   VecPlaceArray(ctx->x1,px);
690:   VecPlaceArray(ctx->x2,px+m);
691:   VecPlaceArray(ctx->y1,py);
692:   VecPlaceArray(ctx->y2,py+m);
693:   /* y2 = C*x1 + K*x2 */
694:   MatMult(ctx->C,ctx->x1,ctx->y1);
695:   MatMult(ctx->K,ctx->x2,ctx->y2);
696:   VecAXPY(ctx->y2,ctx->sfactor,ctx->y1);
697:   /* y1 = K*x1 */
698:   MatMult(ctx->K,ctx->x1,ctx->y1);
699:   VecResetArray(ctx->x1);
700:   VecResetArray(ctx->x2);
701:   VecResetArray(ctx->y1);
702:   VecResetArray(ctx->y2);
703:   VecRestoreArrayRead(x,&px);
704:   VecRestoreArray(y,&py);
705:   return(0);
706: }

710: PetscErrorCode MatMult_Linear_H1B(Mat B,Vec x,Vec y)
711: {
712:   PetscErrorCode    ierr;
713:   PEP_LINEAR        *ctx;
714:   const PetscScalar *px;
715:   PetscScalar       *py;
716:   PetscInt          m;

719:   MatShellGetContext(B,(void**)&ctx);
720:   MatGetLocalSize(ctx->M,&m,NULL);
721:   VecGetArrayRead(x,&px);
722:   VecGetArray(y,&py);
723:   VecPlaceArray(ctx->x1,px);
724:   VecPlaceArray(ctx->x2,px+m);
725:   VecPlaceArray(ctx->y1,py);
726:   VecPlaceArray(ctx->y2,py+m);
727:   /* y1 = K*x2 */
728:   MatMult(ctx->K,ctx->x2,ctx->y1);
729:   /* y2 = -M*x1 */
730:   MatMult(ctx->M,ctx->x1,ctx->y2);
731:   VecScale(ctx->y2,-ctx->sfactor*ctx->sfactor);
732:   VecResetArray(ctx->x1);
733:   VecResetArray(ctx->x2);
734:   VecResetArray(ctx->y1);
735:   VecResetArray(ctx->y2);
736:   VecRestoreArrayRead(x,&px);
737:   VecRestoreArray(y,&py);
738:   return(0);
739: }

743: PetscErrorCode MatGetDiagonal_Linear_H1A(Mat A,Vec diag)
744: {
746:   PEP_LINEAR     *ctx;
747:   PetscScalar    *pd;
748:   PetscInt       m;

751:   MatShellGetContext(A,(void**)&ctx);
752:   MatGetLocalSize(ctx->M,&m,NULL);
753:   VecGetArray(diag,&pd);
754:   VecPlaceArray(ctx->x1,pd);
755:   VecPlaceArray(ctx->x2,pd+m);
756:   MatGetDiagonal(ctx->K,ctx->x1);
757:   VecCopy(ctx->x1,ctx->x2);
758:   VecResetArray(ctx->x1);
759:   VecResetArray(ctx->x2);
760:   VecRestoreArray(diag,&pd);
761:   return(0);
762: }

766: PetscErrorCode MatGetDiagonal_Linear_H1B(Mat B,Vec diag)
767: {

771:   VecSet(diag,0.0);
772:   return(0);
773: }

777: PetscErrorCode MatCreateExplicit_Linear_H1A(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
778: {

782:   SlepcMatTile(1.0,ctx->K,0.0,ctx->K,ctx->sfactor,ctx->C,1.0,ctx->K,A);
783:   return(0);
784: }

788: PetscErrorCode MatCreateExplicit_Linear_H1B(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
789: {

793:   SlepcMatTile(0.0,ctx->K,1.0,ctx->K,-ctx->sfactor*ctx->sfactor,ctx->M,0.0,ctx->K,B);
794:   return(0);
795: }

797: /* - - - H2 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

801: PetscErrorCode MatMult_Linear_H2A(Mat A,Vec x,Vec y)
802: {
803:   PetscErrorCode    ierr;
804:   PEP_LINEAR        *ctx;
805:   const PetscScalar *px;
806:   PetscScalar       *py;
807:   PetscInt          m;

810:   MatShellGetContext(A,(void**)&ctx);
811:   MatGetLocalSize(ctx->M,&m,NULL);
812:   VecGetArrayRead(x,&px);
813:   VecGetArray(y,&py);
814:   VecPlaceArray(ctx->x1,px);
815:   VecPlaceArray(ctx->x2,px+m);
816:   VecPlaceArray(ctx->y1,py);
817:   VecPlaceArray(ctx->y2,py+m);
818:   /* y1 = -K*x2 */
819:   MatMult(ctx->K,ctx->x2,ctx->y1);
820:   VecScale(ctx->y1,-1.0);
821:   /* y2 = M*x1 */
822:   MatMult(ctx->M,ctx->x1,ctx->y2);
823:   VecScale(ctx->y2,ctx->sfactor*ctx->sfactor);
824:   VecResetArray(ctx->x1);
825:   VecResetArray(ctx->x2);
826:   VecResetArray(ctx->y1);
827:   VecResetArray(ctx->y2);
828:   VecRestoreArrayRead(x,&px);
829:   VecRestoreArray(y,&py);
830:   return(0);
831: }

835: PetscErrorCode MatMult_Linear_H2B(Mat B,Vec x,Vec y)
836: {
837:   PetscErrorCode    ierr;
838:   PEP_LINEAR        *ctx;
839:   const PetscScalar *px;
840:   PetscScalar       *py;
841:   PetscInt          m;

844:   MatShellGetContext(B,(void**)&ctx);
845:   MatGetLocalSize(ctx->M,&m,NULL);
846:   VecGetArrayRead(x,&px);
847:   VecGetArray(y,&py);
848:   VecPlaceArray(ctx->x1,px);
849:   VecPlaceArray(ctx->x2,px+m);
850:   VecPlaceArray(ctx->y1,py);
851:   VecPlaceArray(ctx->y2,py+m);
852:   /* y1 = M*x1 + C*x2 */
853:   MatMult(ctx->M,ctx->x1,ctx->y1);
854:   VecScale(ctx->y1,ctx->sfactor*ctx->sfactor);
855:   MatMult(ctx->C,ctx->x2,ctx->y2);
856:   VecAXPY(ctx->y1,ctx->sfactor,ctx->y2);
857:   /* y2 = M*x2 */
858:   MatMult(ctx->M,ctx->x2,ctx->y2);
859:   VecScale(ctx->y2,ctx->sfactor*ctx->sfactor);
860:   VecResetArray(ctx->x1);
861:   VecResetArray(ctx->x2);
862:   VecResetArray(ctx->y1);
863:   VecResetArray(ctx->y2);
864:   VecRestoreArrayRead(x,&px);
865:   VecRestoreArray(y,&py);
866:   return(0);
867: }

871: PetscErrorCode MatGetDiagonal_Linear_H2A(Mat A,Vec diag)
872: {

876:   VecSet(diag,0.0);
877:   return(0);
878: }

882: PetscErrorCode MatGetDiagonal_Linear_H2B(Mat B,Vec diag)
883: {
885:   PEP_LINEAR     *ctx;
886:   PetscScalar    *pd;
887:   PetscInt       m;

890:   MatShellGetContext(B,(void**)&ctx);
891:   MatGetLocalSize(ctx->M,&m,NULL);
892:   VecGetArray(diag,&pd);
893:   VecPlaceArray(ctx->x1,pd);
894:   VecPlaceArray(ctx->x2,pd+m);
895:   MatGetDiagonal(ctx->M,ctx->x1);
896:   VecScale(ctx->x1,ctx->sfactor*ctx->sfactor);
897:   VecCopy(ctx->x1,ctx->x2);
898:   VecResetArray(ctx->x1);
899:   VecResetArray(ctx->x2);
900:   VecRestoreArray(diag,&pd);
901:   return(0);
902: }

906: PetscErrorCode MatCreateExplicit_Linear_H2A(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
907: {

911:   SlepcMatTile(0.0,ctx->K,-1.0,ctx->K,ctx->sfactor*ctx->sfactor,ctx->M,0.0,ctx->K,A);
912:   return(0);
913: }

917: PetscErrorCode MatCreateExplicit_Linear_H2B(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
918: {

922:   SlepcMatTile(ctx->sfactor*ctx->sfactor,ctx->M,ctx->sfactor,ctx->C,0.0,ctx->C,ctx->sfactor*ctx->sfactor,ctx->M,B);
923:   return(0);
924: }