Actual source code: test2.c

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

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Test exponential function.\n\n";

 13: #include <slepcfn.h>

 15: int main(int argc,char **argv)
 16: {
 17:   FN             fn,fncopy;
 18:   PetscScalar    x,y,yp,tau,eta,alpha,beta;
 19:   char           strx[50],str[50];

 21:   SlepcInitialize(&argc,&argv,(char*)0,help);
 22:   FNCreate(PETSC_COMM_WORLD,&fn);
 23:   FNSetFromOptions(fn);

 25:   /* plain exponential exp(x) */
 26:   FNSetType(fn,FNEXP);
 27:   FNView(fn,NULL);
 28:   x = 2.2;
 29:   SlepcSNPrintfScalar(strx,sizeof(strx),x,PETSC_FALSE);
 30:   FNEvaluateFunction(fn,x,&y);
 31:   FNEvaluateDerivative(fn,x,&yp);
 32:   SlepcSNPrintfScalar(str,sizeof(str),y,PETSC_FALSE);
 33:   PetscPrintf(PETSC_COMM_WORLD,"  f(%s)=%s\n",strx,str);
 34:   SlepcSNPrintfScalar(str,sizeof(str),yp,PETSC_FALSE);
 35:   PetscPrintf(PETSC_COMM_WORLD,"  f'(%s)=%s\n",strx,str);

 37:   /* exponential with scaling factors eta*exp(tau*x) */
 38:   FNSetType(fn,FNEXP);
 39:   tau = -0.2;
 40:   eta = 1.3;
 41:   FNSetScale(fn,tau,eta);
 42:   FNView(fn,NULL);
 43:   x = 2.2;
 44:   SlepcSNPrintfScalar(strx,sizeof(strx),x,PETSC_FALSE);
 45:   FNEvaluateFunction(fn,x,&y);
 46:   FNEvaluateDerivative(fn,x,&yp);
 47:   SlepcSNPrintfScalar(str,sizeof(str),y,PETSC_FALSE);
 48:   PetscPrintf(PETSC_COMM_WORLD,"  f(%s)=%s\n",strx,str);
 49:   SlepcSNPrintfScalar(str,sizeof(str),yp,PETSC_FALSE);
 50:   PetscPrintf(PETSC_COMM_WORLD,"  f'(%s)=%s\n",strx,str);

 52:   /* test FNDuplicate */
 53:   FNDuplicate(fn,PetscObjectComm((PetscObject)fn),&fncopy);

 55:   /* test FNGetScale */
 56:   FNGetScale(fncopy,&alpha,&beta);
 57:   PetscPrintf(PETSC_COMM_WORLD,"Parameters:\n - alpha: ");
 58:   SlepcSNPrintfScalar(str,sizeof(str),alpha,PETSC_FALSE);
 59:   PetscPrintf(PETSC_COMM_WORLD,"%s ",str);
 60:   PetscPrintf(PETSC_COMM_WORLD,"\n - beta: ");
 61:   SlepcSNPrintfScalar(str,sizeof(str),beta,PETSC_FALSE);
 62:   PetscPrintf(PETSC_COMM_WORLD,"%s ",str);
 63:   PetscPrintf(PETSC_COMM_WORLD,"\n");

 65:   FNDestroy(&fn);
 66:   FNDestroy(&fncopy);
 67:   SlepcFinalize();
 68:   return 0;
 69: }

 71: /*TEST

 73:    test:
 74:       suffix: 1
 75:       nsize: 1
 76:       filter: grep -v "computing matrix functions"

 78: TEST*/