Actual source code: qeplin_n2.c
1: /*
3: Linearization for general 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 I ] [ I 0 ] [ l*x ]
34: */
38: PetscErrorCode MatMult_Linear_N2A(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 = x2 */
59: VecCopy(ctx->x2,ctx->y2);
60: VecResetArray(ctx->x1);
61: VecResetArray(ctx->x2);
62: VecResetArray(ctx->y1);
63: VecResetArray(ctx->y2);
64: VecRestoreArrayRead(x,&px);
65: VecRestoreArray(y,&py);
66: return(0);
67: }
71: PetscErrorCode MatMult_Linear_N2B(Mat B,Vec x,Vec y)
72: {
73: PetscErrorCode ierr;
74: QEP_LINEAR *ctx;
75: const PetscScalar *px;
76: PetscScalar *py;
77: PetscInt m;
78:
80: MatShellGetContext(B,(void**)&ctx);
81: MatGetLocalSize(ctx->M,&m,PETSC_NULL);
82: VecGetArrayRead(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 = x1 */
94: VecCopy(ctx->x1,ctx->y2);
95: VecResetArray(ctx->x1);
96: VecResetArray(ctx->x2);
97: VecResetArray(ctx->y1);
98: VecResetArray(ctx->y2);
99: VecRestoreArrayRead(x,&px);
100: VecRestoreArray(y,&py);
101: return(0);
102: }
106: PetscErrorCode MatGetDiagonal_Linear_N2A(Mat A,Vec diag)
107: {
109: QEP_LINEAR *ctx;
110: PetscScalar *pd;
111: PetscInt m;
112:
114: MatShellGetContext(A,(void**)&ctx);
115: MatGetLocalSize(ctx->M,&m,PETSC_NULL);
116: VecGetArray(diag,&pd);
117: VecPlaceArray(ctx->x1,pd);
118: VecPlaceArray(ctx->x2,pd+m);
119: MatGetDiagonal(ctx->K,ctx->x1);
120: VecScale(ctx->x1,-1.0);
121: VecSet(ctx->x2,1.0);
122: VecResetArray(ctx->x1);
123: VecResetArray(ctx->x2);
124: VecRestoreArray(diag,&pd);
125: return(0);
126: }
130: PetscErrorCode MatGetDiagonal_Linear_N2B(Mat B,Vec diag)
131: {
133: QEP_LINEAR *ctx;
134: PetscScalar *pd;
135: PetscInt m;
136:
138: MatShellGetContext(B,(void**)&ctx);
139: MatGetLocalSize(ctx->M,&m,PETSC_NULL);
140: VecGetArray(diag,&pd);
141: VecPlaceArray(ctx->x1,pd);
142: VecPlaceArray(ctx->x2,pd+m);
143: MatGetDiagonal(ctx->C,ctx->x1);
144: VecScale(ctx->x1,ctx->sfactor);
145: VecSet(ctx->x2,0.0);
146: VecResetArray(ctx->x1);
147: VecResetArray(ctx->x2);
148: VecRestoreArray(diag,&pd);
149: return(0);
150: }
154: PetscErrorCode MatCreateExplicit_Linear_N2A(MPI_Comm comm,QEP_LINEAR *ctx,Mat *A)
155: {
157: PetscInt M,N,m,n;
158: Mat Id;
159:
161: MatGetSize(ctx->M,&M,&N);
162: MatGetLocalSize(ctx->M,&m,&n);
163: MatCreate(((PetscObject)ctx->M)->comm,&Id);
164: MatSetSizes(Id,m,n,M,N);
165: MatSetFromOptions(Id);
166: MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
167: MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
168: MatShift(Id,1.0);
169: SlepcMatTile(-1.0,ctx->K,0.0,Id,0.0,Id,1.0,Id,A);
170: MatDestroy(&Id);
171: return(0);
172: }
176: PetscErrorCode MatCreateExplicit_Linear_N2B(MPI_Comm comm,QEP_LINEAR *ctx,Mat *B)
177: {
179: PetscInt M,N,m,n;
180: Mat Id;
181:
183: MatGetSize(ctx->M,&M,&N);
184: MatGetLocalSize(ctx->M,&m,&n);
185: MatCreate(((PetscObject)ctx->M)->comm,&Id);
186: MatSetSizes(Id,m,n,M,N);
187: MatSetFromOptions(Id);
188: MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
189: MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
190: MatShift(Id,1.0);
191: SlepcMatTile(ctx->sfactor,ctx->C,ctx->sfactor*ctx->sfactor,ctx->M,1.0,Id,0.0,Id,B);
192: MatDestroy(&Id);
193: return(0);
194: }