-
Notifications
You must be signed in to change notification settings - Fork 0
/
main.cpp
2419 lines (2047 loc) · 94.5 KB
/
main.cpp
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
#include <iostream>
// # define MADNESS_HAS_LIBXC 0
// #define USE_GENTENSOR 0 // only needed if madness was configured with `-D ENABLE_GENTENSOR=1
// #include<madness/chem/nemo.h>
// #include<madness.h>
// #include<madness/chem.h>
// #include<madness/misc/info.h>
// #include<madness/mra/nonlinsol.h>
// #include<madness/world/timing_utilities.h>
#include<madchem.h>
using namespace madness;
enum Uplo {upper, lower};
static bool transform_c=false;
static bool debug=false;
static double alpha1=constants::fine_structure_constant;
static double shift=0.0;
static double epsilon=1.e-12;
static bool use_ble=false;
static const double bohr_rad=52917.7211;
static const double lo=1.e-7;
double compute_gamma(const double nuclear_charge) {
return sqrt(1-nuclear_charge*nuclear_charge*alpha1*alpha1);
}
struct stepfunction : public FunctionFunctorInterface<double,3> {
int axis=-1;
stepfunction(const int axis) : axis(axis) {
MADNESS_CHECK(axis>=0 && axis<3);
}
double operator()(const coord_3d& r) const override { return r[axis]/(r.normf()+epsilon); }
Level special_level() override {return 20;};
std::vector<Vector<double, 3UL>> special_points() const override {
coord_3d o={0.0,0.0,0.0};
return {o};
}
};
struct ncf_cusp {
double a=1.3;
double Z=-1;
ncf_cusp(double a, double Z) :a(a),Z(Z) {}
double operator()(const double& r) const {
if (a<0.0) return 1.0;
return 1.0+(1.0/(a-1.0))*exp(-a*Z*r);
}
};
struct Sigma_ncf_cusp {
double a=1.3;
double Z=-1;
Sigma_ncf_cusp(double a, double Z) :a(a),Z(Z) {}
double operator()(const double& r) const {
if (a<0.0) return 0.0;
const double denominator=1.0 + 1.0 / (a - 1.0) * exp(-a * Z * r);
const double numerator=-a*Z/(a-1.0) * exp(-a * Z * r);
return numerator/denominator;
}
};
struct ncf_singularity {
double gamma=0.0;
ncf_singularity(double gamma) : gamma(gamma) {}
double operator()(const double& r) const {
if (gamma<0.0) return 1.0;
return 1.0/(std::pow(r,1.0-gamma)+epsilon);
}
};
struct ncf : public FunctionFunctorInterface<double,3> {
ncf_singularity ncfs;
ncf_cusp ncfc;
long power=1;
ncf(double gamma, double a, double Z) : ncfs(gamma), ncfc(a,Z) {}
double operator()(const coord_3d& r) const override {
double rr=r.normf();
return std::pow(ncfs(rr) * ncfc(rr),power);
}
Level special_level() override {return 20;};
std::vector<Vector<double, 3UL>> special_points() const override {
coord_3d o={0.0,0.0,0.0};
return {o};
}
};
double generalized_laguerre(const double alpha, const long n, const double r) {
if (n<0) return 0.0;
else if (n==0) return 1.0;
else if (n==1) return (1.0+alpha-r);
else if (n==2) return 0.5*r*r - (alpha + 2.0)*r + 0.5*(alpha+1)*(alpha+2);
else
MADNESS_EXCEPTION("generalized Laguerre polynomial implemented only up to order n=2",1);
}
//Functor to make the fermi nuclear charge distribution (NOT the potential, and NOT normalized) for a given center
//This is based on Visscher and Dyall's 1997 paper on nuclear charge distributions.
class FermiNucDistFunctor : public FunctionFunctorInterface<double,3> {
private:
int m_A;
double m_T;
double m_C;
std::vector<coord_3d> m_R;
public:
// Constructor
FermiNucDistFunctor(int& Z, coord_3d R, double bohr_rad){
//m_T = 0.000043463700858425666; //2.3 fm in bohr
m_T = 2.3/bohr_rad;
m_R.push_back(R);
//find atomic mass numbers for each atom. This list matches that of Visccher and Dyall (1997)
int Alist[116] = {1,4,7,9,11,12,14,16,19,20,23,24,27,28,31,32,35,40,39,40,45,48,51,52,55,56,59,58,63,64,69,74,75,80,79,84,85,88,89,90,93,98,98,102,103,106,107,114,115,120,121,130,127,132,133,138,139,140,141,144,145,152,153,158,159,162,162,168,169,174,175,180,181,184,187,192,193,195,197,202,205,208,209,209,210,222,223,226,227,232,231,238,237,244,243,247,247,251,252,257,258,259,262,261,262,263,262,265,266,264,272,277,284,289,288,292};
m_A = Alist[Z-1];
double PI = constants::pi;
if(m_A < 5){
m_C = 0.000022291*pow(m_A, 1.0/3.0) - 0.0000090676;
}
else{
m_C = sqrt(5.0/3.0*pow((0.836*pow(m_A,1.0/3.0)+0.570)/bohr_rad,2) - 7.0/3.0*pow(PI*m_T/4.0/log(3.0),2));
}
}
//overload () operator
double operator() (const coord_3d&r) const {
double x = r[0] - m_R[0][0];
double y = r[1] - m_R[0][1];
double z = r[2] - m_R[0][2];
double rr = sqrt(x*x+y*y+z*z);
double result = 1.0/(1.0+exp(4.0*log(3.0)*(rr-m_C)/m_T));
return result;
}
//Because the distribution is only nonzero in a small window around the center, need to create a special point
std::vector<coord_3d> special_points() const {
return m_R;
}
madness::Level special_level() {
return 18;
}
//Print the parameters of the Fermi nuclear charge distribution
void print_details(World& world){
//Constants necessary to print the details. Technically need to use bohr_rad parameter here
int Alist[116] = {1,4,7,9,11,12,14,16,19,20,23,24,27,28,31,32,35,40,39,40,45,48,51,52,55,56,59,58,63,64,69,74,75,80,79,84,85,88,89,90,93,98,98,102,103,106,107,114,115,120,121,130,127,132,133,138,139,140,141,144,145,152,153,158,159,162,162,168,169,174,175,180,181,184,187,192,193,195,197,202,205,208,209,209,210,222,223,226,227,232,231,238,237,244,243,247,247,251,252,257,258,259,262,261,262,263,262,265,266,264,272,277,284,289,288,292};
double T = 2.3/52917.72490083583;
double PI = constants::pi;
if(world.rank()==0){
for(int i = 0; i < 116; i++){
double RMS = (0.836*pow(Alist[i],1.0/3.0)+0.570)/52917.72490083583;
double C;
if(Alist[i] < 5){
C = 0.000022291*pow(Alist[i], 1.0/3.0) - 0.0000090676;
}
else{
C = sqrt(5.0/3.0*pow(RMS,2)-7.0/3.0*pow(PI*T/4.0/log(3.0),2));
}
double xi = 3.0/2.0/pow(RMS,2);
printf("Z: %3i, A: %3i, RMS: %.10e, C: %.10e, xi: %.10e\n", i+1, Alist[i], RMS, C, xi);
}
}
}
};
//Creates the fermi nuclear potential from the charge distribution. Also calculates the nuclear repulsion energy
real_function_3d make_fermi_potential(World& world, double& nuclear_repulsion_energy, const double nuclear_charge){
real_convolution_3d op=CoulombOperator(world,lo,FunctionDefaults<3>::get_thresh());
if(world.rank()==0) print("\n***Making a Fermi Potential***");
//Get list of atom coordinates
// std::vector<coord_3d> Rlist = Init_params.molecule.get_all_coords_vec();
std::vector<coord_3d> Rlist = {{0.0,0.0,0.0}};
std::vector<int> Zlist(Rlist.size());
unsigned int num_atoms = Rlist.size();
//variables for upcoming loop
real_function_3d temp, potential;
double tempnorm;
//Go through the atoms in the molecule and construct the total charge distribution due to all nuclei
for(unsigned int i = 0; i < num_atoms; i++){
// Zlist[i] = Init_params.molecule.get_atomic_number(i);
Zlist[i] = nuclear_charge;
FermiNucDistFunctor rho(Zlist[i], Rlist[i],bohr_rad);
temp = real_factory_3d(world).functor(rho).truncate_mode(0);
tempnorm = temp.trace();
temp.scale(-Zlist[i]/tempnorm);
if(i == 0){
potential = temp;
//rho.print_details(world);
}
else{
potential += temp;
}
}
//Potential is found by application of the coulomb operator to the charge distribution
potential = apply(op,potential);
//Calculate the nuclear repulsion energy
//It doesn't change iteration to iteration, so we want to calculate it once and store the result
//We calculate it inside this function because here we already have access to the nuclear charges and coordinates
nuclear_repulsion_energy = 0.0;
double rr;
for(unsigned int m = 0; m < num_atoms; m++){
for(unsigned int n = m+1; n < num_atoms; n++){
coord_3d dist = Rlist[m] - Rlist[n];
rr = std::sqrt(dist[0]*dist[0]+dist[1]*dist[1]+dist[2]*dist[2]);
nuclear_repulsion_energy += Zlist[m]*Zlist[n]/rr;
}
}
return potential;
}
/// returns the complex value of a given spherical harmonic
struct SphericalHarmonics{
long l, m;
bool zero=false;
SphericalHarmonics(const long l, const long m) : l(l), m(m) {
if (l<0) this->l=-l-1;
if (abs(this->m)>this->l) zero=true;
}
double_complex operator()(const coord_3d& xyz) const {
if (zero) return {0.0,0.0};
const double r=xyz.normf();
const double r2=r*r;
const double x=xyz[0];
const double y=xyz[1];
const double z=xyz[2];
const double_complex i=double_complex(0.0,1.0);
auto stepx=stepfunction(0);
auto stepy=stepfunction(1);
auto stepz=stepfunction(2);
double_complex step_x_m_iy=stepx(xyz) - i* stepy(xyz);
double_complex step_x_p_iy=stepx(xyz) + i* stepy(xyz);
if ((l==0) and (m== 0)) return 0.5*sqrt(1.0/constants::pi);
// if ((l==1) and (m==-1)) return 0.5*sqrt(1.5/constants::pi) * (x - i*y)/r;
// if ((l==1) and (m== 0)) return 0.5*sqrt(3.0/constants::pi) * z/r;
// if ((l==1) and (m== 1)) return -0.5*sqrt(1.5/constants::pi) * (x + i*y)/r;
if ((l==1) and (m==-1)) return 0.5*sqrt(1.5/constants::pi) * (stepx(xyz) - i *stepy(xyz));
if ((l==1) and (m== 0)) return 0.5*sqrt(3.0/constants::pi) * stepz(xyz);
if ((l==1) and (m== 1)) return -0.5*sqrt(1.5/constants::pi) * (stepx(xyz) + i *stepy(xyz));
if ((l==2) and (m==-2)) return 0.25*sqrt(7.5/constants::pi) * step_x_m_iy * step_x_m_iy;
if ((l==2) and (m==-1)) return 0.5 *sqrt(7.5/constants::pi) * step_x_m_iy * stepz(xyz);
if ((l==2) and (m== 0)) return 0.25*sqrt(5.0/constants::pi) * (3.0*stepz(xyz)*stepz(xyz)- 1.0);
if ((l==2) and (m== 1)) return -0.5 *sqrt(7.5/constants::pi) * step_x_p_iy * stepz(xyz);
if ((l==2) and (m== 2)) return 0.25*sqrt(7.5/constants::pi) * step_x_p_iy * step_x_p_iy;
if ((l==3) and (m==-3)) return 0.125*sqrt(35/constants::pi) * step_x_m_iy * step_x_m_iy * step_x_m_iy;
if ((l==3) and (m==-2)) return 0.25*sqrt(105/constants::pi) * step_x_m_iy * step_x_m_iy * stepz(xyz);
if ((l==3) and (m==-1)) return 0.125*sqrt(21/constants::pi) * step_x_m_iy * (5*stepz(xyz)*stepz(xyz) - 1.0);
if ((l==3) and (m== 0)) return 0.25*sqrt(7.0/constants::pi) * (5*stepz(xyz)*stepz(xyz)*stepz(xyz) - 3.0*stepz(xyz));
if ((l==3) and (m== 1)) return -0.125*sqrt(21/constants::pi) * step_x_p_iy * (5*stepz(xyz)*stepz(xyz) - 1.0);
if ((l==3) and (m== 2)) return 0.25*sqrt(105/constants::pi) * step_x_p_iy * step_x_p_iy * stepz(xyz);
if ((l==3) and (m== 3)) return -0.125*sqrt(35/constants::pi) * step_x_p_iy * step_x_p_iy * step_x_p_iy;
MADNESS_EXCEPTION("out of range in SphericalHarmonics",1);
return double_complex(0.0,0.0);
}
};
struct sgl_guess {
long n,l,m;
double Z;
sgl_guess(const long n, const long l, const long m, const double Z) : n(n), l(l), m(m), Z(Z) {}
double_complex operator()(const coord_3d& coord) const {
double_complex Y=SphericalHarmonics(l,m)(coord);
double r=coord.normf();
double rho=2.0*Z*r/n;
double R= generalized_laguerre(2.*l+1.,n-l-1,rho);
double e=exp(-0.5*rho);
return e*R*std::pow(rho,double(l))*Y;
}
double energy() const {
return 0.5*Z*Z/(n*n);
}
complex_function_3d get_wf(World& world) const {
std::vector<coord_3d> special_points(1,coord_3d({0.0,0.0,0.0}));
complex_function_3d wf = complex_factory_3d(world).functor(*this).special_points(special_points);
double norm=wf.norm2();
wf.scale(1.0/norm);
return wf;
}
};
/// following Dyall: p 104f
struct Xi {
long k, component, a;
double m, l;
Xi(const long k, const double m, const double component, const double j, const long l)
: k(k), m(m), component(component), l(l) {
MADNESS_CHECK(component==0 or component==1);
a=compute_a(j,l);
}
static long compute_a(const double j, const long l) {
if (j>l) return 1;
if (j<l) return -1;
throw;
}
double_complex operator()(const coord_3d& xyz) const {
double_complex result;
double denominator2=2.0*l+1;
if (component==0) {
double numerator=a*sqrt((l+0.5+a*m)/denominator2);
auto Y=SphericalHarmonics(k,lround(m-0.5));
result=numerator*Y(xyz);
} else if (component==1) {
double numerator=sqrt((l+0.5-a*m)/denominator2);
auto Y=SphericalHarmonics(k,lround(m+0.5));
result=numerator*Y(xyz);
}
return result;
}
};
struct Omega {
long k, component;
double m;
Omega(const long k, const double m, const double component) : k(k), m(m), component(component) {
MADNESS_CHECK(component==0 or component==1);
}
double_complex operator()(const coord_3d& xyz) const {
double sgnk= (k>0) ? 1.0 : -1.0;
double_complex result;
if (component==0) {
double nn=sqrt((k+0.5-m)/(2*k+1));
auto Y=SphericalHarmonics(k,lround(m-0.5));
result=nn*Y(xyz);
} else if (component==1) {
double nn=sqrt((k+0.5+m)/(2*k+1));
auto Y=SphericalHarmonics(k,lround(m+0.5));
result=-sgnk*nn*Y(xyz);
}
return result;
}
};
void set_version(const int v, const double nuclear_charge) {
if (v==1) {
transform_c=false;
shift=0.0;
} else if (v==2) {
transform_c=false;
shift=0.0;
} else if (v==3) {
transform_c=true;
shift=-compute_gamma(nuclear_charge)/(alpha1*alpha1);
} else {
MADNESS_EXCEPTION("unknown version",1);
}
}
double compute_electronic_energy(const double energy) {
double c=1.0/alpha1;
return energy-c*c-shift;
}
template<typename T, std::size_t NDIM>
class MyDerivativeOperator : public SCFOperatorBase<T,NDIM> {
typedef Function<T,NDIM> functionT;
typedef std::vector<functionT> vecfuncT;
typedef Tensor<T> tensorT;
public:
MyDerivativeOperator(World& world, const int axis1) : world(world), axis(axis1) {}
void set_ble1() {ble=true;};
std::string info() const {return "D"+std::to_string(axis);}
functionT operator()(const functionT& ket) const {
vecfuncT vket(1,ket);
return this->operator()(vket)[0];
}
vecfuncT operator()(const vecfuncT& vket) const {
auto gradop = free_space_derivative<T,NDIM>(world, axis);
if (ble) gradop.set_ble1();
vecfuncT dvket=apply(world, gradop, vket, false);
world.gop.fence();
return dvket;
}
T operator()(const functionT& bra, const functionT& ket) const {
vecfuncT vbra(1,bra), vket(1,ket);
Tensor<T> tmat=this->operator()(vbra,vket);
return tmat(0l,0l);
}
tensorT operator()(const vecfuncT& vbra, const vecfuncT& vket) const {
const auto bra_equiv_ket = &vbra == &vket;
vecfuncT dvket=this->operator()(vket);
return matrix_inner(world,vbra,dvket, bra_equiv_ket);
}
private:
World& world;
int axis;
bool ble=false;
};
template<typename T, std::size_t NDIM>
class n_times_n_dot_p_minus_p : public SCFOperatorBase<T,NDIM> {
typedef Function<T,NDIM> functionT;
typedef std::vector<functionT> vecfuncT;
typedef Tensor<T> tensorT;
public:
n_times_n_dot_p_minus_p(World& world, const int axis1) : world(world), axis(axis1) {}
std::string info() const override {return "n"+std::to_string(axis)+" * n.p";}
functionT operator()(const functionT& ket) const override {
vecfuncT vket(1,ket);
return this->operator()(vket)[0];
}
vecfuncT operator()(const vecfuncT& vket) const override {
auto D0 = free_space_derivative<T,NDIM>(world,0);
auto D1 = free_space_derivative<T,NDIM>(world,1);
auto D2 = free_space_derivative<T,NDIM>(world,2);
auto Daxis = free_space_derivative<T,NDIM>(world,axis);
if (use_ble) {
D0.set_ble1();
D1.set_ble1();
D2.set_ble1();
}
world.gop.fence();
const long axis1=axis;
auto stepx=stepfunction(0);
auto stepy=stepfunction(1);
auto stepz=stepfunction(2);
auto step_axis=stepfunction(axis1);
real_function_3d n0=real_factory_3d(world).functor(stepx);
real_function_3d n1=real_factory_3d(world).functor(stepy);
real_function_3d n2=real_factory_3d(world).functor(stepz);
real_function_3d naxis=real_factory_3d(world).functor(step_axis);
auto result=naxis * (n0 * apply(world,D0,vket) + n1 * apply(world,D1,vket) + n2*apply(world,D2,vket)) - apply(world,Daxis,vket);
return truncate(result);
}
T operator()(const functionT& bra, const functionT& ket) const override {
vecfuncT vbra(1,bra), vket(1,ket);
Tensor<T> tmat=this->operator()(vbra,vket);
return tmat(0l,0l);
}
tensorT operator()(const vecfuncT& vbra, const vecfuncT& vket) const override {
const auto bra_equiv_ket = &vbra == &vket;
vecfuncT dvket=this->operator()(vket);
return matrix_inner(world,vbra,dvket, bra_equiv_ket);
}
private:
World& world;
int axis;
};
template<typename T, std::size_t NDIM>
class AngularMomentum_with_step : public SCFOperatorBase<T,NDIM> {
typedef Function<T,NDIM> functionT;
typedef std::vector<functionT> vecfuncT;
typedef Tensor<T> tensorT;
public:
AngularMomentum_with_step(World& world, const int axis1) : world(world), axis(axis1) {}
std::string info() const {return "1/r L"+std::to_string(axis);}
functionT operator()(const functionT& ket) const {
vecfuncT vket(1,ket);
return this->operator()(vket)[0];
}
vecfuncT operator()(const vecfuncT& vket) const {
// if axis==0 -> lx = y pz - z py etc
int index1=(axis+1)%3;
int index2=(axis+2)%3;
auto gradop1 = free_space_derivative<T,NDIM>(world,index1);
auto gradop2 = free_space_derivative<T,NDIM>(world,index2);
if (use_ble) {
gradop1.set_ble1();
gradop2.set_ble1();
}
const double_complex ii={0.0,1.0};
vecfuncT p1=-ii*apply(world, gradop1, vket);
vecfuncT p2=-ii*apply(world, gradop2, vket);
world.gop.fence();
auto step1=stepfunction(index1);
auto step2=stepfunction(index2);
real_function_3d r1=real_factory_3d(world).functor(step1);
real_function_3d r2=real_factory_3d(world).functor(step2);
auto result=r1*p2 - r2*p1;
return result;
}
T operator()(const functionT& bra, const functionT& ket) const {
vecfuncT vbra(1,bra), vket(1,ket);
Tensor<T> tmat=this->operator()(vbra,vket);
return tmat(0l,0l);
}
tensorT operator()(const vecfuncT& vbra, const vecfuncT& vket) const {
const auto bra_equiv_ket = &vbra == &vket;
vecfuncT dvket=this->operator()(vket);
return matrix_inner(world,vbra,dvket, bra_equiv_ket);
}
private:
World& world;
int axis;
};
/// defines a 4-spinor
class Spinor {
public:
vector_complex_function_3d components;
World& world() const {return components.front().world();}
Spinor() {
components.resize(4);
}
Spinor(World& world) {
components=zero_functions_compressed<double_complex,3>(world,4);
}
Spinor(const vector_complex_function_3d& components) : components(components){}
Spinor& operator+=(const Spinor& other) {
components+=other.components;
return *this;
}
Spinor& operator-=(const Spinor& other) {
components-=other.components;
return *this;
}
Spinor operator+(const Spinor& other) const {
Spinor result;
result.components=copy(world(),components); // deep copy
result.components+=other.components;
return result;
}
Spinor operator-(const Spinor& other) const {
Spinor result;
result.components=copy(world(),components); // deep copy
result.components-=other.components;
return result;
}
Spinor& truncate() {
madness::truncate(components);
return *this;
}
void print_norms(std::string text) {
auto norms = norm2s(world(), components);
auto realnorms = norm2s(world(), real(components));
auto imagnorms = norm2s(world(), imag(components));
// real_function_3d x=real_factory_3d(world()).functor([](const coord_3d& r) {return r[0];});
// real_function_3d y=real_factory_3d(world()).functor([](const coord_3d& r) {return r[1];});
// real_function_3d z=real_factory_3d(world()).functor([](const coord_3d& r) {return r[2];});
// auto x1=inner(this->components,x*(this->components));
// auto y1=inner(this->components,y*(this->components));
// auto z1=inner(this->components,z*(this->components));
// auto x2=inner(this->components,x*x*(this->components));
// auto y2=inner(this->components,y*y*(this->components));
// auto z2=inner(this->components,z*z*(this->components));
print(text,norms);
print(" -- real,imag",realnorms,imagnorms);
// print(" -- moments ",x1,y1,z1,x2,y2,z2);
// plot(text);
}
void plot(const std::string filename) const {
plot_plane(world(), real(components), "Re_"+filename);
plot_plane(world(), imag(components), "Im_"+filename);
}
friend double_complex inner(const Spinor& bra, const Spinor& ket) {
return inner(bra.components,ket.components);
}
};
template<typename T>
Spinor operator*(const T fac, const Spinor& arg) {
return Spinor(fac*arg.components);
}
template<typename T>
Spinor operator*(const Spinor& arg, const T fac) {
return Spinor(fac*arg.components);
}
Spinor copy(const Spinor& other) {
return Spinor(copy(other.world(),other.components));
}
std::vector<Spinor> copy(const std::vector<Spinor>& other) {
std::vector<Spinor> result;
for (auto& oo : other) result.push_back(copy(oo.world(),oo.components));
return result;
}
template<typename T>
std::vector<Spinor> operator*(const std::vector<Spinor>& arg, const T fac) {
std::vector<Spinor> result=copy(arg);
for (auto& r : result) r=r*fac;
return result;
}
double_complex inner(const std::vector<Spinor>& bra, const std::vector<Spinor> ket) {
double_complex result=0.0;
for (int i=0; i<bra.size(); ++i) result+=inner(bra[i],ket[i]);
return result;
}
Tensor<double_complex> matrix_inner(const std::vector<Spinor>& bra, const std::vector<Spinor> ket) {
Tensor<double_complex> result(bra.size(),ket.size());
for (int i=0; i<bra.size(); ++i) {
for (int j=0; j<ket.size(); ++j) {
result(i,j)=inner(bra[i],ket[j]);
}
}
return result;
}
std::vector<Spinor> operator-=(std::vector<Spinor>& left, const std::vector<Spinor>& right) {
for (int i=0; i<right.size(); ++i) left[i]-=right[i];
return left;
}
std::vector<Spinor> operator+=(std::vector<Spinor>& left, const std::vector<Spinor>& right) {
for (int i=0; i<right.size(); ++i) left[i]+=right[i];
return left;
}
std::vector<Spinor> operator-(const std::vector<Spinor>& left, const std::vector<Spinor>& right) {
std::vector<Spinor> result=copy(left);
result-=right;
return result;
}
std::vector<Spinor> truncate(std::vector<Spinor> arg) {
for (auto& a : arg) a.truncate();
return arg;
}
struct LProjector {
long lmax=3;
World& world;
std::map<std::pair<long,long>,complex_function_3d> Ylm;
LProjector(World& world) : world(world) {
for (long l=0; l<lmax; ++l) {
for (long m=-l; m<l+1; ++m) {
SphericalHarmonics Y(l,m);
auto lm=std::make_pair(l,m);
auto functor=[&Y](const coord_3d& r){return Y(r)*100.0*exp(-10.0*r.normf());};
Ylm[lm]=complex_factory_3d(world).functor(functor);
}
}
}
void analyze(const Spinor& f, const std::string text="") const {
madness::print(text);
for (int c=0; c<4; ++c) {
auto bla=f.components[c];
LProjector lproj(world);
lproj(bla,"component "+std::to_string(c));
}
}
void operator()(const complex_function_3d& f, const std::string text="") const {
madness::print(" "+text);
for (long l=0; l<lmax; ++l) {
for (long m = -l; m < l+1; ++m) {
auto lm=std::make_pair(l,m);
const complex_function_3d& Y=Ylm.find(lm)->second;
double_complex ovlp=inner(Y,f);
if (std::abs(ovlp)>1.e-7) madness::print("< lm | f> ",l,m,ovlp);
}
}
}
};
// The default constructor for functions does not initialize
// them to any value, but the solver needs functions initialized
// to zero for which we also need the world object.
struct spinorallocator {
World& world;
const int n;
/// @param[in] world the world
/// @param[in] nn the number of functions in a given vector
spinorallocator(World& world, const int n) : world(world), n(n) {
}
/// allocate a vector of n empty functions
std::vector<Spinor> operator()() {
std::vector<Spinor> r(n);
for (auto & rr : r) rr=Spinor(world);
return r;
}
};
/// class defining an operator in matrix form, fixed to size (4,4)
class MatrixOperator {
public:
typedef std::vector<std::pair<double_complex,std::shared_ptr<SCFOperatorBase<double_complex,3>>>> opT;
explicit MatrixOperator(const int i=4, const int j=4) {
elements.resize(i);
for (auto& e : elements) e.resize(j);
}
std::size_t nrow() const {return elements.size();}
std::size_t ncol() const {
if (elements.size()) return elements[0].size();
return 0;
}
virtual Spinor operator()(const Spinor& arg) const {
World& world=arg.components[0].world();
double_complex norm1=inner(arg,arg);
Spinor result;
for (auto& c : result.components) c=complex_factory_3d(world).compressed();
for (int i=0; i<ncol(); ++i) {
for (int j=0; j<nrow(); ++j) {
const opT& ops=elements[i][j];
for (const auto& op : ops) {
auto fac=op.first;
result.components[i]+=fac * (*op.second)(arg.components[j]);
}
}
}
double_complex norm2=inner(arg,arg);
result.truncate();
if (std::abs(norm1-norm2)/std::abs(norm1)>1.e-10) throw;
return result;
}
virtual std::vector<Spinor> operator()(const std::vector<Spinor>& arg) const {
std::vector<Spinor> result;
for (auto& a : arg) result.push_back((*this)(a));
return result;
}
/// add a submatrix to this
/// @param[in] istart row where to add the submatrix
/// @param[in] jstart column where to add the submatrix
void add_submatrix(int istart, int jstart, const MatrixOperator& submatrix) {
if (istart+submatrix.ncol()>ncol()) throw std::runtime_error("submatrix too large: too many columns");
if (jstart+submatrix.nrow()>nrow()) throw std::runtime_error("submatrix too large: too many rows");
for (int i=istart; i<istart+submatrix.ncol(); ++i) {
for (int j=jstart; j<jstart+submatrix.nrow(); ++j) {
for (auto& elem : submatrix.elements[i-istart][j-jstart]) elements[i][j].push_back(elem);
}
}
}
MatrixOperator& operator+=(const MatrixOperator& other) {
for (int i=0; i<ncol(); ++i) {
for (int j=0; j<nrow(); ++j) {
for (auto& elem : other.elements[i][j]) elements[i][j].push_back(elem);
}
}
return *this;
}
MatrixOperator operator+(const MatrixOperator& other) const {
MatrixOperator result;
for (int i=0; i<ncol(); ++i) {
for (int j=0; j<nrow(); ++j) {
for (auto& elem : this->elements[i][j]) result.elements[i][j].push_back(elem);
for (auto& elem : other.elements[i][j]) result.elements[i][j].push_back(elem);
}
}
return result;
}
void add_operator(const int i, const int j, const double_complex& fac, const std::shared_ptr<SCFOperatorBase<double_complex,3>> op) {
elements[i][j].push_back(std::make_pair(fac,op));
}
void print(std::string name="") const {
madness::print(name);
for (int i=0; i<ncol(); i++) {
for (int j=0; j<nrow(); j++) {
const opT& ops=elements[i][j];
for (auto op : ops) {
auto fac=op.first;
std::cout << " " << fac << " " << (op.second)->info();
}
std::cout << " ||| ";
}
std::cout << std::endl;
}
}
/// matrix containing prefactor and operator
std::vector<std::vector<opT>> elements;
};
class Metric : public MatrixOperator{
public:
double c0=1.0, c1=1.0, c2=1.0, c3=1.0;
virtual Spinor operator()(const Spinor& arg) const {
auto result=copy(arg);
result.components[0].scale(c0);
result.components[1].scale(c1);
result.components[2].scale(c2);
result.components[3].scale(c3);
return result;
}
void print() const {
madness::print("metric ",c0,c1,c2,c3);
}
};
Metric N_metric() {
Metric result;
result.c0=alpha1*alpha1;
result.c1=alpha1*alpha1;
return result;
}
Metric M_metric() {
Metric result;
result.c2=alpha1*alpha1;
result.c3=alpha1*alpha1;
return result;
}
void show_norms(const Spinor& bra, const Spinor& ket, const std::string& name) {
Metric m;
if (transform_c) m=M_metric();
auto norms = inner(ket.world(), m(bra).components, ket.components);
madness::print("norms of ",name,": ",norms);
}
MatrixOperator make_Hdiag(World& world, const LocalPotentialOperator<double_complex,3>& V1) {
MatrixOperator Hv;
Hv.add_operator(0,0, 1.0,std::make_shared<LocalPotentialOperator<double_complex,3>>(V1));
Hv.add_operator(1,1, 1.0,std::make_shared<LocalPotentialOperator<double_complex,3>>(V1));
Hv.add_operator(2,2, 1.0,std::make_shared<LocalPotentialOperator<double_complex,3>>(V1));
Hv.add_operator(3,3, 1.0,std::make_shared<LocalPotentialOperator<double_complex,3>>(V1));
return Hv;
}
/// this is the nuclear potential on the diagonal
MatrixOperator make_Hv(World& world, const double nuclear_charge) {
complex_function_3d V=complex_factory_3d(world)
.functor([&nuclear_charge](const coord_3d& r){return double_complex(-nuclear_charge/(r.normf()+1.e-13));});
auto V1=LocalPotentialOperator<double_complex,3>(world,"V",V);
return make_Hdiag(world,V1);
}
MatrixOperator make_sigma_sn(World& world, const Uplo& uplo,
const LocalPotentialOperator<double_complex, 3>& x,
const LocalPotentialOperator<double_complex, 3>& y,
const LocalPotentialOperator<double_complex, 3>& z,
const double factor) {
MatrixOperator H(4,4);
const double_complex ii=double_complex(0.0,1.0);
const double_complex one=double_complex(1.0,0.0);
if (uplo==upper) {
H.add_operator(0, 2, one * factor, std::make_shared<LocalPotentialOperator<double_complex, 3>>(z));
H.add_operator(0, 3, one * factor, std::make_shared<LocalPotentialOperator<double_complex, 3>>(x));
H.add_operator(0, 3, -ii * factor, std::make_shared<LocalPotentialOperator<double_complex, 3>>(y));
H.add_operator(1, 2, one * factor, std::make_shared<LocalPotentialOperator<double_complex, 3>>(x));
H.add_operator(1, 2, ii * factor, std::make_shared<LocalPotentialOperator<double_complex, 3>>(y));
H.add_operator(1, 3,-one * factor, std::make_shared<LocalPotentialOperator<double_complex, 3>>(z));
} else if (uplo==lower) {
H.add_operator(2, 0, one * factor, std::make_shared<LocalPotentialOperator<double_complex, 3>>(z));
H.add_operator(2, 1, one * factor, std::make_shared<LocalPotentialOperator<double_complex, 3>>(x));
H.add_operator(2, 1, -ii * factor, std::make_shared<LocalPotentialOperator<double_complex, 3>>(y));
H.add_operator(3, 0, one * factor, std::make_shared<LocalPotentialOperator<double_complex, 3>>(x));
H.add_operator(3, 0, ii * factor, std::make_shared<LocalPotentialOperator<double_complex, 3>>(y));
H.add_operator(3, 1,-one * factor, std::make_shared<LocalPotentialOperator<double_complex, 3>>(z));
}
return H;
}
/// returns a (4,4) hermition matrix with H=\vec \alpha \vec {xyz}
MatrixOperator make_alpha_sn(World& world,
const LocalPotentialOperator<double_complex, 3>& x,
const LocalPotentialOperator<double_complex, 3>& y,
const LocalPotentialOperator<double_complex, 3>& z,
const bool ll_negative) {
double fac=ll_negative ? -1.0 : 1.0;
MatrixOperator H(4,4);
H+=make_sigma_sn(world,upper,x,y,z,1.0);
H+=make_sigma_sn(world,lower,x,y,z,fac);
return H;
}
/// returns a (2,2) matrix: i(gamma_1 - 1) c 1/r \vec\sigma\vec n \vec \sigma l
MatrixOperator make_snsl(World& world,const double Z) {
MatrixOperator sp(2,2);
const double_complex ii=double_complex(0.0,1.0);
const double gamma1=compute_gamma(Z);
double_complex prefac=(gamma1-1.0)*ii/alpha1;
auto ox=n_times_n_dot_p_minus_p<double_complex,3>(world,0);
auto oy=n_times_n_dot_p_minus_p<double_complex,3>(world,1);
auto oz=n_times_n_dot_p_minus_p<double_complex,3>(world,2);
sp.add_operator(0,0, prefac,std::make_shared<n_times_n_dot_p_minus_p<double_complex,3>>(oz));
sp.add_operator(0,1, prefac,std::make_shared<n_times_n_dot_p_minus_p<double_complex,3>>(ox));
sp.add_operator(0,1,-ii*prefac,std::make_shared<n_times_n_dot_p_minus_p<double_complex,3>>(oy));
sp.add_operator(1,0, prefac,std::make_shared<n_times_n_dot_p_minus_p<double_complex,3>>(ox));
sp.add_operator(1,0, ii*prefac,std::make_shared<n_times_n_dot_p_minus_p<double_complex,3>>(oy));
sp.add_operator(1,1, -prefac,std::make_shared<n_times_n_dot_p_minus_p<double_complex,3>>(oz));
return sp;
}
/// returns a (2,2) matrix: Z/r \vec\sigma\vec l
MatrixOperator make_Zrsl(World& world, const double Z) {
MatrixOperator sp(2,2);
const double_complex ii=double_complex(0.0,1.0);
auto lx=AngularMomentum_with_step<double_complex,3>(world,0);
auto ly=AngularMomentum_with_step<double_complex,3>(world,1);
auto lz=AngularMomentum_with_step<double_complex,3>(world,2);
sp.add_operator(0,0, Z,std::make_shared<AngularMomentum_with_step<double_complex,3>>(lz));
sp.add_operator(0,1, Z,std::make_shared<AngularMomentum_with_step<double_complex,3>>(lx));
sp.add_operator(0,1,-ii*Z,std::make_shared<AngularMomentum_with_step<double_complex,3>>(ly));
sp.add_operator(1,0, Z,std::make_shared<AngularMomentum_with_step<double_complex,3>>(lx));
sp.add_operator(1,0, ii*Z,std::make_shared<AngularMomentum_with_step<double_complex,3>>(ly));
sp.add_operator(1,1, -Z,std::make_shared<AngularMomentum_with_step<double_complex,3>>(lz));
return sp;
}
/// Hv for ansatz 1
MatrixOperator make_Hv_reg1(World& world, const double nuclear_charge) {
const double_complex ii=double_complex(0.0,1.0);
const double c=1.0/alpha1;
double gamma=compute_gamma(nuclear_charge);
complex_function_3d x_div_r2_f=complex_factory_3d(world).functor([&gamma,&c,&ii](const coord_3d& r){return -c*ii*(gamma-1)*r[0]/(inner(r,r)+epsilon);});
complex_function_3d y_div_r2_f=complex_factory_3d(world).functor([&gamma,&c,&ii](const coord_3d& r){return -c*ii*(gamma-1)*r[1]/(inner(r,r)+epsilon);});
complex_function_3d z_div_r2_f=complex_factory_3d(world).functor([&gamma,&c,&ii](const coord_3d& r){return -c*ii*(gamma-1)*r[2]/(inner(r,r)+epsilon);});
auto x_div_r2=LocalPotentialOperator<double_complex,3>(world,"x/r2",x_div_r2_f);
auto y_div_r2=LocalPotentialOperator<double_complex,3>(world,"y/r2",y_div_r2_f);
auto z_div_r2=LocalPotentialOperator<double_complex,3>(world,"z/r2",z_div_r2_f);
return make_alpha_sn(world,x_div_r2,y_div_r2,z_div_r2,false);
}
/// Hv for ansatz 2