diff --git a/Src/BalanceTree/FCostCell.hpp b/Src/BalanceTree/FCostCell.hpp
index 771b1d3ef95f26d5da9e42e5f81390a701771b91..3cf9cb862b75e5fb593737b9351689f0e4d4c0bf 100644
--- a/Src/BalanceTree/FCostCell.hpp
+++ b/Src/BalanceTree/FCostCell.hpp
@@ -7,6 +7,17 @@
 
 #include <type_traits>
 
+
+/** 
+ * \brief Empty trait class.
+ * \author Quentin Khan
+ *
+ * This class is used to check whether a cell class has FCostCell in its
+ * inheritance tree.
+ */
+class FCostCellTypeTrait {};
+
+
 /**
  * \brief Cell with a cost memory for balance computations.
  * \author Quentin Khan
@@ -18,7 +29,7 @@
  * \tparam CostType The type to use in order to store the cost. Defaults to FSize.
  */
 template<typename BaseClass, typename CostType = FSize>
-class FCostCell : public BaseClass {
+class FCostCell : public BaseClass, virtual public FCostCellTypeTrait {
     static_assert(std::is_arithmetic<CostType>::value,
                   "The cell cost type must be an arithmetic type.");
 
diff --git a/Tests/noDist/AlgoLoaderCostZones.hpp b/Tests/noDist/AlgoLoaderCostZones.hpp
new file mode 100644
index 0000000000000000000000000000000000000000..a633f3ccdd0695e8beeed9f6acd11a91c3542f3b
--- /dev/null
+++ b/Tests/noDist/AlgoLoaderCostZones.hpp
@@ -0,0 +1,82 @@
+#ifndef _ALGOLOADERCOSTZONES_HPP_
+#define _ALGOLOADERCOSTZONES_HPP_
+
+#include "PerfTestUtils.hpp"
+
+#include "Core/FFmmAlgorithm.hpp"
+
+#include "BalanceTree/FFmmAlgorithmThreadBalanced.hpp"
+#include "BalanceTree/FCostCell.hpp"
+#include "BalanceTree/FCostZones.hpp"
+
+/**
+ * \brief Algorithm loader for the CostZones algorithm.
+ *
+ * \warning : This loader requires that the KernelLoader supply a type definition
+ * for a `CostKernelClass`
+ */
+template <class _TreeLoader, template<typename> class _KernelLoader>
+class AlgoLoaderCostZones : public FAlgoLoader<_TreeLoader, _KernelLoader> {
+public:
+    using TreeLoader     = _TreeLoader;
+    using KernelLoader   = _KernelLoader<TreeLoader>;
+
+    using FReal          = typename TreeLoader::FReal;
+    using CellClass      = typename TreeLoader::CellClass;
+    using ContainerClass = typename TreeLoader::ContainerClass;
+    using LeafClass      = typename TreeLoader::LeafClass;
+    using OctreeClass    = typename TreeLoader::OctreeClass;
+    using KernelClass    = typename KernelLoader::KernelClass;
+    using CostKernelClass= typename KernelLoader::CostKernelClass;
+
+    static_assert(std::is_base_of<FCostCellTypeTrait, CellClass>::value,
+        "The tree cells must derive from FCostCell.");
+
+    using FMMClass = FFmmAlgorithmThreadBalanced
+        <OctreeClass, CellClass, ContainerClass, KernelClass, LeafClass>;    
+    using CostFmmClass = FFmmAlgorithm
+        <OctreeClass, CellClass, ContainerClass, CostKernelClass, LeafClass>;
+
+    TreeLoader& _treeLoader;
+    KernelLoader& _kernelLoader;
+
+    /// Builds the loader
+    AlgoLoaderCostZones(FPerfTestParams& /*params*/,
+                        TreeLoader& treeLoader,
+                        KernelLoader& kernelLoader) :
+        _treeLoader(treeLoader),
+        _kernelLoader(kernelLoader) {
+        
+    }
+
+    /// Computes the tree cells costs then runs the costzones and FMM algorithms.
+    void run() {
+
+        OctreeClass* p_tree = &(_treeLoader._tree);
+
+        // Compute tree cells costs
+        CostFmmClass costAlgo(p_tree, &(_kernelLoader._costKernel));
+
+        this->time.tic();
+        costAlgo.execute();
+        this->time.tac();
+        std::cout << "Generating tree cost: " << this->time.elapsed() << "s.\n";
+
+        FCostZones<OctreeClass, CellClass> costzones(p_tree, omp_get_max_threads());
+
+        this->time.tic();
+        costzones.run();
+        this->time.tac();
+        std::cout << "Generating cost zones: " << this->time.elapsed() << "s.\n";
+
+
+        this->time.tic();
+        FMMClass algo(p_tree, &(_kernelLoader._kernel), costzones.getZoneBounds(), costzones.getLeafZoneBounds());
+        algo.execute();
+        this->time.tac();
+    }
+};
+
+
+
+#endif
diff --git a/Tests/noDist/AlgoLoaderTask.hpp b/Tests/noDist/AlgoLoaderTask.hpp
new file mode 100644
index 0000000000000000000000000000000000000000..42d9dd9ec90cff5d38ef6a75404141870b84bdc3
--- /dev/null
+++ b/Tests/noDist/AlgoLoaderTask.hpp
@@ -0,0 +1,45 @@
+#ifndef _ALGOLOADERTASK_HPP_
+#define _ALGOLOADERTASK_HPP_
+
+#include "PerfTestUtils.hpp"
+
+#include "Core/FFmmAlgorithmTask.hpp"
+
+
+template <class _TreeLoader, template<typename> class _KernelLoader>
+class AlgoLoaderTask : public FAlgoLoader<_TreeLoader, _KernelLoader> {
+public:
+    using TreeLoader     = _TreeLoader;
+    using KernelLoader   = _KernelLoader<TreeLoader>;
+
+
+    using FReal          = typename TreeLoader::FReal;
+    using CellClass      = typename TreeLoader::CellClass;
+    using ContainerClass = typename TreeLoader::ContainerClass;
+    using LeafClass      = typename TreeLoader::LeafClass;
+    using OctreeClass    = typename TreeLoader::OctreeClass;
+    using KernelClass    = typename KernelLoader::KernelClass;
+
+    using FMMClass = FFmmAlgorithmTask<OctreeClass, CellClass, ContainerClass, KernelClass, LeafClass>;
+    
+    TreeLoader& _treeLoader;
+    KernelLoader& _kernelLoader;
+
+    AlgoLoaderTask(FPerfTestParams& /*params*/,
+                   TreeLoader& treeLoader,
+                   KernelLoader& kernelLoader) :
+        _treeLoader(treeLoader),
+        _kernelLoader(kernelLoader) {
+        
+    }
+
+
+    void run() {
+        FMMClass algo(&(_treeLoader._tree), &(_kernelLoader._kernel));
+        algo.execute();
+    }
+};
+
+
+
+#endif
diff --git a/Tests/noDist/AlgoLoaderThread.hpp b/Tests/noDist/AlgoLoaderThread.hpp
new file mode 100644
index 0000000000000000000000000000000000000000..7ccec611b719e7026bbfd314fe1623eae954020c
--- /dev/null
+++ b/Tests/noDist/AlgoLoaderThread.hpp
@@ -0,0 +1,42 @@
+#ifndef _ALGOLOADERTHREAD_HPP_
+#define _ALGOLOADERTHREAD_HPP_
+
+#include "PerfTestUtils.hpp"
+
+#include "Core/FFmmAlgorithmThread.hpp"
+
+
+template <class _TreeLoader, template<typename> class _KernelLoader>
+class AlgoLoaderThread : public FAlgoLoader<_TreeLoader, _KernelLoader> {
+public:
+    using TreeLoader     = _TreeLoader;
+    using KernelLoader   = _KernelLoader<TreeLoader>;
+
+    using FReal          = typename TreeLoader::FReal;
+    using CellClass      = typename TreeLoader::CellClass;
+    using ContainerClass = typename TreeLoader::ContainerClass;
+    using LeafClass      = typename TreeLoader::LeafClass;
+    using OctreeClass    = typename TreeLoader::OctreeClass;
+    using KernelClass    = typename KernelLoader::KernelClass;
+
+    using FMMClass = FFmmAlgorithmThread<OctreeClass, CellClass, ContainerClass, KernelClass, LeafClass>;
+    
+    TreeLoader& _treeLoader;
+    KernelLoader& _kernelLoader;
+
+    AlgoLoaderThread(FPerfTestParams& /*params*/,
+                     TreeLoader& treeLoader,
+                     KernelLoader& kernelLoader) :
+        _treeLoader(treeLoader),
+        _kernelLoader(kernelLoader) {
+        
+    }
+
+
+    void run() {
+        FMMClass algo(&(_treeLoader._tree), &(_kernelLoader._kernel), false);
+        algo.execute();
+    }
+};
+
+#endif
diff --git a/Tests/noDist/KernelLoaderFChebSym.hpp b/Tests/noDist/KernelLoaderFChebSym.hpp
new file mode 100644
index 0000000000000000000000000000000000000000..e1cfd6e0e2f70a8e85a59efa46b6231061237016
--- /dev/null
+++ b/Tests/noDist/KernelLoaderFChebSym.hpp
@@ -0,0 +1,67 @@
+#ifndef _KERNELLOADERFCHEBSYM_HPP_
+#define _KERNELLOADERFCHEBSYM_HPP_
+
+#include "PerfTestUtils.hpp"
+
+#include "Kernels/Interpolation/FInterpMatrixKernel.hpp"
+#include "Kernels/Chebyshev/FChebSymKernel.hpp"
+
+#include "BalanceTree/FChebSymCostKernel.hpp"
+
+/**
+ * \brief Kernel loader for the symetric Chebyshev kernel.
+ *
+ * \warning This loader requires that TreeLoader::CellClass inherits from
+ * FChebCell.
+ *
+ * \note This loader also provides the typedef CostKernelClass and a member
+ * _costKernel that cam be used by the AlgoLoaderCostZones.
+ */
+template <typename _TreeLoader>
+class KernelLoaderFChebSym : public FKernelLoader<_TreeLoader> {
+
+public:
+    // Required type definitions
+    using TreeLoader     = _TreeLoader;
+    using FReal          = typename TreeLoader::FReal;
+    /// Must derive from FChebCell
+    using CellClass      = typename TreeLoader::CellClass;
+    using ContainerClass = typename TreeLoader::ContainerClass;
+    using OctreeClass    = typename TreeLoader::OctreeClass;
+
+    using MatrixKernelClass = FInterpMatrixKernelR<FReal>;
+    using KernelClass       = FChebSymKernel <FReal, CellClass, ContainerClass,
+                                              MatrixKernelClass, TreeLoader::ORDER>;
+    /// Kernel class used to compute the tree cell costs.
+    using CostKernelClass = FChebSymCostKernel<FReal, CellClass, ContainerClass, 
+                                               MatrixKernelClass, TreeLoader::ORDER,
+                                               OctreeClass>;
+
+    // Meaningfull (?) error message.
+    static_assert(std::is_base_of<FChebCell<FReal,TreeLoader::ORDER>,CellClass>::value,
+                  "TreeLoader::CellClass must derive from FChebCell");    
+
+
+    const FReal epsilon = 1e-4;
+    const MatrixKernelClass _matrixKernel;
+    /// kernel used to compute the tree cells interactions.
+    KernelClass _kernel;
+    /// Kernel used to compute the tree cells costs.
+    CostKernelClass _costKernel;
+
+    /// Builds the loader and loads the kernel.
+    KernelLoaderFChebSym(FPerfTestParams& /*params*/, TreeLoader& treeLoader) :
+        _matrixKernel(),
+        _kernel(treeLoader._tree.getHeight(),
+                treeLoader._tree.getBoxWidth(),
+                treeLoader._tree.getBoxCenter(),
+                &_matrixKernel),
+        _costKernel(&(treeLoader._tree), epsilon){
+
+    }
+
+
+};
+
+
+#endif
diff --git a/Tests/noDist/PerfTest.cpp b/Tests/noDist/PerfTest.cpp
new file mode 100644
index 0000000000000000000000000000000000000000..9f694d905f3599a216ef275448d8839b732334f6
--- /dev/null
+++ b/Tests/noDist/PerfTest.cpp
@@ -0,0 +1,85 @@
+
+/**
+ * \file
+ * \author Quentin Khan
+ *
+ * This program is used to run different performance tests for the various
+ * algorithms that have been implemented for ScalFMM.
+ *
+ * See the PerfUtils.hpp file classes for some more in depth information. Run
+ * with argument --help for usage information.
+ */
+
+
+#include <iostream>
+#include <string>
+
+#include "Utils/FParameters.hpp"
+#include "Utils/FParameterNames.hpp"
+
+#include "PerfTestUtils.hpp"
+
+#include "TreeLoaderFCheb.hpp"
+
+#include "KernelLoaderFChebSym.hpp"
+
+#include "AlgoLoaderThread.hpp"
+#include "AlgoLoaderTask.hpp"
+#include "AlgoLoaderCostZones.hpp"
+
+/**
+ * \brief Runs a generic sequence of actions to use an algorithm.
+ *
+ * This function runs the basic steps that are needed to run an FMM algorithm
+ * over a set of particles. It does the following steps :
+ *
+ *    - Load a tree using the class defined as a TreeLoader
+ *    - Prepares the needed kernels using the KernelLoader
+ *    - Prepares and runs the algorithm using the AlgorithmLoader
+ *
+ * See documentation of FTreeLoader, FKernelLoader, FAlgoLoader.
+ */
+template <class TreeLoader,
+          template <typename TL> class KernelLoader,
+          template <typename TL, template <typename TL> class KL> class AlgoLoader>
+void runperf(FPerfTestParams& params)
+{
+    TreeLoader treeLoader(params);
+    KernelLoader<TreeLoader> kernelLoader(params, treeLoader);
+    AlgoLoader<TreeLoader, KernelLoader> algoLoader(params, treeLoader, kernelLoader);
+    algoLoader.run();
+}
+
+int main (int argc, char** argv)
+{
+    FHelpDescribeAndExit(argc, argv,
+                         "Driver for Chebyshev interpolation kernel  (1/r kernel).",
+                         FParameterDefinitions::InputFile,
+                         FParameterDefinitions::OctreeHeight,
+                         FParameterDefinitions::OctreeSubHeight,
+                         FParameterDefinitions::NbThreads,
+                         {{"--algo"},"Algorithm to run (costzones, basic-static,"
+                                         " basic-dynamic, task)"});
+    FPerfTestParams params;
+    {
+        using namespace FParameterDefinitions;
+        using namespace FParameters;
+        params.filename = getStr(argc,argv,InputFile.options,
+                                 "../Data/unitCubeXYZQ100.bfma");
+        params.treeHeight = getValue(argc, argv, OctreeHeight.options, 5);
+        params.subTreeHeight = getValue(argc, argv, OctreeSubHeight.options, 2);
+        params.nbThreads = getValue(argc, argv, NbThreads.options, 1);
+        params.algo = getStr(argc,argv,{"--algo"},"task");
+    }
+
+
+    if( "basic-dynamic" == params.algo ) {
+        runperf<TreeLoaderFCheb<>, KernelLoaderFChebSym, AlgoLoaderThread>(params);
+    } else if( "task" == params.algo ) {
+        runperf<TreeLoaderFCheb<>, KernelLoaderFChebSym, AlgoLoaderTask>(params);
+    } else if ( "costzones" == params.algo ) {
+        runperf<TreeLoaderFCheb<>, KernelLoaderFChebSym, AlgoLoaderCostZones>(params);
+    }
+    
+
+}
diff --git a/Tests/noDist/PerfTestUtils.hpp b/Tests/noDist/PerfTestUtils.hpp
new file mode 100644
index 0000000000000000000000000000000000000000..a3b632250af3e3167abf26c32bb2265bcea711d6
--- /dev/null
+++ b/Tests/noDist/PerfTestUtils.hpp
@@ -0,0 +1,175 @@
+#ifndef _PERFTESTUTILS_HPP_
+#define _PERFTESTUTILS_HPP_
+
+#include <string>
+
+#include "Utils/FTic.hpp"
+#include "Files/FFmaGenericLoader.hpp"
+
+#include "Containers/FOctree.hpp"
+
+/**
+ * \brief Store the PerfTest program parameters.
+ */
+struct FPerfTestParams {
+    int subTreeHeight = 2; ///< Subtree height.
+    int treeHeight = 5;    ///< Tree height.
+    int nbThreads = 1;     ///< Maximum number of threads (when used). 
+    std::string filename = ""; ///< Particles file.
+    std::string algo = "task"; ///< Algorithm to run.
+};
+
+
+/**
+ * \brief Base class for tree loaders.
+ *
+ * This class itself does not provide anything but a base on which to build tree
+ * loaders. A tree loader should satisfy the following rules.
+ *
+ *    - Define the public typedefs : CellClass, ContainerClass, LeafClass,
+ *      OctreeClass.
+ *    - Provide public acces to a member of type OctreeClass _tree as the tree
+ *      that is loaded.
+ *    - Tree loading must happen at construction.
+ *    - It may provide any other members or typdefs required by a special
+ *      FKernelLoader or FAlgoLoader.
+ *
+ * For convenience, this class provides a timer and a basic loadTree method that
+ * should be enough to load a tree from and FMA file.
+ *
+ * \note It is not mandatory that a loader inherit from this class. It must
+ * however follow the aforementioned rules.
+ */
+class FTreeLoader {
+public:
+    /// A timer
+    /** Is used to time the loadTree method.
+     */
+    FTic time;
+protected:
+
+    /**
+     * \brief A template which type is always false.
+     *
+     * This template is only expanded by the compiler when it is requested
+     * (ie. the compiler will not try to optimize out its value.). Must be used
+     * to create false static_assert to catch unintended use of a template.
+     */
+    template <typename... Args>
+    struct false_type {
+        bool value = false;
+    };  
+
+    /**
+     * \brief Failure method for unimplemented loadTree templates.
+     *
+     * This template will catch unspecialised call to the loadTree method and
+     * will cause the compilation to fail with a (somewhat) meaningfull message.
+     */
+    template <typename...Args>
+    void loadTree(Args...) {
+        static_assert(false_type<Args...>::value,
+                      "I don't know how to load this tree with this loader...");
+    }
+
+
+    /**
+     * \brief Simple method to load a tree from a FMA file.
+     *
+     * The template parameters are usualy guessed by the compiler.
+     *
+     * \tparam OctreeClass The class of the tree to fill.
+     * \tparam FReal The floating point type.
+     *
+     * \param loader The file loader to read from the file.
+     * \param tree The tree to be filled.
+     */
+    template <class OctreeClass, typename FReal>
+    void loadTree(FFmaGenericLoader<FReal>& loader, OctreeClass& tree) {
+        std::cout << "Creating & inserting particles" << std::flush;
+
+        time.tic();
+
+        FPoint<FReal> position;
+        FReal physicalValue = 0.0;
+        for(FSize idxPart = 0 ; idxPart < loader.getNumberOfParticles() ; ++idxPart) {
+            // Read particle per particle from file
+            loader.fillParticle(&position,&physicalValue);
+            // put particle in octree
+            tree.insert(position, idxPart, physicalValue);
+        }
+
+        time.tac();
+        std::cout << " Done  (" << time.elapsed() << "s)." << std::endl;
+    }
+
+};
+
+/**
+ * \brief Base class for kernel loaders.
+ *
+ * This class itself does not provide anything but a base on which to build
+ * kernel loaders. A kernel loader should satisfy the following rules.
+ *
+ *    - Define the public typedefs : TreeLoader, KernelClass.
+ *    - Provide public acces to a member of type Kernelclass _kernel as the
+ *      kernel that is loaded.
+ *    - Kernel loading must happen at construction.
+ *    - It may provide any other members or typdefs required by a special
+ *      FAlgoLoader.
+ *
+ * For convenience, this class provides a timer.
+ *
+ * \tparam _TreeLoader The tree loader that was used.
+ *
+ * \note It is not mandatory that a loader inherit from this class. It must
+ * however follow the aforementioned rules.
+ */
+template<class _TreeLoader>
+class FKernelLoader {
+    /// The tree loader that was used (see FTreeLoader).
+    using TreeLoader = _TreeLoader;
+public:
+        FTic time;
+};
+
+/**
+ * \brief Base class for algorithm loaders.
+ *
+ * This class itself does not provide anything but a base on which to build
+ * algorithm loaders. A kernel loader should satisfy the following rules.
+ *
+ *    - Define the public typedefs : TreeLoader, KernelLoader.
+ *    - Provide public acces to a member of type
+ *      \link TreeLoader Treeloader::OctreeClass* \endlink` _algo`
+ *      as the algorithm that is loaded. This pointer should be valid from the
+ *      end of the ::run method to the destruction of the loader.
+ *    - It may provide any other members or typdefs.
+ *
+ * For convenience, this class provides a timer.
+ *
+ * \tparam _TreeLoader The tree loader that was used.
+ * \tparam _KernelLoader The kernel loader *template* that was used, the
+ *         KernelLoader type will then be _KernelLoader<_TreeLoader>.
+ *
+ * \note It is not mandatory that a loader inherit from this class. It must
+ * however follow the aforementioned rules.
+ */
+template <class _TreeLoader, template<typename> class _KernelLoader>
+class FAlgoLoader {
+    /// The tree loader that was used (see FTreeLoader).
+    using TreeLoader = _TreeLoader;
+    /// The kernel loader that was used (see FKernelLoader).
+    using KernelLoader = _KernelLoader<TreeLoader>;
+public:
+    /// A timer.
+    FTic time;
+    /// Method that runs the algorithm.
+    virtual void run() = 0;
+};
+
+
+
+
+
+#endif
diff --git a/Tests/noDist/TreeLoaderFCheb.hpp b/Tests/noDist/TreeLoaderFCheb.hpp
new file mode 100644
index 0000000000000000000000000000000000000000..6b1a8a9011251057673ee5dd55288fc0dca1d5bb
--- /dev/null
+++ b/Tests/noDist/TreeLoaderFCheb.hpp
@@ -0,0 +1,48 @@
+#ifndef _TREELOADERFCHEB_HPP_
+#define _TREELOADERFCHEB_HPP_
+
+#include "PerfTestUtils.hpp"
+
+#include "Kernels/Chebyshev/FChebCell.hpp"
+#include "Containers/FOctree.hpp"
+#include "Components/FSimpleLeaf.hpp"
+#include "Kernels/P2P/FP2PParticleContainerIndexed.hpp"
+
+#include "BalanceTree/FCostCell.hpp"
+
+/**
+ * \brief Tree loader for a Chebyshev cell type tree.
+ *
+ * See FTreeLoader documentation.
+ */
+template <typename _FReal = double>
+class TreeLoaderFCheb : public FTreeLoader {
+public:
+    using FReal = _FReal;
+    enum {ORDER = 7}; ///< Chebyshev interpolation order.
+
+    // Required type definitions.
+    using CellClass          = FCostCell<FChebCell<FReal, ORDER>>;
+    using ContainerClass     = FP2PParticleContainerIndexed<FReal>;
+    using LeafClass          = FSimpleLeaf<FReal, ContainerClass >;
+    using OctreeClass        = FOctree<FReal, CellClass, ContainerClass, LeafClass>;
+
+    /// File loader.
+    FFmaGenericLoader<FReal> _loader;
+    /// Required tree member.
+    OctreeClass _tree;
+
+    /// Constructs the loader and loads the tree.
+    TreeLoaderFCheb(FPerfTestParams& params):
+        _loader(params.filename),
+        _tree(params.treeHeight,
+              params.subTreeHeight,
+              _loader.getBoxWidth(),
+              _loader.getCenterOfBox()) {
+        this->loadTree(_loader, _tree);
+    }
+
+
+};
+
+#endif