Subversion Repositories slepc-dev

Rev

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

Rev Author Line No. Line
1376 slepc 1
/*
2
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3
      SLEPc - Scalable Library for Eigenvalue Problem Computations
4
      Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
1333 slepc 5
 
1376 slepc 6
      This file is part of SLEPc. See the README file for conditions of use
7
      and additional information.
8
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9
*/
10
 
1333 slepc 11
#include "zpetsc.h"
12
#include "slepcsvd.h"
13
 
14
#ifdef PETSC_HAVE_FORTRAN_CAPS
15
#define svdmonitordefault_           SVDMONITORDEFAULT
16
#define svdmonitorlg_                SVDMONITORLG
17
#define svdview_                     SVDVIEW
18
#define svdcreate_                   SVDCREATE
19
#define svdsettype_                  SVDSETTYPE
20
#define svdgettype_                  SVDGETTYPE
1345 slepc 21
#define svdgetip_                    SVDGETIP
1333 slepc 22
#define svdmonitorset_               SVDMONITORSET
23
#define svdgettransposemode_         SVDGETTRANSPOSEMODE
24
#define svdgetwhichsingulartriplets_ SVDGETWHICHSINGULARTRIPLETS
25
#define svdsetoptionsprefix_         SVDSETOPTIONSPREFIX
26
#define svdappendoptionsprefix_      SVDAPPENDOPTIONSPREFIX
27
#define svdgetoptionsprefix_         SVDGETOPTIONSPREFIX
28
#define svdgetconvergedreason_       SVDGETCONVERGEDREASON
29
#elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
30
#define svdmonitordefault_           svdmonitordefault
31
#define svdmonitorlg_                svdmonitorlg
32
#define svdview_                     svdview
33
#define svdcreate_                   svdcreate
34
#define svdsettype_                  svdsettype
35
#define svdgettype_                  svdgettype
1345 slepc 36
#define svdgetip_                    svdgetip
1333 slepc 37
#define svdmonitorset_               svdmonitorset
38
#define svdgettransposemode_         svdgettransposemode
39
#define svdgetwhichsingulartriplets_ svdgetwhichsingulartriplets
40
#define svdsetoptionsprefix_         svdsetoptionsprefix
41
#define svdappendoptionsprefix_      svdappendoptionsprefix
42
#define svdgetoptionsprefix_         svdgetoptionsprefix
43
#define svdgetconvergedreason_       svdgetconvergedreason
44
#endif
45
 
46
EXTERN_C_BEGIN
47
static void (PETSC_STDCALL *f1)(SVD*,int*,int*,PetscReal*,PetscReal*,int*,void*,PetscErrorCode*);
48
static void (PETSC_STDCALL *f2)(void*,PetscErrorCode*);
49
 
50
/*
51
   These are not usually called from Fortran but allow Fortran users
52
   to transparently set these monitors from .F code, hence no STDCALL
53
*/
54
void svdmonitordefault_(SVD *svd,int *it,int *nconv,PetscReal *sigma,PetscReal *errest,int *nest,void *ctx,PetscErrorCode *ierr)
55
{
56
  *ierr = SVDMonitorDefault(*svd,*it,*nconv,sigma,errest,*nest,ctx);
57
}
58
 
59
void svdmonitorlg_(SVD *svd,int *it,int *nconv,PetscReal *sigma,PetscReal *errest,int *nest,void *ctx,PetscErrorCode *ierr)
60
{
61
  *ierr = SVDMonitorLG(*svd,*it,*nconv,sigma,errest,*nest,ctx);
62
}
63
EXTERN_C_END
64
 
65
/* These are not extern C because they are passed into non-extern C user level functions */
66
static PetscErrorCode ourmonitor(SVD svd,int i,int nc,PetscReal *sigma,PetscReal *d,int l,void* ctx)
67
{
68
  PetscErrorCode ierr = 0;
69
  (*f1)(&svd,&i,&nc,sigma,d,&l,ctx,&ierr);CHKERRQ(ierr);
70
  return 0;
71
}
72
 
73
static PetscErrorCode ourdestroy(void* ctx)
74
{
75
  PetscErrorCode ierr = 0;
76
  (*f2)(ctx,&ierr);CHKERRQ(ierr);
77
  return 0;
78
}
79
 
80
EXTERN_C_BEGIN
81
 
82
void PETSC_STDCALL svdview_(SVD *svd,PetscViewer *viewer, PetscErrorCode *ierr)
83
{
84
  PetscViewer v;
85
  PetscPatchDefaultViewers_Fortran(viewer,v);
86
  *ierr = SVDView(*svd,v);
87
}
88
 
