Actual source code: butterfly.c

slepc-3.18.2 2023-01-26
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    This example implements one of the problems found at
 12:        NLEVP: A Collection of Nonlinear Eigenvalue Problems,
 13:        The University of Manchester.
 14:    The details of the collection can be found at:
 15:        [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
 16:            Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.

 18:    The butterfly problem is a quartic PEP with T-even structure.
 19: */

 21: static char help[] = "Quartic polynomial eigenproblem with T-even structure.\n\n"
 22:   "The command line options are:\n"
 23:   "  -m <m>, grid size, the dimension of the matrices is n=m*m.\n"
 24:   "  -c <array>, problem parameters, must be 10 comma-separated real values.\n\n";

 26: #include <slepcpep.h>

 28: #define NMAT 5

 30: int main(int argc,char **argv)
 31: {
 32:   Mat            A[NMAT];         /* problem matrices */
 33:   PEP            pep;             /* polynomial eigenproblem solver context */
 34:   PetscInt       n,m=8,k,II,Istart,Iend,i,j;
 35:   PetscReal      c[10] = { 0.6, 1.3, 1.3, 0.1, 0.1, 1.2, 1.0, 1.0, 1.2, 1.0 };
 36:   PetscBool      flg;
 37:   PetscBool      terse;

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

 42:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 43:   n = m*m;
 44:   k = 10;
 45:   PetscOptionsGetRealArray(NULL,NULL,"-c",c,&k,&flg);
 47:   PetscPrintf(PETSC_COMM_WORLD,"\nButterfly problem, n=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n\n",n,m);

 49:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 50:                      Compute the polynomial matrices
 51:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 53:   /* initialize matrices */
 54:   for (i=0;i<NMAT;i++) {
 55:     MatCreate(PETSC_COMM_WORLD,&A[i]);
 56:     MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n);
 57:     MatSetFromOptions(A[i]);
 58:     MatSetUp(A[i]);
 59:   }
 60:   MatGetOwnershipRange(A[0],&Istart,&Iend);

 62:   /* A0 */
 63:   for (II=Istart;II<Iend;II++) {
 64:     i = II/m; j = II-i*m;
 65:     MatSetValue(A[0],II,II,4.0*c[0]/6.0+4.0*c[1]/6.0,INSERT_VALUES);
 66:     if (j>0) MatSetValue(A[0],II,II-1,c[0]/6.0,INSERT_VALUES);
 67:     if (j<m-1) MatSetValue(A[0],II,II+1,c[0]/6.0,INSERT_VALUES);
 68:     if (i>0) MatSetValue(A[0],II,II-m,c[1]/6.0,INSERT_VALUES);
 69:     if (i<m-1) MatSetValue(A[0],II,II+m,c[1]/6.0,INSERT_VALUES);
 70:   }

 72:   /* A1 */
 73:   for (II=Istart;II<Iend;II++) {
 74:     i = II/m; j = II-i*m;
 75:     if (j>0) MatSetValue(A[1],II,II-1,c[2],INSERT_VALUES);
 76:     if (j<m-1) MatSetValue(A[1],II,II+1,-c[2],INSERT_VALUES);
 77:     if (i>0) MatSetValue(A[1],II,II-m,c[3],INSERT_VALUES);
 78:     if (i<m-1) MatSetValue(A[1],II,II+m,-c[3],INSERT_VALUES);
 79:   }

 81:   /* A2 */
 82:   for (II=Istart;II<Iend;II++) {
 83:     i = II/m; j = II-i*m;
 84:     MatSetValue(A[2],II,II,-2.0*c[4]-2.0*c[5],INSERT_VALUES);
 85:     if (j>0) MatSetValue(A[2],II,II-1,c[4],INSERT_VALUES);
 86:     if (j<m-1) MatSetValue(A[2],II,II+1,c[4],INSERT_VALUES);
 87:     if (i>0) MatSetValue(A[2],II,II-m,c[5],INSERT_VALUES);
 88:     if (i<m-1) MatSetValue(A[2],II,II+m,c[5],INSERT_VALUES);
 89:   }

 91:   /* A3 */
 92:   for (II=Istart;II<Iend;II++) {
 93:     i = II/m; j = II-i*m;
 94:     if (j>0) MatSetValue(A[3],II,II-1,c[6],INSERT_VALUES);
 95:     if (j<m-1) MatSetValue(A[3],II,II+1,-c[6],INSERT_VALUES);
 96:     if (i>0) MatSetValue(A[3],II,II-m,c[7],INSERT_VALUES);
 97:     if (i<m-1) MatSetValue(A[3],II,II+m,-c[7],INSERT_VALUES);
 98:   }

100:   /* A4 */
101:   for (II=Istart;II<Iend;II++) {
102:     i = II/m; j = II-i*m;
103:     MatSetValue(A[4],II,II,2.0*c[8]+2.0*c[9],INSERT_VALUES);
104:     if (j>0) MatSetValue(A[4],II,II-1,-c[8],INSERT_VALUES);
105:     if (j<m-1) MatSetValue(A[4],II,II+1,-c[8],INSERT_VALUES);
106:     if (i>0) MatSetValue(A[4],II,II-m,-c[9],INSERT_VALUES);
107:     if (i<m-1) MatSetValue(A[4],II,II+m,-c[9],INSERT_VALUES);
108:   }

110:   /* assemble matrices */
111:   for (i=0;i<NMAT;i++) MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY);
112:   for (i=0;i<NMAT;i++) MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY);

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:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

