Actual source code: chebyshev.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: A polynomial filter based on a truncated Chebyshev series.
12: */
14: #include <slepc/private/stimpl.h>
15: #include "filter.h"
17: #define DAMPING_LANCZOS_POW 2
19: /*
20: Gateway to Chebyshev for evaluating y=p(A)*x
21: (Clenshaw with a vector)
22: */
23: static PetscErrorCode MatMult_Chebyshev(Mat A,Vec x,Vec y)
24: {
25: ST st;
26: ST_FILTER *ctx;
27: CHEBYSHEV_CTX cheby;
28: PetscInt degree,i,w0=0,w1=1;
29: PetscReal s1,s2;
31: PetscFunctionBegin;
32: PetscCall(MatShellGetContext(A,&st));
33: ctx = (ST_FILTER*)st->data;
34: cheby = (CHEBYSHEV_CTX)ctx->data;
35: degree = ctx->polyDegree;
36: s1 = 4/(ctx->right - ctx->left);
37: s2 = -2*(ctx->right+ctx->left)/(ctx->right-ctx->left);
39: PetscCall(VecSet(st->work[w0],0.0));
40: PetscCall(VecCopy(x,y));
41: PetscCall(VecScale(y,cheby->damping_coeffs[degree]*cheby->coeffs[degree]));
43: for (i=0;i<degree;i++) {
44: /* swap work vectors, save a copy of current y */
45: w0 = 1-w0;
46: w1 = 1-w1;
47: PetscCall(VecCopy(y,st->work[w0]));
48: /* compute next y */
49: if (i == degree-1) {
50: s1 = s1/2;
51: s2 = s2/2;
52: }
53: PetscCall(VecAXPBYPCZ(y,cheby->damping_coeffs[degree-i-1]*cheby->coeffs[degree-i-1],-1,s2,x,st->work[w1]));
54: PetscCall(MatMult(ctx->T,st->work[w0],st->work[w1]));
55: PetscCall(VecAXPY(y,s1,st->work[w1]));
56: }
57: PetscFunctionReturn(PETSC_SUCCESS);
58: }
60: /*
61: Gateway to Chebyshev for evaluating C=p(A)*B
62: (Clenshaw with a matrix)
63: */
64: static PetscErrorCode MatMatMult_Chebyshev(Mat A,Mat B,Mat C,void *pctx)
65: {
66: ST st;
67: ST_FILTER *ctx;
68: CHEBYSHEV_CTX cheby;
69: PetscInt degree,i,m1,m2,w0=0,w1=1;
70: PetscReal s1,s2;
72: PetscFunctionBegin;
73: PetscCall(MatShellGetContext(A,&st));
74: ctx = (ST_FILTER*)st->data;
75: cheby = (CHEBYSHEV_CTX)ctx->data;
76: degree = ctx->polyDegree;
77: s1 = 4/(ctx->right - ctx->left);
78: s2 = -2*(ctx->right+ctx->left)/(ctx->right-ctx->left);
80: if (ctx->nW) { /* check if work matrices must be resized */
81: PetscCall(MatGetSize(B,NULL,&m1));
82: PetscCall(MatGetSize(ctx->W[0],NULL,&m2));
83: if (m1!=m2) {
84: PetscCall(MatDestroyMatrices(ctx->nW,&ctx->W));
85: ctx->nW = 0;
86: }
87: }
88: if (!ctx->nW) { /* allocate work matrices */
89: ctx->nW = 2;
90: PetscCall(PetscMalloc1(ctx->nW,&ctx->W));
91: for (i=0;i<ctx->nW;i++) PetscCall(MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ctx->W[i]));
92: }
94: PetscCall(MatScale(ctx->W[w0],0));
95: PetscCall(MatCopy(B,C,SAME_NONZERO_PATTERN));
96: PetscCall(MatScale(C,cheby->damping_coeffs[degree]*cheby->coeffs[degree]));
98: for (i=0;i<degree;i++) {
99: /* swap work matrices, save a copy of current C */
100: w0 = 1-w0;
101: w1 = 1-w1;
102: PetscCall(MatCopy(C,ctx->W[w0],SAME_NONZERO_PATTERN));
103: if (i == degree-1) {
104: s1 = s1/2;
105: s2 = s2/2;
106: }
107: PetscCall(MatScale(C,s2));
108: PetscCall(MatAXPY(C,-1,ctx->W[w1],SAME_NONZERO_PATTERN));
109: PetscCall(MatAXPY(C,cheby->damping_coeffs[degree-i-1]*cheby->coeffs[degree-i-1],B,SAME_NONZERO_PATTERN));
110: PetscCall(MatMatMult(ctx->T,ctx->W[w0],MAT_REUSE_MATRIX,PETSC_CURRENT,&ctx->W[w1]));
111: PetscCall(MatAXPY(C,s1,ctx->W[w1],SAME_NONZERO_PATTERN));
112: }
113: PetscFunctionReturn(PETSC_SUCCESS);
114: }
116: static PetscErrorCode Chebyshev_Compute_Coefficients(ST st)
117: {
118: ST_FILTER *ctx = (ST_FILTER*)st->data;
119: CHEBYSHEV_CTX cheby = (CHEBYSHEV_CTX)ctx->data;
120: PetscInt i,degree;
121: PetscReal s_inta,s_intb,acosa,acosb,t1,t2,s;
123: PetscFunctionBegin;
124: degree = ctx->polyDegree;
125: PetscCall(PetscFree(cheby->coeffs));
126: PetscCall(PetscMalloc1(degree+1,&cheby->coeffs));
127: s_inta = PetscMax((2.0*ctx->inta - (ctx->right + ctx->left))/(ctx->right - ctx->left),-1.0);
128: s_intb = PetscMin((2.0*ctx->intb - (ctx->right + ctx->left))/(ctx->right - ctx->left),1.0);
129: acosa = PetscAcosReal(s_inta);
130: acosb = PetscAcosReal(s_intb);
131: t1 = (acosa - acosb)/PETSC_PI;
132: t2 = t1;
133: cheby->coeffs[0] = t1;
134: for (i=0;i<degree;i++) {
135: cheby->coeffs[i+1] = 2.0*(PetscSinReal((i+1)*acosa) - PetscSinReal((i+1)*acosb)) / ((i+1)*PETSC_PI);
136: t1 = t1 + cheby->coeffs[i+1]*PetscCosReal((i+1)*acosa);
137: t2 = t2 + cheby->coeffs[i+1]*PetscCosReal((i+1)*acosb);
138: }
139: if (s_inta == -1.0) s = 0.5/PetscAbs(t1);
140: else if (s_intb == 1.0) s = 0.5/PetscMin(PetscAbs(t1),PetscAbs(t2));
141: else s = 0.5/PetscMin(PetscAbs(t1),PetscAbs(t2));
142: for (i=0;i<degree+1;i++) cheby->coeffs[i] = s*cheby->coeffs[i];
143: PetscFunctionReturn(PETSC_SUCCESS);
144: }
146: static PetscErrorCode Chebyshev_Compute_Damping_Coefficients(ST st)
147: {
148: ST_FILTER *ctx = (ST_FILTER*)st->data;
149: CHEBYSHEV_CTX cheby = (CHEBYSHEV_CTX)ctx->data;
150: PetscInt i,degree;
151: PetscReal invdeg;
153: PetscFunctionBegin;
154: degree = ctx->polyDegree;
155: PetscCall(PetscFree(cheby->damping_coeffs));
156: PetscCall(PetscMalloc1(degree+1,&cheby->damping_coeffs));
157: cheby->damping_coeffs[0] = 1;
158: switch (ctx->damping) {
159: case ST_FILTER_DAMPING_NONE:
160: for (i=1;i<=degree;i++) cheby->damping_coeffs[i] = 1;
161: break;
162: case ST_FILTER_DAMPING_LANCZOS:
163: for (i=1;i<=degree;i++) cheby->damping_coeffs[i] = PetscPowRealInt(PetscSinReal(i*PETSC_PI/(degree+1))/(i*PETSC_PI/(degree+1)),DAMPING_LANCZOS_POW);
164: break;
165: case ST_FILTER_DAMPING_JACKSON:
166: invdeg = 1/((PetscReal)degree+2);
167: for (i=1;i<=degree;i++) cheby->damping_coeffs[i] = invdeg*((degree+2-i)*PetscCosReal(PETSC_PI*i*invdeg)+PetscSinReal(PETSC_PI*i*invdeg)/PetscTanReal(PETSC_PI*invdeg));
168: break;
169: case ST_FILTER_DAMPING_FEJER:
170: for (i=1;i<=degree;i++) cheby->damping_coeffs[i] = 1-i/(PetscReal)(degree+1);
171: break;
172: default:
173: SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONG,"Invalid Chebyshev filter damping type");
174: }
175: PetscFunctionReturn(PETSC_SUCCESS);
176: }
178: /*
179: Computes the coefficients and damping coefficients for the Chebyshev filter
181: Creates the shifted (and scaled) matrix and the filter P(z).
182: G is a shell matrix whose MatMult() and MatMatMult() applies the filter.
183: */
184: static PetscErrorCode STComputeOperator_Filter_Chebyshev(ST st,Mat *G)
185: {
186: ST_FILTER *ctx = (ST_FILTER*)st->data;
187: PetscInt n,m,N,M;
189: PetscFunctionBegin;
190: PetscCall(STSetWorkVecs(st,2));
191: PetscCall(Chebyshev_Compute_Coefficients(st));
192: PetscCall(Chebyshev_Compute_Damping_Coefficients(st));
193: /* create shell matrix*/
194: PetscCall(MatDestroy(&ctx->T));
195: PetscCall(MatDuplicate(st->A[0],MAT_COPY_VALUES,&ctx->T));
196: if (!*G) {
197: PetscCall(MatGetSize(ctx->T,&N,&M));
198: PetscCall(MatGetLocalSize(ctx->T,&n,&m));
199: PetscCall(MatCreateShell(PetscObjectComm((PetscObject)st),n,m,N,M,st,G));
200: PetscCall(MatShellSetOperation(*G,MATOP_MULT,(PetscErrorCodeFn*)MatMult_Chebyshev));
201: PetscCall(MatShellSetMatProductOperation(*G,MATPRODUCT_AB,NULL,MatMatMult_Chebyshev,NULL,MATDENSE,MATDENSE));
202: PetscCall(MatShellSetMatProductOperation(*G,MATPRODUCT_AB,NULL,MatMatMult_Chebyshev,NULL,MATDENSECUDA,MATDENSECUDA));
203: PetscCall(MatShellSetMatProductOperation(*G,MATPRODUCT_AB,NULL,MatMatMult_Chebyshev,NULL,MATDENSEHIP,MATDENSEHIP));
204: }
205: PetscFunctionReturn(PETSC_SUCCESS);
206: }
208: static PetscErrorCode Chebyshev_Clenshaw_Damping_Real(ST st,PetscReal x,PetscReal *px)
209: {
210: ST_FILTER *ctx = (ST_FILTER*)st->data;
211: CHEBYSHEV_CTX cheby = (CHEBYSHEV_CTX)ctx->data;
212: PetscReal t1=0,t2=0;
213: PetscInt i,degree;
215: PetscFunctionBegin;
216: degree = ctx->polyDegree;
217: PetscCheck(cheby->coeffs,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_NULL,"Chebyshev coefficients are not computed");
218: PetscCheck(cheby->damping_coeffs,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_NULL,"Chebyshev damping coefficients are not computed");
219: x = 2*x;
220: for (i=degree;i>1;i=i-2) {
221: t2 = cheby->damping_coeffs[i]*cheby->coeffs[i] + x*t1 - t2;
222: t1 = cheby->damping_coeffs[i-1]*cheby->coeffs[i-1] + x*t2 - t1;
223: }
224: if (degree % 2) {
225: t2 = cheby->damping_coeffs[1]*cheby->coeffs[1] + x*t1 - t2;
226: t1 = t2;
227: }
228: *px = cheby->damping_coeffs[0]*cheby->coeffs[0] + .5*x*t1 - t2;
229: PetscFunctionReturn(PETSC_SUCCESS);
230: }
232: /*
233: Computes the threshold for the Chebyshev filter
235: The threshold is P(scale(ctx->inta)), where ctx->inta is the lower bound of the interval.
236: It should be the same value if the upper bound was used instead.
237: */
238: static PetscErrorCode STFilterGetThreshold_Filter_Chebyshev(ST st,PetscReal *gamma)
239: {
240: ST_FILTER *ctx = (ST_FILTER*)st->data;
241: PetscReal c,e;
243: PetscFunctionBegin;
244: /* scale ctx->inta */
245: c = (ctx->right + ctx->left) / 2.0;
246: e = (ctx->right - ctx->left) / 2.0;
247: PetscCall(Chebyshev_Clenshaw_Damping_Real(st,(ctx->inta-c)/e,gamma));
248: PetscFunctionReturn(PETSC_SUCCESS);
249: }
251: static PetscErrorCode STDestroy_Filter_Chebyshev(ST st)
252: {
253: ST_FILTER *ctx = (ST_FILTER*)st->data;
254: CHEBYSHEV_CTX cheby = (CHEBYSHEV_CTX)ctx->data;
256: PetscFunctionBegin;
257: PetscCall(PetscFree(cheby->coeffs));
258: PetscCall(PetscFree(cheby->damping_coeffs));
259: PetscCall(PetscFree(cheby));
260: PetscFunctionReturn(PETSC_SUCCESS);
261: }
263: PetscErrorCode STCreate_Filter_Chebyshev(ST st)
264: {
265: ST_FILTER *ctx = (ST_FILTER*)st->data;
266: CHEBYSHEV_CTX cheby;
268: PetscFunctionBegin;
269: PetscCall(PetscNew(&cheby));
270: ctx->data = (void*)cheby;
272: ctx->computeoperator = STComputeOperator_Filter_Chebyshev;
273: ctx->getthreshold = STFilterGetThreshold_Filter_Chebyshev;
274: ctx->destroy = STDestroy_Filter_Chebyshev;
275: PetscFunctionReturn(PETSC_SUCCESS);
276: }