Actual source code: test2.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 interface functions.\n\n";

 13: #include <slepcnep.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat                  A[3],B;      /* problem matrices */
 18:   FN                   f[3],g;      /* problem functions */
 19:   NEP                  nep;         /* eigenproblem solver context */
 20:   DS                   ds;
 21:   RG                   rg;
 22:   PetscReal            tol;
 23:   PetscScalar          coeffs[2],target;
 24:   PetscInt             n=20,i,its,nev,ncv,mpd,Istart,Iend,nterm;
 25:   PetscBool            twoside;
 26:   NEPWhich             which;
 27:   NEPConvergedReason   reason;
 28:   NEPType              type;
 29:   NEPRefine            refine;
 30:   NEPRefineScheme      rscheme;
 31:   NEPConv              conv;
 32:   NEPStop              stop;
 33:   NEPProblemType       ptype;
 34:   MatStructure         mstr;
 35:   PetscViewerAndFormat *vf;

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

 41:   /*
 42:      Matrices
 43:   */
 44:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A[0]));
 45:   PetscCall(MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n));
 46:   PetscCall(MatSetFromOptions(A[0]));
 47:   PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
 48:   for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(A[0],i,i,i+1,INSERT_VALUES));
 49:   PetscCall(MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY));
 50:   PetscCall(MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY));

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

 60:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A[2]));
 61:   PetscCall(MatSetSizes(A[2],PETSC_DECIDE,PETSC_DECIDE,n,n));
 62:   PetscCall(MatSetFromOptions(A[2]));
 63:   PetscCall(MatGetOwnershipRange(A[1],&Istart,&Iend));
 64:   for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(A[2],i,i,n/(PetscReal)(i+1),INSERT_VALUES));
 65:   PetscCall(MatAssemblyBegin(A[2],MAT_FINAL_ASSEMBLY));
 66:   PetscCall(MatAssemblyEnd(A[2],MAT_FINAL_ASSEMBLY));

 68:   /*
 69:      Functions: f0=-lambda, f1=1.0, f2=sqrt(lambda)
 70:   */
 71:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[0]));
 72:   PetscCall(FNSetType(f[0],FNRATIONAL));
 73:   coeffs[0] = -1.0; coeffs[1] = 0.0;
 74:   PetscCall(FNRationalSetNumerator(f[0],2,coeffs));

 76:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[1]));
 77:   PetscCall(FNSetType(f[1],FNRATIONAL));
 78:   coeffs[0] = 1.0;
 79:   PetscCall(FNRationalSetNumerator(f[1],1,coeffs));

 81:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[2]));
 82:   PetscCall(FNSetType(f[2],FNSQRT));

 84:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 85:              Create eigensolver and test interface functions
 86:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 87:   PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
 88:   PetscCall(NEPSetSplitOperator(nep,3,A,f,SAME_NONZERO_PATTERN));
 89:   PetscCall(NEPGetSplitOperatorInfo(nep,&nterm,&mstr));
 90:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Nonlinear function with %" PetscInt_FMT " terms, with %s nonzero pattern\n",nterm,MatStructures[mstr]));
 91:   PetscCall(NEPGetSplitOperatorTerm(nep,0,&B,&g));
 92:   PetscCall(MatView(B,NULL));
 93:   PetscCall(FNView(g,NULL));

 95:   PetscCall(NEPSetType(nep,NEPRII));
 96:   PetscCall(NEPGetType(nep,&type));
 97:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Type set to %s\n",type));
 98:   PetscCall(NEPGetTwoSided(nep,&twoside));
 99:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Two-sided flag = %s\n",twoside?"true":"false"));

101:   PetscCall(NEPGetProblemType(nep,&ptype));
102:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Problem type before changing = %d",(int)ptype));
103:   PetscCall(NEPSetProblemType(nep,NEP_RATIONAL));
104:   PetscCall(NEPGetProblemType(nep,&ptype));
105:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d.\n",(int)ptype));

107:   PetscCall(NEPSetRefine(nep,NEP_REFINE_SIMPLE,1,1e-9,2,NEP_REFINE_SCHEME_EXPLICIT));
108:   PetscCall(NEPGetRefine(nep,&refine,NULL,&tol,&its,&rscheme));
109:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Refinement: %s, tol=%g, its=%" PetscInt_FMT ", scheme=%s\n",NEPRefineTypes[refine],(double)tol,its,NEPRefineSchemes[rscheme]));

111:   PetscCall(NEPSetTarget(nep,1.1));
112:   PetscCall(NEPGetTarget(nep,&target));
113:   PetscCall(NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE));
114:   PetscCall(NEPGetWhichEigenpairs(nep,&which));
115:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Which = %d, target = %g\n",(int)which,(double)PetscRealPart(target)));

117:   PetscCall(NEPSetDimensions(nep,1,12,PETSC_DEFAULT));
118:   PetscCall(NEPGetDimensions(nep,&nev,&ncv,&mpd));
119:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Dimensions: nev=%" PetscInt_FMT ", ncv=%" PetscInt_FMT ", mpd=%" PetscInt_FMT "\n",nev,ncv,mpd));

121:   PetscCall(NEPSetTolerances(nep,1.0e-6,200));
122:   PetscCall(NEPGetTolerances(nep,&tol,&its));
123:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Tolerance = %.6f, max_its = %" PetscInt_FMT "\n",(double)tol,its));

125:   PetscCall(NEPSetConvergenceTest(nep,NEP_CONV_ABS));
126:   PetscCall(NEPGetConvergenceTest(nep,&conv));
127:   PetscCall(NEPSetStoppingTest(nep,NEP_STOP_BASIC));
128:   PetscCall(NEPGetStoppingTest(nep,&stop));
129:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Convergence test = %d, stopping test = %d\n",(int)conv,(int)stop));

131:   PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf));
132:   PetscCall(NEPMonitorSet(nep,(PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))NEPMonitorFirst,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy));
133:   PetscCall(NEPMonitorCancel(nep));

135:   PetscCall(NEPGetDS(nep,&ds));
136:   PetscCall(DSView(ds,NULL));
137:   PetscCall(NEPSetFromOptions(nep));

139:   PetscCall(NEPGetRG(nep,&rg));
140:   PetscCall(RGView(rg,NULL));

142:   PetscCall(NEPSolve(nep));
143:   PetscCall(NEPGetConvergedReason(nep,&reason));
144:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d\n",(int)reason));

146:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147:                     Display solution and clean up
148:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149:   PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
150:   PetscCall(NEPDestroy(&nep));
151:   PetscCall(MatDestroy(&A[0]));
152:   PetscCall(MatDestroy(&A[1]));
153:   PetscCall(MatDestroy(&A[2]));
154:   PetscCall(FNDestroy(&f[0]));
155:   PetscCall(FNDestroy(&f[1]));
156:   PetscCall(FNDestroy(&f[2]));
157:   PetscCall(SlepcFinalize());
158:   return 0;
159: }

161: /*TEST

163:    test:
164:       suffix: 1
165:       args: -nep_view
166:       filter: grep -v tolerance

168: TEST*/