Actual source code: fnrational.c
slepc-3.6.1 2015-09-03
1: /*
2: Rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc/private/fnimpl.h> /*I "slepcfn.h" I*/
25: #include <slepcblaslapack.h>
27: typedef struct {
28: PetscScalar *pcoeff; /* numerator coefficients */
29: PetscInt np; /* length of array pcoeff, p(x) has degree np-1 */
30: PetscScalar *qcoeff; /* denominator coefficients */
31: PetscInt nq; /* length of array qcoeff, q(x) has degree nq-1 */
32: } FN_RATIONAL;
36: PetscErrorCode FNEvaluateFunction_Rational(FN fn,PetscScalar x,PetscScalar *y)
37: {
38: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
39: PetscInt i;
40: PetscScalar p,q;
43: if (!ctx->np) p = 1.0;
44: else {
45: p = ctx->pcoeff[0];
46: for (i=1;i<ctx->np;i++)
47: p = ctx->pcoeff[i]+x*p;
48: }
49: if (!ctx->nq) *y = p;
50: else {
51: q = ctx->qcoeff[0];
52: for (i=1;i<ctx->nq;i++)
53: q = ctx->qcoeff[i]+x*q;
54: *y = p/q;
55: }
56: return(0);
57: }
61: PetscErrorCode FNEvaluateFunctionMat_Rational(FN fn,Mat A,Mat B)
62: {
63: #if defined(PETSC_MISSING_LAPACK_GESV)
65: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESV - Lapack routines are unavailable");
66: #else
68: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
69: PetscBLASInt n,ld,*ipiv,info;
70: PetscInt i,j,m;
71: PetscScalar *Aa,*Ba,*W,*P,*Q,one=1.0,zero=0.0;
74: MatDenseGetArray(A,&Aa);
75: MatDenseGetArray(B,&Ba);
76: MatGetSize(A,&m,NULL);
77: PetscBLASIntCast(m,&n);
78: ld = n;
79: P = Ba;
80: PetscMalloc3(m*m,&Q,m*m,&W,ld,&ipiv);
81: PetscMemzero(P,m*m*sizeof(PetscScalar));
82: if (!ctx->np) {
83: for (i=0;i<m;i++) P[i+i*ld] = 1.0;
84: } else {
85: for (i=0;i<m;i++) P[i+i*ld] = ctx->pcoeff[0];
86: for (j=1;j<ctx->np;j++) {
87: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&one,P,&ld,Aa,&ld,&zero,W,&ld));
88: PetscMemcpy(P,W,m*m*sizeof(PetscScalar));
89: for (i=0;i<m;i++) P[i+i*ld] += ctx->pcoeff[j];
90: }
91: }
92: if (ctx->nq) {
93: PetscMemzero(Q,m*m*sizeof(PetscScalar));
94: for (i=0;i<m;i++) Q[i+i*ld] = ctx->qcoeff[0];
95: for (j=1;j<ctx->nq;j++) {
96: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&one,Q,&ld,Aa,&ld,&zero,W,&ld));
97: PetscMemcpy(Q,W,m*m*sizeof(PetscScalar));
98: for (i=0;i<m;i++) Q[i+i*ld] += ctx->qcoeff[j];
99: }
100: PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&n,&n,Q,&ld,ipiv,P,&ld,&info));
101: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESV %d",info);
102: }
103: PetscFree3(Q,W,ipiv);
104: MatDenseRestoreArray(A,&Aa);
105: MatDenseRestoreArray(B,&Ba);
106: return(0);
107: #endif
108: }
112: PetscErrorCode FNEvaluateDerivative_Rational(FN fn,PetscScalar x,PetscScalar *yp)
113: {
114: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
115: PetscInt i;
116: PetscScalar p,q,pp,qp;
119: if (!ctx->np) {
120: p = 1.0;
121: pp = 0.0;
122: } else {
123: p = ctx->pcoeff[0];
124: pp = 0.0;
125: for (i=1;i<ctx->np;i++) {
126: pp = p+x*pp;
127: p = ctx->pcoeff[i]+x*p;
128: }
129: }
130: if (!ctx->nq) *yp = pp;
131: else {
132: q = ctx->qcoeff[0];
133: qp = 0.0;
134: for (i=1;i<ctx->nq;i++) {
135: qp = q+x*qp;
136: q = ctx->qcoeff[i]+x*q;
137: }
138: *yp = (pp*q-p*qp)/(q*q);
139: }
140: return(0);
141: }
145: PetscErrorCode FNView_Rational(FN fn,PetscViewer viewer)
146: {
148: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
149: PetscBool isascii;
150: PetscInt i;
151: char str[50];
154: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
155: if (isascii) {
156: if (fn->alpha!=(PetscScalar)1.0 || fn->beta!=(PetscScalar)1.0) {
157: SlepcSNPrintfScalar(str,50,fn->alpha,PETSC_FALSE);
158: PetscViewerASCIIPrintf(viewer," Scale factors: alpha=%s,",str);
159: SlepcSNPrintfScalar(str,50,fn->beta,PETSC_FALSE);
160: PetscViewerASCIIPrintf(viewer," beta=%s\n",str);
161: }
162: if (!ctx->nq) {
163: if (!ctx->np) {
164: PetscViewerASCIIPrintf(viewer," Constant: 1.0\n");
165: } else if (ctx->np==1) {
166: SlepcSNPrintfScalar(str,50,ctx->pcoeff[0],PETSC_FALSE);
167: PetscViewerASCIIPrintf(viewer," Constant: %s\n",str);
168: } else {
169: PetscViewerASCIIPrintf(viewer," Polynomial: ");
170: for (i=0;i<ctx->np-1;i++) {
171: SlepcSNPrintfScalar(str,50,ctx->pcoeff[i],PETSC_TRUE);
172: PetscViewerASCIIPrintf(viewer,"%s*x^%1D",str,ctx->np-i-1);
173: }
174: SlepcSNPrintfScalar(str,50,ctx->pcoeff[ctx->np-1],PETSC_TRUE);
175: PetscViewerASCIIPrintf(viewer,"%s\n",str);
176: }
177: } else if (!ctx->np) {
178: PetscViewerASCIIPrintf(viewer," Inverse polinomial: 1 / (");
179: for (i=0;i<ctx->nq-1;i++) {
180: SlepcSNPrintfScalar(str,50,ctx->qcoeff[i],PETSC_TRUE);
181: PetscViewerASCIIPrintf(viewer,"%s*x^%1D",str,ctx->nq-i-1);
182: }
183: SlepcSNPrintfScalar(str,50,ctx->qcoeff[ctx->nq-1],PETSC_TRUE);
184: PetscViewerASCIIPrintf(viewer,"%s)\n",str);
185: } else {
186: PetscViewerASCIIPrintf(viewer," Rational function: (");
187: for (i=0;i<ctx->np-1;i++) {
188: SlepcSNPrintfScalar(str,50,ctx->pcoeff[i],PETSC_TRUE);
189: PetscViewerASCIIPrintf(viewer,"%s*x^%1D",str,ctx->np-i-1);
190: }
191: SlepcSNPrintfScalar(str,50,ctx->pcoeff[ctx->np-1],PETSC_TRUE);
192: PetscViewerASCIIPrintf(viewer,"%s) / (",str);
193: for (i=0;i<ctx->nq-1;i++) {
194: SlepcSNPrintfScalar(str,50,ctx->qcoeff[i],PETSC_TRUE);
195: PetscViewerASCIIPrintf(viewer,"%s*x^%1D",str,ctx->nq-i-1);
196: }
197: SlepcSNPrintfScalar(str,50,ctx->qcoeff[ctx->nq-1],PETSC_TRUE);
198: PetscViewerASCIIPrintf(viewer,"%s)\n",str);
199: }
200: }
201: return(0);
202: }
206: static PetscErrorCode FNRationalSetNumerator_Rational(FN fn,PetscInt np,PetscScalar *pcoeff)
207: {
209: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
210: PetscInt i;
213: ctx->np = np;
214: PetscFree(ctx->pcoeff);
215: if (np) {
216: PetscMalloc1(np,&ctx->pcoeff);
217: PetscLogObjectMemory((PetscObject)fn,np*sizeof(PetscScalar));
218: for (i=0;i<np;i++) ctx->pcoeff[i] = pcoeff[i];
219: }
220: return(0);
221: }
225: /*@
226: FNRationalSetNumerator - Sets the parameters defining the numerator of the
227: rational function.
229: Logically Collective on FN
231: Input Parameters:
232: + fn - the math function context
233: . np - number of coefficients
234: - pcoeff - coefficients (array of scalar values)
236: Notes:
237: Let the rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials.
238: This function provides the coefficients of the numerator p(x).
239: Hence, p(x) is of degree np-1.
240: If np is zero, then the numerator is assumed to be p(x)=1.
242: In polynomials, high order coefficients are stored in the first positions
243: of the array, e.g. to represent x^2-3 use {1,0,-3}.
245: Level: intermediate
247: .seealso: FNRationalSetDenominator(), FNRationalGetNumerator()
248: @*/
249: PetscErrorCode FNRationalSetNumerator(FN fn,PetscInt np,PetscScalar *pcoeff)
250: {
256: if (np<0) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Argument np cannot be negative");
258: PetscTryMethod(fn,"FNRationalSetNumerator_C",(FN,PetscInt,PetscScalar*),(fn,np,pcoeff));
259: return(0);
260: }
264: static PetscErrorCode FNRationalGetNumerator_Rational(FN fn,PetscInt *np,PetscScalar *pcoeff[])
265: {
267: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
268: PetscInt i;
271: if (np) *np = ctx->np;
272: if (pcoeff) {
273: if (!ctx->np) *pcoeff = NULL;
274: else {
275: PetscMalloc1(ctx->np,pcoeff);
276: for (i=0;i<ctx->np;i++) (*pcoeff)[i] = ctx->pcoeff[i];
277: }
278: }
279: return(0);
280: }
284: /*@
285: FNRationalGetNumerator - Gets the parameters that define the numerator of the
286: rational function.
288: Not Collective
290: Input Parameter:
291: . fn - the math function context
293: Output Parameters:
294: + np - number of coefficients
295: - pcoeff - coefficients (array of scalar values, length nq)
297: Notes:
298: The values passed by user with FNRationalSetNumerator() are returned (or null
299: pointers otherwise).
300: The pcoeff array should be freed by the user when no longer needed.
302: Level: intermediate
304: .seealso: FNRationalSetNumerator()
305: @*/
306: PetscErrorCode FNRationalGetNumerator(FN fn,PetscInt *np,PetscScalar *pcoeff[])
307: {
312: PetscTryMethod(fn,"FNRationalGetNumerator_C",(FN,PetscInt*,PetscScalar**),(fn,np,pcoeff));
313: return(0);
314: }
318: static PetscErrorCode FNRationalSetDenominator_Rational(FN fn,PetscInt nq,PetscScalar *qcoeff)
319: {
321: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
322: PetscInt i;
325: ctx->nq = nq;
326: PetscFree(ctx->qcoeff);
327: if (nq) {
328: PetscMalloc1(nq,&ctx->qcoeff);
329: PetscLogObjectMemory((PetscObject)fn,nq*sizeof(PetscScalar));
330: for (i=0;i<nq;i++) ctx->qcoeff[i] = qcoeff[i];
331: }
332: return(0);
333: }
337: /*@
338: FNRationalSetDenominator - Sets the parameters defining the denominator of the
339: rational function.
341: Logically Collective on FN
343: Input Parameters:
344: + fn - the math function context
345: . nq - number of coefficients
346: - qcoeff - coefficients (array of scalar values)
348: Notes:
349: Let the rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials.
350: This function provides the coefficients of the denominator q(x).
351: Hence, q(x) is of degree nq-1.
352: If nq is zero, then the function is assumed to be polynomial, r(x) = p(x).
354: In polynomials, high order coefficients are stored in the first positions
355: of the array, e.g. to represent x^2-3 use {1,0,-3}.
357: Level: intermediate
359: .seealso: FNRationalSetNumerator(), FNRationalGetDenominator()
360: @*/
361: PetscErrorCode FNRationalSetDenominator(FN fn,PetscInt nq,PetscScalar *qcoeff)
362: {
368: if (nq<0) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Argument nq cannot be negative");
370: PetscTryMethod(fn,"FNRationalSetDenominator_C",(FN,PetscInt,PetscScalar*),(fn,nq,qcoeff));
371: return(0);
372: }
376: static PetscErrorCode FNRationalGetDenominator_Rational(FN fn,PetscInt *nq,PetscScalar *qcoeff[])
377: {
379: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
380: PetscInt i;
383: if (nq) *nq = ctx->nq;
384: if (qcoeff) {
385: if (!ctx->nq) *qcoeff = NULL;
386: else {
387: PetscMalloc1(ctx->nq,qcoeff);
388: for (i=0;i<ctx->nq;i++) (*qcoeff)[i] = ctx->qcoeff[i];
389: }
390: }
391: return(0);
392: }
396: /*@
397: FNRationalGetDenominator - Gets the parameters that define the denominator of the
398: rational function.
400: Not Collective
402: Input Parameter:
403: . fn - the math function context
405: Output Parameters:
406: + nq - number of coefficients
407: - qcoeff - coefficients (array of scalar values, length nq)
409: Notes:
410: The values passed by user with FNRationalSetDenominator() are returned (or null
411: pointers otherwise).
412: The qcoeff array should be freed by the user when no longer needed.
414: Level: intermediate
416: .seealso: FNRationalSetDenominator()
417: @*/
418: PetscErrorCode FNRationalGetDenominator(FN fn,PetscInt *nq,PetscScalar *qcoeff[])
419: {
424: PetscTryMethod(fn,"FNRationalGetDenominator_C",(FN,PetscInt*,PetscScalar**),(fn,nq,qcoeff));
425: return(0);
426: }
430: PetscErrorCode FNSetFromOptions_Rational(PetscOptions *PetscOptionsObject,FN fn)
431: {
433: #define PARMAX 10
434: PetscScalar array[PARMAX];
435: PetscInt i,k;
436: PetscBool flg;
439: PetscOptionsHead(PetscOptionsObject,"FN Rational Options");
441: k = PARMAX;
442: for (i=0;i<k;i++) array[i] = 0;
443: PetscOptionsScalarArray("-fn_rational_numerator","Numerator coefficients (one or more scalar values separated with a comma without spaces)","FNRationalSetNumerator",array,&k,&flg);
444: if (flg) {
445: FNRationalSetNumerator(fn,k,array);
446: }
448: k = PARMAX;
449: for (i=0;i<k;i++) array[i] = 0;
450: PetscOptionsScalarArray("-fn_rational_denominator","Denominator coefficients (one or more scalar values separated with a comma without spaces)","FNRationalSetDenominator",array,&k,&flg);
451: if (flg) {
452: FNRationalSetDenominator(fn,k,array);
453: }
455: PetscOptionsTail();
456: return(0);
457: }
461: PetscErrorCode FNDuplicate_Rational(FN fn,MPI_Comm comm,FN *newfn)
462: {
464: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data,*ctx2;
465: PetscInt i;
468: PetscNewLog(*newfn,&ctx2);
469: (*newfn)->data = (void*)ctx2;
470: ctx2->np = ctx->np;
471: if (ctx->np) {
472: PetscMalloc1(ctx->np,&ctx2->pcoeff);
473: PetscLogObjectMemory((PetscObject)(*newfn),ctx->np*sizeof(PetscScalar));
474: for (i=0;i<ctx->np;i++) ctx2->pcoeff[i] = ctx->pcoeff[i];
475: }
476: ctx2->nq = ctx->nq;
477: if (ctx->nq) {
478: PetscMalloc1(ctx->nq,&ctx2->qcoeff);
479: PetscLogObjectMemory((PetscObject)(*newfn),ctx->nq*sizeof(PetscScalar));
480: for (i=0;i<ctx->nq;i++) ctx2->qcoeff[i] = ctx->qcoeff[i];
481: }
482: return(0);
483: }
487: PetscErrorCode FNDestroy_Rational(FN fn)
488: {
490: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
493: PetscFree(ctx->pcoeff);
494: PetscFree(ctx->qcoeff);
495: PetscFree(fn->data);
496: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetNumerator_C",NULL);
497: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetNumerator_C",NULL);
498: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetDenominator_C",NULL);
499: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetDenominator_C",NULL);
500: return(0);
501: }
505: PETSC_EXTERN PetscErrorCode FNCreate_Rational(FN fn)
506: {
508: FN_RATIONAL *ctx;
511: PetscNewLog(fn,&ctx);
512: fn->data = (void*)ctx;
514: fn->ops->evaluatefunction = FNEvaluateFunction_Rational;
515: fn->ops->evaluatederivative = FNEvaluateDerivative_Rational;
516: fn->ops->evaluatefunctionmat = FNEvaluateFunctionMat_Rational;
517: fn->ops->setfromoptions = FNSetFromOptions_Rational;
518: fn->ops->view = FNView_Rational;
519: fn->ops->duplicate = FNDuplicate_Rational;
520: fn->ops->destroy = FNDestroy_Rational;
521: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetNumerator_C",FNRationalSetNumerator_Rational);
522: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetNumerator_C",FNRationalGetNumerator_Rational);
523: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetDenominator_C",FNRationalSetDenominator_Rational);
524: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetDenominator_C",FNRationalGetDenominator_Rational);
525: return(0);
526: }