Anasazi Version of the Day
Loading...
Searching...
No Matches
FortranRoutines.h
1// @HEADER
2// *****************************************************************************
3// Anasazi: Block Eigensolvers Package
4//
5// Copyright 2004 NTESS and the Anasazi contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10// This software is a result of the research described in the report
11//
12// "A comparison of algorithms for modal analysis in the absence
13// of a sparse direct method", P. Arbenz, R. Lehoucq, and U. Hetmaniuk,
14// Sandia National Laboratories, Technical report SAND2003-1028J.
15//
16// It is based on the Epetra, AztecOO, and ML packages defined in the Trilinos
17// framework ( http://trilinos.org/ ).
18
19#ifndef FORTRAN_ROUTINES
20#define FORTRAN_ROUTINES
21
22#include "Epetra_BLAS_wrappers.h"
23#include "Epetra_LAPACK_wrappers.h"
24
25#ifdef EPETRA_MPI
26#include <mpi.h>
27#endif
28
29// LOGICAL as 4 bytes
30typedef int LOGICAL;
31
32#ifdef __cplusplus
33extern "C" {
34#endif
35
36// Double precision BLAS 1 //
37void PREFIX F77_FUNC(dswap,DSWAP)(int *n, double x[], int* incx, double y[], int* incy);
38
39// Double precision LAPACK //
40void PREFIX F77_FUNC(dgeqrf,DGEQRF)(int *M, int *N, double *A, int *lda, double *tau,
41 double *work, int *lwork, int *info);
42void PREFIX F77_FUNC(dormqr,DORMQR)(Epetra_fcd, Epetra_fcd, int *M, int *N, int *K, double *A,
43 int *lda, double *tau, double *C, int *ldc, double *work, int *lwork,
44 int *info);
45void PREFIX F77_FUNC(dsteqr,DSTEQR)(Epetra_fcd, int *N, double *D, double *E, double *Z,
46 int *ldz, double *work, int *info);
47
48#if defined (INTEL_CXML)
49int PREFIX F77_FUNC(ilaenv,ILAENV)(int *ispec, char *NAME, unsigned int len_name, char *OPTS,
50 unsigned int len_opts, int *N1, int *N2, int *N3, int *N4);
51#else
52int PREFIX F77_FUNC(ilaenv,ILAENV)(int *ispec, char *NAME, char *OPTS, int *N1, int *N2,
53 int *N3, int *N4, int len_name, int len_opts);
54#endif
55
56// Double precision customized ARPACK routines //
57#if defined (INTEL_CXML)
58void PREFIX F77_FUNC(mydsaupd,MYDSAUPD)(int *, char *, unsigned int, int *, char *,
59 unsigned int, int *, double *, double *, int *, double *, int *, int *,
60 int *, double *, double *, int *, int *, int *);
61void PREFIX F77_FUNC(dseupd,DSEUPD)(LOGICAL *rvec, char *HOWMNY, unsigned int len_howny,
62 LOGICAL *select, double *D, double *Z, int *ldz, double *sigma, char *BMAT,
63 unsigned int len_bmat, int *N, char *which, unsigned int len_which,
64 int *nev, double *tol, double *resid, int *ncv, double *V, int *ldv,
65 int *iparam, int *ipntr, double *workd, double *workl, int *lworkl,
66 int *info);
67#else
68void PREFIX F77_FUNC(mydsaupd,MYDSAUPD)(int *, char *, int *, char *, int *,
69 double *, double *, int *, double *, int *, int *, int *, double *,
70 double *, int *, int *, int *, int, int);
71void PREFIX F77_FUNC(dseupd,DSEUPD)(LOGICAL *rvec, char *HOWMNY, LOGICAL *select, double *D,
72 double *Z, int *ldz, double *sigma, char *BMAT, int *N, char *which,
73 int *nev, double *tol, double *resid, int *ncv, double *V, int *ldv,
74 int *iparam, int *ipntr, double *workd, double *workl, int *lworkl,
75 int *info, int len_howmny, int len_bmat, int len_which);
76#endif
77
78#ifdef EPETRA_MPI
79
80// Double precision customized PARPACK routines //
81#if defined (INTEL_CXML)
82void PREFIX F77_FUNC(mypdsaupd,MYPDSAUPD)(MPI_Comm *, int *, char *, unsigned int, int *,
83 char *, unsigned int, int *, double *, double *, int *, double *,
84 int *, int *, int *, double *, double *, int *, int *, int *);
85void PREFIX F77_FUNC(pdseupd,PDSEUPD)(MPI_Comm *MyComm, LOGICAL *rvec, char *HOWMNY,
86 unsigned int len_howmny, LOGICAL *select, double *D, double *Z,
87 int *ldz, double *sigma, char *BMAT, unsigned int len_bmat, int *N,
88 char *which, unsigned int len_which, int *nev, double *tol, double *resid,
89 int *ncv, double *V, int *ldv, int *iparam, int *ipntr, double *workd,
90 double *workl, int *lworkl, int *info);
91#else
92void PREFIX F77_FUNC(mypdsaupd,MYPDSAUPD)(MPI_Comm *, int *, char *, int *, char *, int *,
93 double *, double *, int *, double *, int *, int *, int *, double *,
94 double *, int *, int *, int *, int, int);
95void PREFIX F77_FUNC(pdseupd,PDSEUPD)(MPI_Comm *MyComm, LOGICAL *rvec, char *HOWMNY,
96 LOGICAL *select, double *D, double *Z, int *ldz, double *sigma,
97 char *BMAT, int *N, char *which, int *nev, double *tol, double *resid,
98 int *ncv, double *V, int *ldv, int *iparam, int *ipntr, double *workd,
99 double *workl, int *lworkl, int *info, int len_howmny,
100 int len_bmat, int len_which);
101#endif
102
103#endif
104
105#ifdef __cplusplus
106}
107#endif
108
109class FortranRoutines {
110
111 public:
112
113 // Double precision BLAS 1 //
114 void SCAL_INCX(int N, double ALPHA, double *X, int incX) const;
115 void SWAP(int N, double *X, int incx, double *Y, int incy) const;
116
117 // Double precision LAPACK //
118 void GEQRF(int M, int N, double *A, int lda, double *tau, double *work, int lwork,
119 int *info) const;
120 void ORMQR(char SIDE, char TRANS, int M, int N, int K, double *A, int lda, double *tau,
121 double *C, int ldc, double *work, int lwork, int *info) const;
122 void SPEV(char JOBZ, char UPLO, int N, double *A, double *W, double *Z, int ldz,
123 double *work, int *info) const;
124 void STEQR(char COMPZ, int N, double *D, double *E, double *Z, int ldz, double *work,
125 int *info) const;
126 void SYEV(char JOBZ, char UPLO, int N, double *A, int lda, double *W, double *work,
127 int lwork, int *info) const;
128 void SYGV(int itype, char JOBZ, char UPLO, int N, double *A, int lda, double *B, int ldb,
129 double *W, double *work, int lwork, int *info) const;
130
131 int LAENV(int ispec, char *NAME, char *OPTS, int N1, int N2, int N3, int N4,
132 int len_name, int len_opts) const;
133
134 // Double precision ARPACK routines
135 void SAUPD(int *ido, char BMAT, int N, char *which, int nev, double tol, double *resid,
136 int ncv, double *V, int ldv, int *iparam, int *ipntr, double *workd, double *workl,
137 int lworkl, int *info, int verbose) const;
138 void SEUPD(LOGICAL rvec, char HOWMNY, LOGICAL *select, double *D, double *Z, int ldz,
139 double sigma, char BMAT, int N, char *which, int nev, double tol, double *resid,
140 int ncv, double *V, int ldv, int *iparam, int *ipntr, double *workd,
141 double *workl, int lworkl, int *info) const;
142
143#ifdef EPETRA_MPI
144 // Double precision PARPACK routines
145 void PSAUPD(MPI_Comm MyComm, int *ido, char BMAT, int N, char *which, int nev, double tol,
146 double *resid, int ncv, double *V, int ldv, int *iparam, int *ipntr, double *workd,
147 double *workl, int lworkl, int *info, int verbose) const;
148 void PSEUPD(MPI_Comm MyComm, LOGICAL rvec, char HOWMNY, LOGICAL *select, double *D, double *Z,
149 int ldz, double sigma, char BMAT, int N, char *which, int nev, double tol,
150 double *resid, int ncv, double *V, int ldv, int *iparam, int *ipntr, double *workd,
151 double *workl, int lworkl, int *info) const;
152
153#endif
154
155};
156
157#endif