Actual source code: test1.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7:
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 QEP without calling QEPSetFromOptions (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 <qep_type> = qep type to test.\n"
27: " -epstype <eps_type> = eps type to test (for linear).\n\n";
29: #include <slepcqep.h>
33: int main(int argc,char **argv)
34: {
35: Mat M,C,K; /* problem matrices */
36: QEP qep; /* quadratic eigenproblem solver context */
37: const QEPType type;
38: PetscInt N,n=10,m,Istart,Iend,II,nev,maxit,i,j;
39: PetscBool flag;
40: char qeptype[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(PETSC_NULL,"-n",&n,PETSC_NULL);
50: PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flag);
51: if(!flag) m=n;
52: N = n*m;
53: PetscOptionsGetString(PETSC_NULL,"-type",qeptype,30,PETSC_NULL);
54: PetscOptionsGetString(PETSC_NULL,"-epstype",epstype,30,&flag);
55: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%D (%Dx%D grid)",N,n,m);
56: PetscPrintf(PETSC_COMM_WORLD,"\nQEP type: %s",qeptype);
57: if (flag) {
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*K+k*X+M)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:
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: }
81: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
82: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
84: /* C is the zero matrix */
85: MatCreate(PETSC_COMM_WORLD,&C);
86: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
87: MatSetFromOptions(C);
88: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
89: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
90:
91: /* M is the identity matrix */
92: MatCreate(PETSC_COMM_WORLD,&M);
93: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
94: MatSetFromOptions(M);
95: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
96: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
97: MatShift(M,1.0);
98:
99: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: Create the eigensolver and set various options
101: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103: /*
104: Create eigensolver context
105: */
106: QEPCreate(PETSC_COMM_WORLD,&qep);
108: /*
109: Set matrices and problem type
110: */
111: QEPSetOperators(qep,M,C,K);
112: QEPSetProblemType(qep,QEP_GENERAL);
113: QEPSetDimensions(qep,4,20,PETSC_DEFAULT);
114: QEPSetTolerances(qep,PETSC_SMALL,PETSC_DEFAULT);
116: /*
117: Set solver type at runtime
118: */
119: QEPSetType(qep,qeptype);
120: if (flag) {
121: PetscTypeCompare((PetscObject)qep,QEPLINEAR,&flag);
122: if (flag) {
123: QEPLinearGetEPS(qep,&eps);
124: EPSSetType(eps,epstype);
125: EPSGetST(eps,&st);
126: STGetKSP(st,&ksp);
127: KSPGetPC(ksp,&pc);
128: PCSetType(pc,PCJACOBI);
129: PetscTypeCompare((PetscObject)eps,EPSJD,&flag);
130: if (flag) {
131: EPSJDSetInitialSize(eps,1);
132: EPSJDSetFix(eps,PetscSqrtReal(PETSC_SMALL));
133: }
134: PetscTypeCompare((PetscObject)eps,EPSGD,&flag);
135: if (flag) {
136: EPSGDSetInitialSize(eps,1);
137: }
138: }
139: }
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Solve the eigensystem
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145: QEPSolve(qep);
147: /*
148: Optional: Get some information from the solver and display it
149: */
150: QEPGetType(qep,&type);
151: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
152: QEPGetDimensions(qep,&nev,PETSC_NULL,PETSC_NULL);
153: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
154: QEPGetTolerances(qep,PETSC_NULL,&maxit);
155: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: maxit=%D\n",maxit);
157: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158: Display solution and clean up
159: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161: QEPPrintSolution(qep,PETSC_NULL);
163: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164: Free work space
165: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
167: QEPDestroy(&qep);
168: MatDestroy(&M);
169: MatDestroy(&C);
170: MatDestroy(&K);
171: SlepcFinalize();
172: return 0;
173: }