Actual source code: test1.c

slepc-3.17.1 2022-04-11
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[] = "Computes exp(t*A)*v for a matrix loaded from file.\n\n"
 12:   "The command line options are:\n"
 13:   "  -t <sval>, where <sval> = scalar value that multiplies the argument.\n"
 14:   "  -file <filename>, where <filename> = matrix file in PETSc binary form.\n\n";

 16: #include <slepcmfn.h>

 18: int main(int argc,char **argv)
 19: {
 20:   Mat                A;           /* problem matrix */
 21:   MFN                mfn;
 22:   FN                 f;
 23:   PetscReal          norm;
 24:   PetscScalar        t=2.0;
 25:   Vec                v,y;
 26:   PetscViewer        viewer;
 27:   PetscBool          flg;
 28:   char               filename[PETSC_MAX_PATH_LEN];
 29:   MFNConvergedReason reason;

 31:   SlepcInitialize(&argc,&argv,(char*)0,help);

 33:   PetscOptionsGetScalar(NULL,NULL,"-t",&t,NULL);
 34:   PetscPrintf(PETSC_COMM_WORLD,"\nMatrix exponential y=exp(t*A)*e, loaded from file\n\n");

 36:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 37:                 Load matrix A from binary file
 38:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 40:   PetscOptionsGetString(NULL,NULL,"-file",filename,sizeof(filename),&flg);

 43: #if defined(PETSC_USE_COMPLEX)
 44:   PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrix from a binary file...\n");
 45: #else
 46:   PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrix from a binary file...\n");
 47: #endif
 48:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
 49:   MatCreate(PETSC_COMM_WORLD,&A);
 50:   MatSetFromOptions(A);
 51:   MatLoad(A,viewer);
 52:   PetscViewerDestroy(&viewer);

 54:   /* set v = ones(n,1) */
 55:   MatCreateVecs(A,NULL,&y);
 56:   MatCreateVecs(A,NULL,&v);
 57:   VecSet(v,1.0);

 59:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 60:                 Create the solver and set various options
 61:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 63:   MFNCreate(PETSC_COMM_WORLD,&mfn);
 64:   MFNSetOperator(mfn,A);
 65:   MFNGetFN(mfn,&f);
 66:   FNSetType(f,FNEXP);
 67:   FNSetScale(f,t,1.0);
 68:   MFNSetFromOptions(mfn);

 70:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 71:                       Solve the problem, y=exp(t*A)*v
 72:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 74:   MFNSolve(mfn,v,y);
 75:   MFNGetConvergedReason(mfn,&reason);
 77:   VecNorm(y,NORM_2,&norm);
 78:   PetscPrintf(PETSC_COMM_WORLD," Computed vector at time t=%.4g has norm %g\n",(double)PetscRealPart(t),(double)norm);

 80:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 81:                       Solve the problem, y=exp(t*A^T)*v
 82:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 84:   MFNSolveTranspose(mfn,v,y);
 85:   MFNGetConvergedReason(mfn,&reason);
 87:   VecNorm(y,NORM_2,&norm);
 88:   PetscPrintf(PETSC_COMM_WORLD," With transpose: computed vector has norm %g\n\n",(double)norm);

 90:   /*
 91:      Free work space
 92:   */
 93:   MFNDestroy(&mfn);
 94:   MatDestroy(&A);
 95:   VecDestroy(&v);
 96:   VecDestroy(&y);
 97:   SlepcFinalize();
 98:   return 0;
 99: }

101: /*TEST

103:    testset:
104:       args: -file ${DATAFILESPATH}/matrices/real/bfw782a.petsc -mfn_type {{krylov expokit}} -t 0.05
105:       requires: double !complex datafilespath !defined(PETSC_USE_64BIT_INDICES)
106:       output_file: output/test1_1.out
107:       test:
108:          suffix: 1
109:       test:
110:          suffix: 1_cuda
111:          args: -mat_type aijcusparse
112:          requires: cuda

114:    testset:
115:       args: -file ${DATAFILESPATH}/matrices/complex/qc324.petsc -mfn_type {{krylov expokit}}
116:       requires: double complex datafilespath !defined(PETSC_USE_64BIT_INDICES)
117:       output_file: output/test1_2.out
118:       test:
119:          suffix: 2
120:       test:
121:          suffix: 2_cuda
122:          args: -mat_type aijcusparse
123:          requires: cuda

125: TEST*/