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
8: #include sinvert.h
12: static int MatSinvert_Mult(Mat A,Vec x,Vec y)
13: {
14: int ierr;
15: CTX_SINV *ctx;
16: PetscScalar alpha;
18: MatShellGetContext(A,(void**)&ctx);
19: alpha = -ctx->sigma;
21: if (ctx->B) { /* y = (A - sB) x */
22: MatMult(ctx->B,x,ctx->w);
23: MatMult(ctx->A,x,y);
24: VecAXPY(&alpha,ctx->w,y);
25: }
26: else { /* y = (A - sI) x */
27: MatMult(ctx->A,x,y);
28: VecAXPY(&alpha,x,y);
29: }
30: return(0);
31: }
35: static int MatSinvert_GetDiagonal(Mat A,Vec diag)
36: {
37: int ierr;
38: CTX_SINV *ctx;
39: PetscScalar alpha;
40: Vec diagb;
42: MatShellGetContext(A,(void**)&ctx);
43: alpha = -ctx->sigma;
45: MatGetDiagonal(ctx->A,diag);
46: if (ctx->B) {
47: VecDuplicate(diag,&diagb);
48: MatGetDiagonal(ctx->B,diagb);
49: VecAXPY(&alpha,diagb,diag);
50: VecDestroy(diagb);
51: }
52: else {
53: VecShift(&alpha,diag);
54: }
55: return(0);
56: }
60: static int MatSinvert_Destroy(Mat A)
61: {
62: CTX_SINV *ctx;
63: int ierr;
65: MatShellGetContext(A,(void**)&ctx);
66: if (ctx->B) { VecDestroy(ctx->w); }
67: PetscFree(ctx);
68: return(0);
69: }
73: int MatCreateMatSinvert(ST st,Mat *mat)
74: {
75: int n, m, N, M, ierr;
76: PetscTruth hasA, hasB;
77: CTX_SINV *ctx;
80: PetscNew(CTX_SINV,&ctx);
81: PetscMemzero(ctx,sizeof(CTX_SINV));
82: PetscLogObjectMemory(st,sizeof(CTX_SINV));
83: ctx->A = st->A;
84: ctx->B = st->B;
85: ctx->sigma = st->sigma;
86: if (st->B) { VecDuplicate(st->vec,&ctx->w); }
87: MatGetSize(st->A,&M,&N);
88: MatGetLocalSize(st->A,&m,&n);
89: MatCreateShell(st->comm,m,n,M,N,(void*)ctx,mat);
90: MatShellSetOperation(*mat,MATOP_MULT,(void(*)(void))MatSinvert_Mult);
91: MatShellSetOperation(*mat,MATOP_DESTROY,(void(*)(void))MatSinvert_Destroy);
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))MatSinvert_GetDiagonal);
97: }
99: return(0);
100: }