diff --git a/cmake_support/FindFFTW.cmake b/cmake_support/FindFFTW.cmake
index 008ba3d93c94808427b14b7106e26c50243d81fd..b1be94c3efaf55600b8ecdb519a85cb57923f520 100644
--- a/cmake_support/FindFFTW.cmake
+++ b/cmake_support/FindFFTW.cmake
@@ -15,11 +15,21 @@ else (FFTW_INCLUDE_PATH)
   find_path (FFTW_INCLUDE_PATH fftw3.h)
   if (_DOUBLE_PREC)
     find_library (FFTW_LIBRARIES NAMES fftw3)
+    find_library (FFTW_THREADS_LIB NAMES fftw3_threads)
   elseif(NOT _DOUBLE_PREC)
     find_library (FFTW_LIBRARIES NAMES fftw3f)
+    find_library (FFTW_THREADS_LIB NAMES fftw3f_threads)
   endif()
+  if(FFTW_THREADS_LIB)
+    SET(FFTW_USE_THREADS TRUE)
+    SET(FFTW_LIBRARIES ${FFTW_LIBRARIES} ${FFTW_THREADS_LIB} )
+  else(FFTW_THREADS_LIB)
+  endif(FFTW_THREADS_LIB)
 endif (FFTW_INCLUDE_PATH)
 
+
+
+
 # handle the QUIETLY and REQUIRED arguments and set FFTW_FOUND to TRUE if 
 # all listed variables are TRUE
 include (FindPackageHandleStandardArgs)
diff --git a/modules/config/CMakeLists.txt b/modules/config/CMakeLists.txt
index 12abad743846debb764683f191c9c363a6a5ba2a..e8f153cff4742e75e969ffc25ab6223c8ed95f4e 100644
--- a/modules/config/CMakeLists.txt
+++ b/modules/config/CMakeLists.txt
@@ -43,6 +43,11 @@ if (ENABLE_SPNAV)
 else()
   set(spnav_enabled 0)
 endif()
+if (FFTW_USE_THREADS)
+  set(fftw_use_threads 1)
+else()
+  set(fftw_use_threads 0)
+endif()
 
 set(config_hh_generator "CMake")
 set(CONFIG_HH_FILE "${CMAKE_CURRENT_SOURCE_DIR}/config.hh")
diff --git a/modules/config/config.hh.in b/modules/config/config.hh.in
index 4d6fc6ff1fd4d58d400161820eece3d9078ac572..20ac7eab47dfa8685e9700f8d60edbf7df84b1e9 100644
--- a/modules/config/config.hh.in
+++ b/modules/config/config.hh.in
@@ -30,4 +30,5 @@
 #define OST_DOUBLE_PRECISION @double_prec@
 #define OST_STATIC_PROPERTY_WORKAROUND @static_props@
 #define OST_SPNAV_ENABLED @spnav_enabled@
+#define OST_FFT_USE_THREADS @fftw_use_threads@
 #endif
diff --git a/modules/img/alg/src/CMakeLists.txt b/modules/img/alg/src/CMakeLists.txt
index da982b40c80144e34b7bb5d34f8206bdf959a309..804cd683ec0bd9fa08fe6e95b72510e92e77e4f7 100644
--- a/modules/img/alg/src/CMakeLists.txt
+++ b/modules/img/alg/src/CMakeLists.txt
@@ -93,9 +93,9 @@ line_average.hh
 rscrosscorr.hh
 )
 
-
+include(${QT_USE_FILE})
 module(NAME img_alg SOURCES "${OST_IMG_ALG_SOURCES}" 
        HEADERS "${OST_IMG_ALG_HEADERS}" 
        HEADER_OUTPUT_DIR ost/img/alg
        DEPENDS_ON img
-       LINK ${FFTW_LIBRARIES})
+       LINK ${FFTW_LIBRARIES} ${QT_QTCORE_LIBRARY})
diff --git a/modules/img/alg/src/fft.cc b/modules/img/alg/src/fft.cc
index f38bd69f0c0d753de11c615a4b1d80c5e7d76fb2..de8480e55d150426655187d8129437c404c89b45 100644
--- a/modules/img/alg/src/fft.cc
+++ b/modules/img/alg/src/fft.cc
@@ -25,6 +25,7 @@
 #include <boost/shared_ptr.hpp>
 
 #include <fftw3.h>
+#include <QThread>
 
 #include <ost/message.hh>
 #include <ost/img/value_util.hh>
@@ -54,8 +55,14 @@ int half_plus_one(int n) {
 
 } // anon ns
 
