Actual source code: test8.c

slepc-3.18.1 2022-11-02
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;

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

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

 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++) MatSetValue(A[0],i,i,i+1,INSERT_VALUES);
 42:   MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY);
 43:   MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY);

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

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

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

 69:   FNCreate(PETSC_COMM_WORLD,&f[1]);
 70:   FNSetType(f[1],FNRATIONAL);
 71:   coeffs[0] = 1.0; coeffs[1] = 0.0;
 72:   FNRationalSetNumerator(f[1],2,coeffs);

 74:   FNCreate(PETSC_COMM_WORLD,&f[2]);
 75:   FNSetType(f[2],FNRATIONAL);
 76:   coeffs[0] = 1.0; coeffs[1] = 0.0; coeffs[2] = 0.0;
 77:   FNRationalSetNumerator(f[2],3,coeffs);

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

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

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

131:   NEPDestroy(&nep);
132:   MatDestroy(&A[0]);
133:   MatDestroy(&A[1]);
134:   MatDestroy(&A[2]);
135:   FNDestroy(&f[0]);
136:   FNDestroy(&f[1]);
137:   FNDestroy(&f[2]);
138:   SlepcFinalize();
139:   return 0;
140: }

142: /*TEST

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

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

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

162:    test:
163:       suffix: 4
164:       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
165:       requires: x double

167: TEST*/