-
Notifications
You must be signed in to change notification settings - Fork 104
/
makemaggriddh.3
151 lines (151 loc) · 5.78 KB
/
makemaggriddh.3
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
.\" Automatically generated by Pandoc 3.1.3
.\"
.\" Define V font for inline verbatim, using C font in formats
.\" that render this, and otherwise B font.
.ie "\f[CB]x\f[]"x" \{\
. ftr V B
. ftr VI BI
. ftr VB B
. ftr VBI BI
.\}
.el \{\
. ftr V CR
. ftr VI CI
. ftr VB CB
. ftr VBI CBI
.\}
.TH "makemaggriddh" "1" "2021-02-15" "Fortran 95" "SHTOOLS 4.11"
.hy
.SH MakeMagGridDH
.PP
Create 2D cylindrical maps on a flattened ellipsoid of all three vector
components of the magnetic field, the magnitude of the magnetic field,
and the magnetic potential.
.SH Usage
.PP
call MakeMagGridDH (\f[V]cilm\f[R], \f[V]lmax\f[R], \f[V]r0\f[R],
\f[V]a\f[R], \f[V]f\f[R], \f[V]rad\f[R], \f[V]theta\f[R], \f[V]phi\f[R],
\f[V]total\f[R], \f[V]n\f[R], \f[V]sampling\f[R], \f[V]lmax_calc\f[R],
\f[V]pot_grid\f[R], \f[V]extend\f[R], \f[V]exitstatus\f[R])
.SH Parameters
.TP
\f[V]cilm\f[R] : input, real(dp), dimension (2, \f[V]lmax\f[R]+1, \f[V]lmax\f[R]+1)
The real Schmidt semi-normalized spherical harmonic coefficients to be
expanded in the space domain.
The coefficients \f[V]C1lm\f[R] and \f[V]C2lm\f[R] refer to the cosine
(\f[V]Clm\f[R]) and sine (\f[V]Slm\f[R]) coefficients, respectively,
with \f[V]Clm=cilm(1,l+1,m+1)\f[R] and \f[V]Slm=cilm(2,l+1,m+1)\f[R].
Alternatively, \f[V]C1lm\f[R] and \f[V]C2lm\f[R] correspond to the
positive and negative order coefficients, respectively.
The coefficients are assumed to have units of nT.
.TP
\f[V]lmax\f[R] : input, integer(int32)
The maximum spherical harmonic degree of the coefficients
\f[V]cilm\f[R].
This determines the number of samples of the output grids,
\f[V]n=2*lmax+2\f[R], and the latitudinal sampling interval,
\f[V]90/(lmax+1)\f[R].
.TP
\f[V]r0\f[R] : input, real(dp)
The reference radius of the spherical harmonic coefficients.
.TP
\f[V]a\f[R] : input, real(dp)
The semi-major axis of the flattened ellipsoid on which the field is
computed.
.TP
\f[V]f\f[R] : input, real(dp)
The flattening of the reference ellipsoid: i.e.,
\f[V]F=(R_equator-R_pole)/R_equator\f[R].
.TP
\f[V]rad\f[R] : output, real(dp), dimension(nlat, nlong)
A 2D map of radial component of the magnetic field that conforms to the
sampling theorem of Driscoll and Healy (1994).
If \f[V]sampling\f[R] is 1, the grid is equally sampled and is
dimensioned as (\f[V]n\f[R] by \f[V]n\f[R]), where \f[V]n\f[R] is
\f[V]2lmax+2\f[R].
If sampling is 2, the grid is equally spaced and is dimensioned as
(\f[V]n\f[R] by 2\f[V]n\f[R]).
The first latitudinal band of the grid corresponds to 90 N, the
latitudinal sampling interval is 180/\f[V]n\f[R] degrees, and the
default behavior is to exclude the latitudinal band for 90 S.
The first longitudinal band of the grid is 0 E, by default the
longitudinal band for 360 E is not included, and the longitudinal
sampling interval is 360/\f[V]n\f[R] for an equally sampled and
180/\f[V]n\f[R] for an equally spaced grid, respectively.
If \f[V]extend\f[R] is 1, the longitudinal band for 360 E and the
latitudinal band for 90 S will be included, which increases each of the
dimensions of the grid by 1.
.TP
\f[V]theta\f[R] : output, real(dp), dimension(nlat, nlong)
A 2D equally sampled or equally spaced grid of the theta component of
the magnetic field.
.TP
\f[V]phi\f[R] : output, real(dp), dimension(nlat, nlong)
A 2D equally sampled or equally spaced grid of the phi component of the
magnetic field.
.TP
\f[V]total\f[R] : output, real(dp), dimension(nlat, nlong)
A 2D equally sampled or equally spaced grid of the total magnetic field
strength.
.TP
\f[V]n\f[R] : output, integer(int32)
The number of samples in latitude of the output grids.
This is equal to \f[V]2lmax+2\f[R].
.TP
\f[V]sampling\f[R] : optional, input, integer(int32), default = 1
If 1 (default) the output grids are equally sampled (\f[V]n\f[R] by
\f[V]n\f[R]).
If 2, the grids are equally spaced (\f[V]n\f[R] by 2\f[V]n\f[R]).
.TP
\f[V]lmaxcalc\f[R] : optional, input, integer(int32), default = \f[V]lmax\f[R]
The maximum spherical harmonic degree used in evaluating the functions.
This must be less than or equal to \f[V]lmax\f[R].
.TP
\f[V]potgrid\f[R] : output, real(dp), dimension(nlat, nlong)
A 2D equally sampled or equaly spaced grid of the magnetic potential.
.TP
\f[V]extend\f[R] : input, optional, integer(int32), default = 0
If 1, compute the longitudinal band for 360 E and the latitudinal band
for 90 S.
This increases each of the dimensions of the grids by 1.
.TP
\f[V]exitstatus\f[R] : output, optional, integer(int32)
If present, instead of executing a STOP when an error is encountered,
the variable exitstatus will be returned describing the error.
0 = No errors; 1 = Improper dimensions of input array; 2 = Improper
bounds for input variable; 3 = Error allocating memory; 4 = File IO
error.
.SH Description
.PP
\f[V]MakeMagGridDH\f[R] will create 2-dimensional cylindrical maps from
the spherical harmonic coefficients \f[V]cilm\f[R] of all three
components of the magnetic field, the total field strength, and the
magnetic potential.
The magnetic potential is given by
.PP
\f[V]V = R0 Sum_{l=1}\[ha]LMAX (R0/r)\[ha]{l+1} Sum_{m=-l}\[ha]l C_{lm} Y_{lm}\f[R]
.PP
and the magnetic field is
.PP
\f[V]B = - Grad V\f[R].
.PP
The coefficients are referenced to a radius \f[V]r0\f[R], and the
function is computed on a flattened ellipsoid with semi-major axis
\f[V]a\f[R] (i.e., the mean equatorial radius) and flattening
\f[V]f\f[R].
.PP
The default is to use an input grid that is equally sampled (\f[V]n\f[R]
by \f[V]n\f[R]), but this can be changed to use an equally spaced grid
(\f[V]n\f[R] by 2\f[V]n\f[R]) by the optional argument
\f[V]sampling\f[R].
The redundant longitudinal band for 360 E and the latitudinal band for
90 S are excluded by default, but these can be computed by specifying
the optional argument \f[V]extend\f[R].
.SH Reference
.PP
Driscoll, J.R.
and D.M.
Healy, Computing Fourier transforms and convolutions on the 2-sphere,
Adv.
Appl.
Math., 15, 202-250, 1994.