Actual source code: qeplin_n2.c

  1: /*                       

  3:    Linearization for general 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 = [ -K   0 ]     B = [ C  M ]     z = [  x  ]
 34:                                 [  0   I ]         [ I  0 ]         [ l*x ]
 35:  */

 39: PetscErrorCode MatMult_QEPLINEAR_N2A(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*x1 */
 56:   MatMult(ctx->K,ctx->x1,ctx->y1);
 57:   VecScale(ctx->y1,-1.0);
 58:   /* y2 = x2 */
 59:   VecCopy(ctx->x2,ctx->y2);
 60:   VecResetArray(ctx->x1);
 61:   VecResetArray(ctx->x2);
 62:   VecResetArray(ctx->y1);
 63:   VecResetArray(ctx->y2);
 64:   VecRestoreArray(x,&px);
 65:   VecRestoreArray(y,&py);
 66:   return(0);
 67: }

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

105: PetscErrorCode MatGetDiagonal_QEPLINEAR_N2A(Mat A,Vec diag)
106: {
108:   QEP_LINEAR     *ctx;
109:   PetscScalar    *pd;
110:   PetscInt       m;
111: 
113:   MatShellGetContext(A,(void**)&ctx);
114:   MatGetLocalSize(ctx->M,&m,PETSC_NULL);
115:   VecGetArray(diag,&pd);
116:   VecPlaceArray(ctx->x1,pd);
117:   VecPlaceArray(ctx->x2,pd+m);
118:   MatGetDiagonal(ctx->K,ctx->x1);
119:   VecScale(ctx->x1,-1.0);
120:   VecSet(ctx->x2,1.0);
121:   VecResetArray(ctx->x1);
122:   VecResetArray(ctx->x2);
123:   VecRestoreArray(diag,&pd);
124:   return(0);
125: }

129: PetscErrorCode MatGetDiagonal_QEPLINEAR_N2B(Mat B,Vec diag)
130: {
132:   QEP_LINEAR     *ctx;
133:   PetscScalar    *pd;
134:   PetscInt       m;
135: 
137:   MatShellGetContext(B,(void**)&ctx);
138:   MatGetLocalSize(ctx->M,&m,PETSC_NULL);
139:   VecGetArray(diag,&pd);
140:   VecPlaceArray(ctx->x1,pd);
141:   VecPlaceArray(ctx->x2,pd+m);
142:   MatGetDiagonal(ctx->C,ctx->x1);
143:   VecScale(ctx->x1,ctx->sfactor);
144:   VecSet(ctx->x2,0.0);
145:   VecResetArray(ctx->x1);
146:   VecResetArray(ctx->x2);
147:   VecRestoreArray(diag,&pd);
148:   return(0);
149: }

153: PetscErrorCode MatCreateExplicit_QEPLINEAR_N2A(MPI_Comm comm,QEP_LINEAR *ctx,Mat *A)
154: {
156:   PetscInt       M,N,m,n,i,start,end,ncols;
157:   const PetscInt    *cols;
158:   const PetscScalar *vals;
159: 
161:   MatGetSize(ctx->M,&M,&N);
162:   MatGetLocalSize(ctx->M,&m,&n);
163:   MatCreate(comm,A);
164:   MatSetSizes(*A,m+n,m+n,M+N,M+N);
165:   MatSetFromOptions(*A);
166:   MatGetOwnershipRange(ctx->M,&start,&end);
167:   for (i=start;i<end;i++) {
168:     MatGetRow(ctx->K,i,&ncols,&cols,&vals);
169:     MatSetValues(*A,1,&i,ncols,cols,vals,INSERT_VALUES);
170:     MatRestoreRow(ctx->K,i,&ncols,&cols,&vals);
171:     MatSetValue(*A,i+M,i+M,-1.0,INSERT_VALUES);
172:   }
173:   MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
174:   MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
175:   MatScale(*A,-1.0);
176:   return(0);
177: }

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