Actual source code: test9.c
slepc-3.17.2 2022-08-09
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;
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(0);
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: SlepcInitialize(&argc,&argv,(char*)0,help);
54: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
55: n = m*m;
56: k = 10;
57: PetscOptionsGetRealArray(NULL,NULL,"-c",c,&k,&flg);
59: PetscPrintf(PETSC_COMM_WORLD,"\nButterfly problem, n=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n\n",n,m);
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Compute the polynomial matrices
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65: /* initialize matrices */
66: for (i=0;i<NMAT;i++) {
67: MatCreate(PETSC_COMM_WORLD,&A[i]);
68: MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n);
69: MatSetFromOptions(A[i]);
70: MatSetUp(A[i]);
71: }
72: MatGetOwnershipRange(A[0],&Istart,&Iend);
74: /* A0 */
75: for (II=Istart;II<Iend;II++) {
76: i = II/m; j = II-i*m;
77: MatSetValue(A[0],II,II,4.0*c[0]/6.0+4.0*c[1]/6.0,INSERT_VALUES);
78: if (j>0) MatSetValue(A[0],II,II-1,c[0]/6.0,INSERT_VALUES);
79: if (j<m-1) MatSetValue(A[0],II,II+1,c[0]/6.0,INSERT_VALUES);
80: if (i>0) MatSetValue(A[0],II,II-m,c[1]/6.0,INSERT_VALUES);
81: if (i<m-1) MatSetValue(A[0],II,II+m,c[1]/6.0,INSERT_VALUES);
82: }
84: /* A1 */
85: for (II=Istart;II<Iend;II++) {
86: i = II/m; j = II-i*m;
87: if (j>0) MatSetValue(A[1],II,II-1,c[2],INSERT_VALUES);
88: if (j<m-1) MatSetValue(A[1],II,II+1,-c[2],INSERT_VALUES);
89: if (i>0) MatSetValue(A[1],II,II-m,c[3],INSERT_VALUES);
90: if (i<m-1) MatSetValue(A[1],II,II+m,-c[3],INSERT_VALUES);
91: }
93: /* A2 */
94: for (II=Istart;II<Iend;II++) {
95: i = II/m; j = II-i*m;
96: MatSetValue(A[2],II,II,-2.0*c[4]-2.0*c[5],INSERT_VALUES);
97: if (j>0) MatSetValue(A[2],II,II-1,c[4],INSERT_VALUES);
98: if (j<m-1) MatSetValue(A[2],II,II+1,c[4],INSERT_VALUES);
99: if (i>0) MatSetValue(A[2],II,II-m,c[5],INSERT_VALUES);
100: if (i<m-1) MatSetValue(A[2],II,II+m,c[5],INSERT_VALUES);
101: }
103: /* A3 */
104: for (II=Istart;II<Iend;II++) {
105: i = II/m; j = II-i*m;
106: if (j>0) MatSetValue(A[3],II,II-1,c[6],INSERT_VALUES);
107: if (j<m-1) MatSetValue(A[3],II,II+1,-c[6],INSERT_VALUES);
108: if (i>0) MatSetValue(A[3],II,II-m,c[7],INSERT_VALUES);
109: if (i<m-1) MatSetValue(A[3],II,II+m,-c[7],INSERT_VALUES);
110: }
112: /* A4 */
113: for (II=Istart;II<Iend;II++) {
114: i = II/m; j = II-i*m;
115: MatSetValue(A[4],II,II,2.0*c[8]+2.0*c[9],INSERT_VALUES);
116: if (j>0) MatSetValue(A[4],II,II-1,-c[8],INSERT_VALUES);
117: if (j<m-1) MatSetValue(A[4],II,II+1,-c[8],INSERT_VALUES);
118: if (i>0) MatSetValue(A[4],II,II-m,-c[9],INSERT_VALUES);
119: if (i<m-1) MatSetValue(A[4],II,II+m,-c[9],INSERT_VALUES);
120: }
122: /* assemble matrices */
123: for (i=0;i<NMAT;i++) MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY);
124: for (i=0;i<NMAT;i++) MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY);
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Create the eigensolver and solve the problem
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: PEPCreate(PETSC_COMM_WORLD,&pep);
131: PEPSetOptionsPrefix(pep,"check_");
132: PEPAppendOptionsPrefix(pep,"myprefix_");
133: PEPGetOptionsPrefix(pep,&prefix);
134: PetscPrintf(PETSC_COMM_WORLD,"PEP prefix is currently: %s\n\n",prefix);
136: PEPSetOperators(pep,NMAT,A);
137: PEPSetEigenvalueComparison(pep,MyEigenSort,NULL);
138: PEPSetFromOptions(pep);
139: PEPSolve(pep);
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Display solution and clean up
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145: /* show detailed info unless -terse option is given by user */
146: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
147: if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
148: else {
149: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
150: PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
151: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
152: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
153: }
154: PEPDestroy(&pep);
155: for (i=0;i<NMAT;i++) MatDestroy(&A[i]);
156: SlepcFinalize();
157: return 0;
158: }
160: /*TEST
162: test:
163: args: -check_myprefix_pep_nev 4 -terse
164: requires: double
166: TEST*/