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.
6: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
7: SLEPc - Scalable Library for Eigenvalue Problem Computations
8: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
10: This file is part of SLEPc.
12: SLEPc is free software: you can redistribute it and/or modify it under the
13: terms of version 3 of the GNU Lesser General Public License as published by
14: the Free Software Foundation.
16: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
17: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
19: more details.
21: You should have received a copy of the GNU Lesser General Public License
22: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
23: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
24: */
26: #include <slepc-private/stimpl.h>
28: typedef struct {
29: PetscScalar alpha;
30: ST st;
31: Vec z;
32: PetscInt nmat;
33: PetscInt *matIdx;
34: } ST_SHELLMAT;
38: PetscErrorCode STMatShellShift(Mat A,PetscScalar alpha)
39: {
41: ST_SHELLMAT *ctx;
44: MatShellGetContext(A,(void**)&ctx);
45: ctx->alpha = alpha;
46: return(0);
47: }
51: static PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
52: {
54: ST_SHELLMAT *ctx;
55: ST st;
56: PetscInt i;
59: MatShellGetContext(A,(void**)&ctx);
60: st = ctx->st;
61: if (ctx->alpha != 0.0) {
62: MatMult(st->A[ctx->matIdx[ctx->nmat-1]],x,y);
63: if (ctx->nmat>1) { /* */
64: for (i=ctx->nmat-2;i>=0;i--) {
65: MatMult(st->A[ctx->matIdx[i]],x,ctx->z);
66: VecAYPX(y,ctx->alpha,ctx->z);
67: }
68: } else { /* y = (A + alpha*I) x */
69: VecAXPY(y,ctx->alpha,x);
70: }
71: } else {
72: MatMult(st->A[ctx->matIdx[0]],x,y);
73: }
74: return(0);
75: }
79: static PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
80: {
82: ST_SHELLMAT *ctx;
83: ST st;
84: PetscInt i;
87: MatShellGetContext(A,(void**)&ctx);
88: st = ctx->st;
89: if (ctx->alpha != 0.0) {
90: MatMultTranspose(st->A[ctx->matIdx[ctx->nmat-1]],x,y);
91: if (st->nmat>1) { /* y = (A + alpha*B) x */
92: for (i=ctx->nmat-2;i>=0;i--) {
93: MatMultTranspose(st->A[ctx->matIdx[i]],x,y);
94: VecAYPX(y,ctx->alpha,ctx->z);
95: }
96: } else { /* y = (A + alpha*I) x */
97: VecAXPY(y,ctx->alpha,x);
98: }
99: } else {
100: MatMultTranspose(st->A[ctx->matIdx[0]],x,y);
101: }
102: return(0);
103: }
107: static PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec diag)
108: {
110: ST_SHELLMAT *ctx;
111: ST st;
112: Vec diagb;
113: PetscInt i;
116: MatShellGetContext(A,(void**)&ctx);
117: st = ctx->st;
118: if (ctx->alpha != 0.0) {
119: MatGetDiagonal(st->A[ctx->matIdx[ctx->nmat-1]],diag);
120: if (st->nmat>1) {
121: VecDuplicate(diag,&diagb);
122: for (i=ctx->nmat-2;i>=0;i--) {
123: MatGetDiagonal(st->A[ctx->matIdx[i]],diagb);
124: VecAYPX(diag,ctx->alpha,diagb);
125: }
126: VecDestroy(&diagb);
127: } else {
128: VecShift(diag,ctx->alpha);
129: }
130: } else {
131: MatGetDiagonal(st->A[ctx->matIdx[0]],diag);
132: }
133: return(0);
134: }
138: static PetscErrorCode MatDestroy_Shell(Mat A)
139: {
141: ST_SHELLMAT *ctx;
144: MatShellGetContext(A,(void**)&ctx);
145: VecDestroy(&ctx->z);
146: PetscFree(ctx->matIdx);
147: PetscFree(ctx);
148: return(0);
149: }
153: PetscErrorCode STMatShellCreate(ST st,PetscScalar alpha,PetscInt nmat,PetscInt *matIdx,Mat *mat)
154: {
156: PetscInt n,m,N,M,i;
157: PetscBool has=PETSC_FALSE,hasA,hasB;
158: ST_SHELLMAT *ctx;
161: MatGetSize(st->A[0],&M,&N);
162: MatGetLocalSize(st->A[0],&m,&n);
163: PetscNew(ST_SHELLMAT,&ctx);
164: ctx->st = st;
165: ctx->alpha = alpha;
166: ctx->nmat = matIdx?nmat:st->nmat;
167: PetscMalloc(ctx->nmat*sizeof(PetscInt),&ctx->matIdx);
168: if (matIdx) {
169: for (i=0;i<ctx->nmat;i++) ctx->matIdx[i] = matIdx[i];
170: } else {
171: ctx->matIdx[0] = 0;
172: if (ctx->nmat>1) ctx->matIdx[1] = 1;
173: }
174: MatGetVecs(st->A[0],&ctx->z,NULL);
175: MatCreateShell(PetscObjectComm((PetscObject)st),m,n,M,N,(void*)ctx,mat);
176: MatShellSetOperation(*mat,MATOP_MULT,(void(*)(void))MatMult_Shell);
177: MatShellSetOperation(*mat,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell);
178: MatShellSetOperation(*mat,MATOP_DESTROY,(void(*)(void))MatDestroy_Shell);
180: MatHasOperation(st->A[ctx->matIdx[0]],MATOP_GET_DIAGONAL,&hasA);
181: if (st->nmat>1) {
182: has = hasA;
183: for (i=1;i<ctx->nmat;i++) {
184: MatHasOperation(st->A[ctx->matIdx[i]],MATOP_GET_DIAGONAL,&hasB);
185: has = (has && hasB)? PETSC_TRUE: PETSC_FALSE;
186: }
187: }
188: if ((hasA && st->nmat==1) || has) {
189: MatShellSetOperation(*mat,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Shell);
190: }
191: return(0);
192: }