SHTns 3.7
SHT_example.f

A Fortran example program that performs backward and forward Spherical Harmonic Transforms using SHTns.

A Fortran example program that performs backward and forward Spherical Harmonic Transforms using SHTns. Compile using :

make SHT_fort_ex
See also
Fortran API.
Basic Fortran 77 interface
1!
2! Copyright (c) 2010-2018 Centre National de la Recherche Scientifique.
3! written by Nathanael Schaeffer (CNRS, ISTerre, Grenoble, France).
4!
5! nathanael.schaeffer@univ-grenoble-alpes.fr
6!
7! This software is governed by the CeCILL license under French law and
8! abiding by the rules of distribution of free software. You can use,
9! modify and/or redistribute the software under the terms of the CeCILL
10! license as circulated by CEA, CNRS and INRIA at the following URL
11! "http://www.cecill.info".
12!
13! The fact that you are presently reading this means that you have had
14! knowledge of the CeCILL license and that you accept its terms.
15!
16
17
18!! A Fortran example program that performs backward and forward Spherical Harmonic Transforms using SHTns
19 PROGRAM sht_example
20
21 IMPLICIT NONE
22! import useful parameters for shtns initialization
23 include 'shtns.f'
24
25 integer*4 lmax, mmax, mres
26 integer*4 nlat, nphi
27 integer*4 nlm
28 integer*4 layout
29 integer*4 norm
30 real*8 eps_polar
31 complex*16, allocatable :: Slm(:), Tlm(:)
32 real*8, allocatable :: sh(:,:), th(:,:)
33
34 integer i,lm, m
35
36! set size of transform
37 lmax = 5
38 mmax = 2
39 mres = 2
40 nphi = 6
41 nlat = 8
42
43! compute sizes required for arrays.
44 call shtns_calc_nlm(nlm, lmax, mmax, mres)
45 print*,'NLM=',nlm
46
47! displays information during initialization
48 call shtns_verbose(1)
49
50! enable multi-threaded transform (OpenMP) if supported.
51 call shtns_use_threads(0)
52
53! init SHT. SHT_PHI_CONTIGUOUS is defined in 'shtns.f'
54 layout = sht_phi_contiguous
55 call shtns_init_sh_gauss(layout, lmax, mmax, mres, nlat, nphi)
56
57! alternatively, you can use the two following calls, giving more control on the initialization process,
58! namely you can choose a normalization with 'norm' and control the polar optimization
59! with 'eps_polar' : from 0 (no optimization) to 1e-6 (agressive optimization).
60! norm = SHT_ORTHONORMAL + SHT_REAL_NORM
61! call shtns_set_size(lmax, mmax, mres, norm)
62! eps_polar = 1.e-10
63! call shtns_precompute(SHT_GAUSS, layout, eps_polar, nlat, nphi)
64
65! display information:
66 call shtns_print_cfg()
67
68! allocate memory for spectral and spatial representation
69 allocate ( slm(1:nlm), tlm(1:nlm) )
70 allocate ( sh(1:nphi,1:nlat), th(1:nphi,1:nlat) )
71 slm(:) = 0.0
72
73! get index for given l and m
74 call shtns_lmidx(lm, 1, 0)
75 slm(lm) = 1.0
76 print*
77 print*,slm(:)
78
79 call shtns_sh_to_spat(slm, sh)
80
81! print all spatial data
82 print*
83 do i=1,nphi
84 print*,'phi index=',i
85 print*,sh(i,:)
86 enddo
87
88! spatial to spectral transform. DOES NOT PRESERVE THE CONTENT OF INPUT 'Sh'
89 call shtns_spat_to_sh(sh, slm)
90
91 print*
92 print*,slm(:)
93
94! print SH data, grouping by m (shows how to access SH coefficients)
95 print*
96 do m=0, mres*mmax, mres
97 print*,'m=', m
98 call shtns_lmidx(lm,m,m)
99 print*,slm(lm : lm+(lmax-m))
100 enddo
101
102! Legendre transform (fixed m, no fft):
103! we need a larger spatial grid for these functions to work (nlat >= 16)
104 nlat = 32
105 call shtns_init_sh_gauss(layout, lmax, mmax, mres, nlat, nphi)
106 m = 1*mres
107 call shtns_lmidx(lm, m, m)
108 slm(lm+1) = 1.0
109 call shtns_sh_to_spat_ml(m/mres, slm(lm), tlm, lmax)
110 print*
111 print*,'spectral m=', m
112 print*,slm(lm : lm+lmax-m+1)
113 print*,'spatial m=', m
114 print*,tlm(1 : nlat)
115
116 stop
117 END
int shtns_use_threads(int num_threads)
Enables multi-thread transform using OpenMP with num_threads (if available). Returns number of thread...
Definition sht_init.c:1605