Actual source code: test8.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 NEP view and monitor functionality.\n\n";

 13: #include <slepcnep.h>

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

 28:   PetscFunctionBeginUser;
 29:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 30:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n));

 32:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 33:         Generate the matrices
 34:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 36:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A[0]));
 37:   PetscCall(MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n));
 38:   PetscCall(MatSetFromOptions(A[0]));
 39:   PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
 40:   for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(A[0],i,i,i+1,INSERT_VALUES));
 41:   PetscCall(MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY));
 42:   PetscCall(MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY));

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

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

 58:   /*
 59:      Functions: f0=1.0, f1=lambda, f2=lambda^2
 60:   */
 61:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[0]));
 62:   PetscCall(FNSetType(f[0],FNRATIONAL));
 63:   coeffs[0] = 1.0;
 64:   PetscCall(FNRationalSetNumerator(f[0],1,coeffs));

 66:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[1]));
 67:   PetscCall(FNSetType(f[1],FNRATIONAL));
 68:   coeffs[0] = 1.0; coeffs[1] = 0.0;
 69:   PetscCall(FNRationalSetNumerator(f[1],2,coeffs));

 71:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[2]));
 72:   PetscCall(FNSetType(f[2],FNRATIONAL));
 73:   coeffs[0] = 1.0; coeffs[1] = 0.0; coeffs[2] = 0.0;
 74:   PetscCall(FNRationalSetNumerator(f[2],3,coeffs));

 76:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 77:                       Create the NEP solver
 78:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 79:   PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
 80:   PetscCall(PetscObjectSetName((PetscObject)nep,"nep"));
 81:   PetscCall(NEPSetSplitOperator(nep,3,A,f,SAME_NONZERO_PATTERN));
 82:   PetscCall(NEPSetTarget(nep,1.1));
 83:   PetscCall(NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE));
 84:   PetscCall(NEPSetFromOptions(nep));

 86:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 87:                 Solve the eigensystem and display solution
 88:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 89:   PetscCall(NEPSolve(nep));
 90:   PetscCall(NEPGetConverged(nep,&nconv));
 91:   PetscCall(NEPGetIterationNumber(nep,&its));
 92:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," %" PetscInt_FMT " converged eigenpairs after %" PetscInt_FMT " iterations\n",nconv,its));
 93:   if (nconv>0) {
 94:     PetscCall(MatCreateVecs(A[0],&xr,&xi));
 95:     PetscCall(NEPGetEigenpair(nep,0,&kr,&ki,xr,xi));
 96:     PetscCall(VecDestroy(&xr));
 97:     PetscCall(VecDestroy(&xi));
 98:     PetscCall(NEPGetErrorEstimate(nep,0,&errest));
 99:   }
100:   PetscCall(NEPErrorView(nep,NEP_ERROR_BACKWARD,NULL));

102:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103:                    Check file containing the eigenvalues
104:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105:   PetscCall(PetscOptionsGetString(NULL,NULL,"-checkfile",filename,sizeof(filename),&checkfile));
106:   if (checkfile) {
107: #if defined(PETSC_HAVE_COMPLEX)
108:     PetscComplex *eigs,eval;
109:     PetscCall(PetscMalloc1(nconv,&eigs));
110:     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer));
111:     PetscCall(PetscViewerBinaryRead(viewer,eigs,nconv,NULL,PETSC_COMPLEX));
112:     PetscCall(PetscViewerDestroy(&viewer));
113:     for (i=0;i<nconv;i++) {
114:       PetscCall(NEPGetEigenpair(nep,i,&kr,&ki,NULL,NULL));
115: #if defined(PETSC_USE_COMPLEX)
116:       eval = kr;
117: #else
118:       eval = PetscCMPLX(kr,ki);
119: #endif
120:       PetscCheck(eval==eigs[i],PETSC_COMM_WORLD,PETSC_ERR_FILE_UNEXPECTED,"Eigenvalues in the file do not match");
121:     }
122:     PetscCall(PetscFree(eigs));
123: #else
124:     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"The -checkfile option requires C99 complex numbers");
125: #endif
126:   }

128:   PetscCall(NEPDestroy(&nep));
129:   PetscCall(MatDestroy(&A[0]));
130:   PetscCall(MatDestroy(&A[1]));
131:   PetscCall(MatDestroy(&A[2]));
132:   PetscCall(FNDestroy(&f[0]));
133:   PetscCall(FNDestroy(&f[1]));
134:   PetscCall(FNDestroy(&f[2]));
135:   PetscCall(SlepcFinalize());
136:   return 0;
137: }

139: /*TEST

141:    test:
142:       suffix: 1
143:       args: -nep_type slp -nep_target -.5 -nep_error_backward ::ascii_info_detail -nep_view_values -nep_error_absolute ::ascii_matlab -nep_monitor_all -nep_converged_reason -nep_view
144:       filter: grep -v "tolerance" | grep -v "problem type" | sed -e "s/[+-]0\.0*i//g" -e "s/+0i//" -e "s/[+-][0-9]\.[0-9]*e-[0-9]*i//g" -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
145:       requires: double

147:    test:
148:       suffix: 2
149:       args: -nep_type rii -nep_target -.5 -nep_rii_hermitian -nep_monitor -nep_view_values ::ascii_matlab
150:       filter: sed -e "s/[+-][0-9]\.[0-9]*e-[0-9]*i//" -e "s/([0-9]\.[0-9]*e[+-]\([0-9]*\))/(removed)/g"
151:       requires: double

153:    test:
154:       suffix: 3
155:       args: -nep_type slp -nep_nev 4 -nep_view_values binary:myvalues.bin -checkfile myvalues.bin -nep_error_relative ::ascii_matlab
156:       filter: sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
157:       requires: double c99_complex

159:    test:
160:       suffix: 4
161:       args: -nep_type slp -nep_nev 4 -nep_monitor draw::draw_lg -nep_monitor_all draw::draw_lg -nep_view_values draw -draw_save myeigen.ppm -draw_virtual
162:       requires: x double

164: TEST*/