Actual source code: ex49.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: */

 11: static char help[] = "User-defined split preconditioner when solving a generalized eigenproblem.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 14:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";

 16: #include <slepceps.h>

 18: int main(int argc,char **argv)
 19: {
 20:   Mat            A,B,A0,B0,mats[2]; /* problem matrices and sparser approximations */
 21:   EPS            eps;               /* eigenproblem solver context */
 22:   ST             st;
 23:   PetscInt       N,n=24,m,Istart,Iend,II,i,j;
 24:   PetscBool      flag,terse;

 26:   PetscFunctionBeginUser;
 27:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));

 29:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 30:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
 31:   if (!flag) m=n;
 32:   N = n*m;
 33:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nGHEP with split preconditioner, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));

 35:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 36:           Compute the problem matrices A and B
 37:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 39:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 40:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
 41:   PetscCall(MatSetFromOptions(A));

 43:   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
 44:   PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N));
 45:   PetscCall(MatSetFromOptions(B));

 47:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 48:   for (II=Istart;II<Iend;II++) {
 49:     i = II/n; j = II-i*n;
 50:     if (i>0) PetscCall(MatSetValue(A,II,II-n,-0.2,INSERT_VALUES));
 51:     if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-0.2,INSERT_VALUES));
 52:     if (j>0) PetscCall(MatSetValue(A,II,II-1,-3.0,INSERT_VALUES));
 53:     if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-3.0,INSERT_VALUES));
 54:     PetscCall(MatSetValue(A,II,II,7.0,INSERT_VALUES));
 55:     PetscCall(MatSetValue(B,II,II,2.0,INSERT_VALUES));
 56:   }
 57:   if (Istart==0) {
 58:     PetscCall(MatSetValue(B,0,0,6.0,INSERT_VALUES));
 59:     PetscCall(MatSetValue(B,0,1,-1.0,INSERT_VALUES));
 60:     PetscCall(MatSetValue(B,1,0,-1.0,INSERT_VALUES));
 61:     PetscCall(MatSetValue(B,1,1,1.0,INSERT_VALUES));
 62:   }
 63:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 64:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
 65:   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
 66:   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));

 68:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 69:           Compute sparser approximations A0 and B0
 70:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 72:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A0));
 73:   PetscCall(MatSetSizes(A0,PETSC_DECIDE,PETSC_DECIDE,N,N));
 74:   PetscCall(MatSetFromOptions(A0));

 76:   PetscCall(MatCreate(PETSC_COMM_WORLD,&B0));
 77:   PetscCall(MatSetSizes(B0,PETSC_DECIDE,PETSC_DECIDE,N,N));
 78:   PetscCall(MatSetFromOptions(B0));

 80:   PetscCall(MatGetOwnershipRange(A0,&Istart,&Iend));
 81:   for (II=Istart;II<Iend;II++) {
 82:     i = II/n; j = II-i*n;
 83:     if (j>0) PetscCall(MatSetValue(A0,II,II-1,-3.0,INSERT_VALUES));
 84:     if (j<n-1) PetscCall(MatSetValue(A0,II,II+1,-3.0,INSERT_VALUES));
 85:     PetscCall(MatSetValue(A0,II,II,7.0,INSERT_VALUES));
 86:     PetscCall(MatSetValue(B0,II,II,2.0,INSERT_VALUES));
 87:   }
 88:   if (Istart==0) {
 89:     PetscCall(MatSetValue(B0,0,0,6.0,INSERT_VALUES));
 90:     PetscCall(MatSetValue(B0,1,1,1.0,INSERT_VALUES));
 91:   }
 92:   PetscCall(MatAssemblyBegin(A0,MAT_FINAL_ASSEMBLY));
 93:   PetscCall(MatAssemblyEnd(A0,MAT_FINAL_ASSEMBLY));
 94:   PetscCall(MatAssemblyBegin(B0,MAT_FINAL_ASSEMBLY));
 95:   PetscCall(MatAssemblyEnd(B0,MAT_FINAL_ASSEMBLY));

 97:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 98:                 Create the eigensolver and set various options
 99:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

101:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
102:   PetscCall(EPSSetOperators(eps,A,B));
103:   PetscCall(EPSSetProblemType(eps,EPS_GHEP));
104:   PetscCall(EPSGetST(eps,&st));
105:   PetscCall(STSetType(st,STSINVERT));
106:   mats[0] = A0; mats[1] = B0;
107:   PetscCall(STSetSplitPreconditioner(st,2,mats,SUBSET_NONZERO_PATTERN));
108:   PetscCall(EPSSetTarget(eps,0.0));
109:   PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
110:   PetscCall(EPSSetFromOptions(eps));

112:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113:                  Solve the eigensystem and display solution
114:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

116:   PetscCall(EPSSolve(eps));

118:   /* show detailed info unless -terse option is given by user */
119:   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
120:   if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
121:   else {
122:     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
123:     PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
124:     PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
125:     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
126:   }
127:   PetscCall(EPSDestroy(&eps));
128:   PetscCall(MatDestroy(&A));
129:   PetscCall(MatDestroy(&B));
130:   PetscCall(MatDestroy(&A0));
131:   PetscCall(MatDestroy(&B0));
132:   PetscCall(SlepcFinalize());
133:   return 0;
134: }

136: /*TEST

138:    testset:
139:       args: -eps_nev 4 -terse
140:       output_file: output/ex49_1.out
141:       requires: !single
142:       test:
143:          suffix: 1
144:       test:
145:          suffix: 1_jd
146:          args: -eps_type jd -st_type precond
147:       test:
148:          suffix: 1_lobpcg
149:          args: -eps_type lobpcg -st_type precond -eps_smallest_real -st_shift 0.2

151:    testset:
152:       args: -eps_type ciss -eps_all -rg_type ellipse -rg_ellipse_center 0 -rg_ellipse_radius 0.34 -rg_ellipse_vscale .2 -st_ksp_type gmres -terse
153:       output_file: output/ex49_2.out
154:       test:
155:          suffix: 2
156:       test:
157:          suffix: 2_nost
158:          args: -eps_ciss_usest 0
159:          requires: !single
160:       test:
161:          suffix: 2_par
162:          nsize: 2
163:          args: -eps_ciss_partitions 2
164:          requires: !single

166: TEST*/