Actual source code: dvd_testconv.c

  1: /*
  2:   SLEPc eigensolver: "davidson"

  4:   Step: test for convergence

  6:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  7:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  8:    Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain

 10:    This file is part of SLEPc.

 12:    SLEPc is free software: you can redistribute it and/or modify it under  the
 13:    terms of version 3 of the GNU Lesser General Public License as published by
 14:    the Free Software Foundation.

 16:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 17:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 18:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 19:    more details.

 21:    You  should have received a copy of the GNU Lesser General  Public  License
 22:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 23:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 24: */

 26:  #include davidson.h

 28: PetscBool dvd_testconv_basic_0(dvdDashboard *d, PetscScalar eigvr,
 29:                                 PetscScalar eigvi, PetscReal r,
 30:                                 PetscReal *err);
 31: PetscBool dvd_testconv_slepc_0(dvdDashboard *d, PetscScalar eigvr,
 32:                                 PetscScalar eigvi, PetscReal r,
 33:                                 PetscReal *err);

 37: PetscErrorCode dvd_testconv_basic(dvdDashboard *d, dvdBlackboard *b)
 38: {
 39:   PetscErrorCode  ierr;

 42:   /* Setup the step */
 43:   if (b->state >= DVD_STATE_CONF) {
 44:     PetscFree(d->testConv_data);
 45:     d->testConv = dvd_testconv_basic_0;
 46:   }
 47:   return(0);
 48: }

 52: PetscBool dvd_testconv_basic_0(dvdDashboard *d, PetscScalar eigvr,
 53:                                 PetscScalar eigvi, PetscReal r,
 54:                                 PetscReal *err)
 55: {
 56:   PetscBool       conv;
 57:   PetscReal       eig_norm, errest;

 60:   eig_norm = SlepcAbsEigenvalue(eigvr, eigvi);
 61:   errest = r/eig_norm;
 62:   conv = (errest <= d->tol) ? PETSC_TRUE : PETSC_FALSE;
 63:   if (err) *err = errest;
 64:   PetscFunctionReturn(conv);
 65: }

 69: PetscErrorCode dvd_testconv_slepc(dvdDashboard *d, dvdBlackboard *b)
 70: {
 71:   PetscErrorCode  ierr;

 74:   /* Setup the step */
 75:   if (b->state >= DVD_STATE_CONF) {
 76:     PetscFree(d->testConv_data);
 77:     d->testConv = dvd_testconv_slepc_0;
 78:   }
 79:   return(0);
 80: }

 84: PetscBool dvd_testconv_slepc_0(dvdDashboard *d, PetscScalar eigvr,
 85:                                 PetscScalar eigvi, PetscReal r,
 86:                                 PetscReal *err)
 87: {
 88:   PetscErrorCode  ierr;

 91:   (*d->eps->converged)(d->eps, eigvr, eigvi, r, err, d->eps->convergedctx);
 92:   CHKERRABORT(PetscObjectComm((PetscObject)d->eps), ierr);
 93:   PetscFunctionReturn(*err<d->eps->tol ? PETSC_TRUE : PETSC_FALSE);
 94: }