Actual source code: test10.c

slepc-3.21.2 2024-09-25
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[] = "Tests a user-defined convergence test in PEP (based on ex16.c).\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: /*
 19:   MyConvergedRel - Convergence test relative to the norm of M (given in ctx).
 20: */
 21: PetscErrorCode MyConvergedRel(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 22: {
 23:   PetscReal norm = *(PetscReal*)ctx;

 25:   PetscFunctionBegin;
 26:   *errest = res/norm;
 27:   PetscFunctionReturn(PETSC_SUCCESS);
 28: }

 30: int main(int argc,char **argv)
 31: {
 32:   Mat            M,C,K,A[3];      /* problem matrices */
 33:   PEP            pep;             /* polynomial eigenproblem solver context */
 34:   PetscInt       N,n=10,m,Istart,Iend,II,nev,i,j;
 35:   PetscBool      flag;
 36:   PetscReal      norm;

 38:   PetscFunctionBeginUser;
 39:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

 41:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 42:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
 43:   if (!flag) m=n;
 44:   N = n*m;
 45:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));

 47:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 48:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 49:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 51:   /* K is the 2-D Laplacian */
 52:   PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
 53:   PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N));
 54:   PetscCall(MatSetFromOptions(K));
 55:   PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
 56:   for (II=Istart;II<Iend;II++) {
 57:     i = II/n; j = II-i*n;
 58:     if (i>0) PetscCall(MatSetValue(K,II,II-n,-1.0,INSERT_VALUES));
 59:     if (i<m-1) PetscCall(MatSetValue(K,II,II+n,-1.0,INSERT_VALUES));
 60:     if (j>0) PetscCall(MatSetValue(K,II,II-1,-1.0,INSERT_VALUES));
 61:     if (j<n-1) PetscCall(MatSetValue(K,II,II+1,-1.0,INSERT_VALUES));
 62:     PetscCall(MatSetValue(K,II,II,4.0,INSERT_VALUES));
 63:   }
 64:   PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
 65:   PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));

 67:   /* C is the 1-D Laplacian on horizontal lines */
 68:   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
 69:   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N));
 70:   PetscCall(MatSetFromOptions(C));
 71:   PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
 72:   for (II=Istart;II<Iend;II++) {
 73:     i = II/n; j = II-i*n;
 74:     if (j>0) PetscCall(MatSetValue(C,II,II-1,-1.0,INSERT_VALUES));
 75:     if (j<n-1) PetscCall(MatSetValue(C,II,II+1,-1.0,INSERT_VALUES));
 76:     PetscCall(MatSetValue(C,II,II,2.0,INSERT_VALUES));
 77:   }
 78:   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
 79:   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));

 81:   /* M is a diagonal matrix */
 82:   PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
 83:   PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N));
 84:   PetscCall(MatSetFromOptions(M));
 85:   PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
 86:   for (II=Istart;II<Iend;II++) PetscCall(MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES));
 87:   PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
 88:   PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));

 90:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 91:                 Create the eigensolver and set various options
 92:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 94:   PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
 95:   A[0] = K; A[1] = C; A[2] = M;
 96:   PetscCall(PEPSetOperators(pep,3,A));
 97:   PetscCall(PEPSetProblemType(pep,PEP_HERMITIAN));
 98:   PetscCall(PEPSetDimensions(pep,4,20,PETSC_DEFAULT));

100:   /* setup convergence test relative to the norm of M */
101:   PetscCall(MatNorm(M,NORM_1,&norm));
102:   PetscCall(PEPSetConvergenceTestFunction(pep,MyConvergedRel,&norm,NULL));
103:   PetscCall(PEPSetFromOptions(pep));

105:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106:                       Solve the eigensystem
107:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

109:   PetscCall(PEPSolve(pep));
110:   PetscCall(PEPGetDimensions(pep,&nev,NULL,NULL));
111:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));

113:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114:                     Display solution and clean up
115:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

117:   PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
118:   PetscCall(PEPDestroy(&pep));
119:   PetscCall(MatDestroy(&M));
120:   PetscCall(MatDestroy(&C));
121:   PetscCall(MatDestroy(&K));
122:   PetscCall(SlepcFinalize());
123:   return 0;
124: }

126: /*TEST

128:    testset:
129:       requires: double
130:       suffix: 1

132: TEST*/