LCOV - code coverage report
Current view: top level - pep/tests - test11.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 86 86 100.0 %
Date: 2024-11-23 00:39:48 Functions: 2 2 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /*
       2             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       3             :    SLEPc - Scalable Library for Eigenvalue Problem Computations
       4             :    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
       5             : 
       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].
      13             : 
      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             : */
      21             : 
      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";
      28             : 
      29             : #include <slepcpep.h>
      30             : 
      31             : /*
      32             :    User-defined routines
      33             : */
      34             : PetscErrorCode MyStoppingTest(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*);
      35             : 
      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;
      40             : 
      41           1 : int main(int argc,char **argv)
      42             : {
      43           1 :   Mat            M,C,K,A[3];      /* problem matrices */
      44           1 :   PEP            pep;             /* polynomial eigenproblem solver context */
      45           1 :   RG             rg;              /* region object */
      46           1 :   ST             st;
      47           1 :   CTX_SPRING     *ctx;
      48           1 :   PetscBool      terse;
      49           1 :   PetscViewer    viewer;
      50           1 :   PetscInt       n=30,Istart,Iend,i,mpd;
      51           1 :   PetscReal      mu=1.0,tau=10.0,kappa=5.0;
      52             : 
      53           1 :   PetscFunctionBeginUser;
      54           1 :   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
      55             : 
      56           1 :   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
      57           1 :   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL));
      58           1 :   PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
      59           1 :   PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL));
      60           1 :   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));
      61             : 
      62             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      63             :      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
      64             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      65             : 
      66             :   /* K is a tridiagonal */
      67           1 :   PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
      68           1 :   PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
      69           1 :   PetscCall(MatSetFromOptions(K));
      70           1 :   PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
      71          31 :   for (i=Istart;i<Iend;i++) {
      72          30 :     if (i>0) PetscCall(MatSetValue(K,i,i-1,-kappa,INSERT_VALUES));
      73          30 :     PetscCall(MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES));
      74          30 :     if (i<n-1) PetscCall(MatSetValue(K,i,i+1,-kappa,INSERT_VALUES));
      75             :   }
      76           1 :   PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
      77           1 :   PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
      78             : 
      79             :   /* C is a tridiagonal */
      80           1 :   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
      81           1 :   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n));
      82           1 :   PetscCall(MatSetFromOptions(C));
      83           1 :   PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
      84          31 :   for (i=Istart;i<Iend;i++) {
      85          30 :     if (i>0) PetscCall(MatSetValue(C,i,i-1,-tau,INSERT_VALUES));
      86          30 :     PetscCall(MatSetValue(C,i,i,tau*3.0,INSERT_VALUES));
      87          30 :     if (i<n-1) PetscCall(MatSetValue(C,i,i+1,-tau,INSERT_VALUES));
      88             :   }
      89           1 :   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
      90           1 :   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
      91             : 
      92             :   /* M is a diagonal matrix */
      93           1 :   PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
      94           1 :   PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
      95           1 :   PetscCall(MatSetFromOptions(M));
      96           1 :   PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
      97          31 :   for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(M,i,i,mu,INSERT_VALUES));
      98           1 :   PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
      99           1 :   PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
     100             : 
     101             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
     102             :                 Create the eigensolver and set various options
     103             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     104             : 
     105           1 :   PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
     106           1 :   A[0] = K; A[1] = C; A[2] = M;
     107           1 :   PetscCall(PEPSetOperators(pep,3,A));
     108           1 :   PetscCall(PEPSetProblemType(pep,PEP_GENERAL));
     109             : 
     110             :   /*
     111             :      Define the region containing the eigenvalues of interest
     112             :   */
     113           1 :   PetscCall(PEPGetRG(pep,&rg));
     114           1 :   PetscCall(RGSetType(rg,RGINTERVAL));
     115           1 :   PetscCall(RGIntervalSetEndpoints(rg,-0.5057,-0.5052,-0.001,0.001));
     116           1 :   PetscCall(PEPSetTarget(pep,-0.43));
     117           1 :   PetscCall(PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE));
     118           1 :   PetscCall(PEPGetST(pep,&st));
     119           1 :   PetscCall(STSetType(st,STSINVERT));
     120             : 
     121             :   /*
     122             :      Set solver options. In particular, we must allocate sufficient
     123             :      storage for all eigenpairs that may converge (ncv). This is
     124             :      application-dependent.
     125             :   */
     126           1 :   mpd = 40;
     127           1 :   PetscCall(PEPSetDimensions(pep,2*mpd,3*mpd,mpd));
     128           1 :   PetscCall(PEPSetTolerances(pep,SLEPC_DEFAULT_TOL,2000));
     129           1 :   PetscCall(PetscNew(&ctx));
     130           1 :   ctx->lastnconv = 0;
     131           1 :   ctx->nreps     = 0;
     132           1 :   PetscCall(PEPSetStoppingTestFunction(pep,MyStoppingTest,(void*)ctx,NULL));
     133             : 
     134           1 :   PetscCall(PEPSetFromOptions(pep));
     135             : 
     136             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
     137             :                       Solve the eigensystem
     138             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     139             : 
     140           1 :   PetscCall(PEPSolve(pep));
     141             : 
     142             :   /* show detailed info unless -terse option is given by user */
     143           1 :   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
     144           1 :   PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
     145           1 :   PetscCall(PEPConvergedReasonView(pep,viewer));
     146           1 :   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
     147           1 :   if (!terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,viewer));
     148           1 :   PetscCall(PetscViewerPopFormat(viewer));
     149             : 
     150           1 :   PetscCall(PEPDestroy(&pep));
     151           1 :   PetscCall(MatDestroy(&M));
     152           1 :   PetscCall(MatDestroy(&C));
     153           1 :   PetscCall(MatDestroy(&K));
     154           1 :   PetscCall(PetscFree(ctx));
     155           1 :   PetscCall(SlepcFinalize());
     156             :   return 0;
     157             : }
     158             : 
     159             : /*
     160             :     Function for user-defined stopping test.
     161             : 
     162             :     Ignores the value of nev. It only takes into account the number of
     163             :     eigenpairs that have converged in recent outer iterations (restarts);
     164             :     if no new eigenvalues have converged in the last few restarts,
     165             :     we stop the iteration, assuming that no more eigenvalues are present
     166             :     inside the region.
     167             : */
     168          13 : PetscErrorCode MyStoppingTest(PEP pep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,PEPConvergedReason *reason,void *ptr)
     169             : {
     170          13 :   CTX_SPRING     *ctx = (CTX_SPRING*)ptr;
     171             : 
     172          13 :   PetscFunctionBeginUser;
     173             :   /* check usual termination conditions, but ignoring the case nconv>=nev */
     174          13 :   PetscCall(PEPStoppingBasic(pep,its,max_it,nconv,PETSC_INT_MAX,reason,NULL));
     175          13 :   if (*reason==PEP_CONVERGED_ITERATING) {
     176             :     /* check if nconv is the same as before */
     177          13 :     if (nconv==ctx->lastnconv) ctx->nreps++;
     178             :     else {
     179           1 :       ctx->lastnconv = nconv;
     180           1 :       ctx->nreps     = 0;
     181             :     }
     182             :     /* check if no eigenvalues converged in last 10 restarts */
     183          13 :     if (nconv && ctx->nreps>10) *reason = PEP_CONVERGED_USER;
     184             :   }
     185          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     186             : }
     187             : 
     188             : /*TEST
     189             : 
     190             :    test:
     191             :       args: -terse
     192             :       requires: !single
     193             :       suffix: 1
     194             : 
     195             : TEST*/

Generated by: LCOV version 1.14