32#ifndef MADNESS_MRA_MRAIMPL_H__INCLUDED
33#define MADNESS_MRA_MRAIMPL_H__INCLUDED
36#error "mraimpl.h should ONLY be included in one of the mraX.cc files (x=1..6)"
53 bool isnan(
const std::complex<T>&
v) {
69 template <
typename T, std::
size_t NDIM>
71 if (!
two_scale_hg(
k, &hg))
throw "failed to get twoscale coefficients";
78 h1 =
copy(hg(sk,sk2));
79 g0 =
copy(hg(sk2,sk));
89 template <
typename T, std::
size_t NDIM>
99 for (
int mu=0;
mu<npt; ++
mu) {
102 for (
int j=0; j<
k; ++j) {
103 quad_phi(
mu,j) = phi[j];
104 quad_phiw(
mu,j) = quad_w(
mu)*phi[j];
110 template <
typename T, std::
size_t NDIM>
118 template <
typename T, std::
size_t NDIM>
122 for (
const auto& [key, node] : coeffs) {
124 if (key.level() > 0) {
127 if (pit == coeffs.end()) {
128 print(world.rank(),
"FunctionImpl: verify: MISSING PARENT",key,parent);
133 const nodeT& pnode = pit->second;
134 if (!pnode.has_children()) {
135 print(world.rank(),
"FunctionImpl: verify: PARENT THINKS IT HAS NO CHILDREN",key,parent);
143 typename dcT::const_iterator cit = coeffs.find(kit.key()).get();
144 if (cit == coeffs.end()) {
145 if (node.has_children()) {
146 print(world.rank(),
"FunctionImpl: verify: MISSING CHILD",key,kit.key());
153 if (! node.has_children()) {
154 print(world.rank(),
"FunctionImpl: verify: UNEXPECTED CHILD",key,kit.key());
168 template<
typename T, std::
size_t NDIM>
174 auto check_internal_coeff_size = [&state, &
k](
const coeffT&
c) {
176 return c.dim(0)==2*
k;
182 auto check_leaf_coeff_size = [&state, &
k](
const coeffT&
c) {
193 for (
const auto& [key, node] : coeffs) {
194 const auto&
c=node.coeff();
195 const bool is_internal=node.has_children();
196 const bool is_leaf=not node.has_children();
197 if (is_internal) good=good and check_internal_coeff_size(
c);
198 if (is_leaf) good=good and check_leaf_coeff_size(
c);
200 print(
"incorrect size of coefficients for key",key,
"state",state,
c.dim(0));;
206 template <
typename T, std::
size_t NDIM>
208 return coeffs.get_pmap();
222 template <
typename T, std::
size_t NDIM>
224 const double beta,
const implT&
g,
const bool fence) {
229 ProcessID owner = coeffs.owner(cdata.key0);
230 if (world.rank() == owner) {
238 apply_opT apply_op(
this);
240 woT::task(world.rank(), &implT:: template forward_traverse<coeff_opT,apply_opT>,
241 coeff_op, apply_op, cdata.key0);
245 if (fence) world.gop.fence();
249 template <
typename T, std::
size_t NDIM>
255 template <
typename T, std::
size_t NDIM>
261 template <
typename T, std::
size_t NDIM>
267 template <
typename T, std::
size_t NDIM>
272 template <
typename T, std::
size_t NDIM>
277 template <
typename T, std::
size_t NDIM>
282 template <
typename T, std::
size_t NDIM>
287 template <
typename T, std::
size_t NDIM>
292 template <
typename T, std::
size_t NDIM>
299 template <
typename T, std::
size_t NDIM>
305 template <
typename T, std::
size_t NDIM>
311 template <
typename T, std::
size_t NDIM>
318 template <
typename T, std::
size_t NDIM>
321 template <
typename T, std::
size_t NDIM>
324 template <
typename T, std::
size_t NDIM>
327 template <
typename T, std::
size_t NDIM>
330 template <
typename T, std::
size_t NDIM>
333 template <
typename T, std::
size_t NDIM>
336 template <
typename T, std::
size_t NDIM>
339 template <
typename T, std::
size_t NDIM>
342 template <
typename T, std::
size_t NDIM>
345 template <
typename T, std::
size_t NDIM>
348 template <
typename T, std::
size_t NDIM>
351 template <
typename T, std::
size_t NDIM>
353 timer_accumulate.accumulate(time);
356 template <
typename T, std::
size_t NDIM>
358 if (world.rank()==0) {
359 timer_accumulate.print(
"accumulate");
360 timer_target_driven.print(
"target_driven");
361 timer_lr_result.print(
"result2low_rank");
365 template <
typename T, std::
size_t NDIM>
367 if (world.rank()==0) {
368 timer_accumulate.reset();
369 timer_target_driven.reset();
370 timer_lr_result.reset();
377 template <
typename T, std::
size_t NDIM>
382 if (world.rank() == coeffs.owner(cdata.key0)) {
383 if (is_compressed()) {
384 truncate_spawn(cdata.key0,tol);
386 truncate_reconstructed_spawn(cdata.key0,tol);
393 template <
typename T, std::
size_t NDIM>
402 template <
typename T, std::
size_t NDIM>
409 std::vector<Tensor<double> > localinfo_vec(1,localinfo);
410 std::vector<Tensor<double> > printinfo=world.gop.concat0(localinfo_vec);
414 if (world.rank()==0) do_print_plane(
filename,printinfo,xaxis,yaxis,el2);
422 template <
typename T, std::
size_t NDIM>
425 user_to_sim<NDIM>(el2,x_sim);
435 const keyT& key = it->first;
436 const nodeT& node = it->second;
444 double scale=std::pow(0.5,
double(n));
445 double xloleft =
scale*l[xaxis];
446 double yloleft =
scale*l[yaxis];
447 double xhiright =
scale*(l[xaxis]+1);
448 double yhiright =
scale*(l[yaxis]+1);
459 if ((user[0]<-5.0) or (user[1]<-5.0) or (user[2]>5.0) or (user[3]>5.0))
continue;
465 const double maxrank=40;
471 const int npt = cdata.npt + 1;
482 plotinfo(counter,0)=color;
483 plotinfo(counter,1)=user[0];
484 plotinfo(counter,2)=user[1];
485 plotinfo(counter,3)=user[2];
486 plotinfo(counter,4)=user[3];
493 else plotinfo=plotinfo(
Slice(0,counter-1),
Slice(
_));
498 template <
typename T, std::
size_t NDIM>
500 const int xaxis,
const int yaxis,
const coordT el2) {
507 pFile = fopen(
filename.c_str(),
"w");
511 fprintf(pFile,
"\\psset{unit=1cm}\n");
512 fprintf(pFile,
"\\begin{pspicture}(%4.2f,%4.2f)(%4.2f,%4.2f)\n",
515 fprintf(pFile,
"\\pslinewidth=0.1pt\n");
517 for (std::vector<
Tensor<double> >::const_iterator it=plotinfo.begin(); it!=plotinfo.end(); ++it) {
522 for (
long i=0; i<localinfo.
dim(0); ++i) {
524 fprintf(pFile,
"\\newhsbcolor{mycolor}{%8.4f 1.0 0.7}\n",localinfo(i,0));
525 fprintf(pFile,
"\\psframe["
528 "(%12.8f,%12.8f)(%12.8f,%12.8f)\n",
529 localinfo(i,1),localinfo(i,2),localinfo(i,3),localinfo(i,4));
535 fprintf(pFile,
"\\end{pspicture}\n");
541 template <
typename T, std::
size_t NDIM>
545 std::vector<keyT> local_keys=local_leaf_keys();
548 std::vector<keyT> all_keys=world.gop.concat0(local_keys);
552 if (world.rank()==0) do_print_grid(
filename,all_keys);
557 template <
typename T, std::
size_t NDIM>
561 std::vector<keyT> keys(coeffs.size());
568 const keyT& key = it->first;
569 const nodeT& node = it->second;
570 if (node.
is_leaf()) keys[i++]=key;
583 template <
typename T, std::
size_t NDIM>
590 const size_t npt = qx.
dim(0);
593 long npoints=power<NDIM>(npt);
595 long nboxes=keys.size();
599 pFile = fopen(
filename.c_str(),
"w");
601 fprintf(pFile,
"%ld\n",npoints*nboxes);
602 fprintf(pFile,
"%ld points per box and %ld boxes \n",npoints,nboxes);
605 typename std::vector<keyT>::const_iterator key_it=keys.begin();
606 for (key_it=keys.begin(); key_it!=keys.end(); ++key_it) {
608 const keyT& key=*key_it;
609 fprintf(pFile,
"# key: %8d",key.
level());
616 const double h = std::pow(0.5,
double(n));
623 for (
size_t i=0; i<npt; ++i) {
624 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
625 for (
size_t j=0; j<npt; ++j) {
626 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
627 for (
size_t k=0;
k<npt; ++
k) {
628 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
634 fprintf(pFile,
"%18.12f %18.12f %18.12f\n",
c[0],
c[1],
c[2]);
648 template <
typename T, std::
size_t NDIM>
654 const int MAXLEVEL1 = 20;
655 const int MAXLEVEL2 = 10;
662 return tol*std::min(1.0,
pow(0.5,
double(std::min(key.
level(),MAXLEVEL1)))*
L);
666 return tol*std::min(1.0,
pow(0.25,
double(std::min(key.
level(),MAXLEVEL2)))*
L*
L);
683 const static double fac=1.0/std::pow(2,
NDIM*0.5);
687 return tol*std::min(1.0,
pow(0.5,
double(std::min(key.
level(),MAXLEVEL1)))*
L);
695 template <
typename T, std::
size_t NDIM>
697 std::vector<Slice> s(
NDIM);
699 for (std::size_t i=0; i<
NDIM; ++i)
700 s[i] = cdata.s[l[i]&1];
706 template <
typename T, std::
size_t NDIM>
708 const keyT& child,
const keyT& parent,
const coeffT& coeff)
const {
718 if (coeff.
dim(0)==2*
f->get_k()) result=coeff;
719 else if (coeff.
dim(0)==
f->get_k()) {
720 result(
f->cdata.s0)+=coeff;
729 const coeffT scoeff=
f->parent_to_child(coeff,parent,child);
730 result(
f->cdata.s0)+=scoeff;
738 template <
typename T, std::
size_t NDIM>
743 for (
typename dcT::iterator it= coeffs.begin(); it!=end; ++it) {
745 nodeT& node=it->second;
746 if (key.
level()>max_level) coeffs.erase(key);
749 this->undo_redundant(
true);
754 template <
typename T, std::
size_t NDIM>
768 template <
typename T, std::
size_t NDIM>
770 const std::vector<tensorT>&
c,
772 if (key == cdata.key0 && coeffs.owner(key)!=world.rank())
return;
775 std::unique_ptr<typename dcT::accessor[]> acc(
new typename dcT::accessor[
v.size()]);
776 for (
unsigned int i=0; i<
c.size(); i++) {
779 bool exists = !
v[i]->coeffs.insert(acc[i],key);
791 for (
unsigned int i=0; i<
v.size(); i++) {
792 done &= acc[i]->second.has_coeff();
797 std::vector<tensorT>
d(
v.size());
798 for (
unsigned int i=0; i<
v.size(); i++) {
799 if (acc[i]->second.has_coeff()) {
802 s(cdata.s0) = acc[i]->second.coeff().full_tensor();
803 acc[i]->second.clear_coeff();
805 acc[i]->second.set_has_children(
true);
811 const keyT& child = kit.key();
812 std::vector<Slice> cp = child_patch(child);
813 std::vector<tensorT> childc(
v.size());
814 for (
unsigned int i=0; i<
v.size(); i++) {
815 if (
d[i].size()) childc[i] =
copy(
d[i](cp));
817 woT::task(coeffs.owner(child), &implT::refine_to_common_level,
v, childc, child);
823 template <
typename T, std::
size_t NDIM>
825 if (world.size()> 1000)
828 box_interior[from] = ni;
832 template <
typename T, std::
size_t NDIM>
834 if (world.size() >= 1000)
836 for (
int i=0; i<world.size(); ++i)
837 box_leaf[i] = box_interior[i] == 0;
839 long nleaf=0, ninterior=0;
842 const nodeT& node = it->second;
848 this->send(0, &implT::put_in_box, world.rank(), nleaf, ninterior);
850 if (world.rank() == 0) {
851 for (
int i=0; i<world.size(); ++i) {
852 printf(
"load: %5d %8ld %8ld\n", i, box_leaf[i], box_interior[i]);
858 template <
typename T, std::
size_t NDIM>
864 template <
typename T, std::
size_t NDIM>
868 double test = 2*
lo*hi + hi*hi;
875 template <
typename T, std::
size_t NDIM>
878 coeffs.insert(acc,key);
879 nodeT& node = acc->second;
901 const keyT& child = kit.key();
902 if (
d.size() > 0) ss =
copy(
d(child_patch(child)));
904 woT::task(coeffs.owner(child), &implT::sum_down_spawn, child, ss);
909 if (
c.size() <= 0)
c =
coeffT(cdata.vk,targs);
914 template <
typename T, std::
size_t NDIM>
917 if (world.rank() == coeffs.owner(cdata.key0)) sum_down_spawn(cdata.key0,
coeffT());
918 if (fence) world.gop.fence();
922 template <
typename T, std::
size_t NDIM>
926 const std::pair<keyT,coeffT>& left,
927 const std::pair<keyT,coeffT>& center,
928 const std::pair<keyT,coeffT>& right) {
929 D->forward_do_diff1(
f,
this,key,left,center,right);
933 template <
typename T, std::
size_t NDIM>
937 const std::pair<keyT,coeffT>& left,
938 const std::pair<keyT,coeffT>& center,
939 const std::pair<keyT,coeffT>& right) {
940 D->do_diff1(
f,
this,key,left,center,right);
945 template <
typename T, std::
size_t NDIM>
947 typedef std::pair<keyT,coeffT> argT;
948 for (
const auto& [key, node]:
f->coeffs) {
949 if (node.has_coeff()) {
951 argT center(key,node.coeff());
959 if (fence) world.gop.fence();
964 template <
typename T, std::
size_t NDIM>
977 template <
typename T, std::
size_t NDIM>
985 std::vector<long> vkhalf=std::vector<long>(
NDIM/2,cdata.vk[0]);
1012 template <
typename T, std::
size_t NDIM>
1019 if (ket_only)
return coeff_ket;
1022 coeffT val_ket=coeffs2values(key,coeff_ket);
1040 coeff_result=
coeffT(values2coeffs(key,val_ket2),this->get_tensor_args());
1045 val_ket=val_ket.
convert(get_tensor_args());
1047 coeff_result=values2coeffs(key,val_result);
1051 return coeff_result;
1056 template <
typename T, std::
size_t NDIM>
1060 const_cast<implT*
>(&
f)->flo_unary_op_node_inplace(
do_mapdim(map,*
this),fence);
1065 template <
typename T, std::
size_t NDIM>
1068 const_cast<implT*
>(&
f)->flo_unary_op_node_inplace(
do_mirror(mirrormap,*
this),fence);
1075 template <
typename T, std::
size_t NDIM>
1077 const std::vector<long>&
mirror,
bool fence) {
1087 template <
typename T, std::
size_t NDIM>
1091 this->scale_inplace(0.5,
true);
1098 template <
typename T, std::
size_t NDIM>
1106 template <
typename T, std::
size_t NDIM>
1114 template <
typename T, std::
size_t NDIM>
1116 std::list<keyT> to_be_erased;
1117 for (
auto it=coeffs.begin(); it!=coeffs.end(); ++it) {
1118 const keyT& key=it->first;
1120 if (key.
level()==n) node.set_is_leaf(
true);
1121 if (key.
level()>n) to_be_erased.push_back(key);
1123 for (
auto& key : to_be_erased) coeffs.erase(key);
1130 template <
typename T, std::
size_t NDIM>
1133 flo_unary_op_node_inplace(
1151 template <
typename T, std::
size_t NDIM>
1159 template <
typename T, std::
size_t NDIM>
1180 template <
typename T, std::
size_t NDIM>
1188 template <
typename T, std::
size_t NDIM>
1200 template <
typename T, std::
size_t NDIM>
1206 const tensorT h[2] = {cdata.h0T, cdata.h1T};
1214 for (
size_t ii=0; ii<
NDIM; ++ii) matrices[ii]=
h[kit.key().translation()[ii]%2];
1230 template <
typename T, std::
size_t NDIM>
1235 const tensorT h[2] = {cdata.h0, cdata.h1};
1248 template <
typename T, std::
size_t NDIM>
1250 long kmin = std::min(cdata.k,old.
cdata.k);
1251 std::vector<Slice> s(
NDIM,
Slice(0,kmin-1));
1254 const keyT& key = it->first;
1255 const nodeT& node = it->second;
1259 coeffs.replace(key,
nodeT(
c,
false));
1262 coeffs.replace(key,nodeT(coeffT(),
true));
1269 template <
typename T, std::
size_t NDIM>
1271 return coeffs.probe(key) && coeffs.find(key).get()->second.has_children();
1274 template <
typename T, std::
size_t NDIM>
1276 return coeffs.probe(key) && (not coeffs.find(key).get()->second.has_children());
1280 template <
typename T, std::
size_t NDIM>
1282 for (
unsigned int i=0; i<
v.size(); ++i) {
1291 template <
typename T, std::
size_t NDIM>
1294 for (
typename dcT::iterator it=coeffs.begin(); it!=end; ++it) {
1295 it->second.set_norm_tree(0.0);
1296 it->second.set_snorm(0.0);
1297 it->second.set_dnorm(0.0);
1302 template <
typename T, std::
size_t NDIM>
1305 for (
typename dcT::iterator it=coeffs.begin(); it!=end; ++it) {
1306 const keyT& key = it->first;
1308 const auto found = coeffs.find(acc,key);
1310 nodeT& node = acc->second;
1318 int ndir =
static_cast<int>(std::pow(
static_cast<double>(3),
static_cast<int>(
NDIM)));
1319 std::vector< Future <bool> >
v = future_vector_factory<bool>(ndir);
1324 for (std::size_t
d=0;
d<
NDIM; ++
d) {
1332 keyT neigh = neighbor(key,
keyT(key.
level(),l), is_periodic);
1335 v[i++] = this->
task(coeffs.owner(neigh), &implT::exists_and_has_children, neigh);
1341 woT::task(world.rank(), &implT::broaden_op, key,
v);
1353 template <
typename T, std::
size_t NDIM>
1356 if (world.rank() == coeffs.owner(cdata.key0))
1357 woT::task(world.rank(), &implT::trickle_down_op, cdata.key0,
coeffT());
1358 if (fence) world.gop.fence();
1364 template <
typename T, std::
size_t NDIM>
1375 if (it == coeffs.end()) {
1377 it = coeffs.find(key).get();
1384 if (node.coeff().has_no_data()) node.coeff()=
coeffT(cdata.vk,targs);
1387 if (node.has_children()) {
1389 if (key.level() > 0)
d += s;
1392 const keyT& child = kit.key();
1396 woT::task(coeffs.owner(child), &implT::trickle_down_op, child, ss);
1401 node.coeff().reduce_rank(
thresh);
1406 template <
typename T, std::
size_t NDIM>
1410 if (current_state==finalstate)
return;
1421 remove_internal_coefficients(fence);
1425 remove_internal_coefficients(fence);
1442 remove_leaf_coefficients(fence);
1460 print(
"could not respect no-fence parameter in change_tree_state");
1467 template <
typename T, std::
size_t NDIM>
1470 if (is_reconstructed())
return;
1472 if (is_redundant() or is_nonstandard_with_leaves()) {
1474 this->remove_internal_coefficients(fence);
1478 if (world.rank() == coeffs.owner(cdata.key0))
1479 woT::task(world.rank(), &implT::reconstruct_op, cdata.key0,coeffT(),
true);
1480 }
else if (is_nonstandard()) {
1483 if (world.rank() == coeffs.owner(cdata.key0))
1484 woT::task(world.rank(), &implT::reconstruct_op, cdata.key0,coeffT(),
false);
1488 if (fence) world.gop.fence();
1499 template <
typename T, std::
size_t NDIM>
1503 set_tree_state(newstate);
1506 bool redundant1=(tree_state==
redundant);
1508 if (world.rank() == coeffs.owner(cdata.key0)) {
1510 compress_spawn(cdata.key0, nonstandard1, keepleaves1, redundant1);
1516 template <
typename T, std::
size_t NDIM>
1518 flo_unary_op_node_inplace(remove_internal_coeffs(),fence);
1521 template <
typename T, std::
size_t NDIM>
1523 flo_unary_op_node_inplace(remove_leaf_coeffs(),fence);
1527 template <
typename T, std::
size_t NDIM>
1531 if (is_redundant())
return;
1532 MADNESS_CHECK_THROW(is_reconstructed(),
"impl::make_redundant() wants a reconstructed tree");
1537 template <
typename T, std::
size_t NDIM>
1546 template <
typename T, std::
size_t NDIM>
1548 if (world.rank() == coeffs.owner(cdata.key0))
1549 norm_tree_spawn(cdata.key0);
1554 template <
typename T, std::
size_t NDIM>
1560 double value =
v[i].get();
1564 coeffs.task(key, &nodeT::set_norm_tree,
sum);
1569 template <
typename T, std::
size_t NDIM>
1571 nodeT& node = coeffs.find(key).get()->second;
1573 std::vector< Future<double> >
v = future_vector_factory<double>(1<<
NDIM);
1576 v[i] = woT::task(coeffs.owner(kit.key()), &implT::norm_tree_spawn, kit.key());
1578 return woT::task(world.rank(),&implT::norm_tree_op, key,
v);
1592 template <
typename T, std::
size_t NDIM>
1595 nodeT& node = coeffs.find(key).get()->second;
1602 std::vector<Future<coeffT> >
v = future_vector_factory<coeffT>(1<<
NDIM);
1605 v[i] = woT::task(coeffs.owner(kit.key()), &implT::truncate_reconstructed_spawn, kit.key(),tol,
TaskAttributes::hipri());
1616 template <
typename T, std::
size_t NDIM>
1623 for (
size_t i=0; i<
v.size(); ++i)
if (
v[i].get().has_no_data())
return coeffT();
1630 const auto found = coeffs.find(acc, key);
1636 d(child_patch(kit.key())) +=
v[i].get().full_tensor();
1642 const double error=
d.normf();
1644 nodeT& node = coeffs.find(key).get()->second;
1646 if (
error < truncate_tol(tol,key)) {
1649 coeffs.erase(kit.key());
1653 acc->second.set_coeff(ss);
1667 template <
typename T, std::
size_t NDIM>
1669 const std::vector<
Future<std::pair<coeffT,double> > >&
v,
bool nonstandard1) {
1677 double norm_tree2=0.0;
1680 d(child_patch(kit.key())) +=
v[i].get().first.full_tensor();
1681 norm_tree2+=
v[i].get().second*
v[i].get().second;
1686 timer_filter.accumulate(cpu1-cpu0);
1690 const auto found = coeffs.find(acc, key);
1692 MADNESS_CHECK_THROW(!acc->second.has_coeff(),
"compress_op: existing coeffs where there should be none");
1700 double snorm=ss.
normf();
1702 if (key.
level()> 0 && !nonstandard1)
d(cdata.s0) = 0.0;
1705 double dnorm=dd.
normf();
1708 acc->second.set_snorm(snorm);
1709 acc->second.set_dnorm(dnorm);
1712 acc->second.set_coeff(dd);
1714 timer_compress_svd.accumulate(cpu1-cpu0);
1717 return std::make_pair(ss,snorm);
1726 template <
typename T, std::
size_t NDIM>
1727 std::pair<typename FunctionImpl<T,NDIM>::coeffT,
double>
1732 double norm_tree2=0.0;
1734 d(child_patch(kit.key())) +=
v[i].get().first.full_tensor();
1735 norm_tree2+=
v[i].get().second*
v[i].get().second;
1747 double dnorm=
d.normf();
1748 double snorm=s.normf();
1751 const auto found = coeffs.find(acc, key);
1754 acc->second.set_coeff(s);
1755 acc->second.set_dnorm(dnorm);
1756 acc->second.set_snorm(snorm);
1764 template <
typename T, std::
size_t NDIM>
1767 if (is_compressed())
return;
1769 flo_unary_op_node_inplace(do_standard(
this),fence);
1777 template <
typename T, std::
size_t NDIM>
1782 tight_args.thresh*=0.01;
1785 flo_unary_op_node_inplace(do_consolidate_buffer(tight_args),
true);
1787 if (printme) printf(
"time in consolidate_buffer %8.4f\n",end1-begin1);
1793 flo_unary_op_node_inplace(do_reduce_rank(targs),
true);
1795 if (printme) printf(
"time in do_reduce_rank %8.4f\n",end1-begin1);
1799 flo_unary_op_node_inplace(do_change_tensor_type(targs,*
this),
true);
1801 if (printme) printf(
"time in do_change_tensor_type %8.4f\n",end1-begin1);
1805 flo_unary_op_node_inplace(do_truncate_NS_leafs(
this),
true);
1807 if (printme) printf(
"time in do_truncate_NS_leafs %8.4f\n",end1-begin1);
1810 double elapsed=end-begin;
1820 template <
typename T, std::
size_t NDIM>
1829 template <
typename T, std::
size_t NDIM>
1841 template <
typename T, std::
size_t NDIM>
1843 std::size_t maxdepth = 0;
1846 std::size_t
N = (std::size_t) it->first.level();
1855 template <
typename T, std::
size_t NDIM>
1857 std::size_t maxdepth = max_local_depth();
1858 world.gop.max(maxdepth);
1863 template <
typename T, std::
size_t NDIM>
1865 std::size_t maxsize = 0;
1866 maxsize = coeffs.size();
1867 world.gop.max(maxsize);
1872 template <
typename T, std::
size_t NDIM>
1874 std::size_t minsize = 0;
1875 minsize = coeffs.size();
1876 world.gop.min(minsize);
1881 template <
typename T, std::
size_t NDIM>
1883 std::size_t
sum = 0;
1884 sum = coeffs.size();
1890 template <
typename T, std::
size_t NDIM>
1892 std::size_t
sum = 0;
1893 for (
const auto& [key,node] : coeffs) {
1894 if (node.has_coeff())
sum+=node.size();
1900 template <
typename T, std::
size_t NDIM>
1902 std::size_t
sum = size_local();
1908 template <
typename T, std::
size_t NDIM>
1910 std::size_t
sum = coeffs.size() * (
sizeof(
keyT) +
sizeof(
nodeT));
1913 const nodeT& node = it->second;
1921 template <
typename T, std::
size_t NDIM>
1924 for (
auto& [key,node] : coeffs) {
1925 if (node.has_coeff())
sum+=node.coeff().nCoeff();
1931 template <
typename T, std::
size_t NDIM>
1933 std::size_t
sum = nCoeff_local();
1940 template <
typename T, std::
size_t NDIM>
1942 const size_t tsize=this->tree_size();
1944 const size_t ncoeff=this->nCoeff();
1946 const double d=
sizeof(
T);
1947 const double fac=1024*1024*1024;
1951 double local = norm2sq_local();
1952 this->world.gop.sum(local);
1953 this->world.gop.fence();
1957 if (this->world.rank()==0) {
1959 constexpr std::size_t
bufsize=128;
1961 snprintf(buf,
bufsize,
"%40s at time %.1fs: norm/tree/#coeff/size: %7.5f %zu, %6.3f m, %6.3f GByte",
1962 (
name.c_str()), wall,
norm, tsize,
double(ncoeff)*1.e-6,
double(ncoeff)/fac*
d);
1963 print(std::string(buf));
1968 template <
typename T, std::
size_t NDIM>
1970 if (this->targs.tt==
TT_FULL)
return;
1973 if (is_compressed())
k0=2*
k;
1978 if (world.rank()==0)
print(
"n.size(),k0,dim",n.
size(),
k0,dim);
1981 const nodeT& node = it->second;
1997 if (world.rank()==0) {
1998 print(
"configurations number of nodes");
1999 print(
" full rank ",n_full);
2000 for (
unsigned int i=0; i<n.
size(); i++) {
2001 print(
" ",i,
" ",n[i]);
2003 print(
" large rank ",n_large);
2008 for (
unsigned int i=0; i<std::min(3l,n.
size()); i++) nlog[0]+=n[i];
2009 for (
unsigned int i=3; i<std::min(10l,n.
size()); i++) nlog[1]+=n[i];
2010 for (
unsigned int i=10; i<std::min(30l,n.
size()); i++) nlog[2]+=n[i];
2011 for (
unsigned int i=30; i<std::min(100l,n.
size()); i++) nlog[3]+=n[i];
2012 for (
unsigned int i=100; i<std::min(300l,n.
size()); i++) nlog[4]+=n[i];
2013 for (
unsigned int i=300; i<std::min(1000l,n.
size()); i++) nlog[5]+=n[i];
2015 std::vector<std::string> slog={
"3",
"10",
"30",
"100",
"300",
"1000"};
2016 for (
unsigned int i=0; i<nlog.
size(); i++) {
2017 print(
" < ",slog[i],
" ",nlog[i]);
2019 print(
" large rank ",n_large);
2024 template <
typename T, std::
size_t NDIM>
2027 const int k = cdata.k;
2034 for (
int p=0;
p<
k; ++
p)
2037 else if (
NDIM == 2) {
2038 for (
int p=0;
p<
k; ++
p)
2039 for (
int q=0;
q<
k; ++
q)
2042 else if (
NDIM == 3) {
2043 for (
int p=0;
p<
k; ++
p)
2044 for (
int q=0;
q<
k; ++
q)
2045 for (
int r=0; r<
k; ++r)
2046 sum +=
c(
p,
q,r)*px[0][
p]*px[1][
q]*px[2][r];
2048 else if (
NDIM == 4) {
2049 for (
int p=0;
p<
k; ++
p)
2050 for (
int q=0;
q<
k; ++
q)
2051 for (
int r=0; r<
k; ++r)
2052 for (
int s=0; s<
k; ++s)
2053 sum +=
c(
p,
q,r,s)*px[0][
p]*px[1][
q]*px[2][r]*px[3][s];
2055 else if (
NDIM == 5) {
2056 for (
int p=0;
p<
k; ++
p)
2057 for (
int q=0;
q<
k; ++
q)
2058 for (
int r=0; r<
k; ++r)
2059 for (
int s=0; s<
k; ++s)
2060 for (
int t=0; t<
k; ++t)
2061 sum +=
c(
p,
q,r,s,t)*px[0][
p]*px[1][
q]*px[2][r]*px[3][s]*px[4][t];
2063 else if (
NDIM == 6) {
2064 for (
int p=0;
p<
k; ++
p)
2065 for (
int q=0;
q<
k; ++
q)
2066 for (
int r=0; r<
k; ++r)
2067 for (
int s=0; s<
k; ++s)
2068 for (
int t=0; t<
k; ++t)
2069 for (
int u=0;
u<
k; ++
u)
2070 sum +=
c(
p,
q,r,s,t,
u)*px[0][
p]*px[1][
q]*px[2][r]*px[3][s]*px[4][t]*px[5][
u];
2078 template <
typename T, std::
size_t NDIM>
2090 if (it == coeffs.end()) {
2092 it = coeffs.find(key).get();
2094 nodeT& node = it->second;
2105 if (!
d.has_data())
d =
coeffT(cdata.v2k,targs);
2106 if (accumulate_NS and (key.
level() > 0))
d(cdata.s0) += s;
2107 if (
d.dim(0)==2*get_k()) {
2112 const keyT& child = kit.key();
2116 woT::task(coeffs.owner(child), &implT::reconstruct_op, child, ss, accumulate_NS);
2126 if (s.has_no_data()) ss=
coeffT(cdata.vk,targs);
2132 template <
typename T, std::
size_t NDIM>
2135 std::vector<long> npt(
NDIM,qx.
dim(0));
2141 template <
typename T, std::
size_t NDIM>
2144 std::vector<long> npt(
NDIM,qx.
dim(0));
2150 template <
typename T, std::
size_t NDIM>
2159 const double h = std::pow(0.5,
double(n));
2161 const int npt = qx.
dim(0);
2168 for (std::size_t i = 0; i <
NDIM; i++) {
2169 c1[i] = cell(i,0) +
h*cell_width[i]*(l[i] + qx((
long)0));
2170 c2[i] = cell(i,0) +
h*cell_width[i]*(l[i] + qx(npt-1));
2172 if (
f.screened(c1, c2)) {
2178 bool vectorized =
f.supports_vectorized();
2180 T* fvptr = fval.
ptr();
2182 double* x1 =
new double[npt];
2184 for (
int i=0; i<npt; ++i, ++idx) {
2185 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2189 f(xvals, fvptr, npt);
2192 else if (
NDIM == 2) {
2193 double* x1 =
new double[npt*npt];
2194 double* x2 =
new double[npt*npt];
2196 for (
int i=0; i<npt; ++i) {
2197 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2198 for (
int j=0; j<npt; ++j, ++idx) {
2199 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2205 f(xvals, fvptr, npt*npt);
2209 else if (
NDIM == 3) {
2210 double* x1 =
new double[npt*npt*npt];
2211 double* x2 =
new double[npt*npt*npt];
2212 double* x3 =
new double[npt*npt*npt];
2214 for (
int i=0; i<npt; ++i) {
2215 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2216 for (
int j=0; j<npt; ++j) {
2217 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2218 for (
int k=0;
k<npt; ++
k, ++idx) {
2219 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2227 f(xvals, fvptr, npt*npt*npt);
2232 else if (
NDIM == 4) {
2233 double* x1 =
new double[npt*npt*npt*npt];
2234 double* x2 =
new double[npt*npt*npt*npt];
2235 double* x3 =
new double[npt*npt*npt*npt];
2236 double* x4 =
new double[npt*npt*npt*npt];
2238 for (
int i=0; i<npt; ++i) {
2239 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2240 for (
int j=0; j<npt; ++j) {
2241 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2242 for (
int k=0;
k<npt; ++
k) {
2243 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2244 for (
int m=0;
m<npt; ++
m, ++idx) {
2245 c[3] = cell(3,0) +
h*cell_width[3]*(l[3] + qx(
m));
2254 Vector<double*,4> xvals {x1, x2, x3, x4};
2255 f(xvals, fvptr, npt*npt*npt*npt);
2261 else if (
NDIM == 5) {
2262 double* x1 =
new double[npt*npt*npt*npt*npt];
2263 double* x2 =
new double[npt*npt*npt*npt*npt];
2264 double* x3 =
new double[npt*npt*npt*npt*npt];
2265 double* x4 =
new double[npt*npt*npt*npt*npt];
2266 double* x5 =
new double[npt*npt*npt*npt*npt];
2268 for (
int i=0; i<npt; ++i) {
2269 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2270 for (
int j=0; j<npt; ++j) {
2271 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2272 for (
int k=0;
k<npt; ++
k) {
2273 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2274 for (
int m=0;
m<npt; ++
m) {
2275 c[3] = cell(3,0) +
h*cell_width[3]*(l[3] + qx(
m));
2276 for (
int n=0; n<npt; ++n, ++idx) {
2277 c[4] = cell(4,0) +
h*cell_width[4]*(l[4] + qx(n));
2288 Vector<double*,5> xvals {x1, x2, x3, x4, x5};
2289 f(xvals, fvptr, npt*npt*npt*npt*npt);
2296 else if (
NDIM == 6) {
2297 double* x1 =
new double[npt*npt*npt*npt*npt*npt];
2298 double* x2 =
new double[npt*npt*npt*npt*npt*npt];
2299 double* x3 =
new double[npt*npt*npt*npt*npt*npt];
2300 double* x4 =
new double[npt*npt*npt*npt*npt*npt];
2301 double* x5 =
new double[npt*npt*npt*npt*npt*npt];
2302 double* x6 =
new double[npt*npt*npt*npt*npt*npt];
2304 for (
int i=0; i<npt; ++i) {
2305 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2306 for (
int j=0; j<npt; ++j) {
2307 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2308 for (
int k=0;
k<npt; ++
k) {
2309 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2310 for (
int m=0;
m<npt; ++
m) {
2311 c[3] = cell(3,0) +
h*cell_width[3]*(l[3] + qx(
m));
2312 for (
int n=0; n<npt; ++n) {
2313 c[4] = cell(4,0) +
h*cell_width[4]*(l[4] + qx(n));
2314 for (
int p=0;
p<npt; ++
p, ++idx) {
2315 c[5] = cell(5,0) +
h*cell_width[5]*(l[5] + qx(
p));
2328 Vector<double*,6> xvals {x1, x2, x3, x4, x5, x6};
2329 f(xvals, fvptr, npt*npt*npt*npt*npt*npt);
2347 for (
int i=0; i<npt; ++i) {
2348 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2353 else if (
NDIM == 2) {
2354 for (
int i=0; i<npt; ++i) {
2355 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2356 for (
int j=0; j<npt; ++j) {
2357 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2363 else if (
NDIM == 3) {
2364 for (
int i=0; i<npt; ++i) {
2365 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2366 for (
int j=0; j<npt; ++j) {
2367 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2368 for (
int k=0;
k<npt; ++
k) {
2369 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2376 else if (
NDIM == 4) {
2377 for (
int i=0; i<npt; ++i) {
2378 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2379 for (
int j=0; j<npt; ++j) {
2380 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2381 for (
int k=0;
k<npt; ++
k) {
2382 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2383 for (
int m=0;
m<npt; ++
m) {
2384 c[3] = cell(3,0) +
h*cell_width[3]*(l[3] + qx(
m));
2385 fval(i,j,
k,
m) =
f(
c);
2392 else if (
NDIM == 5) {
2393 for (
int i=0; i<npt; ++i) {
2394 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2395 for (
int j=0; j<npt; ++j) {
2396 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2397 for (
int k=0;
k<npt; ++
k) {
2398 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2399 for (
int m=0;
m<npt; ++
m) {
2400 c[3] = cell(3,0) +
h*cell_width[3]*(l[3] + qx(
m));
2401 for (
int n=0; n<npt; ++n) {
2402 c[4] = cell(4,0) +
h*cell_width[4]*(l[4] + qx(n));
2403 fval(i,j,
k,
m,n) =
f(
c);
2411 else if (
NDIM == 6) {
2412 for (
int i=0; i<npt; ++i) {
2413 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2414 for (
int j=0; j<npt; ++j) {
2415 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2416 for (
int k=0;
k<npt; ++
k) {
2417 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2418 for (
int m=0;
m<npt; ++
m) {
2419 c[3] = cell(3,0) +
h*cell_width[3]*(l[3] + qx(
m));
2420 for (
int n=0; n<npt; ++n) {
2421 c[4] = cell(4,0) +
h*cell_width[4]*(l[4] + qx(n));
2422 for (
int p=0;
p<npt; ++
p) {
2423 c[5] = cell(5,0) +
h*cell_width[5]*(l[5] + qx(
p));
2424 fval(i,j,
k,
m,n,
p) =
f(
c);
2439 template <
typename T, std::
size_t NDIM>
2445 template <
typename T, std::
size_t NDIM>
2457 template <
typename T, std::
size_t NDIM>
2462 if (do_refine && key.
level() < max_refine_level) {
2465 std::vector<Vector<double,NDIM> > newspecialpts;
2466 if (key.
level() < special_level && specialpts.size() > 0) {
2468 const auto bperiodic = bc.is_periodic();
2469 for (
unsigned int i = 0; i < specialpts.size(); ++i) {
2474 newspecialpts.push_back(specialpts[i]);
2488 const keyT& child = it.key();
2489 r(child_patch(child)) =
project(child);
2493 if (truncate_on_project)
s0 =
copy(
d(cdata.s0));
2500 if (newspecialpts.size() > 0 || dnorm >=truncate_tol(
thresh,key.
level())) {
2503 const keyT& child = it.key();
2506 p = world.random_proc();
2509 p = coeffs.owner(child);
2512 woT::task(
p, &implT::project_refine_op, child, do_refine, newspecialpts);
2516 if (truncate_on_project) {
2518 coeffs.replace(key,
nodeT(s,
false));
2523 const keyT& child = it.key();
2525 coeffs.replace(child,
nodeT(s,
false));
2535 template <
typename T, std::
size_t NDIM>
2537 std::vector<long> v0(
NDIM,0
L);
2538 std::vector<long> v1(
NDIM,1L);
2541 if (is_compressed()) {
2542 if (world.rank() == coeffs.owner(cdata.key0)) {
2545 nodeT& node = it->second;
2556 for (
typename dcT::iterator it=coeffs.begin(); it!=coeffs.end(); ++it) {
2557 Level n = it->first.level();
2558 nodeT& node = it->second;
2566 node.
coeff()(s) += tt;
2573 if (fence) world.gop.fence();
2576 template <
typename T, std::
size_t NDIM>
2579 if (is_compressed()) initial_level = std::max(initial_level,1);
2580 if (coeffs.is_local(key)) {
2581 if (is_compressed()) {
2582 if (key.
level() == initial_level) {
2586 coeffs.replace(key,
nodeT(
coeffT(cdata.v2k,targs),
true));
2590 if (key.
level()<initial_level) {
2594 coeffs.replace(key,
nodeT(
coeffT(cdata.vk,targs),
false));
2598 if (key.
level() < initial_level) {
2600 insert_zero_down_to_initial_level(kit.key());
2607 template <
typename T, std::
size_t NDIM>
2611 if (it == coeffs.end()) {
2615 coeffs.replace(key,
nodeT());
2616 it = coeffs.find(key).get();
2618 nodeT& node = it->second;
2620 std::vector< Future<bool> >
v = future_vector_factory<bool>(1<<
NDIM);
2625 return woT::task(world.rank(),&implT::truncate_op, key, tol,
v);
2634 if (dnorm < truncate_tol(tol,key)) {
2643 template <
typename T, std::
size_t NDIM>
2647 for (
int i=0; i<(1<<
NDIM); ++i)
if (
v[i].get())
return true;
2648 nodeT& node = coeffs.find(key).get()->second;
2655 if (key.
level() > 1) {
2657 if (dnorm < truncate_tol(tol,key)) {
2662 coeffs.erase(kit.key());
2671 template <
typename T, std::
size_t NDIM>
2673 if (world.rank() == 0) do_print_tree(cdata.key0, os, maxlevel);
2675 if (world.rank() == 0) os.flush();
2680 template <
typename T, std::
size_t NDIM>
2683 if (it == coeffs.end()) {
2685 for (
int i=0; i<key.
level(); ++i) os <<
" ";
2686 os << key <<
" missing --> " << coeffs.owner(key) <<
"\n";
2689 const nodeT& node = it->second;
2690 for (
int i=0; i<key.
level(); ++i) os <<
" ";
2691 os << key <<
" " << node <<
" --> " << coeffs.owner(key) <<
"\n";
2694 do_print_tree(kit.key(),os,maxlevel);
2700 template <
typename T, std::
size_t NDIM>
2702 std::multimap<Level, std::tuple<tranT, std::string>>
data;
2703 if (world.rank() == 0) do_print_tree_json(cdata.key0,
data, maxlevel);
2705 if (world.rank() == 0) {
2706 for (
Level level = 0; level != maxlevel; ++level) {
2707 if (
data.count(level) == 0)
2712 os <<
"\"" << level <<
"\":{";
2713 os <<
"\"level\": " << level <<
",";
2714 os <<
"\"nodes\":{";
2715 auto range =
data.equal_range(level);
2716 for (
auto it = range.first; it != range.second; ++it) {
2717 os <<
"\"" << std::get<0>(it->second) <<
"\":"
2718 << std::get<1>(it->second);
2719 if (std::next(it) != range.second)
2731 template <
typename T, std::
size_t NDIM>
2734 if (it == coeffs.end()) {
2738 const nodeT& node = it->second;
2739 std::ostringstream oss;
2742 oss <<
",\"owner\": " << coeffs.owner(key) <<
"}";
2743 auto node_json_str = oss.str();
2747 do_print_tree_json(kit.key(),
data, maxlevel);
2753 template <
typename T, std::
size_t NDIM>
2756 if (world.rank() == 0) do_print_tree_graphviz(cdata.key0, os, maxlevel);
2758 if (world.rank() == 0) os.flush();
2762 template <
typename T, std::
size_t NDIM>
2766 static int64_t value(
const keyT& key) {
2768 for (int64_t j = 0; j <= key.
level()-1; ++j) {
2769 result += (1 << j*
NDIM);
2777 if (it != coeffs.end()) {
2778 const nodeT& node = it->second;
2781 os << uniqhash::value(key) <<
" -> " << uniqhash::value(kit.key()) <<
"\n";
2782 do_print_tree_graphviz(kit.key(),os,maxlevel);
2788 template <
typename T, std::
size_t NDIM>
2792 if (not functor)
MADNESS_EXCEPTION(
"FunctionImpl: project: confusion about function?",0);
2795 if (functor->provides_coeff())
return functor->coeff(key).full_tensor_copy();
2800 tensorT workq(cdata.vq,
false);
2809 template <
typename T, std::
size_t NDIM>
2811 if (coeffs.probe(key)) {
2812 return Future<double>(coeffs.find(key).get()->second.get_norm_tree());
2816 return woT::task(coeffs.owner(parent), &implT::get_norm_tree_recursive, parent,
TaskAttributes::hipri());
2820 template <
typename T, std::
size_t NDIM>
2824 if (coeffs.probe(key)) {
2825 const nodeT& node = coeffs.find(key).get()->second;
2829 result.
set(std::pair<keyT,coeffT>(key,node.
coeff()));
2833 result.
set(std::pair<keyT,coeffT>(key,
coeffT()));
2840 if (coeffs.is_local(parent))
2848 template <
typename T, std::
size_t NDIM>
2852 if (coeffs.probe(key)) {
2853 const nodeT& node = coeffs.find(key).get()->second;
2856 result.
set(std::pair<keyT,coeffT>(key,node.
coeff()));
2870 template <
typename T, std::
size_t NDIM>
2892 nodeT& node = it->second;
2898 for (std::size_t i=0; i<
NDIM; ++i) {
2899 double xi = x[i]*2.0;
2901 if (li == 2) li = 1;
2913 template <
typename T, std::
size_t NDIM>
2920 while (key.
level() <= maxlevel) {
2921 if (coeffs.owner(key) ==
me) {
2924 if (it != coeffs.end()) {
2925 nodeT& node = it->second;
2931 for (std::size_t i=0; i<
NDIM; ++i) {
2932 double xi = x[i]*2.0;
2934 if (li == 2) li = 1;
2940 return std::pair<bool,T>(
false,0.0);
2943 template <
typename T, std::
size_t NDIM>
2965 nodeT& node = it->second;
2971 for (std::size_t i=0; i<
NDIM; ++i) {
2972 double xi = x[i]*2.0;
2974 if (li == 2) li = 1;
2985 template <
typename T, std::
size_t NDIM>
3007 nodeT& node = it->second;
3013 for (std::size_t i=0; i<
NDIM; ++i) {
3014 double xi = x[i]*2.0;
3016 if (li == 2) li = 1;
3028 template <
typename T, std::
size_t NDIM>
3039 template <
typename T, std::
size_t NDIM>
3042 coeffT shalf=t(cdata.sh);
3045 sfull(cdata.sh)-=shalf;
3049 template <
typename T, std::
size_t NDIM>
3055 if (t.
rank()==0)
return;
3057 for (
long i=0; i<t.
rank(); ++i) {
3060 tnorm(
c, &lo1, &hi1);
3068 template <
typename A,
typename B>
3075 template <
typename T, std::
size_t NDIM>
3092 template <
typename T, std::
size_t NDIM>
3100 template <
typename T, std::
size_t NDIM>
3106 template <
typename T, std::
size_t NDIM>
3114template <
typename T, std::
size_t NDIM>
3120 template <
typename T, std::
size_t NDIM>
3126 template <
typename T, std::
size_t NDIM>
3131 template <
typename T, std::
size_t NDIM>
3136 template <
typename T, std::
size_t NDIM>
3141 for (
int mu=0;
mu<cdata.npt; ++
mu) {
3142 double xmu =
scale*(cdata.quad_x(
mu)+lc) - lp;
3145 for (
int i=0; i<
k; ++i) phi(i,
mu) =
p[i];
3150 template <
typename T, std::
size_t NDIM>
3160 coeffT result = fcube_for_mul<T>(child, parent, s);
3162 result =
transform(result,cdata.quad_phiw);
3168 template <
typename T, std::
size_t NDIM>
3171 std::vector<long> v0(
NDIM,0);
3173 if (is_compressed()) {
3174 if (world.rank() == coeffs.owner(cdata.key0)) {
3176 if (it != coeffs.end()) {
3177 const nodeT& node = it->second;
3184 const keyT& key = it->first;
3185 const nodeT& node = it->second;
3202 }
else if (l >= two2n) {
3206 }
while (l >= two2n);
3215 return l >= 0 && l < two2n;
3218 template <
typename T, std::
size_t NDIM>
3227 return keyT::invalid();
3233 template <
typename T, std::
size_t NDIM>
3241 return keyT::invalid();
3247 template <
typename T, std::
size_t NDIM>
3251 typedef std::pair< Key<NDIM>,
coeffT > argT;
3261 template <
typename T, std::
size_t NDIM>
3263 bool nonstandard1,
bool keepleaves,
bool redundant1) {
3264 if (!coeffs.probe(key))
print(
"missing node",key);
3268 nodeT& node = coeffs.find(key).get()->second;
3272 std::vector< Future<std::pair<coeffT,double> > >
v = future_vector_factory<std::pair<coeffT,double> >(1<<
NDIM);
3277 v[i] = woT::task(coeffs.owner(kit.key()), &implT::compress_spawn, kit.key(),
3280 if (redundant1)
return woT::task(world.rank(),&implT::make_redundant_op, key,
v);
3281 return woT::task(world.rank(),&implT::compress_op, key,
v, nonstandard1);
3288 if (key.
level()==0) {
3300 coeffT sdcoeff(cdata.v2k,this->get_tensor_type());
3301 sdcoeff(cdata.s0)+=node.
coeff();
3302 node.
coeff()=sdcoeff;
3314 auto snorm=(keepleaves) ? node.
coeff().
normf() : 0.0;
3324 template <
typename T, std::
size_t NDIM>
3327 const coordT& plotlo,
const coordT& plothi,
const std::vector<long>& npt,
3328 bool eval_refine)
const {
3333 for (std::size_t i=0; i<
NDIM; ++i) {
3335 h[i] = (plothi[i]-plotlo[i])/(npt[i]-1);
3344 const Vector<Translation,NDIM>& l = key.
translation();
3345 const double twon =
pow(2.0,
double(n));
3346 const tensorT& coeff = coeffs.find(key).get()->second.coeff().full_tensor();
3354 for (std::size_t
d=0;
d<
NDIM; ++
d) {
3357 boxhi[
d] = boxlo[
d]+fac;
3359 if (boxlo[
d] > plothi[
d] || boxhi[
d] < plotlo[
d]) {
3361 npttotal = boxnpt[
d] = 0;
3365 else if (npt[
d] == 1) {
3367 boxlo[
d] = boxhi[
d] = plotlo[
d];
3372 boxlo[
d] = std::max(boxlo[
d],plotlo[
d]);
3373 boxhi[
d] = std::min(boxhi[
d],plothi[
d]);
3376 double xlo = long((boxlo[
d]-plotlo[
d])/
h[
d])*
h[
d] + plotlo[
d];
3377 if (xlo < boxlo[
d]) xlo +=
h[
d];
3379 double xhi = long((boxhi[
d]-plotlo[
d])/
h[
d])*
h[
d] + plotlo[
d];
3380 if (xhi > boxhi[
d]) xhi -=
h[
d];
3383 boxnpt[
d] = long(round((boxhi[
d] - boxlo[
d])/
h[
d])) + 1;
3385 npttotal *= boxnpt[
d];
3390 for (std::size_t
d=0;
d<
NDIM; ++
d) {
3391 double xd = boxlo[
d] + it[
d]*
h[
d];
3392 x[
d] = twon*xd - l[
d];
3395 ind[
d] = long(round((xd-plotlo[
d])/
h[
d]));
3406 T tmp = eval_cube(n, x, coeff);
3416 template <
typename T, std::
size_t NDIM>
3419 const std::vector<long>& npt,
3420 const bool eval_refine)
const {
3428 const nodeT& node = it->second;
3429 if (node.has_coeff()) {
3430 woT::task(world.rank(), &implT::plot_cube_kernel,
3437 world.taskq.fence();
3445 fprintf(
f,
"%.6e\n",t);
3449 fprintf(
f,
"%.6e %.6e\n", t.real(), t.imag());
3452 template <
typename T, std::
size_t NDIM>
3456 const std::vector<long>& npt,
3460 const char* element[6] = {
"lines",
"quads",
"cubes",
"cubes4D",
"cubes5D",
"cubes6D"};
3465 if (world.
rank() == 0) {
3469 fprintf(
f,
"object 1 class gridpositions counts ");
3470 for (std::size_t
d=0;
d<
NDIM; ++
d) fprintf(
f,
" %ld",npt[
d]);
3473 fprintf(
f,
"origin ");
3474 for (std::size_t
d=0;
d<
NDIM; ++
d) fprintf(
f,
" %.6e", cell(
d,0));
3478 fprintf(
f,
"delta ");
3479 for (std::size_t
c=0;
c<
d; ++
c) fprintf(
f,
" 0");
3481 if (npt[
d]>1)
h = (cell(
d,1)-cell(
d,0))/(npt[
d]-1);
3483 for (std::size_t
c=
d+1;
c<
NDIM; ++
c) fprintf(
f,
" 0");
3488 fprintf(
f,
"object 2 class gridconnections counts ");
3489 for (std::size_t
d=0;
d<
NDIM; ++
d) fprintf(
f,
" %ld",npt[
d]);
3491 fprintf(
f,
"attribute \"element type\" string \"%s\"\n", element[
NDIM-1]);
3492 fprintf(
f,
"attribute \"ref\" string \"positions\"\n");
3496 for (std::size_t
d=0;
d<
NDIM; ++
d) npoint *= npt[
d];
3497 const char* iscomplex =
"";
3499 const char* isbinary =
"";
3500 if (binary) isbinary =
"binary";
3501 fprintf(
f,
"object 3 class array type double %s rank 0 items %d %s data follows\n",
3502 iscomplex, npoint, isbinary);
3506 Tensor<T> r =
function.eval_cube(cell, npt);
3508 if (world.
rank() == 0) {
3512 fwrite((
void *) r.ptr(),
sizeof(
T), r.size(),
f);
3516 for (IndexIterator it(npt); it; ++it) {
3523 fprintf(
f,
"object \"%s\" class field\n",
filename);
3524 fprintf(
f,
"component \"positions\" value 1\n");
3525 fprintf(
f,
"component \"connections\" value 2\n");
3526 fprintf(
f,
"component \"data\" value 3\n");
3527 fprintf(
f,
"\nend\n");
3533 template <std::
size_t NDIM>
3539 max_refine_level = 30;
3544 truncate_on_project =
true;
3545 apply_randomize =
false;
3546 project_randomize =
false;
3549 cell = make_default_cell();
3550 recompute_cell_info();
3551 set_default_pmap(world);
3554 template <std::
size_t NDIM>
3557 return std::make_shared<LevelPmap< Key<NDIM> >>(world);
3561 template <std::
size_t NDIM>
3563 pmap = make_default_pmap(world);
3564 pmap_nproc = world.
nproc();
3568 template <std::
size_t NDIM>
3570 std::cout <<
"Function Defaults:" << std::endl;
3571 std::cout <<
" Dimension " <<
": " <<
NDIM << std::endl;
3572 std::cout <<
" k" <<
": " <<
k << std::endl;
3573 std::cout <<
" thresh" <<
": " <<
thresh << std::endl;
3574 std::cout <<
" initial_level" <<
": " << initial_level << std::endl;
3575 std::cout <<
" special_level" <<
": " << special_level << std::endl;
3576 std::cout <<
" max_refine_level" <<
": " << max_refine_level << std::endl;
3577 std::cout <<
" truncate_mode" <<
": " <<
truncate_mode << std::endl;
3578 std::cout <<
" refine" <<
": " <<
refine << std::endl;
3579 std::cout <<
" autorefine" <<
": " << autorefine << std::endl;
3580 std::cout <<
" debug" <<
": " <<
debug << std::endl;
3581 std::cout <<
" truncate_on_project" <<
": " << truncate_on_project << std::endl;
3582 std::cout <<
" apply_randomize" <<
": " << apply_randomize << std::endl;
3583 std::cout <<
" project_randomize" <<
": " << project_randomize << std::endl;
3584 std::cout <<
" bc" <<
": " << get_bc() << std::endl;
3585 std::cout <<
" tt" <<
": " << tt << std::endl;
3586 std::cout <<
" cell" <<
": " << cell << std::endl;
3589 template <
typename T, std::
size_t NDIM>
3590 const FunctionCommonData<T,NDIM>*
FunctionCommonData<T,NDIM>::data[
MAXK] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
double w(double t, double eps)
Definition DKops.h:22
double q(double t)
Definition DKops.h:18
std::complex< double > double_complex
Definition cfft.h:14
Definition test_ar.cc:118
Definition test_ar.cc:141
Definition test_tree.cc:78
long dim(int i) const
Returns the size of dimension i.
Definition basetensor.h:147
long size() const
Returns the number of elements in the tensor.
Definition basetensor.h:138
This class is used to specify boundary conditions for all operators.
Definition bc.h:72
a class to track where relevant (parent) coeffs are
Definition funcimpl.h:791
Tri-diagonal operator traversing tree primarily for derivative operator.
Definition derivative.h:73
static std::vector< Key< NDIM > > disp
standard displacements to be used with standard kernels (range-unrestricted, no lattice sum)
Definition displacements.h:53
static array_of_bools< NDIM > periodic_axes
along which axes lattice summation is performed?
Definition displacements.h:54
static std::vector< Key< NDIM > > disp_periodic[64]
displacements to be used with lattice-summed kernels
Definition displacements.h:55
ElementaryInterface (formerly FunctorInterfaceWrapper) interfaces a c-function.
Definition function_interface.h:273
FunctionCommonData holds all Function data common for given k.
Definition function_common_data.h:52
static const FunctionCommonData< T, NDIM > & get(int k)
Definition function_common_data.h:111
static void _init_quadrature(int k, int npt, Tensor< double > &quad_x, Tensor< double > &quad_w, Tensor< double > &quad_phi, Tensor< double > &quad_phiw, Tensor< double > &quad_phit)
Initialize the quadrature information.
Definition mraimpl.h:91
void _init_twoscale()
Private. Initialize the twoscale coefficients.
Definition mraimpl.h:70
FunctionDefaults holds default paramaters as static class members.
Definition funcdefaults.h:100
static Tensor< double > make_default_cell_width()
Definition funcdefaults.h:136
static std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > make_default_pmap(World &world)
Makes a default process map for the given world.
Definition mraimpl.h:3555
static bool truncate_on_project
If true initial projection inserts at n-1 not n.
Definition funcdefaults.h:114
static bool apply_randomize
If true use randomization for load balancing in apply integral operator.
Definition funcdefaults.h:115
static double cell_volume
Volume of simulation cell.
Definition funcdefaults.h:121
static int k
Wavelet order.
Definition funcdefaults.h:105
static int truncate_mode
Truncation method.
Definition funcdefaults.h:110
static void set_default_pmap(World &world)
Definition mraimpl.h:3562
static int pmap_nproc
Number of processes assumed by pmap, -1 indicates uninitialized pmap.
Definition funcdefaults.h:125
static Tensor< double > cell_width
Width of simulation cell in each dimension.
Definition funcdefaults.h:119
static double thresh
Truncation threshold.
Definition funcdefaults.h:106
static bool debug
Controls output of debug info.
Definition funcdefaults.h:113
static int special_level
Minimum level for fine scale projection of special boxes.
Definition funcdefaults.h:108
static double cell_min_width
Size of smallest dimension.
Definition funcdefaults.h:122
static std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > pmap
Default mapping of keys to processes.
Definition funcdefaults.h:124
static const Tensor< double > & get_cell_width()
Returns the width of each user cell dimension.
Definition funcdefaults.h:370
static const BoundaryConditions< NDIM > & get_bc()
Returns the default boundary conditions.
Definition funcdefaults.h:311
static Tensor< double > cell
cell[NDIM][2] Simulation cell, cell(0,0)=xlo, cell(0,1)=xhi, ...
Definition funcdefaults.h:118
static int initial_level
Initial level for fine scale projection.
Definition funcdefaults.h:107
static bool autorefine
Whether to autorefine in multiplication, etc.
Definition funcdefaults.h:112
static Tensor< double > make_default_cell()
Definition funcdefaults.h:127
static std::optional< BoundaryConditions< NDIM > > bc
Default boundary conditions, not initialized by default and must be set explicitly before use.
Definition funcdefaults.h:117
static void print()
Definition mraimpl.h:3569
static Tensor< double > rcell_width
Reciprocal of width.
Definition funcdefaults.h:120
static double get_cell_min_width()
Returns the minimum width of any user cell dimension.
Definition funcdefaults.h:380
static bool project_randomize
If true use randomization for load balancing in project/refine.
Definition funcdefaults.h:116
static const Tensor< double > & get_cell()
Gets the user cell for the simulation.
Definition funcdefaults.h:348
static void set_defaults(World &world)
Definition mraimpl.h:3534
static int max_refine_level
Level at which to stop refinement.
Definition funcdefaults.h:109
static TensorType tt
structure of the tensor in FunctionNode
Definition funcdefaults.h:123
static bool refine
Whether to refine new functions.
Definition funcdefaults.h:111
Abstract base class interface required for functors used as input to Functions.
Definition function_interface.h:68
Definition funcimpl.h:5536
FunctionImpl holds all Function state to facilitate shallow copy semantics.
Definition funcimpl.h:945
bool is_nonstandard() const
Definition mraimpl.h:273
T eval_cube(Level n, coordT &x, const tensorT &c) const
Definition mraimpl.h:2025
void do_print_tree_graphviz(const keyT &key, std::ostream &os, Level maxlevel) const
Functor for the do_print_tree method (using GraphViz)
Definition mraimpl.h:2763
void change_tensor_type1(const TensorArgs &targs, bool fence)
change the tensor type of the coefficients in the FunctionNode
Definition mraimpl.h:1099
void evaldepthpt(const Vector< double, NDIM > &xin, const keyT &keyin, const typename Future< Level >::remote_refT &ref)
Get the depth of the tree at a point in simulation coordinates.
Definition mraimpl.h:2944
void scale_inplace(const T q, bool fence)
In-place scale by a constant.
Definition mraimpl.h:3115
void gaxpy_oop_reconstructed(const double alpha, const implT &f, const double beta, const implT &g, const bool fence)
perform: this= alpha*f + beta*g, invoked by result
Definition mraimpl.h:223
bool is_redundant() const
Returns true if the function is redundant.
Definition mraimpl.h:262
FunctionNode< Q, NDIM > nodeT
Type of node.
Definition funcimpl.h:955
std::size_t nCoeff_local() const
Returns the number of coefficients in the function for this MPI rank.
Definition mraimpl.h:1922
void print_size(const std::string name) const
print tree size and size
Definition mraimpl.h:1941
void print_info() const
Prints summary of data distribution.
Definition mraimpl.h:833
void abs_inplace(bool fence)
Definition mraimpl.h:3127
void print_timer() const
Definition mraimpl.h:357
void evalR(const Vector< double, NDIM > &xin, const keyT &keyin, const typename Future< long >::remote_refT &ref)
Get the rank of leaf box of the tree at a point in simulation coordinates.
Definition mraimpl.h:2986
const FunctionCommonData< T, NDIM > & cdata
Definition funcimpl.h:983
void do_print_grid(const std::string filename, const std::vector< keyT > &keys) const
print the grid in xyz format
Definition mraimpl.h:584
std::size_t nCoeff() const
Returns the number of coefficients in the function ... collective global sum.
Definition mraimpl.h:1932
keyT neighbor_in_volume(const keyT &key, const keyT &disp) const
Returns key of general neighbor that resides in-volume.
Definition mraimpl.h:3234
void compress(const TreeState newstate, bool fence)
compress the wave function
Definition mraimpl.h:1500
std::pair< coeffT, double > compress_op(const keyT &key, const std::vector< Future< std::pair< coeffT, double > > > &v, bool nonstandard)
calculate the wavelet coefficients using the sum coefficients of all child nodes
Definition mraimpl.h:1668
Future< bool > truncate_spawn(const keyT &key, double tol)
Returns true if after truncation this node has coefficients.
Definition mraimpl.h:2608
Future< double > norm_tree_spawn(const keyT &key)
Definition mraimpl.h:1570
std::vector< keyT > local_leaf_keys() const
return the keys of the local leaf boxes
Definition mraimpl.h:558
void do_print_tree(const keyT &key, std::ostream &os, Level maxlevel) const
Functor for the do_print_tree method.
Definition mraimpl.h:2681
void unset_functor()
Definition mraimpl.h:312
void do_print_plane(const std::string filename, std::vector< Tensor< double > > plotinfo, const int xaxis, const int yaxis, const coordT el2)
print the MRA structure
Definition mraimpl.h:499
std::pair< Key< NDIM >, ShallowNode< T, NDIM > > find_datum(keyT key) const
return the a std::pair<key, node>, which MUST exist
Definition mraimpl.h:965
void set_functor(const std::shared_ptr< FunctionFunctorInterface< T, NDIM > > functor1)
Definition mraimpl.h:293
bool verify_tree_state_local() const
check that the tree state and the coeffs are consistent
Definition mraimpl.h:169
const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > & get_pmap() const
Definition mraimpl.h:207
void sock_it_to_me(const keyT &key, const RemoteReference< FutureImpl< std::pair< keyT, coeffT > > > &ref) const
Walk up the tree returning pair(key,node) for first node with coefficients.
Definition mraimpl.h:2821
double get_thresh() const
Definition mraimpl.h:328
void trickle_down(bool fence)
sum all the contributions from all scales after applying an operator in mod-NS form
Definition mraimpl.h:1354
std::pair< coeffT, double > make_redundant_op(const keyT &key, const std::vector< Future< std::pair< coeffT, double > > > &v)
similar to compress_op, but insert only the sum coefficients in the tree
Definition mraimpl.h:1728
void set_autorefine(bool value)
Definition mraimpl.h:337
tensorT filter(const tensorT &s) const
Transform sum coefficients at level n to sums+differences at level n-1.
Definition mraimpl.h:1152
void chop_at_level(const int n, const bool fence=true)
remove all nodes with level higher than n
Definition mraimpl.h:1115
void print_tree_json(std::ostream &os=std::cout, Level maxlevel=10000) const
Definition mraimpl.h:2701
coeffT parent_to_child_NS(const keyT &child, const keyT &parent, const coeffT &coeff) const
Directly project parent NS coeffs to child NS coeffs.
Definition mraimpl.h:707
void mapdim(const implT &f, const std::vector< long > &map, bool fence)
Permute the dimensions of f according to map, result on this.
Definition mraimpl.h:1057
bool is_compressed() const
Returns true if the function is compressed.
Definition mraimpl.h:250
void mirror(const implT &f, const std::vector< long > &mirror, bool fence)
mirror the dimensions of f according to map, result on this
Definition mraimpl.h:1066
void print_stats() const
print the number of configurations per node
Definition mraimpl.h:1969
void broaden(const array_of_bools< NDIM > &is_periodic, bool fence)
Definition mraimpl.h:1303
coeffT truncate_reconstructed_op(const keyT &key, const std::vector< Future< coeffT > > &v, const double tol)
given the sum coefficients of all children, truncate or not
Definition mraimpl.h:1617
void fcube(const keyT &key, const FunctionFunctorInterface< T, NDIM > &f, const Tensor< double > &qx, tensorT &fval) const
Evaluate function at quadrature points in the specified box.
Definition mraimpl.h:2446
void forward_do_diff1(const DerivativeBase< T, NDIM > *D, const implT *f, const keyT &key, const std::pair< keyT, coeffT > &left, const std::pair< keyT, coeffT > ¢er, const std::pair< keyT, coeffT > &right)
Definition mraimpl.h:923
std::vector< Slice > child_patch(const keyT &child) const
Returns patch referring to coeffs of child in parent box.
Definition mraimpl.h:696
void print_tree_graphviz(std::ostream &os=std::cout, Level maxlevel=10000) const
Definition mraimpl.h:2754
std::size_t min_nodes() const
Returns the min number of nodes on a processor.
Definition mraimpl.h:1873
void make_redundant(const bool fence)
convert this to redundant, i.e. have sum coefficients on all levels
Definition mraimpl.h:1528
std::size_t max_nodes() const
Returns the max number of nodes on a processor.
Definition mraimpl.h:1864
coeffT upsample(const keyT &key, const coeffT &coeff) const
upsample the sum coefficients of level 1 to sum coeffs on level n+1
Definition mraimpl.h:1231
void flo_unary_op_node_inplace(const opT &op, bool fence)
Definition funcimpl.h:2229
std::size_t size_local() const
Returns the number of coefficients in the function for each rank.
Definition mraimpl.h:1891
void plot_cube_kernel(archive::archive_ptr< Tensor< T > > ptr, const keyT &key, const coordT &plotlo, const coordT &plothi, const std::vector< long > &npt, bool eval_refine) const
Definition mraimpl.h:3325
T trace_local() const
Returns int(f(x),x) in local volume.
Definition mraimpl.h:3169
void print_grid(const std::string filename) const
Definition mraimpl.h:542
Future< std::pair< coeffT, double > > compress_spawn(const keyT &key, bool nonstandard, bool keepleaves, bool redundant1)
Invoked on node where key is local.
Definition mraimpl.h:3262
bool get_autorefine() const
Definition mraimpl.h:334
void phi_for_mul(Level np, Translation lp, Level nc, Translation lc, Tensor< double > &phi) const
Compute the Legendre scaling functions for multiplication.
Definition mraimpl.h:3137
Future< std::pair< keyT, coeffT > > find_me(const keyT &key) const
find_me. Called by diff_bdry to get coefficients of boundary function
Definition mraimpl.h:3249
TensorType get_tensor_type() const
Definition mraimpl.h:319
void remove_leaf_coefficients(const bool fence)
Definition mraimpl.h:1522
void insert_zero_down_to_initial_level(const keyT &key)
Initialize nodes to zero function at initial_level of refinement.
Definition mraimpl.h:2577
void do_diff1(const DerivativeBase< T, NDIM > *D, const implT *f, const keyT &key, const std::pair< keyT, coeffT > &left, const std::pair< keyT, coeffT > ¢er, const std::pair< keyT, coeffT > &right)
Definition mraimpl.h:934
void do_print_tree_json(const keyT &key, std::multimap< Level, std::tuple< tranT, std::string > > &data, Level maxlevel) const
Functor for the do_print_tree_json method.
Definition mraimpl.h:2732
void finalize_sum()
after summing up we need to do some cleanup;
Definition mraimpl.h:1821
dcT coeffs
The coefficients.
Definition funcimpl.h:988
bool exists_and_is_leaf(const keyT &key) const
Definition mraimpl.h:1275
void verify_tree() const
Verify tree is properly constructed ... global synchronization involved.
Definition mraimpl.h:111
void set_tensor_args(const TensorArgs &t)
Definition mraimpl.h:325
Range< typename dcT::const_iterator > rangeT
Definition funcimpl.h:5665
std::size_t real_size() const
Returns the number of coefficients in the function ... collective global sum.
Definition mraimpl.h:1909
bool exists_and_has_children(const keyT &key) const
Definition mraimpl.h:1270
void sum_down_spawn(const keyT &key, const coeffT &s)
is this the same as trickle_down() ?
Definition mraimpl.h:876
keyT neighbor(const keyT &key, const keyT &disp, const array_of_bools< NDIM > &is_periodic) const
Returns key of general neighbor enforcing BC.
Definition mraimpl.h:3219
void norm_tree(bool fence)
compute for each FunctionNode the norm of the function inside that node
Definition mraimpl.h:1547
bool has_leaves() const
Definition mraimpl.h:288
bool verify_parents_and_children() const
check that parents and children are consistent
Definition mraimpl.h:119
void reconstruct_op(const keyT &key, const coeffT &s, const bool accumulate_NS=true)
Definition mraimpl.h:2079
const coeffT parent_to_child(const coeffT &s, const keyT &parent, const keyT &child) const
Directly project parent coeffs to child coeffs.
Definition mraimpl.h:3152
void undo_redundant(const bool fence)
convert this from redundant to standard reconstructed form
Definition mraimpl.h:1538
GenTensor< Q > coeffT
Type of tensor used to hold coeffs.
Definition funcimpl.h:956
const keyT & key0() const
Returns cdata.key0.
Definition mraimpl.h:394
double finalize_apply()
after apply we need to do some cleanup;
Definition mraimpl.h:1778
const dcT & get_coeffs() const
Definition mraimpl.h:343
double norm2sq_local() const
Returns the square of the local norm ... no comms.
Definition mraimpl.h:1830
const FunctionCommonData< T, NDIM > & get_cdata() const
Definition mraimpl.h:349
void sum_down(bool fence)
After 1d push operator must sum coeffs down the tree to restore correct scaling function coefficients...
Definition mraimpl.h:915
bool noautorefine(const keyT &key, const tensorT &t) const
Always returns false (for when autorefine is not wanted)
Definition mraimpl.h:859
double truncate_tol(double tol, const keyT &key) const
Returns the truncation threshold according to truncate_method.
Definition mraimpl.h:649
bool autorefine_square_test(const keyT &key, const nodeT &t) const
Returns true if this block of coeffs needs autorefining.
Definition mraimpl.h:865
void erase(const Level &max_level)
truncate tree at a certain level
Definition mraimpl.h:739
void reconstruct(bool fence)
reconstruct this tree – respects fence
Definition mraimpl.h:1468
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:3659
coeffT assemble_coefficients(const keyT &key, const coeffT &coeff_ket, const coeffT &vpotential1, const coeffT &vpotential2, const tensorT &veri) const
given several coefficient tensors, assemble a result tensor
Definition mraimpl.h:1013
static void tnorm(const tensorT &t, double *lo, double *hi)
Computes norm of low/high-order polyn. coeffs for autorefinement test.
Definition mraimpl.h:3029
std::pair< bool, T > eval_local_only(const Vector< double, NDIM > &xin, Level maxlevel)
Evaluate function only if point is local returning (true,value); otherwise return (false,...
Definition mraimpl.h:2915
std::size_t max_depth() const
Returns the maximum depth of the tree ... collective ... global sum/broadcast.
Definition mraimpl.h:1856
std::size_t size() const
Returns the number of coefficients in the function ... collective global sum.
Definition mraimpl.h:1901
void reduce_rank(const double thresh, bool fence)
reduce the rank of the coefficients tensors
Definition mraimpl.h:1107
std::shared_ptr< FunctionFunctorInterface< T, NDIM > > get_functor()
Definition mraimpl.h:300
tensorT unfilter(const tensorT &s) const
Transform sums+differences at level n to sum coefficients at level n+1.
Definition mraimpl.h:1181
Tensor< T > eval_plot_cube(const coordT &plotlo, const coordT &plothi, const std::vector< long > &npt, const bool eval_refine=false) const
Definition mraimpl.h:3417
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:1407
Future< coeffT > truncate_reconstructed_spawn(const keyT &key, const double tol)
truncate using a tree in reconstructed form
Definition mraimpl.h:1593
void map_and_mirror(const implT &f, const std::vector< long > &map, const std::vector< long > &mirror, bool fence)
map and mirror the translation index and the coefficients, result on this
Definition mraimpl.h:1076
void truncate(double tol, bool fence)
Truncate according to the threshold with optional global fence.
Definition mraimpl.h:378
bool is_reconstructed() const
Returns true if the function is compressed.
Definition mraimpl.h:256
double norm_tree_op(const keyT &key, const std::vector< Future< double > > &v)
Definition mraimpl.h:1555
void reset_timer()
Definition mraimpl.h:366
void refine_to_common_level(const std::vector< FunctionImpl< T, NDIM > * > &v, const std::vector< tensorT > &c, const keyT key)
Refine multiple functions down to the same finest level.
Definition mraimpl.h:769
int get_k() const
Definition mraimpl.h:340
void eval(const Vector< double, NDIM > &xin, const keyT &keyin, const typename Future< T >::remote_refT &ref)
Evaluate the function at a point in simulation coordinates.
Definition mraimpl.h:2871
bool truncate_op(const keyT &key, double tol, const std::vector< Future< bool > > &v)
Definition mraimpl.h:2644
void zero_norm_tree()
Definition mraimpl.h:1292
std::size_t max_local_depth() const
Returns the maximum local depth of the tree ... no communications.
Definition mraimpl.h:1842
tensorT project(const keyT &key) const
Definition mraimpl.h:2789
double check_symmetry_local() const
Returns some asymmetry measure ... no comms.
Definition mraimpl.h:755
Future< double > get_norm_tree_recursive(const keyT &key) const
Definition mraimpl.h:2810
bool is_redundant_after_merge() const
Returns true if the function is redundant_after_merge.
Definition mraimpl.h:268
Key< NDIM > keyT
Type of key.
Definition funcimpl.h:954
bool is_on_demand() const
Definition mraimpl.h:283
void accumulate_timer(const double time) const
Definition mraimpl.h:352
void trickle_down_op(const keyT &key, const coeffT &s)
sum all the contributions from all scales after applying an operator in mod-NS form
Definition mraimpl.h:1365
void set_thresh(double value)
Definition mraimpl.h:331
Tensor< double > print_plane_local(const int xaxis, const int yaxis, const coordT &el2)
collect the data for a plot of the MRA structure locally on each node
Definition mraimpl.h:423
void sock_it_to_me_too(const keyT &key, const RemoteReference< FutureImpl< std::pair< keyT, coeffT > > > &ref) const
Definition mraimpl.h:2849
void broaden_op(const keyT &key, const std::vector< Future< bool > > &v)
Definition mraimpl.h:1281
void print_plane(const std::string filename, const int xaxis, const int yaxis, const coordT &el2)
Print a plane ("xy", "xz", or "yz") containing the point x to file.
Definition mraimpl.h:403
void print_tree(std::ostream &os=std::cout, Level maxlevel=10000) const
Definition mraimpl.h:2672
void project_refine_op(const keyT &key, bool do_refine, const std::vector< Vector< double, NDIM > > &specialpts)
Definition mraimpl.h:2458
std::size_t tree_size() const
Returns the size of the tree structure of the function ... collective global sum.
Definition mraimpl.h:1882
void add_scalar_inplace(T t, bool fence)
Adds a constant to the function. Local operation, optional fence.
Definition mraimpl.h:2536
tensorT downsample(const keyT &key, const std::vector< Future< coeffT > > &v) const
downsample the sum coefficients of level n+1 to sum coeffs on level n
Definition mraimpl.h:1201
void abs_square_inplace(bool fence)
Definition mraimpl.h:3132
void put_in_box(ProcessID from, long nl, long ni) const
Definition mraimpl.h:824
TensorArgs get_tensor_args() const
Definition mraimpl.h:322
void average(const implT &rhs)
take the average of two functions, similar to: this=0.5*(this+rhs)
Definition mraimpl.h:1088
void diff(const DerivativeBase< T, NDIM > *D, const implT *f, bool fence)
Definition mraimpl.h:946
void square_inplace(bool fence)
Pointwise squaring of function with optional global fence.
Definition mraimpl.h:3121
void remove_internal_coefficients(const bool fence)
Definition mraimpl.h:1517
void compute_snorm_and_dnorm(bool fence=true)
compute norm of s and d coefficients for all nodes
Definition mraimpl.h:1131
void standard(bool fence)
Changes non-standard compressed form to standard compressed form.
Definition mraimpl.h:1765
bool is_nonstandard_with_leaves() const
Definition mraimpl.h:278
FunctionNode holds the coefficients, etc., at each node of the 2^NDIM-tree.
Definition funcimpl.h:127
bool has_coeff() const
Returns true if there are coefficients in this node.
Definition funcimpl.h:200
void clear_coeff()
Clears the coefficients (has_coeff() will subsequently return false)
Definition funcimpl.h:295
bool is_leaf() const
Returns true if this does not have children.
Definition funcimpl.h:213
void set_has_children(bool flag)
Sets has_children attribute to value of flag.
Definition funcimpl.h:254
double get_norm_tree() const
Gets the value of norm_tree.
Definition funcimpl.h:311
void set_snorm(const double sn)
set the precomputed norm of the (virtual) s coefficients
Definition funcimpl.h:321
void print_json(std::ostream &s) const
Definition funcimpl.h:466
bool has_children() const
Returns true if this node has children.
Definition funcimpl.h:207
void set_coeff(const coeffT &coeffs)
Takes a shallow copy of the coeff — same as this->coeff()=coeff.
Definition funcimpl.h:285
void set_dnorm(const double dn)
set the precomputed norm of the (virtual) d coefficients
Definition funcimpl.h:326
coeffT & coeff()
Returns a non-const reference to the tensor containing the coeffs.
Definition funcimpl.h:227
void set_norm_tree(double norm_tree)
Sets the value of norm_tree.
Definition funcimpl.h:306
A multiresolution adaptive numerical function.
Definition mra.h:139
Implements the functionality of futures.
Definition future.h:74
A future is a possibly yet unevaluated value.
Definition future.h:369
T & get(bool dowork=true) &
Gets the value, waiting if necessary.
Definition future.h:570
remote_refT remote_ref(World &world) const
Returns a structure used to pass references to another process.
Definition future.h:671
void set(const Future< T > &other)
A.set(B), where A and B are futures ensures A has/will have the same value as B.
Definition future.h:504
RemoteReference< FutureImpl< T > > remote_refT
Definition future.h:394
Definition lowranktensor.h:59
GenTensor convert(const TensorArgs &targs) const
Definition gentensor.h:198
long dim(const int i) const
return the number of entries in dimension i
Definition lowranktensor.h:391
Tensor< T > full_tensor_copy() const
Definition gentensor.h:206
long ndim() const
Definition lowranktensor.h:386
constexpr bool is_full_tensor() const
Definition gentensor.h:224
size_t real_size() const
Definition gentensor.h:214
GenTensor< T > & emul(const GenTensor< T > &other)
Inplace multiply by corresponding elements of argument Tensor.
Definition lowranktensor.h:631
float_scalar_type normf() const
Definition lowranktensor.h:406
void reduce_rank(const double &eps)
Definition gentensor.h:217
long rank() const
Definition gentensor.h:212
const Tensor< T > & full_tensor() const
Definition gentensor.h:200
TensorType tensor_type() const
Definition gentensor.h:221
bool has_data() const
Definition gentensor.h:210
const long * dims() const
return the number of entries in dimension i
Definition lowranktensor.h:397
GenTensor & gaxpy(const T alpha, const GenTensor &other, const T beta)
Definition lowranktensor.h:580
IsSupported< TensorTypeData< Q >, GenTensor< T > & >::type scale(Q fac)
Inplace multiplication by scalar of supported type (legacy name)
Definition lowranktensor.h:426
constexpr bool is_svd_tensor() const
Definition gentensor.h:222
Iterates in lexical order thru all children of a key.
Definition key.h:467
Key is the index for a node of the 2^NDIM-tree.
Definition key.h:69
Level level() const
Definition key.h:168
bool is_neighbor_of(const Key &key, const array_of_bools< NDIM > &bperiodic) const
Assuming keys are at the same level, returns true if displaced by no more than 1 in any direction.
Definition key.h:287
bool is_valid() const
Checks if a key is valid.
Definition key.h:123
bool thisKeyContains(const Vector< double, NDIM > &x, const unsigned int &dim0, const unsigned int &dim1) const
check if this MultiIndex contains point x, disregarding these two dimensions
Definition key.h:319
bool is_invalid() const
Checks if a key is invalid.
Definition key.h:118
Key parent(int generation=1) const
Returns the key of the parent.
Definition key.h:252
const Vector< Translation, NDIM > & translation() const
Definition key.h:173
Range, vaguely a la Intel TBB, to encapsulate a random-access, STL-like start and end iterator with c...
Definition range.h:64
Simple structure used to manage references/pointers to remote instances.
Definition worldref.h:395
double weights(const unsigned int &i) const
return the weight
Definition srconf.h:671
const Tensor< T > flat_vector(const unsigned int &idim) const
return shallow copy of a slice of one of the vectors, flattened to (r,kVec)
Definition srconf.h:545
Definition SVDTensor.h:42
long rank() const
Definition SVDTensor.h:77
A slice defines a sub-range or patch of a dimension.
Definition slice.h:103
static TaskAttributes hipri()
Definition thread.h:456
static TaskAttributes generator()
Definition thread.h:448
Traits class to specify support of numeric types.
Definition type_data.h:56
A tensor is a multidimensional array.
Definition tensor.h:317
float_scalar_type normf() const
Returns the Frobenius norm of the tensor.
Definition tensor.h:1726
Tensor< T > & fill(T x)
Inplace fill with a scalar (legacy name)
Definition tensor.h:562
T sum() const
Returns the sum of all elements of the tensor.
Definition tensor.h:1662
Tensor< T > reshape(int ndimnew, const long *d)
Returns new view/tensor reshaping size/number of dimensions to conforming tensor.
Definition tensor.h:1384
T * ptr()
Returns a pointer to the internal data.
Definition tensor.h:1825
IsSupported< TensorTypeData< Q >, Tensor< T > & >::type scale(Q x)
Inplace multiplication by scalar of supported type (legacy name)
Definition tensor.h:686
Tensor< T > & emul(const Tensor< T > &t)
Inplace multiply by corresponding elements of argument Tensor.
Definition tensor.h:1799
bool has_data() const
Definition tensor.h:1887
iterator begin()
Returns an iterator to the beginning of the local data (no communication)
Definition worlddc.h:1357
implT::const_iterator const_iterator
Definition worlddc.h:1135
iterator end()
Returns an iterator past the end of the local data (no communication)
Definition worlddc.h:1371
Future< iterator > futureT
Definition worlddc.h:1138
implT::iterator iterator
Definition worlddc.h:1134
implT::accessor accessor
Definition worlddc.h:1136
void fence(bool debug=false)
Synchronizes all processes in communicator AND globally ensures no pending AM or tasks.
Definition worldgop.cc:161
A parallel world class.
Definition world.h:132
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
ProcessID nproc() const
Returns the number of processes in this World (same as MPI_Comm_size()).
Definition world.h:325
Wrapper for an opaque pointer for serialization purposes.
Definition archive.h:851
syntactic sugar for std::array<bool, N>
Definition array_of_bools.h:19
char * p(char *buf, const char *name, int k, int initial_level, double thresh, int order)
Definition derivatives.cc:72
const std::size_t bufsize
Definition derivatives.cc:16
static double lo
Definition dirac-hatom.cc:23
static bool debug
Definition dirac-hatom.cc:16
Provides FunctionCommonData, FunctionImpl and FunctionFactory.
static double function(const coord_3d &r)
Normalized gaussian.
Definition functionio.cc:100
auto T(World &world, response_space &f) -> response_space
Definition global_functions.cc:28
Tensor< TENSOR_RESULT_TYPE(T, Q) > & fast_transform(const Tensor< T > &t, const Tensor< Q > &c, Tensor< TENSOR_RESULT_TYPE(T, Q) > &result, Tensor< TENSOR_RESULT_TYPE(T, Q) > &workspace)
Restricted but heavily optimized form of transform()
Definition tensor.h:2444
const double beta
Definition gygi_soltion.cc:62
static const double v
Definition hatom_sf_dirac.cc:20
static double u(double r, double c)
Definition he.cc:20
Tensor< double > op(const Tensor< double > &x)
Definition kain.cc:508
static double pow(const double *a, const double *b)
Definition lda.h:74
Macros and tools pertaining to the configuration of MADNESS.
#define MADNESS_PRAGMA_CLANG(x)
Definition madness_config.h:200
#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_EXCEPTION(msg, value)
Macro for throwing a MADNESS exception.
Definition madness_exception.h:119
#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
static const bool VERIFY_TREE
Definition mra.h:57
Definition potentialmanager.cc:41
Namespace for all elements and tools of MADNESS.
Definition DFParameters.h:10
bool two_scale_hg(int k, Tensor< double > *hg)
Definition twoscale.cc:151
@ BC_FREE
Definition bc.h:53
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
static const char * filename
Definition legendre.cc:96
static const std::vector< Slice > ___
Entire dimension.
Definition slice.h:128
static bool enforce_in_volume(Level n, const Translation &l)
Definition mraimpl.h:3213
double abs(double x)
Definition complexfun.h:48
static double cpu_time()
Returns the cpu time in seconds relative to an arbitrary origin.
Definition timers.h:127
Vector< double, 3 > coordT
Definition corepotential.cc:54
GenTensor< TENSOR_RESULT_TYPE(R, Q)> general_transform(const GenTensor< R > &t, const Tensor< Q > c[])
Definition gentensor.h:274
response_space scale(response_space a, double b)
void legendre_scaling_functions(double x, long k, double *p)
Evaluate the first k Legendre scaling functions.
Definition legendre.cc:85
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
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
static Key< NDIM > simpt2key(const Vector< T, NDIM > &pt, Level n)
Definition funcdefaults.h:453
TreeState
Definition funcdefaults.h:59
@ nonstandard_after_apply
s and d coeffs, state after operator application
Definition funcdefaults.h:64
@ on_demand
no coeffs anywhere, but a functor providing if necessary
Definition funcdefaults.h:67
@ 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
static void user_to_sim(const Vector< double, NDIM > &xuser, Vector< double, NDIM > &xsim)
Convert user coords (cell[][]) to simulation coords ([0,1]^ndim)
Definition funcdefaults.h:444
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
Tensor< double > tensorT
Definition distpm.cc:21
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
int64_t Translation
Definition key.h:57
void plotdx(const Function< T, NDIM > &f, const char *filename, const Tensor< double > &cell=FunctionDefaults< NDIM >::get_cell(), const std::vector< long > &npt=std::vector< long >(NDIM, 201L), bool binary=true)
Writes an OpenDX format file with a cube/slice of points on a uniform grid.
Definition mraimpl.h:3453
Function< T, NDIM > mirror(const Function< T, NDIM > &f, const std::vector< long > &mirrormap, bool fence=true)
Generate a new function by mirroring within the dimensions .. optional fence.
Definition mra.h:2383
static void verify_tree(World &world, const std::vector< Function< T, NDIM > > &v)
Definition SCF.cc:74
static const Slice _(0,-1, 1)
void change_tensor_type(GenTensor< T > &t, const TensorArgs &targs)
change representation to targ.tt
Definition gentensor.h:284
int Level
Definition key.h:58
std::enable_if< std::is_base_of< ProjectorBase, projT >::value, OuterProjector< projT, projQ > >::type outer(const projT &p0, const projQ &p1)
Definition projector.h:457
TreeState get_tree_state(const Function< T, NDIM > &f)
get tree state of a function
Definition mra.h:2841
bool gauss_legendre(int n, double xlo, double xhi, double *x, double *w)
Definition legendre.cc:226
static double pop(std::vector< double > &v)
Definition SCF.cc:115
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:226
Tensor< T > fcube(const Key< NDIM > &, T(*f)(const Vector< double, NDIM > &), const Tensor< double > &)
Definition mraimpl.h:2133
TensorType
low rank representations of tensors (see gentensor.h)
Definition gentensor.h:120
@ TT_2D
Definition gentensor.h:120
@ TT_FULL
Definition gentensor.h:120
static void dxprintvalue(FILE *f, const double t)
Definition mraimpl.h:3444
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
NDIM & f
Definition mra.h:2528
void error(const char *msg)
Definition world.cc:142
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:2854
NDIM const Function< R, NDIM > & g
Definition mra.h:2528
double wall_time()
Returns the wall time in seconds relative to an arbitrary origin.
Definition timers.cc:48
static bool print_timings
Definition SCF.cc:106
constexpr Vector< T, sizeof...(Ts)+1 > vec(T t, Ts... ts)
Factory function for creating a madness::Vector.
Definition vector.h:750
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:2484
static std::vector< double > ttt
Definition SCF.cc:107
Function< T, NDIM > project(const Function< T, NDIM > &other, int k=FunctionDefaults< NDIM >::get_k(), double thresh=FunctionDefaults< NDIM >::get_thresh(), bool fence=true)
Definition mra.h:2511
std::string name(const FuncType &type, const int ex=-1)
Definition ccpairfunction.h:28
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:2096
static const int MAXK
The maximum wavelet order presently supported.
Definition funcdefaults.h:54
static bool enforce_bc(bool is_periodic, Level n, Translation &l)
Definition mraimpl.h:3193
bool isnan(const std::complex< T > &v)
Definition mraimpl.h:53
const double mu
Definition navstokes_cosines.cc:95
static const double b
Definition nonlinschro.cc:119
static const double d
Definition nonlinschro.cc:121
static const double a
Definition nonlinschro.cc:118
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
static const long k
Definition rk.cc:44
const double xi
Exponent for delta function approx.
Definition siam_example.cc:60
Definition test_ar.cc:204
Definition test_ccpairfunction.cc:22
add two functions f and g: result=alpha * f + beta * g
Definition funcimpl.h:3598
"put" this on g
Definition funcimpl.h:2654
change representation of nodes' coeffs to low rank, optional fence
Definition funcimpl.h:2687
check symmetry wrt particle exchange
Definition funcimpl.h:2360
compute the norm of the wavelet coefficients
Definition funcimpl.h:4495
Definition funcimpl.h:2714
Definition funcimpl.h:1474
mirror dimensions of this, write result on f
Definition funcimpl.h:2588
map this on f
Definition funcimpl.h:2508
mirror dimensions of this, write result on f
Definition funcimpl.h:2538
Definition funcimpl.h:5590
reduce the rank of the nodes, optional fence
Definition funcimpl.h:2334
remove all coefficients of internal nodes
Definition funcimpl.h:2280
Definition funcimpl.h:4567
shallow-copy, pared-down version of FunctionNode, for special purpose only
Definition funcimpl.h:749
TensorArgs holds the arguments for creating a LowRankTensor.
Definition gentensor.h:134
double thresh
Definition gentensor.h:135
Definition mraimpl.h:3101
void operator()(const Key< NDIM > &key, Tensor< T > &t) const
Definition mraimpl.h:3102
void serialize(Archive &ar)
Definition mraimpl.h:3103
Definition mraimpl.h:3107
void operator()(const Key< NDIM > &key, Tensor< T > &t) const
Definition mraimpl.h:3108
void serialize(Archive &ar)
Definition mraimpl.h:3109
Definition mraimpl.h:3069
void operator()(const A &a, const B &b) const
Definition mraimpl.h:3070
void serialize(Archive &ar)
Definition mraimpl.h:3072
Definition mraimpl.h:3076
void operator()(const Key< NDIM > &key, FunctionNode< T, NDIM > &node) const
Definition mraimpl.h:3084
void serialize(Archive &ar)
Definition mraimpl.h:3087
T q
Definition mraimpl.h:3077
scaleinplace()
Definition mraimpl.h:3078
scaleinplace(T q)
Definition mraimpl.h:3080
void operator()(const Key< NDIM > &key, Tensor< T > &t) const
Definition mraimpl.h:3081
Definition mraimpl.h:3093
void serialize(Archive &ar)
Definition mraimpl.h:3097
void operator()(const Key< NDIM > &key, Tensor< T > &t) const
Definition mraimpl.h:3094
insert/replaces the coefficients into the function
Definition funcimpl.h:692
Definition lowrankfunction.h:336
int np
Definition tdse1d.cc:165
static const double s0
Definition tdse4.cc:83
AtomicInt sum
Definition test_atomicint.cc:46
int me
Definition test_binsorter.cc:10
double norm(const T i1)
Definition test_cloud.cc:85
int task(int i)
Definition test_runtime.cpp:4
void e()
Definition test_sig.cc:75
#define N
Definition testconv.cc:37
static const double alpha
Definition testcosine.cc:10
static const int truncate_mode
Definition testcosine.cc:14
constexpr std::size_t NDIM
Definition testgconv.cc:54
double h(const coord_1d &r)
Definition testgconv.cc:175
double g1(const coord_t &r)
Definition testgconv.cc:122
std::size_t axis
Definition testpdiff.cc:59
double k0
Definition testperiodic.cc:66
Defines and implements WorldObject.
Implements WorldContainer.
Defines and implements a concurrent hashmap.
#define PROFILE_FUNC
Definition worldprofile.h:209
#define PROFILE_MEMBER_FUNC(classname)
Definition worldprofile.h:210
#define PROFILE_BLOCK(name)
Definition worldprofile.h:208
int ProcessID
Used to clearly identify process number/rank.
Definition worldtypes.h:43
Key< D > keyT
Definition writecoeff2.cc:11
void test()
Definition y.cc:696