Actual source code: test11.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: */
10: /*
11: Example based on spring problem in NLEVP collection [1]. See the parameters
12: meaning at Example 2 in [2].
14: [1] T. Betcke, N. J. Higham, V. Mehrmann, C. Schroder, and F. Tisseur,
15: NLEVP: A Collection of Nonlinear Eigenvalue Problems, MIMS EPrint
16: 2010.98, November 2010.
17: [2] F. Tisseur, Backward error and condition of polynomial eigenvalue
18: problems, Linear Algebra and its Applications, 309 (2000), pp. 339--361,
19: April 2000.
20: */
22: static char help[] = "Illustrates the use of a user-defined stopping test.\n\n"
23: "The command line options are:\n"
24: " -n <n> ... number of grid subdivisions.\n"
25: " -mu <value> ... mass (default 1).\n"
26: " -tau <value> ... damping constant of the dampers (default 10).\n"
27: " -kappa <value> ... damping constant of the springs (default 5).\n\n";
29: #include <slepcpep.h>
31: /*
32: User-defined routines
33: */
34: PetscErrorCode MyStoppingTest(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*);
36: typedef struct {
37: PetscInt lastnconv; /* last value of nconv; used in stopping test */
38: PetscInt nreps; /* number of repetitions of nconv; used in stopping test */
39: } CTX_SPRING;
41: int main(int argc,char **argv)
42: {
43: Mat M,C,K,A[3]; /* problem matrices */
44: PEP pep; /* polynomial eigenproblem solver context */
45: RG rg; /* region object */
46: ST st;
47: CTX_SPRING *ctx;
48: PetscBool terse;
49: PetscViewer viewer;
50: PetscInt n=30,Istart,Iend,i,mpd;
51: PetscReal mu=1.0,tau=10.0,kappa=5.0;
53: SlepcInitialize(&argc,&argv,(char*)0,help);
55: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
56: PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL);
57: PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
58: PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL);
59: PetscPrintf(PETSC_COMM_WORLD,"\nDamped mass-spring system, n=%" PetscInt_FMT " mu=%g tau=%g kappa=%g\n\n",n,(double)mu,(double)tau,(double)kappa);
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65: /* K is a tridiagonal */
66: MatCreate(PETSC_COMM_WORLD,&K);
67: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
68: MatSetFromOptions(K);
69: MatSetUp(K);
70: MatGetOwnershipRange(K,&Istart,&Iend);
71: for (i=Istart;i<Iend;i++) {
72: if (i>0) MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
73: MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
74: if (i<n-1) MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
75: }
76: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
77: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
79: /* C is a tridiagonal */
80: MatCreate(PETSC_COMM_WORLD,&C);
81: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
82: MatSetFromOptions(C);
83: MatSetUp(C);
84: MatGetOwnershipRange(C,&Istart,&Iend);
85: for (i=Istart;i<Iend;i++) {
86: if (i>0) MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
87: MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
88: if (i<n-1) MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
89: }
90: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
91: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
93: /* M is a diagonal matrix */
94: MatCreate(PETSC_COMM_WORLD,&M);
95: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
96: MatSetFromOptions(M);
97: MatSetUp(M);
98: MatGetOwnershipRange(M,&Istart,&Iend);
99: for (i=Istart;i<Iend;i++) MatSetValue(M,i,i,mu,INSERT_VALUES);
100: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
101: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
103: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: Create the eigensolver and set various options
105: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107: PEPCreate(PETSC_COMM_WORLD,&pep);
108: A[0] = K; A[1] = C; A[2] = M;
109: PEPSetOperators(pep,3,A);
110: PEPSetProblemType(pep,PEP_GENERAL);
111: PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT);
113: /*
114: Define the region containing the eigenvalues of interest
115: */
116: PEPGetRG(pep,&rg);
117: RGSetType(rg,RGINTERVAL);
118: RGIntervalSetEndpoints(rg,-0.5057,-0.5052,-0.001,0.001);
119: PEPSetTarget(pep,-0.43);
120: PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE);
121: PEPGetST(pep,&st);
122: STSetType(st,STSINVERT);
124: /*
125: Set solver options. In particular, we must allocate sufficient
126: storage for all eigenpairs that may converge (ncv). This is
127: application-dependent.
128: */
129: mpd = 40;
130: PEPSetDimensions(pep,2*mpd,3*mpd,mpd);
131: PEPSetTolerances(pep,PETSC_DEFAULT,2000);
132: PetscNew(&ctx);
133: ctx->lastnconv = 0;
134: ctx->nreps = 0;
135: PEPSetStoppingTestFunction(pep,MyStoppingTest,(void*)ctx,NULL);
137: PEPSetFromOptions(pep);
139: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140: Solve the eigensystem
141: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143: PEPSolve(pep);
145: /* show detailed info unless -terse option is given by user */
146: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
147: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
148: PEPConvergedReasonView(pep,viewer);
149: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
150: if (!terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,viewer);
151: PetscViewerPopFormat(viewer);
153: PEPDestroy(&pep);
154: MatDestroy(&M);
155: MatDestroy(&C);
156: MatDestroy(&K);
157: PetscFree(ctx);
158: SlepcFinalize();
159: return 0;
160: }
162: /*
163: Function for user-defined stopping test.
165: Ignores the value of nev. It only takes into account the number of
166: eigenpairs that have converged in recent outer iterations (restarts);
167: if no new eigenvalues have converged in the last few restarts,
168: we stop the iteration, assuming that no more eigenvalues are present
169: inside the region.
170: */
171: PetscErrorCode MyStoppingTest(PEP pep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,PEPConvergedReason *reason,void *ptr)
172: {
173: CTX_SPRING *ctx = (CTX_SPRING*)ptr;
176: /* check usual termination conditions, but ignoring the case nconv>=nev */
177: PEPStoppingBasic(pep,its,max_it,nconv,PETSC_MAX_INT,reason,NULL);
178: if (*reason==PEP_CONVERGED_ITERATING) {
179: /* check if nconv is the same as before */
180: if (nconv==ctx->lastnconv) ctx->nreps++;
181: else {
182: ctx->lastnconv = nconv;
183: ctx->nreps = 0;
184: }
185: /* check if no eigenvalues converged in last 10 restarts */
186: if (nconv && ctx->nreps>10) *reason = PEP_CONVERGED_USER;
187: }
188: PetscFunctionReturn(0);
189: }
191: /*TEST
193: test:
194: args: -terse
195: suffix: 1
197: TEST*/