Actual source code: test2.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, 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 exponential function.\n\n";

 24: #include <slepcfn.h>

 28: int main(int argc,char **argv)
 29: {
 31:   FN             fn;
 32:   PetscScalar    x,y,yp,tau,eta,alpha,beta;
 33:   char           strx[50],str[50];

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

 38:   /* plain exponential exp(x) */
 39:   FNSetType(fn,FNEXP);
 40:   FNView(fn,NULL);
 41:   x = 2.2;
 42:   SlepcSNPrintfScalar(strx,50,x,PETSC_FALSE);
 43:   FNEvaluateFunction(fn,x,&y);
 44:   FNEvaluateDerivative(fn,x,&yp);
 45:   SlepcSNPrintfScalar(str,50,y,PETSC_FALSE);
 46:   PetscPrintf(PETSC_COMM_WORLD,"  f(%s)=%s\n",strx,str);
 47:   SlepcSNPrintfScalar(str,50,yp,PETSC_FALSE);
 48:   PetscPrintf(PETSC_COMM_WORLD,"  f'(%s)=%s\n",strx,str);

 50:   /* exponential with scaling factors eta*exp(tau*x) */
 51:   FNSetType(fn,FNEXP);
 52:   tau = -0.2;
 53:   eta = 1.3;
 54:   FNSetScale(fn,tau,eta);
 55:   FNView(fn,NULL);
 56:   x = 2.2;
 57:   SlepcSNPrintfScalar(strx,50,x,PETSC_FALSE);
 58:   FNEvaluateFunction(fn,x,&y);
 59:   FNEvaluateDerivative(fn,x,&yp);
 60:   SlepcSNPrintfScalar(str,50,y,PETSC_FALSE);
 61:   PetscPrintf(PETSC_COMM_WORLD,"  f(%s)=%s\n",strx,str);
 62:   SlepcSNPrintfScalar(str,50,yp,PETSC_FALSE);
 63:   PetscPrintf(PETSC_COMM_WORLD,"  f'(%s)=%s\n",strx,str);

 65:   /* test FNGetScale */
 66:   FNGetScale(fn,&alpha,&beta);
 67:   PetscPrintf(PETSC_COMM_WORLD,"Parameters:\n - alpha: ");
 68:   SlepcSNPrintfScalar(str,50,alpha,PETSC_FALSE);
 69:   PetscPrintf(PETSC_COMM_WORLD,"%s ",str);
 70:   PetscPrintf(PETSC_COMM_WORLD,"\n - beta: ");
 71:   SlepcSNPrintfScalar(str,50,beta,PETSC_FALSE);
 72:   PetscPrintf(PETSC_COMM_WORLD,"%s ",str);
 73:   PetscPrintf(PETSC_COMM_WORLD,"\n");

 75:   FNDestroy(&fn);
 76:   SlepcFinalize();
 77:   return ierr;
 78: }