Actual source code: test11.c
slepc-3.16.3 2022-04-11
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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;
50: PetscViewer viewer;
51: PetscInt n=30,Istart,Iend,i,mpd;
52: PetscReal mu=1.0,tau=10.0,kappa=5.0;
54: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
56: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
57: PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL);
58: PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
59: PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL);
60: PetscPrintf(PETSC_COMM_WORLD,"\nDamped mass-spring system, n=%D mu=%g tau=%g kappa=%g\n\n",n,(double)mu,(double)tau,(double)kappa);
62: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
64: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66: /* K is a tridiagonal */
67: MatCreate(PETSC_COMM_WORLD,&K);
68: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
69: MatSetFromOptions(K);
70: MatSetUp(K);
71: MatGetOwnershipRange(K,&Istart,&Iend);
72: for (i=Istart;i<Iend;i++) {
73: if (i>0) {
74: MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
75: }
76: MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
77: if (i<n-1) {
78: MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
79: }
80: }
81: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
82: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
84: /* C is a tridiagonal */
85: MatCreate(PETSC_COMM_WORLD,&C);
86: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
87: MatSetFromOptions(C);
88: MatSetUp(C);
89: MatGetOwnershipRange(C,&Istart,&Iend);
90: for (i=Istart;i<Iend;i++) {
91: if (i>0) {
92: MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
93: }
94: MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
95: if (i<n-1) {
96: MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
97: }
98: }
99: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
100: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
102: /* M is a diagonal matrix */
103: MatCreate(PETSC_COMM_WORLD,&M);
104: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
105: MatSetFromOptions(M);
106: MatSetUp(M);
107: MatGetOwnershipRange(M,&Istart,&Iend);
108: for (i=Istart;i<Iend;i++) {
109: MatSetValue(M,i,i,mu,INSERT_VALUES);
110: }
111: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
112: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
114: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115: Create the eigensolver and set various options
116: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118: PEPCreate(PETSC_COMM_WORLD,&pep);
119: A[0] = K; A[1] = C; A[2] = M;
120: PEPSetOperators(pep,3,A);
121: PEPSetProblemType(pep,PEP_GENERAL);
122: PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT);
124: /*
125: Define the region containing the eigenvalues of interest
126: */
127: PEPGetRG(pep,&rg);
128: RGSetType(rg,RGINTERVAL);
129: RGIntervalSetEndpoints(rg,-0.5057,-0.5052,-0.001,0.001);
130: PEPSetTarget(pep,-0.43);
131: PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE);
132: PEPGetST(pep,&st);
133: STSetType(st,STSINVERT);
135: /*
136: Set solver options. In particular, we must allocate sufficient
137: storage for all eigenpairs that may converge (ncv). This is
138: application-dependent.
139: */
140: mpd = 40;
141: PEPSetDimensions(pep,2*mpd,3*mpd,mpd);
142: PEPSetTolerances(pep,PETSC_DEFAULT,2000);
143: PetscNew(&ctx);
144: ctx->lastnconv = 0;
145: ctx->nreps = 0;
146: PEPSetStoppingTestFunction(pep,MyStoppingTest,(void*)ctx,NULL);
148: PEPSetFromOptions(pep);
150: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151: Solve the eigensystem
152: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154: PEPSolve(pep);
156: /* show detailed info unless -terse option is given by user */
157: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
158: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
159: PEPConvergedReasonView(pep,viewer);
160: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
161: if (!terse) {
162: PEPErrorView(pep,PEP_ERROR_BACKWARD,viewer);
163: }
164: PetscViewerPopFormat(viewer);
166: PEPDestroy(&pep);
167: MatDestroy(&M);
168: MatDestroy(&C);
169: MatDestroy(&K);
170: PetscFree(ctx);
171: SlepcFinalize();
172: return ierr;
173: }
175: /*
176: Function for user-defined stopping test.
178: Ignores the value of nev. It only takes into account the number of
179: eigenpairs that have converged in recent outer iterations (restarts);
180: if no new eigenvalues have converged in the last few restarts,
181: we stop the iteration, assuming that no more eigenvalues are present
182: inside the region.
183: */
184: PetscErrorCode MyStoppingTest(PEP pep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,PEPConvergedReason *reason,void *ptr)
185: {
187: CTX_SPRING *ctx = (CTX_SPRING*)ptr;
190: /* check usual termination conditions, but ignoring the case nconv>=nev */
191: PEPStoppingBasic(pep,its,max_it,nconv,PETSC_MAX_INT,reason,NULL);
192: if (*reason==PEP_CONVERGED_ITERATING) {
193: /* check if nconv is the same as before */
194: if (nconv==ctx->lastnconv) ctx->nreps++;
195: else {
196: ctx->lastnconv = nconv;
197: ctx->nreps = 0;
198: }
199: /* check if no eigenvalues converged in last 10 restarts */
200: if (nconv && ctx->nreps>10) *reason = PEP_CONVERGED_USER;
201: }
202: return(0);
203: }
205: /*TEST
207: test:
208: args: -terse
209: suffix: 1
211: TEST*/