Actual source code: ex49.c

slepc-3.20.1 2023-11-27
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,(char*)0,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));
 42:   PetscCall(MatSetUp(A));

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

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

 70:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 71:           Compute sparser approximations A0 and B0
 72:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 74:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A0));
 75:   PetscCall(MatSetSizes(A0,PETSC_DECIDE,PETSC_DECIDE,N,N));
 76:   PetscCall(MatSetFromOptions(A0));
 77:   PetscCall(MatSetUp(A0));

 79:   PetscCall(MatCreate(PETSC_COMM_WORLD,&B0));
 80:   PetscCall(MatSetSizes(B0,PETSC_DECIDE,PETSC_DECIDE,N,N));
 81:   PetscCall(MatSetFromOptions(B0));
 82:   PetscCall(MatSetUp(B0));

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

101:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102:                 Create the eigensolver and set various options
103:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

116:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117:                  Solve the eigensystem and display solution
118:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

120:   PetscCall(EPSSolve(eps));

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

140: /*TEST

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

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

170: TEST*/