Actual source code: test1.c
slepc-3.7.2 2016-07-19
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Test the solution of a PEP without calling PEPSetFromOptions (based on ex16.c).\n\n"
23: "The command line options are:\n"
24: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
25: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n"
26: " -type <pep_type> = pep type to test.\n"
27: " -epstype <eps_type> = eps type to test (for linear).\n\n";
29: #include <slepcpep.h>
33: int main(int argc,char **argv)
34: {
35: Mat M,C,K,A[3]; /* problem matrices */
36: PEP pep; /* polynomial eigenproblem solver context */
37: PEPType type;
38: PetscInt N,n=10,m,Istart,Iend,II,nev,maxit,i,j;
39: PetscBool flag,isgd2,epsgiven;
40: char peptype[30] = "linear",epstype[30] = "";
41: EPS eps;
42: ST st;
43: KSP ksp;
44: PC pc;
47: SlepcInitialize(&argc,&argv,(char*)0,help);
49: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
50: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
51: if (!flag) m=n;
52: N = n*m;
53: PetscOptionsGetString(NULL,NULL,"-type",peptype,30,NULL);
54: PetscOptionsGetString(NULL,NULL,"-epstype",epstype,30,&epsgiven);
55: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%D (%Dx%D grid)",N,n,m);
56: PetscPrintf(PETSC_COMM_WORLD,"\nPEP type: %s",peptype);
57: if (epsgiven) {
58: PetscPrintf(PETSC_COMM_WORLD,"\nEPS type: %s",epstype);
59: }
60: PetscPrintf(PETSC_COMM_WORLD,"\n\n");
62: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
64: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66: /* K is the 2-D Laplacian */
67: MatCreate(PETSC_COMM_WORLD,&K);
68: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N);
69: MatSetFromOptions(K);
70: MatSetUp(K);
71: MatGetOwnershipRange(K,&Istart,&Iend);
72: for (II=Istart;II<Iend;II++) {
73: i = II/n; j = II-i*n;
74: if (i>0) { MatSetValue(K,II,II-n,-1.0,INSERT_VALUES); }
75: if (i<m-1) { MatSetValue(K,II,II+n,-1.0,INSERT_VALUES); }
76: if (j>0) { MatSetValue(K,II,II-1,-1.0,INSERT_VALUES); }
77: if (j<n-1) { MatSetValue(K,II,II+1,-1.0,INSERT_VALUES); }
78: MatSetValue(K,II,II,4.0,INSERT_VALUES);
79: }
80: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
81: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
83: /* C is the 1-D Laplacian on horizontal lines */
84: MatCreate(PETSC_COMM_WORLD,&C);
85: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
86: MatSetFromOptions(C);
87: MatSetUp(C);
88: MatGetOwnershipRange(C,&Istart,&Iend);
89: for (II=Istart;II<Iend;II++) {
90: i = II/n; j = II-i*n;
91: if (j>0) { MatSetValue(C,II,II-1,-1.0,INSERT_VALUES); }
92: if (j<n-1) { MatSetValue(C,II,II+1,-1.0,INSERT_VALUES); }
93: MatSetValue(C,II,II,2.0,INSERT_VALUES);
94: }
95: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
96: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
98: /* M is a diagonal matrix */
99: MatCreate(PETSC_COMM_WORLD,&M);
100: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
101: MatSetFromOptions(M);
102: MatSetUp(M);
103: MatGetOwnershipRange(M,&Istart,&Iend);
104: for (II=Istart;II<Iend;II++) {
105: MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES);
106: }
107: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
108: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
110: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111: Create the eigensolver and set various options
112: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: PEPCreate(PETSC_COMM_WORLD,&pep);
115: A[0] = K; A[1] = C; A[2] = M;
116: PEPSetOperators(pep,3,A);
117: PEPSetProblemType(pep,PEP_GENERAL);
118: PEPSetDimensions(pep,4,20,PETSC_DEFAULT);
119: PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT);
121: /*
122: Set solver type at runtime
123: */
124: PEPSetType(pep,peptype);
125: if (epsgiven) {
126: PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&flag);
127: if (flag) {
128: PEPLinearGetEPS(pep,&eps);
129: PetscStrcmp(epstype,"gd2",&isgd2);
130: if (isgd2) {
131: EPSSetType(eps,EPSGD);
132: EPSGDSetDoubleExpansion(eps,PETSC_TRUE);
133: } else {
134: EPSSetType(eps,epstype);
135: }
136: EPSGetST(eps,&st);
137: STGetKSP(st,&ksp);
138: KSPGetPC(ksp,&pc);
139: PCSetType(pc,PCJACOBI);
140: PetscObjectTypeCompare((PetscObject)eps,EPSGD,&flag);
141: }
142: PEPLinearSetExplicitMatrix(pep,PETSC_TRUE);
143: }
144: PetscObjectTypeCompare((PetscObject)pep,PEPQARNOLDI,&flag);
145: if (flag) {
146: PEPGetST(pep,&st);
147: STSetTransform(st,PETSC_TRUE);
148: }
150: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151: Solve the eigensystem
152: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154: PEPSolve(pep);
156: /*
157: Optional: Get some information from the solver and display it
158: */
159: PEPGetType(pep,&type);
160: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
161: PEPGetDimensions(pep,&nev,NULL,NULL);
162: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
163: PEPGetTolerances(pep,NULL,&maxit);
164: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: maxit=%D\n",maxit);
166: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167: Display solution and clean up
168: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
171: PEPDestroy(&pep);
172: MatDestroy(&M);
173: MatDestroy(&C);
174: MatDestroy(&K);
175: SlepcFinalize();
176: return ierr;
177: }