Actual source code: test5.c

slepc-3.21.1 2024-04-26
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:   PetscFunctionBeginUser;
 28:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 29:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 30:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nPEP of diagonal problem, n=%" PetscInt_FMT "\n\n",n));

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

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

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

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

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

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

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

115: /*TEST

117:    test:
118:       suffix: 1
119:       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
120:       requires: !single
121:       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"

123:    test:
124:       suffix: 2
125:       args: -n 12 -pep_largest_real -pep_monitor -pep_view_values ::ascii_matlab
126:       requires: double
127:       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/"

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

134:    test:
135:       suffix: 4
136:       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
137:       requires: double !complex
138:       filter: grep -v "tolerance" | sed -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g"

140:    test:
141:       suffix: 5
142:       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
143:       requires: x double

145: TEST*/