diff --git a/jacobi/openmp/not_opt/Makefile b/jacobi/openmp/not_opt/Makefile
index fc4833da0e6826ad3d36680ad49de5393c24039e..9ec273e11a1b937239434170c612555136db02f6 100644
--- a/jacobi/openmp/not_opt/Makefile
+++ b/jacobi/openmp/not_opt/Makefile
@@ -51,7 +51,7 @@ valgrind_callgrind: $(PROG_CALLGRIND)
 
 valgrind_cachegrind: $(PROG_CACHEGRIND)
 	@echo 'oooOOO... valgrind_cachegrind ...OOOooo'
-	valgrind --tool=cachegrind --log-file=valgrind_cachegrind_log_.%p.txt ./$< 128 128 2
+	valgrind --tool=cachegrind --cache-sim=yes --log-file=valgrind_cachegrind_log_.%p.txt ./$< 128 128 2
 	@echo '$$ cg_annotate --auto=yes cachegrind.out.<pid> | less'
 	@echo '(kcachegrind is required in order to visualize the output using the GUI)'
 	@echo 'oooOOO... valgrind_cachegrind ...OOOooo'
diff --git a/jacobi/openmp/not_opt/src/jacobi_2D_omp_not_opt.c b/jacobi/openmp/not_opt/src/jacobi_2D_omp_not_opt.c
index d61d7fd456d76679efd65d95f616d019038c2b07..5768ff47d1139db9dd6f5f6fde3d5fff2a090ba9 100644
--- a/jacobi/openmp/not_opt/src/jacobi_2D_omp_not_opt.c
+++ b/jacobi/openmp/not_opt/src/jacobi_2D_omp_not_opt.c
@@ -234,9 +234,6 @@ void JacobiAlgorithm(MyData      **const restrict Phi,
     #pragma omp barrier
 #endif /* DEBUG */
 
-
-    /***** MISSING LOOP(s) PARALLELIZATION *****/
-    
     for (int j=jbeg ; j<=jend ; j++)
       {
 	for (int i=ibeg ; i<=iend ; i++)
@@ -244,10 +241,7 @@ void JacobiAlgorithm(MyData      **const restrict Phi,
 	    Phi[j][i] = 0.25 * (Phi0[j][i-1] + Phi0[j][i+1] +
 				Phi0[j-1][i] + Phi0[j+1][i]);
 
-	    #pragma omp critical
-	    {
-	      *error += delta[X] * delta[Y] * fabs(Phi[j][i] - Phi0[j][i]);
-	    }
+	    *error += delta[X] * delta[Y] * fabs(Phi[j][i] - Phi0[j][i]);
 	  } /* loop over columns */
       } /* loop over rows */
   } /* omp parallel */
diff --git a/jacobi/serial/not_opt/Makefile b/jacobi/serial/not_opt/Makefile
index 191f891f79c02871d09e8960b75ea9a71d369c47..fef753e16ce11f9fffc7082ac3b1dfafd85ec1c7 100644
--- a/jacobi/serial/not_opt/Makefile
+++ b/jacobi/serial/not_opt/Makefile
@@ -48,7 +48,7 @@ valgrind_callgrind: $(PROG_CALLGRIND)
 
 valgrind_cachegrind: $(PROG_CACHEGRIND)
 	@echo 'oooOOO... valgrind_cachegrind ...OOOooo'
-	valgrind --tool=cachegrind --log-file=valgrind_cachegrind_log_.%p.txt ./$< 128 128
+	valgrind --tool=cachegrind --cache-sim=yes --log-file=valgrind_cachegrind_log_.%p.txt ./$< 128 128
 	@echo '$$ cg_annotate --auto=yes cachegrind.out.<pid> | less'
 	@echo '(kcachegrind is required in order to visualize the output using the GUI)'
 	@echo 'oooOOO... valgrind_cachegrind ...OOOooo'
diff --git a/jacobi_solutions/serial/not_opt/Makefile b/jacobi_solutions/serial/not_opt/Makefile
index 191f891f79c02871d09e8960b75ea9a71d369c47..fef753e16ce11f9fffc7082ac3b1dfafd85ec1c7 100644
--- a/jacobi_solutions/serial/not_opt/Makefile
+++ b/jacobi_solutions/serial/not_opt/Makefile
@@ -48,7 +48,7 @@ valgrind_callgrind: $(PROG_CALLGRIND)
 
 valgrind_cachegrind: $(PROG_CACHEGRIND)
 	@echo 'oooOOO... valgrind_cachegrind ...OOOooo'
-	valgrind --tool=cachegrind --log-file=valgrind_cachegrind_log_.%p.txt ./$< 128 128
+	valgrind --tool=cachegrind --cache-sim=yes --log-file=valgrind_cachegrind_log_.%p.txt ./$< 128 128
 	@echo '$$ cg_annotate --auto=yes cachegrind.out.<pid> | less'
 	@echo '(kcachegrind is required in order to visualize the output using the GUI)'
 	@echo 'oooOOO... valgrind_cachegrind ...OOOooo'
