Actual source code: test1.c
slepc-3.20.1 2023-11-27
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 the solution of a PEP without calling PEPSetFromOptions (based on ex16.c).\n\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
14: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n"
15: " -type <pep_type> = pep type to test.\n"
16: " -epstype <eps_type> = eps type to test (for linear).\n\n";
18: #include <slepcpep.h>
20: int main(int argc,char **argv)
21: {
22: Mat M,C,K,A[3]; /* problem matrices */
23: PEP pep; /* polynomial eigenproblem solver context */
24: PetscInt N,n=10,m,Istart,Iend,II,nev,i,j;
25: PetscReal keep;
26: PetscBool flag,isgd2,epsgiven,lock;
27: char peptype[30] = "linear",epstype[30] = "";
28: EPS eps;
29: ST st;
30: KSP ksp;
31: PC pc;
33: PetscFunctionBeginUser;
34: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
36: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
37: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
38: if (!flag) m=n;
39: N = n*m;
40: PetscCall(PetscOptionsGetString(NULL,NULL,"-type",peptype,sizeof(peptype),NULL));
41: PetscCall(PetscOptionsGetString(NULL,NULL,"-epstype",epstype,sizeof(epstype),&epsgiven));
42: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
44: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
46: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
48: /* K is the 2-D Laplacian */
49: PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
50: PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N));
51: PetscCall(MatSetFromOptions(K));
52: PetscCall(MatSetUp(K));
53: PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
54: for (II=Istart;II<Iend;II++) {
55: i = II/n; j = II-i*n;
56: if (i>0) PetscCall(MatSetValue(K,II,II-n,-1.0,INSERT_VALUES));
57: if (i<m-1) PetscCall(MatSetValue(K,II,II+n,-1.0,INSERT_VALUES));
58: if (j>0) PetscCall(MatSetValue(K,II,II-1,-1.0,INSERT_VALUES));
59: if (j<n-1) PetscCall(MatSetValue(K,II,II+1,-1.0,INSERT_VALUES));
60: PetscCall(MatSetValue(K,II,II,4.0,INSERT_VALUES));
61: }
62: PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
63: PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
65: /* C is the 1-D Laplacian on horizontal lines */
66: PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
67: PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N));
68: PetscCall(MatSetFromOptions(C));
69: PetscCall(MatSetUp(C));
70: PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
71: for (II=Istart;II<Iend;II++) {
72: i = II/n; j = II-i*n;
73: if (j>0) PetscCall(MatSetValue(C,II,II-1,-1.0,INSERT_VALUES));
74: if (j<n-1) PetscCall(MatSetValue(C,II,II+1,-1.0,INSERT_VALUES));
75: PetscCall(MatSetValue(C,II,II,2.0,INSERT_VALUES));
76: }
77: PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
78: PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
80: /* M is a diagonal matrix */
81: PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
82: PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N));
83: PetscCall(MatSetFromOptions(M));
84: PetscCall(MatSetUp(M));
85: PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
86: for (II=Istart;II<Iend;II++) PetscCall(MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES));
87: PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
88: PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Create the eigensolver and set various options
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
95: A[0] = K; A[1] = C; A[2] = M;
96: PetscCall(PEPSetOperators(pep,3,A));
97: PetscCall(PEPSetProblemType(pep,PEP_GENERAL));
98: PetscCall(PEPSetDimensions(pep,4,20,PETSC_DEFAULT));
99: PetscCall(PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT));
101: /*
102: Set solver type at runtime
103: */
104: PetscCall(PEPSetType(pep,peptype));
105: if (epsgiven) {
106: PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&flag));
107: if (flag) {
108: PetscCall(PEPLinearGetEPS(pep,&eps));
109: PetscCall(PetscStrcmp(epstype,"gd2",&isgd2));
110: if (isgd2) {
111: PetscCall(EPSSetType(eps,EPSGD));
112: PetscCall(EPSGDSetDoubleExpansion(eps,PETSC_TRUE));
113: } else PetscCall(EPSSetType(eps,epstype));
114: PetscCall(EPSGetST(eps,&st));
115: PetscCall(STGetKSP(st,&ksp));
116: PetscCall(KSPGetPC(ksp,&pc));
117: PetscCall(PCSetType(pc,PCJACOBI));
118: PetscCall(PetscObjectTypeCompare((PetscObject)eps,EPSGD,&flag));
119: }
120: PetscCall(PEPLinearSetExplicitMatrix(pep,PETSC_TRUE));
121: }
122: PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPQARNOLDI,&flag));
123: if (flag) {
124: PetscCall(STCreate(PETSC_COMM_WORLD,&st));
125: PetscCall(STSetTransform(st,PETSC_TRUE));
126: PetscCall(PEPSetST(pep,st));
127: PetscCall(STDestroy(&st));
128: PetscCall(PEPQArnoldiGetRestart(pep,&keep));
129: PetscCall(PEPQArnoldiGetLocking(pep,&lock));
130: if (!lock && keep<0.6) PetscCall(PEPQArnoldiSetRestart(pep,0.6));
131: }
132: PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPTOAR,&flag));
133: if (flag) {
134: PetscCall(PEPTOARGetRestart(pep,&keep));
135: PetscCall(PEPTOARGetLocking(pep,&lock));
136: if (!lock && keep<0.6) PetscCall(PEPTOARSetRestart(pep,0.6));
137: }
139: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140: Solve the eigensystem
141: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143: PetscCall(PEPSolve(pep));
144: PetscCall(PEPGetDimensions(pep,&nev,NULL,NULL));
145: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
147: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148: Display solution and clean up
149: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151: PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
152: PetscCall(PEPDestroy(&pep));
153: PetscCall(MatDestroy(&M));
154: PetscCall(MatDestroy(&C));
155: PetscCall(MatDestroy(&K));
156: PetscCall(SlepcFinalize());
157: return 0;
158: }
160: /*TEST
162: testset:
163: args: -m 11
164: output_file: output/test1_1.out
165: filter: sed -e "s/1.16403/1.16404/g" | sed -e "s/1.65362i/1.65363i/g" | sed -e "s/-1.16404-1.65363i, -1.16404+1.65363i/-1.16404+1.65363i, -1.16404-1.65363i/" | sed -e "s/-0.51784-1.31039i, -0.51784+1.31039i/-0.51784+1.31039i, -0.51784-1.31039i/"
166: requires: !single
167: test:
168: suffix: 1
169: args: -type {{toar qarnoldi linear}}
170: test:
171: suffix: 1_linear_gd
172: args: -type linear -epstype gd
173: requires: !__float128
175: TEST*/