qsym2/angmom/
sh_conversion.rs

1//! Conversions between spherical/solid harmonics and Cartesian functions.
2
3use std::cmp::Ordering;
4
5use factorial::Factorial;
6use ndarray::{Array2, Axis};
7use num::{BigUint, Complex};
8use num_traits::{cast::ToPrimitive, Zero};
9
10use crate::basis::ao::{CartOrder, PureOrder};
11use crate::permutation::PermutableCollection;
12
13#[cfg(test)]
14#[path = "sh_conversion_tests.rs"]
15mod sh_conversion_tests;
16
17/// Calculates the number of combinations of `n` things taken `r` at a time (signed arguments).
18///
19/// If $`n < 0`$ or $`r < 0`$ or $`r > n`$, `0` is returned.
20///
21/// # Arguments
22///
23/// * `n` - Number of things.
24/// * `r` - Number of elements taken.
25///
26/// # Returns
27///
28/// The number of combinations.
29fn comb(n: i32, r: i32) -> BigUint {
30    if n < 0 || r < 0 || r > n {
31        BigUint::zero()
32    } else {
33        let nu = u32::try_from(n).expect("Unable to convert `n` to `u32`.");
34        let ru = u32::try_from(r).expect("Unable to convert `r` to `u32`.");
35        (nu - ru + 1..=nu).product::<BigUint>()
36            / BigUint::from(ru)
37                .checked_factorial()
38                .unwrap_or_else(|| panic!("Unable to compute the factorial of {ru}."))
39    }
40}
41
42/// Calculates the number of combinations of `n` things taken `r` at a time (unsigned arguments).
43///
44/// If $`r > n`$, `0` is returned.
45///
46/// # Arguments
47///
48/// * `n` - Number of things.
49/// * `r` - Number of elements taken.
50///
51/// # Returns
52///
53/// The number of combinations.
54fn combu(nu: u32, ru: u32) -> BigUint {
55    if ru > nu {
56        BigUint::zero()
57    } else {
58        (nu - ru + 1..=nu).product::<BigUint>()
59            / BigUint::from(ru)
60                .checked_factorial()
61                .unwrap_or_else(|| panic!("Unable to compute the factorial of {ru}."))
62    }
63}
64
65/// Calculates the number of permutations of `n` things taken `r` at a time (signed arguments).
66///
67/// If $`n < 0`$ or $`r < 0`$ or $`r > n`$, `0` is returned.
68///
69/// # Arguments
70///
71/// * `n` - Number of things.
72/// * `r` - Number of elements taken.
73///
74/// # Returns
75///
76/// The number of permutations.
77#[allow(dead_code)]
78fn perm(n: i32, r: i32) -> BigUint {
79    if n < 0 || r < 0 || r > n {
80        BigUint::zero()
81    } else {
82        let nu = u32::try_from(n).expect("Unable to convert `n` to `u32`.");
83        let ru = u32::try_from(r).expect("Unable to convert `r` to `u32`.");
84        (nu - ru + 1..=nu).product::<BigUint>()
85    }
86}
87
88/// Calculates the number of permutations of `n` things taken `r` at a time (unsigned arguments).
89///
90/// If $`r > n`$, `0` is returned.
91///
92/// # Arguments
93///
94/// * `n` - Number of things.
95/// * `r` - Number of elements taken.
96///
97/// # Returns
98///
99/// The number of permutations.
100#[allow(dead_code)]
101fn permu(n: u32, r: u32) -> BigUint {
102    if r > n {
103        BigUint::zero()
104    } else {
105        (n - r + 1..=n).product::<BigUint>()
106    }
107}
108
109/// Obtains the normalisation constant for a solid harmonic Gaussian, as given in Equation 8 of
110/// Schlegel, H. B. & Frisch, M. J. Transformation between Cartesian and pure spherical harmonic
111/// Gaussians. *International Journal of Quantum Chemistry* **54**, 83–87 (1995),
112/// [DOI](https://doi.org/10.1002/qua.560540202).
113///
114/// The complex solid harmonic Gaussian is defined in Equation 1 of the above reference as
115///
116/// ```math
117/// \tilde{g}(\alpha, l, m, n, \mathbf{r})
118///     = \tilde{N}(n, \alpha) Y_l^m r^n e^{-\alpha r^2},
119/// ```
120///
121/// where $`Y_l^m`$ is a complex spherical harmonic of degree $`l`$ and order $`m`$.
122///
123/// # Arguments
124///
125/// * `n` - The non-negative exponent of the radial part of the solid harmonic Gaussian.
126/// * `alpha` - The coefficient on the exponent of the Gaussian term.
127///
128/// # Returns
129///
130/// The normalisation constant $`\tilde{N}(n, \alpha)`$.
131fn norm_sph_gaussian(n: u32, alpha: f64) -> f64 {
132    let num = (BigUint::from(2u64).pow(2 * n + 3)
133        * BigUint::from(u64::from(n) + 1)
134            .checked_factorial()
135            .unwrap_or_else(|| panic!("Unable to compute the factorial of {}.", u64::from(n) + 1)))
136    .to_f64()
137    .expect("Unable to convert a `BigUint` value to `f64`.")
138        * alpha.powf(f64::from(n) + 1.5);
139    let den = BigUint::from(2 * u64::from(n) + 2)
140        .checked_factorial()
141        .unwrap_or_else(|| {
142            panic!(
143                "Unable to compute the factorial of {}.",
144                2 * u64::from(n) + 2
145            )
146        })
147        .to_f64()
148        .expect("Unable to convert a `BigUint` value to `f64`.")
149        * std::f64::consts::PI.sqrt();
150    (num / den).sqrt()
151}
152
153/// Obtains the normalisation constant for a Cartesian Gaussian, as given in Equation 9 of
154/// Schlegel, H. B. & Frisch, M. J. Transformation between Cartesian and pure spherical harmonic
155/// Gaussians. *International Journal of Quantum Chemistry* **54**, 83–87 (1995),
156/// [DOI](https://doi.org/10.1002/qua.560540202).
157///
158/// The Cartesian Gaussian is defined in Equation 2 of the above reference as
159///
160/// ```math
161/// g(\alpha, l_x, l_y, l_z, \mathbf{r})
162///     = N(l_x, l_y, l_z, \alpha) x^{l_x} y^{l_y} z^{l_z} e^{-\alpha r^2}.
163/// ```
164///
165/// # Arguments
166///
167/// * `lcartqns` - A tuple of $`(l_x, l_y, l_z)`$ specifying the non-negative exponents of
168///   the Cartesian components of the Cartesian Gaussian.
169/// * `alpha` - The coefficient on the exponent of the Gaussian term.
170///
171/// # Returns
172///
173/// The normalisation constant $`N(l_x, l_y, l_z, \alpha)`$.
174fn norm_cart_gaussian(lcartqns: (u32, u32, u32), alpha: f64) -> f64 {
175    let (lx, ly, lz) = lcartqns;
176    let lcart = lx + ly + lz;
177    let num = (BigUint::from(2u32).pow(2 * lcart)
178        * BigUint::from(lx)
179            .checked_factorial()
180            .unwrap_or_else(|| panic!("Unable to compute the factorial of {lx}."))
181        * BigUint::from(ly)
182            .checked_factorial()
183            .unwrap_or_else(|| panic!("Unable to compute the factorial of {ly}."))
184        * BigUint::from(lz)
185            .checked_factorial()
186            .unwrap_or_else(|| panic!("Unable to compute the factorial of {lz}.")))
187    .to_f64()
188    .expect("Unable to convert a `BigUint` value to `f64`.")
189        * alpha.powf(f64::from(lcart) + 1.5);
190    let den = (BigUint::from(2 * lx)
191        .checked_factorial()
192        .unwrap_or_else(|| panic!("Unable to compute the factorial of {}.", 2 * lx))
193        * BigUint::from(2 * ly)
194            .checked_factorial()
195            .unwrap_or_else(|| panic!("Unable to compute the factorial of {}.", 2 * ly))
196        * BigUint::from(2 * lz)
197            .checked_factorial()
198            .unwrap_or_else(|| panic!("Unable to compute the factorial of {}.", 2 * lz)))
199    .to_f64()
200    .expect("Unable to convert a `BigUint` value to `f64`.")
201        * std::f64::consts::PI.powi(3).sqrt();
202    (num / den).sqrt()
203}
204
205/// Obtain the complex coefficients $`c(l, m_l, n, l_x, l_y, l_z)`$ based on Equation 15 of
206/// Schlegel, H. B. & Frisch, M. J. Transformation between Cartesian and pure spherical harmonic
207/// Gaussians. *International Journal of Quantum Chemistry* **54**, 83–87 (1995),
208/// [DOI](https://doi.org/10.1002/qua.560540202), but more generalised
209/// for $`l \leq l_{\mathrm{cart}} = l_x + l_y + l_z`$.
210///
211/// Let $`\tilde{g}(\alpha, l, m_l, l_{\mathrm{cart}}, \mathbf{r})`$ be a complex solid
212/// harmonic Gaussian as defined in Equation 1 of
213/// the above reference with $`n = l_{\mathrm{cart}}`$, and let
214/// $`g(\alpha, l_x, l_y, l_z, \mathbf{r})`$ be a Cartesian Gaussian as defined in
215/// Equation 2 of the above reference.
216/// The complex coefficients $`c(l, m_l, n, l_x, l_y, l_z)`$ effect the transformation
217///
218/// ```math
219/// \tilde{g}(\alpha, l, m_l, l_{\mathrm{cart}}, \mathbf{r})
220/// = \sum_{l_x+l_y+l_z=l_{\mathrm{cart}}}
221///     c(l, m_l, l_{\mathrm{cart}}, l_x, l_y, l_z)
222///     g(\alpha, l_x, l_y, l_z, \mathbf{r})
223/// ```
224///
225/// and are given by
226///
227/// ```math
228/// c(l, m_l, l_{\mathrm{cart}}, l_x, l_y, l_z)
229/// = \frac{\tilde{N}(l_{\mathrm{cart}}, \alpha)}{N(l_x, l_y, l_z, \alpha)}
230///     \tilde{c}(l, m_l, l_{\mathrm{cart}}, l_x, l_y, l_z).
231/// ```
232///
233/// The normalisation constants $`\tilde{N}(l_{\mathrm{cart}}, \alpha)`$
234/// and $`N(l_x, l_y, l_z, \alpha)`$ are given in Equations 8 and 9 of
235/// the above reference, and for $`n = l_{\mathrm{cart}}`$, this ratio turns out to be
236/// independent of $`\alpha`$.
237/// The more general form of $`\tilde{c}(l, m_l, l_{\mathrm{cart}}, l_x, l_y, l_z)`$ has been
238/// derived to be
239///
240/// ```math
241/// \tilde{c}(l, m_l, l_{\mathrm{cart}}, l_x, l_y, l_z)
242/// = \frac{\lambda_{\mathrm{cs}}}{2^l l!}
243///     \sqrt{\frac{(2l+1)(l-\lvert m_l \rvert)!}{4\pi(l+\lvert m_l \rvert)!}}
244///     \sum_{i=0}^{(l-\lvert m_l \rvert)/2}
245///         {l\choose i} \frac{(-1)^i(2l-2i)!}{(l-\lvert m_l \rvert -2i)!}\\
246///     \sum_{p=0}^{\lvert m_l \rvert} {{\lvert m_l \rvert} \choose p}
247///         (\pm \mathbb{i})^{\lvert m_l \rvert-p}
248///     \sum_{q=0}^{\Delta l/2} {{\Delta l/2} \choose q} {i \choose j_q}
249///     \sum_{k=0}^{j_q} {q \choose t_{pk}} {j_q \choose k}
250/// ```
251///
252/// where $`+\mathbb{i}`$ applies for $`m_l > 0`$, $`-\mathbb{i}`$
253/// for $`m_l \le 0`$, $`\lambda_{\mathrm{cs}}`$ is the Condon--Shortley
254/// phase given by
255///
256/// ```math
257/// \lambda_{\mathrm{cs}} =
258///     \begin{cases}
259///         (-1)^{m_l} & m_l > 0 \\
260///         1          & m_l \leq 0
261///     \end{cases}
262/// ```
263///
264/// and
265///
266/// ```math
267/// t_{pk} = \frac{l_x-p-2k}{2} \quad \textrm{and} \quad
268/// j_q = \frac{l_x+l_y-\lvert m_l \rvert-2q}{2}.
269/// ```
270///
271/// If $`\Delta l`$ is odd, $`\tilde{c}(l, m_l, l_{\mathrm{cart}}, l_x, l_y, l_z)`$ must vanish.
272/// When  $`t_{pk}`$ or $`j_q`$ is a half-integer, the inner sum in which it is involved
273/// evaluates to zero.
274///
275/// # Arguments
276///
277/// * `lpureqns` - A tuple of $`(l, m_l)`$ specifying the quantum numbers for the spherical
278///   harmonic component of the solid harmonic Gaussian.
279/// * `lcartqns` - A tuple of $`(l_x, l_y, l_z)`$ specifying the exponents of the Cartesian
280///   components of the Cartesian Gaussian.
281/// * `csphase` - If `true`, the Condon--Shortley phase will be used as defined above.
282///   If `false`, this phase will be set to unity.
283///
284/// # Returns
285///
286/// The complex factor $`c(l, m_l, l_{\mathrm{cart}}, l_x, l_y, l_z)`$.
287///
288/// # Panics
289///
290/// Panics when any required factorials cannot be computed.
291#[allow(clippy::too_many_lines)]
292pub fn complexc(lpureqns: (u32, i32), lcartqns: (u32, u32, u32), csphase: bool) -> Complex<f64> {
293    let (l, m) = lpureqns;
294    let li32 = i32::try_from(l).unwrap_or_else(|_| panic!("Cannot convert `{l}` to `i32`."));
295    assert!(
296        m.unsigned_abs() <= l,
297        "m must be between -l and l (inclusive)."
298    );
299    let (lx, ly, lz) = lcartqns;
300    let lxi32 = i32::try_from(lx).unwrap_or_else(|_| panic!("Cannot convert `{lx}` to `i32`."));
301    let lyi32 = i32::try_from(ly).unwrap_or_else(|_| panic!("Cannot convert `{ly}` to `i32`."));
302    let lzi32 = i32::try_from(lz).unwrap_or_else(|_| panic!("Cannot convert `{lz}` to `i32`."));
303    let lcart = lx + ly + lz;
304    let lcarti32 = lxi32 + lyi32 + lzi32;
305    let dl = lcarti32 - li32;
306    if dl % 2 != 0 {
307        return Complex::<f64>::zero();
308    }
309
310    let num = f64::from(
311        (2 * l + 1)
312            * (l - m.unsigned_abs())
313                .checked_factorial()
314                .unwrap_or_else(|| {
315                    panic!(
316                        "Unable to compute the factorial of {}.",
317                        l - m.unsigned_abs()
318                    )
319                }),
320    );
321    let den = 4.0
322        * std::f64::consts::PI
323        * f64::from(
324            (l + m.unsigned_abs())
325                .checked_factorial()
326                .unwrap_or_else(|| {
327                    panic!(
328                        "Unable to compute the factorial of {}.",
329                        l + m.unsigned_abs()
330                    )
331                }),
332        );
333    let mut prefactor =
334        1.0 / f64::from(
335            2u32.pow(l)
336                * l.checked_factorial()
337                    .unwrap_or_else(|| panic!("Unable to compute the factorial of {l}.")),
338        ) * (num / den).sqrt();
339    if csphase && m > 0 {
340        prefactor *=
341            f64::from((-1i32).pow(u32::try_from(m).expect("Unable to convert `m` to `u32`.")));
342    }
343    let ntilde = norm_sph_gaussian(lcart, 1.0);
344    let n = norm_cart_gaussian(lcartqns, 1.0);
345
346    let si =
347        (0..=((l - m.unsigned_abs()).div_euclid(2))).fold(Complex::<f64>::zero(), |acc_si, i| {
348            // i <= (l - |m|) / 2
349            let ii32 =
350                i32::try_from(i).unwrap_or_else(|_| panic!("Cannot convert `{i}` to `i32`."));
351            let mut ifactor = combu(l, i)
352                .to_f64()
353                .expect("Unable to convert a `BigUint` value to `f64`.")
354                * BigUint::from(2 * l - 2 * i)
355                    .checked_factorial()
356                    .unwrap_or_else(|| {
357                        panic!("Unable to compute the factorial of {}.", 2 * l - 2 * i)
358                    })
359                    .to_f64()
360                    .unwrap_or_else(|| {
361                        panic!(
362                            "Unable to convert the factorial of {} to `f64`.",
363                            2 * l - 2 * i
364                        )
365                    })
366                / BigUint::from(l - m.unsigned_abs() - 2 * i)
367                    .checked_factorial()
368                    .unwrap_or_else(|| {
369                        panic!(
370                            "Unable to compute the factorial of {}.",
371                            l - m.unsigned_abs() - 2 * i
372                        )
373                    })
374                    .to_f64()
375                    .unwrap_or_else(|| {
376                        panic!(
377                            "Unable to convert the factorial of {} to `f64`.",
378                            l - m.unsigned_abs() - 2 * i
379                        )
380                    });
381            if i % 2 == 1 {
382                ifactor *= -1.0;
383            };
384            let sp = (0..=(m.unsigned_abs())).fold(Complex::<f64>::zero(), |acc_sp, p| {
385                let pi32 =
386                    i32::try_from(p).unwrap_or_else(|_| panic!("Cannot convert `{p}` to `i32`."));
387                let pfactor = if m > 0 {
388                    combu(m.unsigned_abs(), p)
389                        .to_f64()
390                        .expect("Unable to convert a `BigUint` value to `f64`.")
391                        * Complex::<f64>::i().powu(m.unsigned_abs() - p)
392                } else {
393                    combu(m.unsigned_abs(), p)
394                        .to_f64()
395                        .expect("Unable to convert a `BigUint` value to `f64`.")
396                        * (-1.0 * Complex::<f64>::i()).powu(m.unsigned_abs() - p)
397                };
398                let sq = (0..=(dl.div_euclid(2))).fold(Complex::<f64>::zero(), |acc_sq, q| {
399                    let jq_num = lxi32 + lyi32 - 2 * q - m.abs();
400                    if jq_num.rem_euclid(2) == 0 {
401                        let jq = jq_num.div_euclid(2);
402                        let qfactor = (comb(dl.div_euclid(2), q) * comb(ii32, jq))
403                            .to_f64()
404                            .expect("Unable to convert a `BigUint` value to `f64`.");
405                        let sk = (0..=jq).fold(Complex::<f64>::zero(), |acc_sk, k| {
406                            let tpk_num = lxi32 - pi32 - 2 * k;
407                            if tpk_num.rem_euclid(2) == 0 {
408                                let tpk = tpk_num.div_euclid(2);
409                                let kfactor = (comb(q, tpk) * comb(jq, k))
410                                    .to_f64()
411                                    .expect("Unable to convert a `BigUint` value to `f64`.");
412                                acc_sk + kfactor
413                            } else {
414                                acc_sk
415                            }
416                        });
417                        acc_sq + qfactor * sk
418                    } else {
419                        acc_sq
420                    }
421                });
422                acc_sp + pfactor * sq
423            });
424            acc_si + ifactor * sp
425        });
426    (ntilde / n) * prefactor * si
427}
428
429/// Calculates the overlap between two normalised Cartesian Gaussians of the same order and radial
430/// width, as given in Equation 19 of Schlegel, H. B. & Frisch, M. J. Transformation between
431/// Cartesian and pure spherical harmonic Gaussians. *International Journal of Quantum Chemistry*
432/// **54**, 83–87 (1995), [DOI](https://doi.org/10.1002/qua.560540202).
433///
434/// # Arguments
435///
436/// * `lcartqns1` - A tuple of $`(l_x, l_y, l_z)`$ specifying the exponents of the Cartesian
437///   components of the first Cartesian Gaussian.
438/// * `lcartqns2` - A tuple of $`(l_x, l_y, l_z)`$ specifying the exponents of the Cartesian
439///   components of the first Cartesian Gaussian.
440///
441/// # Returns
442///
443/// The overlap between the two specified normalised Cartesian Gaussians.
444fn cartov(lcartqns1: (u32, u32, u32), lcartqns2: (u32, u32, u32)) -> f64 {
445    let (lx1, ly1, lz1) = lcartqns1;
446    let (lx2, ly2, lz2) = lcartqns2;
447    let lcart1 = lx1 + ly1 + lz1;
448    let lcart2 = lx2 + ly2 + lz2;
449    assert_eq!(
450        lcart1, lcart2,
451        "Only Cartesian Gaussians of the same order are supported."
452    );
453
454    if (lx1 + lx2).rem_euclid(2) == 0
455        && (ly1 + ly2).rem_euclid(2) == 0
456        && (lz1 + lz2).rem_euclid(2) == 0
457    {
458        let num1 = (BigUint::from(lx1 + lx2)
459            .checked_factorial()
460            .unwrap_or_else(|| panic!("Unable to compute the factorial of {}.", lx1 + lx2))
461            * BigUint::from(ly1 + ly2)
462                .checked_factorial()
463                .unwrap_or_else(|| panic!("Unable to compute the factorial of {}.", ly1 + ly2))
464            * BigUint::from(lz1 + lz2)
465                .checked_factorial()
466                .unwrap_or_else(|| panic!("Unable to compute the factorial of {}.", lz1 + lz2)))
467        .to_f64()
468        .expect("Unable to convert a `BigUint` value to `f64`.");
469
470        let den1 = (BigUint::from((lx1 + lx2).div_euclid(2))
471            .checked_factorial()
472            .unwrap_or_else(|| {
473                panic!(
474                    "Unable to compute the factorial of {}.",
475                    (lx1 + lx2).div_euclid(2)
476                )
477            })
478            * BigUint::from((ly1 + ly2).div_euclid(2))
479                .checked_factorial()
480                .unwrap_or_else(|| {
481                    panic!(
482                        "Unable to compute the factorial of {}.",
483                        (ly1 + ly2).div_euclid(2)
484                    )
485                })
486            * BigUint::from((lz1 + lz2).div_euclid(2))
487                .checked_factorial()
488                .unwrap_or_else(|| {
489                    panic!(
490                        "Unable to compute the factorial of {}.",
491                        (lz1 + lz2).div_euclid(2)
492                    )
493                }))
494        .to_f64()
495        .expect("Unable to convert a `BigUint` value to `f64`.");
496
497        let num2 = (BigUint::from(lx1)
498            .checked_factorial()
499            .unwrap_or_else(|| panic!("Unable to compute the factorial of {lx1}."))
500            * BigUint::from(ly1)
501                .checked_factorial()
502                .unwrap_or_else(|| panic!("Unable to compute the factorial of {ly1}."))
503            * BigUint::from(lz1)
504                .checked_factorial()
505                .unwrap_or_else(|| panic!("Unable to compute the factorial of {lz1}."))
506            * BigUint::from(lx2)
507                .checked_factorial()
508                .unwrap_or_else(|| panic!("Unable to compute the factorial of {lx2}."))
509            * BigUint::from(ly2)
510                .checked_factorial()
511                .unwrap_or_else(|| panic!("Unable to compute the factorial of {ly2}."))
512            * BigUint::from(lz2)
513                .checked_factorial()
514                .unwrap_or_else(|| panic!("Unable to compute the factorial of {lz2}.")))
515        .to_f64()
516        .expect("Unable to convert a `BigUint` value to `f64`.");
517
518        let den2 = (BigUint::from(2 * lx1)
519            .checked_factorial()
520            .unwrap_or_else(|| panic!("Unable to compute the factorial of {}.", 2 * lx1))
521            * BigUint::from(2 * ly1)
522                .checked_factorial()
523                .unwrap_or_else(|| panic!("Unable to compute the factorial of {}.", 2 * ly1))
524            * BigUint::from(2 * lz1)
525                .checked_factorial()
526                .unwrap_or_else(|| panic!("Unable to compute the factorial of {}.", 2 * lz1))
527            * BigUint::from(2 * lx2)
528                .checked_factorial()
529                .unwrap_or_else(|| panic!("Unable to compute the factorial of {}.", 2 * lx2))
530            * BigUint::from(2 * ly2)
531                .checked_factorial()
532                .unwrap_or_else(|| panic!("Unable to compute the factorial of {}.", 2 * ly2))
533            * BigUint::from(2 * lz2)
534                .checked_factorial()
535                .unwrap_or_else(|| panic!("Unable to compute the factorial of {}.", 2 * lz2)))
536        .to_f64()
537        .expect("Unable to convert a `BigUint` value to `f64`.");
538
539        (num1 / den1) * (num2 / den2).sqrt()
540    } else {
541        0.0
542    }
543}
544
545/// Computes the inverse complex coefficients $`c^{-1}(l_x, l_y, l_z, l, m_l, l_{\mathrm{cart}})`$
546/// based on Equation 18 of Schlegel, H. B. & Frisch, M. J. Transformation between
547/// Cartesian and pure spherical harmonic Gaussians. *International Journal of Quantum Chemistry*
548/// **54**, 83–87 (1995), [DOI](https://doi.org/10.1002/qua.560540202), but more generalised for
549/// $`l \leq l_{\mathrm{cart}} = l_x + l_y + l_z`$.
550///
551/// Let $`\tilde{g}(\alpha, l, m_l, l_{\mathrm{cart}}, \mathbf{r})`$ be a complex solid
552/// harmonic Gaussian as defined in Equation 1 of the above reference with
553/// $`n = l_{\mathrm{cart}}`$, and let $`g(\alpha, l_x, l_y, l_z, \mathbf{r})`$ be a Cartesian
554/// Gaussian as defined in Equation 2 of the above reference. The inverse complex coefficients
555/// $`c^{-1}(l_x, l_y, l_z, l, m_l, l_{\mathrm{cart}})`$ effect the inverse transformation
556///
557/// ```math
558/// g(\alpha, l_x, l_y, l_z, \mathbf{r})
559/// = \sum_{l \le l_{\mathrm{cart}} = l_x+l_y+l_z} \sum_{m_l = -l}^{l}
560///     c^{-1}(l_x, l_y, l_z, l, m_l, l_{\mathrm{cart}})
561///     \tilde{g}(\alpha, l, m_l, l_{\mathrm{cart}}, \mathbf{r}).
562/// ```
563///
564/// # Arguments
565///
566/// * `lcartqns` - A tuple of $`(l_x, l_y, l_z)`$ specifying the exponents of the Cartesian
567///   components of the Cartesian Gaussian.
568/// * `lpureqns` - A tuple of $`(l, m_l)`$ specifying the quantum numbers for the spherical
569///   harmonic component of the solid harmonic Gaussian.
570/// * `csphase` - If `true`, the Condon--Shortley phase will be used as defined in
571///   [`complexc`]. If `false`, this phase will be set to unity.
572///
573/// # Returns
574///
575/// $`c^{-1}(l_x, l_y, l_z, l, m_l, l_{\mathrm{cart}})`$.
576pub fn complexcinv(lcartqns: (u32, u32, u32), lpureqns: (u32, i32), csphase: bool) -> Complex<f64> {
577    let (lx, ly, lz) = lcartqns;
578    let lcart = lx + ly + lz;
579    let mut cinv = Complex::<f64>::zero();
580    for lx2 in 0..=lcart {
581        for ly2 in 0..=(lcart - lx2) {
582            let lz2 = lcart - lx2 - ly2;
583            cinv += cartov(lcartqns, (lx2, ly2, lz2))
584                * complexc(lpureqns, (lx2, ly2, lz2), csphase).conj();
585        }
586    }
587    cinv
588}
589
590/// Obtains the transformation matrix $`\boldsymbol{\Upsilon}^{(l)}`$ allowing complex spherical
591/// harmonics to be expressed as linear combinations of real spherical harmonics.
592///
593/// Let $`Y_{lm}`$ be a real spherical harmonic of degree $`l`$. Then, a complex spherical
594/// harmonic of degree $`l`$ and order $`m`$ is given by
595///
596/// ```math
597/// Y_l^m =
598///     \begin{cases}
599///         \frac{\lambda_{\mathrm{cs}}}{\sqrt{2}}
600///         \left(Y_{l\lvert m \rvert}
601///               - \mathbb{i} Y_{l,-\lvert m \rvert}\right)
602///         & \mathrm{if}\ m < 0 \\
603///         Y_{l0} & \mathrm{if}\ m = 0 \\
604///         \frac{\lambda_{\mathrm{cs}}}{\sqrt{2}}
605///         \left(Y_{l\lvert m \rvert}
606///               + \mathbb{i} Y_{l,-\lvert m \rvert}\right)
607///         & \mathrm{if}\ m > 0 \\
608///     \end{cases}
609/// ```
610///
611/// where $`\lambda_{\mathrm{cs}}`$ is the Condon--Shortley phase as defined in [`complexc`].
612/// The linear combination coefficients can then be gathered into a square matrix
613/// $`\boldsymbol{\Upsilon}^{(l)}`$ of dimensions $`(2l+1)\times(2l+1)`$ such that
614///
615/// ```math
616///     Y_l^m = \sum_{m'} Y_{lm'} \Upsilon^{(l)}_{m'm}.
617/// ```
618///
619/// # Arguments
620///
621/// * `l` - The spherical harmonic degree.
622/// * `csphase` - If `true`, $`\lambda_{\mathrm{cs}}`$ is as defined in [`complexc`]. If `false`,
623///   $`\lambda_{\mathrm{cs}} = 1`$.
624/// * `pureorder` - A [`PureOrder`] struct giving the ordering of the components of the pure
625///   Gaussians.
626///
627/// # Returns
628///
629/// The $`\boldsymbol{\Upsilon}^{(l)}`$ matrix.
630pub fn sh_c2r_mat(l: u32, csphase: bool, pureorder: &PureOrder) -> Array2<Complex<f64>> {
631    assert_eq!(pureorder.lpure, l, "Mismatched pure ranks.");
632    let lusize = l as usize;
633    let mut upmat = Array2::<Complex<f64>>::zeros((2 * lusize + 1, 2 * lusize + 1));
634    let po_il = PureOrder::increasingm(l);
635    for &mcomplex in po_il.iter() {
636        let absmreal = mcomplex.unsigned_abs() as usize;
637        match mcomplex.cmp(&0) {
638            Ordering::Less => {
639                // Python-equivalent:
640                // upmat[-absmreal + l, mcomplex + l] = -1.0j / np.sqrt(2)
641                // upmat[+absmreal + l, mcomplex + l] = 1.0 / np.sqrt(2)
642                // mcomplex = -absmreal
643                upmat[(lusize - absmreal, lusize - absmreal)] =
644                    Complex::<f64>::new(0.0, -1.0 / 2.0f64.sqrt());
645                upmat[(lusize + absmreal, lusize - absmreal)] =
646                    Complex::<f64>::new(1.0 / 2.0f64.sqrt(), 0.0);
647            }
648            Ordering::Equal => {
649                upmat[(lusize, lusize)] = Complex::<f64>::from(1.0);
650            }
651            Ordering::Greater => {
652                let lcs = if csphase {
653                    f64::from((-1i32).pow(
654                        u32::try_from(mcomplex).expect("Unable to convert `mcomplex` to `u32`."),
655                    ))
656                } else {
657                    1.0
658                };
659                // Python-equivalent:
660                // upmat[-absmreal + l, mcomplex + l] = lcs * 1.0j / np.sqrt(2)
661                // upmat[+absmreal + l, mcomplex + l] = lcs * 1.0 / np.sqrt(2)
662                // mcomplex = absmreal
663                upmat[(lusize - absmreal, lusize + absmreal)] =
664                    lcs * Complex::<f64>::new(0.0, 1.0 / 2.0f64.sqrt());
665                upmat[(lusize + absmreal, lusize + absmreal)] =
666                    lcs * Complex::<f64>::new(1.0 / 2.0f64.sqrt(), 0.0);
667            }
668        }
669    }
670
671    // upmat is always in increasing-m order. We now permute, if required.
672    if *pureorder != po_il {
673        let perm = pureorder.get_perm_of(&po_il).expect(
674            "Permutation to obtain `pureorder` from the increasing-m order could not be found.",
675        );
676        let image = perm.image();
677        upmat.select(Axis(0), image).select(Axis(1), image)
678    } else {
679        upmat
680    }
681}
682
683/// Obtains the matrix $`\boldsymbol{\Upsilon}^{(l)\dagger}`$ allowing real spherical harmonics
684/// to be expressed as linear combinations of complex spherical harmonics.
685///
686/// Let $`Y_l^m`$ be a complex spherical harmonic of degree $`l`$ and order $`m`$.
687/// Then, a real degree-$`l`$ spherical harmonic $`Y_{lm}`$ can be defined as
688///
689/// ```math
690/// Y_{lm} =
691///     \begin{cases}
692///         \frac{\mathbb{i}}{\sqrt{2}}
693///         \left(Y_l^{-\lvert m \rvert}
694///               - \lambda'_{\mathrm{cs}} Y_l^{\lvert m \rvert}\right)
695///         & \mathrm{if}\ m < 0 \\
696///         Y_l^0 & \mathrm{if}\ m = 0 \\
697///         \frac{1}{\sqrt{2}}
698///         \left(Y_l^{-\lvert m \rvert}
699///               + \lambda'_{\mathrm{cs}} Y_l^{\lvert m \rvert}\right)
700///         & \mathrm{if}\ m > 0 \\
701///     \end{cases}
702/// ```
703///
704/// where $`\lambda'_{\mathrm{cs}} = (-1)^{\lvert m \rvert}`$ if the Condon--Shortley phase as
705/// defined in [`complexc`] is employed for the complex spherical harmonics, and
706/// $`\lambda'_{\mathrm{cs}} = 1`$ otherwise. The linear combination coefficients turn out to be
707/// given by the elements of matrix $`\boldsymbol{\Upsilon}^{(l)\dagger}`$ of dimensions
708/// $`(2l+1)\times(2l+1)`$ such that
709///
710/// ```math
711///     Y_{lm} = \sum_{m'} Y_l^{m'} [\Upsilon^{(l)\dagger}]_{m'm}.
712/// ```
713///
714/// It is obvious from the orthonormality of $`Y_{lm}`$ and $`Y_l^m`$ that
715/// $`\boldsymbol{\Upsilon}^{(l)\dagger} = [\boldsymbol{\Upsilon}^{(l)}]^{-1}`$ where
716/// $`\boldsymbol{\Upsilon}^{(l)}`$ is defined in [`sh_c2r_mat`].
717///
718/// # Arguments
719///
720/// * `l` - The spherical harmonic degree.
721/// * `csphase` - If `true`, $`\lambda_{\mathrm{cs}}`$ is as defined in [`complexc`]. If `false`,
722///   $`\lambda_{\mathrm{cs}} = 1`$.
723/// * `pureorder` - A [`PureOrder`] struct giving the ordering of the components of the pure
724///   Gaussians.
725///
726/// # Returns
727///
728/// The $`\boldsymbol{\Upsilon}^{(l)\dagger}`$ matrix.
729pub fn sh_r2c_mat(l: u32, csphase: bool, pureorder: &PureOrder) -> Array2<Complex<f64>> {
730    let mut mat = sh_c2r_mat(l, csphase, pureorder).t().to_owned();
731    mat.par_mapv_inplace(|x| x.conj());
732    mat
733}
734
735/// Obtains the matrix $`\mathbf{U}^{(l_{\mathrm{cart}}, l)}`$ containing linear combination
736/// coefficients of Cartesian Gaussians in the expansion of a complex solid harmonic Gaussian,
737/// *i.e.*, briefly,
738///
739/// ```math
740/// \tilde{\mathbf{g}}^{\mathsf{T}}(l)
741///     = \mathbf{g}^{\mathsf{T}}(l_{\mathrm{cart}})
742///     \ \mathbf{U}^{(l_{\mathrm{cart}}, l)}.
743/// ```
744///
745/// Let $`\tilde{g}(\alpha, \lambda, l_{\mathrm{cart}}, \mathbf{r})`$ be a complex solid harmonic
746/// Gaussian as defined in Equation 1 of Schlegel, H. B. & Frisch, M. J. Transformation between
747/// Cartesian and pure spherical harmonic Gaussians. *International Journal of Quantum Chemistry*
748/// **54**, 83–87 (1995), [DOI](https://doi.org/10.1002/qua.560540202) with
749/// $`n = l_{\mathrm{cart}}`$, and let $`g(\alpha, \lambda_{\mathrm{cart}}, \mathbf{r})`$ be a
750/// Cartesian Gaussian as defined in Equation 2 of the above reference.
751/// Here, $`\lambda`$ is a single index labelling a complex solid harmonic Gaussian of spherical
752/// harmonic degree $`l`$ and order $`m_l`$, and $`\lambda_{\mathrm{cart}}`$ a single index
753/// labelling a Cartesian Gaussian of degrees $`(l_x, l_y, l_z)`$ such that
754/// $`l_x + l_y + l_z = l_{\mathrm{cart}}`$. We can then write
755///
756/// ```math
757/// \tilde{g}(\alpha, \lambda, l_{\mathrm{cart}}, \mathbf{r})
758/// = \sum_{\lambda_{\mathrm{cart}}}
759///     g(\alpha, \lambda_{\mathrm{cart}}, \mathbf{r})
760///     U^{(l_{\mathrm{cart}}, l)}_{\lambda_{\mathrm{cart}}\lambda}
761/// ```
762///
763/// where $`U^{(l_{\mathrm{cart}}, l)}_{\lambda_{\mathrm{cart}}\lambda}`$
764/// is given by the complex coefficients
765///
766/// ```math
767/// U^{(l_{\mathrm{cart}}, l)}_{\lambda_{\mathrm{cart}}\lambda} =
768///     c(l, m_l, l_{\mathrm{cart}}, l_x, l_y, l_z)
769/// ```
770///
771/// defined in [`complexc`].
772///
773/// $`\mathbf{U}^{(l_{\mathrm{cart}}, l)}`$ has dimensions
774/// $`\frac{1}{2}(l_{\mathrm{cart}}+1)(l_{\mathrm{cart}}+2) \times (2l+1)`$ and contains only
775/// zero elements if $`l`$ and $`l_{\mathrm{cart}}`$ have different parities.
776/// It can be verified that
777/// $`\mathbf{V}^{(l,l_{\mathrm{cart}})}
778/// \ \mathbf{U}^{(l_{\mathrm{cart}}, l)} = \boldsymbol{I}_{2l+1}`$, where
779/// $`\mathbf{V}^{(l,l_{\mathrm{cart}})}`$ is given in [`sh_cart2cl_mat`].
780///
781/// # Arguments
782///
783/// * `lcart` - The total Cartesian degree for the Cartesian Gaussians and
784///   also for the radial part of the solid harmonic Gaussian.
785/// * `l` - The degree of the complex spherical harmonic factor in the solid
786///   harmonic Gaussian.
787/// * `cartorder` - A [`CartOrder`] struct giving the ordering of the components of the Cartesian
788///   Gaussians.
789/// * `csphase` - Set to `true` to use the Condon--Shortley phase in the calculations of the $`c`$
790///   coefficients. See [`complexc`] for more details.
791/// * `pureorder` - A [`PureOrder`] struct giving the ordering of the components of the pure
792///   Gaussians.
793///
794/// # Returns
795///
796/// The $`\mathbf{U}^{(l_{\mathrm{cart}}, l)}`$ matrix.
797pub fn sh_cl2cart_mat(
798    lcart: u32,
799    l: u32,
800    cartorder: &CartOrder,
801    csphase: bool,
802    pureorder: &PureOrder,
803) -> Array2<Complex<f64>> {
804    assert_eq!(cartorder.lcart, lcart, "Mismatched Cartesian ranks.");
805    assert_eq!(pureorder.lpure, l, "Mismatched pure ranks.");
806    let mut umat = Array2::<Complex<f64>>::zeros((
807        ((lcart + 1) * (lcart + 2)).div_euclid(2) as usize,
808        2 * l as usize + 1,
809    ));
810    for (i, &m) in pureorder.iter().enumerate() {
811        for (icart, &lcartqns) in cartorder.iter().enumerate() {
812            umat[(icart, i)] = complexc((l, m), lcartqns, csphase);
813        }
814    }
815    umat
816}
817
818/// Obtains the matrix $`\mathbf{V}^{(l, l_{\mathrm{cart}})}`$ containing linear combination
819/// coefficients of complex solid harmonic Gaussians of a specific degree in the expansion of
820/// Cartesian Gaussians, *i.e.*, briefly,
821///
822/// ```math
823/// \mathbf{g}^{\mathsf{T}}(l_{\mathrm{cart}})
824///     = \tilde{\mathbf{g}}^{\mathsf{T}}(l)
825///     \ \mathbf{V}^{(l, l_{\mathrm{cart}})}.
826/// ```
827///
828/// Let $`\tilde{g}(\alpha, \lambda, l_{\mathrm{cart}}, \mathbf{r})`$ be a complex solid harmonic
829/// Gaussian as defined in Equation 1 of Schlegel, H. B. & Frisch, M. J. Transformation between
830/// Cartesian and pure spherical harmonic Gaussians. *International Journal of Quantum Chemistry*
831/// **54**, 83–87 (1995), [DOI](https://doi.org/10.1002/qua.560540202) with
832/// $`n = l_{\mathrm{cart}}`$, and let $`g(\alpha, \lambda_{\mathrm{cart}}, \mathbf{r})`$ be a
833/// Cartesian Gaussian as defined in Equation 2 of the above reference.  Here, $`\lambda`$ is a
834/// single index labelling a complex solid harmonic Gaussian of spherical harmonic degree $`l`$
835/// and order $`m_l`$, and $`\lambda_{\mathrm{cart}}`$ a single index labelling a Cartesian
836/// Gaussian of degrees $`(l_x, l_y, l_z)`$ such that $`l_x + l_y + l_z = l_{\mathrm{cart}}`$.
837/// We can then write
838///
839/// ```math
840/// g(\alpha, \lambda_{\mathrm{cart}}, \mathbf{r})
841/// = \sum_{\substack{\lambda\\ l \leq l_{\mathrm{cart}}}}
842///     \tilde{g}(\alpha, \lambda, l_{\mathrm{cart}}, \mathbf{r})
843///     V^{(l_{\mathrm{cart}})}_{\lambda\lambda_{\mathrm{cart}}}
844/// ```
845///
846/// where $`V^{(l_{\mathrm{cart}})}_{\lambda\lambda_{\mathrm{cart}}}`$ is given by the inverse
847/// complex coefficients
848///
849/// ```math
850/// V^{(l_{\mathrm{cart}})}_{\lambda\lambda_{\mathrm{cart}}} =
851///     c^{-1}(l_x, l_y, l_z, l, m_l, l_{\mathrm{cart}})
852/// ```
853///
854/// defined in [`complexcinv`].
855///
856/// We can order the rows $`\lambda`$ of $`\mathbf{V}^{(l_{\mathrm{cart}})}`$ that have the same
857/// $`l`$ into rectangular blocks of dimensions
858/// $`(2l+1) \times \frac{1}{2}(l_{\mathrm{cart}}+1)(l_{\mathrm{cart}}+2)`$
859/// which give contributions from complex solid harmonic Gaussians of a particular degree $`l`$.
860/// We denote these blocks $`\mathbf{V}^{(l, l_{\mathrm{cart}})}`$.
861/// They contain only zero elements if $`l`$ and $`l_{\mathrm{cart}}`$ have different parities.
862///
863/// # Arguments
864///
865/// * `l` - The degree of the complex spherical harmonic factor in the solid
866///   harmonic Gaussian.
867/// * `lcart` - The total Cartesian degree for the Cartesian Gaussians and
868///   also for the radial part of the solid harmonic Gaussian.
869/// * `cartorder` - A [`CartOrder`] struct giving the ordering of the components of the Cartesian
870///   Gaussians.
871/// * `csphase` - Set to `true` to use the Condon--Shortley phase in the calculations of the
872///   $`c^{-1}`$ coefficients. See [`complexc`] and [`complexcinv`] for more details.
873/// * `pureorder` - A [`PureOrder`] struct giving the ordering of the components of the pure
874///   Gaussians.
875///
876/// # Returns
877///
878/// The $`\mathbf{V}^{(l, l_{\mathrm{cart}})}`$ block.
879pub fn sh_cart2cl_mat(
880    l: u32,
881    lcart: u32,
882    cartorder: &CartOrder,
883    csphase: bool,
884    pureorder: &PureOrder,
885) -> Array2<Complex<f64>> {
886    assert_eq!(pureorder.lpure, l, "Mismatched pure ranks.");
887    assert_eq!(cartorder.lcart, lcart, "Mismatched Cartesian ranks.");
888    let mut vmat = Array2::<Complex<f64>>::zeros((
889        2 * l as usize + 1,
890        ((lcart + 1) * (lcart + 2)).div_euclid(2) as usize,
891    ));
892    for (icart, &lcartqns) in cartorder.iter().enumerate() {
893        for (i, &m) in pureorder.iter().enumerate() {
894            vmat[(i, icart)] = complexcinv(lcartqns, (l, m), csphase);
895        }
896    }
897    vmat
898}
899
900/// Obtain the matrix $`\mathbf{W}^{(l_{\mathrm{cart}}, l)}`$ containing linear combination
901/// coefficients of Cartesian Gaussians in the expansion of a real solid harmonic Gaussian, *i.e.*,
902/// briefly,
903///
904/// ```math
905/// \bar{\mathbf{g}}^{\mathsf{T}}(l)
906///     = \mathbf{g}^{\mathsf{T}}(l_{\mathrm{cart}})
907///     \ \mathbf{W}^{(l_{\mathrm{cart}}, l)}.
908/// ```
909///
910/// Let $`\bar{g}(\alpha, \lambda, l_{\mathrm{cart}}, \mathbf{r})`$ be
911/// a real solid harmonic Gaussian defined in a similar manner to Equation 1 of Schlegel, H. B.
912/// & Frisch, M. J. Transformation between Cartesian and pure spherical harmonic Gaussians.
913/// *International Journal of Quantum Chemistry* **54**, 83–87 (1995),
914/// [DOI](https://doi.org/10.1002/qua.560540202) with $`n = l_{\mathrm{cart}}`$ but with real
915/// rather than complex spherical harmonic factors, and let
916/// $`g(\alpha, \lambda_{\mathrm{cart}}, \mathbf{r})`$ be a Cartesian Gaussian as defined in
917/// Equation 2 of the above reference. Here, $`\lambda`$ is a single index labelling a complex
918/// solid harmonic Gaussian of spherical harmonic degree $`l`$ and order $`m_l`$, and
919/// $`\lambda_{\mathrm{cart}}`$ a single index labelling a Cartesian Gaussian of degrees
920/// $`(l_x, l_y, l_z)`$ such that $`l_x + l_y + l_z = l_{\mathrm{cart}}`$. We can then write
921///
922/// ```math
923/// \bar{g}(\alpha, \lambda, l_{\mathrm{cart}}, \mathbf{r})
924/// = \sum_{\lambda_{\mathrm{cart}}}
925///     g(\alpha, \lambda_{\mathrm{cart}}, \mathbf{r})
926///     W^{(l_{\mathrm{cart}}, l)}_{\lambda_{\mathrm{cart}}\lambda}.
927/// ```
928///
929/// $`\mathbf{W}^{(l_{\mathrm{cart}}, l)}`$ is given by
930///
931/// ```math
932/// \mathbf{W}^{(l_{\mathrm{cart}}, l)}
933/// = \mathbf{U}^{(l_{\mathrm{cart}}, l)}
934///   \boldsymbol{\Upsilon}^{(l)\dagger},
935/// ```
936///
937/// where $`\boldsymbol{\Upsilon}^{(l)\dagger}`$ is defined in [`sh_r2c_mat`] and
938/// $`\mathbf{U}^{(l_{\mathrm{cart}}, l)}`$ in [`sh_cl2cart_mat`].
939/// $`\mathbf{W}^{(l_{\mathrm{cart}}, l)}`$ must be real.
940/// $`\mathbf{W}^{(l_{\mathrm{cart}}, l)}`$ has dimensions
941/// $`\frac{1}{2}(l_{\mathrm{cart}}+1)(l_{\mathrm{cart}}+2) \times (2l+1)`$ and contains only zero
942/// elements if $`l`$ and $`l_{\mathrm{cart}}`$ have different parities. It can be verified that
943/// $`\mathbf{X}^{(l,l_{\mathrm{cart}})}
944/// \ \mathbf{W}^{(l_{\mathrm{cart}}, l)} = \boldsymbol{I}_{2l+1}`$, where
945/// $`\mathbf{X}^{(l,l_{\mathrm{cart}})}`$ is given in
946/// [`sh_cart2rl_mat`].
947///
948/// # Arguments
949///
950/// * `lcart` - The total Cartesian degree for the Cartesian Gaussians and
951///   also for the radial part of the solid harmonic Gaussian.
952/// * `l` - The degree of the complex spherical harmonic factor in the solid
953///   harmonic Gaussian.
954/// * `cartorder` - A [`CartOrder`] struct giving the ordering of the components of the Cartesian
955///   Gaussians.
956/// * `csphase` - Set to `true` to use the Condon--Shortley phase in the calculations of the $`c`$
957///   coefficients. See [`complexc`] for more details.
958/// * `pureorder` - A [`PureOrder`] struct giving the ordering of the components of the pure
959///   Gaussians.
960///
961/// # Returns
962///
963/// The $`\mathbf{W}^{(l_{\mathrm{cart}}, l)}`$ matrix.
964pub fn sh_rl2cart_mat(
965    lcart: u32,
966    l: u32,
967    cartorder: &CartOrder,
968    csphase: bool,
969    pureorder: &PureOrder,
970) -> Array2<f64> {
971    assert_eq!(cartorder.lcart, lcart, "Mismatched Cartesian ranks.");
972    assert_eq!(pureorder.lpure, l, "Mismatched pure ranks.");
973    let upmatdagger = sh_r2c_mat(l, csphase, pureorder);
974    let umat = sh_cl2cart_mat(lcart, l, cartorder, csphase, pureorder);
975    let wmat = umat.dot(&upmatdagger);
976    assert!(
977        wmat.iter()
978            .all(|x| approx::relative_eq!(x.im, 0.0, max_relative = 1e-7, epsilon = 1e-7)),
979        "wmat is not entirely real."
980    );
981    wmat.map(|x| x.re)
982}
983
984/// Obtains the real matrix $`\mathbf{X}^{(l, l_{\mathrm{cart}})}`$ containing linear combination
985/// coefficients of real solid harmonic Gaussians of a specific degree in the expansion of
986/// Cartesian Gaussians, *i.e.*, briefly,
987///
988/// ```math
989/// \mathbf{g}^{\mathsf{T}}(l_{\mathrm{cart}})
990///     = \bar{\mathbf{g}}^{\mathsf{T}}(l)
991///     \ \mathbf{X}^{(l, l_{\mathrm{cart}})}.
992/// ```
993///
994/// Let $`\bar{g}(\alpha, \lambda, l_{\mathrm{cart}}, \mathbf{r})`$ be a real solid harmonic
995/// Gaussian defined in a similar manner to Equation 1 of Schlegel, H. B. & Frisch, M. J.
996/// Transformation between Cartesian and pure spherical harmonic Gaussians. *International
997/// Journal of Quantum Chemistry* **54**, 83–87 (1995),
998/// [DOI](https://doi.org/10.1002/qua.560540202)
999/// with $`n = l_{\mathrm{cart}}`$, but with real rather than complex spherical harmonic factors,
1000/// and let $`g(\alpha, \lambda_{\mathrm{cart}}, \mathbf{r})`$ be a Cartesian Gaussian as defined
1001/// in Equation 2 of the above reference.  Here, $`\lambda`$ is a single index labelling a real
1002/// solid harmonic Gaussian of spherical harmonic degree $`l`$ and real order $`m_l`$, and
1003/// $`\lambda_{\mathrm{cart}}`$ a single index labelling a Cartesian Gaussian of degrees
1004/// $`(l_x, l_y, l_z)`$ such that $`l_x + l_y + l_z = l_{\mathrm{cart}}`$.
1005/// We can then write
1006///
1007/// ```math
1008/// g(\alpha, \lambda_{\mathrm{cart}}, \mathbf{r})
1009/// = \sum_{\substack{\lambda\\ l \leq l_{\mathrm{cart}}}}
1010///     \bar{g}(\alpha, \lambda, l_{\mathrm{cart}}, \mathbf{r})
1011///     X^{(l_{\mathrm{cart}})}_{\lambda\lambda_{\mathrm{cart}}}.
1012/// ```
1013///
1014/// We can order the rows $`\lambda`$ of $`\mathbf{X}^{(l_{\mathrm{cart}})}`$ that have the same
1015/// $`l`$ into rectangular blocks of dimensions
1016/// $`(2l+1) \times \frac{1}{2}(l_{\mathrm{cart}}+1)(l_{\mathrm{cart}}+2)`$.
1017/// We denote these blocks $`\mathbf{X}^{(l, l_{\mathrm{cart}})}`$ which are given by
1018///
1019/// ```math
1020/// \mathbf{X}^{(l, l_{\mathrm{cart}})}
1021/// = \boldsymbol{\Upsilon}^{(l)} \mathbf{V}^{(l, l_{\mathrm{cart}})},
1022/// ```
1023///
1024/// where $`\boldsymbol{\Upsilon}^{(l)}`$ is defined in
1025/// [`sh_c2r_mat`] and $`\boldsymbol{V}^{(l, l_{\mathrm{cart}})}`$ in [`sh_cart2cl_mat`].
1026/// $`\mathbf{X}^{(l, l_{\mathrm{cart}})}`$ must be real.
1027///
1028/// # Arguments
1029///
1030/// * `l` - The degree of the complex spherical harmonic factor in the solid
1031///   harmonic Gaussian.
1032/// * `lcart` - The total Cartesian degree for the Cartesian Gaussians and
1033///   also for the radial part of the solid harmonic Gaussian.
1034/// * `cartorder` - A [`CartOrder`] struct giving the ordering of the components of the Cartesian
1035///   Gaussians.
1036/// * `csphase` - Set to `true` to use the Condon--Shortley phase in the calculations of the
1037///   $`c^{-1}`$ coefficients. See [`complexc`] and [`complexcinv`] for more details.
1038/// * `pureorder` - A [`PureOrder`] struct giving the ordering of the components of the pure
1039///   Gaussians.
1040///
1041/// # Returns
1042///
1043/// The $`\mathbf{X}^{(l, l_{\mathrm{cart}})}`$ block.
1044pub fn sh_cart2rl_mat(
1045    l: u32,
1046    lcart: u32,
1047    cartorder: &CartOrder,
1048    csphase: bool,
1049    pureorder: &PureOrder,
1050) -> Array2<f64> {
1051    assert_eq!(cartorder.lcart, lcart, "Mismatched Cartesian ranks.");
1052    assert_eq!(pureorder.lpure, l, "Mismatched pure ranks.");
1053    let upmat = sh_c2r_mat(l, csphase, pureorder);
1054    let vmat = sh_cart2cl_mat(l, lcart, cartorder, csphase, pureorder);
1055    let xmat = upmat.dot(&vmat);
1056    assert!(
1057        xmat.iter()
1058            .all(|x| approx::relative_eq!(x.im, 0.0, max_relative = 1e-7, epsilon = 1e-7)),
1059        "xmat is not entirely real."
1060    );
1061    xmat.map(|x| x.re)
1062}
1063
1064/// Returns a list of $`\mathbf{W}^{(l_{\mathrm{cart}}, l)}`$ for
1065/// $`l_{\mathrm{cart}} \ge l \ge 0`$ and $`l \equiv l_{\mathrm{cart}} \mod 2`$.
1066///
1067/// $`\mathbf{W}^{(l_{\mathrm{cart}}, l)}`$ is defined in [`sh_rl2cart_mat`].
1068///
1069/// # Arguments
1070///
1071/// * `lcart` - The total Cartesian degree for the Cartesian Gaussians and
1072///   also for the radial part of the solid harmonic Gaussian.
1073/// * `cartorder` - A [`CartOrder`] struct giving the ordering of the components of the Cartesian
1074///   Gaussians.
1075/// * `csphase` - Set to `true` to use the Condon--Shortley phase in the calculations of the
1076///   $`c`$ coefficients. See [`complexc`] for more details.
1077/// * `pureorder` - A closure to generate a [`PureOrder`] struct giving the ordering of the
1078///   components of the pure Gaussians for a particular value of `l`.
1079///
1080/// # Returns
1081///
1082/// A vector of $`\mathbf{W}^{(l_{\mathrm{cart}}, l)}`$ matrices with
1083/// $`l_{\mathrm{cart}} \ge l \ge 0`$ and $`l \equiv l_{\mathrm{cart}} \mod 2`$ in decreasing
1084/// $`l`$ order.
1085pub fn sh_r2cart(
1086    lcart: u32,
1087    cartorder: &CartOrder,
1088    csphase: bool,
1089    pureorder: fn(u32) -> PureOrder,
1090) -> Vec<Array2<f64>> {
1091    assert_eq!(cartorder.lcart, lcart, "Mismatched Cartesian ranks.");
1092    let lrange = if lcart.rem_euclid(2) == 0 {
1093        #[allow(clippy::range_plus_one)]
1094        (0..lcart + 1).step_by(2).rev()
1095    } else {
1096        #[allow(clippy::range_plus_one)]
1097        (1..lcart + 1).step_by(2).rev()
1098    };
1099    lrange
1100        .map(|l| sh_rl2cart_mat(lcart, l, cartorder, csphase, &pureorder(l)))
1101        .collect()
1102}
1103
1104/// Returns a list of $`\mathbf{X}^{(l, l_{\mathrm{cart}})}`$ for
1105/// $`l_{\mathrm{cart}} \ge l \ge 0`$ and $`l \equiv l_{\mathrm{cart}} \mod 2`$.
1106///
1107/// $`\mathbf{X}^{(l, l_{\mathrm{cart}})}`$ is defined in [`sh_cart2rl_mat`].
1108///
1109/// # Arguments
1110///
1111/// * `lcart` - The total Cartesian degree for the Cartesian Gaussians and
1112///   also for the radial part of the solid harmonic Gaussian.
1113/// * `cartorder` - A [`CartOrder`] struct giving the ordering of the components of the Cartesian
1114///   Gaussians.
1115/// * `csphase` - Set to `true` to use the Condon--Shortley phase in the calculations of the
1116///   $`c^{-1}`$ coefficients. See [`complexc`] and [`complexcinv`] for more details.
1117/// * `pureorder` - A closure to generate a [`PureOrder`] struct giving the ordering of the
1118///   components of the pure Gaussians for a particular value of `l`.
1119///
1120/// # Returns
1121///
1122/// A vector of $`\mathbf{X}^{(l, l_{\mathrm{cart}})}`$ matrices with
1123/// $`l_{\mathrm{cart}} \ge l \ge 0`$ and $`l \equiv l_{\mathrm{cart}} \mod 2`$ in decreasing
1124/// $`l`$ order.
1125pub fn sh_cart2r(
1126    lcart: u32,
1127    cartorder: &CartOrder,
1128    csphase: bool,
1129    pureorder: fn(u32) -> PureOrder,
1130) -> Vec<Array2<f64>> {
1131    assert_eq!(cartorder.lcart, lcart, "Mismatched Cartesian ranks.");
1132    let lrange = if lcart.rem_euclid(2) == 0 {
1133        #[allow(clippy::range_plus_one)]
1134        (0..lcart + 1).step_by(2).rev()
1135    } else {
1136        #[allow(clippy::range_plus_one)]
1137        (1..lcart + 1).step_by(2).rev()
1138    };
1139    lrange
1140        .map(|l| sh_cart2rl_mat(l, lcart, cartorder, csphase, &pureorder(l)))
1141        .collect()
1142}