diff --git a/.gitattributes b/.gitattributes
index a4c6b72429291661e6dae6929b4e898daa125ae3..ddfa7937b849cf7d567d10f873f4ca50321c5e3f 100644
--- a/.gitattributes
+++ b/.gitattributes
@@ -3941,6 +3941,7 @@ RTCP/Cobalt/GPUProc/test/t_cpu_utils.cc -text
 RTCP/Cobalt/GPUProc/test/t_cpu_utils.in_parset -text
 RTCP/Cobalt/GPUProc/test/t_cpu_utils.run -text
 RTCP/Cobalt/GPUProc/test/t_cpu_utils.sh -text
+RTCP/Cobalt/GPUProc/test/tcmpfloat.sh eol=lf
 RTCP/Cobalt/GPUProc/test/testParset.sh eol=lf
 RTCP/Cobalt/GPUProc/test/tstartBGL.in_parset -text
 RTCP/Cobalt/GPUProc/test/tstartBGL.run -text
diff --git a/RTCP/Cobalt/GPUProc/test/CMakeLists.txt b/RTCP/Cobalt/GPUProc/test/CMakeLists.txt
index 7344cd0eb170bdaa3914a853fd9bd133a66f9cc8..9ec33813aff43ae8cf5685a0e6d0a62f9e292d48 100644
--- a/RTCP/Cobalt/GPUProc/test/CMakeLists.txt
+++ b/RTCP/Cobalt/GPUProc/test/CMakeLists.txt
@@ -15,6 +15,7 @@ lofar_add_test(tProductionParsets DEPENDS rtcp)
 
 # cmpfloat is started by scripts for a fuzzy compare of 2 output files with raw floats
 lofar_add_executable(cmpfloat cmpfloat.cc)
+lofar_add_test(tcmpfloat DEPENDS cmpfloat)
 
 # tests that use testParset.sh
 
diff --git a/RTCP/Cobalt/GPUProc/test/cmpfloat.cc b/RTCP/Cobalt/GPUProc/test/cmpfloat.cc
index e63688f9922a45357ec92d4001e7b92afd5c2a75..136d2a8e03604ad20810ffdcba2c7b4f3979512c 100644
--- a/RTCP/Cobalt/GPUProc/test/cmpfloat.cc
+++ b/RTCP/Cobalt/GPUProc/test/cmpfloat.cc
@@ -1,4 +1,4 @@
-//# cmpfloat.cc
+//# cmpfloat.cc: compare floating point values between two binary files
 //# Copyright (C) 2013  ASTRON (Netherlands Institute for Radio Astronomy)
 //# P.O. Box 2, 7990 AA Dwingeloo, The Netherlands
 //#
@@ -21,152 +21,337 @@
 //#include <lofar_config.h>
 
 #include <cstdlib>
+#include <complex>
+#include <string>
+#include <vector>
 #include <iostream>
 #include <fstream>
+#include <limits>
+#include <boost/lexical_cast.hpp>
 
 #include "fpequals.h"
 
