Chris@10: Chris@10: Chris@10: Usage of Multi-threaded FFTW - FFTW 3.3.3 Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10: Chris@10:
Chris@10: Chris@10: Chris@10:

Chris@10: Next: , Chris@10: Previous: Installation and Supported Hardware/Software, Chris@10: Up: Multi-threaded FFTW Chris@10:


Chris@10:
Chris@10: Chris@10:

5.2 Usage of Multi-threaded FFTW

Chris@10: Chris@10:

Here, it is assumed that the reader is already familiar with the usage Chris@10: of the uniprocessor FFTW routines, described elsewhere in this manual. Chris@10: We only describe what one has to change in order to use the Chris@10: multi-threaded routines. Chris@10: Chris@10:

First, programs using the parallel complex transforms should be linked Chris@10: with -lfftw3_threads -lfftw3 -lm on Unix, or -lfftw3_omp Chris@10: -lfftw3 -lm if you compiled with OpenMP. You will also need to link Chris@10: with whatever library is responsible for threads on your system Chris@10: (e.g. -lpthread on GNU/Linux) or include whatever compiler flag Chris@10: enables OpenMP (e.g. -fopenmp with gcc). Chris@10: Chris@10: Chris@10:

Second, before calling any FFTW routines, you should call the Chris@10: function: Chris@10: Chris@10:

     int fftw_init_threads(void);
Chris@10: 
Chris@10:

Chris@10: This function, which need only be called once, performs any one-time Chris@10: initialization required to use threads on your system. It returns zero Chris@10: if there was some error (which should not happen under normal Chris@10: circumstances) and a non-zero value otherwise. Chris@10: Chris@10:

Third, before creating a plan that you want to parallelize, you should Chris@10: call: Chris@10: Chris@10:

     void fftw_plan_with_nthreads(int nthreads);
Chris@10: 
Chris@10:

Chris@10: The nthreads argument indicates the number of threads you want Chris@10: FFTW to use (or actually, the maximum number). All plans subsequently Chris@10: created with any planner routine will use that many threads. You can Chris@10: call fftw_plan_with_nthreads, create some plans, call Chris@10: fftw_plan_with_nthreads again with a different argument, and Chris@10: create some more plans for a new number of threads. Plans already created Chris@10: before a call to fftw_plan_with_nthreads are unaffected. If you Chris@10: pass an nthreads argument of 1 (the default), threads are Chris@10: disabled for subsequent plans. Chris@10: Chris@10:

With OpenMP, to configure FFTW to use all of the currently running Chris@10: OpenMP threads (set by omp_set_num_threads(nthreads) or by the Chris@10: OMP_NUM_THREADS environment variable), you can do: Chris@10: fftw_plan_with_nthreads(omp_get_max_threads()). (The ‘omp_’ Chris@10: OpenMP functions are declared via #include <omp.h>.) Chris@10: Chris@10:

Given a plan, you then execute it as usual with Chris@10: fftw_execute(plan), and the execution will use the number of Chris@10: threads specified when the plan was created. When done, you destroy Chris@10: it as usual with fftw_destroy_plan. As described in Chris@10: Thread safety, plan execution is thread-safe, but plan Chris@10: creation and destruction are not: you should create/destroy Chris@10: plans only from a single thread, but can safely execute multiple plans Chris@10: in parallel. Chris@10: Chris@10:

There is one additional routine: if you want to get rid of all memory Chris@10: and other resources allocated internally by FFTW, you can call: Chris@10: Chris@10:

     void fftw_cleanup_threads(void);
Chris@10: 
Chris@10:

Chris@10: which is much like the fftw_cleanup() function except that it Chris@10: also gets rid of threads-related data. You must not execute any Chris@10: previously created plans after calling this function. Chris@10: Chris@10:

We should also mention one other restriction: if you save wisdom from a Chris@10: program using the multi-threaded FFTW, that wisdom cannot be used Chris@10: by a program using only the single-threaded FFTW (i.e. not calling Chris@10: fftw_init_threads). See Words of Wisdom-Saving Plans. Chris@10: Chris@10: Chris@10: Chris@10: