Actual source code: butterfly.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  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 butterfly problem is a quartic PEP with T-even structure.
 30: */

 32: static char help[] = "NLEVP problem: butterfly.\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;

 52:   SlepcInitialize(&argc,&argv,(char*)0,help);

 54:   PetscOptionsGetInt(NULL,"-m",&m,NULL);
 55:   n = m*m;
 56:   k = 10;
 57:   PetscOptionsGetRealArray(NULL,"-c",c,&k,&flg);
 58:   if (flg && k!=10) SETERRQ1(PETSC_COMM_WORLD,1,"The number of parameters -c should be 10, you provided %D",k); 
 59:   PetscPrintf(PETSC_COMM_WORLD,"\nButterfly problem, n=%D (m=%D)\n\n",n,m);

 61:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 62:                      Compute the polynomial matrices 
 63:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 65:   /* initialize matrices */
 66:   for (i=0;i<NMAT;i++) {
 67:     MatCreate(PETSC_COMM_WORLD,&A[i]);
 68:     MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n);
 69:     MatSetFromOptions(A[i]);
 70:     MatSetUp(A[i]);
 71:   }
 72:   MatGetOwnershipRange(A[0],&Istart,&Iend);

 74:   /* A0 */
 75:   for (II=Istart;II<Iend;II++) {
 76:     i = II/m; j = II-i*m;
 77:     MatSetValue(A[0],II,II,4.0*c[0]/6.0+4.0*c[1]/6.0,INSERT_VALUES);
 78:     if (j>0) { MatSetValue(A[0],II,II-1,c[0]/6.0,INSERT_VALUES); }
 79:     if (j<m-1) { MatSetValue(A[0],II,II+1,c[0]/6.0,INSERT_VALUES); }
 80:     if (i>0) { MatSetValue(A[0],II,II-m,c[1]/6.0,INSERT_VALUES); }
 81:     if (i<m-1) { MatSetValue(A[0],II,II+m,c[1]/6.0,INSERT_VALUES); }
 82:   }

 84:   /* A1 */
 85:   for (II=Istart;II<Iend;II++) {
 86:     i = II/m; j = II-i*m;
 87:     if (j>0) { MatSetValue(A[1],II,II-1,c[2],INSERT_VALUES); }
 88:     if (j<m-1) { MatSetValue(A[1],II,II+1,-c[2],INSERT_VALUES); }
 89:     if (i>0) { MatSetValue(A[1],II,II-m,c[3],INSERT_VALUES); }
 90:     if (i<m-1) { MatSetValue(A[1],II,II+m,-c[3],INSERT_VALUES); }
 91:   }

 93:   /* A2 */
 94:   for (II=Istart;II<Iend;II++) {
 95:     i = II/m; j = II-i*m;
 96:     MatSetValue(A[2],II,II,-2.0*c[4]-2.0*c[5],INSERT_VALUES);
 97:     if (j>0) { MatSetValue(A[2],II,II-1,c[4],INSERT_VALUES); }
 98:     if (j<m-1) { MatSetValue(A[2],II,II+1,c[4],INSERT_VALUES); }
 99:     if (i>0) { MatSetValue(A[2],II,II-m,c[5],INSERT_VALUES); }
100:     if (i<m-1) { MatSetValue(A[2],II,II+m,c[5],INSERT_VALUES); }
101:   }

103:   /* A3 */
104:   for (II=Istart;II<Iend;II++) {
105:     i = II/m; j = II-i*m;
106:     if (j>0) { MatSetValue(A[3],II,II-1,c[6],INSERT_VALUES); }
107:     if (j<m-1) { MatSetValue(A[3],II,II+1,-c[6],INSERT_VALUES); }
108:     if (i>0) { MatSetValue(A[3],II,II-m,c[7],INSERT_VALUES); }
109:     if (i<m-1) { MatSetValue(A[3],II,II+m,-c[7],INSERT_VALUES); }
110:   }

112:   /* A4 */
113:   for (II=Istart;II<Iend;II++) {
114:     i = II/m; j = II-i*m;
115:     MatSetValue(A[4],II,II,2.0*c[8]+2.0*c[9],INSERT_VALUES);
116:     if (j>0) { MatSetValue(A[4],II,II-1,-c[8],INSERT_VALUES); }
117:     if (j<m-1) { MatSetValue(A[4],II,II+1,-c[8],INSERT_VALUES); }
118:     if (i>0) { MatSetValue(A[4],II,II-m,-c[9],INSERT_VALUES); }
119:     if (i<m-1) { MatSetValue(A[4],II,II+m,-c[9],INSERT_VALUES); }
120:   }

122:   /* assemble matrices */
123:   for (i=0;i<NMAT;i++) {
124:     MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY);
125:   }
126:   for (i=0;i<NMAT;i++) {
127:     MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY);
128:   }

130:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
131:                 Create the eigensolver and solve the problem
132:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

134:   PEPCreate(PETSC_COMM_WORLD,&pep);
135:   PEPSetOperators(pep,NMAT,A);
136:   PEPSetFromOptions(pep);
137:   PEPSolve(pep);

139:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140:                     Display solution and clean up
141:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142:   
143:   PEPPrintSolution(pep,NULL);
144:   PEPDestroy(&pep);
145:   for (i=0;i<NMAT;i++) {
146:     MatDestroy(&A[i]);
147:   }
148:   SlepcFinalize();
149:   return 0;
150: }