Actual source code: qeplin_s1.c
1: /*
3: Linearization for Hermitian QEP, companion form 1.
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 = [-K 0 ] z = [ x ]
34: [ -K -C ] [ 0 M ] [ l*x ]
35: */
39: PetscErrorCode MatMult_QEPLINEAR_S1A(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: /* y2 = -(K*x1 + C*x2) */
56: MatMult(ctx->K,ctx->x1,ctx->y2);
57: VecScale(ctx->y2,-1.0);
58: MatMult(ctx->C,ctx->x2,ctx->y1);
59: VecAXPY(ctx->y2,-ctx->sfactor,ctx->y1);
60: /* y1 = -K*x2 */
61: MatMult(ctx->K,ctx->x2,ctx->y1);
62: VecScale(ctx->y1,-1.0);
63: VecResetArray(ctx->x1);
64: VecResetArray(ctx->x2);
65: VecResetArray(ctx->y1);
66: VecResetArray(ctx->y2);
67: VecRestoreArray(x,&px);
68: VecRestoreArray(y,&py);
69: return(0);
70: }
74: PetscErrorCode MatMult_QEPLINEAR_S1B(Mat B,Vec x,Vec y)
75: {
77: QEP_LINEAR *ctx;
78: PetscScalar *px,*py;
79: PetscInt m;
80:
82: MatShellGetContext(B,(void**)&ctx);
83: MatGetLocalSize(ctx->M,&m,PETSC_NULL);
84: VecGetArray(x,&px);
85: VecGetArray(y,&py);
86: VecPlaceArray(ctx->x1,px);
87: VecPlaceArray(ctx->x2,px+m);
88: VecPlaceArray(ctx->y1,py);
89: VecPlaceArray(ctx->y2,py+m);
90: /* y1 = -K*x1 */
91: MatMult(ctx->K,ctx->x1,ctx->y1);
92: VecScale(ctx->y1,-1.0);
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_S1A(Mat A,Vec diag)
108: {
110: QEP_LINEAR *ctx;
111: PetscScalar *pd;
112: PetscInt m;
113:
115: MatShellGetContext(A,(void**)&ctx);
116: MatGetLocalSize(ctx->M,&m,PETSC_NULL);
117: VecGetArray(diag,&pd);
118: VecPlaceArray(ctx->x1,pd);
119: VecPlaceArray(ctx->x2,pd+m);
120: VecSet(ctx->x1,0.0);
121: MatGetDiagonal(ctx->C,ctx->x2);
122: VecScale(ctx->x2,-ctx->sfactor);
123: VecResetArray(ctx->x1);
124: VecResetArray(ctx->x2);
125: VecRestoreArray(diag,&pd);
126: return(0);
127: }
131: PetscErrorCode MatGetDiagonal_QEPLINEAR_S1B(Mat B,Vec diag)
132: {
134: QEP_LINEAR *ctx;
135: PetscScalar *pd;
136: PetscInt m;
137:
139: MatShellGetContext(B,(void**)&ctx);
140: MatGetLocalSize(ctx->M,&m,PETSC_NULL);
141: VecGetArray(diag,&pd);
142: VecPlaceArray(ctx->x1,pd);
143: VecPlaceArray(ctx->x2,pd+m);
144: MatGetDiagonal(ctx->K,ctx->x1);
145: VecScale(ctx->x1,-1.0);
146: MatGetDiagonal(ctx->M,ctx->x2);
147: VecScale(ctx->x2,ctx->sfactor*ctx->sfactor);
148: VecResetArray(ctx->x1);
149: VecResetArray(ctx->x2);
150: VecRestoreArray(diag,&pd);
151: return(0);
152: }
156: PetscErrorCode MatCreateExplicit_QEPLINEAR_S1A(MPI_Comm comm,QEP_LINEAR *ctx,Mat *A)
157: {
159: PetscInt M,N,m,n,i,j,row,start,end,ncols,*pos;
160: PetscScalar *svals;
161: const PetscInt *cols;
162: const PetscScalar *vals;
163:
165: MatGetSize(ctx->M,&M,&N);
166: MatGetLocalSize(ctx->M,&m,&n);
167: MatCreate(comm,A);
168: MatSetSizes(*A,m+n,m+n,M+N,M+N);
169: MatSetFromOptions(*A);
170: PetscMalloc(sizeof(PetscInt)*n,&pos);
171: PetscMalloc(sizeof(PetscScalar)*n,&svals);
172: MatGetOwnershipRange(ctx->M,&start,&end);
173: for (i=start;i<end;i++) {
174: row = i + M;
175: MatGetRow(ctx->K,i,&ncols,&cols,&vals);
176: MatSetValues(*A,1,&row,ncols,cols,vals,INSERT_VALUES);
177: for (j=0;j<ncols;j++)
178: pos[j] = cols[j] + M;
179: MatSetValues(*A,1,&i,ncols,pos,vals,INSERT_VALUES);
180: MatRestoreRow(ctx->K,i,&ncols,&cols,&vals);
181: MatGetRow(ctx->C,i,&ncols,&cols,&vals);
182: for (j=0;j<ncols;j++) {
183: pos[j] = cols[j] + M;
184: svals[j] = vals[j]*ctx->sfactor;
185: }
186: MatSetValues(*A,1,&row,ncols,pos,svals,INSERT_VALUES);
187: MatRestoreRow(ctx->C,i,&ncols,&cols,&vals);
188: }
189: PetscFree(pos);
190: PetscFree(svals);
191: MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
192: MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
193: MatScale(*A,-1.0);
194: return(0);
195: }
199: PetscErrorCode MatCreateExplicit_QEPLINEAR_S1B(MPI_Comm comm,QEP_LINEAR *ctx,Mat *B)
200: {
202: PetscInt M,N,m,n,i,j,row,start,end,ncols,*pos;
203: PetscScalar *svals;
204: const PetscInt *cols;
205: const PetscScalar *vals;
206:
208: MatGetSize(ctx->M,&M,&N);
209: MatGetLocalSize(ctx->M,&m,&n);
210: MatCreate(comm,B);
211: MatSetSizes(*B,m+n,m+n,M+N,M+N);
212: MatSetFromOptions(*B);
213: PetscMalloc(sizeof(PetscInt)*n,&pos);
214: PetscMalloc(sizeof(PetscScalar)*n,&svals);
215: MatGetOwnershipRange(ctx->M,&start,&end);
216: for (i=start;i<end;i++) {
217: MatGetRow(ctx->K,i,&ncols,&cols,&vals);
218: MatSetValues(*B,1,&i,ncols,cols,vals,INSERT_VALUES);
219: MatRestoreRow(ctx->K,i,&ncols,&cols,&vals);
220: }
221: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
222: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
223: MatScale(*B,-1.0);
224: for (i=start;i<end;i++) {
225: row = i + M;
226: MatGetRow(ctx->M,i,&ncols,&cols,&vals);
227: for (j=0;j<ncols;j++) {
228: pos[j] = cols[j] + M;
229: svals[j] = vals[j]*ctx->sfactor*ctx->sfactor;
230: }
231: MatSetValues(*B,1,&row,ncols,pos,svals,INSERT_VALUES);
232: MatRestoreRow(ctx->M,i,&ncols,&cols,&vals);
233: }
234: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
235: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
236: PetscFree(pos);
237: PetscFree(svals);
238: return(0);
239: }