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-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 = [ 0 K ] z = [ x ]
34: [ C K ] [-M 0 ] [ l*x ]
35: */
39: PetscErrorCode MatMult_QEPLINEAR_H1A(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: /* 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: VecRestoreArray(x,&px);
66: VecRestoreArray(y,&py);
67: return(0);
68: }
72: PetscErrorCode MatMult_QEPLINEAR_H1B(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 = K*x2 */
89: MatMult(ctx->K,ctx->x2,ctx->y1);
90: /* y2 = -M*x1 */
91: MatMult(ctx->M,ctx->x1,ctx->y2);
92: VecScale(ctx->y2,-ctx->sfactor*ctx->sfactor);
93: VecResetArray(ctx->x1);
94: VecResetArray(ctx->x2);
95: VecResetArray(ctx->y1);
96: VecResetArray(ctx->y2);
97: VecRestoreArray(x,&px);
98: VecRestoreArray(y,&py);
99: return(0);
100: }
104: PetscErrorCode MatGetDiagonal_QEPLINEAR_H1A(Mat A,Vec diag)
105: {
107: QEP_LINEAR *ctx;
108: PetscScalar *pd;
109: PetscInt m;
110:
112: MatShellGetContext(A,(void**)&ctx);
113: MatGetLocalSize(ctx->M,&m,PETSC_NULL);
114: VecGetArray(diag,&pd);
115: VecPlaceArray(ctx->x1,pd);
116: VecPlaceArray(ctx->x2,pd+m);
117: MatGetDiagonal(ctx->K,ctx->x1);
118: VecCopy(ctx->x1,ctx->x2);
119: VecResetArray(ctx->x1);
120: VecResetArray(ctx->x2);
121: VecRestoreArray(diag,&pd);
122: return(0);
123: }
127: PetscErrorCode MatGetDiagonal_QEPLINEAR_H1B(Mat B,Vec diag)
128: {
130:
132: VecSet(diag,0.0);
133: return(0);
134: }
138: PetscErrorCode MatCreateExplicit_QEPLINEAR_H1A(MPI_Comm comm,QEP_LINEAR *ctx,Mat *A)
139: {
141: PetscInt M,N,m,n,i,j,row,start,end,ncols,*pos;
142: PetscScalar *svals;
143: const PetscInt *cols;
144: const PetscScalar *vals;
145:
147: MatGetSize(ctx->M,&M,&N);
148: MatGetLocalSize(ctx->M,&m,&n);
149: MatCreate(comm,A);
150: MatSetSizes(*A,m+n,m+n,M+N,M+N);
151: MatSetFromOptions(*A);
152: PetscMalloc(sizeof(PetscInt)*n,&pos);
153: PetscMalloc(sizeof(PetscScalar)*n,&svals);
154: MatGetOwnershipRange(ctx->M,&start,&end);
155: for (i=start;i<end;i++) {
156: row = i + M;
157: MatGetRow(ctx->K,i,&ncols,&cols,&vals);
158: MatSetValues(*A,1,&i,ncols,cols,vals,INSERT_VALUES);
159: for (j=0;j<ncols;j++)
160: pos[j] = cols[j] + M;
161: MatSetValues(*A,1,&row,ncols,pos,vals,INSERT_VALUES);
162: MatRestoreRow(ctx->K,i,&ncols,&cols,&vals);
163: MatGetRow(ctx->C,i,&ncols,&cols,&vals);
164: for (j=0;j<ncols;j++)
165: svals[j] = vals[j]*ctx->sfactor;
166: MatSetValues(*A,1,&row,ncols,cols,svals,INSERT_VALUES);
167: MatRestoreRow(ctx->C,i,&ncols,&cols,&vals);
168: }
169: PetscFree(pos);
170: PetscFree(svals);
171: MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
172: MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
173: return(0);
174: }
178: PetscErrorCode MatCreateExplicit_QEPLINEAR_H1B(MPI_Comm comm,QEP_LINEAR *ctx,Mat *B)
179: {
181: PetscInt M,N,m,n,i,j,row,start,end,ncols,*pos;
182: PetscScalar *svals;
183: const PetscInt *cols;
184: const PetscScalar *vals;
185:
187: MatGetSize(ctx->M,&M,&N);
188: MatGetLocalSize(ctx->M,&m,&n);
189: MatCreate(comm,B);
190: MatSetSizes(*B,m+n,m+n,M+N,M+N);
191: MatSetFromOptions(*B);
192: PetscMalloc(sizeof(PetscInt)*n,&pos);
193: PetscMalloc(sizeof(PetscScalar)*n,&svals);
194: MatGetOwnershipRange(ctx->M,&start,&end);
195: for (i=start;i<end;i++) {
196: row = i + M;
197: MatGetRow(ctx->M,i,&ncols,&cols,&vals);
198: for (j=0;j<ncols;j++)
199: svals[j] = vals[j]*ctx->sfactor*ctx->sfactor;
200: MatSetValues(*B,1,&row,ncols,cols,svals,INSERT_VALUES);
201: MatRestoreRow(ctx->M,i,&ncols,&cols,&vals);
202: }
203: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
204: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
205: MatScale(*B,-1.0);
206: for (i=start;i<end;i++) {
207: MatGetRow(ctx->K,i,&ncols,&cols,&vals);
208: for (j=0;j<ncols;j++)
209: pos[j] = cols[j] + M;
210: MatSetValues(*B,1,&i,ncols,pos,vals,INSERT_VALUES);
211: MatRestoreRow(ctx->K,i,&ncols,&cols,&vals);
212: }
213: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
214: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
215: PetscFree(pos);
216: PetscFree(svals);
217: return(0);
218: }