Actual source code: pdde_stability.c
slepc-3.6.3 2016-03-29
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: */
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: /*
46: Function for user-defined eigenvalue ordering criterion.
48: Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
49: one of them as the preferred one according to the criterion.
50: In this example, the preferred value is the one with absolute value closest to 1.
51: */
52: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
53: {
54: PetscReal aa,ab;
57: aa = PetscAbsReal(SlepcAbsEigenvalue(ar,ai)-1.0);
58: ab = PetscAbsReal(SlepcAbsEigenvalue(br,bi)-1.0);
59: *r = aa > ab ? 1 : (aa < ab ? -1 : 0);
60: return(0);
61: }
65: int main(int argc,char **argv)
66: {
67: Mat A[NMAT]; /* problem matrices */
68: PEP pep; /* polynomial eigenproblem solver context */
69: PetscInt m=15,n,II,Istart,Iend,i,j,k;
70: PetscReal h,xi,xj,c[7] = { 2, .3, -2, .2, -2, -.3, -PETSC_PI/2 };
71: PetscScalar alpha,beta,gamma;
72: PetscBool flg,terse;
75: SlepcInitialize(&argc,&argv,(char*)0,help);
76: #if !defined(PETSC_USE_COMPLEX)
77: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires complex scalars");
78: #endif
80: PetscOptionsGetInt(NULL,"-m",&m,NULL);
81: n = m*m;
82: h = PETSC_PI/(m+1);
83: gamma = PetscExpScalar(PETSC_i*c[6]);
84: gamma = gamma/PetscAbsScalar(gamma);
85: k = 7;
86: PetscOptionsGetRealArray(NULL,"-c",c,&k,&flg);
87: if (flg && k!=7) SETERRQ1(PETSC_COMM_WORLD,1,"The number of parameters -c should be 7, you provided %D",k);
88: PetscPrintf(PETSC_COMM_WORLD,"\nPDDE stability, n=%D (m=%D)\n\n",n,m);
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Compute the polynomial matrices
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: /* initialize matrices */
95: for (i=0;i<NMAT;i++) {
96: MatCreate(PETSC_COMM_WORLD,&A[i]);
97: MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n);
98: MatSetFromOptions(A[i]);
99: MatSetUp(A[i]);
100: }
101: MatGetOwnershipRange(A[0],&Istart,&Iend);
103: /* A[1] has a pattern similar to the 2D Laplacian */
104: for (II=Istart;II<Iend;II++) {
105: i = II/m; j = II-i*m;
106: xi = (i+1)*h; xj = (j+1)*h;
107: alpha = c[0]+c[1]*PetscSinReal(xi)+gamma*(c[2]+c[3]*xi*(1.0-PetscExpReal(xi-PETSC_PI)));
108: beta = c[0]+c[1]*PetscSinReal(xj)-gamma*(c[2]+c[3]*xj*(1.0-PetscExpReal(xj-PETSC_PI)));
109: MatSetValue(A[1],II,II,alpha+beta-4.0/(h*h),INSERT_VALUES);
110: if (j>0) { MatSetValue(A[1],II,II-1,1.0/(h*h),INSERT_VALUES); }
111: if (j<m-1) { MatSetValue(A[1],II,II+1,1.0/(h*h),INSERT_VALUES); }
112: if (i>0) { MatSetValue(A[1],II,II-m,1.0/(h*h),INSERT_VALUES); }
113: if (i<m-1) { MatSetValue(A[1],II,II+m,1.0/(h*h),INSERT_VALUES); }
114: }
116: /* A[0] and A[2] are diagonal */
117: for (II=Istart;II<Iend;II++) {
118: i = II/m; j = II-i*m;
119: xi = (i+1)*h; xj = (j+1)*h;
120: alpha = c[4]+c[5]*xi*(PETSC_PI-xi);
121: beta = c[4]+c[5]*xj*(PETSC_PI-xj);
122: MatSetValue(A[0],II,II,alpha,INSERT_VALUES);
123: MatSetValue(A[2],II,II,beta,INSERT_VALUES);
124: }
125:
126: /* assemble matrices */
127: for (i=0;i<NMAT;i++) {
128: MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY);
129: }
130: for (i=0;i<NMAT;i++) {
131: MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY);
132: }
134: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135: Create the eigensolver and solve the problem
136: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138: PEPCreate(PETSC_COMM_WORLD,&pep);
139: PEPSetOperators(pep,NMAT,A);
140: PEPSetEigenvalueComparison(pep,MyEigenSort,NULL);
141: PEPSetDimensions(pep,4,PETSC_DEFAULT,PETSC_DEFAULT);
142: PEPSetFromOptions(pep);
143: PEPSolve(pep);
145: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146: Display solution and clean up
147: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148:
149: /* show detailed info unless -terse option is given by user */
150: PetscOptionsHasName(NULL,"-terse",&terse);
151: if (terse) {
152: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
153: } else {
154: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
155: PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
156: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
157: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
158: }
159: PEPDestroy(&pep);
160: for (i=0;i<NMAT;i++) {
161: MatDestroy(&A[i]);
162: }
163: SlepcFinalize();
164: return 0;
165: }