Actual source code: zkrylovschurf.c

slepc-3.22.2 2024-12-02
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: */

 11: #include <petsc/private/fortranimpl.h>
 12: #include <slepceps.h>

 14: #if defined(PETSC_HAVE_FORTRAN_CAPS)
 15: #define epskrylovschurgetsubintervals_    EPSKRYLOVSCHURGETSUBINTERVALS
 16: #define epskrylovschurgetinertias_        EPSKRYLOVSCHURGETINERTIAS
 17: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 18: #define epskrylovschurgetsubintervals_    epskrylovschurgetsubintervals
 19: #define epskrylovschurgetinertias_        epskrylovschurgetinertias
 20: #endif

 22: SLEPC_EXTERN void epskrylovschurgetsubintervals_(EPS *eps,PetscReal *subint,PetscErrorCode *ierr)
 23: {
 24:   PetscReal *osubint;
 25:   PetscInt  npart;

 27:   CHKFORTRANNULLREAL(subint);
 28:   *ierr = EPSKrylovSchurGetSubintervals(*eps,&osubint); if (*ierr) return;
 29:   *ierr = EPSKrylovSchurGetPartitions(*eps,&npart); if (*ierr) return;
 30:   *ierr = PetscArraycpy(subint,osubint,npart+1); if (*ierr) return;
 31:   *ierr = PetscFree(osubint);
 32: }

 34: SLEPC_EXTERN void epskrylovschurgetinertias_(EPS *eps,PetscInt *nshift,PetscReal *shifts,PetscInt *inertias,PetscErrorCode *ierr)
 35: {
 36:   PetscReal *oshifts;
 37:   PetscInt  *oinertias;
 38:   PetscInt  n;

 40:   CHKFORTRANNULLREAL(shifts);
 41:   CHKFORTRANNULLINTEGER(inertias);
 42:   *ierr = EPSKrylovSchurGetInertias(*eps,&n,&oshifts,&oinertias); if (*ierr) return;
 43:   if (shifts) { *ierr = PetscArraycpy(shifts,oshifts,n); if (*ierr) return; }
 44:   if (inertias) { *ierr = PetscArraycpy(inertias,oinertias,n); if (*ierr) return; }
 45:   *nshift = n;
 46:   *ierr = PetscFree(oshifts);if (*ierr) return;
 47:   *ierr = PetscFree(oinertias);
 48: }