127:   /* show detailed info unless -terse option is given by user */
128:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
129:   if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
130:   else {
131:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
132:     PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
133:     PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
134:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
135:   }
136:   PEPDestroy(&pep);
137:   for (i=0;i<NMAT;i++) MatDestroy(&A[i]);
138:   SlepcFinalize();
139:   return 0;
140: }

142: /*TEST

144:    testset:
145:       args: -pep_nev 4 -st_type sinvert -pep_target 0.01 -terse
146:       output_file: output/butterfly_1.out
147:       test:
148:          suffix: 1_toar
149:          args: -pep_type toar -pep_toar_restart 0.3
150:       test:
151:          suffix: 1_linear
152:          args: -pep_type linear

154:    test:
155:       suffix: 2
156:       args: -pep_type {{toar linear}} -pep_nev 4 -terse
157:       requires: double

159:    testset:
160:       args: -pep_type ciss -rg_type ellipse -rg_ellipse_center 1+.5i -rg_ellipse_radius .15 -terse
161:       requires: complex
162:       filter: sed -e "s/95386/95385/" | sed -e "s/91010/91009/" | sed -e "s/93092/93091/" | sed -e "s/96723/96724/" | sed -e "s/43015/43016/" | sed -e "s/53513/53514/"
163:       output_file: output/butterfly_ciss.out
164:       timeoutfactor: 2
165:       test:
166:          suffix: ciss_hankel
167:          args: -pep_ciss_extraction hankel -pep_ciss_integration_points 40
168:          requires: !single
169:       test:
170:          suffix: ciss_ritz
171:          args: -pep_ciss_extraction ritz
172:       test:
173:          suffix: ciss_caa
174:          args: -pep_ciss_extraction caa -pep_ciss_moments 4
175:       test:
176:          suffix: ciss_part
177:          nsize: 2
178:          args: -pep_ciss_partitions 2
179:       test:
180:          suffix: ciss_refine
181:          args: -pep_ciss_refine_inner 1 -pep_ciss_refine_blocksize 1

183:    testset:
184:       args: -pep_type ciss -rg_type ellipse -rg_ellipse_center .5+.5i -rg_ellipse_radius .25 -pep_ciss_moments 4 -pep_ciss_blocksize 5 -pep_ciss_refine_blocksize 2 -terse
185:       requires: complex double
186:       filter: sed -e "s/46483/46484/" | sed -e "s/54946/54945/" | sed -e "s/48456/48457/" | sed -e "s/74117/74116/" | sed -e "s/37240/37241/"
187:       output_file: output/butterfly_4.out
188:       test:
189:          suffix: 4
190:       test:
191:          suffix: 4_hankel
192:          args: -pep_ciss_extraction hankel

194: TEST*/