Actual source code: qeplin_h2.c

  1: /*                       

  3:    Linearization for gyroscopic QEP, companion form 2.

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

  9:    This file is part of SLEPc.
 10:       
 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 <private/qepimpl.h>         /*I "slepcqep.h" I*/
 26:  #include linearp.h

 28: /*
 29:     Given the quadratic problem (l^2*M + l*C + K)*x = 0 the following
 30:     linearization is employed:

 32:       A*z = l*B*z   where   A = [  0  -K ]     B = [ M  C ]     z = [  x  ]
 33:                                 [  M   0 ]         [ 0  M ]         [ l*x ]
 34:  */

 38: PetscErrorCode MatMult_Linear_H2A(Mat A,Vec x,Vec y)
 39: {
 40:   PetscErrorCode    ierr;
 41:   QEP_LINEAR        *ctx;
 42:   const PetscScalar *px;
 43:   PetscScalar       *py;
 44:   PetscInt          m;
 45: 
 47:   MatShellGetContext(A,(void**)&ctx);
 48:   MatGetLocalSize(ctx->M,&m,PETSC_NULL);
 49:   VecGetArrayRead(x,&px);
 50:   VecGetArray(y,&py);
 51:   VecPlaceArray(ctx->x1,px);
 52:   VecPlaceArray(ctx->x2,px+m);
 53:   VecPlaceArray(ctx->y1,py);
 54:   VecPlaceArray(ctx->y2,py+m);
 55:   /* y1 = -K*x2 */
 56:   MatMult(ctx->K,ctx->x2,ctx->y1);
 57:   VecScale(ctx->y1,-1.0);
 58:   /* y2 = M*x1 */
 59:   MatMult(ctx->M,ctx->x1,ctx->y2);
 60:   VecScale(ctx->y2,ctx->sfactor*ctx->sfactor);
 61:   VecResetArray(ctx->x1);
 62:   VecResetArray(ctx->x2);
 63:   VecResetArray(ctx->y1);
 64:   VecResetArray(ctx->y2);
 65:   VecRestoreArrayRead(x,&px);
 66:   VecRestoreArray(y,&py);
 67:   return(0);
 68: }

 72: PetscErrorCode MatMult_Linear_H2B(Mat B,Vec x,Vec y)
 73: {
 74:   PetscErrorCode    ierr;
 75:   QEP_LINEAR        *ctx;
 76:   const PetscScalar *px;
 77:   PetscScalar       *py;
 78:   PetscInt          m;
 79: 
 81:   MatShellGetContext(B,(void**)&ctx);
 82:   MatGetLocalSize(ctx->M,&m,PETSC_NULL);
 83:   VecGetArrayRead(x,&px);
 84:   VecGetArray(y,&py);
 85:   VecPlaceArray(ctx->x1,px);
 86:   VecPlaceArray(ctx->x2,px+m);
 87:   VecPlaceArray(ctx->y1,py);
 88:   VecPlaceArray(ctx->y2,py+m);
 89:   /* y1 = M*x1 + C*x2 */
 90:   MatMult(ctx->M,ctx->x1,ctx->y1);
 91:   VecScale(ctx->y1,ctx->sfactor*ctx->sfactor);
 92:   MatMult(ctx->C,ctx->x2,ctx->y2);
 93:   VecAXPY(ctx->y1,ctx->sfactor,ctx->y2);
 94:   /* y2 = M*x2 */
 95:   MatMult(ctx->M,ctx->x2,ctx->y2);
 96:   VecScale(ctx->y2,ctx->sfactor*ctx->sfactor);
 97:   VecResetArray(ctx->x1);
 98:   VecResetArray(ctx->x2);
 99:   VecResetArray(ctx->y1);
100:   VecResetArray(ctx->y2);
101:   VecRestoreArrayRead(x,&px);
102:   VecRestoreArray(y,&py);
103:   return(0);
104: }

108: PetscErrorCode MatGetDiagonal_Linear_H2A(Mat A,Vec diag)
109: {
111: 
113:   VecSet(diag,0.0);
114:   return(0);
115: }

119: PetscErrorCode MatGetDiagonal_Linear_H2B(Mat B,Vec diag)
120: {
122:   QEP_LINEAR     *ctx;
123:   PetscScalar    *pd;
124:   PetscInt       m;
125: 
127:   MatShellGetContext(B,(void**)&ctx);
128:   MatGetLocalSize(ctx->M,&m,PETSC_NULL);
129:   VecGetArray(diag,&pd);
130:   VecPlaceArray(ctx->x1,pd);
131:   VecPlaceArray(ctx->x2,pd+m);
132:   MatGetDiagonal(ctx->M,ctx->x1);
133:   VecScale(ctx->x1,ctx->sfactor*ctx->sfactor);
134:   VecCopy(ctx->x1,ctx->x2);
135:   VecResetArray(ctx->x1);
136:   VecResetArray(ctx->x2);
137:   VecRestoreArray(diag,&pd);
138:   return(0);
139: }

143: PetscErrorCode MatCreateExplicit_Linear_H2A(MPI_Comm comm,QEP_LINEAR *ctx,Mat *A)
144: {
146: 
148:   SlepcMatTile(0.0,ctx->K,-1.0,ctx->K,ctx->sfactor*ctx->sfactor,ctx->M,0.0,ctx->K,A);
149:   return(0);
150: }

154: PetscErrorCode MatCreateExplicit_Linear_H2B(MPI_Comm comm,QEP_LINEAR *ctx,Mat *B)
155: {
157: 
159:   SlepcMatTile(ctx->sfactor*ctx->sfactor,ctx->M,ctx->sfactor,ctx->C,0.0,ctx->C,ctx->sfactor*ctx->sfactor,ctx->M,B);
160:   return(0);
161: }