Actual source code: test35.c
slepc-3.19.0 2023-03-31
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 to external package BLOPEX.\n\n"
12: "This is based on ex12.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 <slepceps.h>
18: int main(int argc,char **argv)
19: {
20: Mat A,B; /* matrices */
21: EPS eps; /* eigenproblem solver context */
22: ST st; /* spectral transformation context */
23: KSP ksp;
24: PC pc;
25: PetscInt N,n=35,m,Istart,Iend,II,i,j,bs;
26: PetscBool flag;
28: PetscFunctionBeginUser;
29: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
31: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
32: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
33: if (!flag) m=n;
34: N = n*m;
35: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized eigenproblem with BLOPEX, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
37: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38: Compute the matrices that define the eigensystem, Ax=kBx
39: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
41: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
42: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
43: PetscCall(MatSetFromOptions(A));
44: PetscCall(MatSetUp(A));
46: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
47: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N));
48: PetscCall(MatSetFromOptions(B));
49: PetscCall(MatSetUp(B));
51: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
52: for (II=Istart;II<Iend;II++) {
53: i = II/n; j = II-i*n;
54: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
55: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
56: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
57: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
58: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
59: PetscCall(MatSetValue(B,II,II,2.0,INSERT_VALUES));
60: }
61: if (Istart==0) {
62: PetscCall(MatSetValue(B,0,0,6.0,INSERT_VALUES));
63: PetscCall(MatSetValue(B,0,1,-1.0,INSERT_VALUES));
64: PetscCall(MatSetValue(B,1,0,-1.0,INSERT_VALUES));
65: PetscCall(MatSetValue(B,1,1,1.0,INSERT_VALUES));
66: }
68: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
69: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
70: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
71: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
73: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74: Create the eigensolver and set various options
75: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
78: PetscCall(EPSSetOperators(eps,A,B));
79: PetscCall(EPSSetProblemType(eps,EPS_GHEP));
80: PetscCall(EPSSetType(eps,EPSBLOPEX));
82: /*
83: Set several options
84: */
85: PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
86: PetscCall(EPSGetST(eps,&st));
87: PetscCall(STSetType(st,STPRECOND));
88: PetscCall(STGetKSP(st,&ksp));
89: PetscCall(KSPGetPC(ksp,&pc));
90: PetscCall(KSPSetType(ksp,KSPPREONLY));
91: PetscCall(PCSetType(pc,PCICC));
93: PetscCall(EPSBLOPEXSetBlockSize(eps,4));
94: PetscCall(EPSSetFromOptions(eps));
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Compute all eigenvalues in interval and display info
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100: PetscCall(EPSSolve(eps));
101: PetscCall(EPSBLOPEXGetBlockSize(eps,&bs));
102: PetscCall(PetscPrintf(PETSC_COMM_WORLD," BLOPEX: using block size %" PetscInt_FMT "\n",bs));
104: PetscCall(EPSErrorView(eps,EPS_ERROR_ABSOLUTE,NULL));
106: PetscCall(EPSDestroy(&eps));
107: PetscCall(MatDestroy(&A));
108: PetscCall(MatDestroy(&B));
109: PetscCall(SlepcFinalize());
110: return 0;
111: }
113: /*TEST
115: build:
116: requires: blopex
118: test:
119: suffix: 1
120: args: -eps_nev 8 -eps_view
121: requires: blopex
122: filter: grep -v tolerance | sed -e "s/hermitian/symmetric/"
124: TEST*/