Actual source code: test3.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 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:   Vec                Dl,Dr;
 22:   DS                 ds;
 23:   PetscReal          tol,alpha;
 24:   PetscScalar        target;
 25:   PetscInt           n=20,i,its,nev,ncv,mpd,Istart,Iend,nmat;
 26:   PEPWhich           which;
 27:   PEPConvergedReason reason;
 28:   PEPType            type;
 29:   PEPExtract         extr;
 30:   PEPBasis           basis;
 31:   PEPScale           scale;
 32:   PEPRefine          refine;
 33:   PEPRefineScheme    rscheme;
 34:   PEPConv            conv;
 35:   PEPStop            stop;
 36:   PEPProblemType     ptype;
 37:   PetscViewerAndFormat *vf;

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

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

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

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

 67:   PetscCall(MatCreateVecs(A[0],&Dr,&Dl));
 68:   PetscCall(VecSet(Dl,1.0));
 69:   PetscCall(VecSet(Dr,0.95));

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

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

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

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

 97:   PetscCall(PEPSetScale(pep,PEP_SCALE_BOTH,.1,Dl,Dr,5,1.0));
 98:   PetscCall(PEPGetScale(pep,&scale,&alpha,NULL,NULL,&its,NULL));
 99:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Scaling: %s, alpha=%g, its=%" PetscInt_FMT "\n",PEPScaleTypes[scale],(double)alpha,its));

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

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

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

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

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

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

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

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

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

145:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146:                     Display solution and clean up
147:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148:   PetscCall(PEPErrorView(pep,PEP_ERROR_RELATIVE,NULL));
149:   PetscCall(PEPDestroy(&pep));
150:   PetscCall(MatDestroy(&A[0]));
151:   PetscCall(MatDestroy(&A[1]));
152:   PetscCall(MatDestroy(&A[2]));
153:   PetscCall(VecDestroy(&Dl));
154:   PetscCall(VecDestroy(&Dr));
155:   PetscCall(SlepcFinalize());
156:   return 0;
157: }

159: /*TEST

161:    test:
162:       suffix: 1
163:       args: -pep_tol 1e-6 -pep_ncv 22
164:       filter: sed -e "s/[+-]0\.0*i//g"

166: TEST*/