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)"
52 bool isnan(
const std::complex<T>&
v) {
66 template <
typename T, std::
size_t NDIM>
68 if (!
two_scale_hg(
k, &hg))
throw "failed to get twoscale coefficients";
75 h1 =
copy(hg(sk,sk2));
76 g0 =
copy(hg(sk2,sk));
86 template <
typename T, std::
size_t NDIM>
96 for (
int mu=0;
mu<npt; ++
mu) {
99 for (
int j=0; j<
k; ++j) {
100 quad_phi(
mu,j) = phi[j];
101 quad_phiw(
mu,j) = quad_w(
mu)*phi[j];
108 template <
typename T, std::
size_t NDIM>
116 const keyT& key = it->first;
117 const nodeT& node = it->second;
120 if (is_compressed()) {
140 print(world.rank(),
"FunctionImpl: verify: INCONSISTENT TREE NODE, key =", key,
", node =", node,
141 ", dim[0] =",node.
coeff().
dim(0),
", compressed =",is_compressed());
148 for (
typename dcT::const_iterator it=coeffs.begin(); it!=coeffs.end(); ++it) {
149 const keyT& key = it->first;
150 const nodeT& node = it->second;
152 if (key.level() > 0) {
153 const keyT parent = key.parent();
154 typename dcT::const_iterator pit = coeffs.find(parent).get();
155 if (pit == coeffs.end()) {
156 print(world.rank(),
"FunctionImpl: verify: MISSING PARENT",key,parent);
160 const nodeT& pnode = pit->second;
161 if (!pnode.has_children()) {
162 print(world.rank(),
"FunctionImpl: verify: PARENT THINKS IT HAS NO CHILDREN",key,parent);
169 typename dcT::const_iterator cit = coeffs.find(kit.key()).get();
170 if (cit == coeffs.end()) {
171 if (node.has_children()) {
172 print(world.rank(),
"FunctionImpl: verify: MISSING CHILD",key,kit.key());
178 if (! node.has_children()) {
179 print(world.rank(),
"FunctionImpl: verify: UNEXPECTED CHILD",key,kit.key());
191 template <
typename T, std::
size_t NDIM>
193 return coeffs.get_pmap();
207 template <
typename T, std::
size_t NDIM>
209 const double beta,
const implT&
g,
const bool fence) {
214 ProcessID owner = coeffs.owner(cdata.key0);
215 if (world.rank() == owner) {
223 apply_opT apply_op(
this);
225 woT::task(world.rank(), &implT:: template forward_traverse<coeff_opT,apply_opT>,
226 coeff_op, apply_op, cdata.key0);
230 if (fence) world.gop.fence();
234 template <
typename T, std::
size_t NDIM>
240 template <
typename T, std::
size_t NDIM>
246 template <
typename T, std::
size_t NDIM>
251 template <
typename T, std::
size_t NDIM>
256 template <
typename T, std::
size_t NDIM>
261 template <
typename T, std::
size_t NDIM>
266 template <
typename T, std::
size_t NDIM>
271 template <
typename T, std::
size_t NDIM>
278 template <
typename T, std::
size_t NDIM>
284 template <
typename T, std::
size_t NDIM>
290 template <
typename T, std::
size_t NDIM>
297 template <
typename T, std::
size_t NDIM>
300 template <
typename T, std::
size_t NDIM>
303 template <
typename T, std::
size_t NDIM>
306 template <
typename T, std::
size_t NDIM>
309 template <
typename T, std::
size_t NDIM>
312 template <
typename T, std::
size_t NDIM>
315 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>
332 timer_accumulate.accumulate(time);
335 template <
typename T, std::
size_t NDIM>
337 if (world.rank()==0) {
338 timer_accumulate.print(
"accumulate");
339 timer_target_driven.print(
"target_driven");
340 timer_lr_result.print(
"result2low_rank");
344 template <
typename T, std::
size_t NDIM>
346 if (world.rank()==0) {
347 timer_accumulate.reset();
348 timer_target_driven.reset();
349 timer_lr_result.reset();
356 template <
typename T, std::
size_t NDIM>
361 if (world.rank() == coeffs.owner(cdata.key0)) {
362 if (is_compressed()) {
363 truncate_spawn(cdata.key0,tol);
365 truncate_reconstructed_spawn(cdata.key0,tol);
372 template <
typename T, std::
size_t NDIM>
381 template <
typename T, std::
size_t NDIM>
388 std::vector<Tensor<double> > localinfo_vec(1,localinfo);
389 std::vector<Tensor<double> > printinfo=world.gop.concat0(localinfo_vec);
393 if (world.rank()==0) do_print_plane(
filename,printinfo,xaxis,yaxis,el2);
401 template <
typename T, std::
size_t NDIM>
404 user_to_sim<NDIM>(el2,x_sim);
414 const keyT& key = it->first;
415 const nodeT& node = it->second;
423 double scale=std::pow(0.5,
double(n));
424 double xloleft =
scale*l[xaxis];
425 double yloleft =
scale*l[yaxis];
426 double xhiright =
scale*(l[xaxis]+1);
427 double yhiright =
scale*(l[yaxis]+1);
438 if ((user[0]<-5.0) or (user[1]<-5.0) or (user[2]>5.0) or (user[3]>5.0))
continue;
444 const double maxrank=40;
450 const int npt = cdata.npt + 1;
461 plotinfo(counter,0)=color;
462 plotinfo(counter,1)=user[0];
463 plotinfo(counter,2)=user[1];
464 plotinfo(counter,3)=user[2];
465 plotinfo(counter,4)=user[3];
472 else plotinfo=plotinfo(
Slice(0,counter-1),
Slice(
_));
477 template <
typename T, std::
size_t NDIM>
479 const int xaxis,
const int yaxis,
const coordT el2) {
486 pFile = fopen(
filename.c_str(),
"w");
490 fprintf(pFile,
"\\psset{unit=1cm}\n");
491 fprintf(pFile,
"\\begin{pspicture}(%4.2f,%4.2f)(%4.2f,%4.2f)\n",
494 fprintf(pFile,
"\\pslinewidth=0.1pt\n");
496 for (std::vector<
Tensor<double> >::const_iterator it=plotinfo.begin(); it!=plotinfo.end(); ++it) {
501 for (
long i=0; i<localinfo.
dim(0); ++i) {
503 fprintf(pFile,
"\\newhsbcolor{mycolor}{%8.4f 1.0 0.7}\n",localinfo(i,0));
504 fprintf(pFile,
"\\psframe["
507 "(%12.8f,%12.8f)(%12.8f,%12.8f)\n",
508 localinfo(i,1),localinfo(i,2),localinfo(i,3),localinfo(i,4));
514 fprintf(pFile,
"\\end{pspicture}\n");
520 template <
typename T, std::
size_t NDIM>
524 std::vector<keyT> local_keys=local_leaf_keys();
527 std::vector<keyT> all_keys=world.gop.concat0(local_keys);
531 if (world.rank()==0) do_print_grid(
filename,all_keys);
536 template <
typename T, std::
size_t NDIM>
540 std::vector<keyT> keys(coeffs.size());
547 const keyT& key = it->first;
548 const nodeT& node = it->second;
549 if (node.
is_leaf()) keys[i++]=key;
562 template <
typename T, std::
size_t NDIM>
569 const size_t npt = qx.
dim(0);
572 long npoints=power<NDIM>(npt);
574 long nboxes=keys.size();
578 pFile = fopen(
filename.c_str(),
"w");
580 fprintf(pFile,
"%ld\n",npoints*nboxes);
581 fprintf(pFile,
"%ld points per box and %ld boxes \n",npoints,nboxes);
584 typename std::vector<keyT>::const_iterator key_it=keys.begin();
585 for (key_it=keys.begin(); key_it!=keys.end(); ++key_it) {
587 const keyT& key=*key_it;
588 fprintf(pFile,
"# key: %8d",key.
level());
595 const double h = std::pow(0.5,
double(n));
602 for (
size_t i=0; i<npt; ++i) {
603 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
604 for (
size_t j=0; j<npt; ++j) {
605 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
606 for (
size_t k=0;
k<npt; ++
k) {
607 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
613 fprintf(pFile,
"%18.12f %18.12f %18.12f\n",
c[0],
c[1],
c[2]);
627 template <
typename T, std::
size_t NDIM>
633 const int MAXLEVEL1 = 20;
634 const int MAXLEVEL2 = 10;
641 return tol*std::min(1.0,
pow(0.5,
double(std::min(key.
level(),MAXLEVEL1)))*
L);
645 return tol*std::min(1.0,
pow(0.25,
double(std::min(key.
level(),MAXLEVEL2)))*
L*
L);
662 const static double fac=1.0/std::pow(2,
NDIM*0.5);
666 return tol*std::min(1.0,
pow(0.5,
double(std::min(key.
level(),MAXLEVEL1)))*
L);
674 template <
typename T, std::
size_t NDIM>
676 std::vector<Slice> s(
NDIM);
678 for (std::size_t i=0; i<
NDIM; ++i)
679 s[i] = cdata.s[l[i]&1];
685 template <
typename T, std::
size_t NDIM>
687 const keyT& child,
const keyT& parent,
const coeffT& coeff)
const {
697 if (coeff.
dim(0)==2*
f->get_k()) result=coeff;
698 else if (coeff.
dim(0)==
f->get_k()) {
699 result(
f->cdata.s0)+=coeff;
708 const coeffT scoeff=
f->parent_to_child(coeff,parent,child);
709 result(
f->cdata.s0)+=scoeff;
717 template <
typename T, std::
size_t NDIM>
722 for (
typename dcT::iterator it= coeffs.begin(); it!=end; ++it) {
724 nodeT& node=it->second;
725 if (key.
level()>max_level) coeffs.erase(key);
728 this->undo_redundant(
true);
733 template <
typename T, std::
size_t NDIM>
747 template <
typename T, std::
size_t NDIM>
749 const std::vector<tensorT>&
c,
751 if (key == cdata.key0 && coeffs.owner(key)!=world.rank())
return;
754 std::unique_ptr<typename dcT::accessor[]> acc(
new typename dcT::accessor[
v.size()]);
755 for (
unsigned int i=0; i<
c.size(); i++) {
758 bool exists = !
v[i]->coeffs.insert(acc[i],key);
770 for (
unsigned int i=0; i<
v.size(); i++) {
771 done &= acc[i]->second.has_coeff();
776 std::vector<tensorT>
d(
v.size());
777 for (
unsigned int i=0; i<
v.size(); i++) {
778 if (acc[i]->second.has_coeff()) {
781 s(cdata.s0) = acc[i]->second.coeff().full_tensor_copy();
782 acc[i]->second.clear_coeff();
784 acc[i]->second.set_has_children(
true);
790 const keyT& child = kit.key();
791 std::vector<Slice> cp = child_patch(child);
792 std::vector<tensorT> childc(
v.size());
793 for (
unsigned int i=0; i<
v.size(); i++) {
794 if (
d[i].size()) childc[i] =
copy(
d[i](cp));
796 woT::task(coeffs.owner(child), &implT::refine_to_common_level,
v, childc, child);
802 template <
typename T, std::
size_t NDIM>
804 if (world.size()> 1000)
807 box_interior[from] = ni;
811 template <
typename T, std::
size_t NDIM>
813 if (world.size() >= 1000)
815 for (
int i=0; i<world.size(); ++i)
816 box_leaf[i] = box_interior[i] == 0;
818 long nleaf=0, ninterior=0;
821 const nodeT& node = it->second;
827 this->send(0, &implT::put_in_box, world.rank(), nleaf, ninterior);
829 if (world.rank() == 0) {
830 for (
int i=0; i<world.size(); ++i) {
831 printf(
"load: %5d %8ld %8ld\n", i, box_leaf[i], box_interior[i]);
837 template <
typename T, std::
size_t NDIM>
843 template <
typename T, std::
size_t NDIM>
847 double test = 2*
lo*hi + hi*hi;
854 template <
typename T, std::
size_t NDIM>
857 coeffs.insert(acc,key);
858 nodeT& node = acc->second;
880 const keyT& child = kit.key();
881 if (
d.size() > 0) ss =
copy(
d(child_patch(child)));
883 woT::task(coeffs.owner(child), &implT::sum_down_spawn, child, ss);
888 if (
c.size() <= 0)
c =
coeffT(cdata.vk,targs);
893 template <
typename T, std::
size_t NDIM>
896 if (world.rank() == coeffs.owner(cdata.key0)) sum_down_spawn(cdata.key0,
coeffT());
897 if (fence) world.gop.fence();
901 template <
typename T, std::
size_t NDIM>
905 const std::pair<keyT,coeffT>& left,
906 const std::pair<keyT,coeffT>& center,
907 const std::pair<keyT,coeffT>& right) {
908 D->forward_do_diff1(
f,
this,key,left,center,right);
912 template <
typename T, std::
size_t NDIM>
916 const std::pair<keyT,coeffT>& left,
917 const std::pair<keyT,coeffT>& center,
918 const std::pair<keyT,coeffT>& right) {
919 D->do_diff1(
f,
this,key,left,center,right);
924 template <
typename T, std::
size_t NDIM>
926 typedef std::pair<keyT,coeffT> argT;
927 for (
const auto& [key, node]:
f->coeffs) {
928 if (node.has_coeff()) {
930 argT center(key,node.coeff());
938 if (fence) world.gop.fence();
943 template <
typename T, std::
size_t NDIM>
956 template <
typename T, std::
size_t NDIM>
964 std::vector<long> vkhalf=std::vector<long>(
NDIM/2,cdata.vk[0]);
991 template <
typename T, std::
size_t NDIM>
998 if (ket_only)
return coeff_ket;
1001 coeffT val_ket=coeffs2values(key,coeff_ket);
1019 coeff_result=
coeffT(values2coeffs(key,val_ket2),this->get_tensor_args());
1024 val_ket=val_ket.
convert(get_tensor_args());
1026 coeff_result=values2coeffs(key,val_result);
1030 return coeff_result;
1035 template <
typename T, std::
size_t NDIM>
1039 const_cast<implT*
>(&
f)->flo_unary_op_node_inplace(
do_mapdim(map,*
this),fence);
1044 template <
typename T, std::
size_t NDIM>
1047 const_cast<implT*
>(&
f)->flo_unary_op_node_inplace(
do_mirror(mirrormap,*
this),fence);
1054 template <
typename T, std::
size_t NDIM>
1056 const std::vector<long>&
mirror,
bool fence) {
1066 template <
typename T, std::
size_t NDIM>
1070 this->scale_inplace(0.5,
true);
1077 template <
typename T, std::
size_t NDIM>
1085 template <
typename T, std::
size_t NDIM>
1093 template <
typename T, std::
size_t NDIM>
1095 std::list<keyT> to_be_erased;
1096 for (
auto it=coeffs.begin(); it!=coeffs.end(); ++it) {
1097 const keyT& key=it->first;
1098 nodeT& node=it->second;
1100 if (key.
level()>n) to_be_erased.push_back(key);
1102 for (
auto& key : to_be_erased) coeffs.erase(key);
1109 template <
typename T, std::
size_t NDIM>
1112 flo_unary_op_node_inplace(
1130 template <
typename T, std::
size_t NDIM>
1138 template <
typename T, std::
size_t NDIM>
1159 template <
typename T, std::
size_t NDIM>
1167 template <
typename T, std::
size_t NDIM>
1179 template <
typename T, std::
size_t NDIM>
1185 const tensorT h[2] = {cdata.h0T, cdata.h1T};
1193 for (
size_t ii=0; ii<
NDIM; ++ii) matrices[ii]=
h[kit.key().translation()[ii]%2];
1209 template <
typename T, std::
size_t NDIM>
1214 const tensorT h[2] = {cdata.h0, cdata.h1};
1227 template <
typename T, std::
size_t NDIM>
1229 long kmin = std::min(cdata.k,old.
cdata.k);
1230 std::vector<Slice> s(
NDIM,
Slice(0,kmin-1));
1233 const keyT& key = it->first;
1234 const nodeT& node = it->second;
1238 coeffs.replace(key,
nodeT(
c,
false));
1248 template <
typename T, std::
size_t NDIM>
1250 return coeffs.probe(key) && coeffs.find(key).get()->second.has_children();
1253 template <
typename T, std::
size_t NDIM>
1255 return coeffs.probe(key) && (not coeffs.find(key).get()->second.has_children());
1259 template <
typename T, std::
size_t NDIM>
1261 for (
unsigned int i=0; i<
v.size(); ++i) {
1263 refine_op(true_refine_test(), key);
1270 template <
typename T, std::
size_t NDIM>
1273 for (
typename dcT::iterator it=coeffs.begin(); it!=end; ++it) {
1274 it->second.set_norm_tree(0.0);
1275 it->second.set_snorm(0.0);
1276 it->second.set_dnorm(0.0);
1281 template <
typename T, std::
size_t NDIM>
1284 for (
typename dcT::iterator it=coeffs.begin(); it!=end; ++it) {
1287 const auto found = coeffs.find(acc,key);
1289 nodeT& node = acc->second;
1290 if (node.has_coeff() &&
1291 node.get_norm_tree() != -1.0 &&
1292 node.coeff().normf() >= truncate_tol(
thresh,key)) {
1294 node.set_norm_tree(-1.0);
1297 int ndir =
static_cast<int>(std::pow(
static_cast<double>(3),
static_cast<int>(
NDIM)));
1298 std::vector< Future <bool> >
v = future_vector_factory<bool>(ndir);
1303 for (std::size_t
d=0;
d<
NDIM; ++
d) {
1304 const int odd = key.translation()[
d] & 0x1L;
1311 keyT neigh = neighbor(key, keyT(key.level(),l), is_periodic);
1313 if (neigh.is_valid()) {
1314 v[i++] = this->
task(coeffs.owner(neigh), &implT::exists_and_has_children, neigh);
1320 woT::task(world.rank(), &implT::broaden_op, key,
v);
1332 template <
typename T, std::
size_t NDIM>
1335 if (world.rank() == coeffs.owner(cdata.key0))
1336 woT::task(world.rank(), &implT::trickle_down_op, cdata.key0,
coeffT());
1337 if (fence) world.gop.fence();
1343 template <
typename T, std::
size_t NDIM>
1354 if (it == coeffs.end()) {
1356 it = coeffs.find(key).get();
1363 if (node.coeff().has_no_data()) node.coeff()=
coeffT(cdata.vk,targs);
1366 if (node.has_children()) {
1368 if (key.
level() > 0)
d += s;
1371 const keyT& child = kit.key();
1375 woT::task(coeffs.owner(child), &implT::trickle_down_op, child, ss);
1380 node.coeff().reduce_rank(
thresh);
1385 template <
typename T, std::
size_t NDIM>
1389 if (current_state==finalstate)
return;
1403 bool must_fence=
false;
1410 if (current_state==
redundant) remove_internal_coefficients(fence);
1418 remove_internal_coefficients(
true);
1434 remove_internal_coefficients(
true);
1455 remove_internal_coefficients(
true);
1475 remove_internal_coefficients(
true);
1485 if (must_fence and world.rank()==0) {
1486 print(
"could not respect fence in change_tree_state");
1494 template <
typename T, std::
size_t NDIM>
1497 if (is_reconstructed())
return;
1499 if (is_redundant() or is_nonstandard_with_leaves()) {
1501 this->remove_internal_coefficients(fence);
1505 if (world.rank() == coeffs.owner(cdata.key0))
1506 woT::task(world.rank(), &implT::reconstruct_op, cdata.key0,
coeffT(),
true);
1507 }
else if (is_nonstandard()) {
1510 if (world.rank() == coeffs.owner(cdata.key0))
1511 woT::task(world.rank(), &implT::reconstruct_op, cdata.key0,
coeffT(),
false);
1515 if (fence) world.gop.fence();
1526 template <
typename T, std::
size_t NDIM>
1530 set_tree_state(newstate);
1533 bool redundant1=(tree_state==
redundant);
1535 if (world.rank() == coeffs.owner(cdata.key0)) {
1537 compress_spawn(cdata.key0, nonstandard1, keepleaves1, redundant1);
1543 template <
typename T, std::
size_t NDIM>
1548 template <
typename T, std::
size_t NDIM>
1554 template <
typename T, std::
size_t NDIM>
1558 if (is_redundant())
return;
1559 MADNESS_CHECK_THROW(is_reconstructed(),
"impl::make_redundant() wants a reconstructed tree");
1564 template <
typename T, std::
size_t NDIM>
1573 template <
typename T, std::
size_t NDIM>
1575 if (world.rank() == coeffs.owner(cdata.key0))
1576 norm_tree_spawn(cdata.key0);
1581 template <
typename T, std::
size_t NDIM>
1587 double value =
v[i].get();
1591 coeffs.task(key, &nodeT::set_norm_tree,
sum);
1596 template <
typename T, std::
size_t NDIM>
1598 nodeT& node = coeffs.find(key).get()->second;
1600 std::vector< Future<double> >
v = future_vector_factory<double>(1<<
NDIM);
1603 v[i] = woT::task(coeffs.owner(kit.key()), &implT::norm_tree_spawn, kit.key());
1605 return woT::task(world.rank(),&implT::norm_tree_op, key,
v);
1619 template <
typename T, std::
size_t NDIM>
1622 nodeT& node = coeffs.find(key).get()->second;
1629 std::vector<Future<coeffT> >
v = future_vector_factory<coeffT>(1<<
NDIM);
1632 v[i] = woT::task(coeffs.owner(kit.key()), &implT::truncate_reconstructed_spawn, kit.key(),tol,
TaskAttributes::hipri());
1643 template <
typename T, std::
size_t NDIM>
1650 for (
size_t i=0; i<
v.size(); ++i)
if (
v[i].get().has_no_data())
return coeffT();
1657 const auto found = coeffs.find(acc, key);
1663 d(child_patch(kit.key())) +=
v[i].get().full_tensor_copy();
1669 const double error=
d.normf();
1671 nodeT& node = coeffs.find(key).get()->second;
1674 node.set_has_children(
false);
1676 coeffs.erase(kit.key());
1680 acc->second.set_coeff(ss);
1694 template <
typename T, std::
size_t NDIM>
1696 const std::vector<
Future<std::pair<coeffT,double> > >&
v,
bool nonstandard1) {
1704 double norm_tree2=0.0;
1707 d(child_patch(kit.key())) +=
v[i].get().first.full_tensor_copy();
1708 norm_tree2+=
v[i].get().second*
v[i].get().second;
1713 timer_filter.accumulate(cpu1-cpu0);
1717 const auto found = coeffs.find(acc, key);
1719 MADNESS_CHECK_THROW(!acc->second.has_coeff(),
"compress_op: existing coeffs where there should be none");
1727 double snorm=ss.normf();
1729 if (key.
level()> 0 && !nonstandard1)
d(cdata.s0) = 0.0;
1732 double dnorm=dd.normf();
1735 acc->second.set_snorm(snorm);
1736 acc->second.set_dnorm(dnorm);
1739 acc->second.set_coeff(dd);
1741 timer_compress_svd.accumulate(cpu1-cpu0);
1744 return std::make_pair(ss,snorm);
1753 template <
typename T, std::
size_t NDIM>
1754 std::pair<typename FunctionImpl<T,NDIM>::coeffT,
double>
1759 double norm_tree2=0.0;
1761 d(child_patch(kit.key())) +=
v[i].get().first.full_tensor_copy();
1762 norm_tree2+=
v[i].get().second*
v[i].get().second;
1774 double dnorm=
d.normf();
1775 double snorm=s.normf();
1778 const auto found = coeffs.find(acc, key);
1781 acc->second.set_coeff(s);
1782 acc->second.set_dnorm(dnorm);
1783 acc->second.set_snorm(snorm);
1791 template <
typename T, std::
size_t NDIM>
1794 if (is_compressed())
return;
1796 flo_unary_op_node_inplace(
do_standard(
this),fence);
1804 template <
typename T, std::
size_t NDIM>
1814 if (printme) printf(
"time in consolidate_buffer %8.4f\n",end1-begin1);
1822 if (printme) printf(
"time in do_reduce_rank %8.4f\n",end1-begin1);
1828 if (printme) printf(
"time in do_change_tensor_type %8.4f\n",end1-begin1);
1834 if (printme) printf(
"time in do_truncate_NS_leafs %8.4f\n",end1-begin1);
1837 double elapsed=end-begin;
1847 template <
typename T, std::
size_t NDIM>
1856 template <
typename T, std::
size_t NDIM>
1868 template <
typename T, std::
size_t NDIM>
1870 std::size_t maxdepth = 0;
1873 std::size_t
N = (std::size_t) it->first.level();
1882 template <
typename T, std::
size_t NDIM>
1884 std::size_t maxdepth = max_local_depth();
1885 world.gop.max(maxdepth);
1890 template <
typename T, std::
size_t NDIM>
1892 std::size_t maxsize = 0;
1893 maxsize = coeffs.size();
1894 world.gop.max(maxsize);
1899 template <
typename T, std::
size_t NDIM>
1901 std::size_t minsize = 0;
1902 minsize = coeffs.size();
1903 world.gop.min(minsize);
1908 template <
typename T, std::
size_t NDIM>
1910 std::size_t
sum = 0;
1911 sum = coeffs.size();
1917 template <
typename T, std::
size_t NDIM>
1919 std::size_t
sum = 0;
1920 for (
const auto& [key,node] : coeffs) {
1921 if (node.has_coeff())
sum+=node.size();
1927 template <
typename T, std::
size_t NDIM>
1929 std::size_t
sum = size_local();
1935 template <
typename T, std::
size_t NDIM>
1937 std::size_t
sum = coeffs.size() * (
sizeof(
keyT) +
sizeof(
nodeT));
1940 const nodeT& node = it->second;
1948 template <
typename T, std::
size_t NDIM>
1950 std::size_t
sum = coeffs.size() * (
sizeof(
keyT) +
sizeof(
nodeT));
1953 const nodeT& node = it->second;
1962 template <
typename T, std::
size_t NDIM>
1964 const size_t tsize=this->tree_size();
1966 const size_t ncoeff=this->nCoeff();
1968 const double d=
sizeof(
T);
1969 const double fac=1024*1024*1024;
1973 double local = norm2sq_local();
1974 this->world.gop.sum(local);
1975 this->world.gop.fence();
1979 if (this->world.rank()==0) {
1981 constexpr std::size_t
bufsize=128;
1983 snprintf(buf,
bufsize,
"%40s at time %.1fs: norm/tree/#coeff/size: %7.5f %zu, %6.3f m, %6.3f GByte",
1984 (
name.c_str()), wall,
norm, tsize,
double(ncoeff)*1.e-6,
double(ncoeff)/fac*
d);
1985 print(std::string(buf));
1990 template <
typename T, std::
size_t NDIM>
1992 if (this->targs.tt==
TT_FULL)
return;
1995 if (is_compressed())
k0=2*
k;
2000 if (world.rank()==0)
print(
"n.size(),k0,dim",n.
size(),
k0,dim);
2003 const nodeT& node = it->second;
2019 if (world.rank()==0) {
2020 print(
"configurations number of nodes");
2021 print(
" full rank ",n_full);
2022 for (
unsigned int i=0; i<n.
size(); i++) {
2023 print(
" ",i,
" ",n[i]);
2025 print(
" large rank ",n_large);
2030 for (
unsigned int i=0; i<std::min(3l,n.
size()); i++) nlog[0]+=n[i];
2031 for (
unsigned int i=3; i<std::min(10l,n.
size()); i++) nlog[1]+=n[i];
2032 for (
unsigned int i=10; i<std::min(30l,n.
size()); i++) nlog[2]+=n[i];
2033 for (
unsigned int i=30; i<std::min(100l,n.
size()); i++) nlog[3]+=n[i];
2034 for (
unsigned int i=100; i<std::min(300l,n.
size()); i++) nlog[4]+=n[i];
2035 for (
unsigned int i=300; i<std::min(1000l,n.
size()); i++) nlog[5]+=n[i];
2037 std::vector<std::string> slog={
"3",
"10",
"30",
"100",
"300",
"1000"};
2038 for (
unsigned int i=0; i<nlog.
size(); i++) {
2039 print(
" < ",slog[i],
" ",nlog[i]);
2041 print(
" large rank ",n_large);
2046 template <
typename T, std::
size_t NDIM>
2049 const int k = cdata.k;
2056 for (
int p=0;
p<
k; ++
p)
2059 else if (
NDIM == 2) {
2060 for (
int p=0;
p<
k; ++
p)
2061 for (
int q=0;
q<
k; ++
q)
2064 else if (
NDIM == 3) {
2065 for (
int p=0;
p<
k; ++
p)
2066 for (
int q=0;
q<
k; ++
q)
2067 for (
int r=0; r<
k; ++r)
2068 sum +=
c(
p,
q,r)*px[0][
p]*px[1][
q]*px[2][r];
2070 else if (
NDIM == 4) {
2071 for (
int p=0;
p<
k; ++
p)
2072 for (
int q=0;
q<
k; ++
q)
2073 for (
int r=0; r<
k; ++r)
2074 for (
int s=0; s<
k; ++s)
2075 sum +=
c(
p,
q,r,s)*px[0][
p]*px[1][
q]*px[2][r]*px[3][s];
2077 else if (
NDIM == 5) {
2078 for (
int p=0;
p<
k; ++
p)
2079 for (
int q=0;
q<
k; ++
q)
2080 for (
int r=0; r<
k; ++r)
2081 for (
int s=0; s<
k; ++s)
2082 for (
int t=0; t<
k; ++t)
2083 sum +=
c(
p,
q,r,s,t)*px[0][
p]*px[1][
q]*px[2][r]*px[3][s]*px[4][t];
2085 else if (
NDIM == 6) {
2086 for (
int p=0;
p<
k; ++
p)
2087 for (
int q=0;
q<
k; ++
q)
2088 for (
int r=0; r<
k; ++r)
2089 for (
int s=0; s<
k; ++s)
2090 for (
int t=0; t<
k; ++t)
2091 for (
int u=0;
u<
k; ++
u)
2092 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];
2100 template <
typename T, std::
size_t NDIM>
2112 if (it == coeffs.end()) {
2114 it = coeffs.find(key).get();
2116 nodeT& node = it->second;
2127 if (!
d.has_data())
d =
coeffT(cdata.v2k,targs);
2128 if (accumulate_NS and (key.
level() > 0))
d(cdata.s0) += s;
2129 if (
d.dim(0)==2*get_k()) {
2134 const keyT& child = kit.key();
2138 woT::task(coeffs.owner(child), &implT::reconstruct_op, child, ss, accumulate_NS);
2148 if (s.has_no_data()) ss=
coeffT(cdata.vk,targs);
2154 template <
typename T, std::
size_t NDIM>
2157 std::vector<long> npt(
NDIM,qx.
dim(0));
2163 template <
typename T, std::
size_t NDIM>
2166 std::vector<long> npt(
NDIM,qx.
dim(0));
2172 template <
typename T, std::
size_t NDIM>
2181 const double h = std::pow(0.5,
double(n));
2183 const int npt = qx.
dim(0);
2190 for (std::size_t i = 0; i <
NDIM; i++) {
2191 c1[i] = cell(i,0) +
h*cell_width[i]*(l[i] + qx((
long)0));
2192 c2[i] = cell(i,0) +
h*cell_width[i]*(l[i] + qx(npt-1));
2194 if (
f.screened(c1, c2)) {
2200 bool vectorized =
f.supports_vectorized();
2202 T* fvptr = fval.
ptr();
2204 double* x1 =
new double[npt];
2206 for (
int i=0; i<npt; ++i, ++idx) {
2207 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2211 f(xvals, fvptr, npt);
2214 else if (
NDIM == 2) {
2215 double* x1 =
new double[npt*npt];
2216 double* x2 =
new double[npt*npt];
2218 for (
int i=0; i<npt; ++i) {
2219 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2220 for (
int j=0; j<npt; ++j, ++idx) {
2221 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2227 f(xvals, fvptr, npt*npt);
2231 else if (
NDIM == 3) {
2232 double* x1 =
new double[npt*npt*npt];
2233 double* x2 =
new double[npt*npt*npt];
2234 double* x3 =
new double[npt*npt*npt];
2236 for (
int i=0; i<npt; ++i) {
2237 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2238 for (
int j=0; j<npt; ++j) {
2239 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2240 for (
int k=0;
k<npt; ++
k, ++idx) {
2241 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2249 f(xvals, fvptr, npt*npt*npt);
2254 else if (
NDIM == 4) {
2255 double* x1 =
new double[npt*npt*npt*npt];
2256 double* x2 =
new double[npt*npt*npt*npt];
2257 double* x3 =
new double[npt*npt*npt*npt];
2258 double* x4 =
new double[npt*npt*npt*npt];
2260 for (
int i=0; i<npt; ++i) {
2261 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2262 for (
int j=0; j<npt; ++j) {
2263 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2264 for (
int k=0;
k<npt; ++
k) {
2265 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2266 for (
int m=0;
m<npt; ++
m, ++idx) {
2267 c[3] = cell(3,0) +
h*cell_width[3]*(l[3] + qx(
m));
2277 f(xvals, fvptr, npt*npt*npt*npt);
2283 else if (
NDIM == 5) {
2284 double* x1 =
new double[npt*npt*npt*npt*npt];
2285 double* x2 =
new double[npt*npt*npt*npt*npt];
2286 double* x3 =
new double[npt*npt*npt*npt*npt];
2287 double* x4 =
new double[npt*npt*npt*npt*npt];
2288 double* x5 =
new double[npt*npt*npt*npt*npt];
2290 for (
int i=0; i<npt; ++i) {
2291 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2292 for (
int j=0; j<npt; ++j) {
2293 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2294 for (
int k=0;
k<npt; ++
k) {
2295 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2296 for (
int m=0;
m<npt; ++
m) {
2297 c[3] = cell(3,0) +
h*cell_width[3]*(l[3] + qx(
m));
2298 for (
int n=0; n<npt; ++n, ++idx) {
2299 c[4] = cell(4,0) +
h*cell_width[4]*(l[4] + qx(n));
2311 f(xvals, fvptr, npt*npt*npt*npt*npt);
2318 else if (
NDIM == 6) {
2319 double* x1 =
new double[npt*npt*npt*npt*npt*npt];
2320 double* x2 =
new double[npt*npt*npt*npt*npt*npt];
2321 double* x3 =
new double[npt*npt*npt*npt*npt*npt];
2322 double* x4 =
new double[npt*npt*npt*npt*npt*npt];
2323 double* x5 =
new double[npt*npt*npt*npt*npt*npt];
2324 double* x6 =
new double[npt*npt*npt*npt*npt*npt];
2326 for (
int i=0; i<npt; ++i) {
2327 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2328 for (
int j=0; j<npt; ++j) {
2329 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2330 for (
int k=0;
k<npt; ++
k) {
2331 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2332 for (
int m=0;
m<npt; ++
m) {
2333 c[3] = cell(3,0) +
h*cell_width[3]*(l[3] + qx(
m));
2334 for (
int n=0; n<npt; ++n) {
2335 c[4] = cell(4,0) +
h*cell_width[4]*(l[4] + qx(n));
2336 for (
int p=0;
p<npt; ++
p, ++idx) {
2337 c[5] = cell(5,0) +
h*cell_width[5]*(l[5] + qx(
p));
2351 f(xvals, fvptr, npt*npt*npt*npt*npt*npt);
2365 for (
int i=0; i<npt; ++i) {
2366 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2371 else if (
NDIM == 2) {
2372 for (
int i=0; i<npt; ++i) {
2373 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2374 for (
int j=0; j<npt; ++j) {
2375 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2381 else if (
NDIM == 3) {
2382 for (
int i=0; i<npt; ++i) {
2383 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2384 for (
int j=0; j<npt; ++j) {
2385 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2386 for (
int k=0;
k<npt; ++
k) {
2387 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2394 else if (
NDIM == 4) {
2395 for (
int i=0; i<npt; ++i) {
2396 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2397 for (
int j=0; j<npt; ++j) {
2398 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2399 for (
int k=0;
k<npt; ++
k) {
2400 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2401 for (
int m=0;
m<npt; ++
m) {
2402 c[3] = cell(3,0) +
h*cell_width[3]*(l[3] + qx(
m));
2403 fval(i,j,
k,
m) =
f(
c);
2410 else if (
NDIM == 5) {
2411 for (
int i=0; i<npt; ++i) {
2412 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2413 for (
int j=0; j<npt; ++j) {
2414 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2415 for (
int k=0;
k<npt; ++
k) {
2416 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2417 for (
int m=0;
m<npt; ++
m) {
2418 c[3] = cell(3,0) +
h*cell_width[3]*(l[3] + qx(
m));
2419 for (
int n=0; n<npt; ++n) {
2420 c[4] = cell(4,0) +
h*cell_width[4]*(l[4] + qx(n));
2421 fval(i,j,
k,
m,n) =
f(
c);
2429 else if (
NDIM == 6) {
2430 for (
int i=0; i<npt; ++i) {
2431 c[0] = cell(0,0) +
h*cell_width[0]*(l[0] + qx(i));
2432 for (
int j=0; j<npt; ++j) {
2433 c[1] = cell(1,0) +
h*cell_width[1]*(l[1] + qx(j));
2434 for (
int k=0;
k<npt; ++
k) {
2435 c[2] = cell(2,0) +
h*cell_width[2]*(l[2] + qx(
k));
2436 for (
int m=0;
m<npt; ++
m) {
2437 c[3] = cell(3,0) +
h*cell_width[3]*(l[3] + qx(
m));
2438 for (
int n=0; n<npt; ++n) {
2439 c[4] = cell(4,0) +
h*cell_width[4]*(l[4] + qx(n));
2440 for (
int p=0;
p<npt; ++
p) {
2441 c[5] = cell(5,0) +
h*cell_width[5]*(l[5] + qx(
p));
2442 fval(i,j,
k,
m,n,
p) =
f(
c);
2457 template <
typename T, std::
size_t NDIM>
2463 template <
typename T, std::
size_t NDIM>
2475 template <
typename T, std::
size_t NDIM>
2480 if (do_refine && key.
level() < max_refine_level) {
2483 std::vector<Vector<double,NDIM> > newspecialpts;
2484 if (key.
level() < special_level && specialpts.size() > 0) {
2486 const auto bperiodic = bc.is_periodic();
2487 for (
unsigned int i = 0; i < specialpts.size(); ++i) {
2492 newspecialpts.push_back(specialpts[i]);
2506 const keyT& child = it.key();
2507 r(child_patch(child)) =
project(child);
2511 if (truncate_on_project)
s0 =
copy(
d(cdata.s0));
2518 if (newspecialpts.size() > 0 || dnorm >=truncate_tol(
thresh,key.
level())) {
2521 const keyT& child = it.key();
2524 p = world.random_proc();
2527 p = coeffs.owner(child);
2530 woT::task(
p, &implT::project_refine_op, child, do_refine, newspecialpts);
2534 if (truncate_on_project) {
2536 coeffs.replace(key,
nodeT(s,
false));
2541 const keyT& child = it.key();
2543 coeffs.replace(child,
nodeT(s,
false));
2553 template <
typename T, std::
size_t NDIM>
2555 std::vector<long> v0(
NDIM,0
L);
2556 std::vector<long> v1(
NDIM,1L);
2559 if (is_compressed()) {
2560 if (world.rank() == coeffs.owner(cdata.key0)) {
2563 nodeT& node = it->second;
2574 for (
typename dcT::iterator it=coeffs.begin(); it!=coeffs.end(); ++it) {
2575 Level n = it->first.level();
2576 nodeT& node = it->second;
2584 node.
coeff()(s) += tt;
2591 if (fence) world.gop.fence();
2594 template <
typename T, std::
size_t NDIM>
2597 if (is_compressed()) initial_level = std::max(initial_level,1);
2598 if (coeffs.is_local(key)) {
2599 if (is_compressed()) {
2600 if (key.
level() == initial_level) {
2604 coeffs.replace(key,
nodeT(
coeffT(cdata.v2k,targs),
true));
2608 if (key.
level()<initial_level) {
2612 coeffs.replace(key,
nodeT(
coeffT(cdata.vk,targs),
false));
2616 if (key.
level() < initial_level) {
2618 insert_zero_down_to_initial_level(kit.key());
2625 template <
typename T, std::
size_t NDIM>
2629 if (it == coeffs.end()) {
2633 coeffs.replace(key,
nodeT());
2634 it = coeffs.find(key).get();
2636 nodeT& node = it->second;
2638 std::vector< Future<bool> >
v = future_vector_factory<bool>(1<<
NDIM);
2643 return woT::task(world.rank(),&implT::truncate_op, key, tol,
v);
2652 if (dnorm < truncate_tol(tol,key)) {
2661 template <
typename T, std::
size_t NDIM>
2665 for (
int i=0; i<(1<<
NDIM); ++i)
if (
v[i].get())
return true;
2666 nodeT& node = coeffs.find(key).get()->second;
2673 if (key.
level() > 1) {
2675 if (dnorm < truncate_tol(tol,key)) {
2680 coeffs.erase(kit.key());
2689 template <
typename T, std::
size_t NDIM>
2691 if (world.rank() == 0) do_print_tree(cdata.key0, os, maxlevel);
2693 if (world.rank() == 0) os.flush();
2698 template <
typename T, std::
size_t NDIM>
2701 if (it == coeffs.end()) {
2703 for (
int i=0; i<key.
level(); ++i) os <<
" ";
2704 os << key <<
" missing --> " << coeffs.owner(key) <<
"\n";
2707 const nodeT& node = it->second;
2708 for (
int i=0; i<key.
level(); ++i) os <<
" ";
2709 os << key <<
" " << node <<
" --> " << coeffs.owner(key) <<
"\n";
2712 do_print_tree(kit.key(),os,maxlevel);
2718 template <
typename T, std::
size_t NDIM>
2720 std::multimap<Level, std::tuple<tranT, std::string>>
data;
2721 if (world.rank() == 0) do_print_tree_json(cdata.key0,
data, maxlevel);
2723 if (world.rank() == 0) {
2724 for (
Level level = 0; level != maxlevel; ++level) {
2725 if (
data.count(level) == 0)
2730 os <<
"\"" << level <<
"\":{";
2731 os <<
"\"level\": " << level <<
",";
2732 os <<
"\"nodes\":{";
2733 auto range =
data.equal_range(level);
2734 for (
auto it = range.first; it != range.second; ++it) {
2735 os <<
"\"" << std::get<0>(it->second) <<
"\":"
2736 << std::get<1>(it->second);
2737 if (std::next(it) != range.second)
2749 template <
typename T, std::
size_t NDIM>
2752 if (it == coeffs.end()) {
2756 const nodeT& node = it->second;
2757 std::ostringstream oss;
2760 oss <<
",\"owner\": " << coeffs.owner(key) <<
"}";
2761 auto node_json_str = oss.str();
2765 do_print_tree_json(kit.key(),
data, maxlevel);
2771 template <
typename T, std::
size_t NDIM>
2774 if (world.rank() == 0) do_print_tree_graphviz(cdata.key0, os, maxlevel);
2776 if (world.rank() == 0) os.flush();
2780 template <
typename T, std::
size_t NDIM>
2784 static int64_t value(
const keyT& key) {
2786 for (int64_t j = 0; j <= key.
level()-1; ++j) {
2787 result += (1 << j*
NDIM);
2795 if (it != coeffs.end()) {
2796 const nodeT& node = it->second;
2799 os << uniqhash::value(key) <<
" -> " << uniqhash::value(kit.key()) <<
"\n";
2800 do_print_tree_graphviz(kit.key(),os,maxlevel);
2806 template <
typename T, std::
size_t NDIM>
2810 if (not functor)
MADNESS_EXCEPTION(
"FunctionImpl: project: confusion about function?",0);
2813 if (functor->provides_coeff())
return functor->coeff(key).full_tensor_copy();
2818 tensorT workq(cdata.vq,
false);
2827 template <
typename T, std::
size_t NDIM>
2829 if (coeffs.probe(key)) {
2830 return Future<double>(coeffs.find(key).get()->second.get_norm_tree());
2834 return woT::task(coeffs.owner(parent), &implT::get_norm_tree_recursive, parent,
TaskAttributes::hipri());
2838 template <
typename T, std::
size_t NDIM>
2842 if (coeffs.probe(key)) {
2843 const nodeT& node = coeffs.find(key).get()->second;
2847 result.
set(std::pair<keyT,coeffT>(key,node.
coeff()));
2851 result.
set(std::pair<keyT,coeffT>(key,
coeffT()));
2858 if (coeffs.is_local(parent))
2866 template <
typename T, std::
size_t NDIM>
2870 if (coeffs.probe(key)) {
2871 const nodeT& node = coeffs.find(key).get()->second;
2874 result.
set(std::pair<keyT,coeffT>(key,node.
coeff()));
2888 template <
typename T, std::
size_t NDIM>
2910 nodeT& node = it->second;
2916 for (std::size_t i=0; i<
NDIM; ++i) {
2917 double xi = x[i]*2.0;
2919 if (li == 2) li = 1;
2931 template <
typename T, std::
size_t NDIM>
2938 while (key.
level() <= maxlevel) {
2939 if (coeffs.owner(key) ==
me) {
2942 if (it != coeffs.end()) {
2943 nodeT& node = it->second;
2949 for (std::size_t i=0; i<
NDIM; ++i) {
2950 double xi = x[i]*2.0;
2952 if (li == 2) li = 1;
2958 return std::pair<bool,T>(
false,0.0);
2961 template <
typename T, std::
size_t NDIM>
2983 nodeT& node = it->second;
2989 for (std::size_t i=0; i<
NDIM; ++i) {
2990 double xi = x[i]*2.0;
2992 if (li == 2) li = 1;
3003 template <
typename T, std::
size_t NDIM>
3025 nodeT& node = it->second;
3031 for (std::size_t i=0; i<
NDIM; ++i) {
3032 double xi = x[i]*2.0;
3034 if (li == 2) li = 1;
3046 template <
typename T, std::
size_t NDIM>
3057 template <
typename T, std::
size_t NDIM>
3060 coeffT shalf=t(cdata.sh);
3063 sfull(cdata.sh)-=shalf;
3067 template <
typename T, std::
size_t NDIM>
3073 if (t.
rank()==0)
return;
3075 for (
long i=0; i<t.
rank(); ++i) {
3078 tnorm(
c, &lo1, &hi1);
3086 template <
typename A,
typename B>
3093 template <
typename T, std::
size_t NDIM>
3110 template <
typename T, std::
size_t NDIM>
3118 template <
typename T, std::
size_t NDIM>
3124 template <
typename T, std::
size_t NDIM>
3132template <
typename T, std::
size_t NDIM>
3138 template <
typename T, std::
size_t NDIM>
3144 template <
typename T, std::
size_t NDIM>
3149 template <
typename T, std::
size_t NDIM>
3154 template <
typename T, std::
size_t NDIM>
3159 for (
int mu=0;
mu<cdata.npt; ++
mu) {
3160 double xmu =
scale*(cdata.quad_x(
mu)+lc) - lp;
3163 for (
int i=0; i<
k; ++i) phi(i,
mu) =
p[i];
3168 template <
typename T, std::
size_t NDIM>
3178 coeffT result = fcube_for_mul<T>(child, parent, s);
3180 result =
transform(result,cdata.quad_phiw);
3186 template <
typename T, std::
size_t NDIM>
3189 std::vector<long> v0(
NDIM,0);
3191 if (is_compressed()) {
3192 if (world.rank() == coeffs.owner(cdata.key0)) {
3194 if (it != coeffs.end()) {
3195 const nodeT& node = it->second;
3202 const keyT& key = it->first;
3203 const nodeT& node = it->second;
3220 }
else if (l >= two2n) {
3224 }
while (l >= two2n);
3233 return l >= 0 && l < two2n;
3236 template <
typename T, std::
size_t NDIM>
3245 return keyT::invalid();
3248 return keyT(key.
level(),l);
3251 template <
typename T, std::
size_t NDIM>
3259 return keyT::invalid();
3262 return keyT(key.
level(), l);
3265 template <
typename T, std::
size_t NDIM>
3266 Future< std::pair< Key<NDIM>, GenTensor<T> > >
3269 typedef std::pair< Key<NDIM>,
coeffT > argT;
3279 template <
typename T, std::
size_t NDIM>
3281 bool nonstandard1,
bool keepleaves,
bool redundant1) {
3282 if (!coeffs.probe(key))
print(
"missing node",key);
3286 nodeT& node = coeffs.find(key).get()->second;
3289 if (node.has_children()) {
3290 std::vector< Future<std::pair<coeffT,double> > >
v = future_vector_factory<std::pair<coeffT,double> >(1<<
NDIM);
3295 v[i] = woT::task(coeffs.owner(kit.key()), &implT::compress_spawn, kit.key(),
3298 if (redundant1)
return woT::task(world.rank(),&implT::make_redundant_op, key,
v);
3299 return woT::task(world.rank(),&implT::compress_op, key,
v, nonstandard1);
3306 if (key.
level()==0) {
3310 double snorm=node.coeff().normf();
3311 node.set_dnorm(0.0);
3312 node.set_snorm(snorm);
3313 node.set_norm_tree(snorm);
3318 coeffT sdcoeff(cdata.v2k,this->get_tensor_type());
3319 sdcoeff(cdata.s0)+=node.coeff();
3320 node.coeff()=sdcoeff;
3321 double snorm=node.coeff().normf();
3322 node.set_dnorm(0.0);
3323 node.set_snorm(snorm);
3324 node.set_norm_tree(snorm);
3330 if (not keepleaves) node.clear_coeff();
3332 auto snorm=(keepleaves) ? node.coeff().normf() : 0.0;
3333 node.set_norm_tree(snorm);
3334 node.set_snorm(snorm);
3335 node.set_dnorm(0.0);
3342 template <
typename T, std::
size_t NDIM>
3345 const coordT& plotlo,
const coordT& plothi,
const std::vector<long>& npt,
3346 bool eval_refine)
const {
3351 for (std::size_t i=0; i<
NDIM; ++i) {
3353 h[i] = (plothi[i]-plotlo[i])/(npt[i]-1);
3363 const double twon =
pow(2.0,
double(n));
3364 const tensorT& coeff = coeffs.find(key).get()->second.coeff().full_tensor_copy();
3371 double fac =
pow(0.5,
double(key.
level()));
3373 for (std::size_t
d=0;
d<
NDIM; ++
d) {
3376 boxhi[
d] = boxlo[
d]+fac;
3378 if (boxlo[
d] > plothi[
d] || boxhi[
d] < plotlo[
d]) {
3380 npttotal = boxnpt[
d] = 0;
3384 else if (npt[
d] == 1) {
3386 boxlo[
d] = boxhi[
d] = plotlo[
d];
3391 boxlo[
d] = std::max(boxlo[
d],plotlo[
d]);
3392 boxhi[
d] = std::min(boxhi[
d],plothi[
d]);
3395 double xlo = long((boxlo[
d]-plotlo[
d])/
h[
d])*
h[
d] + plotlo[
d];
3396 if (xlo < boxlo[
d]) xlo +=
h[
d];
3398 double xhi = long((boxhi[
d]-plotlo[
d])/
h[
d])*
h[
d] + plotlo[
d];
3399 if (xhi > boxhi[
d]) xhi -=
h[
d];
3402 boxnpt[
d] = long(round((boxhi[
d] - boxlo[
d])/
h[
d])) + 1;
3404 npttotal *= boxnpt[
d];
3408 for (IndexIterator it(boxnpt); it; ++it) {
3409 for (std::size_t
d=0;
d<
NDIM; ++
d) {
3410 double xd = boxlo[
d] + it[
d]*
h[
d];
3411 x[
d] = twon*xd - l[
d];
3414 ind[
d] = long(round((xd-plotlo[
d])/
h[
d]));
3425 T tmp = eval_cube(n, x, coeff);
3435 template <
typename T, std::
size_t NDIM>
3438 const std::vector<long>& npt,
3439 const bool eval_refine)
const {
3446 const keyT& key = it->first;
3447 const nodeT& node = it->second;
3449 woT::task(world.rank(), &implT::plot_cube_kernel,
3456 world.taskq.fence();
3457 world.gop.sum(r.
ptr(), r.
size());
3464 fprintf(
f,
"%.6e\n",t);
3468 fprintf(
f,
"%.6e %.6e\n", t.real(), t.imag());
3471 template <
typename T, std::
size_t NDIM>
3475 const std::vector<long>& npt,
3479 const char* element[6] = {
"lines",
"quads",
"cubes",
"cubes4D",
"cubes5D",
"cubes6D"};
3484 if (world.
rank() == 0) {
3488 fprintf(
f,
"object 1 class gridpositions counts ");
3489 for (std::size_t
d=0;
d<
NDIM; ++
d) fprintf(
f,
" %ld",npt[
d]);
3492 fprintf(
f,
"origin ");
3493 for (std::size_t
d=0;
d<
NDIM; ++
d) fprintf(
f,
" %.6e", cell(
d,0));
3496 for (std::size_t
d=0;
d<
NDIM; ++
d) {
3497 fprintf(
f,
"delta ");
3498 for (std::size_t
c=0;
c<
d; ++
c) fprintf(
f,
" 0");
3500 if (npt[
d]>1)
h = (cell(
d,1)-cell(
d,0))/(npt[
d]-1);
3501 fprintf(
f,
" %.6e",
h);
3502 for (std::size_t
c=
d+1;
c<
NDIM; ++
c) fprintf(
f,
" 0");
3507 fprintf(
f,
"object 2 class gridconnections counts ");
3508 for (std::size_t
d=0;
d<
NDIM; ++
d) fprintf(
f,
" %ld",npt[
d]);
3510 fprintf(
f,
"attribute \"element type\" string \"%s\"\n", element[
NDIM-1]);
3511 fprintf(
f,
"attribute \"ref\" string \"positions\"\n");
3515 for (std::size_t
d=0;
d<
NDIM; ++
d) npoint *= npt[
d];
3516 const char* iscomplex =
"";
3518 const char* isbinary =
"";
3519 if (binary) isbinary =
"binary";
3520 fprintf(
f,
"object 3 class array type double %s rank 0 items %d %s data follows\n",
3521 iscomplex, npoint, isbinary);
3527 if (world.
rank() == 0) {
3531 fwrite((
void *) r.
ptr(),
sizeof(
T), r.
size(),
f);
3542 fprintf(
f,
"object \"%s\" class field\n",
filename);
3543 fprintf(
f,
"component \"positions\" value 1\n");
3544 fprintf(
f,
"component \"connections\" value 2\n");
3545 fprintf(
f,
"component \"data\" value 3\n");
3546 fprintf(
f,
"\nend\n");
3552 template <std::
size_t NDIM>
3558 max_refine_level = 30;
3563 truncate_on_project =
true;
3564 apply_randomize =
false;
3565 project_randomize =
false;
3568 cell = make_default_cell();
3569 recompute_cell_info();
3570 set_default_pmap(world);
3573 template <std::
size_t NDIM>
3581 template <std::
size_t NDIM>
3583 std::cout <<
"Function Defaults:" << std::endl;
3584 std::cout <<
" Dimension " <<
": " <<
NDIM << std::endl;
3585 std::cout <<
" k" <<
": " <<
k << std::endl;
3586 std::cout <<
" thresh" <<
": " <<
thresh << std::endl;
3587 std::cout <<
" initial_level" <<
": " << initial_level << std::endl;
3588 std::cout <<
" special_level" <<
": " << special_level << std::endl;
3589 std::cout <<
" max_refine_level" <<
": " << max_refine_level << std::endl;
3590 std::cout <<
" truncate_mode" <<
": " <<
truncate_mode << std::endl;
3591 std::cout <<
" refine" <<
": " <<
refine << std::endl;
3592 std::cout <<
" autorefine" <<
": " << autorefine << std::endl;
3593 std::cout <<
" debug" <<
": " <<
debug << std::endl;
3594 std::cout <<
" truncate_on_project" <<
": " << truncate_on_project << std::endl;
3595 std::cout <<
" apply_randomize" <<
": " << apply_randomize << std::endl;
3596 std::cout <<
" project_randomize" <<
": " << project_randomize << std::endl;
3597 std::cout <<
" bc" <<
": " << get_bc() << std::endl;
3598 std::cout <<
" tt" <<
": " << tt << std::endl;
3599 std::cout <<
" cell" <<
": " << cell << std::endl;
3602 template <
typename T, std::
size_t NDIM>
3603 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:266
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:88
void _init_twoscale()
Private. Initialize the twoscale coefficients.
Definition mraimpl.h:67
FunctionDefaults holds default paramaters as static class members.
Definition funcdefaults.h:100
static Tensor< double > make_default_cell_width()
Definition funcdefaults.h:135
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)
Sets the default process map.
Definition mraimpl.h:3574
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:369
static const BoundaryConditions< NDIM > & get_bc()
Returns the default boundary conditions.
Definition funcdefaults.h:310
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:126
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:3582
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:379
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:347
static void set_defaults(World &world)
Definition mraimpl.h:3553
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:5441
FunctionImpl holds all Function state to facilitate shallow copy semantics.
Definition funcimpl.h:945
bool is_nonstandard() const
Definition mraimpl.h:252
T eval_cube(Level n, coordT &x, const tensorT &c) const
Definition mraimpl.h:2047
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:2781
void change_tensor_type1(const TensorArgs &targs, bool fence)
change the tensor type of the coefficients in the FunctionNode
Definition mraimpl.h:1078
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:2962
void scale_inplace(const T q, bool fence)
In-place scale by a constant.
Definition mraimpl.h:3133
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:208
bool is_redundant() const
Returns true if the function is redundant.
Definition mraimpl.h:247
FunctionNode< Q, NDIM > nodeT
Type of node.
Definition funcimpl.h:955
void print_size(const std::string name) const
print tree size and size
Definition mraimpl.h:1963
void print_info() const
Prints summary of data distribution.
Definition mraimpl.h:812
void abs_inplace(bool fence)
Definition mraimpl.h:3145
void print_timer() const
Definition mraimpl.h:336
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:3004
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:563
std::size_t nCoeff() const
Returns the number of coefficients in the function ... collective global sum.
Definition mraimpl.h:1949
keyT neighbor_in_volume(const keyT &key, const keyT &disp) const
Returns key of general neighbor that resides in-volume.
Definition mraimpl.h:3252
void compress(const TreeState newstate, bool fence)
compress the wave function
Definition mraimpl.h:1527
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:1695
Future< bool > truncate_spawn(const keyT &key, double tol)
Returns true if after truncation this node has coefficients.
Definition mraimpl.h:2626
Future< double > norm_tree_spawn(const keyT &key)
Definition mraimpl.h:1597
std::vector< keyT > local_leaf_keys() const
return the keys of the local leaf boxes
Definition mraimpl.h:537
void do_print_tree(const keyT &key, std::ostream &os, Level maxlevel) const
Functor for the do_print_tree method.
Definition mraimpl.h:2699
void unset_functor()
Definition mraimpl.h:291
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:478
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:944
void set_functor(const std::shared_ptr< FunctionFunctorInterface< T, NDIM > > functor1)
Definition mraimpl.h:272
const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > & get_pmap() const
Definition mraimpl.h:192
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:2839
double get_thresh() const
Definition mraimpl.h:307
void trickle_down(bool fence)
sum all the contributions from all scales after applying an operator in mod-NS form
Definition mraimpl.h:1333
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:1755
void set_autorefine(bool value)
Definition mraimpl.h:316
tensorT filter(const tensorT &s) const
Transform sum coefficients at level n to sums+differences at level n-1.
Definition mraimpl.h:1131
void chop_at_level(const int n, const bool fence=true)
remove all nodes with level higher than n
Definition mraimpl.h:1094
void print_tree_json(std::ostream &os=std::cout, Level maxlevel=10000) const
Definition mraimpl.h:2719
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:686
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:1036
bool is_compressed() const
Returns true if the function is compressed.
Definition mraimpl.h:235
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:1045
void print_stats() const
print the number of configurations per node
Definition mraimpl.h:1991
void broaden(const array_of_bools< NDIM > &is_periodic, bool fence)
Definition mraimpl.h:1282
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:1644
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:2464
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:902
std::vector< Slice > child_patch(const keyT &child) const
Returns patch referring to coeffs of child in parent box.
Definition mraimpl.h:675
void print_tree_graphviz(std::ostream &os=std::cout, Level maxlevel=10000) const
Definition mraimpl.h:2772
std::size_t min_nodes() const
Returns the min number of nodes on a processor.
Definition mraimpl.h:1900
void make_redundant(const bool fence)
convert this to redundant, i.e. have sum coefficients on all levels
Definition mraimpl.h:1555
std::size_t max_nodes() const
Returns the max number of nodes on a processor.
Definition mraimpl.h:1891
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:1210
void flo_unary_op_node_inplace(const opT &op, bool fence)
Definition funcimpl.h:2142
std::size_t size_local() const
Returns the number of coefficients in the function for each rank.
Definition mraimpl.h:1918
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:3343
T trace_local() const
Returns int(f(x),x) in local volume.
Definition mraimpl.h:3187
void print_grid(const std::string filename) const
Definition mraimpl.h:521
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:3280
bool get_autorefine() const
Definition mraimpl.h:313
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:3155
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:3267
TensorType get_tensor_type() const
Definition mraimpl.h:298
void remove_leaf_coefficients(const bool fence)
Definition mraimpl.h:1549
void insert_zero_down_to_initial_level(const keyT &key)
Initialize nodes to zero function at initial_level of refinement.
Definition mraimpl.h:2595
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:913
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:2750
void finalize_sum()
after summing up we need to do some cleanup;
Definition mraimpl.h:1848
dcT coeffs
The coefficients.
Definition funcimpl.h:988
bool exists_and_is_leaf(const keyT &key) const
Definition mraimpl.h:1254
void verify_tree() const
Verify tree is properly constructed ... global synchronization involved.
Definition mraimpl.h:109
void set_tensor_args(const TensorArgs &t)
Definition mraimpl.h:304
Range< typename dcT::const_iterator > rangeT
Definition funcimpl.h:5570
std::size_t real_size() const
Returns the number of coefficients in the function ... collective global sum.
Definition mraimpl.h:1936
bool exists_and_has_children(const keyT &key) const
Definition mraimpl.h:1249
void sum_down_spawn(const keyT &key, const coeffT &s)
is this the same as trickle_down() ?
Definition mraimpl.h:855
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:3237
void norm_tree(bool fence)
compute for each FunctionNode the norm of the function inside that node
Definition mraimpl.h:1574
bool has_leaves() const
Definition mraimpl.h:267
void reconstruct_op(const keyT &key, const coeffT &s, const bool accumulate_NS=true)
Definition mraimpl.h:2101
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:3170
void undo_redundant(const bool fence)
convert this from redundant to standard reconstructed form
Definition mraimpl.h:1565
GenTensor< Q > coeffT
Type of tensor used to hold coeffs.
Definition funcimpl.h:956
const keyT & key0() const
Returns cdata.key0.
Definition mraimpl.h:373
double finalize_apply()
after apply we need to do some cleanup;
Definition mraimpl.h:1805
const dcT & get_coeffs() const
Definition mraimpl.h:322
double norm2sq_local() const
Returns the square of the local norm ... no comms.
Definition mraimpl.h:1857
const FunctionCommonData< T, NDIM > & get_cdata() const
Definition mraimpl.h:328
void sum_down(bool fence)
After 1d push operator must sum coeffs down the tree to restore correct scaling function coefficients...
Definition mraimpl.h:894
bool noautorefine(const keyT &key, const tensorT &t) const
Always returns false (for when autorefine is not wanted)
Definition mraimpl.h:838
double truncate_tol(double tol, const keyT &key) const
Returns the truncation threshold according to truncate_method.
Definition mraimpl.h:628
bool autorefine_square_test(const keyT &key, const nodeT &t) const
Returns true if this block of coeffs needs autorefining.
Definition mraimpl.h:844
void erase(const Level &max_level)
truncate tree at a certain level
Definition mraimpl.h:718
void reconstruct(bool fence)
reconstruct this tree – respects fence
Definition mraimpl.h:1495
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:3560
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:992
static void tnorm(const tensorT &t, double *lo, double *hi)
Computes norm of low/high-order polyn. coeffs for autorefinement test.
Definition mraimpl.h:3047
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:2933
std::size_t max_depth() const
Returns the maximum depth of the tree ... collective ... global sum/broadcast.
Definition mraimpl.h:1883
std::size_t size() const
Returns the number of coefficients in the function ... collective global sum.
Definition mraimpl.h:1928
void reduce_rank(const double thresh, bool fence)
reduce the rank of the coefficients tensors
Definition mraimpl.h:1086
std::shared_ptr< FunctionFunctorInterface< T, NDIM > > get_functor()
Definition mraimpl.h:279
tensorT unfilter(const tensorT &s) const
Transform sums+differences at level n to sum coefficients at level n+1.
Definition mraimpl.h:1160
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:3436
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:1386
Future< coeffT > truncate_reconstructed_spawn(const keyT &key, const double tol)
truncate using a tree in reconstructed form
Definition mraimpl.h:1620
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:1055
void truncate(double tol, bool fence)
Truncate according to the threshold with optional global fence.
Definition mraimpl.h:357
bool is_reconstructed() const
Returns true if the function is compressed.
Definition mraimpl.h:241
double norm_tree_op(const keyT &key, const std::vector< Future< double > > &v)
Definition mraimpl.h:1582
void reset_timer()
Definition mraimpl.h:345
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:748
int get_k() const
Definition mraimpl.h:319
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:2889
bool truncate_op(const keyT &key, double tol, const std::vector< Future< bool > > &v)
Definition mraimpl.h:2662
void zero_norm_tree()
Definition mraimpl.h:1271
std::size_t max_local_depth() const
Returns the maximum local depth of the tree ... no communications.
Definition mraimpl.h:1869
tensorT project(const keyT &key) const
Definition mraimpl.h:2807
double check_symmetry_local() const
Returns some asymmetry measure ... no comms.
Definition mraimpl.h:734
Future< double > get_norm_tree_recursive(const keyT &key) const
Definition mraimpl.h:2828
Key< NDIM > keyT
Type of key.
Definition funcimpl.h:954
bool is_on_demand() const
Definition mraimpl.h:262
void accumulate_timer(const double time) const
Definition mraimpl.h:331
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:1344
void set_thresh(double value)
Definition mraimpl.h:310
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:402
void sock_it_to_me_too(const keyT &key, const RemoteReference< FutureImpl< std::pair< keyT, coeffT > > > &ref) const
Definition mraimpl.h:2867
void broaden_op(const keyT &key, const std::vector< Future< bool > > &v)
Definition mraimpl.h:1260
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:382
void print_tree(std::ostream &os=std::cout, Level maxlevel=10000) const
Definition mraimpl.h:2690
void project_refine_op(const keyT &key, bool do_refine, const std::vector< Vector< double, NDIM > > &specialpts)
Definition mraimpl.h:2476
std::size_t tree_size() const
Returns the size of the tree structure of the function ... collective global sum.
Definition mraimpl.h:1909
void add_scalar_inplace(T t, bool fence)
Adds a constant to the function. Local operation, optional fence.
Definition mraimpl.h:2554
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:1180
void abs_square_inplace(bool fence)
Definition mraimpl.h:3150
void put_in_box(ProcessID from, long nl, long ni) const
Definition mraimpl.h:803
TensorArgs get_tensor_args() const
Definition mraimpl.h:301
void average(const implT &rhs)
take the average of two functions, similar to: this=0.5*(this+rhs)
Definition mraimpl.h:1067
void diff(const DerivativeBase< T, NDIM > *D, const implT *f, bool fence)
Definition mraimpl.h:925
void square_inplace(bool fence)
Pointwise squaring of function with optional global fence.
Definition mraimpl.h:3139
void remove_internal_coefficients(const bool fence)
Definition mraimpl.h:1544
void compute_snorm_and_dnorm(bool fence=true)
compute norm of s and d coefficients for all nodes
Definition mraimpl.h:1110
void standard(bool fence)
Changes non-standard compressed form to standard compressed form.
Definition mraimpl.h:1792
bool is_nonstandard_with_leaves() const
Definition mraimpl.h:257
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
void set_is_leaf(bool flag)
Sets has_children attribute to value of !flag.
Definition funcimpl.h:280
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
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:373
T & get(bool dowork=true) &
Gets the value, waiting if necessary.
Definition future.h:574
remote_refT remote_ref(World &world) const
Returns a structure used to pass references to another process.
Definition future.h:675
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:508
RemoteReference< FutureImpl< T > > remote_refT
Definition future.h:398
Definition lowranktensor.h:59
GenTensor convert(const TensorArgs &targs) const
Definition gentensor.h:198
GenTensor full_tensor() const
Definition gentensor.h:200
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
size_t nCoeff() const
Definition gentensor.h:215
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:459
Key is the index for a node of the 2^NDIM-tree.
Definition key.h:68
Level level() const
Definition key.h:161
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:280
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:312
bool is_invalid() const
Checks if a key is invalid.
Definition key.h:111
Key parent(int generation=1) const
Returns the key of the parent.
Definition key.h:245
const Vector< Translation, NDIM > & translation() const
Definition key.h:166
A pmap that locates children on odd levels with their even level parents.
Definition funcimpl.h:105
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:450
static TaskAttributes generator()
Definition thread.h:442
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:1824
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:1798
bool has_data() const
Definition tensor.h:1886
iterator begin()
Returns an iterator to the beginning of the local data (no communication)
Definition worlddc.h:1070
implT::const_iterator const_iterator
Definition worlddc.h:872
iterator end()
Returns an iterator past the end of the local data (no communication)
Definition worlddc.h:1084
Future< iterator > futureT
Definition worlddc.h:875
implT::iterator iterator
Definition worlddc.h:871
implT::accessor accessor
Definition worlddc.h:873
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
Wrapper for an opaque pointer for serialization purposes.
Definition archive.h:850
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:34
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:2443
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
#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
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:182
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:3231
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:1125
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:664
static Key< NDIM > simpt2key(const Vector< T, NDIM > &pt, Level n)
Definition funcdefaults.h:434
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:425
void standard(World &world, std::vector< Function< T, NDIM > > &v, bool fence=true)
Generates standard form of a vector of functions.
Definition vmra.h:239
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:145
const std::vector< Function< T, NDIM > > & reconstruct(const std::vector< Function< T, NDIM > > &v)
reconstruct a vector of functions
Definition vmra.h:158
response_space transpose(response_space &f)
Definition basic_operators.cc:10
int64_t Translation
Definition key.h:56
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:3472
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:2306
static void verify_tree(World &world, const std::vector< Function< T, NDIM > > &v)
Definition SCF.cc:72
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:57
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
bool gauss_legendre(int n, double xlo, double xhi, double *x, double *w)
Definition legendre.cc:226
TreeState get_tree_state(const std::vector< Function< T, NDIM > > &v)
get tree state of a vector of functions
Definition vmra.h:134
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
Tensor< T > fcube(const Key< NDIM > &, T(*f)(const Vector< double, NDIM > &), const Tensor< double > &)
Definition mraimpl.h:2155
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:3463
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:192
NDIM & f
Definition mra.h:2451
void error(const char *msg)
Definition world.cc:139
NDIM const Function< R, NDIM > & g
Definition mra.h:2451
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:104
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:2407
static std::vector< double > ttt
Definition SCF.cc:105
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:2434
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:2037
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:3211
bool isnan(const std::complex< T > &v)
Definition mraimpl.h:52
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:3499
"put" this on g
Definition funcimpl.h:2567
change representation of nodes' coeffs to low rank, optional fence
Definition funcimpl.h:2600
check symmetry wrt particle exchange
Definition funcimpl.h:2273
compute the norm of the wavelet coefficients
Definition funcimpl.h:4396
Definition funcimpl.h:2627
Definition funcimpl.h:1389
mirror dimensions of this, write result on f
Definition funcimpl.h:2501
map this on f
Definition funcimpl.h:2421
mirror dimensions of this, write result on f
Definition funcimpl.h:2451
Definition funcimpl.h:5495
reduce the rank of the nodes, optional fence
Definition funcimpl.h:2247
Changes non-standard compressed form to standard compressed form.
Definition funcimpl.h:4617
given an NS tree resulting from a convolution, truncate leafs if appropriate
Definition funcimpl.h:2168
remove all coefficients of internal nodes
Definition funcimpl.h:2193
remove all coefficients of leaf nodes
Definition funcimpl.h:2210
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:3119
void operator()(const Key< NDIM > &key, Tensor< T > &t) const
Definition mraimpl.h:3120
void serialize(Archive &ar)
Definition mraimpl.h:3121
Definition mraimpl.h:3125
void operator()(const Key< NDIM > &key, Tensor< T > &t) const
Definition mraimpl.h:3126
void serialize(Archive &ar)
Definition mraimpl.h:3127
Definition mraimpl.h:3087
void operator()(const A &a, const B &b) const
Definition mraimpl.h:3088
void serialize(Archive &ar)
Definition mraimpl.h:3090
Definition mraimpl.h:3094
void operator()(const Key< NDIM > &key, FunctionNode< T, NDIM > &node) const
Definition mraimpl.h:3102
void serialize(Archive &ar)
Definition mraimpl.h:3105
T q
Definition mraimpl.h:3095
scaleinplace()
Definition mraimpl.h:3096
scaleinplace(T q)
Definition mraimpl.h:3098
void operator()(const Key< NDIM > &key, Tensor< T > &t) const
Definition mraimpl.h:3099
Definition mraimpl.h:3111
void serialize(Archive &ar)
Definition mraimpl.h:3115
void operator()(const Key< NDIM > &key, Tensor< T > &t) const
Definition mraimpl.h:3112
insert/replaces the coefficients into the function
Definition funcimpl.h:692
Definition lowrankfunction.h:332
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:72
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
void test()
Definition y.cc:696