Actual source code: test3.c

slepc-3.18.2 2023-01-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 PEP interface functions.\n\n";

 13: #include <slepcpep.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat                A[3],B;      /* problem matrices */
 18:   PEP                pep;         /* eigenproblem solver context */
 19:   ST                 st;
 20:   KSP                ksp;
 21:   DS                 ds;
 22:   PetscReal          tol,alpha;
 23:   PetscScalar        target;
 24:   PetscInt           n=20,i,its,nev,ncv,mpd,Istart,Iend,nmat;
 25:   PEPWhich           which;
 26:   PEPConvergedReason reason;
 27:   PEPType            type;
 28:   PEPExtract         extr;
 29:   PEPBasis           basis;
 30:   PEPScale           scale;
 31:   PEPRefine          refine;
 32:   PEPRefineScheme    rscheme;
 33:   PEPConv            conv;
 34:   PEPStop            stop;
 35:   PEPProblemType     ptype;
 36:   PetscViewerAndFormat *vf;

 39:   SlepcInitialize(&argc,&argv,(char*)0,help);
 40:   PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Quadratic Eigenproblem, n=%" PetscInt_FMT "\n\n",n);

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

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

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

 69:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 70:              Create eigensolver and test interface functions
 71:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 72:   PEPCreate(PETSC_COMM_WORLD,&pep);
 73:   PEPSetOperators(pep,3,A);
 74:   PEPGetNumMatrices(pep,&nmat);
 75:   PetscPrintf(PETSC_COMM_WORLD," Polynomial of degree %" PetscInt_FMT "\n",nmat-1);
 76:   PEPGetOperators(pep,0,&B);
 77:   MatView(B,NULL);

 79:   PEPSetType(pep,PEPTOAR);
 80:   PEPGetType(pep,&type);
 81:   PetscPrintf(PETSC_COMM_WORLD," Type set to %s\n",type);

 83:   PEPGetProblemType(pep,&ptype);
 84:   PetscPrintf(PETSC_COMM_WORLD," Problem type before changing = %d",(int)ptype);
 85:   PEPSetProblemType(pep,PEP_HERMITIAN);
 86:   PEPGetProblemType(pep,&ptype);
 87:   PetscPrintf(PETSC_COMM_WORLD," ... changed to %d.\n",(int)ptype);

 89:   PEPGetExtract(pep,&extr);
 90:   PetscPrintf(PETSC_COMM_WORLD," Extraction before changing = %d",(int)extr);
 91:   PEPSetExtract(pep,PEP_EXTRACT_STRUCTURED);
 92:   PEPGetExtract(pep,&extr);
 93:   PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)extr);

 95:   PEPSetScale(pep,PEP_SCALE_SCALAR,.1,NULL,NULL,5,1.0);
 96:   PEPGetScale(pep,&scale,&alpha,NULL,NULL,&its,NULL);
 97:   PetscPrintf(PETSC_COMM_WORLD," Scaling: %s, alpha=%g, its=%" PetscInt_FMT "\n",PEPScaleTypes[scale],(double)alpha,its);

 99:   PEPSetBasis(pep,PEP_BASIS_CHEBYSHEV1);
100:   PEPGetBasis(pep,&basis);
101:   PetscPrintf(PETSC_COMM_WORLD," Polynomial basis: %s\n",PEPBasisTypes[basis]);

103:   PEPSetRefine(pep,PEP_REFINE_SIMPLE,1,1e-9,2,PEP_REFINE_SCHEME_SCHUR);
104:   PEPGetRefine(pep,&refine,NULL,&tol,&its,&rscheme);
105:   PetscPrintf(PETSC_COMM_WORLD," Refinement: %s, tol=%g, its=%" PetscInt_FMT ", scheme=%s\n",PEPRefineTypes[refine],(double)tol,its,PEPRefineSchemes[rscheme]);

107:   PEPSetTarget(pep,4.8);
108:   PEPGetTarget(pep,&target);
109:   PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE);
110:   PEPGetWhichEigenpairs(pep,&which);
111:   PetscPrintf(PETSC_COMM_WORLD," Which = %d, target = %g\n",(int)which,(double)PetscRealPart(target));

113:   PEPSetDimensions(pep,4,PETSC_DEFAULT,PETSC_DEFAULT);
114:   PEPGetDimensions(pep,&nev,&ncv,&mpd);
115:   PetscPrintf(PETSC_COMM_WORLD," Dimensions: nev=%" PetscInt_FMT ", ncv=%" PetscInt_FMT ", mpd=%" PetscInt_FMT "\n",nev,ncv,mpd);

117:   PEPSetTolerances(pep,2.2e-4,200);
118:   PEPGetTolerances(pep,&tol,&its);
119:   PetscPrintf(PETSC_COMM_WORLD," Tolerance = %.5f, max_its = %" PetscInt_FMT "\n",(double)tol,its);

121:   PEPSetConvergenceTest(pep,PEP_CONV_ABS);
122:   PEPGetConvergenceTest(pep,&conv);
123:   PEPSetStoppingTest(pep,PEP_STOP_BASIC);
124:   PEPGetStoppingTest(pep,&stop);
125:   PetscPrintf(PETSC_COMM_WORLD," Convergence test = %d, stopping test = %d\n",(int)conv,(int)stop);

127:   PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);
128:   PEPMonitorSet(pep,(PetscErrorCode (*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))PEPMonitorFirst,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
129:   PEPMonitorCancel(pep);

131:   PEPGetST(pep,&st);
132:   STGetKSP(st,&ksp);
133:   KSPSetTolerances(ksp,1e-8,1e-35,PETSC_DEFAULT,PETSC_DEFAULT);
134:   STView(st,NULL);
135:   PEPGetDS(pep,&ds);
136:   DSView(ds,NULL);

138:   PEPSetFromOptions(pep);
139:   PEPSolve(pep);
140:   PEPGetConvergedReason(pep,&reason);
141:   PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d\n",(int)reason);

143:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144:                     Display solution and clean up
145:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146:   PEPErrorView(pep,PEP_ERROR_RELATIVE,NULL);
147:   PEPDestroy(&pep);
148:   MatDestroy(&A[0]);
149:   MatDestroy(&A[1]);
150:   MatDestroy(&A[2]);
151:   SlepcFinalize();
152:   return 0;
153: }

155: /*TEST

157:    test:
158:       suffix: 1
159:       args: -pep_tol 1e-6 -pep_ncv 22
160:       filter: sed -e "s/[+-]0\.0*i//g"

162: TEST*/