89
void PETSC_STDCALL svdcreate_(MPI_Fint *comm,SVD *svd,PetscErrorCode *ierr)
90
{
91
  *ierr = SVDCreate(MPI_Comm_f2c(*(comm)),svd);
92
}
93
 
94
void PETSC_STDCALL svdsettype_(SVD *svd,CHAR type PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
95
{
96
  char *t;
97
 
98
  FIXCHAR(type,len,t);
99
  *ierr = SVDSetType(*svd,t);
100
  FREECHAR(type,t);
101
}
102
 
103
void PETSC_STDCALL svdgettype_(SVD *svd,CHAR name PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
104
{
105
  const char *tname;
106
 
107
  *ierr = SVDGetType(*svd,&tname);if (*ierr) return;
108
#if defined(PETSC_USES_CPTOFCD)
109
  {
110
    char *t = _fcdtocp(name); int len1 = _fcdlen(name);
111
    *ierr = PetscStrncpy(t,tname,len1);
112
  }
113
#else
114
  *ierr = PetscStrncpy(name,tname,len);
115
#endif
116
  FIXRETURNCHAR(name,len);
117
}
118
 
1345 slepc 119
void PETSC_STDCALL svdgetip_(SVD *svd,IP *ip,int *ierr)
120
{
121
  *ierr = SVDGetIP(*svd,ip);
122
}
123
 
1333 slepc 124
void PETSC_STDCALL svdmonitorset_(SVD *svd,void (PETSC_STDCALL *monitor)(SVD*,int*,int*,PetscReal*,PetscReal*,int*,void*,PetscErrorCode*),
125
                                  void *mctx,void (PETSC_STDCALL *monitordestroy)(void *,PetscErrorCode *),PetscErrorCode *ierr)
126
{
127
  if ((void(*)())monitor == (void(*)())svdmonitordefault_) {
128
    *ierr = SVDMonitorSet(*svd,SVDMonitorDefault,0,0);
129
  } else if ((void(*)())monitor == (void(*)())svdmonitorlg_) {
130
    *ierr = SVDMonitorSet(*svd,SVDMonitorLG,0,0);
131
  } else {
132
    f1  = monitor;
133
    if (FORTRANNULLFUNCTION(monitordestroy)) {
134
      *ierr = SVDMonitorSet(*svd,ourmonitor,mctx,0);
135
    } else {
136
      f2 = monitordestroy;
137
      *ierr = SVDMonitorSet(*svd,ourmonitor,mctx,ourdestroy);
138
    }
139
  }
140
}
141
 
142
void PETSC_STDCALL svdgettransposemode_(SVD *svd,SVDTransposeMode *mode, int *ierr)
143
{
144
  *ierr = SVDGetTransposeMode(*svd,mode);
145
}
146
 
147
void PETSC_STDCALL svdgetwhichsingulartriplets_(SVD *svd,SVDWhich *which, int *ierr)
148
{
149
  *ierr = SVDGetWhichSingularTriplets(*svd,which);
150
}
151
 
152
void PETSC_STDCALL svdsetoptionsprefix_(SVD *svd,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
153
{
154
  char *t;
155
 
156
  FIXCHAR(prefix,len,t);
157
  *ierr = SVDSetOptionsPrefix(*svd,t);
158
  FREECHAR(prefix,t);
159
}
160
 
161
void PETSC_STDCALL svdappendoptionsprefix_(SVD *svd,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
162
{
163
  char *t;
164
 
165
  FIXCHAR(prefix,len,t);
166
  *ierr = SVDAppendOptionsPrefix(*svd,t);
167
  FREECHAR(prefix,t);
168
}
169
 
170
void PETSC_STDCALL svdgetoptionsprefix_(SVD *svd,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
171
{
172
  const char *tname;
173
 
174
  *ierr = SVDGetOptionsPrefix(*svd,&tname);
175
#if defined(PETSC_USES_CPTOFCD)
176
  {
177
    char *t = _fcdtocp(prefix); int len1 = _fcdlen(prefix);
178
    *ierr = PetscStrncpy(t,tname,len1); if (*ierr) return;
179
  }
180
#else
181
  *ierr = PetscStrncpy(prefix,tname,len); if (*ierr) return;
182
#endif
183
  FIXRETURNCHAR(prefix,len);
184
}
185
 
186
void PETSC_STDCALL svdgetconvergedreason_(SVD *svd,SVDConvergedReason *reason,PetscErrorCode *ierr)
187
{
188
  *ierr = SVDGetConvergedReason(*svd,reason);
189
}
190
 
191
EXTERN_C_END