Actual source code: pdde_stability.c
slepc-3.5.2 2014-10-10
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2014, 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: */
21: /*
22: This example implements one of the problems found at
23: NLEVP: A Collection of Nonlinear Eigenvalue Problems,
24: The University of Manchester.
25: The details of the collection can be found at:
26: [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
27: Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.
29: The pdde_stability problem is a complex-symmetric QEP from the stability
30: analysis of a discretized partial delay-differential equation. It requires
31: complex scalars.
32: */
34: static char help[] = "NLEVP problem: pdde_stability.\n\n"
35: "The command line options are:\n"
36: " -m <m>, grid size, the matrices have dimension n=m*m.\n"
37: " -c <a0,b0,a1,b1,a2,b2,phi1>, comma-separated list of 7 real parameters.\n\n";
39: #include <slepcpep.h>
41: #define NMAT 3
45: int main(int argc,char **argv)
46: {
47: Mat A[NMAT]; /* problem matrices */
48: PEP pep; /* polynomial eigenproblem solver context */
49: PetscInt m=15,n,II,Istart,Iend,i,j,k;
50: PetscReal h,xi,xj,c[7] = { 2, .3, -2, .2, -2, -.3, -PETSC_PI/2 };
51: PetscScalar alpha,beta,gamma;
52: PetscBool flg;
55: SlepcInitialize(&argc,&argv,(char*)0,help);
56: #if !defined(PETSC_USE_COMPLEX)
57: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires complex scalars");
58: #endif
60: PetscOptionsGetInt(NULL,"-m",&m,NULL);
61: n = m*m;
62: h = PETSC_PI/(m+1);
63: gamma = PetscExpScalar(PETSC_i*c[6]);
64: gamma = gamma/PetscAbsScalar(gamma);
65: k = 7;
66: PetscOptionsGetRealArray(NULL,"-c",c,&k,&flg);
67: if (flg && k!=7) SETERRQ1(PETSC_COMM_WORLD,1,"The number of parameters -c should be 7, you provided %D",k);
68: PetscPrintf(PETSC_COMM_WORLD,"\nPDDE stability, n=%D (m=%D)\n\n",n,m);
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Compute the polynomial matrices
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74: /* initialize matrices */
75: for (i=0;i<NMAT;i++) {
76: MatCreate(PETSC_COMM_WORLD,&A[i]);
77: MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n);
78: MatSetFromOptions(A[i]);
79: MatSetUp(A[i]);
80: }
81: MatGetOwnershipRange(A[0],&Istart,&Iend);
83: /* A[1] has a pattern similar to the 2D Laplacian */
84: for (II=Istart;II<Iend;II++) {
85: i = II/m; j = II-i*m;
86: xi = (i+1)*h; xj = (j+1)*h;
87: alpha = c[0]+c[1]*PetscSinReal(xi)+gamma*(c[2]+c[3]*xi*(1.0-PetscExpReal(xi-PETSC_PI)));
88: beta = c[0]+c[1]*PetscSinReal(xj)-gamma*(c[2]+c[3]*xj*(1.0-PetscExpReal(xj-PETSC_PI)));
89: MatSetValue(A[1],II,II,alpha+beta-4.0/(h*h),INSERT_VALUES);
90: if (j>0) { MatSetValue(A[1],II,II-1,1.0/(h*h),INSERT_VALUES); }
91: if (j<m-1) { MatSetValue(A[1],II,II+1,1.0/(h*h),INSERT_VALUES); }
92: if (i>0) { MatSetValue(A[1],II,II-m,1.0/(h*h),INSERT_VALUES); }
93: if (i<m-1) { MatSetValue(A[1],II,II+m,1.0/(h*h),INSERT_VALUES); }
94: }
96: /* A[0] and A[2] are diagonal */
97: for (II=Istart;II<Iend;II++) {
98: i = II/m; j = II-i*m;
99: xi = (i+1)*h; xj = (j+1)*h;
100: alpha = c[4]+c[5]*xi*(PETSC_PI-xi);
101: beta = c[4]+c[5]*xj*(PETSC_PI-xj);
102: MatSetValue(A[0],II,II,alpha,INSERT_VALUES);
103: MatSetValue(A[2],II,II,beta,INSERT_VALUES);
104: }
105:
106: /* assemble matrices */
107: for (i=0;i<NMAT;i++) {
108: MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY);
109: }
110: for (i=0;i<NMAT;i++) {
111: MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY);
112: }
114: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115: Create the eigensolver and solve the problem
116: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118: PEPCreate(PETSC_COMM_WORLD,&pep);
119: PEPSetOperators(pep,NMAT,A);
120: PEPSetFromOptions(pep);
121: PEPSolve(pep);
123: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124: Display solution and clean up
125: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126:
127: PEPPrintSolution(pep,NULL);
128: PEPDestroy(&pep);
129: for (i=0;i<NMAT;i++) {
130: MatDestroy(&A[i]);
131: }
132: SlepcFinalize();
133: return 0;
134: }