-using std::cout;
-using std::cerr;
-using std::endl;
-using std::ifstream;
+namespace {
 
+using namespace std;
 using LOFAR::Cobalt::fpEquals;
 
+struct Args {
+  string filename1;
+  string filename2;
+  enum type {FLOAT, DOUBLE, CFLOAT, CDOUBLE} type;
+  size_t skip;
+  size_t nvals;
+  double epsilon;
+  bool verbose;
+};
 
-int main(int argc, char *argv[])
-{
-  cerr.precision(8); // print full float precision (7 + the 0.).
+void pr_usage(const char* progname) {
+  cerr << "Usage: " << progname << " [--type=float|double|cfloat|cdouble]"
+                                   " [--skip=bytes] [--size=nvals]"
+                                   " [--epsilon=fpval] [--verbose] file1 file2" << endl;
+  cerr << "  --type     interpret input data as array of specified type" << endl;
+  cerr << "             cfloat means complex float. Default: double" << endl;
+  cerr << "  --skip     number of bytes to skip before comparison starts. Default: 0" << endl;
+  cerr << "  --size     must compare number of values of type. Default: until EOF" << endl;
+  cerr << "  --epsilon  maximum absolute difference tolerance. Default: std::numeric_limits<T>::epsilon()" << endl;
+  cerr << "  --verbose  print some info to stdout, regardless of exit status" << endl;
+}
 
-  // Default epsilon.
-  double epsilon = std::numeric_limits<float>::epsilon();
+bool parseArgs(int argc, char *argv[], Args &args) {
+  // defaults
+  args.type = args.DOUBLE;
+  args.skip = 0;
+  args.nvals = std::numeric_limits<size_t>::max();
+  args.epsilon = std::numeric_limits<double>::epsilon();
+  args.verbose = false;
 
-  if (argc < 3 || argc > 4)
-  {
-    cerr << "Usage: " << argv[0] << " [" << epsilon << "] <file1> <file2>" << endl;
-    cerr << "  where the optional floating point argument overrides the comparison epsilon" << endl;
-    return 1;
-  }
+  const string typePrefix("--type=");
+  const string skipPrefix("--skip=");
+  const string sizePrefix("--size=");
+  const string epsilonPrefix("--epsilon=");
+  const string verbosePrefix("--verbose");
 
-  char *filename1;
-  char *filename2;
-
-  if (argc == 3)
-  {
-    filename1 = argv[1];
-    filename2 = argv[2];
-  } else { // argc == 4
-    filename1 = argv[2];
-    filename2 = argv[3];
-
-    epsilon = std::atof(argv[1]);
-    if (epsilon <= 0.0 || epsilon > 1.0)
-    {
-      cerr << "Epsilon command line argument is out of range" << endl;
-      return 1;
-    } else {
-      cout << "Using an epsilon of " << epsilon << endl;
+  bool ok = true;
+  bool epsSet = false;
+  unsigned nfiles = 0;
+
+  for (int i = 1; i < argc; i++) {
+    string opt(argv[i]);
+    string val;
+
+    if (opt.compare(0, typePrefix.size(), typePrefix) == 0) {
+      val = opt.erase(0, typePrefix.size());
+
+      if (val == "float") {
+        args.type = args.FLOAT;
+        if (!epsSet)
+          args.epsilon = std::numeric_limits<float>::epsilon();
+      } else if (val == "double") {
+        args.type = args.DOUBLE;
+        if (!epsSet)
+          args.epsilon = std::numeric_limits<double>::epsilon();
+      } else if (val == "cfloat") {
+        args.type = args.CFLOAT;
+        if (!epsSet)
+          args.epsilon = std::numeric_limits<float>::epsilon();
+      } else if (val == "cdouble") {
+        args.type = args.CDOUBLE;
+        if (!epsSet)
+          args.epsilon = std::numeric_limits<double>::epsilon();
+      } else {
+        cerr << "Error: invalid value in --type argument: " << val << endl;
+        ok = false;
+      }
+    } else if (opt.compare(0, skipPrefix.size(), skipPrefix) == 0) {
+      val = opt.erase(0, skipPrefix.size());
+      try {
+        args.skip = boost::lexical_cast<size_t>(val);
+        if ((ssize_t)args.skip < 0)
+          throw boost::bad_lexical_cast();
+      } catch (boost::bad_lexical_cast& exc) {
+        cerr << "Error: invalid value in --skip argument: " << val << endl;
+        ok = false;
+      }
+    } else if (opt.compare(0, sizePrefix.size(), sizePrefix) == 0) {
+      val = opt.erase(0, sizePrefix.size());
+      try {
+        args.nvals = boost::lexical_cast<size_t>(val);
+        if ((ssize_t)args.nvals < 0)
+          throw boost::bad_lexical_cast();
+      } catch (boost::bad_lexical_cast& exc) {
+        cerr << "Error: invalid value in --size argument: " << val << endl;
+        ok = false;
+      }
+    } else if (opt.compare(0, epsilonPrefix.size(), epsilonPrefix) == 0) {
+      val = opt.erase(0, epsilonPrefix.size());
+      try {
+        args.epsilon = boost::lexical_cast<double>(val);
+        args.epsilon = std::abs(args.epsilon);
+        epsSet = true;
+      } catch (boost::bad_lexical_cast& exc) {
+        cerr << "Error: invalid value in --epsilon argument: " << val << endl;
+        ok = false;
+      }
+    } else if (opt == verbosePrefix) {
+      args.verbose = true;
+    } else { // filename
+      if (nfiles == 0) {
+        args.filename1 = opt;
+      } else if (nfiles == 1) {
+        args.filename2 = opt;
+      }
+      nfiles += 1;
     }
   }
-  float eps = (float)epsilon; // atm, we only cmp single precision floats
 
-  ifstream ifs1(filename1, std::ios::binary);
-  if (!ifs1)
-  {
-    cerr << "Failed to open file " << filename1 << endl;
-    return 1;
+  if (nfiles != 2) {
+    cerr << "Error: need 2 file arguments, got " << nfiles << endl;
+    ok = false;
   }
 
-  ifstream ifs2(filename2, std::ios::binary);
-  if (!ifs2)
-  {
-    cerr << "Failed to open file " << filename2 << endl;
-    return 1;
+  return ok;
+}
+
+template <typename T>
+bool compareValues(T v1, T v2, double epsilon, size_t pos,
+                   T& maxFactor, T& minFactor) {
+  if (!fpEquals(v1, v2, (T)epsilon)) {
+    cerr << "Error: value diff beyond epsilon at compared value " << pos << ": "
+         << v1 << " " << v2 << endl;
+
+    T factor = v2 / v1; // inf is fine, NaN if eps was set to 0 or odd data
+    if (maxFactor == T(1.0)) {
+      // first unequal val, so 1.0 must be as initialized (not a factor)
+      maxFactor = minFactor = factor;
+    } else if (factor > maxFactor) {
+      maxFactor = factor;
+    } else if (factor < minFactor) {
+      minFactor = factor;
+    }
+
+    return false;
+  }
+
+  return true;
+}
+
+// Note the plural form of the complex identifiers: both factors are in the cval
+template <typename T>
+bool compareValues(complex<T> v1, complex<T> v2, double epsilon, size_t pos,
+                   complex<T>& maxFactors, complex<T>& minFactors) {
+  if (!fpEquals(v1, v2, (T)epsilon)) {
+    cerr << "Error: value diff beyond epsilon at compared value " << pos << ": "
+         << v1 << " " << v2 << endl;
+
+    T realFactor = v2.real() / v1.real(); // idem as above
+    T imagFactor = v2.imag() / v1.imag(); // idem
+    if (maxFactors == T(1.0)) {
+      // first unequal val, so 1.0 must be as initialized (not a factor)
+      maxFactors.real() = minFactors.real() = realFactor;
+      maxFactors.imag() = minFactors.imag() = imagFactor;
+    } else {
+      if (realFactor > maxFactors.real()) {
+        maxFactors.real(realFactor);
+      } else if (realFactor < minFactors.real()) {
+        minFactors.real(realFactor);
+      }
+      if (imagFactor > maxFactors.imag()) {
+        maxFactors.imag(imagFactor);
+      } else if (realFactor < minFactors.imag()) {
+        minFactors.imag(imagFactor);
+      }
+    }
+
+    return false;
   }
 
-  const size_t bufLen = 2048;
-  float *buf1 = new float[bufLen];
-  float *buf2 = new float[bufLen];
+  return true;
+}
 
-  int status = 0;
-  size_t total = 0;
+template <typename T>
+void printCommonFactorMessage(const T& maxFactor, const T& minFactor) {
+  // If maxFactor and minFactor are near, then the diff is probably scale only.
+  const T facEps = (T)1e-1; // very loose as values can be of any magnitude (but too loose if no scale, yet min and max < 1e-1. The program then still exits non-zero, but possibly with the wrong message.)
+  bool eq = fpEquals(maxFactor, minFactor, facEps);
+  T avgFac;
+  if (eq) {
+    avgFac = (T)0.5 * (maxFactor + minFactor);
+    cerr << "All errors of vals for this pair of files are within "
+         << facEps << " to a factor " << avgFac << " (inverse="
+         << (T)1.0 / avgFac << ')' << endl;
+  } else
+    cerr << "No clear common factor among all errors: maxFactor=" <<
+            maxFactor << "; minFactor=" << minFactor << endl;
+}
 
-  float maxFactorOff = 0.0f;
-  float minFactorOff = 0.0f;
+// Note the plural form of the complex identifiers: both factors are in the cval
+template <typename T>
+void printCommonFactorMessage(const complex<T>& maxFactors,
+                              const complex<T>& minFactors) {
+  // If maxFactor and minFactor are near, then the diff is probably scale only.
+  // For complex types, also see if real and imag are *-1 of each other (conj).
+  const T facEps = (T)1e-1; // very loose as values can be of any magnitude
+  bool realEq = fpEquals(maxFactors.real(), minFactors.real(), facEps);
+  bool imagEq = fpEquals(maxFactors.imag(), minFactors.imag(), facEps);
+  T avgRealFac, avgImagFac;
+  if (realEq) {
+    avgRealFac = (T)0.5 * (maxFactors.real() + minFactors.real());
+    cerr << "All errors of real vals for this pair of files are within "
+         << facEps << " to a factor " << avgRealFac << " (inverse="
+         << (T)1.0 / avgRealFac << ')' << endl;
+  }
+  if (imagEq) {
+    avgImagFac = (T)0.5 * (maxFactors.imag() + minFactors.imag());
+    cerr << "All errors of imag vals for this pair of files are within "
+         << facEps << " to a factor " << avgImagFac << " (inverse="
+         << (T)1.0 / avgImagFac << ')' << endl;
+  }
+  if (realEq && imagEq && fpEquals(avgRealFac, -avgImagFac, facEps))
+      cerr << "Common real and imag factors appear to (also) differ roughly by "
+           << "a factor -1.0 (likely conjugation error)" << endl;
+  if (!(realEq && imagEq))
+    cerr << "No clear common factor among all errors: maxFactors="
+         << maxFactors << "; minFactors=" << minFactors << endl;
+}
 
-  while (ifs1.good() && ifs2.good()) {
-    size_t len = bufLen;
-    size_t nbytes1, nbytes2;
+template <typename T>
+bool compareStreams(ifstream& ifs1, ifstream& ifs2, size_t skipped,
+             size_t nvals, double epsilon, bool verbose) {
+  bool ok = true;
 
-    ifs1.read(reinterpret_cast<char *>(buf1), bufLen);
-    nbytes1 = ifs1.gcount();
+  cerr.precision(17); // print full double precision on errors
+  T maxFactor = T(1.0); // cmp needs init; 1.0 is not a valid unequality factor
+  T minFactor = T(1.0); // symmetry / good practice, but not used
 
-    ifs2.read(reinterpret_cast<char *>(buf2), bufLen);
-    nbytes2 = ifs2.gcount();
+  size_t i;
+  for (i = 0; i < nvals; i++) {
+    T v1, v2;
+    ifs1.read(reinterpret_cast<char *>(&v1), sizeof(T));
+    ifs2.read(reinterpret_cast<char *>(&v2), sizeof(T));
 
-    if (nbytes1 == 0 || nbytes1 != nbytes2 ||
-        nbytes1 % sizeof(float) != 0 || nbytes2 % sizeof(float) != 0)
-    {
-      cerr << "Failed to read an equal amount of bytes of at least a float from both input streams" << endl;
-      status = 1;
+    // Simultaneous EOF is ok iff nvals wasn't set as prog arg (default is max).
+    bool eof1 = ifs1.eof();
+    bool eof2 = ifs2.eof();
+    if (eof1 && eof2 && nvals == std::numeric_limits<size_t>::max()) {
+      break;
+    } else if (eof1 || eof2) {
+      cerr << "Error: Unexpected EOF in (at least) one stream after comparing "
+           << i << " values" << endl;
+      ok = false;
       break;
     }
 
-    if (nbytes1 < len * sizeof(float))
-      len = nbytes1 / sizeof(float);
-
-    for (size_t i = 0; i < len; i++)
-    {
-      if (!fpEquals(buf1[i], buf2[i], eps))
-      {
-        cerr << "Error: value diff beyond eps at pos " << total + i << ": " << buf1[i] << " " << buf2[i] << endl;
-        status = 2;
-
-        // Try to detect a max and min factor to print at the end.
-        // If near, then the diff is probably scale only. Conj if -1 every odd...
-        float factor = buf1[i] / buf2[i];
-        if (maxFactorOff == 0.0f) {
-          // init
-          maxFactorOff = factor;
-          minFactorOff = factor;
-        } else if (factor > maxFactorOff) {
-          maxFactorOff = factor;
-        } else if (factor < minFactorOff) {
-          minFactorOff = factor;
-        }
-      }
+    size_t nread1 = ifs1.gcount();
+    size_t nread2 = ifs2.gcount();
+    if (nread1 != nread2 || nread1 < sizeof(T)) {
+      cerr << "Failed to read enough data from both streams for another "
+           << "comparison after comparing " << i << " values" << endl;
+      ok = false;
+      break;
     }
 
-    total += len;
+    ok &= compareValues(v1, v2, epsilon, i, maxFactor, minFactor);
   }
 
-  // If cmp error, see if we can easily detect a scale-only error.
-  if (status == 2)
-  {
-    const float facEps = 1e-1; // needs to be very loose as value pairs can be of any magnitude
-    if (std::abs(maxFactorOff - minFactorOff) < facEps)
-      cerr << "All errors of vals for this pair of files are within " <<
-              facEps << " to a factor " << 0.5f * (maxFactorOff + minFactorOff) << endl;
-    else
-      cerr << "No clear common factor among all errors: maxFactor=" <<
-              maxFactorOff << " minFactor=" << minFactorOff << endl;
-  }
+  if (verbose)
+    cout << "Compared " << i << " values (after skipping " << skipped
+         << " bytes)" << endl;
+
+  if (!ok && i > 0)
+    printCommonFactorMessage(maxFactor, minFactor);
+
+  return ok;
+}
+
+} // anon namespace
 
-  if (!ifs1.eof())
-  {
-    cerr << "Error occurred while reading from file " << filename1 << endl;
-    status = 1;
+int main(int argc, char *argv[]) {
+  Args args;
+  if (!parseArgs(argc, argv, args)) {
+    cerr << endl;
+    pr_usage(argv[0]);
+    return 2;
   }
 
-  if (!ifs2.eof())
-  {
-    cerr << "Error occurred while reading from file " << filename2 << endl;
-    status = 1;
+  // open files
+  ifstream ifs1(args.filename1.c_str(), std::ios::binary);
+  if (!ifs1) {
+    cerr << "Failed to open file " << args.filename1 << endl;
+    return 2;
+  }
+  ifstream ifs2(args.filename2.c_str(), std::ios::binary);
+  if (!ifs2) {
+    cerr << "Failed to open file " << args.filename2 << endl;
+    return 2;
   }
 
-  delete[] buf2;
-  delete[] buf1;
+  // skip bytes (e.g. file header)
+  // Don't check, as it turns out that this does not fail if skip > file size.
+  ifs1.seekg(args.skip);
+  ifs2.seekg(args.skip);
+
+  // compare
+  if (args.verbose)
+    cout << "Comparing using an epsilon of " << args.epsilon << endl;
+  bool cmpOk;
+  if (args.type == args.FLOAT)
+    cmpOk = compareStreams<float>(ifs1, ifs2, args.skip, args.nvals,
+                                  args.epsilon, args.verbose);
+  else if (args.type == args.DOUBLE)
+    cmpOk = compareStreams<double>(ifs1, ifs2, args.skip, args.nvals,
+                                   args.epsilon, args.verbose);
+  else if (args.type == args.CFLOAT)
+    cmpOk = compareStreams<complex<float> >(ifs1, ifs2, args.skip, args.nvals,
+                                            args.epsilon, args.verbose);
+  else if (args.type == args.CDOUBLE)
+    cmpOk = compareStreams<complex<double> >(ifs1, ifs2, args.skip, args.nvals,
+                                             args.epsilon, args.verbose);
+  else {
+      cerr << "Internal error: unknown data type" << endl;
+      return 2;
+  }
 
-  return status;
+  return cmpOk ? 0 : 1;
 }
 
diff --git a/RTCP/Cobalt/GPUProc/test/tcmpfloat.py b/RTCP/Cobalt/GPUProc/test/tcmpfloat.py
new file mode 100755
index 0000000000000000000000000000000000000000..778e557d251b8a44c06445b8b99ed20502969c0f
--- /dev/null
+++ b/RTCP/Cobalt/GPUProc/test/tcmpfloat.py
@@ -0,0 +1,91 @@
+#!/usr/bin/python
+# tcmpfloat.py: generate binary input files for tcmpfloat.sh test
+# Copyright (C) 2013  ASTRON (Netherlands Institute for Radio Astronomy)
+# P.O. Box 2, 7990 AA Dwingeloo, The Netherlands
+#
+# This file is part of the LOFAR software suite.
+# The LOFAR software suite is free software: you can redistribute it and/or
+# modify it under the terms of the GNU General Public License as published
+# by the Free Software Foundation, either version 3 of the License, or
+# (at your option) any later version.
+#
+# The LOFAR software suite is distributed in the hope that it will be useful,
+# but WITHOUT ANY WARRANTY; without even the implied warranty of
+# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
+# GNU General Public License for more details.
+#
+# You should have received a copy of the GNU General Public License along
+# with the LOFAR software suite. If not, see <http://www.gnu.org/licenses/>.
+#
+# $Id$
+
+import numpy as np
+
+def main():
+  # Generate all binary input files
+  # to inspect raw floats:  od -t fF file.bin
+  # to inspect raw doubles: od -t fD file.bin
+
+  # Test 1
+  a = [1.0 + 2.0j, 3.0 + 4.0j, -1.1 + -1.2j, -2.1 + -2.2j]
+  ar = np.array(a, dtype=np.complex64)
+  ar.tofile('tcmpfloat-1.1.bin')
+  ar.tofile('tcmpfloat-1.2.bin')
+
+  # Test 2
+  a = [1.0 + 2.0j, 3.0 + 4.0j,  -1.1 + -1.2j, -2.1 + -2.2j, 10.0 + 20.0j]
+  b = [1.1 + 2.1j, 3.1 + 4.1j,  -1.1 + -1.2j, -2.1 + -2.2j, 11.0 + 21.0j]
+  ar = np.array(a, dtype=np.complex128)
+  br = np.array(b, dtype=np.complex128)
+  ar.tofile('tcmpfloat-2.1.bin')
+  br.tofile('tcmpfloat-2.2.bin')
+
+  # Test 3
+  a = [-5.50001 + 6.60001j, -7.70001 + 8.80001j]
+  b = [ 5.50002 + 6.60000j, -7.70000 - 8.80002j]
+  ar = np.array(a, dtype=np.complex128)
+  br = np.array(b, dtype=np.complex128)
+  ar.tofile('tcmpfloat-3.1.bin')
+  br.tofile('tcmpfloat-3.2.bin')
+
+  # Test 4
+  a = [ 1.000000,  2.000000,  3.000000,  4.000000,  5.000000,   6.000000,   7.000000,   8.000000,   9.000000]
+  b = [-2.000001, -4.000001, -6.000001, -8.000001, -9.999998, -12.000001, -13.999999, -16.000000, -18.000001]
+  ar = np.array(a, dtype=np.float32)
+  br = np.array(b, dtype=np.float32)
+  ar.tofile('tcmpfloat-4.1.bin')
+  br.tofile('tcmpfloat-4.2.bin')
+
+  # Test 5
+  a = [-3.000000 - 3.000000j, -2.000000 - 2.000000j, -1.000000 - 1.000000j, 0.000000 + 0.000000j, 1.000000 + 1.000000j, 2.000000 + 2.000000j, 3.000000 + 3.000000j]
+  b = [-3.000000 + 3.000001j, -2.000000 + 2.000001j, -1.000000 + 0.999999j, 0.000000 + 0.000000j, 1.000000 - 1.000001j, 2.000000 - 1.999999j, 3.000000 - 2.999998j]
+  ar = np.array(a, dtype=np.complex64)
+  br = np.array(b, dtype=np.complex64)
+  ar.tofile('tcmpfloat-5.1.bin')
+  br.tofile('tcmpfloat-5.2.bin')
+
+  # Test 6
+  a = [ -3.000000 -  3.000000j, -2.000000 - 2.000000j, -1.000000 - 1.000000j, 0.000000 + 0.000000j, 1.000000 + 1.000000j, 2.000000 + 2.000000j,  3.000000 +  3.000000j]
+  b = [-12.000000 + 12.000001j, -8.000000 + 8.000001j, -4.000000 + 3.999999j, 0.000000 + 0.000000j, 4.000000 - 4.000001j, 8.000000 - 7.999999j, 12.000000 - 11.999998j]
+  ar = np.array(a, dtype=np.complex128)
+  br = np.array(b, dtype=np.complex128)
+  ar.tofile('tcmpfloat-6.1.bin')
+  br.tofile('tcmpfloat-6.2.bin')
+
+  # Test 7 needs no input
+
+  # Test 8 needs no 2nd input
+  ar.tofile('tcmpfloat-8.1.bin')
+
+  # Test 9
+  a = [1.0, 1.0, 1.0, 1.0, 1.0]
+  b = [2.0, 1.0, 1.0, 1.0]
+  ar = np.array(a, dtype=np.float64)
+  br = np.array(b, dtype=np.float64)
+  ar.tofile('tcmpfloat-9.1.bin')
+  br.tofile('tcmpfloat-9.2.bin')
+
+
+if __name__ == "__main__":
+  main()
+
diff --git a/RTCP/Cobalt/GPUProc/test/tcmpfloat.sh b/RTCP/Cobalt/GPUProc/test/tcmpfloat.sh
new file mode 100755
index 0000000000000000000000000000000000000000..19430dd8d984f696462e002db7f2995bfa172935
--- /dev/null
+++ b/RTCP/Cobalt/GPUProc/test/tcmpfloat.sh
@@ -0,0 +1,71 @@
+#!/bin/sh
+# Tests for the cmpfloat.cc program.
+#
+# $Id$
+
+# generate binary input files through tcmpfloat.py
+./runctest.sh tcmpfloat
+if [ $? -ne 0 ]; then echo "Failed to generate test input files"; exit 1; fi
+
+
+status=0
+
+# Test 1: compare a few complex floats
+echo "Test 1"
+./cmpfloat --type=cfloat --verbose tcmpfloat-1.1.bin tcmpfloat-1.2.bin
+if [ $? -ne 0 ]; then echo "TEST ERROR: Test 1 (compare complex floats) failed"; status=1; fi
+
+# Test 2: compare a few complex doubles, skip some bytes and limit nr compared values
+echo "Test 2"
+./cmpfloat --type=cdouble --skip=32 --size=2 --verbose tcmpfloat-2.1.bin tcmpfloat-2.2.bin
+if [ $? -ne 0 ]; then echo "TEST ERROR: Test 2 (compare complex double with skip, size) failed"; status=1; fi
+
+# Test 3: simple double comparison that fails, no scale factor
+echo "Test 3"
+./cmpfloat --verbose tcmpfloat-3.1.bin tcmpfloat-3.2.bin
+if [ $? -ne 1 ]; then echo "TEST ERROR: Test 3 (compare double fails, no scale) failed"; status=1; fi
+
+# Test 4: simple float comparison that fails with scale factor
+echo "Test 4"
+./cmpfloat --type=float --verbose tcmpfloat-4.1.bin tcmpfloat-4.2.bin > tcmpfloat-4.out 2>&1
+if [ $? -ne 1 ]; then echo "TEST ERROR: Test 4 (compare float fails, scale factor) failed"; status=1; fi
+cat tcmpfloat-4.out
+grep inverse tcmpfloat-4.out > /dev/null
+if [ $? -ne 0 ]; then echo "TEST ERROR: Test 4: scale factor (and its inverse) not found"; status=1; fi
+
+# Test 5: simple complex float comparison that fails, because of conjugation
+echo "Test 5"
+./cmpfloat --type=cfloat --verbose tcmpfloat-5.1.bin tcmpfloat-5.2.bin > tcmpfloat-5.out 2>&1
+if [ $? -ne 1 ]; then echo "TEST ERROR: Test 5 (compare complex float fails, conjugation) failed"; status=1; fi
+cat tcmpfloat-5.out
+grep conjugation tcmpfloat-5.out > /dev/null
+if [ $? -ne 0 ]; then echo "TEST ERROR: Test 5: conjugation error not found as such"; status=1; fi
+
+# Test 6: simple complex double comparison that fails with scale factor and conjugation
+echo "Test 6"
+./cmpfloat --type=cdouble --verbose tcmpfloat-6.1.bin tcmpfloat-6.2.bin > tcmpfloat-6.out 2>&1
+if [ $? -ne 1 ]; then echo "TEST ERROR: Test 6 (compare complex double fails, scale factor and conjugation) failed"; status=1; fi
+cat tcmpfloat-6.out
+grep inverse tcmpfloat-6.out > /dev/null
+if [ $? -ne 0 ]; then echo "TEST ERROR: Test 6: scale factor (and its inverse) not found"; status=1; fi
+grep conjugation tcmpfloat-6.out > /dev/null
+if [ $? -ne 0 ]; then echo "TEST ERROR: Test 6: conjugation error not found as such"; status=1; fi
+
+# Test 7: missing filenames
+echo "Test 7"
+./cmpfloat --type=cdouble --verbose
+if [ $? -ne 2 ]; then echo "TEST ERROR: Test 7 (missing filenames) failed"; status=1; fi
+
+# Test 8: non-existing file
+echo "Test 8"
+./cmpfloat tcmpfloat-8.1.bin tcmpfloat-8.2-non-existing.bin
+if [ $? -ne 2 ]; then echo "TEST ERROR: Test 8 (non-existing file) failed"; status=1; fi
+
+# Test 9: after skipping 32 bytes, one file has <4 doubles left to compare
+echo "Test 9"
+./cmpfloat --skip=8 --size=4 tcmpfloat-9.1.bin tcmpfloat-9.2.bin
+if [ $? -ne 1 ]; then echo "TEST ERROR: Test 9 (file too short after skip) failed"; status=1; fi
+
+if [ $status -eq 0 ]; then echo -e "\nAll tcmpfloat tests PASSED"; fi
+exit $status
+
diff --git a/RTCP/Cobalt/GPUProc/test/testParset.sh b/RTCP/Cobalt/GPUProc/test/testParset.sh
index 59b9448fa1a4a9ad3df31d9d34f4dba70be4a8fd..35aed81f1f81720aad83ca5e97ae57906800d0b2 100755
--- a/RTCP/Cobalt/GPUProc/test/testParset.sh
+++ b/RTCP/Cobalt/GPUProc/test/testParset.sh
@@ -120,7 +120,7 @@ function parse_logs
 
       for f in *.MS
       do
-        $testdir/cmpfloat $EPSILON `pwd`/$f $REFDIR/$f || error "Output does not match reference for eps_factor=$eps_factor"
+        $testdir/cmpfloat --type=cfloat --epsilon=$EPSILON --verbose `pwd`/$f $REFDIR/$f || error "Output does not match reference for eps_factor=$eps_factor"
       done
     done
   fi