Actual source code: qeplin_s2.c
1: /*
3: Linearization for Hermitian 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 M ] [ M 0 ] [ l*x ]
35: */
39: PetscErrorCode MatMult_QEPLINEAR_S2A(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 = M*x2 */
59: MatMult(ctx->M,ctx->x2,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_S2B(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 = C*x1 + M*x2 */
89: MatMult(ctx->C,ctx->x1,ctx->y1);
90: VecScale(ctx->y1,ctx->sfactor);
91: MatMult(ctx->M,ctx->x2,ctx->y2);
92: VecAXPY(ctx->y1,ctx->sfactor*ctx->sfactor,ctx->y2);
93: /* y2 = M*x1 */
94: MatMult(ctx->M,ctx->x1,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_S2A(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: MatGetDiagonal(ctx->K,ctx->x1);
121: VecScale(ctx->x1,-1.0);
122: MatGetDiagonal(ctx->M,ctx->x2);
123: VecScale(ctx->x2,ctx->sfactor*ctx->sfactor);
124: VecResetArray(ctx->x1);
125: VecResetArray(ctx->x2);
126: VecRestoreArray(diag,&pd);
127: return(0);
128: }
132: PetscErrorCode MatGetDiagonal_QEPLINEAR_S2B(Mat B,Vec diag)
133: {
135: QEP_LINEAR *ctx;
136: PetscScalar *pd;
137: PetscInt m;
138:
140: MatShellGetContext(B,(void**)&ctx);
141: MatGetLocalSize(ctx->M,&m,PETSC_NULL);
142: VecGetArray(diag,&pd);
143: VecPlaceArray(ctx->x1,pd);
144: VecPlaceArray(ctx->x2,pd+m);
145: MatGetDiagonal(ctx->C,ctx->x1);
146: VecScale(ctx->x1,ctx->sfactor);
147: VecSet(ctx->x2,0.0);
148: VecResetArray(ctx->x1);
149: VecResetArray(ctx->x2);
150: VecRestoreArray(diag,&pd);
151: return(0);
152: }
156: PetscErrorCode MatCreateExplicit_QEPLINEAR_S2A(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: MatGetRow(ctx->K,i,&ncols,&cols,&vals);
175: MatSetValues(*A,1,&i,ncols,cols,vals,INSERT_VALUES);
176: MatRestoreRow(ctx->K,i,&ncols,&cols,&vals);
177: }
178: MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
179: MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
180: MatScale(*A,-1.0);
181: for (i=start;i<end;i++) {
182: row = i + M;
183: MatGetRow(ctx->M,i,&ncols,&cols,&vals);
184: for (j=0;j<ncols;j++) {
185: pos[j] = cols[j] + M;
186: svals[j] = vals[j]*ctx->sfactor*ctx->sfactor;
187: }
188: MatSetValues(*A,1,&row,ncols,pos,svals,INSERT_VALUES);
189: MatRestoreRow(ctx->M,i,&ncols,&cols,&vals);
190: }
191: MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
192: MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
193: PetscFree(pos);
194: PetscFree(svals);
195: return(0);
196: }
200: PetscErrorCode MatCreateExplicit_QEPLINEAR_S2B(MPI_Comm comm,QEP_LINEAR *ctx,Mat *B)
201: {
203: PetscInt M,N,m,n,i,j,row,start,end,ncols,*pos;
204: PetscScalar *svals;
205: const PetscInt *cols;
206: const PetscScalar *vals;
207:
209: MatGetSize(ctx->M,&M,&N);
210: MatGetLocalSize(ctx->M,&m,&n);
211: MatCreate(comm,B);
212: MatSetSizes(*B,m+n,m+n,M+N,M+N);
213: MatSetFromOptions(*B);
214: PetscMalloc(sizeof(PetscInt)*n,&pos);
215: PetscMalloc(sizeof(PetscScalar)*n,&svals);
216: MatGetOwnershipRange(ctx->M,&start,&end);
217: for (i=start;i<end;i++) {
218: row = i + M;
219: MatGetRow(ctx->C,i,&ncols,&cols,&vals);
220: for (j=0;j<ncols;j++) {
221: svals[j] = vals[j]*ctx->sfactor;
222: }
223: MatSetValues(*B,1,&i,ncols,cols,svals,INSERT_VALUES);
224: MatRestoreRow(ctx->C,i,&ncols,&cols,&vals);
225: MatGetRow(ctx->M,i,&ncols,&cols,&vals);
226: for (j=0;j<ncols;j++) {
227: pos[j] = cols[j] + M;
228: svals[j] = vals[j]*ctx->sfactor*ctx->sfactor;
229: }
230: MatSetValues(*B,1,&i,ncols,pos,svals,INSERT_VALUES);
231: MatSetValues(*B,1,&row,ncols,cols,svals,INSERT_VALUES);
232: MatRestoreRow(ctx->M,i,&ncols,&cols,&vals);
233: }
234: PetscFree(pos);
235: PetscFree(svals);
236: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
237: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
238: return(0);
239: }