about summary refs log tree commit diff
path: root/src/mvlmm.cpp
diff options
context:
space:
mode:
authorPjotr Prins2018-09-28 15:22:11 +0000
committerPjotr Prins2018-09-28 15:22:11 +0000
commitc21dfe50df1c1399190a7c5c100c5d0b14ab7ef8 (patch)
treec81d820dd504f9956afb3f1778289e1ee9c21918 /src/mvlmm.cpp
parent8597ff255020980a7e79603331b5d1dd6a7a42b9 (diff)
downloadpangemma-c21dfe50df1c1399190a7c5c100c5d0b14ab7ef8.tar.gz
Fixes regression #179
Diffstat (limited to 'src/mvlmm.cpp')
-rw-r--r--src/mvlmm.cpp10
1 files changed, 9 insertions, 1 deletions
diff --git a/src/mvlmm.cpp b/src/mvlmm.cpp
index a3c7c85..4b8db05 100644
--- a/src/mvlmm.cpp
+++ b/src/mvlmm.cpp
@@ -256,7 +256,15 @@ double EigenProc(const gsl_matrix *V_g, const gsl_matrix *V_e, gsl_vector *D_l,
   gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, V_e_hi, VgVehi, 0.0, Lambda);
 
   // Eigen decomposition of Lambda.
-  EigenDecomp_Zeroed(Lambda, U_l, D_l, 0);
+  EigenDecomp(Lambda, U_l, D_l, 0);
+
+  for (size_t i = 0; i < d_size; i++) {
+    d = gsl_vector_get(D_l, i);
+    if (d < 0) {
+      gsl_vector_set(D_l, i, 0);
+    }
+  }
+
 
   // Calculate UltVeh and UltVehi.
   gsl_blas_dgemm(CblasTrans, CblasNoTrans, 1.0, U_l, V_e_h, 0.0, UltVeh);