/** \mainpage High performance Spherical Harmonic Transform for numerical simulations. * * SHTns is a high performance library for Spherical Harmonic Transform written in C, aimed at numerical simulation (fluid flows, mhd, ...) in spherical geometries.\n * * Main features : * - \link bench blazingly fast\endlink * - distributed under the open source \link license CeCILL License\endlink (GPL compatible) * - both scalar and \link vsh vector transforms\endlink * - backward and forward (synthesis and analysis) functions * - flexible truncation * - spatial data can be stored in latitude-major or longitude-major order arrays. * - various conventions (normalization and Condon-Shortley phase) * - can be used from Fortran, c/c++, and Python programs * - a highly efficient Gauss algorithm working with Gauss nodes (based on Gauss-Legendre quadrature) * - an algorithm using DCT for regular nodes (based on a generalized Fejer quadrature) * - support for SSE2, SSE3 and AVX vectorization with gcc. * - parallel transforms with OpenMP (for Gauss grid only). * - synthesis (inverse transform) at any coordinate (not constrained to a grid) useful for rendering purposes. * - ability to choose the optimal spatial sizes for a given spherical harmonic truncation. * - on-the-fly transforms : saving memory and bandwidth, and can even be faster depending on architecture. * - accurate up to spherical harmonic degree l=16383 (at least). * - \link rotation rotation functions\endlink to rotate spherical harmonics (beta). * - \link operators special spectral operator\endlink functions that do not require a transform (multiply by cos(theta)...) . * - scalar transforms for complex spatial fields. * * Using DCT acceleration and other classical \link opt optimizations\endlink, it is intended to be \link bench very fast\endlink. * * It requires the FFTW library for Fast Fourier Transforms. * * If you use SHTns for research work, please cite the paper: Efficient Spherical Harmonic Transforms aimed at pseudo-spectral numerical simulations, * also available from arXiv. * * If you accept the open source \link license CeCILL License\endlink (GPL compatible french License), you can download SHTns. * * Please report bugs and feature request on the issue tracker. * * \see shtns.h for the definitions of variables, macros and functions. * \see The example programs (in \link SHT_example.f Fortran \endlink, \link SHT_example.c C \endlink and \link SHT_example.py Python \endlink) to get started. * \see The organisation of data used by SHTns is described in \ref spat. * \see The description of \ref opt. * * \author SHTns is written by Nathanael Schaeffer (CNRS). Email: nathanael.schaeffer@ujf-grenoble.fr * * \image html logo_cnrs_small.png * */ /** \page compil Compiling and installing SHTns To compile SHTns, you need a C compiler and the FFTW library version 3.0 (or better 3.3) or later, carefully optimized and installed. On modern x86 processors, we recommend using GCC 4.0 or later as it produces the fastest vectorized code. First run \code ./configure --help \endcode to see the available configure options. Then, for example run \code ./configure --prefix=$HOME --enable-openmp \endcode if you want to install the multi-threaded transform library to $HOME/lib and the headers to $HOME/include. To install the Python extension, see \ref python. The configure script will detect your FFTW library capabilities and other system dependant features. You may optionally adapt the resulting \c Makefile and \c sht_config.h to your architecture and compiler. Then compile with \c make and install with \c make \c install. SHTns requires FFTW, which can be installed as a compiled package on many operating systems. However, we recommend to compile it and install it yourself, as it can lead to much better performance. If FFTW has not been compiled and tuned carefully for your machine, you may get bad performance with SHTns (down to 50% slower !). Furthermore, if FFTW is more recent than 3.3, or has not been configured with \c --enable-openmp (see below), SHTns cannot do the Fourier transform part in parallel. \section fftw Compiling and installing FFTW You can download the latest FFTW package from http://www.fftw.org, compile it and install it. In order to use the parallel transforms, make sure to pass \c --enable-openmp to the configure script. We also recommend passing \c --enable-sse2 or \c --enable-avx to the configure script on x86 architectures. You may also take a look at the FFTW documentation for compilation advice. For example: \code ./configure --enable-openmp --enable-shared --enable-avx --prefix=$HOME/usr make make install \endcode That's it, FFTW is now installed and ready for use. \subsection mkl Using MKL instead of FFTW SHTns can use MKL instead of FFTW. To do this, pass \c --enable-mkl to the configure script. Performance can be slightly better, but there is no DCT support, so only gauss grid is available. \note Currently the call of MKL through the fftw wrappers is not thread safe. The user is responsible for working around this MKL problem as indicated by intel. If you want to call SHT function from multiple threads, you must tell MKL how many threads you will use BEFORE initializing shtns : \code #include "fftw3_mkl.h" fftw3_mkl.number_of_user_threads = 4; sht = shtns_create(...); ... \endcode \section makefile Preparing the Makefile Run \c ./configure in the SHTns directory. You can use \c --enable-openmp to enable multi-threaded transforms, and \c --enable-long-double to (maybe) increase accuracy during initialization (not recommended). You can then edit the resulting Makefile: \li set \c PREFIX= to the desired install path. \li further editing is needed if you don't use gcc \section config Changing default options Simple tunings can be done by passing options to \c ./configure (see \c ./configure \c --help) \section comp Compiling \li type \code make \endcode to compile SHTns. it will produce the library \c libshtns.a \li type \code make PREFIX=/usr/local install \endcode to install the library in your system (set PREFIX to the desired location, default is $HOME). \li type \code make docs \endcode to generate this documentation, placed in the doc/html/ subdirectory. \li type \code make time_SHT \endcode to build the test and timing program. */ /** \page spat Spatial data layouts and grids used by SHTns The angular coordinates on a spherical shell are the co-latitude \f$ \theta \f$ and the longitude \f$ \phi \f$. The spatial discretization is defined by a call to \ref shtns_init; or to \ref shtns_set_grid or \ref shtns_set_grid_auto after \ref shtns_create. A field \e A (or one of its component in case of a vector field) is discretized on an ordered grid, which consists of - NPHI equally spaced nodes in longitude, spanning the range of angle between 0 (included) and \f$2\pi\f$/MRES (excluded), - NLAT gauss nodes or equally spaced nodes (depending on \ref shtns_type) in latitude, spanning angles from 0 to \f$ \pi \f$. There are constraints on the sizes NLAT and NPHI. The sampling theorem requires NPHI > 2*MMAX, with MMAX the highest Fourier mode to be resolved. Similarly NLAT > LMAX is required for a Gauss-Legendre quadrature. When dealing with non-linear equation, or more generally when analysing data that is not band-limited to LMAX and/or MMAX, anti-aliasing constraints exist. For a given non-linear order of the equation N (typically 2), the anti-aliasing conditions are : - NPHI > (N+1)*MMAX - NLAT > (N+1)*LMAX/2 for Gauss-Legendre quadrature - NLAT > N*LMAX for our DCT quadrature Furthermore, FFT is more efficient for certain values of NPHI and/or NLAT. If you do not require specific NLAT or NPHI, you can let SHTns choose them for you, by using \ref shtns_set_grid_auto with nphi and/or nlat set to zero. Identifying a grid point with its indices \c it and \c ip in latitudinal and longitudinal direction respectively (indices start at 0), one has : - \f$ \theta \f$ = \c acos(ct[it]) where \c ct = \ref shtns_info.ct is initialized by the call to \ref shtns_init - \f$ \phi \f$ = \c ip * 2\f$\pi\f$/(NPHI*MRES) and \ref PHI_RAD and \ref PHI_DEG can be used to get the phi angle corresponding to \c ip. - How to access \f$ A(\theta,\phi) \f$ depends on the data layout, but it is always a contiguous array of double precision floating point values (64 bit) Currently, three data layouts are supported. \ref phi_fast, \ref theta_fast, and \ref native. \section phi_fast Contiguous longitudes In this layout, increasing longitudes are stored next to each other for each latitude. That is \f$ A(\theta,\phi) \f$ = \c A[it*NPHI + ip] in C or \c A(ip,it) in Fortran. Use \ref SHT_PHI_CONTIGUOUS to instruct \ref shtns_init to use this spatial data layout : \code shtns_init ( sht_gauss | SHT_PHI_CONTIGUOUS, ... ) \endcode will tell shtns to precompute everything for a gauss grid, and spatial data stored with longitude varying fastest. \section theta_fast Contiguous latitudes In this layout, increasing latitude are stored next to each other for each longitude. That is \f$ A(\theta,\phi) \f$ = \c A[ip*NLAT + it] in C or \c A(it,ip) in Fortran. Use \ref SHT_THETA_CONTIGUOUS to instruct \ref shtns_init to use this spatial data layout. \section native Native layout The native way of storing spatial field data (which will help you achieve the best performance with SHTns) is the same as the \ref theta_fast layout, except that it requires you to allocate slightly more memory for a field than the NLAT*NPHI double values. Namely NLAT*(NPHI/2+1)*2 doubles are required (instead of NLAT*NPHI) to be able tu use the in-place FFT of FFTW. In Fortran this means you will allocate data as if the phi direction had (NPHI/2+1)*2 points instead fo NPHI. The additional space, is located at the end of the useful data, and you don't need to worry about it. To instruct \ref shtns_init that your spatial data has been set up using this layout, use \ref SHT_NATIVE_LAYOUT. \note One must be careful when allocating spatial data, as the Fourier transform requires some additional space.\n Precisely, the number of required \c double is \ref NSPAT_ALLOC = NLAT*(NPHI/2+1)*2. The first NLAT*NPHI values are the spatial data. \note In addition, spatial data must be allocated using the \c fftw_malloc function, for example using : \code A = fftw_malloc( NSPAT_ALLOC * sizeof(double) );\endcode \note It is not recommended to use the native layout with Fortran or Python. */ /** \page spec Spherical Harmonics storage and normalization \section spec_data Spherical Harmonic coefficients layout The collection of Alm (spherical harmonics coefficient of degree \e l and order \e m) is stored in an array of \c complex \c double floating point values. The spherical harmonic truncation is defined by a call to \ref shtns_init or to \ref shtns_create which returns an \c shtns_cfg object (pointer to \ref shtns_info) from which the the size of the array, \ref shtns_info.nlm, can be read (see examples). The field \em A is decomposed on the basis of spherical harmonics Ylm (degree \e l, order \e m) : \f[ A(\theta,\phi) = \sum_{l,m} A_l^m Y_l^m(\theta,\phi)\f] The series is truncated at degree LMAX and order MMAX*MRES, and only order that are multiple of MRES are described. The NLM coefficients Alm are stored in a one-dimensional array, grouped by \e m.\n The details of the storage may depend on the algorithm, and to access a particular coefficient, you should always use the following macros : - \ref LM(l,m) gives the index in the array for coefficient of degree \e l and order \e m. - \ref LiM(l,im) gives the index in the array for coefficient of degree \e l and order \e im*MRES. Hence \c Alm[LM(l,m)] is the (complex) SH coefficient of degree \e l and order \e m. \note Since SHTns only deals with real spatial data, the Spherical Harmonic coefficients with negative \e m are related to the complex conjugate of the positive \e m coefficients : \f$ A_l^{-m} = (-1)^{m} (A_l^m)^*\f$. For efficiency purposes, only the coefficients with positive \e m are stored. \note There is no test on bounds, so that LM(l,m) with \e l or \e m that is not described, will give undefined result (in particular for an \e m that is not a multiple of \c MRES) To perform loops on all coefficients, use the macros \ref LM_LOOP and \ref LM_L_LOOP in combination of the arrays \ref shtns_info.li which give the degree \e l and functions of it for any array index. \section Allocation Allocation for a SH description is simply done with : \code Alm = fftw_malloc( NLM * sizeof(complex double) );\endcode The use of fftw_malloc ensures proper alignment of data for SSE vectorization (if fftw has been properly compiled with --enable-sse2 or --enable-avx). \section norm Normalization Several normalizations for the spherical harmonics exist (details on wikipedia). SHTns lets you choose one of the following : \ref shtns_init uses the default (defined by \ref SHT_DEFAULT_NORM) but you can choose another normalization by calling \ref shtns_create instead. The Legendre associated functions are defined by : \f[ P_l^m (x) = (-1)^m\ (1-x^2)^{m/2}\ \frac{d^m}{dx^m}P_l(x) \f] which includes the Condon-Shortley phase \f$ (-1)^m \f$ by default. For reference, function \ref SH_to_point gives a simple explicit implementation of the inverse SH transform (synthesis). \subsection norm_opt Normalization variations Use \ref shtns_create with \ref SHT_NO_CS_PHASE to disable the Condon-Shortley phase. For example, to initialize a Schmidt semi-normalized transform of maximum degree 16 without Condon-Shortrley phase, use \code shtns_create(16, 16, 1, sht_schmidt | SHT_NO_CS_PHASE); \endcode By default, SHTns uses "complex" spherical harmonic normalization. However, as it deals only with real functions, the m<0 coefficients are not stored, since they are complex conjugate of the m>0 ones. Thus, one must pay attention, that the m>0 coefficient have to be counted twice in a total energy calculation. An alternative to this "complex" normalization is the "real" normalization, for which the energy is simply the sum of the coefficients (m>=0) squared. Use \ref shtns_create with \ref SHT_REAL_NORM to use a "real" spherical harmonic normalization. This is the usual "real" spherical harmonics, if one takes the complex conjugate of the coefficients. A few useful examples, for orthonormal spherical harmonics : - a constant unit value on the sphere is represented by the coefficient \f$ c_0^0 = \sqrt{4\pi} \f$. - \f$ \cos \theta \f$ is represented by the coefficient \f$ c_1^0 = \sqrt{4\pi/3} \f$. - \f$ \sin \theta \, \cos \phi \f$ is represented by \f$ c_1^1 = -\sqrt{2\pi/3} \f$ (with Condon-Shortley phase) - with a "real" normalization instead, one would have \f$ c_1^1 = -\sqrt{4\pi/3} = -c_1^0 \f$. The functions \ref sh00_1(), \ref sh10_ct(), \ref sh11_st() \ref shlm_e1() will help to build simple spectral fields, no matter what normalization you choose. \subsection nrj Energy \note The energy of a scalar field can be retrieved from it's spherical harmonic representation as: \f[ \int ||q||^2 = \sum_{l,m \geq 0} C_m N_l \: |Q_l^m|^2 \f] where \f$ C_m \f$ and \f$ N_l \f$ are the following normalization factor: - if \ref SHT_REAL_NORM has been used, \f$ C_m = 1 \f$. Otherwise \f$ C_m = 2 - \delta_{m0} \f$ (ie \f$ C_m=1 \f$ for \f$ m=0 \f$ and \f$ C_m = 2 \f$ for \f$ m>0 \f$). - for \ref sht_orthonormal : \f$ N_l = 1 \f$ - for \ref sht_fourpi : \f$ N_l = 4\pi \f$ - for \ref sht_schmidt : \f$ N_l = 4\pi/(2l+1) \f$ \subsection Troubleshooting If you have problems in getting your spherical harmonic coefficient right, you may check the following : - Schmidt semi-normalization ? try multiplying or dividing by \f$ \sqrt{2l+1} \f$ - Condon-Shortley phase ? try multiplying by \f$ (-1)^m \f$ - Real or Complex ? try multiplying or dividing the m>0 coefficients by \f$ \sqrt{2} \f$ - Relative sign ? try to take the complex conjugate of the coefficients. - Does the array include l=0 ? */ /** \page vsh Vector Spherical Harmonics as implemented in SHTns SHTns also provides vector transforms. \section vsh_def Radial - Spheroidal - Toroidal decomposition The components of a vector field on a sphere \f$\mathbf{v}(\theta,\phi) = (v_r, v_\theta, v_\phi)\f$ in canonical spherical coordinates cannot be transformed as scalars, because they do not behave like scalar fields under rotation for example (actually only the radial component does). To overcome this, we use the following decomposition: \f[ \mathbf{v} = Q(\theta,\phi) \mathbf{e_r}+ r \nabla S(\theta,\phi) - \mathbf{r} \times \nabla T(\theta,\phi) \f] where Q, S and T are respectively the radial, spheroidal and toroidal scalar fields. We have explicitely: \f[ \mathbf{v} = \left( \begin{array}{rcl} Q & & \\ \partial_\theta S & + & \frac{1}{\sin\theta} \partial_\phi T \\ \frac{1}{\sin\theta} \partial_\phi S & - &\partial_\theta T \end{array} \right) \f] Expanding Q, S and T (which are scalar fields) in spherical harmonics gives: \f[ \mathbf{v} = \sum_{l,m} \left( Q_l^m \, Y_l^m(\theta,\phi) \mathbf{e_r}+ S_l^m \, r \nabla Y_l^m(\theta,\phi) - T_l^m \, \mathbf{r} \times \nabla Y_l^m(\theta,\phi) \right) \f] If this only involves a simple scalar transform for the radial component, it does involve derivatives of spherical harmonics for the tangential components. The functions \ref SHqst_to_spat and \ref spat_to_SHqst transform vector fields from their spherical coordinate spatial representation \c Vr, \c Vt, \c Vp to their spectral representation using \c Qlm, \c Slm and \c Tlm. The functions \ref SHsphtor_to_spat and \ref spat_to_SHsphtor do the same for a tangential vector on the sphere (without radial component). \see A closely related definition of vector spherical harmonics on wikipedia. \note The energy of the vector field can be retrieved from it's orthonormal spherical harmonic representation as: \f[ \int ||\mathbf{v}||^2 = \sum_{l,m \geq 0} C_m N_l \: \left( |Q_l^m|^2 \, + \, l(l+1)\,\left( |S_l^m|^2 + |T_l^m|^2 \right) \right) \f] where \f$ C_m \f$ and \f$ N_l \f$ are defined in the \ref nrj section. \section vsh_3d Special case of 3D divergence-free fields For a divergenceless 3D vector, the radial scalar Q and the spheroidal scalar S can be derived from the same poloidal scalar P : \f[ Q = \frac{l(l+1)}{r} P \f] \f[ S = \frac{1}{r} \frac{\partial \, rP}{\partial r} \f] which corresponds to the poloidal/toroidal decomposition that ensures the vector field to be divergence-free. \f[ \mathbf{v} = \nabla \times (T \mathbf{r}) + \nabla \times \nabla \times (P \mathbf{r}) \f] */ /** \page opt Optimizations implemented in SHTns SHTns is an implementation of the Spherical Harmonic Transform which aims at being accurate and fast, with direct numerical simulations in mind. As such, the following optimizations are implemented : \section opt_fft Use the Fast-Fourier Transform Any serious SHT should use the FFT, as it improves accuracy and speed. SHTns uses the FFTW library for the fast Fourier transform, a portable, flexible and blazingly fast FFT implementation. \section opt_mirror Take advantage of mirror symmetry This is also a classical optimization. Due to the defined symmetry of spherical harmonics with respect to a reflection about the equator, one can reduce by a factor 2 the operation count of both direct and reverse transforms. \section opt_polar Polar optimization A less common, but still classic optimization : high m's spherical harmonics have their magnitude decrease exponentially when approaching the poles. SHTns use a threshold below which the SH is taken to be zero. The default value for this threshold is defined by \ref SHT_DEFAULT_POLAR_OPT, but you can also choose it with the \c eps parameter of the \ref shtns_set_grid function, trading some accuracy for more speed. Around 5% to 15% speed increase are typical values for a SHT with a large MMAX. \section opt_cache Cache optimizations Cache optimizations have been carried out throughout the code. This means ordering coefficients and stripping out systematic zeros. \section opt_size Spatial Size optimization When using \ref shtns_set_grid_auto, you can let SHTns choose the optimal spatial sizes for the spherical harmonic truncation you specified with \ref shtns_create. The spatial sizes are chosen so that no aliasing occurs, and ensuring FFTW will perform as fast as possible. \section opt_vect Explicit Vectorization Most x86 CPUs have support for Single-Instruction-Multiple-Data (SIMD) operations in double precision, allowing to perform the same double precision arithmetic operations on a vector of 2 (SSE2) or 4 (AVX) double precision numbers, effectively multiplying by 2 or 4 the computing power. Modern compilers attempt to vectorize the code, but they cannot do miracles. SHTns uses explicit vectorization extensively (thanks to the GCC vector extensions) to use the full computing power of your CPU. \section opt_dct Discrete Cosine Transform optimization This is a non-trivial optimization mostly efficient on small m's, which can significantly decrease the operation count. It can be noticed that Spherical Harmonics are polynomials of order \f$ l \f$ in \f$ \cos \theta \f$ (multiplied by \f$ \sin \theta \f$ for odd m's). Thus, going to DCT space and summing there can help to further reduce the number of operations. Using a DCT requires a regular grid instead of the usual Gauss grid, and this means doubling the number of grid points to keep an equally accurate quadrature (Tchebychev quadrature), reducing the performance gain from the DCT to zero, when not performing worse. However, SHTns uses a quadrature inspired by the Fejer quadrature (or Clenshaw-Curtis), which allows exact result using only one more grid-point than the classical Gauss-Legendre quadrature. This allows SHTns to efficiently use the DCT optimization, but also to use regular grids with almost half the grid points as usual. The theoretical result is that for a function band-limited by L, L+2 equispaced grid-points allow us to fully recover the spectral content. There is, however, a drawback for non-linear terms : for a product of two functions band-limited to L, we would require L*2+2 grid points, whereas the Gauss-Legendre quadrature would only require (L+1)*3/2 grid-points if we are interested in only the first L modes of the product. There is a slight loss of accuracy due to round-off errors (as an example, relative errors with DCT are around \c 1.e-11 instead of \c 1.e-12 without for lmax around 500). If you can live with it (you really should), the DCT-based algorithm performs much faster on small m's : twice to four times faster for axisymmetric transforms, and between 1.5 and 1.2 times faster on non-axisymmetric ones. For large m's however, the spatial-domain polar optimization takes the lead. Note also that the DCT algorithm, while being usually faster, requires significantly larger initialization time. It is then suitable for performing many transforms (as in a direct numerical simulations). \section opt_runtime Runtime tuning of algorithm The default mode used by SHTns is to measure performance of the different algorithms, and choose the one that performs best (it will also check that the accuracy is good enough). However, there are situation where either the Gauss-Legendre algorithm or a regular grid is required, and you can choose to do so using the adequate \ref shtns_type when calling \ref shtns_init or \ref shtns_set_grid. */ /** \page bench Speed and Accuracy \section Speed SHTns does not implement any "fast" algorithm. However, timings with other Spherical Harmonic Transform tools (including a fast algorithm) show that SHTns performs much faster than any other. Furthermore, even at large sizes, the fast algorithm we tested does not seem to be willing to take the lead.
\f$\ell_{max}\f$ shtools 2.8 (Gauss) libpsht (1 thread) SpharmonicKit2 2.7 (fast) SHTns 2.1 (1 thread, Gauss) SpharmonicKit2/SHTns
63 1.14 ms 1.05 ms 1.1 ms 0.09 ms 12.2
127 3.5 ms 4.7 ms 5.5 ms 0.60 ms 9.2
255 28 ms 27 ms 21 ms 4.2 ms 5.0
511 200 ms 162 ms 110 ms 28 ms 3.9
1023 1.8 s 850 ms 600 ms 216 ms 2.8
2047 13.0 s 4.4 s NA (out of memory) 1.6 s NA
4095 NA (seg fault) 30.5 s 11.8 s NA
Average times for forward or backward scalar transform on an Intel Xeon X5650 (2.67GHz), with gcc 4.4.5 and "-O3 -march=native -ffast-math" compilation options.
\section par_speed Parallel speed SHTns has parallel algorithms since version 2.2. When compared to libpsht (parallelized with OpenMP too), SHTns is faster especially for relatively small sizes.
\f$\ell_{max}\f$ libpsht 20110131 SHTns 2.2.1 libpsht/SHTns
63 5.0 ms 0.05 ms 100
127 5.4 ms 0.22 ms 24.5
255 8.5 ms 1.4 ms 6.1
511 23.5 ms 6.5 ms 3.6
1023 125 ms 43 ms 2.9
2047 700 ms 331 ms 2.1
4095 3.0 s 2.0 s 1.5
Average wall time for forward or backward scalar transform using 12 parallel threads on an Intel Xeon X5650 (2.67GHz), with gcc 4.4.5 and "-O3 -march=native -ffast-math -fopenmp" compilation options.
\section Accuracy We claim that the accuracy of SHTns is as good as it can be with double precision floating point math. Rescaling is performed for large transform where the recurrence relation would otherwise underflow the double precision numbers. SHTns has been tested on x86 architecture with SSE2 double precision floating point math (64 bit) to be accurate up to l=16383 at least. The measured error for a back and forth scalar transform using a Gauss-Legendre algorithm for various truncation levels \c lmax is plotted below. \image html accuracy.png \image latex accuracy.png "Accuracy of our scalar Gauss-Legendre transform" width=10cm Don't trust our word, these results are obtained by running the time_SHT program, shipped with SHTns. For example : \code make time_SHT ./time_SHT 511 -iter=1 -quickinit \endcode */ /** \page usage Using SHTns in a C program To make SHTns known by your C program, add \code #include #include \endcode in the preamble of your source file. First of all you need to initialize SHTns. There are two ways to do this : \li calling \ref shtns_init function, \li or calling \ref shtns_create followed by \ref shtns_set_grid or \ref shtns_set_grid_auto. This way lets you choose normalization and optimization. Multi-threaded transforms can be enabled (if available, see \ref compil) by a call to \ref shtns_use_threads before shtns_init or shtns_create. Then you must allocate some memory (with fftw_malloc so that the memory is properly aligned for vectorized code), and finally you can perform some spherical harmonic transforms. \see Look at shtns.h for the function definitions and documentation. \see See the \link SHT_example.c C example \endlink for a simple usage of SHTns in a C program. When your program is ready, compile it adding these options to the compiler (gcc) : \code -L/path/to/lib/ -I/path/to/include/ -lshtns -lfftw3 -lm \endcode \example SHT_example.c \brief A Fortran example program that performs backward and forward Spherical Harmonic Transforms using SHTns. Compile using: \code make SHT_example \endcode */ /** \page f77 Using SHTns with Fortran 77 SHTns provides an interface to Fortran language (compatible with gfortran). First of all you need to initialize SHTns. There are two ways to do this : \li calling one of the \c shtns_init_* subroutines, \li or calling \c shtns_set_size followed by \c shtns_precompute or \c shtns_precompute_auto. This way lets you choose normalization and optimization. Multi-threaded transforms can be enabled (if available, see \ref compil) by a call to \ref shtns_use_threads before shtns_init or shtns_create. Note that you can call initialization function(s) only once, which means that only one size and grid can be used with the Fortran interface (C and Python do not have this limitation). Then you must allocate some memory, and finaly you can perform some spherical harmonic transforms. \see See The full reference of the \ref fortapi \see See the \link SHT_example.f Fortran example \endlink for a simple usage of SHTns from Fortran language. When your program is ready, compile it adding these options to the compiler (gfortran) : \code -L/path/to/lib/ -I/path/to/include/ -lshtns -lfftw3 -lm -lc \endcode \example SHT_example.f \brief A Fortran example program that performs backward and forward Spherical Harmonic Transforms using SHTns. Compile using : \code make SHT_fort_ex \endcode \see fortapi \see \ref f77 */ /** \page loadsave Loading and saving config (or plans) SHTns provides a simple way to load/save a configuration from/to a file. It allows: - to save initialization time if you often use the same transform sizes. - better reproducibility at bit-level between several runs. - ensure that all MPI processes use the same config/plan. All the calling program has to do is add \ref SHT_LOAD_SAVE_CFG to the flags or layout parameter when calling \ref shtns_set_grid / \ref shtns_set_grid_auto (C) or shtns_precompute / shtns_precompute_auto (Fortran). SHTns will first look for the files \c shtns_cfg and \c shtns_cfg_fftw and use the data stored therein. If no configuration stored in the files correspond to the required one, SHTns will proceed as usual and then store the configuration (unless \ref sht_quick_init mode is used). \section mpi_safe Ensuring same config is used across MPI processes: \code flags |= SHT_LOAD_SAVE_CFG; if (rank > 0) MPI_Barrier(); // all processes except 0 wait here until process 0 has initialized shtns. shtns_set_grid_auto(... flags ...); // this will save the config to a file on process 0 and load from file on other processes. if (rank == 0) MPI_Barrier(); // other processes now resume and load the config from file. \endcode */ /** \page python Using SHTns with Python SHTns provides a Python interface that uses NumPy arrays to perform Spherical Harmonic Transforms. The Python interface has been generated using SWIG (but you don't need it to compile the the Python package). \li Compile using (see also \ref compil) \code ./configure --enable-python make make install \endcode Optionally you can indicate an alternate pyhton interpreter to \c ./configure using PYTHON= (see ./configure --help).
If you cannot install as root, you can alternatively intall the python package as user, replacing \code make install \endcode with \code python setup.py install --user \endcode \li Use within a python script or shell \code import shtns \endcode To start using SHTns you must then create an sht object, from which you can call the transform functions and more. The arrays should be created using NumPy. \see See the \link SHT_example.py Python example \endlink to get started, and the full working example \link shallow_water.py \endlink that solves the non-linear shallow water equations. In case of segfaults, you may have alignment problems (should not occur on 64 bit systems). Please file a bugreport or contact the author. \example SHT_example.py \brief A Python example using NumPy and SHTns that performs backward and forward Spherical Harmonic Transforms. \example shallow_water.py \brief Solves the non-linear barotropically unstable shallow water test case in Python. */ /** \page license License \brief CeCILL Free Software License Agreement \htmlonly

CeCILL FREE SOFTWARE LICENSE AGREEMENT

Version 2.1 dated 2013-06-21

Notice

This Agreement is a Free Software license agreement that is the result of discussions between its authors in order to ensure compliance with the two main principles guiding its drafting:

  • firstly, compliance with the principles governing the distribution of Free Software: access to source code, broad rights granted to users,
  • secondly, the election of a governing law, French law, with which it is conformant, both as regards the law of torts and intellectual property law, and the protection that it offers to both authors and holders of the economic rights over software.

The authors of the CeCILL1 license are:

Commissariat à l'énergie atomique et aux énergies alternatives - CEA, a public scientific, technical and industrial research establishment, having its principal place of business at 25 rue Leblanc, immeuble Le Ponant D, 75015 Paris, France.

Centre National de la Recherche Scientifique - CNRS, a public scientific and technological establishment, having its principal place of business at 3 rue Michel-Ange, 75794 Paris cedex 16, France.

Institut National de Recherche en Informatique et en Automatique - Inria, a public scientific and technological establishment, having its principal place of business at Domaine de Voluceau, Rocquencourt, BP 105, 78153 Le Chesnay cedex, France.

Preamble

The purpose of this Free Software license agreement is to grant users the right to modify and redistribute the software governed by this license within the framework of an open source distribution model.

The exercising of this right is conditional upon certain obligations for users so as to preserve this status for all subsequent redistributions.

In consideration of access to the source code and the rights to copy, modify and redistribute granted by the license, users are provided only with a limited warranty and the software's author, the holder of the economic rights, and the successive licensors only have limited liability.

In this respect, the risks associated with loading, using, modifying and/or developing or reproducing the software by the user are brought to the user's attention, given its Free Software status, which may make it complicated to use, with the result that its use is reserved for developers and experienced professionals having in-depth computer knowledge. Users are therefore encouraged to load and test the suitability of the software as regards their requirements in conditions enabling the security of their systems and/or data to be ensured and, more generally, to use and operate it in the same conditions of security. This Agreement may be freely reproduced and published, provided it is not altered, and that no provisions are either added or removed herefrom.

This Agreement may apply to any or all software for which the holder of the economic rights decides to submit the use thereof to its provisions.

Frequently asked questions can be found on the official website of the CeCILL licenses family (http://www.cecill.info/index.en.html) for any necessary clarification.

Article 1 - DEFINITIONS

For the purpose of this Agreement, when the following expressions commence with a capital letter, they shall have the following meaning:

Agreement: means this license agreement, and its possible subsequent versions and annexes.

Software: means the software in its Object Code and/or Source Code form and, where applicable, its documentation, "as is" when the Licensee accepts the Agreement.

Initial Software: means the Software in its Source Code and possibly its Object Code form and, where applicable, its documentation, "as is" when it is first distributed under the terms and conditions of the Agreement.

Modified Software: means the Software modified by at least one Contribution.

Source Code: means all the Software's instructions and program lines to which access is required so as to modify the Software.

Object Code: means the binary files originating from the compilation of the Source Code.

Holder: means the holder(s) of the economic rights over the Initial Software.

Licensee: means the Software user(s) having accepted the Agreement.

Contributor: means a Licensee having made at least one Contribution.

Licensor: means the Holder, or any other individual or legal entity, who distributes the Software under the Agreement.

Contribution: means any or all modifications, corrections, translations, adaptations and/or new functions integrated into the Software by any or all Contributors, as well as any or all Internal Modules.

Module: means a set of sources files including their documentation that enables supplementary functions or services in addition to those offered by the Software.

External Module: means any or all Modules, not derived from the Software, so that this Module and the Software run in separate address spaces, with one calling the other when they are run.

Internal Module: means any or all Module, connected to the Software so that they both execute in the same address space.

GNU GPL: means the GNU General Public License version 2 or any subsequent version, as published by the Free Software Foundation Inc.

GNU Affero GPL: means the GNU Affero General Public License version 3 or any subsequent version, as published by the Free Software Foundation Inc.

EUPL: means the European Union Public License version 1.1 or any subsequent version, as published by the European Commission.

Parties: mean both the Licensee and the Licensor.

These expressions may be used both in singular and plural form.

Article 2 - PURPOSE

The purpose of the Agreement is the grant by the Licensor to the Licensee of a non-exclusive, transferable and worldwide license for the Software as set forth in Article 5 hereinafter for the whole term of the protection granted by the rights over said Software.

Article 3 - ACCEPTANCE

3.1 The Licensee shall be deemed as having accepted the terms and conditions of this Agreement upon the occurrence of the first of the following events:

  • (i) loading the Software by any or all means, notably, by downloading from a remote server, or by loading from a physical medium;
  • (ii) the first time the Licensee exercises any of the rights granted hereunder.

3.2 One copy of the Agreement, containing a notice relating to the characteristics of the Software, to the limited warranty, and to the fact that its use is restricted to experienced users has been provided to the Licensee prior to its acceptance as set forth in Article 3.1 hereinabove, and the Licensee hereby acknowledges that it has read and understood it.

Article 4 - EFFECTIVE DATE AND TERM

4.1 EFFECTIVE DATE

The Agreement shall become effective on the date when it is accepted by the Licensee as set forth in Article 3.1.

4.2 TERM

The Agreement shall remain in force for the entire legal term of protection of the economic rights over the Software.

Article 5 - SCOPE OF RIGHTS GRANTED

The Licensor hereby grants to the Licensee, who accepts, the following rights over the Software for any or all use, and for the term of the Agreement, on the basis of the terms and conditions set forth hereinafter.

Besides, if the Licensor owns or comes to own one or more patents protecting all or part of the functions of the Software or of its components, the Licensor undertakes not to enforce the rights granted by these patents against successive Licensees using, exploiting or modifying the Software. If these patents are transferred, the Licensor undertakes to have the transferees subscribe to the obligations set forth in this paragraph.

5.1 RIGHT OF USE

The Licensee is authorized to use the Software, without any limitation as to its fields of application, with it being hereinafter specified that this comprises:

  1. permanent or temporary reproduction of all or part of the Software by any or all means and in any or all form.

  2. loading, displaying, running, or storing the Software on any or all medium.

  3. entitlement to observe, study or test its operation so as to determine the ideas and principles behind any or all constituent elements of said Software. This shall apply when the Licensee carries out any or all loading, displaying, running, transmission or storage operation as regards the Software, that it is entitled to carry out hereunder.

5.2 ENTITLEMENT TO MAKE CONTRIBUTIONS

The right to make Contributions includes the right to translate, adapt, arrange, or make any or all modifications to the Software, and the right to reproduce the resulting software.

The Licensee is authorized to make any or all Contributions to the Software provided that it includes an explicit notice that it is the author of said Contribution and indicates the date of the creation thereof.

5.3 RIGHT OF DISTRIBUTION

In particular, the right of distribution includes the right to publish, transmit and communicate the Software to the general public on any or all medium, and by any or all means, and the right to market, either in consideration of a fee, or free of charge, one or more copies of the Software by any means.

The Licensee is further authorized to distribute copies of the modified or unmodified Software to third parties according to the terms and conditions set forth hereinafter.

5.3.1 DISTRIBUTION OF SOFTWARE WITHOUT MODIFICATION

The Licensee is authorized to distribute true copies of the Software in Source Code or Object Code form, provided that said distribution complies with all the provisions of the Agreement and is accompanied by:

  1. a copy of the Agreement,

  2. a notice relating to the limitation of both the Licensor's warranty and liability as set forth in Articles 8 and 9,

and that, in the event that only the Object Code of the Software is redistributed, the Licensee allows effective access to the full Source Code of the Software for a period of at least three years from the distribution of the Software, it being understood that the additional acquisition cost of the Source Code shall not exceed the cost of the data transfer.

5.3.2 DISTRIBUTION OF MODIFIED SOFTWARE

When the Licensee makes a Contribution to the Software, the terms and conditions for the distribution of the resulting Modified Software become subject to all the provisions of this Agreement.

The Licensee is authorized to distribute the Modified Software, in source code or object code form, provided that said distribution complies with all the provisions of the Agreement and is accompanied by:

  1. a copy of the Agreement,

  2. a notice relating to the limitation of both the Licensor's warranty and liability as set forth in Articles 8 and 9,

and, in the event that only the object code of the Modified Software is redistributed,

  1. a note stating the conditions of effective access to the full source code of the Modified Software for a period of at least three years from the distribution of the Modified Software, it being understood that the additional acquisition cost of the source code shall not exceed the cost of the data transfer.

5.3.3 DISTRIBUTION OF EXTERNAL MODULES

When the Licensee has developed an External Module, the terms and conditions of this Agreement do not apply to said External Module, that may be distributed under a separate license agreement.

5.3.4 COMPATIBILITY WITH OTHER LICENSES

The Licensee can include a code that is subject to the provisions of one of the versions of the GNU GPL, GNU Affero GPL and/or EUPL in the Modified or unmodified Software, and distribute that entire code under the terms of the same version of the GNU GPL, GNU Affero GPL and/or EUPL.

The Licensee can include the Modified or unmodified Software in a code that is subject to the provisions of one of the versions of the GNU GPL, GNU Affero GPL and/or EUPL and distribute that entire code under the terms of the same version of the GNU GPL, GNU Affero GPL and/or EUPL.

Article 6 - INTELLECTUAL PROPERTY

6.1 OVER THE INITIAL SOFTWARE

The Holder owns the economic rights over the Initial Software. Any or all use of the Initial Software is subject to compliance with the terms and conditions under which the Holder has elected to distribute its work and no one shall be entitled to modify the terms and conditions for the distribution of said Initial Software.

The Holder undertakes that the Initial Software will remain ruled at least by this Agreement, for the duration set forth in Article 4.2.

6.2 OVER THE CONTRIBUTIONS

The Licensee who develops a Contribution is the owner of the intellectual property rights over this Contribution as defined by applicable law.

6.3 OVER THE EXTERNAL MODULES

The Licensee who develops an External Module is the owner of the intellectual property rights over this External Module as defined by applicable law and is free to choose the type of agreement that shall govern its distribution.

6.4 JOINT PROVISIONS

The Licensee expressly undertakes:

  1. not to remove, or modify, in any manner, the intellectual property notices attached to the Software;

  2. to reproduce said notices, in an identical manner, in the copies of the Software modified or not.

The Licensee undertakes not to directly or indirectly infringe the intellectual property rights on the Software of the Holder and/or Contributors, and to take, where applicable, vis-à-vis its staff, any and all measures required to ensure respect of said intellectual property rights of the Holder and/or Contributors.

Article 7 - RELATED SERVICES

7.1 Under no circumstances shall the Agreement oblige the Licensor to provide technical assistance or maintenance services for the Software.

However, the Licensor is entitled to offer this type of services. The terms and conditions of such technical assistance, and/or such maintenance, shall be set forth in a separate instrument. Only the Licensor offering said maintenance and/or technical assistance services shall incur liability therefor.

7.2 Similarly, any Licensor is entitled to offer to its licensees, under its sole responsibility, a warranty, that shall only be binding upon itself, for the redistribution of the Software and/or the Modified Software, under terms and conditions that it is free to decide. Said warranty, and the financial terms and conditions of its application, shall be subject of a separate instrument executed between the Licensor and the Licensee.

Article 8 - LIABILITY

8.1 Subject to the provisions of Article 8.2, the Licensee shall be entitled to claim compensation for any direct loss it may have suffered from the Software as a result of a fault on the part of the relevant Licensor, subject to providing evidence thereof.

8.2 The Licensor's liability is limited to the commitments made under this Agreement and shall not be incurred as a result of in particular: (i) loss due the Licensee's total or partial failure to fulfill its obligations, (ii) direct or consequential loss that is suffered by the Licensee due to the use or performance of the Software, and (iii) more generally, any consequential loss. In particular the Parties expressly agree that any or all pecuniary or business loss (i.e. loss of data, loss of profits, operating loss, loss of customers or orders, opportunity cost, any disturbance to business activities) or any or all legal proceedings instituted against the Licensee by a third party, shall constitute consequential loss and shall not provide entitlement to any or all compensation from the Licensor.

Article 9 - WARRANTY

9.1 The Licensee acknowledges that the scientific and technical state-of-the-art when the Software was distributed did not enable all possible uses to be tested and verified, nor for the presence of possible defects to be detected. In this respect, the Licensee's attention has been drawn to the risks associated with loading, using, modifying and/or developing and reproducing the Software which are reserved for experienced users.

The Licensee shall be responsible for verifying, by any or all means, the suitability of the product for its requirements, its good working order, and for ensuring that it shall not cause damage to either persons or properties.

9.2 The Licensor hereby represents, in good faith, that it is entitled to grant all the rights over the Software (including in particular the rights set forth in Article 5).

9.3 The Licensee acknowledges that the Software is supplied "as is" by the Licensor without any other express or tacit warranty, other than that provided for in Article 9.2 and, in particular, without any warranty as to its commercial value, its secured, safe, innovative or relevant nature.

Specifically, the Licensor does not warrant that the Software is free from any error, that it will operate without interruption, that it will be compatible with the Licensee's own equipment and software configuration, nor that it will meet the Licensee's requirements.

9.4 The Licensor does not either expressly or tacitly warrant that the Software does not infringe any third party intellectual property right relating to a patent, software or any other property right. Therefore, the Licensor disclaims any and all liability towards the Licensee arising out of any or all proceedings for infringement that may be instituted in respect of the use, modification and redistribution of the Software. Nevertheless, should such proceedings be instituted against the Licensee, the Licensor shall provide it with technical and legal expertise for its defense. Such technical and legal expertise shall be decided on a case-by-case basis between the relevant Licensor and the Licensee pursuant to a memorandum of understanding. The Licensor disclaims any and all liability as regards the Licensee's use of the name of the Software. No warranty is given as regards the existence of prior rights over the name of the Software or as regards the existence of a trademark.

Article 10 - TERMINATION

10.1 In the event of a breach by the Licensee of its obligations hereunder, the Licensor may automatically terminate this Agreement thirty (30) days after notice has been sent to the Licensee and has remained ineffective.

10.2 A Licensee whose Agreement is terminated shall no longer be authorized to use, modify or distribute the Software. However, any licenses that it may have granted prior to termination of the Agreement shall remain valid subject to their having been granted in compliance with the terms and conditions hereof.

Article 11 - MISCELLANEOUS

11.1 EXCUSABLE EVENTS

Neither Party shall be liable for any or all delay, or failure to perform the Agreement, that may be attributable to an event of force majeure, an act of God or an outside cause, such as defective functioning or interruptions of the electricity or telecommunications networks, network paralysis following a virus attack, intervention by government authorities, natural disasters, water damage, earthquakes, fire, explosions, strikes and labor unrest, war, etc.

11.2 Any failure by either Party, on one or more occasions, to invoke one or more of the provisions hereof, shall under no circumstances be interpreted as being a waiver by the interested Party of its right to invoke said provision(s) subsequently.

11.3 The Agreement cancels and replaces any or all previous agreements, whether written or oral, between the Parties and having the same purpose, and constitutes the entirety of the agreement between said Parties concerning said purpose. No supplement or modification to the terms and conditions hereof shall be effective as between the Parties unless it is made in writing and signed by their duly authorized representatives.

11.4 In the event that one or more of the provisions hereof were to conflict with a current or future applicable act or legislative text, said act or legislative text shall prevail, and the Parties shall make the necessary amendments so as to comply with said act or legislative text. All other provisions shall remain effective. Similarly, invalidity of a provision of the Agreement, for any reason whatsoever, shall not cause the Agreement as a whole to be invalid.

11.5 LANGUAGE

The Agreement is drafted in both French and English and both versions are deemed authentic.

Article 12 - NEW VERSIONS OF THE AGREEMENT

12.1 Any person is authorized to duplicate and distribute copies of this Agreement.

12.2 So as to ensure coherence, the wording of this Agreement is protected and may only be modified by the authors of the License, who reserve the right to periodically publish updates or new versions of the Agreement, each with a separate number. These subsequent versions may address new issues encountered by Free Software.

12.3 Any Software distributed under a given version of the Agreement may only be subsequently distributed under the same version of the Agreement or a subsequent version, subject to the provisions of Article 5.3.4.

Article 13 - GOVERNING LAW AND JURISDICTION

13.1 The Agreement is governed by French law. The Parties agree to endeavor to seek an amicable solution to any disagreements or disputes that may arise during the performance of the Agreement.

13.2 Failing an amicable solution within two (2) months as from their occurrence, and unless emergency proceedings are necessary, the disagreements or disputes shall be referred to the Paris Courts having jurisdiction, by the more diligent Party.

\endhtmlonly */