Actual source code: fnrational.c
slepc-3.17.1 2022-04-11
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>
15: #include <slepcblaslapack.h>
17: typedef struct {
18: PetscScalar *pcoeff; /* numerator coefficients */
19: PetscInt np; /* length of array pcoeff, p(x) has degree np-1 */
20: PetscScalar *qcoeff; /* denominator coefficients */
21: PetscInt nq; /* length of array qcoeff, q(x) has degree nq-1 */
22: } FN_RATIONAL;
24: PetscErrorCode FNEvaluateFunction_Rational(FN fn,PetscScalar x,PetscScalar *y)
25: {
26: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
27: PetscInt i;
28: PetscScalar p,q;
30: if (!ctx->np) p = 1.0;
31: else {
32: p = ctx->pcoeff[0];
33: for (i=1;i<ctx->np;i++)
34: p = ctx->pcoeff[i]+x*p;
35: }
36: if (!ctx->nq) *y = p;
37: else {
38: q = ctx->qcoeff[0];
39: for (i=1;i<ctx->nq;i++)
40: q = ctx->qcoeff[i]+x*q;
42: *y = p/q;
43: }
44: PetscFunctionReturn(0);
45: }
47: static PetscErrorCode FNEvaluateFunctionMat_Rational_Private(FN fn,const PetscScalar *Aa,PetscScalar *Ba,PetscInt m,PetscBool firstonly)
48: {
49: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
50: PetscBLASInt n,k,ld,*ipiv,info;
51: PetscInt i,j;
52: PetscScalar *W,*P,*Q,one=1.0,zero=0.0;
54: PetscBLASIntCast(m,&n);
55: ld = n;
56: k = firstonly? 1: n;
57: if (Aa==Ba) PetscMalloc4(m*m,&P,m*m,&Q,m*m,&W,ld,&ipiv);
58: else {
59: P = Ba;
60: PetscMalloc3(m*m,&Q,m*m,&W,ld,&ipiv);
61: }
62: PetscArrayzero(P,m*m);
63: if (!ctx->np) {
64: for (i=0;i<m;i++) P[i+i*ld] = 1.0;
65: } else {
66: for (i=0;i<m;i++) P[i+i*ld] = ctx->pcoeff[0];
67: for (j=1;j<ctx->np;j++) {
68: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&one,P,&ld,Aa,&ld,&zero,W,&ld));
69: PetscArraycpy(P,W,m*m);
70: for (i=0;i<m;i++) P[i+i*ld] += ctx->pcoeff[j];
71: }
72: PetscLogFlops(2.0*n*n*n*(ctx->np-1));
73: }
74: if (ctx->nq) {
75: PetscArrayzero(Q,m*m);
76: for (i=0;i<m;i++) Q[i+i*ld] = ctx->qcoeff[0];
77: for (j=1;j<ctx->nq;j++) {
78: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&one,Q,&ld,Aa,&ld,&zero,W,&ld));
79: PetscArraycpy(Q,W,m*m);
80: for (i=0;i<m;i++) Q[i+i*ld] += ctx->qcoeff[j];
81: }
82: PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&n,&k,Q,&ld,ipiv,P,&ld,&info));
83: SlepcCheckLapackInfo("gesv",info);
84: PetscLogFlops(2.0*n*n*n*(ctx->nq-1)+2.0*n*n*n/3.0+2.0*n*n*k);
85: }
86: if (Aa==Ba) {
87: PetscArraycpy(Ba,P,m*k);
88: PetscFree4(P,Q,W,ipiv);
89: } else PetscFree3(Q,W,ipiv);
90: PetscFunctionReturn(0);
91: }
93: PetscErrorCode FNEvaluateFunctionMat_Rational(FN fn,Mat A,Mat B)
94: {
95: PetscInt m;
96: const PetscScalar *Aa;
97: PetscScalar *Ba;
99: MatDenseGetArrayRead(A,&Aa);
100: MatDenseGetArray(B,&Ba);
101: MatGetSize(A,&m,NULL);
102: FNEvaluateFunctionMat_Rational_Private(fn,Aa,Ba,m,PETSC_FALSE);
103: MatDenseRestoreArrayRead(A,&Aa);
104: MatDenseRestoreArray(B,&Ba);
105: PetscFunctionReturn(0);
106: }
108: PetscErrorCode FNEvaluateFunctionMatVec_Rational(FN fn,Mat A,Vec v)
109: {
110: PetscInt m;
111: const PetscScalar *Aa;
112: PetscScalar *Ba;
113: Mat B;
115: FN_AllocateWorkMat(fn,A,&B);
116: MatDenseGetArrayRead(A,&Aa);
117: MatDenseGetArray(B,&Ba);
118: MatGetSize(A,&m,NULL);
119: FNEvaluateFunctionMat_Rational_Private(fn,Aa,Ba,m,PETSC_TRUE);
120: MatDenseRestoreArrayRead(A,&Aa);
121: MatDenseRestoreArray(B,&Ba);
122: MatGetColumnVector(B,v,0);
123: FN_FreeWorkMat(fn,&B);
124: PetscFunctionReturn(0);
125: }
127: PetscErrorCode FNEvaluateDerivative_Rational(FN fn,PetscScalar x,PetscScalar *yp)
128: {
129: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
130: PetscInt i;
131: PetscScalar p,q,pp,qp;
133: if (!ctx->np) {
134: p = 1.0;
135: pp = 0.0;
136: } else {
137: p = ctx->pcoeff[0];
138: pp = 0.0;
139: for (i=1;i<ctx->np;i++) {
140: pp = p+x*pp;
141: p = ctx->pcoeff[i]+x*p;
142: }
143: }
144: if (!ctx->nq) *yp = pp;
145: else {
146: q = ctx->qcoeff[0];
147: qp = 0.0;
148: for (i=1;i<ctx->nq;i++) {
149: qp = q+x*qp;
150: q = ctx->qcoeff[i]+x*q;
151: }
153: *yp = (pp*q-p*qp)/(q*q);
154: }
155: PetscFunctionReturn(0);
156: }
158: PetscErrorCode FNView_Rational(FN fn,PetscViewer viewer)
159: {
160: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
161: PetscBool isascii;
162: PetscInt i;
163: char str[50];
165: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
166: if (isascii) {
167: if (fn->alpha!=(PetscScalar)1.0 || fn->beta!=(PetscScalar)1.0) {
168: SlepcSNPrintfScalar(str,sizeof(str),fn->alpha,PETSC_FALSE);
169: PetscViewerASCIIPrintf(viewer," Scale factors: alpha=%s,",str);
170: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
171: SlepcSNPrintfScalar(str,sizeof(str),fn->beta,PETSC_FALSE);
172: PetscViewerASCIIPrintf(viewer," beta=%s\n",str);
173: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
174: }
175: if (!ctx->nq) {
176: if (!ctx->np) PetscViewerASCIIPrintf(viewer," Constant: 1.0\n");
177: else if (ctx->np==1) {
178: SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[0],PETSC_FALSE);
179: PetscViewerASCIIPrintf(viewer," Constant: %s\n",str);
180: } else {
181: PetscViewerASCIIPrintf(viewer," Polynomial: ");
182: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
183: for (i=0;i<ctx->np-1;i++) {
184: SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[i],PETSC_TRUE);
185: PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->np-i-1);
186: }
187: SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[ctx->np-1],PETSC_TRUE);
188: PetscViewerASCIIPrintf(viewer,"%s\n",str);
189: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
190: }
191: } else if (!ctx->np) {
192: PetscViewerASCIIPrintf(viewer," Inverse polinomial: 1 / (");
193: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
194: for (i=0;i<ctx->nq-1;i++) {
195: SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[i],PETSC_TRUE);
196: PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->nq-i-1);
197: }
198: SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[ctx->nq-1],PETSC_TRUE);
199: PetscViewerASCIIPrintf(viewer,"%s)\n",str);
200: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
201: } else {
202: PetscViewerASCIIPrintf(viewer," Rational function: (");
203: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
204: for (i=0;i<ctx->np-1;i++) {
205: SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[i],PETSC_TRUE);
206: PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->np-i-1);
207: }
208: SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[ctx->np-1],PETSC_TRUE);
209: PetscViewerASCIIPrintf(viewer,"%s) / (",str);
210: for (i=0;i<ctx->nq-1;i++) {
211: SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[i],PETSC_TRUE);
212: PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->nq-i-1);
213: }
214: SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[ctx->nq-1],PETSC_TRUE);
215: PetscViewerASCIIPrintf(viewer,"%s)\n",str);
216: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
217: }
218: }
219: PetscFunctionReturn(0);
220: }
222: static PetscErrorCode FNRationalSetNumerator_Rational(FN fn,PetscInt np,PetscScalar *pcoeff)
223: {
224: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
225: PetscInt i;
228: ctx->np = np;
229: PetscFree(ctx->pcoeff);
230: if (np) {
231: PetscMalloc1(np,&ctx->pcoeff);
232: PetscLogObjectMemory((PetscObject)fn,np*sizeof(PetscScalar));
233: for (i=0;i<np;i++) ctx->pcoeff[i] = pcoeff[i];
234: }
235: PetscFunctionReturn(0);
236: }
238: /*@C
239: FNRationalSetNumerator - Sets the parameters defining the numerator of the
240: rational function.
242: Logically Collective on fn
244: Input Parameters:
245: + fn - the math function context
246: . np - number of coefficients
247: - pcoeff - coefficients (array of scalar values)
249: Notes:
250: Let the rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials.
251: This function provides the coefficients of the numerator p(x).
252: Hence, p(x) is of degree np-1.
253: If np is zero, then the numerator is assumed to be p(x)=1.
255: In polynomials, high order coefficients are stored in the first positions
256: of the array, e.g. to represent x^2-3 use {1,0,-3}.
258: Level: intermediate
260: .seealso: FNRationalSetDenominator(), FNRationalGetNumerator()
261: @*/
262: PetscErrorCode FNRationalSetNumerator(FN fn,PetscInt np,PetscScalar *pcoeff)
263: {
267: PetscTryMethod(fn,"FNRationalSetNumerator_C",(FN,PetscInt,PetscScalar*),(fn,np,pcoeff));
268: PetscFunctionReturn(0);
269: }
271: static PetscErrorCode FNRationalGetNumerator_Rational(FN fn,PetscInt *np,PetscScalar *pcoeff[])
272: {
273: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
274: PetscInt i;
276: if (np) *np = ctx->np;
277: if (pcoeff) {
278: if (!ctx->np) *pcoeff = NULL;
279: else {
280: PetscMalloc1(ctx->np,pcoeff);
281: for (i=0;i<ctx->np;i++) (*pcoeff)[i] = ctx->pcoeff[i];
282: }
283: }
284: PetscFunctionReturn(0);
285: }
287: /*@C
288: FNRationalGetNumerator - Gets the parameters that define the numerator of the
289: rational function.
291: Not Collective
293: Input Parameter:
294: . fn - the math function context
296: Output Parameters:
297: + np - number of coefficients
298: - pcoeff - coefficients (array of scalar values, length nq)
300: Notes:
301: The values passed by user with FNRationalSetNumerator() are returned (or null
302: pointers otherwise).
303: The pcoeff array should be freed by the user when no longer needed.
305: Level: intermediate
307: .seealso: FNRationalSetNumerator()
308: @*/
309: PetscErrorCode FNRationalGetNumerator(FN fn,PetscInt *np,PetscScalar *pcoeff[])
310: {
312: PetscUseMethod(fn,"FNRationalGetNumerator_C",(FN,PetscInt*,PetscScalar**),(fn,np,pcoeff));
313: PetscFunctionReturn(0);
314: }
316: static PetscErrorCode FNRationalSetDenominator_Rational(FN fn,PetscInt nq,PetscScalar *qcoeff)
317: {
318: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
319: PetscInt i;
322: ctx->nq = nq;
323: PetscFree(ctx->qcoeff);
324: if (nq) {
325: PetscMalloc1(nq,&ctx->qcoeff);
326: PetscLogObjectMemory((PetscObject)fn,nq*sizeof(PetscScalar));
327: for (i=0;i<nq;i++) ctx->qcoeff[i] = qcoeff[i];
328: }
329: PetscFunctionReturn(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: PetscFunctionReturn(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: PetscFunctionReturn(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: PetscFunctionReturn(0);
408: }
410: PetscErrorCode FNSetFromOptions_Rational(PetscOptionItems *PetscOptionsObject,FN fn)
411: {
412: #define PARMAX 10
413: PetscScalar array[PARMAX];
414: PetscInt i,k;
415: PetscBool flg;
417: PetscOptionsHead(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: PetscOptionsTail();
430: PetscFunctionReturn(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: PetscLogObjectMemory((PetscObject)(*newfn),ctx->np*sizeof(PetscScalar));
442: for (i=0;i<ctx->np;i++) ctx2->pcoeff[i] = ctx->pcoeff[i];
443: }
444: ctx2->nq = ctx->nq;
445: if (ctx->nq) {
446: PetscMalloc1(ctx->nq,&ctx2->qcoeff);
447: PetscLogObjectMemory((PetscObject)(*newfn),ctx->nq*sizeof(PetscScalar));
448: for (i=0;i<ctx->nq;i++) ctx2->qcoeff[i] = ctx->qcoeff[i];
449: }
450: PetscFunctionReturn(0);
451: }
453: PetscErrorCode FNDestroy_Rational(FN fn)
454: {
455: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
457: PetscFree(ctx->pcoeff);
458: PetscFree(ctx->qcoeff);
459: PetscFree(fn->data);
460: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetNumerator_C",NULL);
461: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetNumerator_C",NULL);
462: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetDenominator_C",NULL);
463: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetDenominator_C",NULL);
464: PetscFunctionReturn(0);
465: }
467: SLEPC_EXTERN PetscErrorCode FNCreate_Rational(FN fn)
468: {
469: FN_RATIONAL *ctx;
471: PetscNewLog(fn,&ctx);
472: fn->data = (void*)ctx;
474: fn->ops->evaluatefunction = FNEvaluateFunction_Rational;
475: fn->ops->evaluatederivative = FNEvaluateDerivative_Rational;
476: fn->ops->evaluatefunctionmat[0] = FNEvaluateFunctionMat_Rational;
477: fn->ops->evaluatefunctionmatvec[0] = FNEvaluateFunctionMatVec_Rational;
478: fn->ops->setfromoptions = FNSetFromOptions_Rational;
479: fn->ops->view = FNView_Rational;
480: fn->ops->duplicate = FNDuplicate_Rational;
481: fn->ops->destroy = FNDestroy_Rational;
482: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetNumerator_C",FNRationalSetNumerator_Rational);
483: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetNumerator_C",FNRationalGetNumerator_Rational);
484: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetDenominator_C",FNRationalSetDenominator_Rational);
485: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetDenominator_C",FNRationalGetDenominator_Rational);
486: PetscFunctionReturn(0);
487: }