Actual source code: ex16.c
slepc-3.6.1 2015-09-03
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, 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[] = "Quadratic eigenproblem for testing the PEP object.\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\n";
27: #include <slepcpep.h>
31: int main(int argc,char **argv)
32: {
33: Mat M,C,K,A[3]; /* problem matrices */
34: PEP pep; /* polynomial eigenproblem solver context */
35: PEPType type;
36: PetscInt N,n=10,m,Istart,Iend,II,nev,i,j;
37: PetscBool flag,terse;
40: SlepcInitialize(&argc,&argv,(char*)0,help);
42: PetscOptionsGetInt(NULL,"-n",&n,NULL);
43: PetscOptionsGetInt(NULL,"-m",&m,&flag);
44: if (!flag) m=n;
45: N = n*m;
46: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%D (%Dx%D grid)\n\n",N,n,m);
48: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
50: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
52: /* K is the 2-D Laplacian */
53: MatCreate(PETSC_COMM_WORLD,&K);
54: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N);
55: MatSetFromOptions(K);
56: MatSetUp(K);
58: MatGetOwnershipRange(K,&Istart,&Iend);
59: for (II=Istart;II<Iend;II++) {
60: i = II/n; j = II-i*n;
61: if (i>0) { MatSetValue(K,II,II-n,-1.0,INSERT_VALUES); }
62: if (i<m-1) { MatSetValue(K,II,II+n,-1.0,INSERT_VALUES); }
63: if (j>0) { MatSetValue(K,II,II-1,-1.0,INSERT_VALUES); }
64: if (j<n-1) { MatSetValue(K,II,II+1,-1.0,INSERT_VALUES); }
65: MatSetValue(K,II,II,4.0,INSERT_VALUES);
66: }
68: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
69: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
71: /* C is the zero matrix */
72: MatCreate(PETSC_COMM_WORLD,&C);
73: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
74: MatSetFromOptions(C);
75: MatSetUp(C);
76: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
77: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
79: /* M is the identity matrix */
80: MatCreate(PETSC_COMM_WORLD,&M);
81: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
82: MatSetFromOptions(M);
83: MatSetUp(M);
84: MatGetOwnershipRange(M,&Istart,&Iend);
85: for (i=Istart;i<Iend;i++) {
86: MatSetValue(M,i,i,1.0,INSERT_VALUES);
87: }
88: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
89: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
91: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92: Create the eigensolver and set various options
93: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95: /*
96: Create eigensolver context
97: */
98: PEPCreate(PETSC_COMM_WORLD,&pep);
100: /*
101: Set matrices and problem type
102: */
103: A[0] = K; A[1] = C; A[2] = M;
104: PEPSetOperators(pep,3,A);
105: PEPSetProblemType(pep,PEP_HERMITIAN);
107: /*
108: Set solver parameters at runtime
109: */
110: PEPSetFromOptions(pep);
112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: Solve the eigensystem
114: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116: PEPSolve(pep);
118: /*
119: Optional: Get some information from the solver and display it
120: */
121: PEPGetType(pep,&type);
122: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
123: PEPGetDimensions(pep,&nev,NULL,NULL);
124: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Display solution and clean up
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: /* show detailed info unless -terse option is given by user */
131: PetscOptionsHasName(NULL,"-terse",&terse);
132: if (terse) {
133: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
134: } else {
135: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
136: PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
137: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
138: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
139: }
140: PEPDestroy(&pep);
141: MatDestroy(&M);
142: MatDestroy(&C);
143: MatDestroy(&K);
144: SlepcFinalize();
145: return 0;
146: }