Actual source code: stshellmat.c
slepc-3.7.1 2016-05-27
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-2016, 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: PetscScalar *coeffs;
31: ST st;
32: Vec z;
33: PetscInt nmat;
34: PetscInt *matIdx;
35: } ST_SHELLMAT;
39: PetscErrorCode STMatShellShift(Mat A,PetscScalar alpha)
40: {
42: ST_SHELLMAT *ctx;
45: MatShellGetContext(A,(void**)&ctx);
46: ctx->alpha = alpha;
47: return(0);
48: }
52: /*
53: For i=0:nmat-1 computes y = (sum_i (coeffs[i]*alpha^i*st->A[idx[i]]))x
54: If null coeffs computes with coeffs[i]=1.0
55: */
56: static PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
57: {
59: ST_SHELLMAT *ctx;
60: ST st;
61: PetscInt i;
62: PetscScalar t=1.0,c;
65: MatShellGetContext(A,(void**)&ctx);
66: st = ctx->st;
67: MatMult(st->A[ctx->matIdx[0]],x,y);
68: if (ctx->coeffs && ctx->coeffs[0]!=1.0) {
69: VecScale(y,ctx->coeffs[0]);
70: }
71: if (ctx->alpha!=0.0) {
72: for (i=1;i<ctx->nmat;i++) {
73: MatMult(st->A[ctx->matIdx[i]],x,ctx->z);
74: t *= ctx->alpha;
75: c = (ctx->coeffs)?t*ctx->coeffs[i]:t;
76: VecAXPY(y,c,ctx->z);
77: }
78: if (ctx->nmat==1) { /* y = (A + alpha*I) x */
79: VecAXPY(y,ctx->alpha,x);
80: }
81: }
82: return(0);
83: }
87: static PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
88: {
90: ST_SHELLMAT *ctx;
91: ST st;
92: PetscInt i;
93: PetscScalar t=1.0,c;
96: MatShellGetContext(A,(void**)&ctx);
97: st = ctx->st;
98: MatMultTranspose(st->A[ctx->matIdx[0]],x,y);
99: if (ctx->coeffs && ctx->coeffs[0]!=1.0) {
100: VecScale(y,ctx->coeffs[0]);
101: }
102: if (ctx->alpha!=0.0) {
103: for (i=1;i<ctx->nmat;i++) {
104: MatMultTranspose(st->A[ctx->matIdx[i]],x,ctx->z);
105: t *= ctx->alpha;
106: c = (ctx->coeffs)?t*ctx->coeffs[i]:t;
107: VecAXPY(y,c,ctx->z);
108: }
109: if (ctx->nmat==1) { /* y = (A + alpha*I) x */
110: VecAXPY(y,ctx->alpha,x);
111: }
112: }
113: return(0);
114: }
118: static PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec diag)
119: {
121: ST_SHELLMAT *ctx;
122: ST st;
123: Vec diagb;
124: PetscInt i;
125: PetscScalar t=1.0,c;
128: MatShellGetContext(A,(void**)&ctx);
129: st = ctx->st;
130: MatGetDiagonal(st->A[ctx->matIdx[0]],diag);
131: if (ctx->coeffs && ctx->coeffs[0]!=1.0) {
132: VecScale(diag,ctx->coeffs[0]);
133: }
134: if (ctx->alpha!=0.0) {
135: if (ctx->nmat==1) { /* y = (A + alpha*I) x */
136: VecShift(diag,ctx->alpha);
137: } else {
138: VecDuplicate(diag,&diagb);
139: for (i=1;i<ctx->nmat;i++) {
140: MatGetDiagonal(st->A[ctx->matIdx[i]],diagb);
141: t *= ctx->alpha;
142: c = (ctx->coeffs)?t*ctx->coeffs[i]:t;
143: VecAYPX(diag,c,diagb);
144: }
145: VecDestroy(&diagb);
146: }
147: }
148: return(0);
149: }
153: static PetscErrorCode MatDestroy_Shell(Mat A)
154: {
156: ST_SHELLMAT *ctx;
159: MatShellGetContext(A,(void**)&ctx);
160: VecDestroy(&ctx->z);
161: PetscFree(ctx->matIdx);
162: PetscFree(ctx->coeffs);
163: PetscFree(ctx);
164: return(0);
165: }
169: PetscErrorCode STMatShellCreate(ST st,PetscScalar alpha,PetscInt nmat,PetscInt *matIdx,PetscScalar *coeffs,Mat *mat)
170: {
172: PetscInt n,m,N,M,i;
173: PetscBool has=PETSC_FALSE,hasA,hasB;
174: ST_SHELLMAT *ctx;
177: MatGetSize(st->A[0],&M,&N);
178: MatGetLocalSize(st->A[0],&m,&n);
179: PetscNew(&ctx);
180: ctx->st = st;
181: ctx->alpha = alpha;
182: ctx->nmat = matIdx?nmat:st->nmat;
183: PetscMalloc1(ctx->nmat,&ctx->matIdx);
184: if (matIdx) {
185: for (i=0;i<ctx->nmat;i++) ctx->matIdx[i] = matIdx[i];
186: } else {
187: ctx->matIdx[0] = 0;
188: if (ctx->nmat>1) ctx->matIdx[1] = 1;
189: }
190: if (coeffs) {
191: PetscMalloc(ctx->nmat*sizeof(PetscScalar),&ctx->coeffs);
192: for (i=0;i<ctx->nmat;i++) ctx->coeffs[i] = coeffs[i];
193: }
194: MatCreateVecs(st->A[0],&ctx->z,NULL);
195: MatCreateShell(PetscObjectComm((PetscObject)st),m,n,M,N,(void*)ctx,mat);
196: MatShellSetOperation(*mat,MATOP_MULT,(void(*)(void))MatMult_Shell);
197: MatShellSetOperation(*mat,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell);
198: MatShellSetOperation(*mat,MATOP_DESTROY,(void(*)(void))MatDestroy_Shell);
200: MatHasOperation(st->A[ctx->matIdx[0]],MATOP_GET_DIAGONAL,&hasA);
201: if (st->nmat>1) {
202: has = hasA;
203: for (i=1;i<ctx->nmat;i++) {
204: MatHasOperation(st->A[ctx->matIdx[i]],MATOP_GET_DIAGONAL,&hasB);
205: has = (has && hasB)? PETSC_TRUE: PETSC_FALSE;
206: }
207: }
208: if ((hasA && st->nmat==1) || has) {
209: MatShellSetOperation(*mat,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Shell);
210: }
211: return(0);
212: }