Actual source code: test1.c

slepc-3.6.3 2016-03-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

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

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

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

 22: static char help[] = "Test rational function.\n\n";

 24: #include <slepcfn.h>

 28: int main(int argc,char **argv)
 29: {
 31:   FN             fn;
 32:   PetscInt       na,nb;
 33:   PetscScalar    x,y,yp,p[10],q[10],five=5.0;
 34:   char           strx[50],str[50];

 36:   SlepcInitialize(&argc,&argv,(char*)0,help);
 37:   FNCreate(PETSC_COMM_WORLD,&fn);

 39:   /* polynomial p(x) */
 40:   na = 5;
 41:   p[0] = -3.1; p[1] = 1.1; p[2] = 1.0; p[3] = -2.0; p[4] = 3.5;
 42:   FNSetType(fn,FNRATIONAL);
 43:   FNRationalSetNumerator(fn,na,p);
 44:   FNView(fn,NULL);
 45:   x = 2.2;
 46:   SlepcSNPrintfScalar(strx,50,x,PETSC_FALSE);
 47:   FNEvaluateFunction(fn,x,&y);
 48:   FNEvaluateDerivative(fn,x,&yp);
 49:   SlepcSNPrintfScalar(str,50,y,PETSC_FALSE);
 50:   PetscPrintf(PETSC_COMM_WORLD,"  f(%s)=%s\n",strx,str);
 51:   SlepcSNPrintfScalar(str,50,yp,PETSC_FALSE);
 52:   PetscPrintf(PETSC_COMM_WORLD,"  f'(%s)=%s\n",strx,str);

 54:   /* inverse of polynomial 1/q(x) */
 55:   nb = 3;
 56:   q[0] = -3.1; q[1] = 1.1; q[2] = 1.0;
 57:   FNSetType(fn,FNRATIONAL);
 58:   FNRationalSetNumerator(fn,0,NULL);  /* reset previous values */
 59:   FNRationalSetDenominator(fn,nb,q);
 60:   FNView(fn,NULL);
 61:   x = 2.2;
 62:   SlepcSNPrintfScalar(strx,50,x,PETSC_FALSE);
 63:   FNEvaluateFunction(fn,x,&y);
 64:   FNEvaluateDerivative(fn,x,&yp);
 65:   SlepcSNPrintfScalar(str,50,y,PETSC_FALSE);
 66:   PetscPrintf(PETSC_COMM_WORLD,"  f(%s)=%s\n",strx,str);
 67:   SlepcSNPrintfScalar(str,50,yp,PETSC_FALSE);
 68:   PetscPrintf(PETSC_COMM_WORLD,"  f'(%s)=%s\n",strx,str);

 70:   /* rational p(x)/q(x) */
 71:   na = 2; nb = 3;
 72:   p[0] = -3.1; p[1] = 1.1;
 73:   q[0] = 1.0; q[1] = -2.0; q[2] = 3.5;
 74:   FNSetType(fn,FNRATIONAL);
 75:   FNRationalSetNumerator(fn,na,p);
 76:   FNRationalSetDenominator(fn,nb,q);
 77:   FNView(fn,NULL);
 78:   x = 2.2;
 79:   SlepcSNPrintfScalar(strx,50,x,PETSC_FALSE);
 80:   FNEvaluateFunction(fn,x,&y);
 81:   FNEvaluateDerivative(fn,x,&yp);
 82:   SlepcSNPrintfScalar(str,50,y,PETSC_FALSE);
 83:   PetscPrintf(PETSC_COMM_WORLD,"  f(%s)=%s\n",strx,str);
 84:   SlepcSNPrintfScalar(str,50,yp,PETSC_FALSE);
 85:   PetscPrintf(PETSC_COMM_WORLD,"  f'(%s)=%s\n",strx,str);

 87:   /* constant */
 88:   FNSetType(fn,FNRATIONAL);
 89:   FNRationalSetNumerator(fn,1,&five);
 90:   FNRationalSetDenominator(fn,0,NULL);  /* reset previous values */
 91:   FNView(fn,NULL);
 92:   x = 2.2;
 93:   SlepcSNPrintfScalar(strx,50,x,PETSC_FALSE);
 94:   FNEvaluateFunction(fn,x,&y);
 95:   FNEvaluateDerivative(fn,x,&yp);
 96:   SlepcSNPrintfScalar(str,50,y,PETSC_FALSE);
 97:   PetscPrintf(PETSC_COMM_WORLD,"  f(%s)=%s\n",strx,str);
 98:   SlepcSNPrintfScalar(str,50,yp,PETSC_FALSE);
 99:   PetscPrintf(PETSC_COMM_WORLD,"  f'(%s)=%s\n",strx,str);

101:   FNDestroy(&fn);
102:   SlepcFinalize();
103:   return 0;
104: }