Actual source code: test1.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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;
34: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
36: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
37: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
38: if (!flag) m=n;
39: N = n*m;
40: PetscOptionsGetString(NULL,NULL,"-type",peptype,30,NULL);
41: PetscOptionsGetString(NULL,NULL,"-epstype",epstype,30,&epsgiven);
42: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%D (%Dx%D 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: MatCreate(PETSC_COMM_WORLD,&K);
50: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N);
51: MatSetFromOptions(K);
52: MatSetUp(K);
53: MatGetOwnershipRange(K,&Istart,&Iend);
54: for (II=Istart;II<Iend;II++) {
55: i = II/n; j = II-i*n;
56: if (i>0) { MatSetValue(K,II,II-n,-1.0,INSERT_VALUES); }
57: if (i<m-1) { MatSetValue(K,II,II+n,-1.0,INSERT_VALUES); }
58: if (j>0) { MatSetValue(K,II,II-1,-1.0,INSERT_VALUES); }
59: if (j<n-1) { MatSetValue(K,II,II+1,-1.0,INSERT_VALUES); }
60: MatSetValue(K,II,II,4.0,INSERT_VALUES);
61: }
62: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
63: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
65: /* C is the 1-D Laplacian on horizontal lines */
66: MatCreate(PETSC_COMM_WORLD,&C);
67: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
68: MatSetFromOptions(C);
69: MatSetUp(C);
70: MatGetOwnershipRange(C,&Istart,&Iend);
71: for (II=Istart;II<Iend;II++) {
72: i = II/n; j = II-i*n;
73: if (j>0) { MatSetValue(C,II,II-1,-1.0,INSERT_VALUES); }
74: if (j<n-1) { MatSetValue(C,II,II+1,-1.0,INSERT_VALUES); }
75: MatSetValue(C,II,II,2.0,INSERT_VALUES);
76: }
77: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
78: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
80: /* M is a diagonal matrix */
81: MatCreate(PETSC_COMM_WORLD,&M);
82: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
83: MatSetFromOptions(M);
84: MatSetUp(M);
85: MatGetOwnershipRange(M,&Istart,&Iend);
86: for (II=Istart;II<Iend;II++) {
87: MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES);
88: }
89: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
90: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Create the eigensolver and set various options
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: PEPCreate(PETSC_COMM_WORLD,&pep);
97: A[0] = K; A[1] = C; A[2] = M;
98: PEPSetOperators(pep,3,A);
99: PEPSetProblemType(pep,PEP_GENERAL);
100: PEPSetDimensions(pep,4,20,PETSC_DEFAULT);
101: PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT);
103: /*
104: Set solver type at runtime
105: */
106: PEPSetType(pep,peptype);
107: if (epsgiven) {
108: PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&flag);
109: if (flag) {
110: PEPLinearGetEPS(pep,&eps);
111: PetscStrcmp(epstype,"gd2",&isgd2);
112: if (isgd2) {
113: EPSSetType(eps,EPSGD);
114: EPSGDSetDoubleExpansion(eps,PETSC_TRUE);
115: } else {
116: EPSSetType(eps,epstype);
117: }
118: EPSGetST(eps,&st);
119: STGetKSP(st,&ksp);
120: KSPGetPC(ksp,&pc);
121: PCSetType(pc,PCJACOBI);
122: PetscObjectTypeCompare((PetscObject)eps,EPSGD,&flag);
123: }
124: PEPLinearSetExplicitMatrix(pep,PETSC_TRUE);
125: }
126: PetscObjectTypeCompare((PetscObject)pep,PEPQARNOLDI,&flag);
127: if (flag) {
128: STCreate(PETSC_COMM_WORLD,&st);
129: STSetTransform(st,PETSC_TRUE);
130: PEPSetST(pep,st);
131: STDestroy(&st);
132: PEPQArnoldiGetRestart(pep,&keep);
133: PEPQArnoldiGetLocking(pep,&lock);
134: if (!lock && keep<0.6) {
135: PEPQArnoldiSetRestart(pep,0.6);
136: }
137: }
138: PetscObjectTypeCompare((PetscObject)pep,PEPTOAR,&flag);
139: if (flag) {
140: PEPTOARGetRestart(pep,&keep);
141: PEPTOARGetLocking(pep,&lock);
142: if (!lock && keep<0.6) {
143: PEPTOARSetRestart(pep,0.6);
144: }
145: }
147: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148: Solve the eigensystem
149: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151: PEPSolve(pep);
152: PEPGetDimensions(pep,&nev,NULL,NULL);
153: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
155: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156: Display solution and clean up
157: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
160: PEPDestroy(&pep);
161: MatDestroy(&M);
162: MatDestroy(&C);
163: MatDestroy(&K);
164: SlepcFinalize();
165: return ierr;
166: }
168: /*TEST
170: testset:
171: args: -m 11
172: requires: !single !complex
173: output_file: output/test1_1.out
174: test:
175: suffix: 1
176: args: -type {{toar qarnoldi linear}}
177: test:
178: suffix: 1_linear_gd
179: args: -type linear -epstype gd
181: TEST*/