/
boys.h
1980 lines (1749 loc) · 80.7 KB
/
boys.h
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
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/*
* Copyright (C) 2004-2020 Edward F. Valeev
*
* This file is part of Libint.
*
* Libint is free software: you can redistribute it and/or modify
* it under the terms of the GNU Lesser General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* Libint is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Lesser General Public License for more details.
*
* You should have received a copy of the GNU Lesser General Public License
* along with Libint. If not, see <http://www.gnu.org/licenses/>.
*
*/
#ifndef _libint2_src_lib_libint_boys_h_
#define _libint2_src_lib_libint_boys_h_
#if defined(__cplusplus)
#include <iostream>
#include <cstdlib>
#include <cmath>
#include <stdexcept>
#include <libint2/util/vector.h>
#include <cassert>
#include <vector>
#include <algorithm>
#include <limits>
#include <mutex>
#include <type_traits>
#include <memory>
// from now on at least C++11 is required by default
#include <libint2/util/cxxstd.h>
#if LIBINT2_CPLUSPLUS_STD < 2011
# error "Libint2 C++ API requires C++11 support"
#endif
#include <libint2/boys_fwd.h>
#include <libint2/numeric.h>
#include <libint2/initialize.h>
#if HAVE_LAPACK // use F77-type interface for now, switch to LAPACKE later
extern "C" void dgesv_(const int* n,
const int* nrhs, double* A, const int* lda,
int* ipiv, double* b, const int* ldb,
int* info);
#endif
namespace libint2 {
/// holds tables of expensive quantities
template<typename Real>
class ExpensiveNumbers {
public:
ExpensiveNumbers(int ifac = -1, int idf = -1, int ibc = -1) {
if (ifac >= 0) {
fac.resize(ifac + 1);
fac[0] = 1.0;
for (int i = 1; i <= ifac; i++) {
fac[i] = i * fac[i - 1];
}
}
if (idf >= 0) {
df.resize(idf + 1);
/* df[i] gives (i-1)!!, so that (-1)!! is defined... */
df[0] = 1.0;
if (idf >= 1)
df[1] = 1.0;
if (idf >= 2)
df[2] = 1.0;
for (int i = 3; i <= idf; i++) {
df[i] = (i - 1) * df[i - 2];
}
}
if (ibc >= 0) {
bc_.resize((ibc+1)*(ibc+1));
std::fill(bc_.begin(), bc_.end(), Real(0));
bc.resize(ibc+1);
bc[0] = &bc_[0];
for(int i=1; i<=ibc; ++i)
bc[i] = bc[i-1] + (ibc+1);
for(int i=0; i<=ibc; i++)
bc[i][0] = 1.0;
for(int i=0; i<=ibc; i++)
for(int j=1; j<=i; ++j)
bc[i][j] = bc[i][j-1] * Real(i-j+1) / Real(j);
}
for (int i = 0; i < 128; i++) {
twoi1[i] = 1.0 / (Real(2.0) * i + Real(1.0));
ihalf[i] = Real(i) - Real(0.5);
}
}
~ExpensiveNumbers() {
}
std::vector<Real> fac;
std::vector<Real> df;
std::vector<Real*> bc;
// these quantitites are needed with indices <= mmax
// 64 is sufficient to handle up to 4 center integrals with up to L=15 basis functions
// but need higher values for Yukawa integrals ...
Real twoi1[128]; /* 1/(2 i + 1); needed for downward recursion */
Real ihalf[128]; /* i - 0.5, needed for upward recursion */
private:
std::vector<Real> bc_;
};
#define _local_min_macro(a,b) ((a) > (b) ? (a) : (b))
/** Computes the Boys function, \f$ F_m (T) = \int_0^1 u^{2m} \exp(-T u^2) \, {\rm d}u \f$,
* using single algorithm (asymptotic expansion). Slow for the sake of precision control.
* Useful in two cases:
* <ul>
* <li> for reference purposes, if \c Real supports high/arbitrary precision, and </li>
* <li> for moderate values of \f$ T \f$, if \c Real is a low-precision floating-point type.
* N.B. FmEval_Reference2 , which can compute for all practical values of \f$ T \f$ and \f$ m \f$, is recommended
* with standard \c Real types (\c double and \c float). </li>
* </ul>
*
* \note Precision is controlled heuristically, i.e. cannot be guaranteed mathematically;
* will stop if absolute precision is reached, or precision of \c Real is exhausted.
* It is important that \c std::numeric_limits<Real> is defined appropriately.
*
* @tparam Real the type to use for all floating-point computations.
* Must be able to compute logarithm and exponential, i.e.
* log(x) and exp(x), where x is Real, must be valid expressions.
*/
template<typename Real>
struct FmEval_Reference {
static std::shared_ptr<const FmEval_Reference> instance(int /* mmax */, Real /* precision */) {
// thread-safe per C++11 standard [6.7.4]
static auto instance_ = std::make_shared<const FmEval_Reference>();
return instance_;
}
/// computes a single value of \f$ F_m(T) \f$ using MacLaurin series to full precision of @c Real
static Real eval(Real T, size_t m) {
assert(m < 100);
static const Real T_crit = std::numeric_limits<Real>::is_bounded == true ? -log( std::numeric_limits<Real>::min() * 100.5 / 2. ) : Real(0) ;
if (std::numeric_limits<Real>::is_bounded && T > T_crit)
throw std::overflow_error("FmEval_Reference<Real>::eval: Real lacks precision for the given value of argument T");
static const Real half = Real(1)/2;
Real denom = (m + half);
using std::exp;
Real term = exp(-T) / (2 * denom);
Real old_term = 0;
Real sum = term;
const Real epsilon = get_epsilon(T);
const Real epsilon_divided_10 = epsilon / 10;
do {
denom += 1;
old_term = term;
term = old_term * T / denom;
sum += term;
//rel_error = term / sum , hence iterate until rel_error = epsilon
// however, must ensure that contributions are decreasing to ensure that omitted contributions are smaller than epsilon
} while (term > sum * epsilon_divided_10 || old_term < term);
return sum;
}
/// fills up an array of Fm(T) for m in [0,mmax]
/// @param[out] Fm array to be filled in with the Boys function values, must be at least mmax+1 elements long
/// @param[in] T the Boys function argument
/// @param[in] mmax the maximum value of m for which Boys function will be computed;
static void eval(Real* Fm, Real T, size_t mmax) {
// evaluate for mmax using MacLaurin series
// it converges fastest for the largest m -> use it to compute Fmmax(T)
// see JPC 94, 5564 (1990).
for(size_t m=0; m<=mmax; ++m)
Fm[m] = eval(T, m);
return;
}
};
/** Computes the Boys function, \$ F_m (T) = \int_0^1 u^{2m} \exp(-T u^2) \, {\rm d}u \$,
* using multi-algorithm approach (upward recursion for T>=117, and asymptotic summation for T<117).
* This is slow and should be used for reference purposes, e.g. computing the interpolation tables.
* Precision is not always guaranteed as it is limited by the precision of \c Real type.
* When \c Real is \c double, can maintain absolute precision of epsilon for up to m=40.
*
* @tparam Real the type to use for all floating-point computations.
* Must be able to compute logarithm, exponential, square root, and error function, i.e.
* log(x), exp(x), sqrt(x), and erf(x), where x is Real, must be valid expressions.
*/
template<typename Real>
struct FmEval_Reference2 {
static std::shared_ptr<const FmEval_Reference2> instance(int /* mmax */, Real /* precision */) {
// thread-safe per C++11 standard [6.7.4]
static auto instance_ = std::make_shared<const FmEval_Reference2>();
return instance_;
}
/// fills up an array of Fm(T) for m in [0,mmax]
/// @param[out] Fm array to be filled in with the Boys function values, must be at least mmax+1 elements long
/// @param[in] t the Boys function argument
/// @param[in] mmax the maximum value of m for which Boys function will be computed;
static void eval(Real* Fm, Real t, size_t mmax) {
if (t < Real(117)) {
FmEval_Reference<Real>::eval(Fm,t,mmax);
}
else {
const Real two_over_sqrt_pi{1.128379167095512573896158903121545171688101258657997713688171443421284936882986828973487320404214727};
const Real K = 1/two_over_sqrt_pi;
auto t2 = 2*t;
auto et = exp(-t);
auto sqrt_t = sqrt(t);
Fm[0] = K*erf(sqrt_t)/sqrt_t;
if (mmax > 0)
for(size_t m=0; m<=mmax-1; m++) {
Fm[m+1] = ((2*m + 1)*Fm[m] - et)/(t2);
}
}
}
};
/** Computes the Boys function, \$ F_m (T) = \int_0^1 u^{2m} \exp(-T u^2) \, {\rm d}u \$,
* using 7-th order Chebyshev interpolation.
*/
template <typename Real = double>
class FmEval_Chebyshev7 {
#include <libint2/boys_cheb7.h>
static_assert(std::is_same<Real,double>::value, "FmEval_Chebyshev7 only supports double as the real type");
static constexpr const int ORDER = interpolation_order; //!, interpolation order
static constexpr const int ORDERp1 = ORDER+1; //!< ORDER + 1
static constexpr const Real T_crit = cheb_table_tmax; //!< critical value of T above which safe to use upward recusion
static constexpr Real delta = cheb_table_delta; //!< interval size
static constexpr Real one_over_delta = 1/delta; //! 1/delta
int mmax; //!< the maximum m that is tabulated
ExpensiveNumbers<double> numbers_;
Real *c; /* the Chebyshev coefficients table, T_crit*one_over_delta by mmax*ORDERp1 */
public:
/// \param m_max maximum value of the Boys function index; set to -1 to skip initialization
/// \param precision the desired relative precision
/// \throw std::invalid_argument if \c m_max is greater than \c cheb_table_mmax (see boys_cheb7.h)
/// \throw std::invalid_argument if \c precision is smaller than std::numeric_limits<double>::epsilon()
FmEval_Chebyshev7(int m_max, double precision = std::numeric_limits<double>::epsilon()) :
mmax(m_max), numbers_(14) {
#if defined(__x86_64__) || defined(_M_X64) || defined(__i386) || defined(_M_IX86)
# if !defined(__AVX__) && defined(NDEBUG)
if (libint2::verbose()) {
static bool printed_performance_warning = false;
if (!printed_performance_warning) {
libint2::verbose_stream()
<< "libint2::FmEval_Chebyshev7 on x86(-64) platforms needs AVX support for best performance"
<< std::endl;
printed_performance_warning = true;
}
}
# endif
#endif
if (precision < std::numeric_limits<double>::epsilon())
throw std::invalid_argument(std::string("FmEval_Chebyshev7 does not support precision smaller than ") + std::to_string(std::numeric_limits<double>::epsilon()));
if (mmax > cheb_table_mmax)
throw std::invalid_argument(
"FmEval_Chebyshev7::init() : requested mmax exceeds the "
"hard-coded mmax");
if (m_max >= 0)
init_table();
}
~FmEval_Chebyshev7() {
if (mmax >= 0) {
free(c);
}
}
/// Singleton interface allows to manage the lone instance; adjusts max m values as needed in thread-safe fashion
static std::shared_ptr<const FmEval_Chebyshev7> instance(int m_max, double = 0.0) {
assert(m_max >= 0);
// thread-safe per C++11 standard [6.7.4]
static auto instance_ = std::make_shared<const FmEval_Chebyshev7>(m_max);
while (instance_->max_m() < m_max) {
static std::mutex mtx;
std::lock_guard<std::mutex> lck(mtx);
if (instance_->max_m() < m_max) {
auto new_instance = std::make_shared<const FmEval_Chebyshev7>(m_max);
instance_ = new_instance;
}
}
return instance_;
}
/// @return the maximum value of m for which the Boys function can be computed with this object
int max_m() const { return mmax; }
/// fills in Fm with computed Boys function values for m in [0,mmax]
/// @param[out] Fm array to be filled in with the Boys function values, must be at least mmax+1 elements long
/// @param[in] x the Boys function argument
/// @param[in] mmax the maximum value of m for which Boys function will be computed; mmax must be <= the value returned by max_m
inline void eval(Real* Fm, Real x, int m_max) const {
// large T => use upward recursion
// cost = 1 div + 1 sqrt + (1 + 2*(m-1)) muls
if (x > T_crit) {
const double one_over_x = 1/x;
Fm[0] = 0.88622692545275801365 * sqrt(one_over_x); // see Eq. (9.8.9) in Helgaker-Jorgensen-Olsen
if (m_max == 0)
return;
// this upward recursion formula omits - e^(-x)/(2x), which for x>T_crit is small enough to guarantee full double precision
for (int i = 1; i <= m_max; i++)
Fm[i] = Fm[i - 1] * numbers_.ihalf[i] * one_over_x; // see Eq. (9.8.13)
return;
}
// ---------------------------------------------
// small and intermediate arguments => interpolate Fm and (optional) downward recursion
// ---------------------------------------------
// which interval does this x fall into?
const Real x_over_delta = x * one_over_delta;
const int iv = int(x_over_delta); // the interval index
const Real xd = x_over_delta - (Real)iv - 0.5; // this ranges from -0.5 to 0.5
const int m_min = 0;
#if defined(__AVX__)
const auto x2 = xd*xd;
const auto x3 = x2*xd;
const auto x4 = x2*x2;
const auto x5 = x2*x3;
const auto x6 = x3*x3;
const auto x7 = x3*x4;
libint2::simd::VectorAVXDouble x0vec(1., xd, x2, x3);
libint2::simd::VectorAVXDouble x1vec(x4, x5, x6, x7);
#endif // AVX
const Real *d = c + (ORDERp1) * (iv * (mmax+1) + m_min); // ptr to the interpolation data for m=mmin
int m = m_min;
#if defined(__AVX__)
if (m_max-m >=3) {
const int unroll_size = 4;
const int m_fence = (m_max + 2 - unroll_size);
for(; m<m_fence; m+=unroll_size, d+=ORDERp1*unroll_size) {
libint2::simd::VectorAVXDouble d00v, d01v, d10v, d11v,
d20v, d21v, d30v, d31v;
d00v.load_aligned(d); d01v.load_aligned((d+4));
d10v.load_aligned(d+ORDERp1); d11v.load_aligned((d+4)+ORDERp1);
d20v.load_aligned(d+2*ORDERp1); d21v.load_aligned((d+4)+2*ORDERp1);
d30v.load_aligned(d+3*ORDERp1); d31v.load_aligned((d+4)+3*ORDERp1);
libint2::simd::VectorAVXDouble fm0 = d00v * x0vec + d01v * x1vec;
libint2::simd::VectorAVXDouble fm1 = d10v * x0vec + d11v * x1vec;
libint2::simd::VectorAVXDouble fm2 = d20v * x0vec + d21v * x1vec;
libint2::simd::VectorAVXDouble fm3 = d30v * x0vec + d31v * x1vec;
libint2::simd::VectorAVXDouble sum0123 = horizontal_add(fm0, fm1, fm2, fm3);
sum0123.convert(&Fm[m]);
}
} // unroll_size=4
if (m_max-m >=1) {
const int unroll_size = 2;
const int m_fence = (m_max + 2 - unroll_size);
for(; m<m_fence; m+=unroll_size, d+=ORDERp1*unroll_size) {
libint2::simd::VectorAVXDouble d00v, d01v, d10v, d11v;
d00v.load_aligned(d);
d01v.load_aligned((d+4));
d10v.load_aligned(d+ORDERp1);
d11v.load_aligned((d+4)+ORDERp1);
libint2::simd::VectorAVXDouble fm0 = d00v * x0vec + d01v * x1vec;
libint2::simd::VectorAVXDouble fm1 = d10v * x0vec + d11v * x1vec;
libint2::simd::VectorSSEDouble sum01 = horizontal_add(fm0, fm1);
sum01.convert(&Fm[m]);
}
} // unroll_size=2
{ // no unrolling
for(; m<=m_max; ++m, d+=ORDERp1) {
libint2::simd::VectorAVXDouble d0v, d1v;
d0v.load_aligned(d);
d1v.load_aligned(d+4);
Fm[m] = horizontal_add(d0v * x0vec + d1v * x1vec);
}
}
#else // AVX not available
for(m=m_min; m<=m_max; ++m, d+=ORDERp1) {
Fm[m] = d[0]
+ xd * (d[1]
+ xd * (d[2]
+ xd * (d[3]
+ xd * (d[4]
+ xd * (d[5]
+ xd * (d[6]
+ xd * (d[7])))))));
// // check against the reference value
// if (false) {
// double refvalue = FmEval_Reference2<double>::eval(x, m); // compute F(T)
// if (abs(refvalue - Fm[m]) > 1e-10) {
// std::cout << "T = " << x << " m = " << m << " cheb = "
// << Fm[m] << " ref = " << refvalue << std::endl;
// }
// }
}
#endif
} // eval()
private:
void init_table() {
// get memory
void* result;
int status = posix_memalign(&result, ORDERp1*sizeof(Real), (mmax + 1) * cheb_table_nintervals * ORDERp1 * sizeof(Real));
if (status != 0) {
if (status == EINVAL)
throw std::logic_error(
"FmEval_Chebyshev7::init() : posix_memalign failed, alignment must be a power of 2 at least as large as sizeof(void *)");
if (status == ENOMEM)
throw std::bad_alloc();
abort(); // should be unreachable
}
c = static_cast<Real*>(result);
// copy contents of static table into c
// need all intervals
for (std::size_t iv = 0; iv < cheb_table_nintervals; ++iv) {
// but only values of m up to mmax
std::copy(cheb_table[iv], cheb_table[iv]+(mmax+1)*ORDERp1, c+(iv*(mmax+1))*ORDERp1);
}
}
}; // FmEval_Chebyshev7
#if LIBINT2_CONSTEXPR_STATICS
template <typename Real>
constexpr
double FmEval_Chebyshev7<Real>::cheb_table[FmEval_Chebyshev7<Real>::cheb_table_nintervals][(FmEval_Chebyshev7<Real>::cheb_table_mmax+1)*(FmEval_Chebyshev7<Real>::interpolation_order+1)];
#else
// clang needs an explicit specifalization declaration to avoid warning
# ifdef __clang__
template <typename Real>
double FmEval_Chebyshev7<Real>::cheb_table[FmEval_Chebyshev7<Real>::cheb_table_nintervals][(FmEval_Chebyshev7<Real>::cheb_table_mmax+1)*(FmEval_Chebyshev7<Real>::interpolation_order+1)];
# endif
#endif
#ifndef STATIC_OON
#define STATIC_OON
namespace {
const double oon[] = {0.0, 1.0, 1.0/2.0, 1.0/3.0, 1.0/4.0, 1.0/5.0, 1.0/6.0, 1.0/7.0, 1.0/8.0, 1.0/9.0, 1.0/10.0, 1.0/11.0};
}
#endif
/** Computes the Boys function, \$ F_m (T) = \int_0^1 u^{2m} \exp(-T u^2) \, {\rm d}u \$,
* using Taylor interpolation of up to 8-th order.
* @tparam Real the type to use for all floating-point computations. Following expressions must be valid: exp(Real), pow(Real), fabs(Real), max(Real), and floor(Real).
* @tparam INTERPOLATION_ORDER the interpolation order. The higher the order the less memory this object will need, but the computational cost will increase (usually very slightly)
*/
template<typename Real = double, int INTERPOLATION_ORDER = 7>
class FmEval_Taylor {
public:
static_assert(std::is_same<Real,double>::value, "FmEval_Taylor only supports double as the real type");
static const int max_interp_order = 8;
static const bool INTERPOLATION_AND_RECURSION = false; // compute F_lmax(T) and then iterate down to F_0(T)? Else use interpolation only
const double soft_zero_;
/// Constructs the object to be able to compute Boys funcion for m in [0,mmax], with relative \c precision
FmEval_Taylor(unsigned int mmax, Real precision = std::numeric_limits<Real>::epsilon()) :
soft_zero_(1e-6), cutoff_(precision), numbers_(
INTERPOLATION_ORDER + 1, 2 * (mmax + INTERPOLATION_ORDER - 1)) {
#if defined(__x86_64__) || defined(_M_X64) || defined(__i386) || defined(_M_IX86)
# if !defined(__AVX__) && defined(NDEBUG)
if (libint2::verbose()) {
static bool printed_performance_warning = false;
if (!printed_performance_warning) {
libint2::verbose_stream()
<< "libint2::FmEval_Taylor on x86(-64) platforms needs AVX support for best performance"
<< std::endl;
printed_performance_warning = true;
}
}
# endif
#endif
assert(mmax <= 63);
const double sqrt_pi = std::sqrt(M_PI);
/*---------------------------------------
We are doing Taylor interpolation with
n=TAYLOR_ORDER terms here:
error <= delT^n/(n+1)!
---------------------------------------*/
delT_ = 2.0
* std::pow(cutoff_ * numbers_.fac[INTERPOLATION_ORDER + 1],
1.0 / INTERPOLATION_ORDER);
oodelT_ = 1.0 / delT_;
max_m_ = mmax + INTERPOLATION_ORDER - 1;
T_crit_ = new Real[max_m_ + 1]; /*--- m=0 is included! ---*/
max_T_ = 0;
/*--- Figure out T_crit for each m and put into the T_crit ---*/
for (int m = max_m_; m >= 0; --m) {
/*------------------------------------------
Damped Newton-Raphson method to solve
T^{m-0.5}*exp(-T) = epsilon*Gamma(m+0.5)
The solution is the max T for which to do
the interpolation
------------------------------------------*/
double T = -log(cutoff_);
const double egamma = cutoff_ * sqrt_pi * numbers_.df[2 * m]
/ std::pow(2.0, m);
double T_new = T;
double func;
do {
const double damping_factor = 0.2;
T = T_new;
/* f(T) = the difference between LHS and RHS of the equation above */
func = std::pow(T, m - 0.5) * std::exp(-T) - egamma;
const double dfuncdT = ((m - 0.5) * std::pow(T, m - 1.5)
- std::pow(T, m - 0.5)) * std::exp(-T);
/* f(T) has 2 roots and has a maximum in between. If f'(T) > 0 we are to the left of the hump. Make a big step to the right. */
if (dfuncdT > 0.0) {
T_new *= 2.0;
} else {
/* damp the step */
double deltaT = -func / dfuncdT;
const double sign_deltaT = (deltaT > 0.0) ? 1.0 : -1.0;
const double max_deltaT = damping_factor * T;
if (std::fabs(deltaT) > max_deltaT)
deltaT = sign_deltaT * max_deltaT;
T_new = T + deltaT;
}
if (T_new <= 0.0) {
T_new = T / 2.0;
}
} while (std::fabs(func / egamma) >= soft_zero_);
T_crit_[m] = T_new;
const int T_idx = (int) std::floor(T_new / delT_);
max_T_ = std::max(max_T_, T_idx);
}
// allocate the grid (see the comments below)
{
const int nrow = max_T_ + 1;
const int ncol = max_m_ + 1;
grid_ = new Real*[nrow];
grid_[0] = new Real[nrow * ncol];
//std::cout << "Allocated interpolation table of " << nrow * ncol << " reals" << std::endl;
for (int r = 1; r < nrow; ++r)
grid_[r] = grid_[r - 1] + ncol;
}
/*-------------------------------------------------------
Tabulate the gamma function from t=delT to T_crit[m]:
1) include T=0 though the table is empty for T=0 since
Fm(0) is simple to compute
-------------------------------------------------------*/
/*--- do the mmax first ---*/
for (int T_idx = max_T_; T_idx >= 0; --T_idx) {
const double T = T_idx * delT_;
libint2::FmEval_Reference2<double>::eval(grid_[T_idx], T, max_m_);
}
}
~FmEval_Taylor() {
delete[] T_crit_;
delete[] grid_[0];
delete[] grid_;
}
/// Singleton interface allows to manage the lone instance;
/// adjusts max m and precision values as needed in thread-safe fashion
static std::shared_ptr<const FmEval_Taylor> instance(unsigned int mmax, Real precision = std::numeric_limits<Real>::epsilon()) {
assert(mmax >= 0);
assert(precision >= 0);
// thread-safe per C++11 standard [6.7.4]
static auto instance_ = std::make_shared<const FmEval_Taylor>(mmax, precision);
while (instance_->max_m() < mmax || instance_->precision() > precision) {
static std::mutex mtx;
std::lock_guard<std::mutex> lck(mtx);
if (instance_->max_m() < mmax || instance_->precision() > precision) {
auto new_instance = std::make_shared<const FmEval_Taylor>(mmax, precision);
instance_ = new_instance;
}
}
return instance_;
}
/// @return the maximum value of m for which this object can compute the Boys function
int max_m() const { return max_m_ - INTERPOLATION_ORDER + 1; }
/// @return the precision with which this object can compute the Boys function
Real precision() const { return cutoff_; }
/// computes Boys function values with m index in range [0,mmax]
/// @param[out] Fm array to be filled in with the Boys function values, must be at least mmax+1 elements long
/// @param[in] x the Boys function argument
/// @param[in] mmax the maximum value of m for which Boys function will be computed;
/// it must be <= the value returned by max_m() (this is not checked)
void eval(Real* Fm, Real T, int mmax) const {
const double sqrt_pio2 = 1.2533141373155002512;
const double two_T = 2.0 * T;
// stop recursion at mmin
const int mmin = INTERPOLATION_AND_RECURSION ? mmax : 0;
/*-------------------------------------
Compute Fm(T) from mmax down to mmin
-------------------------------------*/
const bool use_upward_recursion = true;
if (use_upward_recursion) {
// if (T > 30.0) {
if (T > T_crit_[0]) {
const double one_over_x = 1.0/T;
Fm[0] = 0.88622692545275801365 * sqrt(one_over_x); // see Eq. (9.8.9) in Helgaker-Jorgensen-Olsen
if (mmax == 0)
return;
// this upward recursion formula omits - e^(-x)/(2x), which for x>T_crit is <1e-15
for (int i = 1; i <= mmax; i++)
Fm[i] = Fm[i - 1] * numbers_.ihalf[i] * one_over_x; // see Eq. (9.8.13)
return;
}
}
// since Tcrit grows with mmax, this condition only needs to be determined once
if (T > T_crit_[mmax]) {
double pow_two_T_to_minusjp05 = std::pow(two_T, -mmax - 0.5);
for (int m = mmax; m >= mmin; --m) {
/*--- Asymptotic formula ---*/
Fm[m] = numbers_.df[2 * m] * sqrt_pio2 * pow_two_T_to_minusjp05;
pow_two_T_to_minusjp05 *= two_T;
}
}
else
{
const int T_ind = (int) (0.5 + T * oodelT_);
const Real h = T_ind * delT_ - T;
const Real* F_row = grid_[T_ind] + mmin;
#if defined (__AVX__)
libint2::simd::VectorAVXDouble h0123, h4567;
if (INTERPOLATION_ORDER == 3 || INTERPOLATION_ORDER == 7) {
const double h2 = h*h*oon[2];
const double h3 = h2*h*oon[3];
h0123 = libint2::simd::VectorAVXDouble (1.0, h, h2, h3);
if (INTERPOLATION_ORDER == 7) {
const double h4 = h3*h*oon[4];
const double h5 = h4*h*oon[5];
const double h6 = h5*h*oon[6];
const double h7 = h6*h*oon[7];
h4567 = libint2::simd::VectorAVXDouble (h4, h5, h6, h7);
}
}
// libint2::simd::VectorAVXDouble h0123(1.0);
// libint2::simd::VectorAVXDouble h4567(1.0);
#endif
int m = mmin;
if (mmax-m >=1) {
const int unroll_size = 2;
const int m_fence = (mmax + 2 - unroll_size);
for(; m<m_fence; m+=unroll_size, F_row+=unroll_size) {
#if defined(__AVX__)
if (INTERPOLATION_ORDER == 3 || INTERPOLATION_ORDER == 7) {
libint2::simd::VectorAVXDouble fr0_0123; fr0_0123.load(F_row);
libint2::simd::VectorAVXDouble fr1_0123; fr1_0123.load(F_row+1);
libint2::simd::VectorSSEDouble fm01 = horizontal_add(fr0_0123*h0123, fr1_0123*h0123);
if (INTERPOLATION_ORDER == 7) {
libint2::simd::VectorAVXDouble fr0_4567; fr0_4567.load(F_row+4);
libint2::simd::VectorAVXDouble fr1_4567; fr1_4567.load(F_row+5);
fm01 += horizontal_add(fr0_4567*h4567, fr1_4567*h4567);
}
fm01.convert(&Fm[m]);
}
else {
#endif
Real total0 = 0.0, total1 = 0.0;
for(int i=INTERPOLATION_ORDER; i>=1; --i) {
total0 = oon[i]*h*(F_row[i] + total0);
total1 = oon[i]*h*(F_row[i+1] + total1);
}
Fm[m] = F_row[0] + total0;
Fm[m+1] = F_row[1] + total1;
#if defined(__AVX__)
}
#endif
}
} // unroll_size = 2
if (m<=mmax) { // unroll_size = 1
#if defined(__AVX__)
if (INTERPOLATION_ORDER == 3 || INTERPOLATION_ORDER == 7) {
libint2::simd::VectorAVXDouble fr0123; fr0123.load(F_row);
if (INTERPOLATION_ORDER == 7) {
libint2::simd::VectorAVXDouble fr4567; fr4567.load(F_row+4);
// libint2::simd::VectorSSEDouble fm = horizontal_add(fr0123*h0123, fr4567*h4567);
// Fm[m] = horizontal_add(fm);
Fm[m] = horizontal_add(fr0123*h0123 + fr4567*h4567);
}
else { // INTERPOLATION_ORDER == 3
Fm[m] = horizontal_add(fr0123*h0123);
}
}
else {
#endif
Real total = 0.0;
for(int i=INTERPOLATION_ORDER; i>=1; --i) {
total = oon[i]*h*(F_row[i] + total);
}
Fm[m] = F_row[0] + total;
#if defined(__AVX__)
}
#endif
} // unroll_size = 1
// check against the reference value
// if (false) {
// double refvalue = FmEval_Reference2<double>::eval(T, mmax); // compute F(T) with m=mmax
// if (abs(refvalue - Fm[mmax]) > 1e-14) {
// std::cout << "T = " << T << " m = " << mmax << " cheb = "
// << Fm[mmax] << " ref = " << refvalue << std::endl;
// }
// }
} // if T < T_crit
/*------------------------------------
And then do downward recursion in j
------------------------------------*/
if (INTERPOLATION_AND_RECURSION && mmin > 0) {
const Real exp_mT = std::exp(-T);
for (int m = mmin - 1; m >= 0; --m) {
Fm[m] = (exp_mT + two_T * Fm[m+1]) * numbers_.twoi1[m];
}
}
}
private:
Real **grid_; /* Table of "exact" Fm(T) values. Row index corresponds to
values of T (max_T+1 rows), column index to values
of m (max_m+1 columns) */
Real delT_; /* The step size for T, depends on cutoff */
Real oodelT_; /* 1.0 / delT_, see above */
Real cutoff_; /* Tolerance cutoff used in all computations of Fm(T) */
int max_m_; /* Maximum value of m in the table, depends on cutoff
and the number of terms in Taylor interpolation */
int max_T_; /* Maximum index of T in the table, depends on cutoff
and m */
Real *T_crit_; /* Maximum T for each row, depends on cutoff;
for a given m and T_idx <= max_T_idx[m] use Taylor interpolation,
for a given m and T_idx > max_T_idx[m] use the asymptotic formula */
ExpensiveNumbers<double> numbers_;
/**
* Power series estimate of the error introduced by replacing
* \f$ F_m(T) = \int_0^1 \exp(-T t^2) t^{2 m} \, \mathrm{d} t \f$ with analytically
* integrable \f$ \int_0^\infty \exp(-T t^2) t^{2 m} \, \mathrm{d} t = \frac{(2m-1)!!}{2^{m+1}} \sqrt{\frac{\pi}{T^{2m+1}}} \f$
* @param m
* @param T
* @return the error estimate
*/
static double truncation_error(unsigned int m, double T) {
const double m2= m * m;
const double m3= m2 * m;
const double m4= m2 * m2;
const double T2= T * T;
const double T3= T2 * T;
const double T4= T2 * T2;
const double T5= T2 * T3;
const double result = exp(-T) * (105 + 16*m4 + 16*m3*(T - 8) - 30*T + 12*T2
- 8*T3 + 16*T4 + 8*m2*(43 - 9*T + 2*T2) +
4*m*(-88 + 23*T - 8*T2 + 4*T3))/(32*T5);
return result;
}
/**
* Leading-order estimate of the error introduced by replacing
* \f$ F_m(T) = \int_0^1 \exp(-T t^2) t^{2 m} \, \mathrm{d} t \f$ with analytically
* integrable \f$ \int_0^\infty \exp(-T t^2) t^{2 m} \, \mathrm{d} t = \frac{(2m-1)!!}{2^{m+1}} \sqrt{\frac{\pi}{T^{2m+1}}} \f$
* @param m
* @param T
* @return the error estimate
*/
static double truncation_error(double T) {
const double result = exp(-T) /(2*T);
return result;
}
};
namespace detail {
constexpr double pow10(long exp) {
return (exp == 0) ? 1. : ((exp > 0) ? 10. * pow10(exp-1) : 0.1 * pow10(exp+1));
}
}
//////////////////////////////////////////////////////////
/// core integral for Yukawa and exponential interactions
//////////////////////////////////////////////////////////
/**
* Evaluates core integral for Gaussian integrals over the Yukawa potential \f$ \exp(- \zeta r) / r \f$ and
* the exponential interaction \f$ \exp(- \zeta r) \f$
* @tparam Real real type
* @warning only @p Real = double is supported
* @warning guarantees absolute precision of only about 1e-14
*/
template<typename Real = double>
struct TennoGmEval {
private:
#include <libint2/tenno_cheb.h>
static_assert(std::is_same<Real,double>::value, "TennoGmEval only supports double as the real type");
static const int mmin_ = -1;
static constexpr const Real Tmax = (1 << cheb_table_tmaxlog2); //!< critical value of T above which use upward recursion
static constexpr const Real Umax = detail::pow10(cheb_table_umaxlog10); //!< max value of U for which to interpolate
static constexpr const Real Umin = detail::pow10(cheb_table_uminlog10); //!< min value of U for which to interpolate
static constexpr const std::size_t ORDERp1 = interpolation_order + 1;
static constexpr const Real maxabsprecision = 1.4e-14; //!< guaranteed abs precision of the interpolation table for m>0
public:
/// \param m_max maximum value of the Gm function index
/// \param precision the desired *absolute* precision (relative precision for most intervals will be below epsilon, but for large T/U values and high m relative precision is low
/// \throw std::invalid_argument if \c m_max is greater than \c cheb_table_mmax (see tenno_cheb.h)
/// \throw std::invalid_argument if \c precision is smaller than \c maxabsprecision
TennoGmEval(unsigned int mmax, Real precision = -1) :
mmax_(mmax), precision_(precision),
numbers_() {
#if defined(__x86_64__) || defined(_M_X64) || defined(__i386) || defined(_M_IX86)
# if !defined(__AVX__) && defined(NDEBUG)
if (libint2::verbose()) {
static bool printed_performance_warning = false;
if (!printed_performance_warning) {
libint2::verbose_stream()
<< "libint2::TennoGmEval on x86(-64) platforms needs AVX support for best performance"
<< std::endl;
printed_performance_warning = true;
}
}
# endif
#endif
// if (precision_ < maxabsprecision)
// throw std::invalid_argument(
// std::string(
// "TennoGmEval does not support precision smaller than ") +
// std::to_string(maxabsprecision));
if (mmax > cheb_table_mmax)
throw std::invalid_argument(
"TennoGmEval::init() : requested mmax exceeds the "
"hard-coded mmax");
init_table();
}
~TennoGmEval() {
if (c_ != nullptr)
free(c_);
}
/// Singleton interface allows to manage the lone instance; adjusts max m values as needed in thread-safe fashion
static std::shared_ptr<const TennoGmEval> instance(int m_max, double = 0) {
assert(m_max >= 0);
// thread-safe per C++11 standard [6.7.4]
static auto instance_ = std::make_shared<const TennoGmEval>(m_max);
while (instance_->max_m() < m_max) {
static std::mutex mtx;
std::lock_guard<std::mutex> lck(mtx);
if (instance_->max_m() < m_max) {
auto new_instance = std::make_shared<const TennoGmEval>(m_max);
instance_ = new_instance;
}
}
return instance_;
}
unsigned int max_m() const { return mmax_; }
/// @return the precision with which this object can compute the result
Real precision() const { return precision_; }
/// @param[in] Gm pointer to array of @c mmax+1 @c Real elements, on
/// return this contains the core integral for Yukawa
/// interaction, \f$ \exp(-zeta r_{12})/r_{12} \f$ ,
/// namely \f$ G_{m}(T,U) = \int_0^1 t^{2m} \exp(U(1-t^{-2}) - Tt^2) dt, m \in [0,m_\mathrm{max}] \f$, where
/// \f$ T = \rho |\vec{P} - \vec{Q}|^2 \f$ and
/// \f$ U = \zeta^2 / (4 \rho) \f$
/// @param[in] one_over_rho \f$ 1/\rho \f$
/// @param[in] T \f$ T \f$
/// @param[in] mmax \f$ m_\mathrm{max} \f$
/// @param[in] zeta \f$ \zeta \f$
void eval_yukawa(Real* Gm, Real one_over_rho, Real T, size_t mmax, Real zeta) const {
assert(mmax <= mmax_);
assert(T >= 0);
const auto U = 0.25 * zeta * zeta * one_over_rho;
assert(U >= 0);
if (T > Tmax || U < Umin) {
eval_Gm_urr(Gm, T, U, mmax, 0); // no need for G_-1
} else {
interpolate_Gm<false>(Gm, T, U, 0, mmax);
}
}
/// @param[in] Gm pointer to array of @c mmax+1 @c Real elements, on
/// return this contains the core integral for Slater-type
/// geminal, \f$ \exp(-zeta r_{12}) \f$ ,
/// namely \f$ \sqrt{U} (G_{m-1}(T,U) - G_m(T,U)), m \in [0,m_\mathrm{max}] \f$ where\
// \f$ G_{m}(T,U) = \int_0^1 t^{2m} \exp(U(1-t^{-2}) - Tt^2) dt \f$, where
/// \f$ T = \rho |\vec{P} - \vec{Q}|^2 \f$ and
/// \f$ U = \zeta^2 / (4 \rho) \f$
/// @param[in] one_over_rho \f$ 1/\rho \f$
/// @param[in] T \f$ T \f$
/// @param[in] mmax \f$ m_\mathrm{max} \f$
/// @param[in] zeta \f$ \zeta \f$
void eval_slater(Real* Gm, Real one_over_rho, Real T, size_t mmax, Real zeta) const {
assert(mmax <= mmax_);
assert(T >= 0);
const auto U = 0.25 * zeta * zeta * one_over_rho;
assert(U > 0); // integral factorizes into 2 overlaps for U = 0
const auto zeta_over_two_rho = 0.5 * zeta * one_over_rho;
if (T > Tmax) {
eval_Gm_urr(Gm, T, U, mmax, -1);
} else {
interpolate_Gm<true>(Gm, T, U, zeta_over_two_rho, mmax);
}
}
private:
/// Computes G_0 and (optionally) G_{-1}
/// @tparam compute_Gm1, if true, compute G_0 and G_{-1}, otherwise G_0 only
/// @param[in] T the value of \f$ T \f$
/// @param[in] U the value of \f$ U \f$
/// @return if @c compute_Gm1==true return {G_0,G_{-1}}, otherwise {G_0,0}
template <bool compute_Gm1>
static inline std::tuple<Real,Real> eval_G0_and_maybe_Gm1(Real T, Real U) {
assert(T >= 0);
assert(U >= 0);
const Real sqrtPi(1.7724538509055160272981674833411451827975494561224);
Real G_m1 = 0;
Real G_0 = 0;
if (U == 0) { // \sqrt{U} G_{-1} is finite, need to handle that case separately
assert(compute_Gm1 == false);
if (T < std::numeric_limits<Real>::epsilon()) {
G_0 = 1;
}
else {
const Real sqrt_T = sqrt(T);
const Real sqrtPi_over_2 = sqrtPi * 0.5;
G_0 = sqrtPi_over_2 * erf(sqrt_T) / sqrt_T;
}
}
else if (T > std::numeric_limits<Real>::epsilon()) { // U > 0
const Real sqrt_U = sqrt(U);
const Real sqrt_T = sqrt(T);
const Real oo_sqrt_T = 1 / sqrt_T;
const Real oo_sqrt_U = compute_Gm1 ? (1 / sqrt_U) : 0;
const Real kappa = sqrt_U - sqrt_T;
const Real lambda = sqrt_U + sqrt_T;
const Real sqrtPi_over_4 = sqrtPi * 0.25;
const Real pfac = sqrtPi_over_4;
const Real erfc_k = exp(kappa * kappa - T) * erfc(kappa);
const Real erfc_l = exp(lambda * lambda - T) * erfc(lambda);
G_m1 = compute_Gm1 ? pfac * (erfc_k + erfc_l) * oo_sqrt_U : 0.0;
G_0 = pfac * (erfc_k - erfc_l) * oo_sqrt_T;
}
else { // T = 0, U > 0
const Real exp_U = exp(U);
const Real sqrt_U = sqrt(U);