Actual source code: butterfly.c
slepc-3.7.3 2016-09-29
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: */
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 butterfly problem is a quartic PEP with T-even structure.
30: */
32: static char help[] = "Quartic polynomial eigenproblem with T-even structure.\n\n"
33: "The command line options are:\n"
34: " -m <m>, grid size, the dimension of the matrices is n=m*m.\n"
35: " -c <array>, problem parameters, must be 10 comma-separated real values.\n\n";
37: #include <slepcpep.h>
39: #define NMAT 5
43: int main(int argc,char **argv)
44: {
45: Mat A[NMAT]; /* problem matrices */
46: PEP pep; /* polynomial eigenproblem solver context */
47: PetscInt n,m=8,k,II,Istart,Iend,i,j;
48: PetscReal c[10] = { 0.6, 1.3, 1.3, 0.1, 0.1, 1.2, 1.0, 1.0, 1.2, 1.0 };
49: PetscBool flg;
50: PetscBool terse;
53: SlepcInitialize(&argc,&argv,(char*)0,help);
55: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
56: n = m*m;
57: k = 10;
58: PetscOptionsGetRealArray(NULL,NULL,"-c",c,&k,&flg);
59: if (flg && k!=10) SETERRQ1(PETSC_COMM_WORLD,1,"The number of parameters -c should be 10, you provided %D",k);
60: PetscPrintf(PETSC_COMM_WORLD,"\nButterfly problem, n=%D (m=%D)\n\n",n,m);
62: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: Compute the polynomial matrices
64: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66: /* initialize matrices */
67: for (i=0;i<NMAT;i++) {
68: MatCreate(PETSC_COMM_WORLD,&A[i]);
69: MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n);
70: MatSetFromOptions(A[i]);
71: MatSetUp(A[i]);
72: }
73: MatGetOwnershipRange(A[0],&Istart,&Iend);
75: /* A0 */
76: for (II=Istart;II<Iend;II++) {
77: i = II/m; j = II-i*m;
78: MatSetValue(A[0],II,II,4.0*c[0]/6.0+4.0*c[1]/6.0,INSERT_VALUES);
79: if (j>0) { MatSetValue(A[0],II,II-1,c[0]/6.0,INSERT_VALUES); }
80: if (j<m-1) { MatSetValue(A[0],II,II+1,c[0]/6.0,INSERT_VALUES); }
81: if (i>0) { MatSetValue(A[0],II,II-m,c[1]/6.0,INSERT_VALUES); }
82: if (i<m-1) { MatSetValue(A[0],II,II+m,c[1]/6.0,INSERT_VALUES); }
83: }
85: /* A1 */
86: for (II=Istart;II<Iend;II++) {
87: i = II/m; j = II-i*m;
88: if (j>0) { MatSetValue(A[1],II,II-1,c[2],INSERT_VALUES); }
89: if (j<m-1) { MatSetValue(A[1],II,II+1,-c[2],INSERT_VALUES); }
90: if (i>0) { MatSetValue(A[1],II,II-m,c[3],INSERT_VALUES); }
91: if (i<m-1) { MatSetValue(A[1],II,II+m,-c[3],INSERT_VALUES); }
92: }
94: /* A2 */
95: for (II=Istart;II<Iend;II++) {
96: i = II/m; j = II-i*m;
97: MatSetValue(A[2],II,II,-2.0*c[4]-2.0*c[5],INSERT_VALUES);
98: if (j>0) { MatSetValue(A[2],II,II-1,c[4],INSERT_VALUES); }
99: if (j<m-1) { MatSetValue(A[2],II,II+1,c[4],INSERT_VALUES); }
100: if (i>0) { MatSetValue(A[2],II,II-m,c[5],INSERT_VALUES); }
101: if (i<m-1) { MatSetValue(A[2],II,II+m,c[5],INSERT_VALUES); }
102: }
104: /* A3 */
105: for (II=Istart;II<Iend;II++) {
106: i = II/m; j = II-i*m;
107: if (j>0) { MatSetValue(A[3],II,II-1,c[6],INSERT_VALUES); }
108: if (j<m-1) { MatSetValue(A[3],II,II+1,-c[6],INSERT_VALUES); }
109: if (i>0) { MatSetValue(A[3],II,II-m,c[7],INSERT_VALUES); }
110: if (i<m-1) { MatSetValue(A[3],II,II+m,-c[7],INSERT_VALUES); }
111: }
113: /* A4 */
114: for (II=Istart;II<Iend;II++) {
115: i = II/m; j = II-i*m;
116: MatSetValue(A[4],II,II,2.0*c[8]+2.0*c[9],INSERT_VALUES);
117: if (j>0) { MatSetValue(A[4],II,II-1,-c[8],INSERT_VALUES); }
118: if (j<m-1) { MatSetValue(A[4],II,II+1,-c[8],INSERT_VALUES); }
119: if (i>0) { MatSetValue(A[4],II,II-m,-c[9],INSERT_VALUES); }
120: if (i<m-1) { MatSetValue(A[4],II,II+m,-c[9],INSERT_VALUES); }
121: }
123: /* assemble matrices */
124: for (i=0;i<NMAT;i++) {
125: MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY);
126: }
127: for (i=0;i<NMAT;i++) {
128: MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY);
129: }
131: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132: Create the eigensolver and solve the problem
133: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135: PEPCreate(PETSC_COMM_WORLD,&pep);
136: PEPSetOperators(pep,NMAT,A);
137: PEPSetFromOptions(pep);
138: PEPSolve(pep);
140: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141: Display solution and clean up
142: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143:
144: /* show detailed info unless -terse option is given by user */
145: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
146: if (terse) {
147: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
148: } else {
149: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
150: PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
151: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
152: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
153: }
154: PEPDestroy(&pep);
155: for (i=0;i<NMAT;i++) {
156: MatDestroy(&A[i]);
157: }
158: SlepcFinalize();
159: return ierr;
160: }