Actual source code: fnrational.c
slepc-3.18.1 2022-11-02
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: Rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials
12: */
14: #include <slepc/private/fnimpl.h>
16: typedef struct {
17: PetscScalar *pcoeff; /* numerator coefficients */
18: PetscInt np; /* length of array pcoeff, p(x) has degree np-1 */
19: PetscScalar *qcoeff; /* denominator coefficients */
20: PetscInt nq; /* length of array qcoeff, q(x) has degree nq-1 */
21: } FN_RATIONAL;
23: PetscErrorCode FNEvaluateFunction_Rational(FN fn,PetscScalar x,PetscScalar *y)
24: {
25: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
26: PetscInt i;
27: PetscScalar p,q;
29: if (!ctx->np) p = 1.0;
30: else {
31: p = ctx->pcoeff[0];
32: for (i=1;i<ctx->np;i++)
33: p = ctx->pcoeff[i]+x*p;
34: }
35: if (!ctx->nq) *y = p;
36: else {
37: q = ctx->qcoeff[0];
38: for (i=1;i<ctx->nq;i++)
39: q = ctx->qcoeff[i]+x*q;
41: *y = p/q;
42: }
43: return 0;
44: }
46: /*
47: Horner evaluation of P=p(A)
48: d = degree of polynomial; coeff = coefficients of polynomial; W = workspace
49: */
50: static PetscErrorCode EvaluatePoly(Mat A,Mat P,Mat W,PetscInt d,PetscScalar *coeff)
51: {
52: PetscInt j;
54: MatZeroEntries(P);
55: if (!d) MatShift(P,1.0);
56: else {
57: MatShift(P,coeff[0]);
58: for (j=1;j<d;j++) {
59: MatMatMult(P,A,MAT_REUSE_MATRIX,PETSC_DEFAULT,&W);
60: MatCopy(W,P,SAME_NONZERO_PATTERN);
61: MatShift(P,coeff[j]);
62: }
63: }
64: return 0;
65: }
67: PetscErrorCode FNEvaluateFunctionMat_Rational(FN fn,Mat A,Mat B)
68: {
69: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
70: Mat P,Q,W,F;
71: PetscBool iscuda;
73: if (A==B) MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&P);
74: else P = B;
75: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&W);
77: EvaluatePoly(A,P,W,ctx->np,ctx->pcoeff);
78: if (ctx->nq) {
79: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q);
80: EvaluatePoly(A,Q,W,ctx->nq,ctx->qcoeff);
81: PetscObjectTypeCompare((PetscObject)A,MATSEQDENSECUDA,&iscuda);
82: MatGetFactor(Q,iscuda?MATSOLVERCUDA:MATSOLVERPETSC,MAT_FACTOR_LU,&F);
83: MatLUFactorSymbolic(F,Q,NULL,NULL,NULL);
84: MatLUFactorNumeric(F,Q,NULL);
85: MatMatSolve(F,P,P);
86: MatDestroy(&F);
87: MatDestroy(&Q);
88: }
90: if (A==B) {
91: MatCopy(P,B,SAME_NONZERO_PATTERN);
92: MatDestroy(&P);
93: }
94: MatDestroy(&W);
95: return 0;
96: }
98: PetscErrorCode FNEvaluateFunctionMatVec_Rational(FN fn,Mat A,Vec v)
99: {
100: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
101: Mat P,Q,W,F;
102: Vec b;
103: PetscBool iscuda;
105: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&P);
106: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&W);
108: EvaluatePoly(A,P,W,ctx->np,ctx->pcoeff);
109: if (ctx->nq) {
110: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q);
111: EvaluatePoly(A,Q,W,ctx->nq,ctx->qcoeff);
112: PetscObjectTypeCompare((PetscObject)A,MATSEQDENSECUDA,&iscuda);
113: MatGetFactor(Q,iscuda?MATSOLVERCUDA:MATSOLVERPETSC,MAT_FACTOR_LU,&F);
114: MatLUFactorSymbolic(F,Q,NULL,NULL,NULL);
115: MatLUFactorNumeric(F,Q,NULL);
116: MatCreateVecs(P,&b,NULL);
117: MatGetColumnVector(P,b,0);
118: MatSolve(F,b,v);
119: VecDestroy(&b);
120: MatDestroy(&F);
121: MatDestroy(&Q);
122: } else MatGetColumnVector(P,v,0);
124: MatDestroy(&P);
125: MatDestroy(&W);
126: return 0;
127: }
129: PetscErrorCode FNEvaluateDerivative_Rational(FN fn,PetscScalar x,PetscScalar *yp)
130: {
131: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
132: PetscInt i;
133: PetscScalar p,q,pp,qp;
135: if (!ctx->np) {
136: p = 1.0;
137: pp = 0.0;
138: } else {
139: p = ctx->pcoeff[0];
140: pp = 0.0;
141: for (i=1;i<ctx->np;i++) {
142: pp = p+x*pp;
143: p = ctx->pcoeff[i]+x*p;
144: }
145: }
146: if (!ctx->nq) *yp = pp;
147: else {
148: q = ctx->qcoeff[0];
149: qp = 0.0;
150: for (i=1;i<ctx->nq;i++) {
151: qp = q+x*qp;
152: q = ctx->qcoeff[i]+x*q;
153: }
155: *yp = (pp*q-p*qp)/(q*q);
156: }
157: return 0;
158: }
160: PetscErrorCode FNView_Rational(FN fn,PetscViewer viewer)
161: {
162: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
163: PetscBool isascii;
164: PetscInt i;
165: char str[50];
167: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
168: if (isascii) {
169: if (fn->alpha!=(PetscScalar)1.0 || fn->beta!=(PetscScalar)1.0) {
170: SlepcSNPrintfScalar(str,sizeof(str),fn->alpha,PETSC_FALSE);
171: PetscViewerASCIIPrintf(viewer," scale factors: alpha=%s,",str);
172: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
173: SlepcSNPrintfScalar(str,sizeof(str),fn->beta,PETSC_FALSE);
174: PetscViewerASCIIPrintf(viewer," beta=%s\n",str);
175: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
176: }
177: if (!ctx->nq) {
178: if (!ctx->np) PetscViewerASCIIPrintf(viewer," constant: 1.0\n");
179: else if (ctx->np==1) {
180: SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[0],PETSC_FALSE);
181: PetscViewerASCIIPrintf(viewer," constant: %s\n",str);
182: } else {
183: PetscViewerASCIIPrintf(viewer," polynomial: ");
184: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
185: for (i=0;i<ctx->np-1;i++) {
186: SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[i],PETSC_TRUE);
187: PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->np-i-1);
188: }
189: SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[ctx->np-1],PETSC_TRUE);
190: PetscViewerASCIIPrintf(viewer,"%s\n",str);
191: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
192: }
193: } else if (!ctx->np) {
194: PetscViewerASCIIPrintf(viewer," inverse polynomial: 1 / (");
195: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
196: for (i=0;i<ctx->nq-1;i++) {
197: SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[i],PETSC_TRUE);
198: PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->nq-i-1);
199: }
200: SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[ctx->nq-1],PETSC_TRUE);
201: PetscViewerASCIIPrintf(viewer,"%s)\n",str);
202: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
203: } else {
204: PetscViewerASCIIPrintf(viewer," rational function: (");
205: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
206: for (i=0;i<ctx->np-1;i++) {
207: SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[i],PETSC_TRUE);
208: PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->np-i-1);
209: }
210: SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[ctx->np-1],PETSC_TRUE);
211: PetscViewerASCIIPrintf(viewer,"%s) / (",str);
212: for (i=0;i<ctx->nq-1;i++) {
213: SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[i],PETSC_TRUE);
214: PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->nq-i-1);
215: }
216: SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[ctx->nq-1],PETSC_TRUE);
217: PetscViewerASCIIPrintf(viewer,"%s)\n",str);
218: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
219: }
220: }
221: return 0;
222: }
224: static PetscErrorCode FNRationalSetNumerator_Rational(FN fn,PetscInt np,PetscScalar *pcoeff)
225: {
226: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
227: PetscInt i;
230: ctx->np = np;
231: PetscFree(ctx->pcoeff);
232: if (np) {
233: PetscMalloc1(np,&ctx->pcoeff);
234: for (i=0;i<np;i++) ctx->pcoeff[i] = pcoeff[i];
235: }
236: return 0;
237: }
239: /*@C
240: FNRationalSetNumerator - Sets the parameters defining the numerator of the
241: rational function.
243: Logically Collective on fn
245: Input Parameters:
246: + fn - the math function context
247: . np - number of coefficients
248: - pcoeff - coefficients (array of scalar values)
250: Notes:
251: Let the rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials.
252: This function provides the coefficients of the numerator p(x).
253: Hence, p(x) is of degree np-1.
254: If np is zero, then the numerator is assumed to be p(x)=1.
256: In polynomials, high order coefficients are stored in the first positions
257: of the array, e.g. to represent x^2-3 use {1,0,-3}.
259: Level: intermediate
261: .seealso: FNRationalSetDenominator(), FNRationalGetNumerator()
262: @*/
263: PetscErrorCode FNRationalSetNumerator(FN fn,PetscInt np,PetscScalar *pcoeff)
264: {
268: PetscTryMethod(fn,"FNRationalSetNumerator_C",(FN,PetscInt,PetscScalar*),(fn,np,pcoeff));
269: return 0;
270: }
272: static PetscErrorCode FNRationalGetNumerator_Rational(FN fn,PetscInt *np,PetscScalar *pcoeff[])
273: {
274: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
275: PetscInt i;
277: if (np) *np = ctx->np;
278: if (pcoeff) {
279: if (!ctx->np) *pcoeff = NULL;
280: else {
281: PetscMalloc1(ctx->np,pcoeff);
282: for (i=0;i<ctx->np;i++) (*pcoeff)[i] = ctx->pcoeff[i];
283: }
284: }
285: return 0;
286: }
288: /*@C
289: FNRationalGetNumerator - Gets the parameters that define the numerator of the
290: rational function.
292: Not Collective
294: Input Parameter:
295: . fn - the math function context
297: Output Parameters:
298: + np - number of coefficients
299: - pcoeff - coefficients (array of scalar values, length nq)
301: Notes:
302: The values passed by user with FNRationalSetNumerator() are returned (or null
303: pointers otherwise).
304: The pcoeff array should be freed by the user when no longer needed.
306: Level: intermediate
308: .seealso: FNRationalSetNumerator()
309: @*/
310: PetscErrorCode FNRationalGetNumerator(FN fn,PetscInt *np,PetscScalar *pcoeff[])
311: {
313: PetscUseMethod(fn,"FNRationalGetNumerator_C",(FN,PetscInt*,PetscScalar**),(fn,np,pcoeff));
314: return 0;
315: }
317: static PetscErrorCode FNRationalSetDenominator_Rational(FN fn,PetscInt nq,PetscScalar *qcoeff)
318: {
319: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
320: PetscInt i;
323: ctx->nq = nq;
324: PetscFree(ctx->qcoeff);
325: if (nq) {
326: PetscMalloc1(nq,&ctx->qcoeff);
327: for (i=0;i<nq;i++) ctx->qcoeff[i] = qcoeff[i];
328: }
329: return 0;
330: }
332: /*@C
333: FNRationalSetDenominator - Sets the parameters defining the denominator of the
334: rational function.
336: Logically Collective on fn
338: Input Parameters:
339: + fn - the math function context
340: . nq - number of coefficients
341: - qcoeff - coefficients (array of scalar values)
343: Notes:
344: Let the rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials.
345: This function provides the coefficients of the denominator q(x).
346: Hence, q(x) is of degree nq-1.
347: If nq is zero, then the function is assumed to be polynomial, r(x) = p(x).
349: In polynomials, high order coefficients are stored in the first positions
350: of the array, e.g. to represent x^2-3 use {1,0,-3}.
352: Level: intermediate
354: .seealso: FNRationalSetNumerator(), FNRationalGetDenominator()
355: @*/
356: PetscErrorCode FNRationalSetDenominator(FN fn,PetscInt nq,PetscScalar *qcoeff)
357: {
361: PetscTryMethod(fn,"FNRationalSetDenominator_C",(FN,PetscInt,PetscScalar*),(fn,nq,qcoeff));
362: return 0;
363: }
365: static PetscErrorCode FNRationalGetDenominator_Rational(FN fn,PetscInt *nq,PetscScalar *qcoeff[])
366: {
367: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
368: PetscInt i;
370: if (nq) *nq = ctx->nq;
371: if (qcoeff) {
372: if (!ctx->nq) *qcoeff = NULL;
373: else {
374: PetscMalloc1(ctx->nq,qcoeff);
375: for (i=0;i<ctx->nq;i++) (*qcoeff)[i] = ctx->qcoeff[i];
376: }
377: }
378: return 0;
379: }
381: /*@C
382: FNRationalGetDenominator - Gets the parameters that define the denominator of the
383: rational function.
385: Not Collective
387: Input Parameter:
388: . fn - the math function context
390: Output Parameters:
391: + nq - number of coefficients
392: - qcoeff - coefficients (array of scalar values, length nq)
394: Notes:
395: The values passed by user with FNRationalSetDenominator() are returned (or a null
396: pointer otherwise).
397: The qcoeff array should be freed by the user when no longer needed.
399: Level: intermediate
401: .seealso: FNRationalSetDenominator()
402: @*/
403: PetscErrorCode FNRationalGetDenominator(FN fn,PetscInt *nq,PetscScalar *qcoeff[])
404: {
406: PetscUseMethod(fn,"FNRationalGetDenominator_C",(FN,PetscInt*,PetscScalar**),(fn,nq,qcoeff));
407: return 0;
408: }
410: PetscErrorCode FNSetFromOptions_Rational(FN fn,PetscOptionItems *PetscOptionsObject)
411: {
412: #define PARMAX 10
413: PetscScalar array[PARMAX];
414: PetscInt i,k;
415: PetscBool flg;
417: PetscOptionsHeadBegin(PetscOptionsObject,"FN Rational Options");
419: k = PARMAX;
420: for (i=0;i<k;i++) array[i] = 0;
421: PetscOptionsScalarArray("-fn_rational_numerator","Numerator coefficients (one or more scalar values separated with a comma without spaces)","FNRationalSetNumerator",array,&k,&flg);
422: if (flg) FNRationalSetNumerator(fn,k,array);
424: k = PARMAX;
425: for (i=0;i<k;i++) array[i] = 0;
426: PetscOptionsScalarArray("-fn_rational_denominator","Denominator coefficients (one or more scalar values separated with a comma without spaces)","FNRationalSetDenominator",array,&k,&flg);
427: if (flg) FNRationalSetDenominator(fn,k,array);
429: PetscOptionsHeadEnd();
430: return 0;
431: }
433: PetscErrorCode FNDuplicate_Rational(FN fn,MPI_Comm comm,FN *newfn)
434: {
435: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data,*ctx2 = (FN_RATIONAL*)(*newfn)->data;
436: PetscInt i;
438: ctx2->np = ctx->np;
439: if (ctx->np) {
440: PetscMalloc1(ctx->np,&ctx2->pcoeff);
441: for (i=0;i<ctx->np;i++) ctx2->pcoeff[i] = ctx->pcoeff[i];
442: }
443: ctx2->nq = ctx->nq;
444: if (ctx->nq) {
445: PetscMalloc1(ctx->nq,&ctx2->qcoeff);
446: for (i=0;i<ctx->nq;i++) ctx2->qcoeff[i] = ctx->qcoeff[i];
447: }
448: return 0;
449: }
451: PetscErrorCode FNDestroy_Rational(FN fn)
452: {
453: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
455: PetscFree(ctx->pcoeff);
456: PetscFree(ctx->qcoeff);
457: PetscFree(fn->data);
458: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetNumerator_C",NULL);
459: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetNumerator_C",NULL);
460: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetDenominator_C",NULL);
461: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetDenominator_C",NULL);
462: return 0;
463: }
465: SLEPC_EXTERN PetscErrorCode FNCreate_Rational(FN fn)
466: {
467: FN_RATIONAL *ctx;
469: PetscNew(&ctx);
470: fn->data = (void*)ctx;
472: fn->ops->evaluatefunction = FNEvaluateFunction_Rational;
473: fn->ops->evaluatederivative = FNEvaluateDerivative_Rational;
474: fn->ops->evaluatefunctionmat[0] = FNEvaluateFunctionMat_Rational;
475: fn->ops->evaluatefunctionmatvec[0] = FNEvaluateFunctionMatVec_Rational;
476: #if defined(PETSC_HAVE_CUDA)
477: fn->ops->evaluatefunctionmatcuda[0] = FNEvaluateFunctionMat_Rational;
478: fn->ops->evaluatefunctionmatveccuda[0] = FNEvaluateFunctionMatVec_Rational;
479: #endif
480: fn->ops->setfromoptions = FNSetFromOptions_Rational;
481: fn->ops->view = FNView_Rational;
482: fn->ops->duplicate = FNDuplicate_Rational;
483: fn->ops->destroy = FNDestroy_Rational;
484: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetNumerator_C",FNRationalSetNumerator_Rational);
485: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetNumerator_C",FNRationalGetNumerator_Rational);
486: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetDenominator_C",FNRationalSetDenominator_Rational);
487: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetDenominator_C",FNRationalGetDenominator_Rational);
488: return 0;
489: }