33#ifndef MADNESS_MRA_VMRA_H__INCLUDED
34#define MADNESS_MRA_VMRA_H__INCLUDED
133 template <
typename T, std::
size_t NDIM>
137 if (std::any_of(
v.begin(),
v.end(), [](
const Function<T,NDIM>&
f) {return not f.is_initialized();})) {
140 TreeState state=
v[0].get_impl()->get_tree_state();
141 for (
const auto&
f :
v) {
148 template <
typename T, std::
size_t NDIM>
161 template <
typename T, std::
size_t NDIM>
170 template <
typename T, std::
size_t NDIM>
176 template <
typename T, std::
size_t NDIM>
185 template <
typename T, std::
size_t NDIM>
195 template <
typename T, std::
size_t NDIM>
198 for (
const auto&
f : vf)
f.refine(
false);
205 template <
typename T, std::
size_t NDIM>
211 std::vector<FunctionImpl<T,NDIM>*> v_ptr;
214 for (
unsigned int i=0; i<vf.size(); ++i) {
215 if (vf[i].is_initialized()) v_ptr.push_back(vf[i].
get_impl().get());
219 std::sort(v_ptr.begin(),v_ptr.end());
220 typename std::vector<FunctionImpl<T, NDIM>*>::iterator it;
221 it = std::unique(v_ptr.begin(), v_ptr.end());
224 std::vector< Tensor<T> >
c(v_ptr.size());
225 v_ptr[0]->refine_to_common_level(v_ptr,
c, key0);
226 if (fence) v_ptr[0]->world.gop.fence();
228 for (
unsigned int i=0; i<vf.size(); i++) vf[i].
verify_tree();
232 template <
typename T, std::
size_t NDIM>
242 template <
typename T, std::
size_t NDIM>
255 template <
typename T, std::
size_t NDIM>
258 const bool fence=
true) {
260 if (
v.size()==0)
return v;
265 for (
const auto&
f :
v)
266 if (
f.is_initialized()) {
275 auto change_initial_to_intermediate =[](
const std::vector<Function<T,NDIM>>&
v,
280 if (
f.is_initialized() and
f.get_impl()->get_tree_state()==initialstate) {
281 f.change_tree_state(intermediatestate,
false);
315 template<
typename T, std::
size_t NDIM>
327 print(
"ensure_tree_state_respecting_fence failed");
328 throw std::runtime_error(
"ensure_tree_state_respecting_fence failed");
334 template <
typename T, std::
size_t NDIM>
346 vv.truncate(tol,
false);
355 template <
typename T, std::
size_t NDIM>
357 double tol=0.0,
bool fence=
true) {
358 if (
v.size()>0)
truncate(
v[0].world(),
v,tol,fence);
365 template <
typename T, std::
size_t NDIM>
367 double thresh=0.0,
bool fence=
true) {
368 if (
v.size()==0)
return v;
369 for (
auto& vv :
v) vv.reduce_rank(
thresh,
false);
370 if (fence)
v[0].world().gop.fence();
376 template <
typename T, std::
size_t NDIM>
377 std::vector< Function<T,NDIM> >
384 std::vector< Function<T,NDIM> > df(
v.size());
385 for (
unsigned int i=0; i<
v.size(); ++i) {
386 df[i] =
D(
v[i],
false);
393 template <
typename T, std::
size_t NDIM>
394 std::vector< Function<T,NDIM> >
396 std::vector< Function<T,NDIM> > r(n);
397 for (
int i=0; i<n; ++i) {
403 print(
"zero_functions_tree_state: unknown tree state");
404 throw std::runtime_error(
"zero_functions_tree_state: unknown tree state");
414 template <
typename T, std::
size_t NDIM>
415 std::vector< Function<T,NDIM> >
417 return zero_functions_tree_state<T,NDIM>(world,n,
reconstructed,fence);
421 template <
typename T, std::
size_t NDIM>
422 std::vector< Function<T,NDIM> >
424 return zero_functions_tree_state<T,NDIM>(world,n,
compressed,fence);
428 template <
typename T, std::
size_t NDIM>
429 std::vector< Function<T,NDIM> >
432 return zero_functions_tree_state<T,NDIM>(world,n,state,fence);
438 template<
typename T, std::
size_t NDIM>
440 if (vf_in.size()==0)
return std::vector<Function<T,NDIM>>();
441 World& world=vf_in.front().world();
442 auto vf=
copy(world,vf_in);
444 if (vf.size()==1)
return copy(world,vf_in);
449 for (
int i=0; i<s.dim(0); ++i)
Q(i,i) += 1.5;
456 for (
int i=0; i<
Q.dim(0); ++i)
457 for (
int j=0; j<i; ++j)
473 template <
typename T, std::
size_t NDIM>
477 if(
v.empty())
return v;
484 for(
size_t i=0;i<
v.size();++i) s(i)=1.0/(sqrt(s(i)));
488 for(
size_t i=0;i<
v.size();++i){
489 for(
size_t j=0;j<
v.size();++j){
496 World& world=
v.front().world();
503 template <
typename T, std::
size_t NDIM>
505 if(
v.empty())
return v;
508 World& world=
v.front().world();
518 template <
typename T, std::
size_t NDIM>
524 if(
v.empty())
return v;
529 lindep*=s(s.size()-1);
534 for(
size_t i=0;i<
v.size();++i) {
536 sqrts(i)=1.0/(sqrt(s(i)));
545 for(
size_t i=0;i<
v.size();++i){
546 for(
size_t j=0;j<
v.size();++j){
547 U(i,j)=U(i,j)*(sqrts(j));
552 World& world=
v.front().world();
560 template <
typename T, std::
size_t NDIM>
562 const double lindep){
563 if(
v.empty())
return v;
565 World& world=
v.front().world();
574 template <
typename T, std::
size_t NDIM>
579 if (
v.empty())
return v;
587 World& world=
v.front().world();
595 template <
typename T, std::
size_t NDIM>
597 if(
v.empty())
return v;
599 World& world=
v.front().world();
611 template <
typename T, std::
size_t NDIM>
626 std::vector<Function<T,NDIM> > pv(rank);
636 World& world=
v.front().world();
644 template <
typename T, std::
size_t NDIM>
654 template <
typename T, std::
size_t NDIM>
660 World& world=
v.front().world();
666 template <
typename T, std::
size_t NDIM>
668 std::vector<Function<T,NDIM> >
v=lhs;
669 for (std::size_t i = 0; i < rhs.size(); ++i)
v.push_back(rhs[i]);
673 template <
typename T, std::
size_t NDIM>
675 std::vector<Function<T,NDIM> >result;
676 for(
const auto& x:vv) result=
append(result,x);
680 template<
typename T, std::
size_t NDIM>
682 std::vector<std::shared_ptr<FunctionImpl<T,NDIM>>> result;
683 for (
auto&
f :
v) result.push_back(
f.get_impl());
687 template<
typename T, std::
size_t NDIM>
690 for (std::size_t i=0; i<vimpl.size(); ++i)
v[i].
set_impl(vimpl[i]);
693 template<
typename T, std::
size_t NDIM>
695 std::vector<Function<T,NDIM>>
v(vimpl.size());
696 for (std::size_t i=0; i<vimpl.size(); ++i)
v[i].
set_impl(vimpl[i]);
705 template <
typename T,
typename R, std::
size_t NDIM>
718 std::vector< Function<resultT,NDIM> > vc = zero_functions_compressed<resultT,NDIM>(world,
m);
721 for (
int i=0; i<
m; ++i) {
722 for (
int j=0; j<n; ++j) {
723 if (
c(j,i) !=
R(0.0)) vc[i].gaxpy(resultT(1.0),
v[j],resultT(
c(j,i)),
false);
734 template <
typename T,
typename R, std::
size_t NDIM>
750 vv.get_impl()->get_tree_state()==
reconstructed,
"trees have to be reconstructed in transform_reconstructed");
752 std::vector< Function<resultT,NDIM> > result = zero_functions<resultT,NDIM>(world,
m);
754 for (
int i=0; i<
m; ++i) {
756 for (
int j=0; j<n; ++j) {
757 if (
c(j,i) !=
R(0.0))
v[j].get_impl()->accumulate_trees(*(result[i].
get_impl()),resultT(
c(j,i)),
true);
765 for (
auto& r : result) r.get_impl()->finalize_sum();
773 template <
typename L,
typename R, std::
size_t NDIM>
781 = zero_functions_compressed<TENSOR_RESULT_TYPE(L,R),NDIM>(world,
c.dim(1));
784 vresult[0].vtransform(
v,
c, vresult, tol, fence);
788 template <
typename T,
typename R, std::
size_t NDIM>
804 c.copy_to_replicated(tmp);
807 std::vector< Function<resultT,NDIM> > vc = zero_functions_compressed<resultT,NDIM>(world,
m);
810 for (
int i=0; i<
m; ++i) {
811 for (
int j=0; j<n; ++j) {
812 if (tmp(j,i) !=
R(0.0)) vc[i].gaxpy(1.0,
v[j],tmp(j,i),
false);
822 template <
typename T,
typename Q, std::
size_t NDIM>
825 const std::vector<Q>& factors,
828 for (
unsigned int i=0; i<
v.size(); ++i)
v[i].
scale(factors[i],
false);
833 template <
typename T,
typename Q, std::
size_t NDIM>
839 for (
unsigned int i=0; i<
v.size(); ++i)
v[i].
scale(factor,
false);
844 template <
typename T, std::
size_t NDIM>
848 std::vector<double> norms(
v.size());
850 for (
unsigned int i=0; i<
v.size(); ++i) norms[i] =
v[i].norm2sq_local();
851 world.
gop.
sum(&norms[0], norms.size());
852 for (
unsigned int i=0; i<
v.size(); ++i) norms[i] = sqrt(norms[i]);
857 template <
typename T, std::
size_t NDIM>
862 for (
unsigned int i = 0; i <
v.size(); ++i) norms[i] =
v[i].norm2sq_local();
864 for (
unsigned int i = 0; i <
v.size(); ++i) norms[i] = sqrt(norms[i]);
870 template <
typename T, std::
size_t NDIM>
873 if (
v.size()==0)
return 0.0;
875 std::vector<double> norms(
v.size());
876 for (
unsigned int i=0; i<
v.size(); ++i) norms[i] =
v[i].norm2sq_local();
877 world.
gop.
sum(&norms[0], norms.size());
878 for (
unsigned int i=1; i<
v.size(); ++i) norms[0] += norms[i];
880 return sqrt(norms[0]);
926 template <
typename T, std::
size_t NDIM>
934 const int64_t n =
A.coldim();
935 const int64_t
m =
A.rowdim();
939 const int ichunk = 1000;
940 const int jchunk = 1000;
941 for (int64_t ilo=0; ilo<n; ilo+=ichunk) {
942 int64_t ihi = std::min(ilo + ichunk, n);
943 std::vector< Function<T,NDIM> > ivec(
f.begin()+ilo,
f.begin()+ihi);
944 for (int64_t jlo=0; jlo<
m; jlo+=jchunk) {
945 int64_t jhi = std::min(jlo + jchunk,
m);
946 std::vector< Function<T,NDIM> > jvec(
g.begin()+jlo,
g.begin()+jhi);
949 A.copy_from_replicated_patch(ilo, ihi - 1, jlo, jhi - 1,
P);
960 template <
typename T,
typename R, std::
size_t NDIM>
967 auto tensor_type = [](
const std::vector<Function<T,NDIM>>&
v) {
968 return v.front().get_impl()->get_tensor_type();
974 std::vector<const FunctionImpl<T,NDIM>*> left(
f.size());
975 std::vector<const FunctionImpl<R,NDIM>*> right(
g.size());
976 for (
unsigned int i=0; i<
f.size(); i++) left[i] =
f[i].
get_impl().get();
977 for (
unsigned int i=0; i<
g.size(); i++) right[i]=
g[i].
get_impl().get();
992 template <
typename T,
typename R, std::
size_t NDIM>
998 long n=
f.size(),
m=
g.size();
1004 if ((
void*)(&
f) != (
void*)(&
g))
compress(world,
g);
1006 for (
long i=0; i<n; ++i) {
1008 if (sym) jtop = i+1;
1009 for (
long j=0; j<jtop; ++j) {
1010 r(i,j) =
f[i].inner_local(
g[j]);
1011 if (sym) r(j,i) =
conj(r(i,j));
1036 template <
typename T,
typename R, std::
size_t NDIM>
1041 long n=
f.size(),
m=
g.size();
1045 auto tensor_type = [](
const std::vector<Function<T,NDIM>>&
v) {
1046 return v.front().get_impl()->get_tensor_type();
1052 for (
long i=0; i<n; ++i) r(i) =
f[i].inner_local(
g[i]);
1064 template <
typename T,
typename R, std::
size_t NDIM>
1072 auto tensor_type = [](
const std::vector<Function<T,NDIM>>&
v) {
1073 return v.front().get_impl()->get_tensor_type();
1076 f.change_tree_state(operating_state,
false);
1080 for (
long i=0; i<n; ++i) {
1081 r(i) =
f.inner_local(
g[i]);
1092 template <
typename T,
typename R, std::
size_t NDIM>
1096 if(
f.empty())
return 0.0;
1097 else return inner(
f[0].world(),
f,
g).sum();
1102 template <
typename T,
typename R, std::
size_t NDIM>
1109 a.reconstruct(
false);
1116 template <
typename T,
typename R, std::
size_t NDIM>
1124 a.reconstruct(
false);
1127 for (
unsigned int i=0; i<
v.size(); ++i) {
1128 v[i].norm_tree(
false);
1147 template <
typename T,
typename R, std::
size_t NDIM>
1154 bool symm =
false) {
1156 bool same=(&
f == &
g);
1160 for (
auto& ff :
f) ff.norm_tree(
false);
1161 if (not same)
for (
auto& gg :
g) gg.norm_tree(
false);
1164 std::vector<std::vector<Function<R,NDIM> > >result(
f.size());
1165 std::vector<Function<R,NDIM>> g_i;
1166 for (int64_t i=
f.size()-1; i>=0; --i) {
1168 result[i]=
vmulXX(
f[i],
g, tol,
false);
1170 if (g_i.empty()) g_i =
g;
1172 result[i]=
vmulXX(
f[i], g_i, tol,
false);
1180 template <
typename T, std::
size_t NDIM>
1186 for (
unsigned int i=0; i<
v.size(); ++i) {
1187 v[i].norm_tree(
false);
1193 template <
typename T,
typename R, std::
size_t NDIM>
1205 for (
unsigned int i=0; i<
a.size(); ++i) {
1206 q[i] =
mul(
a[i],
b[i],
false);
1219 template<
typename T, std::
size_t NDIM, std::
size_t LDIM>
1224 std::vector<Function<T,NDIM> > result(
g.size());
1225 for (
auto& r : result) r.set_impl(
f,
false);
1233 for (std::size_t i=0; i<result.size(); ++i) {
1240 for (
auto& ig :
g) ig.get_impl()->undo_redundant(
false);
1245 template<
typename T, std::
size_t NDIM, std::
size_t LDIM>
1247 const std::tuple<int,int,int>
v) {
1248 return partial_mul<T,NDIM,LDIM>(
f,
g,std::array<int,3>({std::get<0>(
v),std::get<1>(
v),std::get<2>(
v)}));
1253 template <
typename T, std::
size_t NDIM>
1254 std::vector< Function<T,NDIM> >
1258 return mul<T,T,NDIM>(world,
v,
v, fence);
1269 template <
typename T, std::
size_t NDIM>
1270 std::vector< Function<typename Tensor<T>::scalar_type,
NDIM> >
1276 std::vector<Function<scalartype,NDIM> > result(
v.size());
1277 for (
size_t i=0; i<
v.size(); ++i) result[i]=
abs_square(
v[i],
false);
1284 template <
typename T, std::
size_t NDIM>
1286 for (
unsigned int j=0; j<
v.size(); ++j) {
1287 v[j].set_thresh(
thresh,
false);
1293 template <
typename T, std::
size_t NDIM>
1294 std::vector< Function<T,NDIM> >
1299 std::vector< Function<T,NDIM> > r =
copy(world,
v);
1300 for (
unsigned int i=0; i<
v.size(); ++i) {
1308 template <
typename T,
typename R, std::
size_t NDIM>
1312 std::vector< Function<R,NDIM> > r(
v.size());
1313 for (
unsigned int i=0; i<
v.size(); ++i) {
1314 r[i] = convert<T,R,NDIM>(
v[i],
false);
1322 template <
typename T, std::
size_t NDIM>
1323 std::vector< Function<T,NDIM> >
1328 std::vector< Function<T,NDIM> > r(
v.size());
1329 for (
unsigned int i=0; i<
v.size(); ++i) {
1330 r[i] =
copy(
v[i],
false);
1338 template <
typename T, std::
size_t NDIM>
1339 std::vector< Function<T,NDIM> >
1342 std::vector< Function<T,NDIM> > r(
v.size());
1343 if (
v.size()>0) r=
copy(
v.front().world(),
v,fence);
1348 template <
typename T, std::
size_t NDIM>
1349 std::vector< Function<T,NDIM> >
1352 const unsigned int n,
1355 std::vector< Function<T,NDIM> > r(n);
1356 for (
unsigned int i=0; i<n; ++i) {
1357 r[i] =
copy(
v,
false);
1372 template <
typename T, std::
size_t NDIM>
1376 bool fence =
true) {
1378 std::vector<Function<T, NDIM>> r(
v.size());
1379 for (
unsigned int i = 0; i <
v.size(); ++i) {
1380 r[i] =
copy(
v[i], pmap,
false);
1387 template <
typename T,
typename R, std::
size_t NDIM>
1399 for (
unsigned int i=0; i<
a.size(); ++i) {
1400 r[i] =
add(
a[i],
b[i],
false);
1407 template <
typename T,
typename R, std::
size_t NDIM>
1418 for (
unsigned int i=0; i<
b.size(); ++i) {
1419 r[i] =
add(
a,
b[i],
false);
1424 template <
typename T,
typename R, std::
size_t NDIM>
1430 return add(world,
a,
b, fence);
1434 template <
typename T,
typename R, std::
size_t NDIM>
1446 for (
unsigned int i=0; i<
a.size(); ++i) {
1447 r[i] =
sub(
a[i],
b[i],
false);
1454 template <
typename T, std::
size_t NDIM>
1461 for (
unsigned int i=0; i<
f.size(); ++i) r.
gaxpy(1.0,
f[i],1.0,
false);
1466 template <
typename T, std::
size_t NDIM>
1474 const int64_t n =
A.coldim();
1475 const int64_t
m =
A.rowdim();
1479 const int ichunk = 1000;
1480 const int jchunk = 1000;
1481 for (int64_t ilo = 0; ilo < n; ilo += ichunk) {
1482 int64_t ihi = std::min(ilo + ichunk, n);
1483 std::vector<Function<T, NDIM>> ivec(
f.begin() + ilo,
f.begin() + ihi);
1484 for (int64_t jlo = 0; jlo <
m; jlo += jchunk) {
1485 int64_t jhi = std::min(jlo + jchunk,
m);
1486 std::vector<Function<T, NDIM>> jvec(
g.begin() + jlo,
g.begin() + jhi);
1489 A.copy_from_replicated_patch(ilo, ihi - 1, jlo, jhi - 1,
P);
1500 template <
typename T,
typename R, std::
size_t NDIM>
1511 std::vector<const FunctionImpl<T, NDIM>*> left(
f.size());
1512 std::vector<const FunctionImpl<R, NDIM>*> right(
g.size());
1513 for (
unsigned int i = 0; i <
f.size(); i++) left[i] =
f[i].
get_impl().get();
1514 for (
unsigned int i = 0; i <
g.size(); i++) right[i] =
g[i].
get_impl().get();
1529 template <
typename T,
typename R, std::
size_t NDIM>
1535 long n=
f.size(),
m=
g.size();
1541 if ((
void*)(&
f) != (
void*)(&
g))
compress(world,
g);
1543 for (
long i=0; i<n; ++i) {
1545 if (sym) jtop = i+1;
1546 for (
long j=0; j<jtop; ++j) {
1548 r(j,i) =
f[i].dot_local(
g[j]);
1550 r(i,j) =
conj(r(j,i));
1552 r(i,j) =
f[i].dot_local(
g[j]);
1563 template <
typename T,
typename R, std::
size_t NDIM>
1570 return sum(world,
mul(world,
a,
b,
true),fence);
1576 template <
typename T,
typename Q,
typename R, std::
size_t NDIM>
1586 if (
a.size()==0)
return std::vector<Function<resultT,NDIM> >();
1588 auto tensor_type = [](
const std::vector<Function<T,NDIM>>&
v) {
1589 return v.front().get_impl()->get_tensor_type();
1593 World& world=
a[0].world();
1594 std::vector<Function<resultT,NDIM> > result(
a.size());
1601 print(
"could not respect fence in gaxpy");
1618 template <
typename T,
typename Q,
typename R, std::
size_t NDIM>
1627 if (
a.size()==0)
return std::vector<Function<resultT,NDIM> >();
1629 World& world=
a[0].world();
1635 print(
"could not respect fence in gaxpy_oop");
1639 std::vector<Function<resultT,NDIM> > result(
a.size());
1640 for (
unsigned int i=0; i<
a.size(); ++i) {
1643 if (fence) world.gop.fence();
1649 template <
typename T,
typename Q,
typename R, std::
size_t NDIM>
1651 if (
a.size() == 0)
return;
1652 World& world=
a.front().world();
1657 template <
typename T,
typename Q,
typename R, std::
size_t NDIM>
1666 if (
a.empty())
return;
1668 auto tensor_type = [](
const std::vector<Function<T,NDIM>>&
v) {
1669 return v.front().get_impl()->get_tensor_type();
1673 bool do_in_reconstructed_state=tensor_type(
a)!=
TT_FULL;
1682 print(
"could not respect fence in gaxpy");
1696 print(
"could not respect fence in gaxpy for a");
1704 print(
"could not respect fence in gaxpy for b");
1711 for (
unsigned int i=0; i<
a.size(); ++i) {
1715 for (
unsigned int i=0; i<
a.size(); ++i)
a[i].
get_impl()->finalize_sum();
1723 template <
typename opT,
typename R, std::
size_t NDIM>
1726 const std::vector< std::shared_ptr<opT> >&
op,
1732 std::vector< Function<R,NDIM> >&
ncf = *
const_cast< std::vector< Function<R,NDIM>
>* >(&
f);
1738 for (
unsigned int i=0; i<
f.size(); ++i) {
1754 template <
typename T,
typename R, std::
size_t NDIM, std::
size_t KDIM>
1763 template <
typename T,
typename R, std::
size_t NDIM, std::
size_t KDIM>
1770 std::vector< Function<R,NDIM> >&
ncf = *
const_cast< std::vector< Function<R,NDIM>
>* >(&
f);
1777 if (
print_timings) printf(
"timer: %20.20s %8.2fs\n",
"make_nonstandard", wall1-wall0);
1780 for (
unsigned int i=0; i<
f.size(); ++i) {
1787 if (
op.destructive()) {
1788 for (
auto& ff :
ncf) ff.clear(
false);
1796 for (
auto& r : result) r.get_impl()->finalize_apply();
1800 for (
auto& r : result) r.get_impl()->print_timer();
1809 template <
typename T, std::
size_t NDIM>
1812 std::vector<double> nn =
norm2s(world,
v);
1813 for (
unsigned int i=0; i<
v.size(); ++i)
v[i].
scale(1.0/nn[i],
false);
1817 template <
typename T, std::
size_t NDIM>
1820 if(world.
rank()==0) std::cout <<
"print_size: " << msg <<
" is empty" << std::endl;
1821 }
else if(
v.size()==1){
1822 v.front().print_size(msg);
1831 template <
typename T, std::
size_t NDIM>
1835 if (x.is_initialized()) size+=x.size_local();
1837 const double d=
sizeof(
T);
1838 const double fac=1024*1024*1024;
1843 template <
typename T, std::
size_t NDIM>
1850 template <
typename T, std::
size_t NDIM>
1853 if (
v.empty())
return 0.0;
1855 const double d=
sizeof(
T);
1856 const double fac=1024*1024*1024;
1859 for(
unsigned int i=0;i<
v.size();i++){
1860 if (
v[i].is_initialized()) size+=
v[i].size();
1868 template <
typename T, std::
size_t NDIM>
1870 const double d=
sizeof(
T);
1871 const double fac=1024*1024*1024;
1872 double size=
f.size();
1881 template <
typename T,
typename opT, std::
size_t NDIM>
1884 const bool fence=
true) {
1887 World& world=vin[0].world();
1890 std::vector<Function<T,NDIM> > vout=zero_functions<T,NDIM>(world,
op.get_result_size());
1891 for (
auto& out : vout) out.set_impl(vin[0],
false);
1902 template <
typename T, std::
size_t NDIM>
1910 template <
typename T, std::
size_t NDIM>
1918 template <
typename T, std::
size_t NDIM>
1926 template <
typename T, std::
size_t NDIM>
1934 template <
typename T, std::
size_t NDIM>
1942 template <
typename T, std::
size_t NDIM>
1950 template <
typename T,
typename R, std::
size_t NDIM>
1954 std::vector<Function<T,NDIM> > tmp=
copy(rhs[0].world(),rhs);
1961 template <
typename T,
typename R, std::
size_t NDIM>
1973 template <
typename T,
typename R, std::
size_t NDIM>
1976 if (
v.size()>0)
return mul(
v[0].world(),
a,
v,
true);
1982 template <
typename T,
typename R, std::
size_t NDIM>
1985 if (
v.size()>0)
return mul(
v[0].world(),
a,
v,
true);
1990 template <
typename T, std::
size_t NDIM>
1993 if (lhs.size() > 0)
gaxpy(lhs.front().world(), 1.0, lhs, 1.0, rhs);
1997 template <
typename T, std::
size_t NDIM>
2001 if (lhs.size() > 0)
gaxpy(lhs.front().world(), 1.0, lhs, -1.0, rhs);
2006 template <
typename T, std::
size_t NDIM>
2007 std::vector<Function<typename Tensor<T>::scalar_type,
NDIM> >
2009 std::vector<Function<typename Tensor<T>::scalar_type,
NDIM> > result(
v.size());
2010 for (std::size_t i=0; i<
v.size(); ++i) result[i]=
real(
v[i],
false);
2011 if (fence and result.size()>0) result[0].world().gop.fence();
2016 template <
typename T, std::
size_t NDIM>
2017 std::vector<Function<typename Tensor<T>::scalar_type,
NDIM> >
2019 std::vector<Function<typename Tensor<T>::scalar_type,
NDIM> > result(
v.size());
2020 for (std::size_t i=0; i<
v.size(); ++i) result[i]=
imag(
v[i],
false);
2021 if (fence and result.size()>0) result[0].world().gop.fence();
2032 template <
typename T, std::
size_t NDIM>
2034 bool refine=
false,
bool fence=
true) {
2040 std::vector< std::shared_ptr< Derivative<T,NDIM> > >
grad=
2041 gradient_operator<T,NDIM>(world);
2043 std::vector<Function<T,NDIM> > result(
NDIM);
2044 for (
size_t i=0; i<
NDIM; ++i) result[i]=
apply(*(
grad[i]),
f,
false);
2050 template <
typename T, std::
size_t NDIM>
2052 bool refine=
false,
bool fence=
true) {
2058 std::vector< std::shared_ptr< Derivative<T,NDIM> > >
grad=
2059 gradient_operator<T,NDIM>(world);
2062 for (
unsigned int i=0; i<
NDIM; ++i) (*
grad[i]).set_ble1();
2064 std::vector<Function<T,NDIM> > result(
NDIM);
2065 for (
unsigned int i=0; i<
NDIM; ++i) result[i]=
apply(*(
grad[i]),
f,
false);
2071 template <
typename T, std::
size_t NDIM>
2073 bool refine=
false,
bool fence=
true) {
2079 std::vector< std::shared_ptr< Derivative<T,NDIM> > >
grad=
2080 gradient_operator<T,NDIM>(world);
2083 for (
unsigned int i=0; i<
NDIM; ++i) (*
grad[i]).set_ble2();
2085 std::vector<Function<T,NDIM> > result(
NDIM);
2086 for (
unsigned int i=0; i<
NDIM; ++i) result[i]=
apply(*(
grad[i]),
f,
false);
2092 template <
typename T, std::
size_t NDIM>
2094 bool refine=
false,
bool fence=
true) {
2100 std::vector< std::shared_ptr< Derivative<T,NDIM> > >
grad=
2101 gradient_operator<T,NDIM>(world);
2104 for (
unsigned int i=0; i<
NDIM; ++i) (*
grad[i]).set_bspline1();
2106 std::vector<Function<T,NDIM> > result(
NDIM);
2107 for (
unsigned int i=0; i<
NDIM; ++i) result[i]=
apply(*(
grad[i]),
f,
false);
2113 template <
typename T, std::
size_t NDIM>
2115 bool refine=
false,
bool fence=
true) {
2121 std::vector< std::shared_ptr< Derivative<T,NDIM> > >
grad=
2122 gradient_operator<T,NDIM>(world);
2125 for (
unsigned int i=0; i<
NDIM; ++i) (*
grad[i]).set_bspline2();
2127 std::vector<Function<T,NDIM> > result(
NDIM);
2128 for (
unsigned int i=0; i<
NDIM; ++i) result[i]=
apply(*(
grad[i]),
f,
false);
2134 template <
typename T, std::
size_t NDIM>
2136 bool refine=
false,
bool fence=
true) {
2142 std::vector< std::shared_ptr< Derivative<T,NDIM> > >
grad=
2143 gradient_operator<T,NDIM>(world);
2146 for (
unsigned int i=0; i<
NDIM; ++i) (*
grad[i]).set_bspline3();
2148 std::vector<Function<T,NDIM> > result(
NDIM);
2149 for (
unsigned int i=0; i<
NDIM; ++i) result[i]=
apply(*(
grad[i]),
f,
false);
2164 template <
typename T, std::
size_t NDIM>
2166 bool do_refine=
false,
bool fence=
true) {
2169 World& world=
v[0].world();
2171 if (do_refine)
refine(world,
v);
2173 std::vector< std::shared_ptr< Derivative<T,NDIM> > >
grad=
2174 gradient_operator<T,NDIM>(world);
2176 std::vector<Function<T,NDIM> > result(
NDIM);
2177 for (
size_t i=0; i<
NDIM; ++i) result[i]=
apply(*(
grad[i]),
v[i],
false);
2179 return sum(world,result,fence);
2190 template <
typename T, std::
size_t NDIM>
2192 bool do_refine=
false,
bool fence=
true) {
2195 World& world=
v[0].world();
2197 if (do_refine)
refine(world,
v);
2199 std::vector< std::shared_ptr< Derivative<T,NDIM> > >
grad=
2200 gradient_operator<T,NDIM>(world);
2202 std::vector<Function<T,NDIM> >
d(
NDIM),dd(
NDIM);
2214 d[0].gaxpy(1.0,dd[0],-1.0,
false);
2215 d[1].gaxpy(1.0,dd[1],-1.0,
false);
2216 d[2].gaxpy(1.0,dd[2],-1.0,
false);
2231 template <
typename T,
typename R, std::
size_t NDIM>
2234 bool do_refine=
false,
bool fence=
true) {
2238 World& world=
f[0].world();
2244 d[0]=
mul(
f[1],
g[2],
false);
2245 d[1]=
mul(
f[2],
g[0],
false);
2246 d[2]=
mul(
f[0],
g[1],
false);
2248 dd[0]=
mul(
f[2],
g[1],
false);
2249 dd[1]=
mul(
f[0],
g[2],
false);
2250 dd[2]=
mul(
f[1],
g[0],
false);
2256 d[0].gaxpy(1.0,dd[0],-1.0,
false);
2257 d[1].gaxpy(1.0,dd[1],-1.0,
false);
2258 d[2].gaxpy(1.0,dd[2],-1.0,
false);
2265 template<
typename T, std::
size_t NDIM>
2282 template<
typename T,
size_t NDIM>
2284 const std::string
name) {
2285 if (world.
rank()==0)
print(
"loading vector of functions",
name);
2287 std::size_t fsize=0;
2290 for (std::size_t i=0; i<fsize; ++i) ar &
f[i];
2294 template<
typename T,
size_t NDIM>
2297 World& world=
f.front().world();
2298 if (world.
rank()==0)
print(
"saving vector of functions",
name);
2300 std::size_t fsize=
f.size();
2302 for (std::size_t i=0; i<fsize; ++i) ar &
f[i];
double q(double t)
Definition DKops.h:18
Definition test_ar.cc:118
long size() const
Returns the number of elements in the tensor.
Definition basetensor.h:138
Implements derivatives operators with variety of boundary conditions on simulation domain.
Definition derivative.h:266
Definition distributed_matrix.h:68
Manages data associated with a row/column/block distributed array.
Definition distributed_matrix.h:388
static void redistribute(World &world, const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > &newpmap)
Sets the default process map and redistributes all functions using the old map.
Definition funcdefaults.h:414
static TensorType get_tensor_type()
Returns the default tensor type.
Definition funcdefaults.h:324
FunctionFactory implements the named-parameter idiom for Function.
Definition function_factory.h:86
FunctionFactory & compressed(bool value=true)
Definition function_factory.h:168
FunctionImpl holds all Function state to facilitate shallow copy semantics.
Definition funcimpl.h:945
static Tensor< TENSOR_RESULT_TYPE(T, R) > inner_local(const std::vector< const FunctionImpl< T, NDIM > * > &left, const std::vector< const FunctionImpl< R, NDIM > * > &right, bool sym)
Definition funcimpl.h:6000
void undo_redundant(const bool fence)
convert this from redundant to standard reconstructed form
Definition mraimpl.h:1534
void multiply(const implT *f, const FunctionImpl< T, LDIM > *g, const int particle)
multiply f (a pair function of NDIM) with an orbital g (LDIM=NDIM/2)
Definition funcimpl.h:3580
void change_tree_state(const TreeState finalstate, bool fence=true)
change the tree state of this function, might or might not respect fence!
Definition mraimpl.h:1403
static Tensor< TENSOR_RESULT_TYPE(T, R)> dot_local(const std::vector< const FunctionImpl< T, NDIM > * > &left, const std::vector< const FunctionImpl< R, NDIM > * > &right, bool sym)
Definition funcimpl.h:6052
FunctionNode holds the coefficients, etc., at each node of the 2^NDIM-tree.
Definition funcimpl.h:127
coeffT & coeff()
Returns a non-const reference to the tensor containing the coeffs.
Definition funcimpl.h:227
A multiresolution adaptive numerical function.
Definition mra.h:139
World & world() const
Returns the world.
Definition mra.h:688
Function< T, NDIM > & gaxpy(const T &alpha, const Function< Q, NDIM > &other, const R &beta, bool fence=true)
Inplace, general bi-linear operation in wavelet basis. No communication except for optional fence.
Definition mra.h:1025
void set_impl(const std::shared_ptr< FunctionImpl< T, NDIM > > &impl)
Replace current FunctionImpl with provided new one.
Definition mra.h:661
void multi_to_multi_op_values(const opT &op, const std::vector< Function< T, NDIM > > &vin, std::vector< Function< T, NDIM > > &vout, const bool fence=true)
apply op on the input vector yielding an output vector of functions
Definition mra.h:1585
bool is_initialized() const
Returns true if the function is initialized.
Definition mra.h:167
long size() const
Definition lowranktensor.h:482
Key is the index for a node of the 2^NDIM-tree.
Definition key.h:69
std::shared_ptr< WorldDCPmapInterface< keyT > > load_balance(double fac=1.0, bool printstuff=false)
Actually does the partitioning of the tree.
Definition lbdeux.h:320
void add_tree(const Function< T, NDIM > &f, const costT &costfn, bool fence=false)
Accumulates cost from a function.
Definition lbdeux.h:289
Convolutions in separated form (including Gaussian)
Definition operator.h:139
A slice defines a sub-range or patch of a dimension.
Definition slice.h:103
A tensor is a multidimensional array.
Definition tensor.h:317
TensorTypeData< T >::scalar_type scalar_type
C++ typename of the real type associated with a complex type.
Definition tensor.h:409
T * ptr()
Returns a pointer to the internal data.
Definition tensor.h:1825
A simple, fixed dimension vector.
Definition vector.h:64
Interface to be provided by any process map.
Definition worlddc.h:82
void fence(bool debug=false)
Synchronizes all processes in communicator AND globally ensures no pending AM or tasks.
Definition worldgop.cc:161
void sum(T *buf, size_t nelem)
Inplace global sum while still processing AM & tasks.
Definition worldgop.h:870
void fence()
Returns after all local tasks have completed.
Definition world_task_queue.h:1384
A parallel world class.
Definition world.h:132
WorldTaskQueue & taskq
Task queue.
Definition world.h:206
ProcessID rank() const
Returns the process rank in this World (same as MPI_Comm_rank()).
Definition world.h:320
WorldGopInterface & gop
Global operations.
Definition world.h:207
An archive for storing local or parallel data wrapping a BinaryFstreamOutputArchive.
Definition parallel_archive.h:321
int integer
Definition crayio.c:25
static const double R
Definition csqrt.cc:46
Declaration and initialization of tree traversal functions and generic derivative.
static double lo
Definition dirac-hatom.cc:23
auto T(World &world, response_space &f) -> response_space
Definition global_functions.cc:34
const double beta
Definition gygi_soltion.cc:62
static const double v
Definition hatom_sf_dirac.cc:20
Tensor< double > op(const Tensor< double > &x)
Definition kain.cc:508
#define rot(x, k)
Definition lookup3.c:72
#define MADNESS_CHECK(condition)
Check a condition — even in a release build the condition is always evaluated so it can have side eff...
Definition madness_exception.h:182
#define MADNESS_ASSERT(condition)
Assert a condition that should be free of side-effects since in release builds this might be a no-op.
Definition madness_exception.h:134
#define MADNESS_CHECK_THROW(condition, msg)
Check a condition — even in a release build the condition is always evaluated so it can have side eff...
Definition madness_exception.h:207
Main include file for MADNESS and defines Function interface.
static const bool VERIFY_TREE
Definition mra.h:57
Namespace for all elements and tools of MADNESS.
Definition DFParameters.h:10
void save_function(const std::vector< Function< T, NDIM > > &f, const std::string name)
save a vector of functions
Definition vmra.h:2295
bool ensure_tree_state_respecting_fence(const std::vector< Function< T, NDIM > > &v, const TreeState state, bool fence)
ensure v has the requested tree state, change the tree state of v if necessary and no fence is given
Definition vmra.h:316
void rr_cholesky(Tensor< T > &A, typename Tensor< T >::scalar_type tol, Tensor< integer > &piv, int &rank)
Compute the rank-revealing Cholesky factorization.
Definition lapack.cc:1203
void make_redundant(World &world, const std::vector< Function< T, NDIM > > &v, bool fence=true)
change tree_state of a vector of functions to redundant
Definition vmra.h:186
std::vector< Function< T, NDIM > > orthonormalize_rrcd(const std::vector< Function< T, NDIM > > &v, Tensor< T > &ovlp, const double tol, Tensor< integer > &piv, int &rank)
Definition vmra.h:612
Function< double, NDIM > abssq(const Function< double_complex, NDIM > &z, bool fence=true)
Returns a new function that is the square of the absolute value of the input.
Definition mra.h:2778
Function< TENSOR_RESULT_TYPE(L, R), NDIM > gaxpy_oop(TENSOR_RESULT_TYPE(L, R) alpha, const Function< L, NDIM > &left, TENSOR_RESULT_TYPE(L, R) beta, const Function< R, NDIM > &right, bool fence=true)
Returns new function alpha*left + beta*right optional fence and no automatic compression.
Definition mra.h:1981
Function< typename TensorTypeData< Q >::scalar_type, NDIM > abs_square(const Function< Q, NDIM > &func)
Definition complexfun.h:121
Function< T, NDIM > square(const Function< T, NDIM > &f, bool fence=true)
Create a new function that is the square of f - global comm only if not reconstructed.
Definition mra.h:2746
response_space scale(response_space a, double b)
Function< TENSOR_RESULT_TYPE(L, R), NDIM > sub(const Function< L, NDIM > &left, const Function< R, NDIM > &right, bool fence=true)
Same as operator- but with optional fence and no automatic compression.
Definition mra.h:2035
std::vector< Function< T, NDIM > > reduce_rank(std::vector< Function< T, NDIM > > v, double thresh=0.0, bool fence=true)
reduces the tensor rank of the coefficient tensor (if applicable)
Definition vmra.h:366
Tensor< double > norm2s_T(World &world, const std::vector< Function< T, NDIM > > &v)
Computes the 2-norms of a vector of functions.
Definition vmra.h:858
std::vector< double > norm2s(World &world, const std::vector< Function< T, NDIM > > &v)
Computes the 2-norms of a vector of functions.
Definition vmra.h:845
std::vector< Function< T, NDIM > > grad_bspline_one(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
Definition vmra.h:2093
void set_impl(std::vector< Function< T, NDIM > > &v, const std::vector< std::shared_ptr< FunctionImpl< T, NDIM > > > vimpl)
Definition vmra.h:688
Function< Q, NDIM > convert(const Function< T, NDIM > &f, bool fence=true)
Type conversion implies a deep copy. No communication except for optional fence.
Definition mra.h:2096
Function< TENSOR_RESULT_TYPE(Q, T), NDIM > mul(const Q alpha, const Function< T, NDIM > &f, bool fence=true)
Returns new function equal to alpha*f(x) with optional fence.
Definition mra.h:1765
void truncate(World &world, response_space &v, double tol, bool fence)
Definition basic_operators.cc:30
std::vector< std::shared_ptr< FunctionImpl< T, NDIM > > > get_impl(const std::vector< Function< T, NDIM > > &v)
Definition vmra.h:681
Function< T, NDIM > div(const std::vector< Function< T, NDIM > > &v, bool do_refine=false, bool fence=true)
shorthand div operator
Definition vmra.h:2165
std::vector< Function< T, NDIM > > orthonormalize_cd(const std::vector< Function< T, NDIM > > &v, Tensor< T > &ovlp)
Definition vmra.h:575
void norm_tree(World &world, const std::vector< Function< T, NDIM > > &v, bool fence=true)
Makes the norm tree for all functions in a vector.
Definition vmra.h:1181
tensorT Q2(const tensorT &s)
Given overlap matrix, return rotation with 2nd order error to orthonormalize the vectors.
Definition SCF.cc:135
std::vector< Function< TENSOR_RESULT_TYPE(T, R), NDIM > > transform(World &world, const std::vector< Function< T, NDIM > > &v, const Tensor< R > &c, bool fence=true)
Transforms a vector of functions according to new[i] = sum[j] old[j]*c[j,i].
Definition vmra.h:707
std::vector< Function< TENSOR_RESULT_TYPE(T, R), NDIM > > cross(const std::vector< Function< T, NDIM > > &f, const std::vector< Function< R, NDIM > > &g, bool do_refine=false, bool fence=true)
shorthand cross operator
Definition vmra.h:2232
TreeState
Definition funcdefaults.h:59
@ nonstandard_after_apply
s and d coeffs, state after operator application
Definition funcdefaults.h:64
@ redundant_after_merge
s coeffs everywhere, must be summed up to yield the result
Definition funcdefaults.h:66
@ reconstructed
s coeffs at the leaves only
Definition funcdefaults.h:60
@ nonstandard
s and d coeffs in internal nodes
Definition funcdefaults.h:62
@ unknown
Definition funcdefaults.h:68
@ compressed
d coeffs in internal nodes, s and d coeffs at the root
Definition funcdefaults.h:61
@ redundant
s coeffs everywhere
Definition funcdefaults.h:65
@ nonstandard_with_leaves
like nonstandard, with s coeffs at the leaves
Definition funcdefaults.h:63
Function< T, NDIM > conj(const Function< T, NDIM > &f, bool fence=true)
Return the complex conjugate of the input function with the same distribution and optional fence.
Definition mra.h:2110
void cholesky(Tensor< T > &A)
Compute the Cholesky factorization.
Definition lapack.cc:1174
std::vector< std::vector< Function< TENSOR_RESULT_TYPE(T, R), NDIM > > > matrix_mul_sparse(World &world, const std::vector< Function< R, NDIM > > &f, const std::vector< Function< R, NDIM > > &g, double tol, bool fence=true, bool symm=false)
Outer product of a vector of functions with a vector of functions using sparsity.
Definition vmra.h:1149
void standard(World &world, std::vector< Function< T, NDIM > > &v, bool fence=true)
Generates standard form of a vector of functions.
Definition vmra.h:243
void compress(World &world, const std::vector< Function< T, NDIM > > &v, bool fence=true)
Compress a vector of functions.
Definition vmra.h:149
const std::vector< Function< T, NDIM > > & reconstruct(const std::vector< Function< T, NDIM > > &v)
reconstruct a vector of functions
Definition vmra.h:162
response_space transpose(response_space &f)
Definition basic_operators.cc:10
std::vector< Function< T, NDIM > > impl2function(const std::vector< std::shared_ptr< FunctionImpl< T, NDIM > > > vimpl)
Definition vmra.h:694
std::vector< Function< T, NDIM > > grad_bpsline_two(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
Definition vmra.h:2114
void set_thresh(World &world, std::vector< Function< T, NDIM > > &v, double thresh, bool fence=true)
Sets the threshold in a vector of functions.
Definition vmra.h:1285
double norm2(World &world, const std::vector< Function< T, NDIM > > &v)
Computes the 2-norm of a vector of functions.
Definition vmra.h:871
std::vector< Function< T, NDIM > > flatten(const std::vector< std::vector< Function< T, NDIM > > > &vv)
Definition vmra.h:674
std::vector< Function< T, NDIM > > orthonormalize_canonical(const std::vector< Function< T, NDIM > > &v, const Tensor< T > &ovlp, double lindep)
Definition vmra.h:519
std::vector< CCPairFunction< T, NDIM > > operator*(const double fac, const std::vector< CCPairFunction< T, NDIM > > &arg)
Definition ccpairfunction.h:1089
static void verify_tree(World &world, const std::vector< Function< T, NDIM > > &v)
Definition SCF.cc:72
std::vector< Function< T, NDIM > > multi_to_multi_op_values(const opT &op, const std::vector< Function< T, NDIM > > &vin, const bool fence=true)
apply op on the input vector yielding an output vector of functions
Definition vmra.h:1882
static const Slice _(0,-1, 1)
Tensor< T > inverse(const Tensor< T > &a_in)
invert general square matrix A
Definition lapack.cc:832
void load_balance(const real_function_6d &f, const bool leaf)
do some load-balancing
Definition madness/chem/mp2.cc:70
TreeState get_tree_state(const Function< T, NDIM > &f)
get tree state of a function
Definition mra.h:2794
std::vector< CCPairFunction< T, NDIM > > operator-(const std::vector< CCPairFunction< T, NDIM > > c1, const std::vector< CCPairFunction< T, NDIM > > &c2)
Definition ccpairfunction.h:1060
Function< TENSOR_RESULT_TYPE(L, R), NDIM > mul_sparse(const Function< L, NDIM > &left, const Function< R, NDIM > &right, double tol, bool fence=true)
Sparse multiplication — left and right must be reconstructed and if tol!=0 have tree of norms already...
Definition mra.h:1806
std::vector< Function< T, NDIM > > partial_mul(const Function< T, NDIM > f, const std::vector< Function< T, LDIM > > g, const int particle)
multiply a high-dimensional function with a low-dimensional function
Definition vmra.h:1220
Function< T, NDIM > gaxpy_oop_reconstructed(const double alpha, const Function< T, NDIM > &left, const double beta, const Function< T, NDIM > &right, const bool fence=true)
Returns new function alpha*left + beta*right optional fence, having both addends reconstructed.
Definition mra.h:1999
std::vector< Function< TENSOR_RESULT_TYPE(T, R), NDIM > > transform_reconstructed(World &world, const std::vector< Function< T, NDIM > > &v, const Tensor< R > &c, bool fence=true)
Transforms a vector of functions according to new[i] = sum[j] old[j]*c[j,i].
Definition vmra.h:736
void print(const T &t, const Ts &... ts)
Print items to std::cout (items separated by spaces) and terminate with a new line.
Definition print.h:225
response_space apply(World &world, std::vector< std::vector< std::shared_ptr< real_convolution_3d > > > &op, response_space &f)
Definition basic_operators.cc:39
std::vector< Function< T, NDIM > > append(const std::vector< Function< T, NDIM > > &lhs, const std::vector< Function< T, NDIM > > &rhs)
combine two vectors
Definition vmra.h:667
@ TT_2D
Definition gentensor.h:120
@ TT_FULL
Definition gentensor.h:120
void refine(World &world, const std::vector< Function< T, NDIM > > &vf, bool fence=true)
refine the functions according to the autorefine criteria
Definition vmra.h:196
void print_size(World &world, const std::vector< Function< T, NDIM > > &v, const std::string &msg="vectorfunction")
Definition vmra.h:1818
NDIM & f
Definition mra.h:2481
Function< TENSOR_RESULT_TYPE(L, R), NDIM > add(const Function< L, NDIM > &left, const Function< R, NDIM > &right, bool fence=true)
Same as operator+ but with optional fence and no automatic compression.
Definition mra.h:1991
const Function< T, NDIM > & change_tree_state(const Function< T, NDIM > &f, const TreeState finalstate, bool fence=true)
change tree state of a function
Definition mra.h:2807
Function< TENSOR_RESULT_TYPE(T, R), NDIM > dot(World &world, const std::vector< Function< T, NDIM > > &a, const std::vector< Function< R, NDIM > > &b, bool fence=true)
Multiplies and sums two vectors of functions r = \sum_i a[i] * b[i].
Definition vmra.h:1565
std::vector< CCPairFunction< T, NDIM > > & operator-=(std::vector< CCPairFunction< T, NDIM > > &rhs, const std::vector< CCPairFunction< T, NDIM > > &lhs)
Definition ccpairfunction.h:1082
std::vector< Function< T, NDIM > > orthonormalize(const std::vector< Function< T, NDIM > > &vf_in)
orthonormalize the vectors
Definition vmra.h:439
NDIM const Function< R, NDIM > & g
Definition mra.h:2481
double wall_time()
Returns the wall time in seconds relative to an arbitrary origin.
Definition timers.cc:48
std::vector< Function< T, NDIM > > grad_ble_one(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
Definition vmra.h:2051
Function< TENSOR_RESULT_TYPE(typename opT::opT, R), NDIM > apply_only(const opT &op, const Function< R, NDIM > &f, bool fence=true)
Apply operator ONLY in non-standard form - required other steps missing !!
Definition mra.h:2184
std::vector< Function< T, NDIM > > grad_bspline_three(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
Definition vmra.h:2135
std::vector< Function< T, NDIM > > zero_functions_compressed(World &world, int n, bool fence=true)
Generates a vector of zero functions (compressed)
Definition vmra.h:423
double inner(response_space &a, response_space &b)
Definition response_functions.h:442
double imag(double x)
Definition complexfun.h:56
void load_function(World &world, std::vector< Function< T, NDIM > > &f, const std::string name)
load a vector of functions
Definition vmra.h:2283
std::vector< Function< TENSOR_RESULT_TYPE(L, R), D > > vmulXX(const Function< L, D > &left, const std::vector< Function< R, D > > &vright, double tol, bool fence=true)
Use the vmra/mul(...) interface instead.
Definition mra.h:1871
void refine_to_common_level(World &world, std::vector< Function< T, NDIM > > &vf, bool fence=true)
refine all functions to a common (finest) level
Definition vmra.h:206
std::vector< Function< T, NDIM > > grad_ble_two(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
Definition vmra.h:2072
static bool print_timings
Definition SCF.cc:104
void normalize(World &world, std::vector< Function< T, NDIM > > &v, bool fence=true)
Normalizes a vector of functions — v[i] = v[i].scale(1.0/v[i].norm2())
Definition vmra.h:1810
std::vector< Function< T, NDIM > > zero_functions(World &world, int n, bool fence=true)
Generates a vector of zero functions (reconstructed)
Definition vmra.h:416
std::vector< CCPairFunction< T, NDIM > > operator+(const std::vector< CCPairFunction< T, NDIM > > c1, const std::vector< CCPairFunction< T, NDIM > > &c2)
Definition ccpairfunction.h:1052
std::vector< Function< T, NDIM > > grad(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
shorthand gradient operator
Definition vmra.h:2033
Function< T, NDIM > multiply(const Function< T, NDIM > f, const Function< T, LDIM > g, const int particle, const bool fence=true)
multiply a high-dimensional function with a low-dimensional function
Definition mra.h:2437
std::vector< Function< T, NDIM > > zero_functions_auto_tree_state(World &world, int n, bool fence=true)
Generates a vector of zero functions, either compressed or reconstructed, depending on tensor type.
Definition vmra.h:430
DistributedMatrix< T > matrix_dot(const DistributedMatrixDistribution &d, const std::vector< Function< T, NDIM > > &f, const std::vector< Function< T, NDIM > > &g, bool sym=false)
Definition vmra.h:1467
std::vector< Function< T, NDIM > > orthonormalize_symmetric(const std::vector< Function< T, NDIM > > &v, const Tensor< T > &ovlp)
symmetric orthonormalization (see e.g. Szabo/Ostlund)
Definition vmra.h:474
std::vector< Function< T, NDIM > > zero_functions_tree_state(World &world, int n, const TreeState state, bool fence=true)
Generates a vector of zero functions with a given tree state.
Definition vmra.h:395
double real(double x)
Definition complexfun.h:52
std::vector< CCPairFunction< T, NDIM > > & operator+=(std::vector< CCPairFunction< T, NDIM > > &lhs, const CCPairFunction< T, NDIM > &rhs)
Definition ccpairfunction.h:1068
Tensor< TENSOR_RESULT_TYPE(T, R) > matrix_dot_old(World &world, const std::vector< Function< T, NDIM > > &f, const std::vector< Function< R, NDIM > > &g, bool sym=false)
Computes the matrix dot product of two function vectors - q(i,j) = dot(f[i],g[j])
Definition vmra.h:1530
std::string name(const FuncType &type, const int ex=-1)
Definition ccpairfunction.h:28
void matrix_inner(DistributedMatrix< T > &A, const std::vector< Function< T, NDIM > > &f, const std::vector< Function< T, NDIM > > &g, bool sym=false)
Definition distpm.cc:46
double get_size(World &world, const std::vector< Function< T, NDIM > > &v)
Definition vmra.h:1851
double get_size_local(World &world, const std::vector< Function< T, NDIM > > &v)
return the size of a vector of functions for each rank
Definition vmra.h:1832
Function< T, NDIM > copy(const Function< T, NDIM > &f, const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > &pmap, bool fence=true)
Create a new copy of the function with different distribution and optional fence.
Definition mra.h:2066
void syev(const Tensor< T > &A, Tensor< T > &V, Tensor< typename Tensor< T >::scalar_type > &e)
Real-symmetric or complex-Hermitian eigenproblem.
Definition lapack.cc:969
Tensor< TENSOR_RESULT_TYPE(T, R) > matrix_inner_old(World &world, const std::vector< Function< T, NDIM > > &f, const std::vector< Function< R, NDIM > > &g, bool sym=false)
Computes the matrix inner product of two function vectors - q(i,j) = inner(f[i],g[j])
Definition vmra.h:993
void make_nonstandard(World &world, std::vector< Function< T, NDIM > > &v, bool fence=true)
Generates non-standard form of a vector of functions.
Definition vmra.h:233
void gaxpy(const double a, ScalarResult< T > &left, const double b, const T &right, const bool fence=true)
the result type of a macrotask must implement gaxpy
Definition macrotaskq.h:245
int distance(const madness::Hash_private::HashIterator< hashT > &it, const madness::Hash_private::HashIterator< hashT > &jt)
Definition worldhashmap.h:616
static long abs(long a)
Definition tensor.h:218
static const double b
Definition nonlinschro.cc:119
static const double d
Definition nonlinschro.cc:121
static const double a
Definition nonlinschro.cc:118
double Q(double a)
Definition relops.cc:20
static const double c
Definition relops.cc:10
static const double m
Definition relops.cc:9
static const double L
Definition rk.cc:46
static const double thresh
Definition rk.cc:45
Definition test_ar.cc:204
double operator()(const Key< 6 > &key, const FunctionNode< double, 6 > &node) const
Definition mp2.h:70
Definition lowrankfunction.h:332
Definition dirac-hatom.cc:108
std::string ok(const bool b)
Definition test6.cc:43
AtomicInt sum
Definition test_atomicint.cc:46
int P
Definition test_binsorter.cc:9
static const double alpha
Definition testcosine.cc:10
constexpr std::size_t NDIM
Definition testgconv.cc:54
#define TENSOR_RESULT_TYPE(L, R)
This macro simplifies access to TensorResultType.
Definition type_data.h:205
#define PROFILE_FUNC
Definition worldprofile.h:209
#define PROFILE_BLOCK(name)
Definition worldprofile.h:208