Actual source code: butterfly.c

slepc-3.22.2 2024-12-02
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;

 39:   PetscFunctionBeginUser;
 40:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));

 42:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
 43:   n = m*m;
 44:   k = 10;
 45:   PetscCall(PetscOptionsGetRealArray(NULL,NULL,"-c",c,&k,&flg));
 46:   PetscCheck(!flg || k==10,PETSC_COMM_WORLD,PETSC_ERR_USER,"The number of parameters -c should be 10, you provided %" PetscInt_FMT,k);
 47:   PetscCall(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:     PetscCall(MatCreate(PETSC_COMM_WORLD,&A[i]));
 56:     PetscCall(MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n));
 57:     PetscCall(MatSetFromOptions(A[i]));
 58:   }
 59:   PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));

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

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

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

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

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

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

113:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114:                 Create the eigensolver and solve the problem
115:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

117:   PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
118:   PetscCall(PEPSetOperators(pep,NMAT,A));
119:   PetscCall(PEPSetFromOptions(pep));
120:   PetscCall(PEPSolve(pep));

122:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123:                     Display solution and clean up
124:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

141: /*TEST

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

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

158:    testset:
159:       args: -pep_type ciss -rg_type ellipse -rg_ellipse_center 1+.5i -rg_ellipse_radius .15 -terse
160:       requires: complex
161:       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/"
162:       output_file: output/butterfly_ciss.out
163:       timeoutfactor: 2
164:       test:
165:          suffix: ciss_hankel
166:          args: -pep_ciss_extraction hankel -pep_ciss_integration_points 60
167:          requires: !single
168:       test:
169:          suffix: ciss_ritz
170:          args: -pep_ciss_extraction ritz
171:       test:
172:          suffix: ciss_caa
173:          args: -pep_ciss_extraction caa -pep_ciss_moments 4
174:       test:
175:          suffix: ciss_part
176:          nsize: 2
177:          args: -pep_ciss_partitions 2
178:       test:
179:          suffix: ciss_refine
180:          args: -pep_ciss_refine_inner 1 -pep_ciss_refine_blocksize 1

182:    testset:
183:       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
184:       requires: complex double
185:       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/"
186:       output_file: output/butterfly_4.out
187:       test:
188:          suffix: 4
189:       test:
190:          suffix: 4_hankel
191:          args: -pep_ciss_extraction hankel

193: TEST*/