Actual source code: ex50.c

slepc-3.17.0 2022-03-31
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 quadratic 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 <slepcpep.h>

 18: int main(int argc,char **argv)
 19: {
 20:   Mat            A[3],P[3];      /* problem matrices and split preconditioner matrices */
 21:   PEP            pep;            /* polynomial eigenproblem solver context */
 22:   ST             st;
 23:   PetscInt       N,n=10,m,Istart,Iend,II,i,j;
 24:   PetscBool      flag,terse;

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

 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 30:   if (!flag) m=n;
 31:   N = n*m;
 32:   PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);

 34:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 35:      Compute the matrices for (k^2*A_2+k*A_1+A_0)x=0, and their approximations
 36:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 38:   /* A[0] is the 2-D Laplacian */
 39:   MatCreate(PETSC_COMM_WORLD,&A[0]);
 40:   MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,N,N);
 41:   MatSetFromOptions(A[0]);
 42:   MatSetUp(A[0]);
 43:   MatCreate(PETSC_COMM_WORLD,&P[0]);
 44:   MatSetSizes(P[0],PETSC_DECIDE,PETSC_DECIDE,N,N);
 45:   MatSetFromOptions(P[0]);
 46:   MatSetUp(P[0]);

 48:   MatGetOwnershipRange(A[0],&Istart,&Iend);
 49:   for (II=Istart;II<Iend;II++) {
 50:     i = II/n; j = II-i*n;
 51:     if (i>0) MatSetValue(A[0],II,II-n,-1.0,INSERT_VALUES);
 52:     if (i<m-1) MatSetValue(A[0],II,II+n,-1.0,INSERT_VALUES);
 53:     if (j>0) MatSetValue(A[0],II,II-1,-1.0,INSERT_VALUES);
 54:     if (j<n-1) MatSetValue(A[0],II,II+1,-1.0,INSERT_VALUES);
 55:     MatSetValue(A[0],II,II,4.0,INSERT_VALUES);
 56:     if (j>0) MatSetValue(P[0],II,II-1,-1.0,INSERT_VALUES);
 57:     if (j<n-1) MatSetValue(P[0],II,II+1,-1.0,INSERT_VALUES);
 58:     MatSetValue(P[0],II,II,4.0,INSERT_VALUES);
 59:   }
 60:   MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY);
 61:   MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY);
 62:   MatAssemblyBegin(P[0],MAT_FINAL_ASSEMBLY);
 63:   MatAssemblyEnd(P[0],MAT_FINAL_ASSEMBLY);

 65:   /* A[1] is the 1-D Laplacian on horizontal lines */
 66:   MatCreate(PETSC_COMM_WORLD,&A[1]);
 67:   MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,N,N);
 68:   MatSetFromOptions(A[1]);
 69:   MatSetUp(A[1]);
 70:   MatCreate(PETSC_COMM_WORLD,&P[1]);
 71:   MatSetSizes(P[1],PETSC_DECIDE,PETSC_DECIDE,N,N);
 72:   MatSetFromOptions(P[1]);
 73:   MatSetUp(P[1]);

 75:   MatGetOwnershipRange(A[1],&Istart,&Iend);
 76:   for (II=Istart;II<Iend;II++) {
 77:     i = II/n; j = II-i*n;
 78:     if (j>0) MatSetValue(A[1],II,II-1,-1.0,INSERT_VALUES);
 79:     if (j<n-1) MatSetValue(A[1],II,II+1,-1.0,INSERT_VALUES);
 80:     MatSetValue(A[1],II,II,2.0,INSERT_VALUES);
 81:     MatSetValue(P[1],II,II,2.0,INSERT_VALUES);
 82:   }
 83:   MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY);
 84:   MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY);
 85:   MatAssemblyBegin(P[1],MAT_FINAL_ASSEMBLY);
 86:   MatAssemblyEnd(P[1],MAT_FINAL_ASSEMBLY);

 88:   /* A[2] is a diagonal matrix */
 89:   MatCreate(PETSC_COMM_WORLD,&A[2]);
 90:   MatSetSizes(A[2],PETSC_DECIDE,PETSC_DECIDE,N,N);
 91:   MatSetFromOptions(A[2]);
 92:   MatSetUp(A[2]);
 93:   MatCreate(PETSC_COMM_WORLD,&P[2]);
 94:   MatSetSizes(P[2],PETSC_DECIDE,PETSC_DECIDE,N,N);
 95:   MatSetFromOptions(P[2]);
 96:   MatSetUp(P[2]);

 98:   MatGetOwnershipRange(A[2],&Istart,&Iend);
 99:   for (II=Istart;II<Iend;II++) {
100:     MatSetValue(A[2],II,II,(PetscReal)(II+1),INSERT_VALUES);
101:     MatSetValue(P[2],II,II,(PetscReal)(II+1),INSERT_VALUES);
102:   }
103:   MatAssemblyBegin(A[2],MAT_FINAL_ASSEMBLY);
104:   MatAssemblyEnd(A[2],MAT_FINAL_ASSEMBLY);
105:   MatAssemblyBegin(P[2],MAT_FINAL_ASSEMBLY);
106:   MatAssemblyEnd(P[2],MAT_FINAL_ASSEMBLY);

108:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109:                 Create the eigensolver and set various options
110:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

112:   PEPCreate(PETSC_COMM_WORLD,&pep);
113:   PEPSetOperators(pep,3,A);
114:   PEPSetProblemType(pep,PEP_HERMITIAN);

116:   PEPGetST(pep,&st);
117:   STSetType(st,STSINVERT);
118:   STSetSplitPreconditioner(st,3,P,SUBSET_NONZERO_PATTERN);

120:   PEPSetTarget(pep,-2.0);
121:   PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE);

123:   /*
124:      Set solver parameters at runtime
125:   */
126:   PEPSetFromOptions(pep);

128:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129:              Solve the eigensystem, display solution and clean up
130:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

132:   PEPSolve(pep);
133:   /* show detailed info unless -terse option is given by user */
134:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
135:   if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
136:   else {
137:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
138:     PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
139:     PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
140:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
141:   }
142:   PEPDestroy(&pep);
143:   MatDestroy(&A[0]);
144:   MatDestroy(&A[1]);
145:   MatDestroy(&A[2]);
146:   MatDestroy(&P[0]);
147:   MatDestroy(&P[1]);
148:   MatDestroy(&P[2]);
149:   SlepcFinalize();
150:   return 0;
151: }

153: /*TEST

155:    testset:
156:       args: -pep_nev 4 -pep_ncv 28 -n 12 -terse
157:       output_file: output/ex50_1.out
158:       requires: double
159:       test:
160:          suffix: 1
161:          args: -pep_type {{toar qarnoldi}}
162:       test:
163:          suffix: 1_linear
164:          args: -pep_type linear -pep_general

166:    testset:
167:       args: -pep_all -n 12 -pep_type ciss -rg_type ellipse -rg_ellipse_center -1+1.5i -rg_ellipse_radius .3 -terse
168:       output_file: output/ex50_2.out
169:       requires: complex double
170:       test:
171:          suffix: 2
172:       test:
173:          suffix: 2_par
174:          nsize: 2
175:          args: -pep_ciss_partitions 2

177: TEST*/