aboutsummaryrefslogtreecommitdiff
path: root/redist
diff options
context:
space:
mode:
authorcnlohr <lohr85@gmail.com>2018-03-17 01:13:28 -0400
committercnlohr <lohr85@gmail.com>2018-03-17 01:13:28 -0400
commit8b9a196232661bd6506a41b7ceca015d684086a8 (patch)
treea37e294ebc9c7944bea3c8e0e1f0bfc67fbc1c6d /redist
parent75a176d2ef65de8450b2054065d5c068a19966a6 (diff)
downloadlibsurvive-8b9a196232661bd6506a41b7ceca015d684086a8.tar.gz
libsurvive-8b9a196232661bd6506a41b7ceca015d684086a8.tar.bz2
Update variable ordering...
Diffstat (limited to 'redist')
-rw-r--r--redist/dclapack.h48
-rw-r--r--redist/dclhelpers.c36
-rw-r--r--redist/dclhelpers.h26
3 files changed, 55 insertions, 55 deletions
diff --git a/redist/dclapack.h b/redist/dclapack.h
index 62f3e9d..5aba214 100644
--- a/redist/dclapack.h
+++ b/redist/dclapack.h
@@ -44,15 +44,15 @@
}
/*
- * B = Transpose(A)
+ * R = Transpose(A)
* A is (n by m)
- * B is (m by n)
+ * R is (m by n)
*/
-#define TRANSP(A,B,n,m) { \
+#define TRANSP(R,A,n,m) { \
int i,j; \
for (i=0; i<n; i++) { \
for (j=0; j<m; j++) { \
- _(B,j,i) = _(A,i,j); \
+ _(R,j,i) = _(A,i,j); \
} \
} \
}
@@ -60,7 +60,7 @@
/*
* Calculate L,U of a matrix A with pivot table
*/
-#define LU(A,L,U,Piv,n) { \
+#define LU(L,U,A,Piv,n) { \
int i,j,k,_tempi; float _tempf; \
for (i=0; i<n; i++) { Piv[i]=i; } \
for (i=0; i<n; i++) { \
@@ -96,14 +96,14 @@
/*
* Pivots a matrix to a different matrix
- * B = Pivot(A) given table 'Piv'
- * A and B are (n by m)
+ * R = Pivot(A) given table 'Piv'
+ * A and R are (n by m)
*/
-#define PIVOT(A,B,Piv,n,m) { \
+#define PIVOT(R,A,Piv,n,m) { \
int i,j; \
for (j=0; j<n; j++) { \
for (i=0; i<m; i++) { \
- _(B,j,i) = _(A,Piv[j],i); \
+ _(R,j,i) = _(A,Piv[j],i); \
} \
} \
}
@@ -113,7 +113,7 @@
* L is n by n (lower triangular)
* B is n by m
*/
-#define L_SUB(L,X,B,n,m) { \
+#define L_SUB(X,L,B,n,m) { \
int i,j,k; \
for (i=0; i<m; i++) { \
for (j=0; j<n; j++) { \
@@ -129,7 +129,7 @@
* U is n by n (upper triangular)
* B is n by m
*/
-#define U_SUB(U,X,B,n,m) { \
+#define U_SUB(X,U,B,n,m) { \
int i,j,k; \
for (i=0; i<m; i++) { \
for (j=n-1; j>=0; j--) { \
@@ -160,14 +160,14 @@
FLOAT C[ORDER][ORDER];
#endif
-#define INV(A,Ainv,n,ORDER) { \
+#define INV(Ainv,A,n,ORDER) { \
INV_SETUP(ORDER) \
int Piv[ORDER]; \
IDENTITY(I,n); \
- LU(A,L,U,Piv,n); \
- PIVOT(I,Ipiv,Piv,n,n); \
- L_SUB(L,C,Ipiv,n,n); \
- U_SUB(U,Ainv,C,n,n); \
+ LU(L,U,A,Piv,n); \
+ PIVOT(Ipiv,I,Piv,n,n); \
+ L_SUB(C,L,Ipiv,n,n); \
+ U_SUB(Ainv,U,C,n,n); \
}
/*
@@ -186,13 +186,13 @@ PRINT(Ainv,n,n); \
* B (m by p)
* C (n by p)
*/
-#define MUL(A,B,C,n,m,p) { \
+#define MUL(R,A,B,n,m,p) { \
int i,j,k; \
for (i=0; i<n; i++) { \
for (j=0; j<p; j++) { \
- _(C,i,j) = 0.0f; \
+ _(R,i,j) = 0.0f; \
for (k=0; k<m; k++) { \
- _(C,i,j) += _(A,i,k) * _(B,k,j); \
+ _(R,i,j) += _(A,i,k) * _(B,k,j); \
} \
} \
} \
@@ -205,13 +205,13 @@ PRINT(Ainv,n,n); \
* C (n by p)
* D (n by p)
*/
-#define MULADD(A,B,C,D,n,m,p) { \
+#define MULADD(R,A,B,C,n,m,p) { \
int i,j,k; \
for (i=0; i<n; i++) { \
for (j=0; j<p; j++) { \
- _(D,i,j) = _(C,i,j); \
+ _(R,i,j) = _(C,i,j); \
for (k=0; k<m; k++) { \
- _(D,i,j) += _(A,i,k) * _(B,k,j); \
+ _(R,i,j) += _(A,i,k) * _(B,k,j); \
} \
} \
} \
@@ -224,7 +224,7 @@ PRINT(Ainv,n,n); \
* C (n by p)
* D (n by p)
*/
-#define GMULADD(A,B,C,D,alpha,beta,n,m,p) { \
+#define GMULADD(R,A,B,C,alpha,beta,n,m,p) { \
int i,j,k; \
float sum; \
for (i=0; i<n; i++) { \
@@ -233,7 +233,7 @@ PRINT(Ainv,n,n); \
for (k=0; k<m; k++) { \
sum += _(A,i,k) * _(B,k,j); \
} \
- _(D,i,j) = alpha * sum + beta * _(C,i,j); \
+ _(R,i,j) = alpha * sum + beta * _(C,i,j); \
} \
} \
}
diff --git a/redist/dclhelpers.c b/redist/dclhelpers.c
index 9d3d0eb..2ee6c4a 100644
--- a/redist/dclhelpers.c
+++ b/redist/dclhelpers.c
@@ -15,49 +15,49 @@ void dclIdentity( DCL_FLOAT * A, int n )
IDENTITY( A, n );
}
-void dclTransp( const DCL_FLOAT * A, DCL_FLOAT * B, int n, int m)
+void dclTransp( DCL_FLOAT * R, const DCL_FLOAT * A, int n, int m )
{
- TRANSP(A,B,n,m);
+ TRANSP(R,A,n,m);
}
-void dclLU( const DCL_FLOAT * A, DCL_FLOAT * L, DCL_FLOAT * U, int * Piv, int n )
+void dclLU( DCL_FLOAT * L, DCL_FLOAT * U, const DCL_FLOAT * A, int * Piv, int n )
{
- LU(A,L,U,Piv,n);
+ LU(L,U,A,Piv,n);
}
-void dclPivot( const DCL_FLOAT * A, DCL_FLOAT * B, int * Piv, int n, int m )
+void dclPivot( DCL_FLOAT * R, const DCL_FLOAT * A, int * Piv, int n, int m )
{
- PIVOT(A,B,Piv,n,m);
+ PIVOT(R,A,Piv,n,m);
}
-void dclLSub( const DCL_FLOAT * L, DCL_FLOAT * X, const DCL_FLOAT * B, int n, int m )
+void dclLSub( DCL_FLOAT * X, const DCL_FLOAT * L, const DCL_FLOAT * B, int n, int m )
{
- L_SUB(L,X,B,n,m);
+ L_SUB(X,L,B,n,m);
}
-void dclUSub( const DCL_FLOAT * U, DCL_FLOAT * X, const DCL_FLOAT * B, int n, int m )
+void dclUSub( DCL_FLOAT * X, const DCL_FLOAT * U, const DCL_FLOAT * B, int n, int m )
{
- U_SUB(U,X,B,n,m);
+ U_SUB(X,U,B,n,m);
}
-void dclInv( const DCL_FLOAT * A, DCL_FLOAT * Ainv, int n )
+void dclInv( DCL_FLOAT * Ainv, const DCL_FLOAT * A, int n )
{
- INV(A,Ainv,n,n);
+ INV(Ainv,A,n,n);
}
-void dclMul( const DCL_FLOAT * A, const DCL_FLOAT * B, DCL_FLOAT * C, int n, int m, int p )
+void dclMul( DCL_FLOAT * R, const DCL_FLOAT * A, const DCL_FLOAT * B, int n, int m, int p )
{
- MUL(A,B,C,n,m,p);
+ MUL(R,A,B,n,m,p);
}
-void dclMulAdd( const DCL_FLOAT * A, const DCL_FLOAT * B, const DCL_FLOAT * C, DCL_FLOAT * D, int n, int m, int p )
+void dclMulAdd( DCL_FLOAT * R, const DCL_FLOAT * A, const DCL_FLOAT * B, const DCL_FLOAT * C, int n, int m, int p )
{
- MULADD(A,B,C,D,n,m,p);
+ MULADD(R,A,B,C,n,m,p);
}
-void dclGMulAdd( const DCL_FLOAT * A, const DCL_FLOAT * B, const DCL_FLOAT * C, DCL_FLOAT * D, DCL_FLOAT alpha, DCL_FLOAT beta, int n, int m, int p )
+void dclGMulAdd( DCL_FLOAT * R, const DCL_FLOAT * A, const DCL_FLOAT * B, const DCL_FLOAT * C, DCL_FLOAT alpha, DCL_FLOAT beta, int n, int m, int p )
{
- GMULADD(A,B,C,D,alpha,beta,n,m,p);
+ GMULADD(R,A,B,C,alpha,beta,n,m,p);
}
diff --git a/redist/dclhelpers.h b/redist/dclhelpers.h
index fbaa7ec..4ada24a 100644
--- a/redist/dclhelpers.h
+++ b/redist/dclhelpers.h
@@ -11,51 +11,51 @@ void dclPrint( const DCL_FLOAT * A, int n, int m );
/* Returns the identity matrix */
void dclIdentity( DCL_FLOAT * A, int n );
-/* B = Transpose(A)
+/* R = Transpose(A)
A is (n by m)
- B is (m by n) */
-void dclTransp( const DCL_FLOAT * A, DCL_FLOAT * B, int n, int m );
+ R is (m by n) */
+void dclTransp( DCL_FLOAT * R, const DCL_FLOAT * A, int n, int m );
/* Calculate L,U of a matrix A with pivot table; the pivot table is output. */
-void dclLU( const DCL_FLOAT * A, DCL_FLOAT * L, DCL_FLOAT * U, int * Piv, int n );
+void dclLU( DCL_FLOAT * L, DCL_FLOAT * U, const DCL_FLOAT * A, int * Piv, int n );
/* Pivots a matrix to a different matrix
- B = Pivot(A) given table 'Piv'
- A and B are (n by m) */
-void dclPivot( const DCL_FLOAT * A, DCL_FLOAT * B, int * Piv, int n, int m );
+ R = Pivot(A) given table 'Piv'
+ A and R are (n by m) */
+void dclPivot( DCL_FLOAT * R, const DCL_FLOAT * A, int * Piv, int n, int m );
/* Solve LX=B for matrix X and B
L is n by n (lower triangular)
B is n by m */
-void dclLSub( const DCL_FLOAT * L, DCL_FLOAT * X, const DCL_FLOAT * B, int n, int m );
+void dclLSub( DCL_FLOAT * X, const DCL_FLOAT * L, const DCL_FLOAT * B, int n, int m );
/* Solve UX=B for matrix X and B
U is n by n (upper triangular)
B is n by m */
-void dclUSub( const DCL_FLOAT * U, DCL_FLOAT * X, const DCL_FLOAT * B, int n, int m );
+void dclUSub( DCL_FLOAT * X, const DCL_FLOAT * U, const DCL_FLOAT * B, int n, int m );
/* Inverts a matrix X (n by n) using the method of LU decomposition */
-void dclInv( const DCL_FLOAT * A, DCL_FLOAT * Ainv, int n );
+void dclInv( DCL_FLOAT * Ainv, const DCL_FLOAT * A, int n );
/* Matrix Multiply C = A * B
A (n by m)
B (m by p)
C (n by p) */
-void dclMul( const DCL_FLOAT * A, const DCL_FLOAT * B, DCL_FLOAT * C, int n, int m, int p );
+void dclMul( DCL_FLOAT * R, const DCL_FLOAT * A, const DCL_FLOAT * B, int n, int m, int p );
/* Matrix Multiply D = A * B + C
A (n by m)
B (m by p)
C (n by p)
D (n by p) */
-void dclMulAdd( const DCL_FLOAT * A, const DCL_FLOAT * B, const DCL_FLOAT * C, DCL_FLOAT * D, int n, int m, int p );
+void dclMulAdd( DCL_FLOAT * R, const DCL_FLOAT * A, const DCL_FLOAT * B, const DCL_FLOAT * C, int n, int m, int p );
/* Matrix Multiply D = alpha * A * B + beta * C
A (n by m)
B (m by p)
C (n by p)
D (n by p) */
-void dclGMulAdd( const DCL_FLOAT * A, const DCL_FLOAT * B, const DCL_FLOAT * C, DCL_FLOAT * D, DCL_FLOAT alpha, DCL_FLOAT beta, int n, int m, int p );
+void dclGMulAdd( DCL_FLOAT * D, const DCL_FLOAT * A, const DCL_FLOAT * B, const DCL_FLOAT * C, DCL_FLOAT alpha, DCL_FLOAT beta, int n, int m, int p );
#endif