Actual source code: ex46.c
slepc-main 2024-12-17
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[] = "Illustrates passing a sparser matrix to build the preconditioner.\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 <slepceps.h>
18: int main(int argc,char **argv)
19: {
20: Mat A,A0; /* operator matrix */
21: EPS eps; /* eigenproblem solver context */
22: ST st;
23: PetscInt N,n=24,m,Istart,Iend,II,i,j;
24: PetscBool flag,terse;
26: PetscFunctionBeginUser;
27: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
29: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
30: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
31: if (!flag) m=n;
32: N = n*m;
33: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nModified 2-D Laplacian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
35: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: Compute the operator matrix A and a sparser approximation A0
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
40: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
41: PetscCall(MatSetFromOptions(A));
42: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
43: for (II=Istart;II<Iend;II++) {
44: i = II/n; j = II-i*n;
45: if (i>0) PetscCall(MatSetValue(A,II,II-n,-0.2,INSERT_VALUES));
46: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-0.2,INSERT_VALUES));
47: if (j>0) PetscCall(MatSetValue(A,II,II-1,-3.0,INSERT_VALUES));
48: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-3.0,INSERT_VALUES));
49: PetscCall(MatSetValue(A,II,II,7.0,INSERT_VALUES));
50: }
51: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
52: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
54: PetscCall(MatCreate(PETSC_COMM_WORLD,&A0));
55: PetscCall(MatSetSizes(A0,PETSC_DECIDE,PETSC_DECIDE,N,N));
56: PetscCall(MatSetFromOptions(A0));
57: PetscCall(MatGetOwnershipRange(A0,&Istart,&Iend));
58: for (II=Istart;II<Iend;II++) {
59: i = II/n; j = II-i*n;
60: if (j>0) PetscCall(MatSetValue(A0,II,II-1,-3.0,INSERT_VALUES));
61: if (j<n-1) PetscCall(MatSetValue(A0,II,II+1,-3.0,INSERT_VALUES));
62: PetscCall(MatSetValue(A0,II,II,7.0,INSERT_VALUES));
63: }
64: PetscCall(MatAssemblyBegin(A0,MAT_FINAL_ASSEMBLY));
65: PetscCall(MatAssemblyEnd(A0,MAT_FINAL_ASSEMBLY));
67: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68: Create the eigensolver and set various options
69: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
72: PetscCall(EPSSetOperators(eps,A,NULL));
73: PetscCall(EPSSetProblemType(eps,EPS_HEP));
74: PetscCall(EPSGetST(eps,&st));
75: PetscCall(STSetType(st,STSINVERT));
76: PetscCall(STSetPreconditionerMat(st,A0));
77: PetscCall(EPSSetTarget(eps,0.0));
78: PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
79: PetscCall(EPSSetFromOptions(eps));
81: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82: Solve the eigensystem and display solution
83: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85: PetscCall(EPSSolve(eps));
87: /* show detailed info unless -terse option is given by user */
88: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
89: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
90: else {
91: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
92: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
93: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
94: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
95: }
96: PetscCall(EPSDestroy(&eps));
97: PetscCall(MatDestroy(&A));
98: PetscCall(MatDestroy(&A0));
99: PetscCall(SlepcFinalize());
100: return 0;
101: }
103: /*TEST
105: testset:
106: args: -eps_nev 4 -terse
107: output_file: output/ex46_1.out
108: requires: !single
109: test:
110: suffix: 1
111: test:
112: suffix: 2
113: args: -st_ksp_type bcgs -st_pc_type bjacobi
115: TEST*/