Actual source code: test5.c

slepc-3.17.2 2022-08-09
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 PEP view and monitor functionality.\n\n";

 13: #include <slepcpep.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat            A[3];
 18:   PEP            pep;
 19:   Vec            xr,xi;
 20:   PetscScalar    kr,ki;
 21:   PetscInt       n=6,Istart,Iend,i,nconv,its;
 22:   PetscReal      errest;
 23:   PetscBool      checkfile;
 24:   char           filename[PETSC_MAX_PATH_LEN];
 25:   PetscViewer    viewer;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);
 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscPrintf(PETSC_COMM_WORLD,"\nPEP of diagonal problem, n=%" PetscInt_FMT "\n\n",n);

 31:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 32:         Generate the matrices
 33:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 34:   MatCreate(PETSC_COMM_WORLD,&A[0]);
 35:   MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n);
 36:   MatSetFromOptions(A[0]);
 37:   MatSetUp(A[0]);
 38:   MatGetOwnershipRange(A[0],&Istart,&Iend);
 39:   for (i=Istart;i<Iend;i++) MatSetValue(A[0],i,i,i+1,INSERT_VALUES);
 40:   MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY);
 41:   MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY);

 43:   MatCreate(PETSC_COMM_WORLD,&A[1]);
 44:   MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,n,n);
 45:   MatSetFromOptions(A[1]);
 46:   MatSetUp(A[1]);
 47:   for (i=Istart;i<Iend;i++) MatSetValue(A[1],i,i,-1.5,INSERT_VALUES);
 48:   MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY);
 49:   MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY);

 51:   MatCreate(PETSC_COMM_WORLD,&A[2]);
 52:   MatSetSizes(A[2],PETSC_DECIDE,PETSC_DECIDE,n,n);
 53:   MatSetFromOptions(A[2]);
 54:   MatSetUp(A[2]);
 55:   for (i=Istart;i<Iend;i++) MatSetValue(A[2],i,i,-1.0/(i+1),INSERT_VALUES);
 56:   MatAssemblyBegin(A[2],MAT_FINAL_ASSEMBLY);
 57:   MatAssemblyEnd(A[2],MAT_FINAL_ASSEMBLY);

 59:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 60:                      Create the PEP solver
 61:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 62:   PEPCreate(PETSC_COMM_WORLD,&pep);
 63:   PetscObjectSetName((PetscObject)pep,"pep");
 64:   PEPSetOperators(pep,3,A);
 65:   PEPSetFromOptions(pep);

 67:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 68:                 Solve the eigensystem and display solution
 69:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 70:   PEPSolve(pep);
 71:   PEPGetConverged(pep,&nconv);
 72:   PEPGetIterationNumber(pep,&its);
 73:   PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT " converged eigenpairs after %" PetscInt_FMT " iterations\n",nconv,its);
 74:   if (nconv>0) {
 75:     MatCreateVecs(A[0],&xr,&xi);
 76:     PEPGetEigenpair(pep,0,&kr,&ki,xr,xi);
 77:     VecDestroy(&xr);
 78:     VecDestroy(&xi);
 79:     PEPGetErrorEstimate(pep,0,&errest);
 80:   }
 81:   PEPErrorView(pep,PEP_ERROR_RELATIVE,NULL);

 83:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 84:                    Check file containing the eigenvalues
 85:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 86:   PetscOptionsGetString(NULL,NULL,"-checkfile",filename,sizeof(filename),&checkfile);
 87:   if (checkfile) {
 88: #if defined(PETSC_HAVE_COMPLEX)
 89:     PetscComplex *eigs,eval;
 90:     PetscMalloc1(nconv,&eigs);
 91:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
 92:     PetscViewerBinaryRead(viewer,eigs,nconv,NULL,PETSC_COMPLEX);
 93:     PetscViewerDestroy(&viewer);
 94:     for (i=0;i<nconv;i++) {
 95:       PEPGetEigenpair(pep,i,&kr,&ki,NULL,NULL);
 96: #if defined(PETSC_USE_COMPLEX)
 97:       eval = kr;
 98: #else
 99:       eval = PetscCMPLX(kr,ki);
100: #endif
102:     }
103:     PetscFree(eigs);
104: #else
105:     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"The -checkfile option requires C99 complex numbers");
106: #endif
107:   }

109:   PEPDestroy(&pep);
110:   MatDestroy(&A[0]);
111:   MatDestroy(&A[1]);
112:   MatDestroy(&A[2]);
113:   SlepcFinalize();
114:   return 0;
115: }

117: /*TEST

119:    test:
120:       suffix: 1
121:       args: -pep_error_backward ::ascii_info_detail -pep_largest_real -pep_view_values -pep_monitor_conv -pep_error_absolute ::ascii_matlab -pep_monitor_all -pep_converged_reason -pep_view
122:       requires: !single
123:       filter: grep -v "tolerance" | grep -v "problem type" | sed -e "s/[+-]0\.0*i//g" -e "s/\([0-9]\.[5]*\)[+-][0-9]\.[0-9]*e-[0-9]*i/\\1/g" -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g"

125:    test:
126:       suffix: 2
127:       args: -n 12 -pep_largest_real -pep_monitor -pep_view_values ::ascii_matlab
128:       requires: double
129:       filter: sed -e "s/[+-][0-9]\.[0-9]*e-[0-9]*i//" -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g" -e "s/5\.\([49]\)999999[0-9]*e+00/5.\\1999999999999999e+00/"

131:    test:
132:       suffix: 3
133:       args: -pep_nev 4 -pep_view_values binary:myvalues.bin -checkfile myvalues.bin
134:       requires: double c99_complex

136:    test:
137:       suffix: 4
138:       args: -pep_nev 4 -pep_ncv 10 -pep_refine -pep_conv_norm -pep_extract none -pep_scale scalar -pep_view -pep_monitor -pep_error_relative ::ascii_info_detail
139:       requires: double !complex
140:       filter: grep -v "tolerance" | sed -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g"

142:    test:
143:       suffix: 5
144:       args: -n 12 -pep_largest_real -pep_monitor draw::draw_lg -pep_monitor_all draw::draw_lg -pep_view_values draw -draw_save myeigen.ppm -draw_virtual
145:       requires: x double

147: TEST*/