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-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 = [-K 0 ] z = [ x ]
33: [ -K -C ] [ 0 M ] [ l*x ]
34: */
38: PetscErrorCode MatMult_Linear_S1A(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: /* 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: VecRestoreArrayRead(x,&px);
68: VecRestoreArray(y,&py);
69: return(0);
70: }
74: PetscErrorCode MatMult_Linear_S1B(Mat B,Vec x,Vec y)
75: {
76: PetscErrorCode ierr;
77: QEP_LINEAR *ctx;
78: const PetscScalar *px;
79: PetscScalar *py;
80: PetscInt m;
81:
83: MatShellGetContext(B,(void**)&ctx);
84: MatGetLocalSize(ctx->M,&m,PETSC_NULL);
85: VecGetArrayRead(x,&px);
86: VecGetArray(y,&py);
87: VecPlaceArray(ctx->x1,px);
88: VecPlaceArray(ctx->x2,px+m);
89: VecPlaceArray(ctx->y1,py);
90: VecPlaceArray(ctx->y2,py+m);
91: /* y1 = -K*x1 */
92: MatMult(ctx->K,ctx->x1,ctx->y1);
93: VecScale(ctx->y1,-1.0);
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_S1A(Mat A,Vec diag)
109: {
111: QEP_LINEAR *ctx;
112: PetscScalar *pd;
113: PetscInt m;
114:
116: MatShellGetContext(A,(void**)&ctx);
117: MatGetLocalSize(ctx->M,&m,PETSC_NULL);
118: VecGetArray(diag,&pd);
119: VecPlaceArray(ctx->x1,pd);
120: VecPlaceArray(ctx->x2,pd+m);
121: VecSet(ctx->x1,0.0);
122: MatGetDiagonal(ctx->C,ctx->x2);
123: VecScale(ctx->x2,-ctx->sfactor);
124: VecResetArray(ctx->x1);
125: VecResetArray(ctx->x2);
126: VecRestoreArray(diag,&pd);
127: return(0);
128: }
132: PetscErrorCode MatGetDiagonal_Linear_S1B(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->K,ctx->x1);
146: VecScale(ctx->x1,-1.0);
147: MatGetDiagonal(ctx->M,ctx->x2);
148: VecScale(ctx->x2,ctx->sfactor*ctx->sfactor);
149: VecResetArray(ctx->x1);
150: VecResetArray(ctx->x2);
151: VecRestoreArray(diag,&pd);
152: return(0);
153: }
157: PetscErrorCode MatCreateExplicit_Linear_S1A(MPI_Comm comm,QEP_LINEAR *ctx,Mat *A)
158: {
160:
162: SlepcMatTile(0.0,ctx->K,-1.0,ctx->K,-1.0,ctx->K,-ctx->sfactor,ctx->C,A);
163: return(0);
164: }
168: PetscErrorCode MatCreateExplicit_Linear_S1B(MPI_Comm comm,QEP_LINEAR *ctx,Mat *B)
169: {
171:
173: SlepcMatTile(-1.0,ctx->K,0.0,ctx->M,0.0,ctx->M,ctx->sfactor*ctx->sfactor,ctx->M,B);
174: return(0);
175: }