Actual source code: test5.c

slepc-3.15.0 2021-03-31
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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:   PetscComplex   *eigs,eval;
 22:   PetscInt       n=6,Istart,Iend,i,nconv,its;
 23:   PetscReal      errest;
 24:   PetscBool      checkfile;
 25:   char           filename[PETSC_MAX_PATH_LEN];
 26:   PetscViewer    viewer;

 29:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 30:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 31:   PetscPrintf(PETSC_COMM_WORLD,"\nPEP of diagonal problem, n=%D\n\n",n);

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

 47:   MatCreate(PETSC_COMM_WORLD,&A[1]);
 48:   MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,n,n);
 49:   MatSetFromOptions(A[1]);
 50:   MatSetUp(A[1]);
 51:   for (i=Istart;i<Iend;i++) {
 52:     MatSetValue(A[1],i,i,-1.5,INSERT_VALUES);
 53:   }
 54:   MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY);
 55:   MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY);

 57:   MatCreate(PETSC_COMM_WORLD,&A[2]);
 58:   MatSetSizes(A[2],PETSC_DECIDE,PETSC_DECIDE,n,n);
 59:   MatSetFromOptions(A[2]);
 60:   MatSetUp(A[2]);
 61:   for (i=Istart;i<Iend;i++) {
 62:     MatSetValue(A[2],i,i,-1.0/(i+1),INSERT_VALUES);
 63:   }
 64:   MatAssemblyBegin(A[2],MAT_FINAL_ASSEMBLY);
 65:   MatAssemblyEnd(A[2],MAT_FINAL_ASSEMBLY);

 67:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 68:                      Create the PEP solver
 69:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 70:   PEPCreate(PETSC_COMM_WORLD,&pep);
 71:   PetscObjectSetName((PetscObject)pep,"pep");
 72:   PEPSetOperators(pep,3,A);
 73:   PEPSetFromOptions(pep);

 75:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 76:                 Solve the eigensystem and display solution
 77:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 78:   PEPSolve(pep);
 79:   PEPGetConverged(pep,&nconv);
 80:   PEPGetIterationNumber(pep,&its);
 81:   PetscPrintf(PETSC_COMM_WORLD," %D converged eigenpairs after %D iterations\n",nconv,its);
 82:   if (nconv>0) {
 83:     MatCreateVecs(A[0],&xr,&xi);
 84:     PEPGetEigenpair(pep,0,&kr,&ki,xr,xi);
 85:     VecDestroy(&xr);
 86:     VecDestroy(&xi);
 87:     PEPGetErrorEstimate(pep,0,&errest);
 88:   }
 89:   PEPErrorView(pep,PEP_ERROR_RELATIVE,NULL);

 91:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 92:                    Check file containing the eigenvalues
 93:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 94:   PetscOptionsGetString(NULL,NULL,"-checkfile",filename,sizeof(filename),&checkfile);
 95:   if (checkfile) {
 96:     PetscMalloc1(nconv,&eigs);
 97:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
 98:     PetscViewerBinaryRead(viewer,eigs,nconv,NULL,PETSC_COMPLEX);
 99:     PetscViewerDestroy(&viewer);
100:     for (i=0;i<nconv;i++) {
101:       PEPGetEigenpair(pep,i,&kr,&ki,NULL,NULL);
102: #if defined(PETSC_USE_COMPLEX)
103:       eval = kr;
104: #else
105:       eval = PetscCMPLX(kr,ki);
106: #endif
107:       if (eval!=eigs[i]) SETERRQ(PETSC_COMM_WORLD,1,"Eigenvalues in the file do not match");
108:     }
109:     PetscFree(eigs);
110:   }

112:   PEPDestroy(&pep);
113:   MatDestroy(&A[0]);
114:   MatDestroy(&A[1]);
115:   MatDestroy(&A[2]);
116:   SlepcFinalize();
117:   return ierr;
118: }

120: /*TEST

122:    test:
123:       suffix: 1
124:       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
125:       requires: !single
126:       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"

128:    test:
129:       suffix: 2
130:       args: -n 12 -pep_largest_real -pep_monitor -pep_view_values ::ascii_matlab
131:       requires: double
132:       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/"

134:    test:
135:       suffix: 3
136:       args: -pep_nev 4 -pep_view_values binary:myvalues.bin -checkfile myvalues.bin
137:       requires: double

139:    test:
140:       suffix: 4
141:       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
142:       requires: double !complex
143:       filter: grep -v "tolerance" | sed -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g"

145:    test:
146:       suffix: 5
147:       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
148:       requires: double

150: TEST*/