Actual source code: test8.c

slepc-3.21.1 2024-04-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: */

 11: static char help[] = "Test interface functions of polynomial JD.\n\n"
 12:   "This is based on ex16.c. 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             M,C,K,A[3];      /* problem matrices */
 21:   PEP             pep;             /* polynomial eigenproblem solver context */
 22:   PetscInt        N,n=10,m,Istart,Iend,II,i,j,midx;
 23:   PetscReal       restart,fix;
 24:   PetscBool       flag,reuse;
 25:   PEPJDProjection proj;

 27:   PetscFunctionBeginUser;
 28:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

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

 36:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 37:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 38:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 40:   /* K is the 2-D Laplacian */
 41:   PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
 42:   PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N));
 43:   PetscCall(MatSetFromOptions(K));
 44:   PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
 45:   for (II=Istart;II<Iend;II++) {
 46:     i = II/n; j = II-i*n;
 47:     if (i>0) PetscCall(MatSetValue(K,II,II-n,-1.0,INSERT_VALUES));
 48:     if (i<m-1) PetscCall(MatSetValue(K,II,II+n,-1.0,INSERT_VALUES));
 49:     if (j>0) PetscCall(MatSetValue(K,II,II-1,-1.0,INSERT_VALUES));
 50:     if (j<n-1) PetscCall(MatSetValue(K,II,II+1,-1.0,INSERT_VALUES));
 51:     PetscCall(MatSetValue(K,II,II,4.0,INSERT_VALUES));
 52:   }
 53:   PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
 54:   PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));

 56:   /* C is the 1-D Laplacian on horizontal lines */
 57:   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
 58:   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N));
 59:   PetscCall(MatSetFromOptions(C));
 60:   PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
 61:   for (II=Istart;II<Iend;II++) {
 62:     i = II/n; j = II-i*n;
 63:     if (j>0) PetscCall(MatSetValue(C,II,II-1,-1.0,INSERT_VALUES));
 64:     if (j<n-1) PetscCall(MatSetValue(C,II,II+1,-1.0,INSERT_VALUES));
 65:     PetscCall(MatSetValue(C,II,II,2.0,INSERT_VALUES));
 66:   }
 67:   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
 68:   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));

 70:   /* M is a diagonal matrix */
 71:   PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
 72:   PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N));
 73:   PetscCall(MatSetFromOptions(M));
 74:   PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
 75:   for (II=Istart;II<Iend;II++) PetscCall(MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES));
 76:   PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
 77:   PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));

 79:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 80:                 Create the eigensolver and set various options
 81:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 83:   PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
 84:   A[0] = K; A[1] = C; A[2] = M;
 85:   PetscCall(PEPSetOperators(pep,3,A));
 86:   PetscCall(PEPSetType(pep,PEPJD));

 88:   /*
 89:      Test interface functions of JD solver
 90:   */
 91:   PetscCall(PEPJDGetRestart(pep,&restart));
 92:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Restart parameter before changing = %g",(double)restart));
 93:   PetscCall(PEPJDSetRestart(pep,0.3));
 94:   PetscCall(PEPJDGetRestart(pep,&restart));
 95:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %g\n",(double)restart));

 97:   PetscCall(PEPJDGetFix(pep,&fix));
 98:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Fix parameter before changing = %g",(double)fix));
 99:   PetscCall(PEPJDSetFix(pep,0.001));
100:   PetscCall(PEPJDGetFix(pep,&fix));
101:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %g\n",(double)fix));

103:   PetscCall(PEPJDGetReusePreconditioner(pep,&reuse));
104:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Reuse preconditioner flag before changing = %d",(int)reuse));
105:   PetscCall(PEPJDSetReusePreconditioner(pep,PETSC_TRUE));
106:   PetscCall(PEPJDGetReusePreconditioner(pep,&reuse));
107:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)reuse));

109:   PetscCall(PEPJDGetProjection(pep,&proj));
110:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Projection type before changing = %d",(int)proj));
111:   PetscCall(PEPJDSetProjection(pep,PEP_JD_PROJECTION_ORTHOGONAL));
112:   PetscCall(PEPJDGetProjection(pep,&proj));
113:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d\n",(int)proj));

115:   PetscCall(PEPJDGetMinimalityIndex(pep,&midx));
116:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Minimality index before changing = %" PetscInt_FMT,midx));
117:   PetscCall(PEPJDSetMinimalityIndex(pep,2));
118:   PetscCall(PEPJDGetMinimalityIndex(pep,&midx));
119:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %" PetscInt_FMT "\n",midx));

121:   PetscCall(PEPSetFromOptions(pep));

123:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124:                       Solve the eigensystem
125:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

127:   PetscCall(PEPSolve(pep));
128:   PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
129:   PetscCall(PEPDestroy(&pep));
130:   PetscCall(MatDestroy(&M));
131:   PetscCall(MatDestroy(&C));
132:   PetscCall(MatDestroy(&K));
133:   PetscCall(SlepcFinalize());
134:   return 0;
135: }

137: /*TEST

139:    test:
140:       args: -n 12 -pep_nev 2 -pep_ncv 21 -pep_conv_abs

142: TEST*/