Actual source code: fnrational.c
slepc-3.7.1 2016-05-27
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-2016, 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: if (q==0.0) SETERRQ(PETSC_COMM_SELF,1,"Function not defined in the requested value");
55: *y = p/q;
56: }
57: return(0);
58: }
62: static PetscErrorCode FNEvaluateFunctionMat_Private(FN fn,PetscScalar *Aa,PetscScalar *Ba,PetscInt m,PetscBool firstonly)
63: {
64: #if defined(PETSC_MISSING_LAPACK_GESV)
66: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESV - Lapack routines are unavailable");
67: #else
69: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
70: PetscBLASInt n,k,ld,*ipiv,info;
71: PetscInt i,j;
72: PetscScalar *W,*P,*Q,one=1.0,zero=0.0;
75: PetscBLASIntCast(m,&n);
76: ld = n;
77: k = firstonly? 1: n;
78: if (Aa==Ba) {
79: PetscMalloc4(m*m,&P,m*m,&Q,m*m,&W,ld,&ipiv);
80: } else {
81: P = Ba;
82: PetscMalloc3(m*m,&Q,m*m,&W,ld,&ipiv);
83: }
84: PetscMemzero(P,m*m*sizeof(PetscScalar));
85: if (!ctx->np) {
86: for (i=0;i<m;i++) P[i+i*ld] = 1.0;
87: } else {
88: for (i=0;i<m;i++) P[i+i*ld] = ctx->pcoeff[0];
89: for (j=1;j<ctx->np;j++) {
90: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&one,P,&ld,Aa,&ld,&zero,W,&ld));
91: PetscMemcpy(P,W,m*m*sizeof(PetscScalar));
92: for (i=0;i<m;i++) P[i+i*ld] += ctx->pcoeff[j];
93: }
94: }
95: if (ctx->nq) {
96: PetscMemzero(Q,m*m*sizeof(PetscScalar));
97: for (i=0;i<m;i++) Q[i+i*ld] = ctx->qcoeff[0];
98: for (j=1;j<ctx->nq;j++) {
99: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&one,Q,&ld,Aa,&ld,&zero,W,&ld));
100: PetscMemcpy(Q,W,m*m*sizeof(PetscScalar));
101: for (i=0;i<m;i++) Q[i+i*ld] += ctx->qcoeff[j];
102: }
103: PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&n,&k,Q,&ld,ipiv,P,&ld,&info));
104: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESV %d",info);
105: }
106: if (Aa==Ba) {
107: PetscMemcpy(Aa,P,m*k*sizeof(PetscScalar));
108: PetscFree4(P,Q,W,ipiv);
109: } else {
110: PetscFree3(Q,W,ipiv);
111: }
112: return(0);
113: #endif
114: }
118: PetscErrorCode FNEvaluateFunctionMat_Rational(FN fn,Mat A,Mat B)
119: {
121: PetscInt m;
122: PetscScalar *Aa,*Ba;
125: MatDenseGetArray(A,&Aa);
126: MatDenseGetArray(B,&Ba);
127: MatGetSize(A,&m,NULL);
128: FNEvaluateFunctionMat_Private(fn,Aa,Ba,m,PETSC_FALSE);
129: MatDenseRestoreArray(A,&Aa);
130: MatDenseRestoreArray(B,&Ba);
131: return(0);
132: }
136: PetscErrorCode FNEvaluateFunctionMatVec_Rational(FN fn,Mat A,Vec v)
137: {
139: PetscInt m;
140: PetscScalar *Aa,*Ba;
141: Mat B;
144: FN_AllocateWorkMat(fn,A,&B);
145: MatDenseGetArray(A,&Aa);
146: MatDenseGetArray(B,&Ba);
147: MatGetSize(A,&m,NULL);
148: FNEvaluateFunctionMat_Private(fn,Aa,Ba,m,PETSC_TRUE);
149: MatDenseRestoreArray(A,&Aa);
150: MatDenseRestoreArray(B,&Ba);
151: MatGetColumnVector(B,v,0);
152: FN_FreeWorkMat(fn,&B);
153: return(0);
154: }
158: PetscErrorCode FNEvaluateDerivative_Rational(FN fn,PetscScalar x,PetscScalar *yp)
159: {
160: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
161: PetscInt i;
162: PetscScalar p,q,pp,qp;
165: if (!ctx->np) {
166: p = 1.0;
167: pp = 0.0;
168: } else {
169: p = ctx->pcoeff[0];
170: pp = 0.0;
171: for (i=1;i<ctx->np;i++) {
172: pp = p+x*pp;
173: p = ctx->pcoeff[i]+x*p;
174: }
175: }
176: if (!ctx->nq) *yp = pp;
177: else {
178: q = ctx->qcoeff[0];
179: qp = 0.0;
180: for (i=1;i<ctx->nq;i++) {
181: qp = q+x*qp;
182: q = ctx->qcoeff[i]+x*q;
183: }
184: if (q==0.0) SETERRQ(PETSC_COMM_SELF,1,"Derivative not defined in the requested value");
185: *yp = (pp*q-p*qp)/(q*q);
186: }
187: return(0);
188: }
192: PetscErrorCode FNView_Rational(FN fn,PetscViewer viewer)
193: {
195: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
196: PetscBool isascii;
197: PetscInt i;
198: char str[50];
201: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
202: if (isascii) {
203: if (fn->alpha!=(PetscScalar)1.0 || fn->beta!=(PetscScalar)1.0) {
204: SlepcSNPrintfScalar(str,50,fn->alpha,PETSC_FALSE);
205: PetscViewerASCIIPrintf(viewer," Scale factors: alpha=%s,",str);
206: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
207: SlepcSNPrintfScalar(str,50,fn->beta,PETSC_FALSE);
208: PetscViewerASCIIPrintf(viewer," beta=%s\n",str);
209: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
210: }
211: if (!ctx->nq) {
212: if (!ctx->np) {
213: PetscViewerASCIIPrintf(viewer," Constant: 1.0\n");
214: } else if (ctx->np==1) {
215: SlepcSNPrintfScalar(str,50,ctx->pcoeff[0],PETSC_FALSE);
216: PetscViewerASCIIPrintf(viewer," Constant: %s\n",str);
217: } else {
218: PetscViewerASCIIPrintf(viewer," Polynomial: ");
219: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
220: for (i=0;i<ctx->np-1;i++) {
221: SlepcSNPrintfScalar(str,50,ctx->pcoeff[i],PETSC_TRUE);
222: PetscViewerASCIIPrintf(viewer,"%s*x^%1D",str,ctx->np-i-1);
223: }
224: SlepcSNPrintfScalar(str,50,ctx->pcoeff[ctx->np-1],PETSC_TRUE);
225: PetscViewerASCIIPrintf(viewer,"%s\n",str);
226: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
227: }
228: } else if (!ctx->np) {
229: PetscViewerASCIIPrintf(viewer," Inverse polinomial: 1 / (");
230: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
231: for (i=0;i<ctx->nq-1;i++) {
232: SlepcSNPrintfScalar(str,50,ctx->qcoeff[i],PETSC_TRUE);
233: PetscViewerASCIIPrintf(viewer,"%s*x^%1D",str,ctx->nq-i-1);
234: }
235: SlepcSNPrintfScalar(str,50,ctx->qcoeff[ctx->nq-1],PETSC_TRUE);
236: PetscViewerASCIIPrintf(viewer,"%s)\n",str);
237: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
238: } else {
239: PetscViewerASCIIPrintf(viewer," Rational function: (");
240: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
241: for (i=0;i<ctx->np-1;i++) {
242: SlepcSNPrintfScalar(str,50,ctx->pcoeff[i],PETSC_TRUE);
243: PetscViewerASCIIPrintf(viewer,"%s*x^%1D",str,ctx->np-i-1);
244: }
245: SlepcSNPrintfScalar(str,50,ctx->pcoeff[ctx->np-1],PETSC_TRUE);
246: PetscViewerASCIIPrintf(viewer,"%s) / (",str);
247: for (i=0;i<ctx->nq-1;i++) {
248: SlepcSNPrintfScalar(str,50,ctx->qcoeff[i],PETSC_TRUE);
249: PetscViewerASCIIPrintf(viewer,"%s*x^%1D",str,ctx->nq-i-1);
250: }
251: SlepcSNPrintfScalar(str,50,ctx->qcoeff[ctx->nq-1],PETSC_TRUE);
252: PetscViewerASCIIPrintf(viewer,"%s)\n",str);
253: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
254: }
255: }
256: return(0);
257: }
261: static PetscErrorCode FNRationalSetNumerator_Rational(FN fn,PetscInt np,PetscScalar *pcoeff)
262: {
264: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
265: PetscInt i;
268: ctx->np = np;
269: PetscFree(ctx->pcoeff);
270: if (np) {
271: PetscMalloc1(np,&ctx->pcoeff);
272: PetscLogObjectMemory((PetscObject)fn,np*sizeof(PetscScalar));
273: for (i=0;i<np;i++) ctx->pcoeff[i] = pcoeff[i];
274: }
275: return(0);
276: }
280: /*@
281: FNRationalSetNumerator - Sets the parameters defining the numerator of the
282: rational function.
284: Logically Collective on FN
286: Input Parameters:
287: + fn - the math function context
288: . np - number of coefficients
289: - pcoeff - coefficients (array of scalar values)
291: Notes:
292: Let the rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials.
293: This function provides the coefficients of the numerator p(x).
294: Hence, p(x) is of degree np-1.
295: If np is zero, then the numerator is assumed to be p(x)=1.
297: In polynomials, high order coefficients are stored in the first positions
298: of the array, e.g. to represent x^2-3 use {1,0,-3}.
300: Level: intermediate
302: .seealso: FNRationalSetDenominator(), FNRationalGetNumerator()
303: @*/
304: PetscErrorCode FNRationalSetNumerator(FN fn,PetscInt np,PetscScalar *pcoeff)
305: {
311: if (np<0) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Argument np cannot be negative");
313: PetscTryMethod(fn,"FNRationalSetNumerator_C",(FN,PetscInt,PetscScalar*),(fn,np,pcoeff));
314: return(0);
315: }
319: static PetscErrorCode FNRationalGetNumerator_Rational(FN fn,PetscInt *np,PetscScalar *pcoeff[])
320: {
322: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
323: PetscInt i;
326: if (np) *np = ctx->np;
327: if (pcoeff) {
328: if (!ctx->np) *pcoeff = NULL;
329: else {
330: PetscMalloc1(ctx->np,pcoeff);
331: for (i=0;i<ctx->np;i++) (*pcoeff)[i] = ctx->pcoeff[i];
332: }
333: }
334: return(0);
335: }
339: /*@
340: FNRationalGetNumerator - Gets the parameters that define the numerator of the
341: rational function.
343: Not Collective
345: Input Parameter:
346: . fn - the math function context
348: Output Parameters:
349: + np - number of coefficients
350: - pcoeff - coefficients (array of scalar values, length nq)
352: Notes:
353: The values passed by user with FNRationalSetNumerator() are returned (or null
354: pointers otherwise).
355: The pcoeff array should be freed by the user when no longer needed.
357: Level: intermediate
359: .seealso: FNRationalSetNumerator()
360: @*/
361: PetscErrorCode FNRationalGetNumerator(FN fn,PetscInt *np,PetscScalar *pcoeff[])
362: {
367: PetscUseMethod(fn,"FNRationalGetNumerator_C",(FN,PetscInt*,PetscScalar**),(fn,np,pcoeff));
368: return(0);
369: }
373: static PetscErrorCode FNRationalSetDenominator_Rational(FN fn,PetscInt nq,PetscScalar *qcoeff)
374: {
376: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
377: PetscInt i;
380: ctx->nq = nq;
381: PetscFree(ctx->qcoeff);
382: if (nq) {
383: PetscMalloc1(nq,&ctx->qcoeff);
384: PetscLogObjectMemory((PetscObject)fn,nq*sizeof(PetscScalar));
385: for (i=0;i<nq;i++) ctx->qcoeff[i] = qcoeff[i];
386: }
387: return(0);
388: }
392: /*@
393: FNRationalSetDenominator - Sets the parameters defining the denominator of the
394: rational function.
396: Logically Collective on FN
398: Input Parameters:
399: + fn - the math function context
400: . nq - number of coefficients
401: - qcoeff - coefficients (array of scalar values)
403: Notes:
404: Let the rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials.
405: This function provides the coefficients of the denominator q(x).
406: Hence, q(x) is of degree nq-1.
407: If nq is zero, then the function is assumed to be polynomial, r(x) = p(x).
409: In polynomials, high order coefficients are stored in the first positions
410: of the array, e.g. to represent x^2-3 use {1,0,-3}.
412: Level: intermediate
414: .seealso: FNRationalSetNumerator(), FNRationalGetDenominator()
415: @*/
416: PetscErrorCode FNRationalSetDenominator(FN fn,PetscInt nq,PetscScalar *qcoeff)
417: {
423: if (nq<0) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Argument nq cannot be negative");
425: PetscTryMethod(fn,"FNRationalSetDenominator_C",(FN,PetscInt,PetscScalar*),(fn,nq,qcoeff));
426: return(0);
427: }
431: static PetscErrorCode FNRationalGetDenominator_Rational(FN fn,PetscInt *nq,PetscScalar *qcoeff[])
432: {
434: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
435: PetscInt i;
438: if (nq) *nq = ctx->nq;
439: if (qcoeff) {
440: if (!ctx->nq) *qcoeff = NULL;
441: else {
442: PetscMalloc1(ctx->nq,qcoeff);
443: for (i=0;i<ctx->nq;i++) (*qcoeff)[i] = ctx->qcoeff[i];
444: }
445: }
446: return(0);
447: }
451: /*@
452: FNRationalGetDenominator - Gets the parameters that define the denominator of the
453: rational function.
455: Not Collective
457: Input Parameter:
458: . fn - the math function context
460: Output Parameters:
461: + nq - number of coefficients
462: - qcoeff - coefficients (array of scalar values, length nq)
464: Notes:
465: The values passed by user with FNRationalSetDenominator() are returned (or null
466: pointers otherwise).
467: The qcoeff array should be freed by the user when no longer needed.
469: Level: intermediate
471: .seealso: FNRationalSetDenominator()
472: @*/
473: PetscErrorCode FNRationalGetDenominator(FN fn,PetscInt *nq,PetscScalar *qcoeff[])
474: {
479: PetscUseMethod(fn,"FNRationalGetDenominator_C",(FN,PetscInt*,PetscScalar**),(fn,nq,qcoeff));
480: return(0);
481: }
485: PetscErrorCode FNSetFromOptions_Rational(PetscOptionItems *PetscOptionsObject,FN fn)
486: {
488: #define PARMAX 10
489: PetscScalar array[PARMAX];
490: PetscInt i,k;
491: PetscBool flg;
494: PetscOptionsHead(PetscOptionsObject,"FN Rational Options");
496: k = PARMAX;
497: for (i=0;i<k;i++) array[i] = 0;
498: PetscOptionsScalarArray("-fn_rational_numerator","Numerator coefficients (one or more scalar values separated with a comma without spaces)","FNRationalSetNumerator",array,&k,&flg);
499: if (flg) {
500: FNRationalSetNumerator(fn,k,array);
501: }
503: k = PARMAX;
504: for (i=0;i<k;i++) array[i] = 0;
505: PetscOptionsScalarArray("-fn_rational_denominator","Denominator coefficients (one or more scalar values separated with a comma without spaces)","FNRationalSetDenominator",array,&k,&flg);
506: if (flg) {
507: FNRationalSetDenominator(fn,k,array);
508: }
510: PetscOptionsTail();
511: return(0);
512: }
516: PetscErrorCode FNDuplicate_Rational(FN fn,MPI_Comm comm,FN *newfn)
517: {
519: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data,*ctx2;
520: PetscInt i;
523: PetscNewLog(*newfn,&ctx2);
524: (*newfn)->data = (void*)ctx2;
525: ctx2->np = ctx->np;
526: if (ctx->np) {
527: PetscMalloc1(ctx->np,&ctx2->pcoeff);
528: PetscLogObjectMemory((PetscObject)(*newfn),ctx->np*sizeof(PetscScalar));
529: for (i=0;i<ctx->np;i++) ctx2->pcoeff[i] = ctx->pcoeff[i];
530: }
531: ctx2->nq = ctx->nq;
532: if (ctx->nq) {
533: PetscMalloc1(ctx->nq,&ctx2->qcoeff);
534: PetscLogObjectMemory((PetscObject)(*newfn),ctx->nq*sizeof(PetscScalar));
535: for (i=0;i<ctx->nq;i++) ctx2->qcoeff[i] = ctx->qcoeff[i];
536: }
537: return(0);
538: }
542: PetscErrorCode FNDestroy_Rational(FN fn)
543: {
545: FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
548: PetscFree(ctx->pcoeff);
549: PetscFree(ctx->qcoeff);
550: PetscFree(fn->data);
551: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetNumerator_C",NULL);
552: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetNumerator_C",NULL);
553: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetDenominator_C",NULL);
554: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetDenominator_C",NULL);
555: return(0);
556: }
560: PETSC_EXTERN PetscErrorCode FNCreate_Rational(FN fn)
561: {
563: FN_RATIONAL *ctx;
566: PetscNewLog(fn,&ctx);
567: fn->data = (void*)ctx;
569: fn->ops->evaluatefunction = FNEvaluateFunction_Rational;
570: fn->ops->evaluatederivative = FNEvaluateDerivative_Rational;
571: fn->ops->evaluatefunctionmat = FNEvaluateFunctionMat_Rational;
572: fn->ops->evaluatefunctionmatvec = FNEvaluateFunctionMatVec_Rational;
573: fn->ops->setfromoptions = FNSetFromOptions_Rational;
574: fn->ops->view = FNView_Rational;
575: fn->ops->duplicate = FNDuplicate_Rational;
576: fn->ops->destroy = FNDestroy_Rational;
577: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetNumerator_C",FNRationalSetNumerator_Rational);
578: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetNumerator_C",FNRationalGetNumerator_Rational);
579: PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetDenominator_C",FNRationalSetDenominator_Rational);
580: PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetDenominator_C",FNRationalGetDenominator_Rational);
581: return(0);
582: }