Actual source code: qeplin_h1.c
1: /*
3: Linearization for gyroscopic 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 = [ K 0 ] B = [ 0 K ] z = [ x ]
33: [ C K ] [-M 0 ] [ l*x ]
34: */
38: PetscErrorCode MatMult_Linear_H1A(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 = C*x1 + K*x2 */
56: MatMult(ctx->C,ctx->x1,ctx->y1);
57: MatMult(ctx->K,ctx->x2,ctx->y2);
58: VecAXPY(ctx->y2,ctx->sfactor,ctx->y1);
59: /* y1 = K*x1 */
60: MatMult(ctx->K,ctx->x1,ctx->y1);
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_H1B(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 = K*x2 */
90: MatMult(ctx->K,ctx->x2,ctx->y1);
91: /* y2 = -M*x1 */
92: MatMult(ctx->M,ctx->x1,ctx->y2);
93: VecScale(ctx->y2,-ctx->sfactor*ctx->sfactor);
94: VecResetArray(ctx->x1);
95: VecResetArray(ctx->x2);
96: VecResetArray(ctx->y1);
97: VecResetArray(ctx->y2);
98: VecRestoreArrayRead(x,&px);
99: VecRestoreArray(y,&py);
100: return(0);
101: }
105: PetscErrorCode MatGetDiagonal_Linear_H1A(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: VecCopy(ctx->x1,ctx->x2);
120: VecResetArray(ctx->x1);
121: VecResetArray(ctx->x2);
122: VecRestoreArray(diag,&pd);
123: return(0);
124: }
128: PetscErrorCode MatGetDiagonal_Linear_H1B(Mat B,Vec diag)
129: {
131:
133: VecSet(diag,0.0);
134: return(0);
135: }
139: PetscErrorCode MatCreateExplicit_Linear_H1A(MPI_Comm comm,QEP_LINEAR *ctx,Mat *A)
140: {
142:
144: SlepcMatTile(1.0,ctx->K,0.0,ctx->K,ctx->sfactor,ctx->C,1.0,ctx->K,A);
145: return(0);
146: }
150: PetscErrorCode MatCreateExplicit_Linear_H1B(MPI_Comm comm,QEP_LINEAR *ctx,Mat *B)
151: {
153:
155: SlepcMatTile(0.0,ctx->K,1.0,ctx->K,-ctx->sfactor*ctx->sfactor,ctx->M,0.0,ctx->K,B);
156: return(0);
157: }