-FFTFnc::FFTFnc(): ori_flag_(false) {}
-FFTFnc::FFTFnc(bool f): ori_flag_(f) {}
+FFTFnc::FFTFnc(): ori_flag_(false)
+{
+  OST_FFTW_fftw_init_threads();
+}
+FFTFnc::FFTFnc(bool f): ori_flag_(f)
+{
+  OST_FFTW_fftw_init_threads();
+}
 
 
 // real spatial -> complex half-frequency
@@ -119,7 +126,7 @@ reinterpret_cast<OST_FFTW_fftw_complex*>(out_state->Data().GetData());
   for(size_t i=0;i<block_count;i++) {
     std::copy(&in_ptr[i*src_size],&in_ptr[(i+1)*src_size],&fftw_in[i*2*dst_size]);
   }
-
+  OST_FFTW_fftw_plan_with_nthreads(std::max<int>(1,QThread::idealThreadCount()));
   OST_FFTW_fftw_plan plan = OST_FFTW_fftw_plan_dft_r2c(rank,n,
 				     fftw_in,fftw_out,
 				     FFTW_ESTIMATE);
@@ -202,6 +209,7 @@ ImageStateBasePtr FFTFnc::VisitState<Complex,HalfFrequencyDomain>(const ComplexH
   Real* fftw_out = reinterpret_cast<Real*>(out_ptr);
   
   assert(sizeof(OST_FFTW_fftw_complex)==sizeof(Complex));
+  OST_FFTW_fftw_plan_with_nthreads(std::max<int>(1,QThread::idealThreadCount()));
   OST_FFTW_fftw_complex* fftw_in =
 reinterpret_cast<OST_FFTW_fftw_complex*>(out_ptr);
 
@@ -247,6 +255,7 @@ ImageStateBasePtr FFTFnc::VisitState<Complex,SpatialDomain>(const ComplexSpatial
 reinterpret_cast<OST_FFTW_fftw_complex*>(out_state->Data().GetData());
 
   // in place transform
+  OST_FFTW_fftw_plan_with_nthreads(std::max<int>(1,QThread::idealThreadCount()));
   OST_FFTW_fftw_plan plan = OST_FFTW_fftw_plan_dft(rank,n,
 				 fftw_out, fftw_out, 
 				 dir, 
@@ -284,6 +293,7 @@ ImageStateBasePtr FFTFnc::VisitState<Complex,FrequencyDomain>(const ComplexFrequ
 reinterpret_cast<OST_FFTW_fftw_complex*>(out_state->Data().GetData());
 
   // in place transform
+  OST_FFTW_fftw_plan_with_nthreads(std::max<int>(1,QThread::idealThreadCount()));
   OST_FFTW_fftw_plan plan = OST_FFTW_fftw_plan_dft(rank,n,
 				 fftw_out, fftw_out, 
 				 dir, 
diff --git a/modules/img/alg/src/fftw_helper.hh b/modules/img/alg/src/fftw_helper.hh
index ab130188bc0d8a004cd32b4b73a5df6712fd17c0..081ff7278e975602bb8c4ce80131e6b18f5a38a4 100644
--- a/modules/img/alg/src/fftw_helper.hh
+++ b/modules/img/alg/src/fftw_helper.hh
@@ -68,4 +68,19 @@ Author: Juergen Haas
 #define OST_FFTW_fftw_plan_dft_c2r fftwf_plan_dft_c2r
 #define OST_FFTW_fftw_plan_many_dft fftwf_plan_many_dft
 #endif
+
+#if OST_FFT_USE_THREADS
+  #if OST_DOUBLE_PRECISION
+    #define OST_FFTW_fftw_init_threads fftw_init_threads
+    #define OST_FFTW_fftw_plan_with_nthreads fftw_plan_with_nthreads
+  #else
+    #define OST_FFTW_fftw_init_threads fftwf_init_threads
+    #define OST_FFTW_fftw_plan_with_nthreads fftwf_plan_with_nthreads
+  #endif
+#else
+  void fftw_noop(unsigned int i=0){}
+  #define OST_FFTW_fftw_init_threads fftw_noop
+  #define OST_FFTW_fftw_plan_with_nthreads fftw_noop
+#endif
+
 #endif
diff --git a/modules/img/alg/src/line_iterator.cc b/modules/img/alg/src/line_iterator.cc
index f303f4b28c5e0be3e341d694264ca5200b413a89..2cd81e5a88584816e8d005b325498d47aee96377 100644
--- a/modules/img/alg/src/line_iterator.cc
+++ b/modules/img/alg/src/line_iterator.cc
@@ -155,7 +155,7 @@ LineIterator LineIterator::operator++(int)
 {
   LineIterator tmp(*this);
   this->operator++();
-  return tmp;;
+  return tmp;
 }
 
 ExtentIterator LineIterator::GetNextLine()