Actual source code: test3.c

slepc-3.13.4 2020-09-02
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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 MFN interface functions.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n\n";

 15:  #include <slepcmfn.h>

 17: int main(int argc,char **argv)
 18: {
 19:   Mat                  A,B;
 20:   MFN                  mfn;
 21:   FN                   f;
 22:   MFNConvergedReason   reason;
 23:   MFNType              type;
 24:   PetscReal            norm,tol;
 25:   Vec                  v,y;
 26:   PetscInt             N,n=4,Istart,Iend,i,j,II,ncv,its,maxit;
 27:   PetscBool            flg,testprefix=PETSC_FALSE;
 28:   const char           *prefix;
 29:   PetscErrorCode       ierr;
 30:   PetscViewerAndFormat *vf;

 32:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 33:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 34:   N = n*n;
 35:   PetscPrintf(PETSC_COMM_WORLD,"\nSquare root of Laplacian y=sqrt(A)*e_1, N=%D (%Dx%D grid)\n\n",N,n,n);
 36:   PetscOptionsGetBool(NULL,NULL,"-test_prefix",&testprefix,NULL);

 38:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 39:                  Compute the discrete 2-D Laplacian, A
 40:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 42:   MatCreate(PETSC_COMM_WORLD,&A);
 43:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 44:   MatSetFromOptions(A);
 45:   MatSetUp(A);

 47:   MatGetOwnershipRange(A,&Istart,&Iend);
 48:   for (II=Istart;II<Iend;II++) {
 49:     i = II/n; j = II-i*n;
 50:     if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
 51:     if (i<n-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
 52:     if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
 53:     if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
 54:     MatSetValue(A,II,II,4.0,INSERT_VALUES);
 55:   }

 57:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 58:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 59:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);

 61:   MatCreateVecs(A,NULL,&v);
 62:   VecSetValue(v,0,1.0,INSERT_VALUES);
 63:   VecAssemblyBegin(v);
 64:   VecAssemblyEnd(v);
 65:   VecDuplicate(v,&y);

 67:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 68:              Create the solver, set the matrix and the function
 69:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 70:   MFNCreate(PETSC_COMM_WORLD,&mfn);
 71:   MFNSetOperator(mfn,A);
 72:   MFNGetFN(mfn,&f);
 73:   FNSetType(f,FNSQRT);

 75:   MFNSetType(mfn,MFNKRYLOV);
 76:   MFNGetType(mfn,&type);
 77:   PetscPrintf(PETSC_COMM_WORLD," Type set to %s\n",type);

 79:   /* test prefix usage */
 80:   if (testprefix) {
 81:     MFNSetOptionsPrefix(mfn,"check_");
 82:     MFNAppendOptionsPrefix(mfn,"myprefix_");
 83:     MFNGetOptionsPrefix(mfn,&prefix);
 84:     PetscPrintf(PETSC_COMM_WORLD," MFN prefix is currently: %s\n",prefix);
 85:   }

 87:   /* test some interface functions */
 88:   MFNGetOperator(mfn,&B);
 89:   MatView(B,PETSC_VIEWER_STDOUT_WORLD);
 90:   MFNSetTolerances(mfn,1e-4,500);
 91:   MFNSetDimensions(mfn,6);
 92:   MFNSetErrorIfNotConverged(mfn,PETSC_TRUE);
 93:   /* test monitors */
 94:   PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);
 95:   MFNMonitorSet(mfn,(PetscErrorCode (*)(MFN,PetscInt,PetscReal,void*))MFNMonitorDefault,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
 96:   /* MFNMonitorCancel(mfn); */
 97:   MFNSetFromOptions(mfn);

 99:   /* query properties and print them */
100:   MFNGetTolerances(mfn,&tol,&maxit);
101:   PetscPrintf(PETSC_COMM_WORLD," Tolerance: %g, max iterations: %D\n",(double)tol,maxit);
102:   MFNGetDimensions(mfn,&ncv);
103:   PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %D\n",ncv);
104:   MFNGetErrorIfNotConverged(mfn,&flg);
105:   if (flg) { PetscPrintf(PETSC_COMM_WORLD," Erroring out if convergence fails\n"); }

107:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108:                            Solve  y=sqrt(A)*v
109:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

111:   MFNSolve(mfn,v,y);
112:   MFNGetConvergedReason(mfn,&reason);
113:   MFNGetIterationNumber(mfn,&its);
114:   PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d\n",(int)reason);
115:   /* PetscPrintf(PETSC_COMM_WORLD," its = %D\n",its); */
116:   VecNorm(y,NORM_2,&norm);
117:   PetscPrintf(PETSC_COMM_WORLD," sqrt(A)*v has norm %g\n",(double)norm);

119:   /*
120:      Free work space
121:   */
122:   MFNDestroy(&mfn);
123:   MatDestroy(&A);
124:   VecDestroy(&v);
125:   VecDestroy(&y);
126:   SlepcFinalize();
127:   return ierr;
128: }

130: /*TEST

132:    test:
133:       suffix: 1
134:       args: -mfn_monitor_cancel -mfn_converged_reason -mfn_view -log_exclude mfn,bv,fn

136:    test:
137:       suffix: 2
138:       args: -test_prefix -check_myprefix_mfn_monitor
139:       filter: sed -e "s/estimate [0-9]\.[0-9]*e[+-]\([0-9]*\)/estimate (removed)/g" | sed -e "s/4.0[0-9]*e-10/4.03e-10/"

141: TEST*/