diff --git a/jacobi_solutions/serial/not_opt/src/jacobi_2D_serial_not_opt.c b/jacobi_solutions/serial/not_opt/src/jacobi_2D_serial_not_opt.c
index ba53f321242432290dca329fae94271189fadf82..f17e977a65b0c6846553a31de4153ad4c956aa24 100644
--- a/jacobi_solutions/serial/not_opt/src/jacobi_2D_serial_not_opt.c
+++ b/jacobi_solutions/serial/not_opt/src/jacobi_2D_serial_not_opt.c
@@ -214,10 +214,13 @@ void JacobiAlgorithm(MyData      **const restrict Phi,
     {
       for (int i=ibeg ; i<=iend ; i++)
 	{
-	  Phi[j][i] = 0.25 * (Phi0[j][i-1] + Phi0[j][i+1] +
-			      Phi0[j-1][i] + Phi0[j+1][i]);
+	  Phi[j][i] = 0.25 * (Phi0[j][i-1] +
+			      Phi0[j][i+1] +
+			      Phi0[j-1][i] +
+			      Phi0[j+1][i]);
                 
-	  *error += delta[X] * delta[Y] * fabs(Phi[j][i] - Phi0[j][i]);
+	  *error += (delta[X] * delta[Y]) *
+	            fabs(Phi[j][i] - Phi0[j][i]);
 	} /* loop over columns */
     } /* loop over rows */
   
diff --git a/jacobi_solutions/serial/not_opt/src/tools.c b/jacobi_solutions/serial/not_opt/src/tools.c
index be2761d170c23fde4d5d4aeda54af00e56abf3ce..5e1f75a00f79bb5f8d32a94334a4f36e21730be9 100644
--- a/jacobi_solutions/serial/not_opt/src/tools.c
+++ b/jacobi_solutions/serial/not_opt/src/tools.c
@@ -51,9 +51,8 @@ void Show_2DdblArray(const MyData **const A,
 
 double seconds()
 {
-  struct timeval tmp;
-  gettimeofday(&tmp, (struct timezone *)0);
-  double sec = tmp.tv_sec + ((double)tmp.tv_usec)/1000000.0;
-
-  return sec;
+  struct timespec ts;
+  return (clock_gettime( CLOCK_PROCESS_CPUTIME_ID, &ts ),
+	   (double)ts.tv_sec +
+	   (double)ts.tv_nsec * 1e-9);
 }
diff --git a/jacobi_solutions/serial/opt/Makefile b/jacobi_solutions/serial/opt/Makefile
index d278dc33018178b113e73ed3f93533d6e0a77861..2a7fd8b407d50f92f26af8f404c312bda712a1fd 100644
--- a/jacobi_solutions/serial/opt/Makefile
+++ b/jacobi_solutions/serial/opt/Makefile
@@ -48,7 +48,7 @@ valgrind_callgrind: $(PROG_CALLGRIND)
 
 valgrind_cachegrind: $(PROG_CACHEGRIND)
 	@echo 'oooOOO... valgrind_cachegrind ...OOOooo'
-	valgrind --tool=cachegrind --log-file=valgrind_cachegrind_log_.%p.txt ./$< 128 128
+	valgrind --tool=cachegrind --cache-sim=yes --log-file=valgrind_cachegrind_log_.%p.txt ./$< 128 128
 	@echo '$$ cg_annotate --auto=yes cachegrind.out.<pid> | less'
 	@echo '(kcachegrind is required in order to visualize the output using the GUI)'
 	@echo 'oooOOO... valgrind_cachegrind ...OOOooo'
diff --git a/jacobi_solutions/serial/opt/src/jacobi_2D_serial_opt.c b/jacobi_solutions/serial/opt/src/jacobi_2D_serial_opt.c
index a14c810679a0d4087e1df3159ec7f1a8d98f238f..7bacb89eac762df1ef6f8d19792c570dd2be746b 100644
--- a/jacobi_solutions/serial/opt/src/jacobi_2D_serial_opt.c
+++ b/jacobi_solutions/serial/opt/src/jacobi_2D_serial_opt.c
@@ -212,13 +212,19 @@ void JacobiAlgorithm(MyData      **const restrict Phi,
     {
       for (int i=ibeg ; i<=iend ; i++)
 	{
-	  Phi[j][i] = 0.25 * (Phi0[j][i-1] + Phi0[j][i+1] +
-			      Phi0[j-1][i] + Phi0[j+1][i]);
-                
-	  *error += delta[X] * delta[Y] * fabs(Phi[j][i] - Phi0[j][i]);
+	  Phi[j][i] = 0.25 * (Phi0[j  ][i-1] +
+			      Phi0[j  ][i+1] +
+			      Phi0[j-1][i  ] +
+			      Phi0[j+1][i  ]);
+
+	  /* avoid fabs from math library */
+	  const MyData diff = (Phi[j][i] - Phi0[j][i]);
+	  *error += ((diff > 0) ? diff : -diff);
 	} /* loop over columns */
     } /* loop over rows */
   
+  *error *= (delta[X] * delta[Y]);
+  
   return;
 }
 
diff --git a/jacobi_solutions/serial/opt/src/tools.c b/jacobi_solutions/serial/opt/src/tools.c
index be2761d170c23fde4d5d4aeda54af00e56abf3ce..5e1f75a00f79bb5f8d32a94334a4f36e21730be9 100644
--- a/jacobi_solutions/serial/opt/src/tools.c
+++ b/jacobi_solutions/serial/opt/src/tools.c
@@ -51,9 +51,8 @@ void Show_2DdblArray(const MyData **const A,
 
 double seconds()
 {
-  struct timeval tmp;
-  gettimeofday(&tmp, (struct timezone *)0);
-  double sec = tmp.tv_sec + ((double)tmp.tv_usec)/1000000.0;
-
-  return sec;
+  struct timespec ts;
+  return (clock_gettime( CLOCK_PROCESS_CPUTIME_ID, &ts ),
+	   (double)ts.tv_sec +
+	   (double)ts.tv_nsec * 1e-9);
 }