Actual source code: test11.c

slepc-3.18.3 2023-03-24
Report Typos and Errors
  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;

 54:   SlepcInitialize(&argc,&argv,(char*)0,help);

 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=%" 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:   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) MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
 74:     MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
 75:     if (i<n-1) MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
 76:   }
 77:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 78:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 80:   /* C is a tridiagonal */
 81:   MatCreate(PETSC_COMM_WORLD,&C);
 82:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 83:   MatSetFromOptions(C);
 84:   MatSetUp(C);
 85:   MatGetOwnershipRange(C,&Istart,&Iend);
 86:   for (i=Istart;i<Iend;i++) {
 87:     if (i>0) MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
 88:     MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
 89:     if (i<n-1) MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
 90:   }
 91:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 92:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 94:   /* M is a diagonal matrix */
 95:   MatCreate(PETSC_COMM_WORLD,&M);
 96:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
 97:   MatSetFromOptions(M);
 98:   MatSetUp(M);
 99:   MatGetOwnershipRange(M,&Istart,&Iend);
100:   for (i=Istart;i<Iend;i++) MatSetValue(M,i,i,mu,INSERT_VALUES);
101:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
102:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

104:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105:                 Create the eigensolver and set various options
106:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

108:   PEPCreate(PETSC_COMM_WORLD,&pep);
109:   A[0] = K; A[1] = C; A[2] = M;
110:   PEPSetOperators(pep,3,A);
111:   PEPSetProblemType(pep,PEP_GENERAL);
112:   PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT);

114:   /*
115:      Define the region containing the eigenvalues of interest
116:   */
117:   PEPGetRG(pep,&rg);
118:   RGSetType(rg,RGINTERVAL);
119:   RGIntervalSetEndpoints(rg,-0.5057,-0.5052,-0.001,0.001);
120:   PEPSetTarget(pep,-0.43);
121:   PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE);
122:   PEPGetST(pep,&st);
123:   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:   PEPSetDimensions(pep,2*mpd,3*mpd,mpd);
132:   PEPSetTolerances(pep,PETSC_DEFAULT,2000);
133:   PetscNew(&ctx);
134:   ctx->lastnconv = 0;
135:   ctx->nreps     = 0;
136:   PEPSetStoppingTestFunction(pep,MyStoppingTest,(void*)ctx,NULL);

138:   PEPSetFromOptions(pep);

140:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141:                       Solve the eigensystem
142:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

144:   PEPSolve(pep);

146:   /* show detailed info unless -terse option is given by user */
147:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
148:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
149:   PEPConvergedReasonView(pep,viewer);
150:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
151:   if (!terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,viewer);
152:   PetscViewerPopFormat(viewer);

154:   PEPDestroy(&pep);
155:   MatDestroy(&M);
156:   MatDestroy(&C);
157:   MatDestroy(&K);
158:   PetscFree(ctx);
159:   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;

177:   /* check usual termination conditions, but ignoring the case nconv>=nev */
178:   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:   return 0;
190: }

192: /*TEST

194:    test:
195:       args: -terse
196:       suffix: 1

198: TEST*/