Mercurial > hg > sv-dependency-builds
view src/fftw-3.3.3/doc/html/Allocating-aligned-memory-in-Fortran.html @ 10:37bf6b4a2645
Add FFTW3
author | Chris Cannam |
---|---|
date | Wed, 20 Mar 2013 15:35:50 +0000 |
parents | |
children |
line wrap: on
line source
<html lang="en"> <head> <title>Allocating aligned memory in Fortran - FFTW 3.3.3</title> <meta http-equiv="Content-Type" content="text/html"> <meta name="description" content="FFTW 3.3.3"> <meta name="generator" content="makeinfo 4.13"> <link title="Top" rel="start" href="index.html#Top"> <link rel="up" href="Calling-FFTW-from-Modern-Fortran.html#Calling-FFTW-from-Modern-Fortran" title="Calling FFTW from Modern Fortran"> <link rel="prev" href="Plan-execution-in-Fortran.html#Plan-execution-in-Fortran" title="Plan execution in Fortran"> <link rel="next" href="Accessing-the-wisdom-API-from-Fortran.html#Accessing-the-wisdom-API-from-Fortran" title="Accessing the wisdom API from Fortran"> <link href="http://www.gnu.org/software/texinfo/" rel="generator-home" title="Texinfo Homepage"> <!-- This manual is for FFTW (version 3.3.3, 25 November 2012). Copyright (C) 2003 Matteo Frigo. Copyright (C) 2003 Massachusetts Institute of Technology. Permission is granted to make and distribute verbatim copies of this manual provided the copyright notice and this permission notice are preserved on all copies. Permission is granted to copy and distribute modified versions of this manual under the conditions for verbatim copying, provided that the entire resulting derived work is distributed under the terms of a permission notice identical to this one. Permission is granted to copy and distribute translations of this manual into another language, under the above conditions for modified versions, except that this permission notice may be stated in a translation approved by the Free Software Foundation. --> <meta http-equiv="Content-Style-Type" content="text/css"> <style type="text/css"><!-- pre.display { font-family:inherit } pre.format { font-family:inherit } pre.smalldisplay { font-family:inherit; font-size:smaller } pre.smallformat { font-family:inherit; font-size:smaller } pre.smallexample { font-size:smaller } pre.smalllisp { font-size:smaller } span.sc { font-variant:small-caps } span.roman { font-family:serif; font-weight:normal; } span.sansserif { font-family:sans-serif; font-weight:normal; } --></style> </head> <body> <div class="node"> <a name="Allocating-aligned-memory-in-Fortran"></a> <p> Next: <a rel="next" accesskey="n" href="Accessing-the-wisdom-API-from-Fortran.html#Accessing-the-wisdom-API-from-Fortran">Accessing the wisdom API from Fortran</a>, Previous: <a rel="previous" accesskey="p" href="Plan-execution-in-Fortran.html#Plan-execution-in-Fortran">Plan execution in Fortran</a>, Up: <a rel="up" accesskey="u" href="Calling-FFTW-from-Modern-Fortran.html#Calling-FFTW-from-Modern-Fortran">Calling FFTW from Modern Fortran</a> <hr> </div> <h3 class="section">7.5 Allocating aligned memory in Fortran</h3> <p><a name="index-alignment-560"></a><a name="index-fftw_005falloc_005freal-561"></a><a name="index-fftw_005falloc_005fcomplex-562"></a>In order to obtain maximum performance in FFTW, you should store your data in arrays that have been specially aligned in memory (see <a href="SIMD-alignment-and-fftw_005fmalloc.html#SIMD-alignment-and-fftw_005fmalloc">SIMD alignment and fftw_malloc</a>). Enforcing alignment also permits you to safely use the new-array execute functions (see <a href="New_002darray-Execute-Functions.html#New_002darray-Execute-Functions">New-array Execute Functions</a>) to apply a given plan to more than one pair of in/out arrays. Unfortunately, standard Fortran arrays do <em>not</em> provide any alignment guarantees. The <em>only</em> way to allocate aligned memory in standard Fortran is to allocate it with an external C function, like the <code>fftw_alloc_real</code> and <code>fftw_alloc_complex</code> functions. Fortunately, Fortran 2003 provides a simple way to associate such allocated memory with a standard Fortran array pointer that you can then use normally. <p>We therefore recommend allocating all your input/output arrays using the following technique: <ol type=1 start=1> <li>Declare a <code>pointer</code>, <code>arr</code>, to your array of the desired type and dimensions. For example, <code>real(C_DOUBLE), pointer :: a(:,:)</code> for a 2d real array, or <code>complex(C_DOUBLE_COMPLEX), pointer :: a(:,:,:)</code> for a 3d complex array. <li>The number of elements to allocate must be an <code>integer(C_SIZE_T)</code>. You can either declare a variable of this type, e.g. <code>integer(C_SIZE_T) :: sz</code>, to store the number of elements to allocate, or you can use the <code>int(..., C_SIZE_T)</code> intrinsic function. e.g. set <code>sz = L * M * N</code> or use <code>int(L * M * N, C_SIZE_T)</code> for an L × M × N array. <li>Declare a <code>type(C_PTR) :: p</code> to hold the return value from FFTW's allocation routine. Set <code>p = fftw_alloc_real(sz)</code> for a real array, or <code>p = fftw_alloc_complex(sz)</code> for a complex array. <li><a name="index-c_005ff_005fpointer-563"></a>Associate your pointer <code>arr</code> with the allocated memory <code>p</code> using the standard <code>c_f_pointer</code> subroutine: <code>call c_f_pointer(p, arr, [...dimensions...])</code>, where <code>[...dimensions...])</code> are an array of the dimensions of the array (in the usual Fortran order). e.g. <code>call c_f_pointer(p, arr, [L,M,N])</code> for an L × M × N array. (Alternatively, you can omit the dimensions argument if you specified the shape explicitly when declaring <code>arr</code>.) You can now use <code>arr</code> as a usual multidimensional array. <li>When you are done using the array, deallocate the memory by <code>call fftw_free(p)</code> on <code>p</code>. </ol> <p>For example, here is how we would allocate an L × M 2d real array: <pre class="example"> real(C_DOUBLE), pointer :: arr(:,:) type(C_PTR) :: p p = fftw_alloc_real(int(L * M, C_SIZE_T)) call c_f_pointer(p, arr, [L,M]) <em>...use arr and arr(i,j) as usual...</em> call fftw_free(p) </pre> <p>and here is an L × M × N 3d complex array: <pre class="example"> complex(C_DOUBLE_COMPLEX), pointer :: arr(:,:,:) type(C_PTR) :: p p = fftw_alloc_complex(int(L * M * N, C_SIZE_T)) call c_f_pointer(p, arr, [L,M,N]) <em>...use arr and arr(i,j,k) as usual...</em> call fftw_free(p) </pre> <p>See <a href="Reversing-array-dimensions.html#Reversing-array-dimensions">Reversing array dimensions</a> for an example allocating a single array and associating both real and complex array pointers with it, for in-place real-to-complex transforms. <!-- --> </body></html>