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-2010, Universidad 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
 26:  #include slepceps.h
 27:  #include linearp.h

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

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

 39: PetscErrorCode MatMult_QEPLINEAR_H2A(Mat A,Vec x,Vec y)
 40: {
 42:   QEP_LINEAR     *ctx;
 43:   PetscScalar    *px,*py;
 44:   PetscInt       m;
 45: 
 47:   MatShellGetContext(A,(void**)&ctx);
 48:   MatGetLocalSize(ctx->M,&m,PETSC_NULL);
 49:   VecGetArray(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:   VecRestoreArray(x,&px);
 66:   VecRestoreArray(y,&py);
 67:   return(0);
 68: }

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

107: PetscErrorCode MatGetDiagonal_QEPLINEAR_H2A(Mat A,Vec diag)
108: {
110: 
112:   VecSet(diag,0.0);
113:   return(0);
114: }

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

142: PetscErrorCode MatCreateExplicit_QEPLINEAR_H2A(MPI_Comm comm,QEP_LINEAR *ctx,Mat *A)
143: {
145:   PetscInt       M,N,m,n,i,j,row,start,end,ncols,*pos;
146:   PetscScalar    *svals;
147:   const PetscInt    *cols;
148:   const PetscScalar *vals;
149: 
151:   MatGetSize(ctx->M,&M,&N);
152:   MatGetLocalSize(ctx->M,&m,&n);
153:   MatCreate(comm,A);
154:   MatSetSizes(*A,m+n,m+n,M+N,M+N);
155:   MatSetFromOptions(*A);
156:   PetscMalloc(sizeof(PetscInt)*n,&pos);
157:   PetscMalloc(sizeof(PetscScalar)*n,&svals);
158:   MatGetOwnershipRange(ctx->M,&start,&end);
159:   for (i=start;i<end;i++) {
160:     MatGetRow(ctx->K,i,&ncols,&cols,&vals);
161:     for (j=0;j<ncols;j++)
162:       pos[j] = cols[j] + M;
163:     MatSetValues(*A,1,&i,ncols,pos,vals,INSERT_VALUES);
164:     MatRestoreRow(ctx->K,i,&ncols,&cols,&vals);
165:   }
166:   MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
167:   MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
168:   MatScale(*A,-1.0);
169:   for (i=start;i<end;i++) {
170:     row = i + M;
171:     MatGetRow(ctx->M,i,&ncols,&cols,&vals);
172:     for (j=0;j<ncols;j++)
173:       svals[j] = vals[j]*ctx->sfactor*ctx->sfactor;
174:     MatSetValues(*A,1,&row,ncols,cols,svals,INSERT_VALUES);
175:     MatRestoreRow(ctx->M,i,&ncols,&cols,&vals);
176:   }
177:   MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
178:   MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
179:   PetscFree(pos);
180:   PetscFree(svals);
181:   return(0);
182: }

186: PetscErrorCode MatCreateExplicit_QEPLINEAR_H2B(MPI_Comm comm,QEP_LINEAR *ctx,Mat *B)
187: {
189:   PetscInt       M,N,m,n,i,j,row,start,end,ncols,*pos;
190:   PetscScalar    *svals;
191:   const PetscInt    *cols;
192:   const PetscScalar *vals;
193: 
195:   MatGetSize(ctx->M,&M,&N);
196:   MatGetLocalSize(ctx->M,&m,&n);
197:   MatCreate(comm,B);
198:   MatSetSizes(*B,m+n,m+n,M+N,M+N);
199:   MatSetFromOptions(*B);
200:   PetscMalloc(sizeof(PetscInt)*n,&pos);
201:   PetscMalloc(sizeof(PetscScalar)*n,&svals);
202:   MatGetOwnershipRange(ctx->M,&start,&end);
203:   for (i=start;i<end;i++) {
204:     row = i + M;
205:     MatGetRow(ctx->C,i,&ncols,&cols,&vals);
206:     for (j=0;j<ncols;j++) {
207:       pos[j] = cols[j] + M;
208:       svals[j] = vals[j]*ctx->sfactor;
209:     }
210:     MatSetValues(*B,1,&i,ncols,pos,svals,INSERT_VALUES);
211:     MatRestoreRow(ctx->C,i,&ncols,&cols,&vals);
212:     MatGetRow(ctx->M,i,&ncols,&cols,&vals);
213:     for (j=0;j<ncols;j++) {
214:       pos[j] = cols[j] + M;
215:       svals[j] = vals[j]*ctx->sfactor*ctx->sfactor;
216:     }
217:     MatSetValues(*B,1,&i,ncols,cols,svals,INSERT_VALUES);
218:     MatSetValues(*B,1,&row,ncols,pos,svals,INSERT_VALUES);
219:     MatRestoreRow(ctx->M,i,&ncols,&cols,&vals);
220:   }
221:   PetscFree(pos);
222:   PetscFree(svals);
223:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
224:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
225:   return(0);
226: }