Actual source code: test9.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[] = "Illustrates use of PEPSetEigenvalueComparison().\n\n"
12: "Based on butterfly.c. The command line options are:\n"
13: " -m <m>, grid size, the dimension of the matrices is n=m*m.\n"
14: " -c <array>, problem parameters, must be 10 comma-separated real values.\n\n";
16: #include <slepcpep.h>
18: #define NMAT 5
20: /*
21: Function for user-defined eigenvalue ordering criterion.
23: Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
24: one of them as the preferred one according to the criterion.
25: In this example, eigenvalues are sorted by magnitude but those with
26: positive real part are preferred.
27: */
28: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
29: {
30: PetscReal rea,reb;
32: PetscFunctionBeginUser;
33: #if defined(PETSC_USE_COMPLEX)
34: rea = PetscRealPart(ar); reb = PetscRealPart(br);
35: #else
36: rea = ar; reb = br;
37: #endif
38: *r = rea<0.0? 1: (reb<0.0? -1: PetscSign(SlepcAbsEigenvalue(br,bi)-SlepcAbsEigenvalue(ar,ai)));
39: PetscFunctionReturn(PETSC_SUCCESS);
40: }
42: int main(int argc,char **argv)
43: {
44: Mat A[NMAT]; /* problem matrices */
45: PEP pep; /* polynomial eigenproblem solver context */
46: PetscInt n,m=8,k,II,Istart,Iend,i,j;
47: PetscReal c[10] = { 0.6, 1.3, 1.3, 0.1, 0.1, 1.2, 1.0, 1.0, 1.2, 1.0 };
48: PetscBool flg;
49: PetscBool terse;
50: const char *prefix;
52: PetscFunctionBeginUser;
53: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
55: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
56: n = m*m;
57: k = 10;
58: PetscCall(PetscOptionsGetRealArray(NULL,NULL,"-c",c,&k,&flg));
59: PetscCheck(!flg || k==10,PETSC_COMM_WORLD,PETSC_ERR_USER,"The number of parameters -c should be 10, you provided %" PetscInt_FMT,k);
60: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nButterfly problem, n=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n\n",n,m));
62: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: Compute the polynomial matrices
64: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66: /* initialize matrices */
67: for (i=0;i<NMAT;i++) {
68: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[i]));
69: PetscCall(MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n));
70: PetscCall(MatSetFromOptions(A[i]));
71: PetscCall(MatSetUp(A[i]));
72: }
73: PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
75: /* A0 */
76: for (II=Istart;II<Iend;II++) {
77: i = II/m; j = II-i*m;
78: PetscCall(MatSetValue(A[0],II,II,4.0*c[0]/6.0+4.0*c[1]/6.0,INSERT_VALUES));
79: if (j>0) PetscCall(MatSetValue(A[0],II,II-1,c[0]/6.0,INSERT_VALUES));
80: if (j<m-1) PetscCall(MatSetValue(A[0],II,II+1,c[0]/6.0,INSERT_VALUES));
81: if (i>0) PetscCall(MatSetValue(A[0],II,II-m,c[1]/6.0,INSERT_VALUES));
82: if (i<m-1) PetscCall(MatSetValue(A[0],II,II+m,c[1]/6.0,INSERT_VALUES));
83: }
85: /* A1 */
86: for (II=Istart;II<Iend;II++) {
87: i = II/m; j = II-i*m;
88: if (j>0) PetscCall(MatSetValue(A[1],II,II-1,c[2],INSERT_VALUES));
89: if (j<m-1) PetscCall(MatSetValue(A[1],II,II+1,-c[2],INSERT_VALUES));
90: if (i>0) PetscCall(MatSetValue(A[1],II,II-m,c[3],INSERT_VALUES));
91: if (i<m-1) PetscCall(MatSetValue(A[1],II,II+m,-c[3],INSERT_VALUES));
92: }
94: /* A2 */
95: for (II=Istart;II<Iend;II++) {
96: i = II/m; j = II-i*m;
97: PetscCall(MatSetValue(A[2],II,II,-2.0*c[4]-2.0*c[5],INSERT_VALUES));
98: if (j>0) PetscCall(MatSetValue(A[2],II,II-1,c[4],INSERT_VALUES));
99: if (j<m-1) PetscCall(MatSetValue(A[2],II,II+1,c[4],INSERT_VALUES));
100: if (i>0) PetscCall(MatSetValue(A[2],II,II-m,c[5],INSERT_VALUES));
101: if (i<m-1) PetscCall(MatSetValue(A[2],II,II+m,c[5],INSERT_VALUES));
102: }
104: /* A3 */
105: for (II=Istart;II<Iend;II++) {
106: i = II/m; j = II-i*m;
107: if (j>0) PetscCall(MatSetValue(A[3],II,II-1,c[6],INSERT_VALUES));
108: if (j<m-1) PetscCall(MatSetValue(A[3],II,II+1,-c[6],INSERT_VALUES));
109: if (i>0) PetscCall(MatSetValue(A[3],II,II-m,c[7],INSERT_VALUES));
110: if (i<m-1) PetscCall(MatSetValue(A[3],II,II+m,-c[7],INSERT_VALUES));
111: }
113: /* A4 */
114: for (II=Istart;II<Iend;II++) {
115: i = II/m; j = II-i*m;
116: PetscCall(MatSetValue(A[4],II,II,2.0*c[8]+2.0*c[9],INSERT_VALUES));
117: if (j>0) PetscCall(MatSetValue(A[4],II,II-1,-c[8],INSERT_VALUES));
118: if (j<m-1) PetscCall(MatSetValue(A[4],II,II+1,-c[8],INSERT_VALUES));
119: if (i>0) PetscCall(MatSetValue(A[4],II,II-m,-c[9],INSERT_VALUES));
120: if (i<m-1) PetscCall(MatSetValue(A[4],II,II+m,-c[9],INSERT_VALUES));
121: }
123: /* assemble matrices */
124: for (i=0;i<NMAT;i++) PetscCall(MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY));
125: for (i=0;i<NMAT;i++) PetscCall(MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY));
127: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: Create the eigensolver and solve the problem
129: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
132: PetscCall(PEPSetOptionsPrefix(pep,"check_"));
133: PetscCall(PEPAppendOptionsPrefix(pep,"myprefix_"));
134: PetscCall(PEPGetOptionsPrefix(pep,&prefix));
135: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"PEP prefix is currently: %s\n\n",prefix));
137: PetscCall(PEPSetOperators(pep,NMAT,A));
138: PetscCall(PEPSetEigenvalueComparison(pep,MyEigenSort,NULL));
139: PetscCall(PEPSetFromOptions(pep));
140: PetscCall(PEPSolve(pep));
142: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143: Display solution and clean up
144: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146: /* show detailed info unless -terse option is given by user */
147: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
148: if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
149: else {
150: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
151: PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
152: PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD));
153: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
154: }
155: PetscCall(PEPDestroy(&pep));
156: for (i=0;i<NMAT;i++) PetscCall(MatDestroy(&A[i]));
157: PetscCall(SlepcFinalize());
158: return 0;
159: }
161: /*TEST
163: test:
164: args: -check_myprefix_pep_nev 4 -terse
165: requires: double
167: TEST*/