Actual source code: test2.c

slepc-3.17.0 2022-03-31
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:   SlepcInitialize(&argc,&argv,(char*)0,help);
 38:   PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n);

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

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

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

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

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

 83:   FNCreate(PETSC_COMM_WORLD,&f[2]);
 84:   FNSetType(f[2],FNSQRT);

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

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

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

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

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

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

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

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

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

137:   NEPGetDS(nep,&ds);
138:   DSView(ds,NULL);
139:   NEPSetFromOptions(nep);

141:   NEPGetRG(nep,&rg);
142:   RGView(rg,NULL);

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

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

163: /*TEST

165:    test:
166:       suffix: 1

168: TEST*/