Actual source code: shellmat.c
1: /*
2: This file contains the subroutines which implement various operations
3: of the matrix associated to the shift-and-invert technique for eigenvalue
4: problems, and also a subroutine to create it.
5: */
7: #include src/st/stimpl.h
11: PetscErrorCode STMatShellMult(Mat A,Vec x,Vec y)
12: {
14: ST ctx;
17: MatShellGetContext(A,(void**)&ctx);
19: MatMult(ctx->A,x,y);
20: if (ctx->sigma != 0.0) {
21: if (ctx->B) { /* y = (A - sB) x */
22: MatMult(ctx->B,x,ctx->w);
23: VecAXPY(y,-ctx->sigma,ctx->w);
24: } else { /* y = (A - sI) x */
25: VecAXPY(y,-ctx->sigma,x);
26: }
27: }
28: return(0);
29: }
33: PetscErrorCode STMatShellMultTranspose(Mat A,Vec x,Vec y)
34: {
36: ST ctx;
39: MatShellGetContext(A,(void**)&ctx);
41: MatMultTranspose(ctx->A,x,y);
42: if (ctx->sigma != 0.0) {
43: if (ctx->B) { /* y = (A - sB) x */
44: MatMultTranspose(ctx->B,x,ctx->w);
45: VecAXPY(y,-ctx->sigma,ctx->w);
46: } else { /* y = (A - sI) x */
47: VecAXPY(y,-ctx->sigma,x);
48: }
49: }
50: return(0);
51: }
55: PetscErrorCode STMatShellGetDiagonal(Mat A,Vec diag)
56: {
58: ST ctx;
59: Vec diagb;
62: MatShellGetContext(A,(void**)&ctx);
64: MatGetDiagonal(ctx->A,diag);
65: if (ctx->sigma != 0.0) {
66: if (ctx->B) {
67: VecDuplicate(diag,&diagb);
68: MatGetDiagonal(ctx->B,diagb);
69: VecAXPY(diag,-ctx->sigma,diagb);
70: VecDestroy(diagb);
71: } else {
72: VecShift(diag,-ctx->sigma);
73: }
74: }
75: return(0);
76: }
80: PetscErrorCode STMatShellCreate(ST st,Mat *mat)
81: {
83: PetscInt n, m, N, M;
84: PetscTruth hasA, hasB;
87: MatGetSize(st->A,&M,&N);
88: MatGetLocalSize(st->A,&m,&n);
89: MatCreateShell(st->comm,m,n,M,N,(void*)st,mat);
90: MatShellSetOperation(*mat,MATOP_MULT,(void(*)(void))STMatShellMult);
91: MatShellSetOperation(*mat,MATOP_MULT_TRANSPOSE,(void(*)(void))STMatShellMultTranspose);
93: MatHasOperation(st->A,MATOP_GET_DIAGONAL,&hasA);
94: if (st->B) { MatHasOperation(st->B,MATOP_GET_DIAGONAL,&hasB); }
95: if ( (hasA && !st->B) || (hasA && hasB) ) {
96: MatShellSetOperation(*mat,MATOP_GET_DIAGONAL,(void(*)(void))STMatShellGetDiagonal);
97: }
99: return(0);
100: }