Actual source code: fnrational.c

  1: /*
  2:    Rational function  r(x) = p(x)/q(x), where p(x) is a polynomial of
  3:    degree na and q(x) is a polynomial of degree nb (can be 0).

  5:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  6:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  7:    Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain

  9:    This file is part of SLEPc.

 11:    SLEPc is free software: you can redistribute it and/or modify it under  the
 12:    terms of version 3 of the GNU Lesser General Public License as published by
 13:    the Free Software Foundation.

 15:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 16:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 17:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 18:    more details.

 20:    You  should have received a copy of the GNU Lesser General  Public  License
 21:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 22:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 23: */

 25: #include <slepc-private/fnimpl.h>      /*I "slepcfn.h" I*/

 29: PetscErrorCode FNEvaluateFunction_Rational(FN fn,PetscScalar x,PetscScalar *y)
 30: {
 31:   PetscInt    i;
 32:   PetscScalar p,q;

 35:   if (!fn->na) p = 1.0;
 36:   else {
 37:     p = fn->alpha[0];
 38:     for (i=1;i<fn->na;i++)
 39:       p = fn->alpha[i]+x*p;
 40:   }
 41:   if (!fn->nb) *y = p;
 42:   else {
 43:     q = fn->beta[0];
 44:     for (i=1;i<fn->nb;i++)
 45:       q = fn->beta[i]+x*q;
 46:     *y = p/q;
 47:   }
 48:   return(0);
 49: }

 53: PetscErrorCode FNEvaluateDerivative_Rational(FN fn,PetscScalar x,PetscScalar *yp)
 54: {
 55:   PetscInt    i;
 56:   PetscScalar p,q,pp,qp;

 59:   if (!fn->na) {
 60:     p = 1.0;
 61:     pp = 0.0;
 62:   } else {
 63:     p = fn->alpha[0];
 64:     pp = 0.0;
 65:     for (i=1;i<fn->na;i++) {
 66:       pp = p+x*pp;
 67:       p = fn->alpha[i]+x*p;
 68:     }
 69:   }
 70:   if (!fn->nb) *yp = pp;
 71:   else {
 72:     q = fn->beta[0];
 73:     qp = 0.0;
 74:     for (i=1;i<fn->nb;i++) {
 75:       qp = q+x*qp;
 76:       q = fn->beta[i]+x*q;
 77:     }
 78:     *yp = (pp*q-p*qp)/(q*q);
 79:   }
 80:   return(0);
 81: }

 85: PetscErrorCode FNView_Rational(FN fn,PetscViewer viewer)
 86: {
 88:   PetscBool      isascii;
 89:   PetscInt       i;
 90:   char           str[50];

 93:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 94:   if (isascii) {
 95:     if (!fn->nb) {
 96:       if (!fn->na) {
 97:         PetscViewerASCIIPrintf(viewer,"  Constant: 1.0\n");
 98:       } else if (fn->na==1) {
 99:         SlepcSNPrintfScalar(str,50,fn->alpha[0],PETSC_FALSE);
100:         PetscViewerASCIIPrintf(viewer,"  Constant: %s\n",str);
101:       } else {
102:         PetscViewerASCIIPrintf(viewer,"  Polynomial: ");
103:         for (i=0;i<fn->na-1;i++) {
104:           SlepcSNPrintfScalar(str,50,fn->alpha[i],PETSC_TRUE);
105:           PetscViewerASCIIPrintf(viewer,"%s*x^%1D",str,fn->na-i-1);
106:         }
107:         SlepcSNPrintfScalar(str,50,fn->alpha[fn->na-1],PETSC_TRUE);
108:         PetscViewerASCIIPrintf(viewer,"%s\n",str);
109:       }
110:     } else if (!fn->na) {
111:       PetscViewerASCIIPrintf(viewer,"  Inverse polinomial: 1 / (");
112:       for (i=0;i<fn->nb-1;i++) {
113:         SlepcSNPrintfScalar(str,50,fn->beta[i],PETSC_TRUE);
114:         PetscViewerASCIIPrintf(viewer,"%s*x^%1D",str,fn->nb-i-1);
115:       }
116:       SlepcSNPrintfScalar(str,50,fn->beta[fn->nb-1],PETSC_TRUE);
117:       PetscViewerASCIIPrintf(viewer,"%s)\n",str);
118:     } else {
119:       PetscViewerASCIIPrintf(viewer,"  Rational function: (");
120:       for (i=0;i<fn->na-1;i++) {
121:         SlepcSNPrintfScalar(str,50,fn->alpha[i],PETSC_TRUE);
122:         PetscViewerASCIIPrintf(viewer,"%s*x^%1D",str,fn->na-i-1);
123:       }
124:         SlepcSNPrintfScalar(str,50,fn->alpha[fn->na-1],PETSC_TRUE);
125:       PetscViewerASCIIPrintf(viewer,"%s) / (",str);
126:       for (i=0;i<fn->nb-1;i++) {
127:         SlepcSNPrintfScalar(str,50,fn->beta[i],PETSC_TRUE);
128:         PetscViewerASCIIPrintf(viewer,"%s*x^%1D",str,fn->nb-i-1);
129:       }
130:       SlepcSNPrintfScalar(str,50,fn->beta[fn->nb-1],PETSC_TRUE);
131:       PetscViewerASCIIPrintf(viewer,"%s)\n",str);
132:     }
133:   }
134:   return(0);
135: }

139: PETSC_EXTERN PetscErrorCode FNCreate_Rational(FN fn)
140: {
142:   fn->ops->evaluatefunction   = FNEvaluateFunction_Rational;
143:   fn->ops->evaluatederivative = FNEvaluateDerivative_Rational;
144:   fn->ops->view               = FNView_Rational;
145:   return(0);
146: }