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-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 = [ -K 0 ] B = [ C M ] z = [ x ]
33: [ 0 M ] [ M 0 ] [ l*x ]
34: */
38: PetscErrorCode MatMult_Linear_S2A(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: /* 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: VecRestoreArrayRead(x,&px);
66: VecRestoreArray(y,&py);
67: return(0);
68: }
72: PetscErrorCode MatMult_Linear_S2B(Mat B,Vec x,Vec y)
73: {
74: PetscErrorCode ierr;
75: QEP_LINEAR *ctx;
76: const PetscScalar *px;
77: PetscScalar *py;
78: PetscInt m;
79:
81: MatShellGetContext(B,(void**)&ctx);
82: MatGetLocalSize(ctx->M,&m,PETSC_NULL);
83: VecGetArrayRead(x,&px);
84: VecGetArray(y,&py);
85: VecPlaceArray(ctx->x1,px);
86: VecPlaceArray(ctx->x2,px+m);
87: VecPlaceArray(ctx->y1,py);
88: VecPlaceArray(ctx->y2,py+m);
89: /* y1 = C*x1 + M*x2 */
90: MatMult(ctx->C,ctx->x1,ctx->y1);
91: VecScale(ctx->y1,ctx->sfactor);
92: MatMult(ctx->M,ctx->x2,ctx->y2);
93: VecAXPY(ctx->y1,ctx->sfactor*ctx->sfactor,ctx->y2);
94: /* y2 = M*x1 */
95: MatMult(ctx->M,ctx->x1,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_S2A(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: MatGetDiagonal(ctx->K,ctx->x1);
122: VecScale(ctx->x1,-1.0);
123: MatGetDiagonal(ctx->M,ctx->x2);
124: VecScale(ctx->x2,ctx->sfactor*ctx->sfactor);
125: VecResetArray(ctx->x1);
126: VecResetArray(ctx->x2);
127: VecRestoreArray(diag,&pd);
128: return(0);
129: }
133: PetscErrorCode MatGetDiagonal_Linear_S2B(Mat B,Vec diag)
134: {
136: QEP_LINEAR *ctx;
137: PetscScalar *pd;
138: PetscInt m;
139:
141: MatShellGetContext(B,(void**)&ctx);
142: MatGetLocalSize(ctx->M,&m,PETSC_NULL);
143: VecGetArray(diag,&pd);
144: VecPlaceArray(ctx->x1,pd);
145: VecPlaceArray(ctx->x2,pd+m);
146: MatGetDiagonal(ctx->C,ctx->x1);
147: VecScale(ctx->x1,ctx->sfactor);
148: VecSet(ctx->x2,0.0);
149: VecResetArray(ctx->x1);
150: VecResetArray(ctx->x2);
151: VecRestoreArray(diag,&pd);
152: return(0);
153: }
157: PetscErrorCode MatCreateExplicit_Linear_S2A(MPI_Comm comm,QEP_LINEAR *ctx,Mat *A)
158: {
160:
162: SlepcMatTile(-1.0,ctx->K,0.0,ctx->M,0.0,ctx->M,ctx->sfactor*ctx->sfactor,ctx->M,A);
163: return(0);
164: }
168: PetscErrorCode MatCreateExplicit_Linear_S2B(MPI_Comm comm,QEP_LINEAR *ctx,Mat *B)
169: {
173: SlepcMatTile(ctx->sfactor,ctx->C,ctx->sfactor*ctx->sfactor,ctx->M,ctx->sfactor*ctx->sfactor,ctx->M,0.0,ctx->M,B);
174: return(0);
175: }