Subversion Repositories slepc-dev

Rev

Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1376 slepc 1
/*
2
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 3
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2575 eromero 4
   Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
1024 slepc 5
 
1672 slepc 6
   This file is part of SLEPc.
7
 
8
   SLEPc is free software: you can redistribute it and/or modify it under  the
9
   terms of version 3 of the GNU Lesser General Public License as published by
10
   the Free Software Foundation.
11
 
12
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
13
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
14
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
15
   more details.
16
 
17
   You  should have received a copy of the GNU Lesser General  Public  License
18
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
1376 slepc 19
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20
*/
21
 
2283 jroman 22
#include <private/fortranimpl.h>
23
#include <slepcst.h>
1024 slepc 24
 
2320 jroman 25
#if defined(PETSC_HAVE_FORTRAN_CAPS)
1781 antodo 26
#define stshellgetcontext_        STSHELLGETCONTEXT
1024 slepc 27
#define stshellsetapply_          STSHELLSETAPPLY
28
#define stshellsetapplytranspose_ STSHELLSETAPPLYTRANSPOSE
29
#define stshellsetbacktransform_  STSHELLSETBACKTRANSFORM
30
#elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
1781 antodo 31
#define stshellgetcontext_        stshellgetcontext
1024 slepc 32
#define stshellsetapply_          stshellsetapply
33
#define stshellsetapplytranspose_ stshellsetapplytranspose
34
#define stshellsetbacktransform_  stshellsetbacktransform
35
#endif
36
 
37
/* These are not extern C because they are passed into non-extern C user level functions */
1780 antodo 38
static PetscErrorCode ourshellapply(ST st,Vec x,Vec y)
1024 slepc 39
{
40
  PetscErrorCode ierr = 0;
1780 antodo 41
  (*(void (PETSC_STDCALL *)(ST*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)st)->fortran_func_pointers[0]))(&st,&x,&y,&ierr);CHKERRQ(ierr);
1024 slepc 42
  return 0;
43
}
44
 
1780 antodo 45
static PetscErrorCode ourshellapplytranspose(ST st,Vec x,Vec y)
1024 slepc 46
{
47
  PetscErrorCode ierr = 0;
1780 antodo 48
  (*(void (PETSC_STDCALL *)(ST*,Vec*,Vec*,PetscErrorCode*))(((PetscObject)st)->fortran_func_pointers[1]))(&st,&x,&y,&ierr);CHKERRQ(ierr);
1024 slepc 49
  return 0;
50
}
51
 
1780 antodo 52
static PetscErrorCode ourshellbacktransform(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
1024 slepc 53
{
54
  PetscErrorCode ierr = 0;
1780 antodo 55
  (*(void (PETSC_STDCALL *)(ST*,PetscInt*,PetscScalar*,PetscScalar*,PetscErrorCode*))(((PetscObject)st)->fortran_func_pointers[2]))(&st,&n,eigr,eigi,&ierr);CHKERRQ(ierr);
1024 slepc 56
  return 0;
57
}
58
 
59
EXTERN_C_BEGIN
60
 
1780 antodo 61
void PETSC_STDCALL stshellgetcontext_(ST *st,void **ctx,PetscErrorCode *ierr)
62
{
63
  *ierr = STShellGetContext(*st,ctx);
64
}
65
 
1024 slepc 66
void PETSC_STDCALL stshellsetapply_(ST *st,void (PETSC_STDCALL *apply)(void*,Vec *,Vec *,PetscErrorCode*),
67
                                    PetscErrorCode *ierr)
68
{
1780 antodo 69
  PetscObjectAllocateFortranPointers(*st,3);
70
  ((PetscObject)*st)->fortran_func_pointers[0] = (PetscVoidFunction)apply;
1024 slepc 71
  *ierr = STShellSetApply(*st,ourshellapply);
72
}
73
 
74
void PETSC_STDCALL stshellsetapplytranspose_(ST *st,void (PETSC_STDCALL *applytranspose)(void*,Vec *,Vec *,PetscErrorCode*),
75
                                             PetscErrorCode *ierr)
76
{
1780 antodo 77
  PetscObjectAllocateFortranPointers(*st,3);
78
  ((PetscObject)*st)->fortran_func_pointers[1] = (PetscVoidFunction)applytranspose;
1024 slepc 79
  *ierr = STShellSetApplyTranspose(*st,ourshellapplytranspose);
80
}
81
 
82
void PETSC_STDCALL stshellsetbacktransform_(ST *st,void (PETSC_STDCALL *backtransform)(void*,PetscScalar*,PetscScalar*,PetscErrorCode*),
83
                                    PetscErrorCode *ierr)
84
{
1781 antodo 85
  PetscObjectAllocateFortranPointers(*st,3);
1780 antodo 86
  ((PetscObject)*st)->fortran_func_pointers[2] = (PetscVoidFunction)backtransform;
1024 slepc 87
  *ierr = STShellSetBackTransform(*st,ourshellbacktransform);
88
}
89
 
90
EXTERN_C_END
91