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-2011, Universitat Politecnica de Valencia, Spain

 10:    This file is part of SLEPc.
 11:       
 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;


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

 49:   return(0);
 50: }

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


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

 68:   PetscFunctionReturn(conv);
 69: }

 73: PetscErrorCode dvd_testconv_slepc(dvdDashboard *d, dvdBlackboard *b)
 74: {
 75:   PetscErrorCode  ierr;


 79:   /* Setup the step */
 80:   if (b->state >= DVD_STATE_CONF) {
 81:     PetscFree(d->testConv_data);
 82:     d->testConv = dvd_testconv_slepc_0;
 83:   }

 85:   return(0);
 86: }

 90: PetscBool dvd_testconv_slepc_0(dvdDashboard *d, PetscScalar eigvr,
 91:                                 PetscScalar eigvi, PetscReal r,
 92:                                 PetscReal *err)
 93: {
 94:   PetscErrorCode  ierr;


 98:   (*d->eps->conv_func)(d->eps, eigvr, eigvi, r, err,
 99:                               d->eps->conv_ctx);
100:   CHKERRABORT(((PetscObject)d->eps)->comm, ierr);

102:   PetscFunctionReturn(*err<d->eps->tol ? PETSC_TRUE : PETSC_FALSE);
103: }