Actual source code: test11.c

slepc-3.21.1 2024-04-26
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;

 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(MatGetOwnershipRange(K,&Istart,&Iend));
 71:   for (i=Istart;i<Iend;i++) {
 72:     if (i>0) PetscCall(MatSetValue(K,i,i-1,-kappa,INSERT_VALUES));
 73:     PetscCall(MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES));
 74:     if (i<n-1) PetscCall(MatSetValue(K,i,i+1,-kappa,INSERT_VALUES));
 75:   }
 76:   PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
 77:   PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));

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

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

101:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102:                 Create the eigensolver and set various options
103:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

111:   /*
112:      Define the region containing the eigenvalues of interest
113:   */
114:   PetscCall(PEPGetRG(pep,&rg));
115:   PetscCall(RGSetType(rg,RGINTERVAL));
116:   PetscCall(RGIntervalSetEndpoints(rg,-0.5057,-0.5052,-0.001,0.001));
117:   PetscCall(PEPSetTarget(pep,-0.43));
118:   PetscCall(PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE));
119:   PetscCall(PEPGetST(pep,&st));
120:   PetscCall(STSetType(st,STSINVERT));

122:   /*
123:      Set solver options. In particular, we must allocate sufficient
124:      storage for all eigenpairs that may converge (ncv). This is
125:      application-dependent.
126:   */
127:   mpd = 40;
128:   PetscCall(PEPSetDimensions(pep,2*mpd,3*mpd,mpd));
129:   PetscCall(PEPSetTolerances(pep,PETSC_DEFAULT,2000));
130:   PetscCall(PetscNew(&ctx));
131:   ctx->lastnconv = 0;
132:   ctx->nreps     = 0;
133:   PetscCall(PEPSetStoppingTestFunction(pep,MyStoppingTest,(void*)ctx,NULL));

135:   PetscCall(PEPSetFromOptions(pep));

137:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138:                       Solve the eigensystem
139:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

141:   PetscCall(PEPSolve(pep));

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

151:   PetscCall(PEPDestroy(&pep));
152:   PetscCall(MatDestroy(&M));
153:   PetscCall(MatDestroy(&C));
154:   PetscCall(MatDestroy(&K));
155:   PetscCall(PetscFree(ctx));
156:   PetscCall(SlepcFinalize());
157:   return 0;
158: }

160: /*
161:     Function for user-defined stopping test.

163:     Ignores the value of nev. It only takes into account the number of
164:     eigenpairs that have converged in recent outer iterations (restarts);
165:     if no new eigenvalues have converged in the last few restarts,
166:     we stop the iteration, assuming that no more eigenvalues are present
167:     inside the region.
168: */
169: PetscErrorCode MyStoppingTest(PEP pep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,PEPConvergedReason *reason,void *ptr)
170: {
171:   CTX_SPRING     *ctx = (CTX_SPRING*)ptr;

173:   PetscFunctionBeginUser;
174:   /* check usual termination conditions, but ignoring the case nconv>=nev */
175:   PetscCall(PEPStoppingBasic(pep,its,max_it,nconv,PETSC_MAX_INT,reason,NULL));
176:   if (*reason==PEP_CONVERGED_ITERATING) {
177:     /* check if nconv is the same as before */
178:     if (nconv==ctx->lastnconv) ctx->nreps++;
179:     else {
180:       ctx->lastnconv = nconv;
181:       ctx->nreps     = 0;
182:     }
183:     /* check if no eigenvalues converged in last 10 restarts */
184:     if (nconv && ctx->nreps>10) *reason = PEP_CONVERGED_USER;
185:   }
186:   PetscFunctionReturn(PETSC_SUCCESS);
187: }

189: /*TEST

191:    test:
192:       args: -terse
193:       requires: !single
194:       suffix: 1

196: TEST*/