ARTS  2.4.0(git:4fb77825)
wigner_functions.cc
Go to the documentation of this file.
1 /* Copyright (C) 2012
2  * Richard Larsson <ric.larsson@gmail.com>
3  *
4  * This program is free software; you can redistribute it and/or modify it
5  * under the terms of the GNU General Public License as published by the
6  * Free Software Foundation; either version 2, or (at your option) any
7  * later version.
8  *
9  * This program is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with this program; if not, write to the Free Software
16  * Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307,
17  * USA. */
18 
27 #include "constants.h"
28 #include "wigner_functions.h"
29 #include <algorithm>
30 #include "arts_omp.h"
31 
32 #if DO_FAST_WIGNER
33 #define WIGNER3 fw3jja6
34 #define WIGNER6 fw6jja
35 #else
36 #define WIGNER3 wig3jj
37 #define WIGNER6 wig6jj
38 #endif
39 
41  const Rational j2,
42  const Rational j3,
43  const Rational m1,
44  const Rational m2,
45  const Rational m3) {
46  const int a = (2 * j1).toInt(), b = (2 * j2).toInt(), c = (2 * j3).toInt(),
47  d = (2 * m1).toInt(), e = (2 * m2).toInt(), f = (2 * m3).toInt();
48  double g;
49  const int j = std::max({std::abs(a),
50  std::abs(b),
51  std::abs(c),
52  std::abs(d),
53  std::abs(e),
54  std::abs(f)}) *
55  3 / 2 +
56  1;
57 
58  wig_temp_init(j);
59  g = WIGNER3(a, b, c, d, e, f);
60  wig_temp_free();
61 
62  return Numeric(g);
63 }
64 
66  const Rational j2,
67  const Rational j3,
68  const Rational l1,
69  const Rational l2,
70  const Rational l3) {
71  const int a = (2 * j1).toInt(), b = (2 * j2).toInt(), c = (2 * j3).toInt(),
72  d = (2 * l1).toInt(), e = (2 * l2).toInt(), f = (2 * l3).toInt();
73  double g;
74  const int j = std::max({std::abs(a),
75  std::abs(b),
76  std::abs(c),
77  std::abs(d),
78  std::abs(e),
79  std::abs(f)});
80 
81  wig_temp_init(j);
82  g = WIGNER6(a, b, c, d, e, f);
83  wig_temp_free();
84 
85  return Numeric(g);
86 }
87 
89  int Ji, int Jf, int Ji_p, int Jf_p, int L, int li, int lf) {
90  return WIGNER3(Ji_p, L, Ji, li, 0, -li) * WIGNER3(Jf_p, L, Jf, -lf, 0, lf) *
91  WIGNER6(Ji, Jf, 2, Jf_p, Ji_p, L) * Numeric(L + 1);
92 }
93 
95  int Nl, int Nk, int Jl, int Jk, int Jl_p, int Jk_p, int L) {
96  return WIGNER3(Nl, Nk, L, 0, 0, 0) * WIGNER6(L, Jk, Jl, 2, Nl, Nk) *
97  WIGNER6(L, Jk_p, Jl_p, 2, Nl, Nk) * WIGNER6(L, Jk, Jl, 2, Jl_p, Jk_p);
98 }
99 
100 
102  int lower;
103  int upper;
104 };
105 
106 
107 
123  const Rational L2)
124 {
125  return {abs(L1-L2).toInt(), (L1+L2).toInt()};
126 }
127 
128 
137  const Wigner3JTriangleLimit lim2)
138 {
139  return
140  {
141  lim1.lower > lim2.lower ?
142  ((lim1.lower % 2) ?
143  lim1.lower + 1 :
144  lim1.lower) :
145  ((lim2.lower % 2) ?
146  lim2.lower + 1 :
147  lim2.lower),
148  lim1.upper < lim2.upper ?
149  ((lim1.upper % 2) ?
150  lim1.upper - 1 :
151  lim1.upper) :
152  ((lim2.upper % 2) ?
153  lim2.upper - 1 :
154  lim2.upper)
155  };
156 }
157 
158 
160 {
161  using Constant::h;
162  using Constant::pow2;
163  using Constant::pow3;
164 
165  constexpr auto B = 43100.44276e6;
166  constexpr auto D = 145.1271e3;
167  constexpr auto H = 49e-3;
168  constexpr auto lB = 59501.3438e6;
169  constexpr auto lD = 58.3680e3;
170  constexpr auto lH = 290.8e-3;
171  constexpr auto gB = -252.58634e6;
172  constexpr auto gD = -243.42;
173  constexpr auto gH = -1.46e-3;
174 
175  const auto lJ = (J * (J + 1)).toNumeric();
176  return h * (
177  B*lJ - D*pow2(lJ) + H*pow3(lJ) -
178  (gB + gD*lJ + gH*pow2(lJ)) +
179  2.0/3.0 * (lB + lD*lJ + lH*pow2(lJ)));
180 }
181 
182 
184 {
185  using Constant::k;
186 
187  auto lambda = 0.39;
188  auto beta = 0.567;
189  auto el = o2_ecs_erot_jn_same(L);
190 
191  return
192  (2*L + 1).toNumeric() / pow(L * (L+1), lambda) * std::exp(-beta * el / (k*T));
193 }
194 
195 
197 {
199  using Constant::h_bar;
200  using Constant::pow2;
201  using Constant::m_u;
202  using Constant::pi;
203  using Constant::k;
204 
205  auto dc = angstrom2meter(0.61);
206 
207  auto en = o2_ecs_erot_jn_same(N);
208  auto enm2 = o2_ecs_erot_jn_same(N-2);
209  auto wnnm2 = (en - enm2) / h_bar;
210 
211  auto v_bar = std::sqrt(8*k*T/(pi * 31.989830 * m_u));
212  auto tauc = dc / v_bar;
213 
214  return 1.0 / pow2(1 + 1.0/24.0 * pow2(wnnm2 * tauc));
215 }
216 
217 
219  const Rational& Ji,
220  const Rational& Jf,
221  const Rational& Ni,
222  const Rational& Nf,
223  const Rational& Si,
224  const Rational& Sf,
225  const Rational& Ji_p,
226  const Rational& Jf_p,
227  const Rational& Ni_p,
228  const Rational& Nf_p,
229  const Rational& n,
230  const Numeric& T)
231 {
233 
234 // std::cout<<"import numpy as np";
235 // std::cout << "wigsym=np.array([";
236 // for (Index L=1; L<=250; L++) {
237 // std::cout<<"["<<L<<", "<< o2_ecs_adiabatic_factor_makarov(L, T)<<","<<o2_ecs_ql_makarov(L, T)<<"],\n";
238 // }
239 // std::cout<<"])\n";
240 // std::exit(1);
241 
242  // Limits above 2 and below the largest index there is
245  find_wigner3j_limits(Nf_p, Nf)));
246 
247  auto f = sqrt(2*Ni+1) * sqrt(2*Ni_p+1) * sqrt(2*Jf+1) * sqrt(2*Jf_p+1) * sqrt(2*Nf+1) * sqrt(2*Nf_p+1) * (2*Ji_p+1).toNumeric();
248  auto g = even(Ji_p + Ji + n) ? 1 : -1; // -1^(Ji_p + Ji + n)
249  auto OmegaNi = o2_ecs_ql_makarov(Ni, T);
250 
251  for (int L=lims.lower; L<=lims.upper; L+=2) {
252  auto OmegaL = o2_ecs_adiabatic_factor_makarov(Rational(L), T);
253  auto QL = o2_ecs_ql_makarov(Rational(L), T);
254  auto a = WIGNER3(Ni_p.toInt(2), Ni.toInt(2), 2*L, 0, 0, 0);
255  auto b = WIGNER3(Nf_p.toInt(2), Nf.toInt(2), 2*L, 0, 0, 0);
256  auto c = WIGNER6(2*L, Ji.toInt(2), Ji_p.toInt(2), Si.toInt(2), Ni_p.toInt(2), Ni.toInt(2));
257  auto d = WIGNER6(2*L, Jf.toInt(2), Jf_p.toInt(2), Sf.toInt(2), Nf_p.toInt(2), Nf.toInt(2));
258  auto e = WIGNER6(2*L, Ji.toInt(2), Ji_p.toInt(2), n.toInt(2), Jf_p.toInt(2), Jf.toInt(2));
259 
260  o2_ecs_wigner_symbol_tran += (a * b * c * d * e * f * g) * (2*L + 1) * QL * (OmegaNi / OmegaL);
261  }
262 
264 }
265 
266 
268 {
269  return
270  (even(Jlo + N) ? 1 : -1) *
271  sqrt(6 * (2*Jlo + 1) * (2*Jup + 1)) *
272  WIGNER6(2, 2, 2, Jup.toInt(2), Jlo.toInt(2), N.toInt(2));
273 }
274 
275 
277  [[maybe_unused]] int fastest,
278  int size)
279 {
280  if (size == 3) {
281  #if DO_FAST_WIGNER
282  fastwigxj_load(FAST_WIGNER_PATH_3J, 3, NULL);
283  #ifdef _OPENMP
284  fastwigxj_thread_dyn_init(3, fastest);
285  #else
286  fastwigxj_dyn_init(3, fastest);
287  #endif
288  #endif
289  wig_table_init(largest, 3);
290 
291  return largest;
292  } else if (size == 6) {
293  #if DO_FAST_WIGNER
294  fastwigxj_load(FAST_WIGNER_PATH_3J, 3, NULL);
295  fastwigxj_load(FAST_WIGNER_PATH_6J, 6, NULL);
296  #ifdef _OPENMP
297  fastwigxj_thread_dyn_init(3, fastest);
298  fastwigxj_thread_dyn_init(6, fastest);
299  #else
300  fastwigxj_dyn_init(3, fastest);
301  fastwigxj_dyn_init(6, fastest);
302  #endif
303  #endif
304  wig_table_init(largest * 2, 6);
305 
306  return largest;
307  } else {
308  return 0;
309  }
310 }
311 
312 
313 bool is_wigner_ready(int j) {
314  extern int wigxjpf_max_prime_decomp;
315  return not(j > wigxjpf_max_prime_decomp);
316 }
317 
318 
319 bool is_wigner3_ready(const Rational& J) {
320  const int test = J.toInt(6) / 2 + 1; // nb. J can be half-valued
321  return is_wigner_ready(test);
322 }
323 
324 
325 bool is_wigner6_ready(const Rational& J) {
326  const int test = J.toInt(4) + 1; // nb. J can be half-valued
327  return is_wigner_ready(test);
328 }
o2_ecs_adiabatic_factor_makarov
Numeric o2_ecs_adiabatic_factor_makarov(Rational N, Numeric T)
Definition: wigner_functions.cc:196
Wigner3JTriangleLimit
Definition: wigner_functions.cc:101
wigner6j
Numeric wigner6j(const Rational j1, const Rational j2, const Rational j3, const Rational l1, const Rational l2, const Rational l3)
Wigner 6J symbol.
Definition: wigner_functions.cc:65
wigner_functions.h
Wigner symbol interactions.
is_wigner3_ready
bool is_wigner3_ready(const Rational &J)
Tells if the function is ready for Wigner 3J calculations.
Definition: wigner_functions.cc:319
o2_ecs_erot_jn_same
Numeric o2_ecs_erot_jn_same(Rational J)
Energy of the J=N line at J.
Definition: wigner_functions.cc:159
abs
#define abs(x)
Definition: legacy_continua.cc:20626
o2_ecs_ql_makarov
Numeric o2_ecs_ql_makarov(Rational L, Numeric T)
Definition: wigner_functions.cc:183
WIGNER6
#define WIGNER6
Definition: wigner_functions.cc:37
wigner3j
Numeric wigner3j(const Rational j1, const Rational j2, const Rational j3, const Rational m1, const Rational m2, const Rational m3)
Wigner 3J symbol.
Definition: wigner_functions.cc:40
find_even_limits
constexpr Wigner3JTriangleLimit find_even_limits(const Wigner3JTriangleLimit lim1, const Wigner3JTriangleLimit lim2)
Combines two limits and return the highest even-numbered low value and the lowest numbered even high ...
Definition: wigner_functions.cc:136
is_wigner6_ready
bool is_wigner6_ready(const Rational &J)
Tells if the function is ready for Wigner 6J calculations.
Definition: wigner_functions.cc:325
make_wigner_ready
Index make_wigner_ready(int largest, [[maybe_unused]] int fastest, int size)
Definition: wigner_functions.cc:276
sqrt
Numeric sqrt(const Rational r)
Square root.
Definition: rational.h:620
Wigner3JTriangleLimit::upper
int upper
Definition: wigner_functions.cc:103
pow
Numeric pow(const Rational base, Numeric exp)
Power of.
Definition: rational.h:628
Conversion::angstrom2meter
constexpr T angstrom2meter(T x)
Definition: constants.h:495
WIGNER3
#define WIGNER3
Definition: wigner_functions.cc:36
even
constexpr bool even(const Rational r)
Returns true if even integer.
Definition: rational.h:762
ARTS::Group::Rational
Rational Rational
Definition: autoarts.h:96
Rational::toInt
constexpr int toInt(int n=1) const
Converts the value to int by n-scaled division in Index form.
Definition: rational.h:157
Numeric
NUMERIC Numeric
The type to use for all floating point numbers.
Definition: matpack.h:33
max
#define max(a, b)
Definition: legacy_continua.cc:20629
co2_ecs_wigner_symbol
Numeric co2_ecs_wigner_symbol(int Ji, int Jf, int Ji_p, int Jf_p, int L, int li, int lf)
Returns the wigner symbol used in Niro etal 2004.
Definition: wigner_functions.cc:88
Wigner3JTriangleLimit::lower
int lower
Definition: wigner_functions.cc:102
N
#define N
Definition: rng.cc:164
o2_ecs_wigner_symbol
Numeric o2_ecs_wigner_symbol(int Nl, int Nk, int Jl, int Jk, int Jl_p, int Jk_p, int L)
Returns the wigner symbol used in Makarov etal 2013.
Definition: wigner_functions.cc:94
Constant::pow3
constexpr T pow3(T x)
power of three
Definition: constants.h:70
constants.h
Constants of physical expressions as constexpr.
find_wigner3j_limits
constexpr Wigner3JTriangleLimit find_wigner3j_limits(const Rational L1, const Rational L2)
Finds the upper and lower limits of L3 given a Wigner-3J symbol:
Definition: wigner_functions.cc:122
o2_makarov2013_reduced_dipole
Numeric o2_makarov2013_reduced_dipole(const Rational &Jup, const Rational &Jlo, const Rational &N)
Returns the reduced dipole moment following Makarov etal 2013.
Definition: wigner_functions.cc:267
Constant::pow2
constexpr T pow2(T x)
power of two
Definition: constants.h:64
Index
INDEX Index
The type to use for all integer numbers and indices.
Definition: matpack.h:39
arts_omp.h
Header file for helper functions for OpenMP.
beta
#define beta
Definition: legacy_continua.cc:21435
o2_ecs_wigner_symbol_tran
Numeric o2_ecs_wigner_symbol_tran(const Rational &Ji, const Rational &Jf, const Rational &Ni, const Rational &Nf, const Rational &Si, const Rational &Sf, const Rational &Ji_p, const Rational &Jf_p, const Rational &Ni_p, const Rational &Nf_p, const Rational &n, const Numeric &T)
Returns the wigner symbol used in Tran etal 2006.
Definition: wigner_functions.cc:218
is_wigner_ready
bool is_wigner_ready(int j)
Tells if the function can deal with the input integer.
Definition: wigner_functions.cc:313
Rational
Implements rational numbers to work with other ARTS types.
Definition: rational.h:54