Actual source code: test3.c
slepc-3.17.0 2022-03-31
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;
38: SlepcInitialize(&argc,&argv,(char*)0,help);
39: PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Quadratic Eigenproblem, n=%" PetscInt_FMT "\n\n",n);
41: MatCreate(PETSC_COMM_WORLD,&A[0]);
42: MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n);
43: MatSetFromOptions(A[0]);
44: MatSetUp(A[0]);
45: MatGetOwnershipRange(A[0],&Istart,&Iend);
46: for (i=Istart;i<Iend;i++) MatSetValue(A[0],i,i,i+1,INSERT_VALUES);
47: MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY);
48: MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY);
50: MatCreate(PETSC_COMM_WORLD,&A[1]);
51: MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,n,n);
52: MatSetFromOptions(A[1]);
53: MatSetUp(A[1]);
54: MatGetOwnershipRange(A[1],&Istart,&Iend);
55: for (i=Istart;i<Iend;i++) MatSetValue(A[1],i,i,1.0,INSERT_VALUES);
56: MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY);
57: MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY);
59: MatCreate(PETSC_COMM_WORLD,&A[2]);
60: MatSetSizes(A[2],PETSC_DECIDE,PETSC_DECIDE,n,n);
61: MatSetFromOptions(A[2]);
62: MatSetUp(A[2]);
63: MatGetOwnershipRange(A[1],&Istart,&Iend);
64: for (i=Istart;i<Iend;i++) MatSetValue(A[2],i,i,n/(PetscReal)(i+1),INSERT_VALUES);
65: MatAssemblyBegin(A[2],MAT_FINAL_ASSEMBLY);
66: MatAssemblyEnd(A[2],MAT_FINAL_ASSEMBLY);
68: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: Create eigensolver and test interface functions
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: PEPCreate(PETSC_COMM_WORLD,&pep);
72: PEPSetOperators(pep,3,A);
73: PEPGetNumMatrices(pep,&nmat);
74: PetscPrintf(PETSC_COMM_WORLD," Polynomial of degree %" PetscInt_FMT "\n",nmat-1);
75: PEPGetOperators(pep,0,&B);
76: MatView(B,NULL);
78: PEPSetType(pep,PEPTOAR);
79: PEPGetType(pep,&type);
80: PetscPrintf(PETSC_COMM_WORLD," Type set to %s\n",type);
82: PEPGetProblemType(pep,&ptype);
83: PetscPrintf(PETSC_COMM_WORLD," Problem type before changing = %d",(int)ptype);
84: PEPSetProblemType(pep,PEP_HERMITIAN);
85: PEPGetProblemType(pep,&ptype);
86: PetscPrintf(PETSC_COMM_WORLD," ... changed to %d.\n",(int)ptype);
88: PEPGetExtract(pep,&extr);
89: PetscPrintf(PETSC_COMM_WORLD," Extraction before changing = %d",(int)extr);
90: PEPSetExtract(pep,PEP_EXTRACT_STRUCTURED);
91: PEPGetExtract(pep,&extr);
92: PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)extr);
94: PEPSetScale(pep,PEP_SCALE_SCALAR,.1,NULL,NULL,5,1.0);
95: PEPGetScale(pep,&scale,&alpha,NULL,NULL,&its,NULL);
96: PetscPrintf(PETSC_COMM_WORLD," Scaling: %s, alpha=%g, its=%" PetscInt_FMT "\n",PEPScaleTypes[scale],(double)alpha,its);
98: PEPSetBasis(pep,PEP_BASIS_CHEBYSHEV1);
99: PEPGetBasis(pep,&basis);
100: PetscPrintf(PETSC_COMM_WORLD," Polynomial basis: %s\n",PEPBasisTypes[basis]);
102: PEPSetRefine(pep,PEP_REFINE_SIMPLE,1,1e-9,2,PEP_REFINE_SCHEME_SCHUR);
103: PEPGetRefine(pep,&refine,NULL,&tol,&its,&rscheme);
104: PetscPrintf(PETSC_COMM_WORLD," Refinement: %s, tol=%g, its=%" PetscInt_FMT ", scheme=%s\n",PEPRefineTypes[refine],(double)tol,its,PEPRefineSchemes[rscheme]);
106: PEPSetTarget(pep,4.8);
107: PEPGetTarget(pep,&target);
108: PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE);
109: PEPGetWhichEigenpairs(pep,&which);
110: PetscPrintf(PETSC_COMM_WORLD," Which = %d, target = %g\n",(int)which,(double)PetscRealPart(target));
112: PEPSetDimensions(pep,4,PETSC_DEFAULT,PETSC_DEFAULT);
113: PEPGetDimensions(pep,&nev,&ncv,&mpd);
114: PetscPrintf(PETSC_COMM_WORLD," Dimensions: nev=%" PetscInt_FMT ", ncv=%" PetscInt_FMT ", mpd=%" PetscInt_FMT "\n",nev,ncv,mpd);
116: PEPSetTolerances(pep,2.2e-4,200);
117: PEPGetTolerances(pep,&tol,&its);
118: PetscPrintf(PETSC_COMM_WORLD," Tolerance = %.5f, max_its = %" PetscInt_FMT "\n",(double)tol,its);
120: PEPSetConvergenceTest(pep,PEP_CONV_ABS);
121: PEPGetConvergenceTest(pep,&conv);
122: PEPSetStoppingTest(pep,PEP_STOP_BASIC);
123: PEPGetStoppingTest(pep,&stop);
124: PetscPrintf(PETSC_COMM_WORLD," Convergence test = %d, stopping test = %d\n",(int)conv,(int)stop);
126: PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);
127: PEPMonitorSet(pep,(PetscErrorCode (*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))PEPMonitorFirst,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
128: PEPMonitorCancel(pep);
130: PEPGetST(pep,&st);
131: STGetKSP(st,&ksp);
132: KSPSetTolerances(ksp,1e-8,1e-35,PETSC_DEFAULT,PETSC_DEFAULT);
133: STView(st,NULL);
134: PEPGetDS(pep,&ds);
135: DSView(ds,NULL);
137: PEPSetFromOptions(pep);
138: PEPSolve(pep);
139: PEPGetConvergedReason(pep,&reason);
140: PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d\n",(int)reason);
142: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143: Display solution and clean up
144: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145: PEPErrorView(pep,PEP_ERROR_RELATIVE,NULL);
146: PEPDestroy(&pep);
147: MatDestroy(&A[0]);
148: MatDestroy(&A[1]);
149: MatDestroy(&A[2]);
150: SlepcFinalize();
151: return 0;
152: }
154: /*TEST
156: test:
157: suffix: 1
158: args: -pep_tol 1e-6 -pep_ncv 22
159: filter: sed -e "s/[+-]0\.0*i//g"
161: TEST*/