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