Actual source code: test5.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[] = "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;

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

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

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

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

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

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

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

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

113: /*TEST

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

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

127:    test:
128:       suffix: 3
129:       args: -pep_nev 4 -pep_view_values binary:myvalues.bin -checkfile myvalues.bin
130:       requires: double

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

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

143: TEST*/