Actual source code: test4.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[] = "Solve a quadratic problem with PEPLINEAR with a user-provided EPS.\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            M,C,K,A[3];
 21:   PEP            pep;
 22:   PetscInt       N,n=10,m,Istart,Iend,II,i,j;
 23:   PetscBool      flag,expmat;
 24:   PetscReal      alpha,beta;
 25:   EPS            eps;
 26:   ST             st;
 27:   KSP            ksp;
 28:   PC             pc;

 30:   PetscFunctionBeginUser;
 31:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 32:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 33:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
 34:   if (!flag) m=n;
 35:   N = n*m;
 36:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));

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

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

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

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

 81:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 82:              Create a standalone EPS with appropriate settings
 83:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 85:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 86:   PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
 87: #if defined(PETSC_USE_COMPLEX)
 88:   PetscCall(EPSSetTarget(eps,0.01*PETSC_i));
 89: #endif
 90:   PetscCall(EPSGetST(eps,&st));
 91:   PetscCall(STSetType(st,STSINVERT));
 92:   PetscCall(STGetKSP(st,&ksp));
 93:   PetscCall(KSPSetType(ksp,KSPBCGS));
 94:   PetscCall(KSPGetPC(ksp,&pc));
 95:   PetscCall(PCSetType(pc,PCJACOBI));
 96:   PetscCall(EPSSetFromOptions(eps));

 98:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 99:              Create the eigensolver and solve the eigensystem
100:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

102:   PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
103:   PetscCall(PetscObjectSetName((PetscObject)pep,"PEP_solver"));
104:   A[0] = K; A[1] = C; A[2] = M;
105:   PetscCall(PEPSetOperators(pep,3,A));
106:   PetscCall(PEPSetType(pep,PEPLINEAR));
107:   PetscCall(PEPSetProblemType(pep,PEP_GENERAL));
108:   PetscCall(PEPLinearSetEPS(pep,eps));
109:   PetscCall(PEPSetFromOptions(pep));
110:   PetscCall(PEPSolve(pep));
111:   PetscCall(PEPLinearGetLinearization(pep,&alpha,&beta));
112:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Linearization with alpha=%g, beta=%g",(double)alpha,(double)beta));
113:   PetscCall(PEPLinearGetExplicitMatrix(pep,&expmat));
114:   if (expmat) PetscCall(PetscPrintf(PETSC_COMM_WORLD," with explicit matrix"));
115:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));

117:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118:                     Display solution and clean up
119:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

121:   PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
122:   PetscCall(PEPDestroy(&pep));
123:   PetscCall(EPSDestroy(&eps));
124:   PetscCall(MatDestroy(&M));
125:   PetscCall(MatDestroy(&C));
126:   PetscCall(MatDestroy(&K));
127:   PetscCall(SlepcFinalize());
128:   return 0;
129: }

131: /*TEST

133:    testset:
134:       args: -pep_linear_explicitmatrix -pep_view_vectors ::ascii_info
135:       test:
136:          suffix: 1_real
137:          requires: !single !complex
138:       test:
139:          suffix: 1
140:          requires: !single complex

142: TEST*/