MADNESS 0.10.1
mra.h
Go to the documentation of this file.
1/*
2 This file is part of MADNESS.
3
4 Copyright (C) 2007,2010 Oak Ridge National Laboratory
5
6 This program is free software; you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation; either version 2 of the License, or
9 (at your option) any later version.
10
11 This program is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with this program; if not, write to the Free Software
18 Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19
20 For more information please contact:
21
22 Robert J. Harrison
23 Oak Ridge National Laboratory
24 One Bethel Valley Road
25 P.O. Box 2008, MS-6367
26
27 email: harrisonrj@ornl.gov
28 tel: 865-241-3937
29 fax: 865-572-0680
30*/
31
32#ifndef MADNESS_MRA_MRA_H__INCLUDED
33#define MADNESS_MRA_MRA_H__INCLUDED
34
35/*!
36 \file mra/mra.h
37 \brief Main include file for MADNESS and defines \c Function interface
38
39 \addtogroup mra
40
41*/
42
43
45#include <madness/misc/misc.h>
47
48#define FUNCTION_INSTANTIATE_1
49#define FUNCTION_INSTANTIATE_2
50#define FUNCTION_INSTANTIATE_3
51#if !defined(HAVE_IBMBGP) || !defined(HAVE_IBMBGQ)
52#define FUNCTION_INSTANTIATE_4
53#define FUNCTION_INSTANTIATE_5
54#define FUNCTION_INSTANTIATE_6
55#endif
56
57static const bool VERIFY_TREE = false; //true
58
59
60namespace madness {
61 /// @brief initialize the internal state of the MADmra library
62 ///
63 /// Reads in (and broadcasts across \p world) the twoscale and autocorrelation coefficients,
64 /// Gauss-Legendre quadrature roots/weights, function defaults and operator displacement lists.
65 /// \warning By default this generates operator displacement lists (see Displacements) for up to 6-d free
66 /// and 3-d periodic boundary conditions. For optimal support for mixed boundary conditions
67 /// (periodic along some axes only) assign the desired boundary conditions
68 /// as default (e.g. `FunctionDefaults<3>::set_bc(BoundaryConditions<3>({BC_FREE, BC_FREE, BC_FREE, BC_FREE, BC_PERIODIC, BC_PERIODIC})`)
69 /// prior to calling this. This will make operator application with such boundary conditions
70 /// as efficient as possible, but will not allow the use of operators with
71 /// other boundary conditions that include periodic axes until Displacements::reset_periodic_axes is invoked.
72 /// By default efficiency is sacrificed for generality.
73 /// \param world broadcast data across this World
74 /// \param argc command-line parameter count
75 /// \param argv command-line parameters array
76 /// \param doprint if true, will log status to std::cout on rank 0 [default=false]
77 /// \param make_stdcout_nice_to_reals if true, will configure std::cout to print reals prettily, according to the MADNESS convention [default=true]
78 void startup(World& world, int argc, char** argv, bool doprint=false, bool make_stdcout_nice_to_reals = true);
79 std::string get_mra_data_dir();
80}
81
82#include <madness/mra/key.h>
85#include <madness/mra/indexit.h>
90#include <madness/mra/lbdeux.h>
92
93// some forward declarations
94namespace madness {
95
96 template<typename T, std::size_t NDIM>
97 class FunctionImpl;
98
99 template<typename T, std::size_t NDIM>
100 class Function;
101
102 template<typename T, std::size_t NDIM>
103 class FunctionNode;
104
105 template<typename T, std::size_t NDIM>
106 class FunctionFactory;
107
108 template<typename T, std::size_t NDIM>
109 class FunctionFunctorInterface;
110
111 template<typename T, std::size_t NDIM>
112 struct leaf_op;
113
114 template<typename T, std::size_t NDIM>
116
117 template<typename T, std::size_t NDIM>
118 struct hartree_leaf_op;
119
120 template<typename T, std::size_t NDIM, std::size_t LDIM, typename opT>
122
123 template<typename T, std::size_t NDIM, typename opT>
124 struct op_leaf_op;
125
126 template<typename T, std::size_t NDIM>
128
129}
130
131
132namespace madness {
133
134 /// \ingroup mra
135 /// \addtogroup function
136
137 /// A multiresolution adaptive numerical function
138 template <typename T, std::size_t NDIM>
140 // We make all of the content of Function and FunctionImpl
141 // public with the intent of avoiding the cumbersome forward
142 // and friend declarations. However, this open access should
143 // not be abused.
144
145 private:
146 std::shared_ptr< FunctionImpl<T,NDIM> > impl;
147
148 public:
149 bool impl_initialized()const{
150 if(impl==NULL) return false;
151 else return true;
152 }
156 typedef Vector<double,NDIM> coordT; ///< Type of vector holding coordinates
157 typedef T typeT;
158 static constexpr std::size_t dimT=NDIM;
159
160
161 /// Asserts that the function is initialized
162 inline void verify() const {
164 }
165
166 /// Returns true if the function is initialized
167 bool is_initialized() const {
168 return impl.get();
169 }
170
171 /// Default constructor makes uninitialized function. No communication.
172
173 /// An uninitialized function can only be assigned to. Any other operation will throw.
174 Function() : impl() {}
175
176
177 /// Constructor from FunctionFactory provides named parameter idiom. Possible non-blocking communication.
182
183
184 /// Copy constructor is \em shallow. No communication, works in either basis.
186 : impl(f.impl) {
187 }
188
189
190 /// Assignment is \em shallow. No communication, works in either basis.
193 if (this != &f) impl = f.impl;
194 return *this;
195 }
196
197 /// Destruction of any underlying implementation is deferred to next global fence.
199
200 /// implements swap algorithm
201 template <typename R, std::size_t MDIM>
203
204
205 /// Evaluates the function at a point in user coordinates. Possible non-blocking comm.
206
207 /// Only the invoking process will receive the result via the future
208 /// though other processes may be involved in the evaluation.
209 ///
210 /// Throws if function is not initialized.
211 Future<T> eval(const coordT& xuser) const {
213 const double eps=1e-15;
214 verify();
216 coordT xsim;
218 // If on the boundary, move the point just inside the
219 // volume so that the evaluation logic does not fail
220 for (std::size_t d=0; d<NDIM; ++d) {
221 if (xsim[d] < -eps) {
222 MADNESS_EXCEPTION("eval: coordinate lower-bound error in dimension", d);
223 }
224 else if (xsim[d] < eps) {
225 xsim[d] = eps;
226 }
227
228 if (xsim[d] > 1.0+eps) {
229 MADNESS_EXCEPTION("eval: coordinate upper-bound error in dimension", d);
230 }
231 else if (xsim[d] > 1.0-eps) {
232 xsim[d] = 1.0-eps;
233 }
234 }
235
236 Future<T> result;
237 impl->eval(xsim, impl->key0(), result.remote_ref(impl->world));
238 return result;
239 }
240
241 /// Evaluate function only if point is local returning (true,value); otherwise return (false,0.0)
242
243 /// maxlevel is the maximum depth to search down to --- the max local depth can be
244 /// computed with max_local_depth();
245 std::pair<bool,T> eval_local_only(const Vector<double,NDIM>& xuser, Level maxlevel) const {
246 const double eps=1e-15;
247 verify();
249 coordT xsim;
251 // If on the boundary, move the point just inside the
252 // volume so that the evaluation logic does not fail
253 for (std::size_t d=0; d<NDIM; ++d) {
254 if (xsim[d] < -eps) {
255 MADNESS_EXCEPTION("eval: coordinate lower-bound error in dimension", d);
256 }
257 else if (xsim[d] < eps) {
258 xsim[d] = eps;
259 }
260
261 if (xsim[d] > 1.0+eps) {
262 MADNESS_EXCEPTION("eval: coordinate upper-bound error in dimension", d);
263 }
264 else if (xsim[d] > 1.0-eps) {
265 xsim[d] = 1.0-eps;
266 }
267 }
268 return impl->eval_local_only(xsim,maxlevel);
269 }
270
271 /// Only the invoking process will receive the result via the future
272 /// though other processes may be involved in the evaluation.
273 ///
274 /// Throws if function is not initialized.
275 ///
276 /// This function is a minimally-modified version of eval()
279 const double eps=1e-15;
280 verify();
282 coordT xsim;
284 // If on the boundary, move the point just inside the
285 // volume so that the evaluation logic does not fail
286 for (std::size_t d=0; d<NDIM; ++d) {
287 if (xsim[d] < -eps) {
288 MADNESS_EXCEPTION("eval: coordinate lower-bound error in dimension", d);
289 }
290 else if (xsim[d] < eps) {
291 xsim[d] = eps;
292 }
293
294 if (xsim[d] > 1.0+eps) {
295 MADNESS_EXCEPTION("eval: coordinate upper-bound error in dimension", d);
296 }
297 else if (xsim[d] > 1.0-eps) {
298 xsim[d] = 1.0-eps;
299 }
300 }
301
302 Future<Level> result;
303 impl->evaldepthpt(xsim, impl->key0(), result.remote_ref(impl->world));
304 return result;
305 }
306
307
308 /// Evaluates the function rank at a point in user coordinates. Possible non-blocking comm.
309
310 /// Only the invoking process will receive the result via the future
311 /// though other processes may be involved in the evaluation.
312 ///
313 /// Throws if function is not initialized.
316 const double eps=1e-15;
317 verify();
319 coordT xsim;
321 // If on the boundary, move the point just inside the
322 // volume so that the evaluation logic does not fail
323 for (std::size_t d=0; d<NDIM; ++d) {
324 if (xsim[d] < -eps) {
325 MADNESS_EXCEPTION("eval: coordinate lower-bound error in dimension", d);
326 }
327 else if (xsim[d] < eps) {
328 xsim[d] = eps;
329 }
330
331 if (xsim[d] > 1.0+eps) {
332 MADNESS_EXCEPTION("eval: coordinate upper-bound error in dimension", d);
333 }
334 else if (xsim[d] > 1.0-eps) {
335 xsim[d] = 1.0-eps;
336 }
337 }
338
339 Future<long> result;
340 impl->evalR(xsim, impl->key0(), result.remote_ref(impl->world));
341 return result;
342 }
343
344 /// Evaluates a cube/slice of points (probably for plotting) ... collective but no fence necessary
345
346 /// All processes receive the entire result (which is a rather severe limit
347 /// on the size of the cube that is possible).
348
349 /// Set eval_refine=true to return the refinment levels of
350 /// the given function.
351
352 /// @param[in] cell A Tensor describe the cube where the function to be evaluated in
353 /// @param[in] npt How many points to evaluate in each dimension
354 /// @param[in] eval_refine Wether to return the refinment levels of the given function
356 const std::vector<long>& npt,
357 bool eval_refine = false) const {
358 MADNESS_ASSERT(static_cast<std::size_t>(cell.dim(0))>=NDIM && cell.dim(1)==2 && npt.size()>=NDIM);
360 const double eps=1e-14;
361 verify();
362 reconstruct();
364 for (std::size_t d=0; d<NDIM; ++d) {
365 simlo[d] = cell(d,0);
366 simhi[d] = cell(d,1);
367 }
370
371 // Move the bounding box infintesimally inside dyadic
372 // points so that the evaluation logic does not fail
373 for (std::size_t d=0; d<NDIM; ++d) {
375 MADNESS_ASSERT(simlo[d] >= 0.0);
376 MADNESS_ASSERT(simhi[d] <= 1.0);
377
378 double delta = eps*(simhi[d]-simlo[d]);
379 simlo[d] += delta;
380 simhi[d] -= 2*delta; // deliberate asymmetry
381 }
382 return impl->eval_plot_cube(simlo, simhi, npt, eval_refine);
383 }
384
385
386 /// Evaluates the function at a point in user coordinates. Collective operation.
387
388 /// Throws if function is not initialized.
389 ///
390 /// This function calls eval, blocks until the result is
391 /// available and then broadcasts the result to everyone.
392 /// Therefore, if you are evaluating many points in parallel
393 /// it is \em vastly less efficient than calling eval
394 /// directly, saving the futures, and then forcing all of the
395 /// results.
396 T operator()(const coordT& xuser) const {
398 verify();
400 T result;
401 if (impl->world.rank() == 0) result = eval(xuser).get();
402 impl->world.gop.broadcast(result);
403 //impl->world.gop.fence();
404 return result;
405 }
406
407 /// Evaluates the function at a point in user coordinates. Collective operation.
408
409 /// See "operator()(const coordT& xuser)" for more info
410 T operator()(double x, double y=0, double z=0, double xx=0, double yy=0, double zz=0) const {
411 coordT r;
412 r[0] = x;
413 if (NDIM>=2) r[1] = y;
414 if (NDIM>=3) r[2] = z;
415 if (NDIM>=4) r[3] = xx;
416 if (NDIM>=5) r[4] = yy;
417 if (NDIM>=6) r[5] = zz;
418 return (*this)(r);
419 }
420
421 /// Throws if function is not initialized.
422 ///
423 /// This function mimics operator() by going through the
424 /// tree looking for the depth of the tree at the point.
425 /// It blocks until the result is
426 /// available and then broadcasts the result to everyone.
427 /// Therefore, if you are evaluating many points in parallel
428 /// it is \em vastly less efficient than calling evaldepthpt
429 /// directly, saving the futures, and then forcing all of the
430 /// results.
431 Level depthpt(const coordT& xuser) const {
433 verify();
435 Level result;
436 if (impl->world.rank() == 0) result = evaldepthpt(xuser).get();
437 impl->world.gop.broadcast(result);
438 //impl->world.gop.fence();
439 return result;
440 }
441
442 /// Returns an estimate of the difference ||this-func||^2 from local data
443
444 /// No communication is performed. If the function is not
445 /// reconstructed, it throws an exception. To get the global
446 /// value either do a global sum of the local values or call
447 /// errsq
448 /// @param[in] func Templated interface to the a user specified function
449 template <typename funcT>
450 double errsq_local(const funcT& func) const {
452 verify();
453 if (!is_reconstructed()) MADNESS_EXCEPTION("Function:errsq_local:not reconstructed",0);
454 return impl->errsq_local(func);
455 }
456
457
458 /// Returns an estimate of the difference ||this-func|| ... global sum performed
459
460 /// If the function is compressed, it is reconstructed first. For efficient use
461 /// especially with many functions, reconstruct them all first, and use errsq_local
462 /// instead so you can perform a global sum on all at the same time.
463 /// @param[in] func Templated interface to the a user specified function
464 template <typename funcT>
465 double err(const funcT& func) const {
467 verify();
471 double local = impl->errsq_local(func);
472 impl->world.gop.sum(local);
473 impl->world.gop.fence();
474 return sqrt(local);
475 }
476
477 /// Verifies the tree data structure ... global sync implied
478 void verify_tree() const {
480 if (impl) impl->verify_tree();
481 }
482
483
484 /// Returns true if compressed, false otherwise. No communication.
485
486 /// If the function is not initialized, returns false.
487 bool is_compressed() const {
489 if (impl)
490 return impl->is_compressed();
491 else
492 return false;
493 }
494
495 /// Returns true if reconstructed, false otherwise. No communication.
496
497 /// If the function is not initialized, returns false.
498 bool is_reconstructed() const {
500 if (impl)
501 return impl->is_reconstructed();
502 else
503 return false;
504 }
505
506 /// Returns true if nonstandard-compressed, false otherwise. No communication.
507
508 /// If the function is not initialized, returns false.
509 bool is_nonstandard() const {
511 return impl ? impl->is_nonstandard() : false;
512 }
513
514 /// Returns true if redundant, false otherwise. No communication.
515
516 /// If the function is not initialized, returns false.
517 bool is_redundant() const {
519 return impl ? impl->is_redundant() : false;
520 }
521
522 /// Returns true if redundant_after_merge, false otherwise. No communication.
523
524 /// If the function is not initialized, returns false.
527 return impl ? impl->is_redundant_after_merge() : false;
528 }
529
530 /// Returns the number of nodes in the function tree ... collective global sum
531 std::size_t tree_size() const {
533 if (!impl) return 0;
534 return impl->tree_size();
535 }
536
537 /// print some info about this
538 void print_size(const std::string name) const {
539 if (!impl) {
540 print("function",name,"not assigned yet");
541 } else {
542 impl->print_size(name);
543 }
544 }
545
546 /// Returns the maximum depth of the function tree ... collective global sum
547 std::size_t max_depth() const {
549 if (!impl) return 0;
550 return impl->max_depth();
551 }
552
553
554 /// Returns the maximum local depth of the function tree ... no communications
555 std::size_t max_local_depth() const {
557 if (!impl) return 0;
558 return impl->max_local_depth();
559 }
560
561
562 /// Returns the max number of nodes on a processor
563 std::size_t max_nodes() const {
565 if (!impl) return 0;
566 return impl->max_nodes();
567 }
568
569 /// Returns the min number of nodes on a processor
570 std::size_t min_nodes() const {
572 if (!impl) return 0;
573 return impl->min_nodes();
574 }
575
576
577 /// Returns the number of coefficients in the function ... collective global sum
578 std::size_t size() const {
580 if (!impl) return 0;
581 return impl->size();
582 }
583
584 /// Return the number of coefficients in the function on this processor
585 std::size_t size_local() const {
587 if (!impl) return 0;
588 return impl->size_local();
589 }
590
591
592 /// Returns value of autorefine flag. No communication.
593 bool autorefine() const {
595 if (!impl) return true;
596 return impl->get_autorefine();
597 }
598
599
600 /// Sets the value of the autorefine flag. Optional global fence.
601
602 /// A fence is required to ensure consistent global state.
603 void set_autorefine(bool value, bool fence = true) {
605 verify();
606 impl->set_autorefine(value);
607 if (fence) impl->world.gop.fence();
608 }
609
610
611 /// Returns value of truncation threshold. No communication.
612 double thresh() const {
614 if (!impl) return 0.0;
615 return impl->get_thresh();
616 }
617
618
619 /// Sets the value of the truncation threshold. Optional global fence.
620
621 /// A fence is required to ensure consistent global state.
622 void set_thresh(double value, bool fence = true) {
624 verify();
625 impl->set_thresh(value);
626 if (fence) impl->world.gop.fence();
627 }
628
629
630 /// Returns the number of multiwavelets (k). No communication.
631 int k() const {
633 verify();
634 return impl->get_k();
635 }
636
637
638 /// Truncate the function with optional fence. Compresses with fence if not compressed.
639
640 /// If the truncation threshold is less than or equal to zero the default value
641 /// specified when the function was created is used.
642 /// If the function is not initialized, it just returns.
643 ///
644 /// Returns this for chaining.
645 /// @param[in] tol Tolerance for truncating the coefficients. Default 0.0 means use the implementation's member value \c thresh instead.
646 /// @param[in] fence Do fence
647 Function<T,NDIM>& truncate(double tol = 0.0, bool fence = true) {
649 if (!impl) return *this;
650 verify();
651// if (!is_compressed()) compress();
652 impl->truncate(tol,fence);
654 return *this;
655 }
656
657
658 /// Returns a shared-pointer to the implementation
659 const std::shared_ptr< FunctionImpl<T,NDIM> >& get_impl() const {
661 verify();
662 return impl;
663 }
664
665 /// Replace current FunctionImpl with provided new one
666 void set_impl(const std::shared_ptr< FunctionImpl<T,NDIM> >& impl) {
668 this->impl = impl;
669 }
670
671
672 /// Replace the current functor with the provided new one
673
674 /// presumably the new functor will be a CompositeFunctor, which will
675 /// change the behavior of the function: multiply the functor with the function
676 void set_functor(const std::shared_ptr<FunctionFunctorInterface<T, NDIM> > functor) {
677 this->impl->set_functor(functor);
678 print("set functor in mra.h");
679 }
680
681 bool is_on_demand() const {return this->impl->is_on_demand();}
682
683 /// Replace current FunctionImpl with a new one using the same parameters & map as f
684
685 /// If zero is true the function is initialized to zero, otherwise it is empty
686 template <typename R>
687 void set_impl(const Function<R,NDIM>& f, bool zero = true) {
688 impl = std::shared_ptr<implT>(new implT(*f.get_impl(), f.get_pmap(), zero));
689 if (zero) world().gop.fence();
690 }
691
692 /// Returns the world
693 World& world() const {
695 verify();
696 return impl->world;
697 }
698
699
700 /// Returns a shared pointer to the process map
701 const std::shared_ptr< WorldDCPmapInterface< Key<NDIM> > >& get_pmap() const {
703 verify();
704 return impl->get_pmap();
705 }
706
707 /// replicate this function according to type
708 ///
709 /// ** note that global operations will return unexpected results **
710 /// Be sure you know what you are doing!
711 void replicate(const DistributionType type, bool fence=true) const {
712 verify();
714 else if (type==DistributionType::NodeReplicated) impl->replicate_on_hosts(fence);
715 else MADNESS_EXCEPTION("Function::replicate: unknown DistributionType",type);
716 }
717
718 /// replicate this function, generating a unique pmap
719
720 /// ** note that global operations will return unexpected results **
721 /// Be sure you know what you are doing!
722 void replicate(bool fence=true) const {
723 verify();
724 impl->replicate(fence);
725 }
726
727 /// replicate this function, one copy per host
728
729 /// map will refer the to first rank on each host to avoid inter-node communication
730 /// ** note that global operations will return unexpected results **
731 /// Be sure you know what you are doing!
732 void replicate_on_hosts(bool fence=true) const {
733 verify();
734 impl->replicate_on_hosts(fence);
735 }
736
737
738 /// distribute this function according to newmap
739 void distribute(std::shared_ptr< WorldDCPmapInterface< Key<NDIM> > > newmap) const {
740 verify();
741 impl->distribute(newmap);
742 }
743
744
745 /// Returns the square of the norm of the local function ... no communication
746
747 /// Works in either basis
748 double norm2sq_local() const {
750 verify();
752 "function must be compressed or reconstructed for norm2sq_local");
753 return impl->norm2sq_local();
754 }
755
756
757 /// Returns the 2-norm of the function ... global sum ... works in either basis
758
759 /// See comments for err() w.r.t. applying to many functions.
760 double norm2() const {
762 verify();
764 double local = impl->norm2sq_local();
765
766 impl->world.gop.sum(local);
767 impl->world.gop.fence();
768 return sqrt(local);
769 }
770
771
772 /// Initializes information about the function norm at all length scales
773 void norm_tree(bool fence = true) const {
775 verify();
778 const_cast<Function<T,NDIM>*>(this)->impl->norm_tree(fence);
779 }
780
781
782 /// Compresses the function, transforming into wavelet basis. Possible non-blocking comm.
783
784 /// By default fence=true meaning that this operation completes before returning,
785 /// otherwise if fence=false it returns without fencing and the user must invoke
786 /// world.gop.fence() to assure global completion before using the function
787 /// for other purposes.
788 ///
789 /// Noop if already compressed or if not initialized.
790 ///
791 /// Since reconstruction/compression do not discard information we define them
792 /// as const ... "logical constness" not "bitwise constness".
793 const Function<T,NDIM>& compress(bool fence = true) const {
795 }
796
797
798 /// Compresses the function retaining scaling function coeffs. Possible non-blocking comm.
799
800 /// By default fence=true meaning that this operation completes before returning,
801 /// otherwise if fence=false it returns without fencing and the user must invoke
802 /// world.gop.fence() to assure global completion before using the function
803 /// for other purposes.
804 ///
805 /// Noop if already compressed or if not initialized.
811
812 /// Converts the function standard compressed form. Possible non-blocking comm.
813
814 /// By default fence=true meaning that this operation completes before returning,
815 /// otherwise if fence=false it returns without fencing and the user must invoke
816 /// world.gop.fence() to assure global completion before using the function
817 /// for other purposes.
818 ///
819 /// Must be already compressed.
820 void standard(bool fence = true) {
822 }
823
824 /// Converts the function to redundant form, i.e. sum coefficients on all levels
825
826 /// By default fence=true meaning that this operation completes before returning,
827 /// otherwise if fence=false it returns without fencing and the user must invoke
828 /// world.gop.fence() to assure global completion before using the function
829 /// for other purposes.
830 ///
831 /// Must be already compressed.
832 void make_redundant(bool fence = true) {
834 }
835
836 /// Reconstructs the function, transforming into scaling function basis. Possible non-blocking comm.
837
838 /// By default fence=true meaning that this operation completes before returning,
839 /// otherwise if fence=false it returns without fencing and the user must invoke
840 /// world.gop.fence() to assure global completion before using the function
841 /// for other purposes.
842 ///
843 /// Noop if already reconstructed or if not initialized.
844 ///
845 /// Since reconstruction/compression do not discard information we define them
846 /// as const ... "logical constness" not "bitwise constness".
847 const Function<T,NDIM>& reconstruct(bool fence = true) const {
849 }
850
851 /// changes tree state to given state
852
853 /// Since reconstruction/compression do not discard information we define them
854 /// as const ... "logical constness" not "bitwise constness".
855 /// @param[in] finalstate The final state of the tree
856 /// @param[in] fence Fence after the operation (might not be respected!!!)
859 if (not impl) return *this;
860 TreeState current_state = impl->get_tree_state();
861 if (finalstate == current_state) return *this;
862 MADNESS_CHECK_THROW(current_state != TreeState::unknown, "unknown tree state");
863
864 impl->change_tree_state(finalstate, fence);
865 if (fence && VERIFY_TREE) verify_tree();
866 return *this;
867 }
868
869 /// Sums scaling coeffs down tree restoring state with coeffs only at leaves. Optional fence. Possible non-blocking comm.
870 void sum_down(bool fence = true) const {
872 verify();
873 MADNESS_CHECK_THROW(impl->get_tree_state()==redundant_after_merge, "sum_down requires a redundant_after_merge state");
874 const_cast<Function<T,NDIM>*>(this)->impl->sum_down(fence);
875 const_cast<Function<T,NDIM>*>(this)->impl->set_tree_state(reconstructed);
876
877 if (fence && VERIFY_TREE) verify_tree(); // Must be after in case nonstandard
878 }
879
880
881 /// Inplace autorefines the function. Optional fence. Possible non-blocking comm.
882 template <typename opT>
883 void refine_general(const opT& op, bool fence = true) const {
885 verify();
887 impl->refine(op, fence);
888 }
889
890
892 bool operator()(implT* impl, const Key<NDIM>& key, const nodeT& t) const {
893 return impl->autorefine_square_test(key, t);
894 }
895
896 template <typename Archive> void serialize (Archive& ar) {}
897 };
898
899 /// Inplace autorefines the function using same test as for squaring.
900
901 /// return this for chaining
902 const Function<T,NDIM>& refine(bool fence = true) const {
904 return *this;
905 }
906
907 /// Inplace broadens support in scaling function basis
909 bool fence = true) const {
910 verify();
911 reconstruct();
912 impl->broaden(bc.is_periodic(), fence);
913 }
914
915
916 /// Clears the function as if constructed uninitialized. Optional fence.
917
918 /// Any underlying data will not be freed until the next global fence.
919 void clear(bool fence = true) {
921 if (impl) {
922 World& world = impl->world;
923 impl.reset();
924 if (fence) world.gop.fence();
925 }
926 }
927
928 /// Process 0 prints a summary of all nodes in the tree (collective)
929 void print_tree(std::ostream& os = std::cout) const {
931 if (impl) impl->print_tree(os);
932 }
933
934 /// same as print_tree() but produces JSON-formatted string
935 /// @warning enclose the result in braces to make it a valid JSON object
936 void print_tree_json(std::ostream& os = std::cout) const {
938 if (impl) impl->print_tree_json(os);
939 }
940
941 /// Process 0 prints a graphviz-formatted output of all nodes in the tree (collective)
942 void print_tree_graphviz(std::ostream& os = std::cout) const {
944 os << "digraph G {" << std::endl;
945 if (impl) impl->print_tree_graphviz(os);
946 os << "}" << std::endl;
947 }
948
949 /// Print a summary of the load balancing info
950
951 /// This is serial and VERY expensive
952 void print_info() const {
954 if (impl) impl->print_info();
955 }
956
958 T (*f)(T);
960 void operator()(const Key<NDIM>& key, Tensor<T>& t) const {
961 UNARY_OPTIMIZED_ITERATOR(T, t, *_p0 = f(*_p0));
962 }
963 template <typename Archive> void serialize(Archive& ar) {}
964 };
965
966 /// Inplace unary operation on function values
967 void unaryop(T (*f)(T)) {
968 // Must fence here due to temporary object on stack
969 // stopping us returning before complete
971 }
972
973
974 /// Inplace unary operation on function values
975 template <typename opT>
976 void unaryop(const opT& op, bool fence=true) {
978 verify();
979 reconstruct();
980 impl->unary_op_value_inplace(op, fence);
981 }
982
983
984 /// Unary operation applied inplace to the coefficients
985 template <typename opT>
986 void unaryop_coeff(const opT& op,
987 bool fence = true) {
989 verify();
990 impl->unary_op_coeff_inplace(op, fence);
991 }
992
993
994 /// Unary operation applied inplace to the nodes
995 template <typename opT>
996 void unaryop_node(const opT& op,
997 bool fence = true) {
999 verify();
1000 impl->unary_op_node_inplace(op, fence);
1001 }
1002
1003
1004
1005
1006 static void doconj(const Key<NDIM>, Tensor<T>& t) {
1008 t.conj();
1009 }
1010
1011 /// Inplace complex conjugate. No communication except for optional fence.
1012
1013 /// Returns this for chaining. Works in either basis.
1017 return *this;
1018 }
1019
1020
1021 /// Inplace, scale the function by a constant. No communication except for optional fence.
1022
1023 /// Works in either basis. Returns reference to this for chaining.
1024 template <typename Q>
1025 Function<T,NDIM>& scale(const Q q, bool fence=true) {
1027 verify();
1028 if (VERIFY_TREE) verify_tree();
1029 impl->scale_inplace(q,fence);
1030 return *this;
1031 }
1032
1033
1034 /// Inplace add scalar. No communication except for optional fence.
1037 verify();
1038 if (VERIFY_TREE) verify_tree();
1039 impl->add_scalar_inplace(t,fence);
1040 return *this;
1041 }
1042
1043
1044 /// Inplace, general bi-linear operation in wavelet basis. No communication except for optional fence.
1045
1046 /// If the functions are not in the wavelet basis an exception is thrown since this routine
1047 /// is intended to be fast and unexpected compression is assumed to be a performance bug.
1048 ///
1049 /// Returns this for chaining, can be in states compressed of redundant_after_merge.
1050 ///
1051 /// this and other may have different distributions and may even live in different worlds
1052 ///
1053 /// this <-- this*alpha + other*beta
1054 template <typename Q, typename R>
1056 const Function<Q,NDIM>& other, const R& beta, bool fence=true) {
1058 verify();
1059 other.verify();
1060
1061 // operation is done either in compressed or reconstructed state
1062 TreeState operating_state=this->get_impl()->get_tensor_type()==TT_FULL ? compressed : reconstructed;
1063
1064 TreeState thisstate=impl->get_tree_state();
1065 TreeState otherstate=other.get_impl()->get_tree_state();
1066
1068 MADNESS_CHECK_THROW(thisstate==compressed, "gaxpy: this must be compressed");
1069 MADNESS_CHECK_THROW(otherstate==compressed, "gaxpy: other must be compressed");
1070 impl->gaxpy_inplace(alpha, *other.get_impl(), beta, fence);
1071
1072 } else if (operating_state==reconstructed) {
1073 // this works both in reconstructed and redundant_after_merge states
1075 "gaxpy: this must be reconstructed or redundant_after_merge");
1077 "gaxpy: other must be reconstructed or redundant_after_merge");
1078
1079 impl->gaxpy_inplace_reconstructed(alpha,*other.get_impl(),beta,fence);
1080 } else {
1081 MADNESS_EXCEPTION("unknown tree state",1);
1082 }
1083 return *this;
1084 }
1085
1086
1087 /// Inplace addition of functions in the wavelet basis
1088
1089 /// Using operator notation forces a global fence after every operation.
1090 /// Functions don't need to be compressed, it's the caller's responsibility
1091 /// to choose an appropriate state with performance, usually compressed for 3d,
1092 /// reconstructed for 6d)
1093 template <typename Q>
1096
1097 // do this in reconstructed or compressed form
1099 this->change_tree_state(operating_state);
1101
1102 MADNESS_ASSERT(impl->get_tree_state() == other.get_impl()->get_tree_state());
1103 if (VERIFY_TREE) verify_tree();
1104 if (VERIFY_TREE) other.verify_tree();
1105 return gaxpy(T(1.0), other, Q(1.0), true);
1106 }
1107
1108
1109 /// Inplace subtraction of functions in the wavelet basis
1110
1111 /// Using operator notation forces a global fence after every operation
1112 template <typename Q>
1115 if (NDIM<=3) {
1116 compress();
1117 other.compress();
1118 } else {
1119 reconstruct();
1120 other.reconstruct();
1121 }
1122 MADNESS_ASSERT(impl->get_tree_state() == other.get_impl()->get_tree_state());
1123 if (VERIFY_TREE) verify_tree();
1124 if (VERIFY_TREE) other.verify_tree();
1125 return gaxpy(T(1.0), other, Q(-1.0), true);
1126 }
1127
1128
1129 /// Inplace scaling by a constant
1130
1131 /// Using operator notation forces a global fence after every operation
1132 template <typename Q>
1134 operator*=(const Q q) {
1136 scale(q,true);
1137 return *this;
1138 }
1139
1140
1141 /// Inplace squaring of function ... global comm only if not reconstructed
1142
1143 /// Returns *this for chaining.
1146 if (!is_reconstructed()) reconstruct();
1147 if (VERIFY_TREE) verify_tree();
1148 impl->square_inplace(fence);
1149 return *this;
1150 }
1151
1152 /// Returns *this for chaining.
1155 if (!is_reconstructed()) reconstruct();
1156 if (VERIFY_TREE) verify_tree();
1157 impl->abs_inplace(fence);
1158 return *this;
1159 }
1160
1161 /// Returns *this for chaining.
1164 if (!is_reconstructed()) reconstruct();
1165 if (VERIFY_TREE) verify_tree();
1166 impl->abs_square_inplace(fence);
1167 return *this;
1168 }
1169
1170 /// Returns local contribution to \c int(f(x),x) ... no communication
1171
1172 /// In the wavelet basis this is just the coefficient of the first scaling
1173 /// function which is a constant. In the scaling function basis we
1174 /// must add up contributions from each box.
1175 T trace_local() const {
1177 if (!impl) return 0.0;
1178 if (VERIFY_TREE) verify_tree();
1179 return impl->trace_local();
1180 }
1181
1182
1183 /// Returns global value of \c int(f(x),x) ... global comm required
1184 T trace() const {
1186 if (!impl) return 0.0;
1187 T sum = impl->trace_local();
1188 impl->world.gop.sum(sum);
1189 impl->world.gop.fence();
1190 return sum;
1191 }
1192
1193
1194 /// Returns local part of inner product ... throws if both not compressed
1195 template <typename R>
1196 TENSOR_RESULT_TYPE(T,R) inner_local(const Function<R,NDIM>& g) const {
1203 return impl->inner_local(*(g.get_impl()));
1204 }
1205
1206 /// Returns local part of dot product ... throws if both not compressed
1208 TENSOR_RESULT_TYPE(T,R) dot_local(const Function<R,NDIM>& g) const {
1211 MADNESS_ASSERT(g.is_compressed());
1213 if (VERIFY_TREE) g.verify_tree();
1214 return impl->dot_local(*(g.get_impl()));
1215 }
1216
1217
1218 /// With this being an on-demand function, fill the MRA tree according to different criteria
1219
1220 /// @param[in] g the function after which the MRA structure is modeled (any basis works)
1221 template<typename R>
1223 MADNESS_ASSERT(g.is_initialized());
1225
1226 // clear what we have
1227 impl->get_coeffs().clear();
1228
1229 //leaf_op<T,NDIM> gnode_is_leaf(g.get_impl().get());
1230 Leaf_op_other<T,NDIM> gnode_is_leaf(g.get_impl().get());
1231 impl->make_Vphi(gnode_is_leaf,fence);
1232 return *this;
1233
1234 }
1235
1236 /// With this being an on-demand function, fill the MRA tree according to different criteria
1237
1238 /// @param[in] op the convolution operator for screening
1239 template<typename opT>
1240 Function<T,NDIM>& fill_tree(const opT& op, bool fence=true) {
1242 // clear what we have
1243 impl->get_coeffs().clear();
1246 impl ->make_Vphi(leaf_op,fence);
1247 return *this;
1248 }
1249
1250 /// With this being an on-demand function, fill the MRA tree according to different criteria
1253 // clear what we have
1254 impl->get_coeffs().clear();
1256 impl->make_Vphi(leaf_op,fence);
1257 return *this;
1258 }
1259
1260 /// Special refinement on 6D boxes where the electrons come close (meet)
1261 /// @param[in] op the convolution operator for screening
1262 template<typename opT>
1263 Function<T,NDIM>& fill_cuspy_tree(const opT& op,const bool fence=true){
1265 // clear what we have
1266 impl->get_coeffs().clear();
1268
1270 impl ->make_Vphi(leaf_op,fence);
1271
1272 return *this;
1273 }
1274
1275 /// Special refinement on 6D boxes where the electrons come close (meet)
1278 // clear what we have
1279 impl->get_coeffs().clear();
1281
1283 impl ->make_Vphi(leaf_op,fence);
1284
1285 return *this;
1286 }
1287
1288 /// Special refinement on 6D boxes for the nuclear potentials (regularized with cusp, non-regularized with singularity)
1289 /// @param[in] op the convolution operator for screening
1290 template<typename opT>
1291 Function<T,NDIM>& fill_nuclear_cuspy_tree(const opT& op,const size_t particle,const bool fence=true){
1293 // clear what we have
1294 impl->get_coeffs().clear();
1296
1298 impl ->make_Vphi(leaf_op,fence);
1299
1300 return *this;
1301 }
1302
1303 /// Special refinement on 6D boxes for the nuclear potentials (regularized with cusp, non-regularized with singularity)
1306 // clear what we have
1307 impl->get_coeffs().clear();
1309
1311 impl ->make_Vphi(leaf_op,fence);
1312
1313 return *this;
1314 }
1315
1316 /// perform the hartree product of f*g, invoked by result
1317 template<size_t LDIM, size_t KDIM, typename opT>
1318 void do_hartree_product(const std::vector<std::shared_ptr<FunctionImpl<T,LDIM>>> left,
1319 const std::vector<std::shared_ptr<FunctionImpl<T,KDIM>>> right,
1320 const opT* op) {
1321
1322 // get the right leaf operator
1324 impl->hartree_product(left,right,leaf_op,true);
1325 impl->finalize_sum();
1326// this->truncate();
1327
1328 }
1329
1330 /// perform the hartree product of f*g, invoked by result
1331 template<size_t LDIM, size_t KDIM>
1332 void do_hartree_product(const std::vector<std::shared_ptr<FunctionImpl<T,LDIM>>> left,
1333 const std::vector<std::shared_ptr<FunctionImpl<T,KDIM>>> right) {
1334
1335// hartree_leaf_op<T,KDIM+LDIM> leaf_op(impl.get(),cdata.s0);
1337 impl->hartree_product(left,right,leaf_op,true);
1338 impl->finalize_sum();
1339// this->truncate();
1340
1341 }
1342
1343 /// Returns the inner product
1344
1345 /// Not efficient for computing multiple inner products
1346 /// @param[in] g Function, optionally on-demand
1347 template <typename R>
1350
1351 // fast return if possible
1352 if (not this->is_initialized()) return 0.0;
1353 if (not g.is_initialized()) return 0.0;
1354
1355 // if this and g are the same, use norm2()
1356 if constexpr (std::is_same_v<T,R>) {
1357 if (this->get_impl() == g.get_impl()) {
1358 TreeState state = this->get_impl()->get_tree_state();
1359 if (not(state == reconstructed or state == compressed))
1361 double norm = this->norm2();
1362 return norm * norm;
1363 }
1364 }
1365
1366 // do it case-by-case
1367 if constexpr (std::is_same_v<R,T>) {
1368 if (this->is_on_demand())
1369 return g.inner_on_demand(*this);
1370 if (g.is_on_demand())
1371 return this->inner_on_demand(g);
1372 }
1373
1375 if (VERIFY_TREE) g.verify_tree();
1376
1377 // compute in compressed form if compression is fast, otherwise in redundant form
1379
1381 g.change_tree_state(operating_state,false);
1382 impl->world.gop.fence();
1383
1384 TENSOR_RESULT_TYPE(T,R) local = impl->inner_local(*g.get_impl());
1385 impl->world.gop.sum(local);
1386 impl->world.gop.fence();
1387
1388 // restore state -- no need for this
1389 // change_tree_state(state,false);
1390 // g.change_tree_state(gstate,false);
1391 // impl->world.gop.fence();
1392
1393 return local;
1394 }
1395
1396 /// Return the local part of inner product with external function ... no communication.
1397 /// If you are going to be doing a bunch of inner_ext calls, set
1398 /// keep_redundant to true and then manually undo_redundant when you
1399 /// are finished.
1400 /// @param[in] f Pointer to function of type T that take coordT arguments. This is the externally provided function
1401 /// @param[in] leaf_refine boolean switch to turn on/off refinement past leaf nodes
1402 /// @param[in] keep_redundant boolean switch to turn on/off undo_redundant
1403 /// @return Returns local part of the inner product, i.e. over the domain of all function nodes on this compute node.
1404 T inner_ext_local(const std::shared_ptr< FunctionFunctorInterface<T,NDIM> > f, const bool leaf_refine=true, const bool keep_redundant=false) const {
1407 T local = impl->inner_ext_local(f, leaf_refine);
1409 return local;
1410 }
1411
1412 /// Return the inner product with external function ... requires communication.
1413 /// If you are going to be doing a bunch of inner_ext calls, set
1414 /// keep_redundant to true and then manually undo_redundant when you
1415 /// are finished.
1416 /// @param[in] f Reference to FunctionFunctorInterface. This is the externally provided function
1417 /// @param[in] leaf_refine boolean switch to turn on/off refinement past leaf nodes
1418 /// @param[in] keep_redundant boolean switch to turn on/off undo_redundant
1419 /// @return Returns the inner product
1420 T inner_ext(const std::shared_ptr< FunctionFunctorInterface<T,NDIM> > f, const bool leaf_refine=true, const bool keep_redundant=false) const {
1423 T local = impl->inner_ext_local(f, leaf_refine);
1424 impl->world.gop.sum(local);
1425 impl->world.gop.fence();
1427 return local;
1428 }
1429
1430 /// Return the inner product with external function ... requires communication.
1431 /// If you are going to be doing a bunch of inner_ext calls, set
1432 /// keep_redundant to true and then manually undo_redundant when you
1433 /// are finished.
1434 /// @param[in] f Reference to FunctionFunctorInterface. This is the externally provided function
1435 /// @param[in] leaf_refine boolean switch to turn on/off refinement past leaf nodes
1436 /// @return Returns the inner product
1438 const bool leaf_refine=true) const {
1440 reconstruct();
1441 T local = impl->inner_adaptive_local(f, leaf_refine);
1442 impl->world.gop.sum(local);
1443 impl->world.gop.fence();
1444 return local;
1445 }
1446
1447 /// Return the local part of gaxpy with external function, this*alpha + f*beta ... no communication.
1448 /// @param[in] alpha prefactor for this Function
1449 /// @param[in] f Pointer to function of type T that take coordT arguments. This is the externally provided function
1450 /// @param[in] beta prefactor for f
1451 template <typename L>
1452 void gaxpy_ext(const Function<L,NDIM>& left, T (*f)(const coordT&), T alpha, T beta, double tol, bool fence=true) const {
1454 if (!left.is_reconstructed()) left.reconstruct();
1455 impl->gaxpy_ext(left.get_impl().get(), f, alpha, beta, tol, fence);
1456 }
1457
1458 /// Returns the inner product for one on-demand function
1459
1460 /// It does work, but it might not give you the precision you expect.
1461 /// The assumption is that the function g returns proper sum
1462 /// coefficients on the MRA tree of this. This might not be the case if
1463 /// g is constructed with an implicit multiplication, e.g.
1464 /// result = <this|g>, with g = 1/r12 | gg>
1465 /// @param[in] g on-demand function
1466 template<typename R>
1468 MADNESS_ASSERT(g.is_on_demand() and (not this->is_on_demand()));
1469
1470 constexpr std::size_t LDIM=std::max(NDIM/2,std::size_t(1));
1471 auto func=dynamic_cast<CompositeFunctorInterface<T,NDIM,LDIM>* >(g.get_impl()->get_functor().get());
1473 func->make_redundant(true);
1474 func->replicate_low_dim_functions(true);
1475 this->reconstruct(); // if this == &g we don't need g to be redundant
1476
1478
1479 TENSOR_RESULT_TYPE(T, R) local = impl->inner_local_on_demand(*g.get_impl());
1480 impl->world.gop.sum(local);
1481 impl->world.gop.fence();
1482
1483 return local;
1484 }
1485
1486 /// project this on the low-dim function g: h(x) = <f(x,y) | g(y)>
1487
1488 /// @param[in] g low-dim function
1489 /// @param[in] dim over which dimensions to be integrated: 0..LDIM-1 or LDIM..NDIM-1
1490 /// @return new function of dimension NDIM-LDIM
1491 template <typename R, size_t LDIM>
1493 if (NDIM<=LDIM) MADNESS_EXCEPTION("confused dimensions in project_out?",1);
1494 MADNESS_CHECK_THROW(dim==0 or dim==1,"dim must be 0 or 1 in project_out");
1495 verify();
1496 typedef TENSOR_RESULT_TYPE(T,R) resultT;
1497 static const size_t KDIM=NDIM-LDIM;
1498
1500 .k(g.k()).thresh(g.thresh());
1501 Function<resultT,KDIM> result=factory; // no empty() here!
1502
1504 g.change_tree_state(redundant,false);
1505 world().gop.fence();
1506 this->get_impl()->project_out(result.get_impl().get(),g.get_impl().get(),dim,true);
1507// result.get_impl()->project_out2(this->get_impl().get(),gimpl,dim);
1508 result.world().gop.fence();
1509 g.change_tree_state(reconstructed,false);
1510 result.get_impl()->trickle_down(false);
1511 result.get_impl()->set_tree_state(reconstructed);
1512 result.world().gop.fence();
1513 return result;
1514 }
1515
1516 Function<T,NDIM/2> dirac_convolution(const bool fence=true) const {
1517 constexpr std::size_t LDIM=NDIM/2;
1518 MADNESS_CHECK_THROW(NDIM==2*LDIM,"NDIM must be even");
1519// // this will be the result function
1522 if(!is_reconstructed()) this->reconstruct();
1523 this->get_impl()->do_dirac_convolution(f.get_impl().get(),fence);
1524 return f;
1525 }
1526
1527 /// Replaces this function with one loaded from an archive using the default processor map
1528
1529 /// Archive can be sequential or parallel.
1530 ///
1531 /// The & operator for serializing will only work with parallel archives.
1532 template <typename Archive>
1533 void load(World& world, Archive& ar) {
1535 // Type checking since we are probably circumventing the archive's own type checking
1536 long magic = 0l, id = 0l, ndim = 0l, k = 0l;
1537 ar & magic & id & ndim & k;
1538 MADNESS_ASSERT(magic == 7776768); // Mellow Mushroom Pizza tel.# in Knoxville
1540 MADNESS_ASSERT(ndim == NDIM);
1541
1542 impl.reset(new implT(FunctionFactory<T,NDIM>(world).k(k).empty()));
1543
1544 impl->load(ar);
1545 }
1546
1547
1548 /// Stores the function to an archive
1549
1550 /// Archive can be sequential or parallel.
1551 ///
1552 /// The & operator for serializing will only work with parallel archives.
1553 template <typename Archive>
1554 void store(Archive& ar) const {
1556 verify();
1557 // For type checking, etc.
1558 ar & long(7776768) & long(TensorTypeData<T>::id) & long(NDIM) & long(k());
1559
1560 impl->store(ar);
1561 }
1562
1563 /// change the tensor type of the coefficients in the FunctionNode
1564
1565 /// @param[in] targs target tensor arguments (threshold and full/low rank)
1566 void change_tensor_type(const TensorArgs& targs, bool fence=true) {
1567 if (not impl) return;
1568 impl->change_tensor_type1(targs,fence);
1569 }
1570
1571
1572 /// This is replaced with left*right ... private
1573 template <typename Q, typename opT>
1575 const opT& op, bool fence) {
1577 func.verify();
1578 MADNESS_ASSERT(func.is_reconstructed());
1579 if (VERIFY_TREE) func.verify_tree();
1580 impl.reset(new implT(*func.get_impl(), func.get_pmap(), false));
1581 impl->unaryXX(func.get_impl().get(), op, fence);
1582 return *this;
1583 }
1584
1585 /// Returns vector of FunctionImpl pointers corresponding to vector of functions
1586 template <typename Q, std::size_t D>
1587 static std::vector< std::shared_ptr< FunctionImpl<Q,D> > > vimpl(const std::vector< Function<Q,D> >& v) {
1589 std::vector< std::shared_ptr< FunctionImpl<Q,D> > > r(v.size());
1590 for (unsigned int i=0; i<v.size(); ++i) r[i] = v[i].get_impl();
1591 return r;
1592 }
1593
1594 /// This is replaced with op(vector of functions) ... private
1595 template <typename opT>
1596 Function<T,NDIM>& multiop_values(const opT& op, const std::vector< Function<T,NDIM> >& vf) {
1597 std::vector<implT*> v(vf.size(),NULL);
1598 for (unsigned int i=0; i<v.size(); ++i) {
1599 if (vf[i].is_initialized()) v[i] = vf[i].get_impl().get();
1600 }
1601 impl->multiop_values(op, v);
1602 world().gop.fence();
1603 if (VERIFY_TREE) verify_tree();
1604
1605 return *this;
1606 }
1607
1608 /// apply op on the input vector yielding an output vector of functions
1609
1610 /// (*this) is just a dummy Function to be able to call internal methods in FuncImpl
1611 /// @param[in] op the operator working on vin
1612 /// @param[in] vin vector of input Functions
1613 /// @param[out] vout vector of output Functions vout = op(vin)
1614 template <typename opT>
1616 const std::vector< Function<T,NDIM> >& vin,
1617 std::vector< Function<T,NDIM> >& vout,
1618 const bool fence=true) {
1619 std::vector<implT*> vimplin(vin.size(),NULL);
1620 for (unsigned int i=0; i<vin.size(); ++i) {
1621 if (vin[i].is_initialized()) vimplin[i] = vin[i].get_impl().get();
1622 }
1623 std::vector<implT*> vimplout(vout.size(),NULL);
1624 for (unsigned int i=0; i<vout.size(); ++i) {
1625 if (vout[i].is_initialized()) vimplout[i] = vout[i].get_impl().get();
1626 }
1627
1628 impl->multi_to_multi_op_values(op, vimplin, vimplout, fence);
1629 if (VERIFY_TREE) verify_tree();
1630
1631 }
1632
1633
1634 /// Multiplication of function * vector of functions using recursive algorithm of mulxx
1635 template <typename L, typename R>
1636 void vmulXX(const Function<L,NDIM>& left,
1637 const std::vector< Function<R,NDIM> >& right,
1638 std::vector< Function<T,NDIM> >& result,
1639 double tol,
1640 bool fence) {
1642
1643 std::vector<FunctionImpl<T,NDIM>*> vresult(right.size());
1644 std::vector<const FunctionImpl<R,NDIM>*> vright(right.size());
1645 for (unsigned int i=0; i<right.size(); ++i) {
1646 result[i].set_impl(left,false);
1647 vresult[i] = result[i].impl.get();
1648 vright[i] = right[i].get_impl().get();
1649 }
1650
1651 left.world().gop.fence(); // Is this still essential? Yes.
1652 vresult[0]->mulXXvec(left.get_impl().get(), vright, vresult, tol, fence);
1653 }
1654
1655 /// Same as \c operator* but with optional fence and no automatic reconstruction
1656
1657 /// f or g are on-demand functions
1658 template<typename L, typename R>
1659 void mul_on_demand(const Function<L,NDIM>& f, const Function<R,NDIM>& g, bool fence=true) {
1660 const FunctionImpl<L,NDIM>* fimpl=f.get_impl().get();
1661 const FunctionImpl<R,NDIM>* gimpl=g.get_impl().get();
1662 if (fimpl->is_on_demand() and gimpl->is_on_demand()) {
1663 MADNESS_EXCEPTION("can't multiply two on-demand functions",1);
1664 }
1665
1666 if (fimpl->is_on_demand()) {
1668 impl->multiply(leaf_op1,gimpl,fimpl,fence);
1669 } else {
1671 impl->multiply(leaf_op1,fimpl,gimpl,fence);
1672 }
1673 }
1674
1675 /// sparse transformation of a vector of functions ... private
1676 template <typename R, typename Q>
1677 void vtransform(const std::vector< Function<R,NDIM> >& v,
1678 const Tensor<Q>& c,
1679 std::vector< Function<T,NDIM> >& vresult,
1680 double tol,
1681 bool fence=true) {
1683 vresult[0].impl->vtransform(vimpl(v), c, vimpl(vresult), tol, fence);
1684 }
1685
1686 /// This is replaced with alpha*left + beta*right ... private
1687 template <typename L, typename R>
1689 T beta, const Function<R,NDIM>& right, bool fence) {
1691 left.verify();
1692 right.verify();
1693 MADNESS_ASSERT(left.is_compressed() && right.is_compressed());
1694 if (VERIFY_TREE) left.verify_tree();
1695 if (VERIFY_TREE) right.verify_tree();
1696 impl.reset(new implT(*left.get_impl(), left.get_pmap(), false));
1697 impl->gaxpy(alpha,*left.get_impl(),beta,*right.get_impl(),fence);
1698 return *this;
1699 }
1700
1701 /// This is replaced with mapdim(f) ... private
1702 Function<T,NDIM>& mapdim(const Function<T,NDIM>& f, const std::vector<long>& map, bool fence) {
1704 f.verify();
1705 if (VERIFY_TREE) f.verify_tree();
1706 for (std::size_t i=0; i<NDIM; ++i) MADNESS_ASSERT(map[i]>=0 && static_cast<std::size_t>(map[i])<NDIM);
1707 impl.reset(new implT(*f.impl, f.get_pmap(), false));
1708 impl->mapdim(*f.impl,map,fence);
1709 return *this;
1710 }
1711
1712 /// This is replaced with mirror(f) ... private
1713
1714 /// similar to mapdim, but maps from x to -x, y to -y, and so on
1715 /// Example: mirror a 3d function on the xy plane: mirror={1,1,-1}
1716 /// @param[in] mirror array of -1 and 1, corresponding to mirror or not
1717 Function<T,NDIM>& mirror(const Function<T,NDIM>& f, const std::vector<long>& mirrormap, bool fence) {
1719 f.verify();
1720 if (VERIFY_TREE) f.verify_tree();
1721 for (std::size_t i=0; i<NDIM; ++i) MADNESS_ASSERT((mirrormap[i]==1) or (mirrormap[i]==-1));
1722 impl.reset(new implT(*f.impl, f.get_pmap(), false));
1723 impl->mirror(*f.impl,mirrormap,fence);
1724 return *this;
1725 }
1726
1727 /// This is replaced with mirror(map(f)) ... private
1728
1729 /// first map then mirror!
1730 /// mirror is similar to mapdim, but maps from x to -x, y to -y, and so on
1731 /// Example: mirror a 3d function on the xy plane: mirror={1,1,-1}
1732 /// Example: c4 rotation of a 3d function around the z axis:
1733 /// x->y, y->-x, z->z: map(1,0,2); mirror(-1,1,1)
1734 /// @param[in] map array holding dimensions
1735 /// @param[in] mirror array of -1 and 1, corresponding to mirror or not
1737 const std::vector<long>& map, const std::vector<long>& mirror,
1738 bool fence) {
1740 f.verify();
1741 if (VERIFY_TREE) f.verify_tree();
1742 for (std::size_t i=0; i<mirror.size(); ++i) MADNESS_ASSERT((mirror[i]==1) or (mirror[i]==-1));
1743 for (std::size_t i=0; i<map.size(); ++i) MADNESS_ASSERT(map[i]>=0 && static_cast<std::size_t>(map[i])<NDIM);
1744
1745 impl.reset(new implT(*f.impl, f.get_pmap(), false));
1746 impl->map_and_mirror(*f.impl,map,mirror,fence);
1747 return *this;
1748 }
1749
1750
1751 /// check symmetry of a function by computing the 2nd derivative
1752 double check_symmetry() const {
1753
1755 if (VERIFY_TREE) verify_tree();
1756 double local = impl->check_symmetry_local();
1757 impl->world.gop.sum(local);
1758 impl->world.gop.fence();
1759 double asy=sqrt(local);
1760 if (this->world().rank()==0) print("asymmetry wrt particle",asy);
1762 return asy;
1763 }
1764
1765 /// reduce the rank of the coefficient tensors
1766 Function<T,NDIM>& reduce_rank(const double thresh=0.0, const bool fence=true) {
1767 verify();
1768 double thresh1= (thresh==0.0) ? impl->get_tensor_args().thresh : thresh;
1769 impl->reduce_rank(thresh1,fence);
1770 return *this;
1771 }
1772
1773 /// remove all nodes with level higher than n
1774 Function<T,NDIM>& chop_at_level(const int n, const bool fence=true) {
1775 verify();
1777 impl->chop_at_level(n,true);
1779 return *this;
1780 }
1781 };
1782
1783// template <typename T, typename opT, std::size_t NDIM>
1784 template <typename T, typename opT, std::size_t NDIM>
1785 Function<T,NDIM> multiop_values(const opT& op, const std::vector< Function<T,NDIM> >& vf) {
1787 r.set_impl(vf[0], false);
1788 r.multiop_values(op, vf);
1789 return r;
1790 }
1791
1792 /// Returns new function equal to alpha*f(x) with optional fence
1793 template <typename Q, typename T, std::size_t NDIM>
1795 mul(const Q alpha, const Function<T,NDIM>& f, bool fence=true) {
1797 f.verify();
1798 if (VERIFY_TREE) f.verify_tree();
1800 result.set_impl(f, false);
1801 result.get_impl()->scale_oop(alpha,*f.get_impl(),fence);
1802 return result;
1803 }
1804
1805
1806 /// Returns new function equal to f(x)*alpha with optional fence
1807 template <typename Q, typename T, std::size_t NDIM>
1809 mul(const Function<T,NDIM>& f, const Q alpha, bool fence=true) {
1811 return mul(alpha,f,fence);
1812 }
1813
1814
1815 /// Returns new function equal to f(x)*alpha
1816
1817 /// Using operator notation forces a global fence after each operation
1818 template <typename Q, typename T, std::size_t NDIM>
1821 return mul(alpha, f, true);
1822 }
1823
1824 /// Returns new function equal to alpha*f(x)
1825
1826 /// Using operator notation forces a global fence after each operation
1827 template <typename Q, typename T, std::size_t NDIM>
1830 return mul(alpha, f, true);
1831 }
1832
1833 /// Sparse multiplication --- left and right must be reconstructed and if tol!=0 have tree of norms already created
1834 template <typename L, typename R,std::size_t NDIM>
1836 mul_sparse(const Function<L,NDIM>& left, const Function<R,NDIM>& right, double tol, bool fence=true) {
1838 left.verify();
1839 right.verify();
1841 if (VERIFY_TREE) left.verify_tree();
1842 if (VERIFY_TREE) right.verify_tree();
1843
1845 result.set_impl(left, false);
1846 result.get_impl()->mulXX(left.get_impl().get(), right.get_impl().get(), tol, fence);
1847 return result;
1848 }
1849
1850 /// Same as \c operator* but with optional fence and no automatic reconstruction
1851 template <typename L, typename R,std::size_t NDIM>
1853 mul(const Function<L,NDIM>& left, const Function<R,NDIM>& right, bool fence=true) {
1854 return mul_sparse(left,right,0.0,fence);
1855 }
1856
1857 /// Generate new function = op(left,right) where op acts on the function values
1858 template <typename L, typename R, typename opT, std::size_t NDIM>
1860 binary_op(const Function<L,NDIM>& left, const Function<R,NDIM>& right, const opT& op, bool fence=true) {
1862 if (!left.is_reconstructed()) left.reconstruct();
1863 if (!right.is_reconstructed()) right.reconstruct();
1864
1866 result.set_impl(left, false);
1867 result.get_impl()->binaryXX(left.get_impl().get(), right.get_impl().get(), op, fence);
1868 return result;
1869 }
1870
1871 /// Out of place application of unary operation to function values with optional fence
1872 template <typename Q, typename opT, std::size_t NDIM>
1874 unary_op(const Function<Q,NDIM>& func, const opT& op, bool fence=true) {
1875 if (!func.is_reconstructed()) func.reconstruct();
1878 result.set_impl(func, false);
1879 result.get_impl()->unaryXXvalues(func.get_impl().get(), op, fence);
1880 return result;
1881 }
1882
1883
1884 /// Out of place application of unary operation to scaling function coefficients with optional fence
1885 template <typename Q, typename opT, std::size_t NDIM>
1887 unary_op_coeffs(const Function<Q,NDIM>& func, const opT& op, bool fence=true) {
1888 if (!func.is_reconstructed()) func.reconstruct();
1890 return result.unary_op_coeffs(func,op,fence);
1891 }
1892
1893 /// Use the vmra/mul(...) interface instead
1894
1895 /// This so that we don't have to have friend functions in a different header.
1896 ///
1897 /// If using sparsity (tol != 0) you must have created the tree of norms
1898 /// already for both left and right.
1899 template <typename L, typename R, std::size_t D>
1900 std::vector< Function<TENSOR_RESULT_TYPE(L,R),D> >
1901 vmulXX(const Function<L,D>& left, const std::vector< Function<R,D> >& vright, double tol, bool fence=true) {
1902 if (vright.size() == 0) return std::vector< Function<TENSOR_RESULT_TYPE(L,R),D> >();
1903 std::vector< Function<TENSOR_RESULT_TYPE(L,R),D> > vresult(vright.size());
1904 vresult[0].vmulXX(left, vright, vresult, tol, fence);
1905 return vresult;
1906 }
1907
1908 /// Multiplies two functions with the new result being of type TensorResultType<L,R>
1909
1910 /// Using operator notation forces a global fence after each operation but also
1911 /// enables us to automatically reconstruct the input functions as required.
1912 template <typename L, typename R, std::size_t NDIM>
1914 operator*(const Function<L,NDIM>& left, const Function<R,NDIM>& right) {
1915 if (!left.is_reconstructed()) left.reconstruct();
1916 if (!right.is_reconstructed()) right.reconstruct();
1917 MADNESS_ASSERT(not (left.is_on_demand() or right.is_on_demand()));
1918 return mul(left,right,true);
1919 }
1920
1921 /// Performs a Hartree/outer product on the two given low-dimensional function vectors
1922
1923 /// @return result(x,y) = \sum_i f_i(x) g_i(y)
1924 template<typename T, std::size_t KDIM, std::size_t LDIM>
1926 hartree_product(const std::vector<Function<T,KDIM>>& left, const std::vector<Function<T,LDIM>>& right) {
1927
1928 MADNESS_CHECK_THROW(left.size()==right.size(), "hartree_product: left and right must have same size");
1929 if (left.size()==0) return Function<T,KDIM+LDIM>();
1930
1932
1934 .k(left.front().k()).thresh(thresh);
1935 Function<T,KDIM+LDIM> result=factory.empty();
1936
1937 // some prep work
1940 std::vector<std::shared_ptr<FunctionImpl<T,KDIM>>> vleft=get_impl(left);
1941 std::vector<std::shared_ptr<FunctionImpl<T,LDIM>>> vright=get_impl(right);
1942
1944
1945 return result;
1946
1947 }
1948
1949 /// Performs a Hartree product on the two given low-dimensional functions
1950 template<typename T, std::size_t KDIM, std::size_t LDIM>
1953 typedef std::vector<Function<T,KDIM>> vector;
1954 return hartree_product(vector({left2}),vector({right2}));
1955 }
1956
1957 /// Performs a Hartree product on the two given low-dimensional functions
1958 template<typename T, std::size_t KDIM, std::size_t LDIM, typename opT>
1961 const opT& op) {
1962
1963 // we need both sum and difference coeffs for error estimation
1964 Function<T,KDIM>& left = const_cast< Function<T,KDIM>& >(left2);
1965 Function<T,LDIM>& right = const_cast< Function<T,LDIM>& >(right2);
1966
1968
1970 .k(left.k()).thresh(thresh);
1971 Function<T,KDIM+LDIM> result=factory.empty();
1972
1973 if (result.world().rank()==0) {
1974 print("incomplete FunctionFactory in Function::hartree_product");
1975 print("thresh: ", thresh);
1976 }
1977 bool same=(left2.get_impl()==right2.get_impl());
1978
1979 // some prep work
1980 left.make_nonstandard(true, true);
1981 right.make_nonstandard(true, true);
1982
1983 std::vector<std::shared_ptr<FunctionImpl<T,KDIM>>> vleft;
1984 std::vector<std::shared_ptr<FunctionImpl<T,LDIM>>> vright;
1985 vleft.push_back(left.get_impl());
1986 vright.push_back(right.get_impl());
1987 result.do_hartree_product(vleft,right,&op);
1988
1989 left.standard(false);
1990 if (not same) right.standard(false);
1991 left2.world().gop.fence();
1992
1993 return result;
1994 }
1995
1996 /// adds beta*right only left: alpha*left + beta*right optional fence and no automatic compression
1997
1998 /// left and right might live in different worlds, the accumulation is non-blocking
1999 template <typename L, typename R,std::size_t NDIM>
2000 void
2002 TENSOR_RESULT_TYPE(L,R) beta, const Function<R,NDIM>& right, bool fence=true) {
2005 left.gaxpy(alpha, right, beta, fence);
2006 }
2007
2008 /// Returns new function alpha*left + beta*right optional fence and no automatic compression
2009 template <typename L, typename R,std::size_t NDIM>
2012 TENSOR_RESULT_TYPE(L,R) beta, const Function<R,NDIM>& right, bool fence=true) {
2015 return result.gaxpy_oop(alpha, left, beta, right, fence);
2016 }
2017
2018 /// Same as \c operator+ but with optional fence and no automatic compression
2019 template <typename L, typename R,std::size_t NDIM>
2021 add(const Function<L,NDIM>& left, const Function<R,NDIM>& right, bool fence=true) {
2022 return gaxpy_oop(TENSOR_RESULT_TYPE(L,R)(1.0), left,
2023 TENSOR_RESULT_TYPE(L,R)(1.0), right, fence);
2024 }
2025
2026
2027 /// Returns new function alpha*left + beta*right optional fence, having both addends reconstructed
2028 template<typename T, std::size_t NDIM>
2030 const double beta, const Function<T,NDIM>& right, const bool fence=true) {
2031 Function<T,NDIM> result;
2032 result.set_impl(right,false);
2033
2036 result.get_impl()->gaxpy_oop_reconstructed(alpha,*left.get_impl(),beta,*right.get_impl(),fence);
2037 return result;
2038
2039 }
2040
2041 /// Adds two functions with the new result being of type TensorResultType<L,R>
2042
2043 /// Using operator notation forces a global fence after each operation
2044 template <typename L, typename R, std::size_t NDIM>
2046 operator+(const Function<L,NDIM>& left, const Function<R,NDIM>& right) {
2047 if (VERIFY_TREE) left.verify_tree();
2048 if (VERIFY_TREE) right.verify_tree();
2049
2050 // no compression for high-dimensional functions
2051 if (NDIM==6) {
2052 left.reconstruct();
2053 right.reconstruct();
2054 return gaxpy_oop_reconstructed(1.0,left,1.0,right,true);
2055 } else {
2056 if (!left.is_compressed()) left.compress();
2057 if (!right.is_compressed()) right.compress();
2058 return add(left,right,true);
2059 }
2060 }
2061
2062 /// Same as \c operator- but with optional fence and no automatic compression
2063 template <typename L, typename R,std::size_t NDIM>
2065 sub(const Function<L,NDIM>& left, const Function<R,NDIM>& right, bool fence=true) {
2066 return gaxpy_oop(TENSOR_RESULT_TYPE(L,R)(1.0), left,
2067 TENSOR_RESULT_TYPE(L,R)(-1.0), right, fence);
2068 }
2069
2070
2071 /// Subtracts two functions with the new result being of type TensorResultType<L,R>
2072
2073 /// Using operator notation forces a global fence after each operation
2074 template <typename L, typename R, std::size_t NDIM>
2076 operator-(const Function<L,NDIM>& left, const Function<R,NDIM>& right) {
2078 // no compression for high-dimensional functions
2079 if (NDIM==6) {
2080 left.reconstruct();
2081 right.reconstruct();
2082 return gaxpy_oop_reconstructed(1.0,left,-1.0,right,true);
2083 } else {
2084 if (!left.is_compressed()) left.compress();
2085 if (!right.is_compressed()) right.compress();
2086 return sub(left,right,true);
2087 }
2088 }
2089
2090 /// Create a new copy of the function with different distribution and optional fence
2091
2092 /// Works in either basis. Different distributions imply
2093 /// asynchronous communication and the optional fence is
2094 /// collective.
2095 template <typename T, std::size_t NDIM>
2097 const std::shared_ptr< WorldDCPmapInterface< Key<NDIM> > >& pmap,
2098 bool fence = true) {
2100 f.verify();
2101 Function<T,NDIM> result;
2103 result.set_impl(std::shared_ptr<implT>(new implT(*f.get_impl(), pmap, false)));
2104 result.get_impl()->copy_coeffs(*f.get_impl(), fence);
2105 if (VERIFY_TREE) result.verify_tree();
2106 return result;
2107 }
2108
2109 /// Create a new copy of the function with the same distribution and optional fence
2110 template <typename T, std::size_t NDIM>
2113 return copy(f, f.get_pmap(), fence);
2114 }
2115
2116 /// Create a new copy of function f living in world (might differ from f.world)
2117
2118 /// uses the default processor map of world
2119 template <typename T, std::size_t NDIM>
2124
2125 // create a new function with pmap distribution, same parameters as f, but no coeffs
2126 Function<T,NDIM> result;
2127 result.set_impl(std::make_shared<implT>(world,*f.get_impl(), pmap, false));
2128 // copy f's coefficients to result
2129 result.get_impl()->copy_coeffs(*f.get_impl(), fence);
2130 return result;
2131 }
2132
2133 /// Type conversion implies a deep copy. No communication except for optional fence.
2134
2135 /// Works in either basis but any loss of precision may result in different errors
2136 /// in applied in a different basis.
2137 ///
2138 /// The new function is formed with the options from the default constructor.
2139 ///
2140 /// There is no automatic type conversion since this is generally a rather dangerous
2141 /// thing and because there would be no way to make the fence optional.
2142 template <typename T, typename Q, std::size_t NDIM>
2145 f.verify();
2146 Function<Q,NDIM> result;
2147 result.set_impl(f, false);
2148 result.get_impl()->copy_coeffs(*f.get_impl(), fence);
2149 return result;
2150 }
2151
2152
2153 /// Return the complex conjugate of the input function with the same distribution and optional fence
2154
2155 /// !!! The fence is actually not optional in the current implementation !!!
2156 template <typename T, std::size_t NDIM>
2159 Function<T,NDIM> result = copy(f,true);
2160 return result.conj(fence);
2161 }
2162
2163 /// Apply operator on a hartree product of two low-dimensional functions
2164
2165 /// Supposed to be something like result= G( f(1)*f(2))
2166 /// the hartree product is never constructed explicitly, but its coeffs are
2167 /// constructed on the fly and processed immediately.
2168 /// @param[in] op the operator
2169 /// @param[in] f1 function of particle 1
2170 /// @param[in] f2 function of particle 2
2171 /// @param[in] fence if we shall fence
2172 /// @return a function of dimension NDIM=LDIM+LDIM
2173 template <typename opT, typename T, std::size_t LDIM>
2174 Function<TENSOR_RESULT_TYPE(typename opT::opT,T), LDIM+LDIM>
2175 apply(const opT& op, const std::vector<Function<T,LDIM>>& f1, const std::vector<Function<T,LDIM>>& f2, bool fence=true) {
2176
2177 World& world=f1.front().world();
2178
2179 typedef TENSOR_RESULT_TYPE(T,typename opT::opT) resultT;
2180 typedef std::vector<Function<T,LDIM>> vecfuncL;
2181
2182 vecfuncL& ff1 = const_cast< vecfuncL& >(f1);
2183 vecfuncL& ff2 = const_cast< vecfuncL& >(f2);
2184
2185 bool same=(ff1[0].get_impl()==ff2[0].get_impl());
2186
2187 reconstruct(world,f1,false);
2188 reconstruct(world,f2,false);
2189 world.gop.fence();
2190 // keep the leaves! They are assumed to be there later
2191 // even for modified op we need NS form for the hartree_leaf_op
2192 for (auto& f : f1) f.make_nonstandard(true,false);
2193 for (auto& f : f2) f.make_nonstandard(true,false);
2194 world.gop.fence();
2195
2196
2199 Function<resultT,LDIM+LDIM> result=factory.empty().fence();
2200
2201 result.get_impl()->reset_timer();
2202 op.reset_timer();
2203
2204 // will fence here
2205 for (size_t i=0; i<f1.size(); ++i)
2206 result.get_impl()->recursive_apply(op, f1[i].get_impl().get(),f2[i].get_impl().get(),false);
2207 world.gop.fence();
2208
2209 if (op.print_timings) {
2210 result.get_impl()->print_timer();
2211 op.print_timer();
2212 }
2213
2214 result.get_impl()->finalize_apply(); // need fence before reconstruct
2215
2216 if (op.modified()) {
2217 result.get_impl()->trickle_down(true);
2218 } else {
2219 result.get_impl()->reconstruct(true);
2220 }
2221 standard(world,ff1,false);
2222 if (not same) standard(world,ff2,false);
2223
2224 return result;
2225 }
2226
2227
2228 /// Apply operator ONLY in non-standard form - required other steps missing !!
2229 template <typename opT, typename R, std::size_t NDIM>
2230 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM>
2231 apply_only(const opT& op, const Function<R,NDIM>& f, bool fence=true) {
2232 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> result;
2233
2234 constexpr std::size_t OPDIM=opT::opdim;
2235 constexpr bool low_dim=(OPDIM*2==NDIM); // apply on some dimensions only
2236
2237 // specialized version for 3D
2238 if (NDIM <= 3 and (not low_dim)) {
2239 result.set_impl(f, false);
2240 result.get_impl()->apply(op, *f.get_impl(), fence);
2241
2242 } else { // general version for higher dimension
2243 //bool print_timings=false;
2244 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> r1;
2245
2246 result.set_impl(f, false);
2247 r1.set_impl(f, false);
2248
2249 result.get_impl()->reset_timer();
2250 op.reset_timer();
2251
2252 result.get_impl()->apply_source_driven(op, *f.get_impl(), fence);
2253
2254 // recursive_apply is about 20% faster than apply_source_driven
2255 //result.get_impl()->recursive_apply(op, f.get_impl().get(),
2256 // r1.get_impl().get(),true); // will fence here
2257
2258 }
2259
2260 return result;
2261 }
2262
2263 /// Apply operator in non-standard form
2264
2265 /// Returns a new function with the same distribution
2266 ///
2267 /// !!! For the moment does NOT respect fence option ... always fences
2268 /// if the operator acts on one particle only the result will be sorted as
2269 /// g.particle=1: g(f) = \int g(x,x') f(x',y) dx' = result(x,y)
2270 /// g.particle=2: g(f) = \int g(y,y') f(x,y') dy' = result(x,y)
2271 /// for the second case it will notably *not* be as it is implemented in the partial inner product!
2272 /// g.particle=2 g(f) = result(x,y)
2273 /// inner(g(y,y'),f(x,y'),1,1) = result(y,x)
2274 /// also note the confusion with the counting of the particles/integration variables
2275 template <typename opT, typename R, std::size_t NDIM>
2276 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM>
2277 apply(const opT& op, const Function<R,NDIM>& f, bool fence=true) {
2278
2279 typedef TENSOR_RESULT_TYPE(typename opT::opT,R) resultT;
2280 Function<R,NDIM>& ff = const_cast< Function<R,NDIM>& >(f);
2282
2283 MADNESS_ASSERT(not f.is_on_demand());
2284 bool print_timings=op.print_timings;
2285
2286 if (VERIFY_TREE) ff.verify_tree();
2287 ff.reconstruct();
2288 if (print_timings) ff.print_size("ff in apply after reconstruct");
2289
2290 if (op.modified()) {
2291
2292 ff.change_tree_state(redundant);
2293// ff.get_impl()->make_redundant(true);
2294 result = apply_only(op, ff, fence);
2295 ff.get_impl()->undo_redundant(false);
2296 result.get_impl()->trickle_down(true);
2297
2298 } else {
2299
2300 // saves the standard() step, which is very expensive in 6D
2301// Function<R,NDIM> fff=copy(ff);
2303 fff.make_nonstandard(op.doleaves, true);
2304 if (print_timings) fff.print_size("ff in apply after make_nonstandard");
2305 if ((print_timings) and (f.world().rank()==0)) {
2306 fff.get_impl()->timer_filter.print("filter");
2307 fff.get_impl()->timer_compress_svd.print("compress_svd");
2308 }
2309 result = apply_only(op, fff, fence);
2310 result.get_impl()->set_tree_state(nonstandard_after_apply);
2311 ff.world().gop.fence();
2312 if (print_timings) result.print_size("result after apply_only");
2313
2314 // svd-tensors need some post-processing
2315 if (result.get_impl()->get_tensor_type()==TT_2D) {
2316 double elapsed=result.get_impl()->finalize_apply();
2317 if (print_timings) printf("time in finalize_apply %8.2f\n",elapsed);
2318 }
2319 if (print_timings) {
2320 result.get_impl()->print_timer();
2321 op.print_timer();
2322 }
2323
2324 result.get_impl()->reconstruct(true);
2325
2326// fff.clear();
2327 if (op.destructive()) {
2328 ff.world().gop.fence();
2329 ff.clear();
2330 } else {
2331 // ff.standard();
2332 ff.reconstruct();
2333 }
2334
2335 }
2336 if (print_timings) result.print_size("result after reconstruction");
2337 return result;
2338 }
2339
2340
2341 template <typename opT, typename R, std::size_t NDIM>
2342 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM>
2343 apply_1d_realspace_push(const opT& op, const Function<R,NDIM>& f, int axis, bool fence=true) {
2345 Function<R,NDIM>& ff = const_cast< Function<R,NDIM>& >(f);
2346 if (VERIFY_TREE) ff.verify_tree();
2347 ff.reconstruct();
2348
2349 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> result;
2350
2351 result.set_impl(ff, false);
2352 result.get_impl()->apply_1d_realspace_push(op, ff.get_impl().get(), axis, fence);
2353 result.get_impl()->set_tree_state(redundant_after_merge);
2354 return result;
2355 }
2356
2357
2358 /// Generate a new function by reordering dimensions ... optional fence
2359
2360 /// You provide an array of dimension NDIM that maps old to new dimensions
2361 /// according to
2362 /// \code
2363 /// newdim = mapdim[olddim]
2364 /// \endcode
2365 /// Works in either scaling function or wavelet basis.
2366 ///
2367 /// Would be easy to modify this to also change the procmap here
2368 /// if desired but presently it uses the same procmap as f.
2369 template <typename T, std::size_t NDIM>
2371 mapdim(const Function<T,NDIM>& f, const std::vector<long>& map, bool fence=true) {
2373 Function<T,NDIM> result;
2374 return result.mapdim(f,map,fence);
2375 }
2376
2377 /// Generate a new function by mirroring within the dimensions .. optional fence
2378
2379 /// similar to mapdim
2380 /// @param[in] mirror array with -1 and 1, corresponding to mirror this dimension or not
2381 template <typename T, std::size_t NDIM>
2383 mirror(const Function<T,NDIM>& f, const std::vector<long>& mirrormap, bool fence=true) {
2385 Function<T,NDIM> result;
2386 return result.mirror(f,mirrormap,fence);
2387 }
2388
2389 /// This is replaced with mirror(map(f)), optional fence
2390
2391 /// first map then mirror!
2392 /// mirror is similar to mapdim, but maps from x to -x, y to -y, and so on
2393 /// Example: mirror a 3d function on the xy plane: mirror={1,1,-1}
2394 /// Example: c4 rotation of a 3d function around the z axis:
2395 /// x->y, y->-x, z->z: map(1,0,2); mirror(-1,1,1)
2396 /// @param[in] map array holding dimensions
2397 /// @param[in] mirror array of -1 and 1, corresponding to mirror or not
2398 template <typename T, std::size_t NDIM>
2400 map_and_mirror(const Function<T,NDIM>& f, const std::vector<long>& map,
2401 const std::vector<long>& mirror, bool fence=true) {
2403 Function<T,NDIM> result;
2404 return result.map_and_mirror(f,map,mirror,fence);
2405 }
2406
2407
2408 /// swap particles 1 and 2
2409
2410 /// param[in] f a function of 2 particles f(1,2)
2411 /// return the input function with particles swapped g(1,2) = f(2,1)
2412 template <typename T, std::size_t NDIM>
2413 typename std::enable_if_t<NDIM%2==0, Function<T,NDIM>>
2415 // this could be done more efficiently for SVD, but it works decently
2416 std::vector<long> map(NDIM);
2417 constexpr std::size_t LDIM=NDIM/2;
2418 static_assert(LDIM*2==NDIM);
2419 for (auto d=0; d<LDIM; ++d) {
2420 map[d]=d+LDIM;
2421 map[d+LDIM]=d;
2422 }
2423// map[0]=3;
2424// map[1]=4;
2425// map[2]=5; // 2 -> 1
2426// map[3]=0;
2427// map[4]=1;
2428// map[5]=2; // 1 -> 2
2429 return mapdim(f,map);
2430 }
2431
2432 /// symmetrize a function
2433
2434 /// @param[in] symmetry possibilities are:
2435 /// (anti-) symmetric particle permutation ("sy_particle", "antisy_particle")
2436 /// symmetric mirror plane ("xy", "xz", "yz")
2437 /// @return a new function symmetrized according to the input parameter
2438 template <typename T, std::size_t NDIM>
2440 symmetrize(const Function<T,NDIM>& f, const std::string symmetry, bool fence=true) {
2441 Function<T,NDIM> result;
2442
2443 MADNESS_ASSERT(NDIM==6); // works only for pair functions
2444 std::vector<long> map(NDIM);
2445
2446 // symmetric particle permutation
2447 if (symmetry=="sy_particle") {
2448 map[0]=3; map[1]=4; map[2]=5;
2449 map[3]=0; map[4]=1; map[5]=2;
2450 } else if (symmetry=="cx") {
2451 map[0]=0; map[1]=2; map[2]=1;
2452 map[3]=3; map[4]=5; map[5]=4;
2453
2454 } else if (symmetry=="cy") {
2455 map[0]=2; map[1]=1; map[2]=0;
2456 map[3]=5; map[4]=4; map[5]=3;
2457
2458 } else if (symmetry=="cz") {
2459 map[0]=1; map[1]=0; map[2]=2;
2460 map[3]=4; map[4]=3; map[5]=5;
2461
2462 } else {
2463 if (f.world().rank()==0) {
2464 print("unknown parameter in symmetrize:",symmetry);
2465 }
2466 MADNESS_EXCEPTION("unknown parameter in symmetrize",1);
2467 }
2468
2469 result.mapdim(f,map,true); // need to fence here
2470 result.get_impl()->average(*f.get_impl());
2471
2472 return result;
2473 }
2474
2475
2476
2477 /// multiply a high-dimensional function with a low-dimensional function
2478
2479 /// @param[in] f NDIM function of 2 particles: f=f(1,2)
2480 /// @param[in] g LDIM function of 1 particle: g=g(1) or g=g(2)
2481 /// @param[in] particle if g=g(1) or g=g(2)
2482 /// @return h(1,2) = f(1,2) * g(p)
2483 template<typename T, std::size_t NDIM, std::size_t LDIM>
2484 Function<T,NDIM> multiply(const Function<T,NDIM> f, const Function<T,LDIM> g, const int particle, const bool fence=true) {
2485
2486 static_assert(LDIM+LDIM==NDIM);
2488
2489 Function<T,NDIM> result;
2490 result.set_impl(f, false);
2491
2492// Function<T,NDIM>& ff = const_cast< Function<T,NDIM>& >(f);
2493// Function<T,LDIM>& gg = const_cast< Function<T,LDIM>& >(g);
2494
2495 f.change_tree_state(redundant,false);
2496 g.change_tree_state(redundant);
2497 FunctionImpl<T,NDIM>* fimpl=f.get_impl().get();
2498 FunctionImpl<T,LDIM>* gimpl=g.get_impl().get();
2499
2500 result.get_impl()->multiply(fimpl,gimpl,particle);
2501 result.world().gop.fence();
2502
2503 f.change_tree_state(reconstructed,false);
2504 g.change_tree_state(reconstructed);
2505 return result;
2506 }
2507
2508
2509 template <typename T, std::size_t NDIM>
2514 bool fence=true)
2515 {
2517 Function<T,NDIM> result = FunctionFactory<T,NDIM>(other.world()).k(k).thresh(thresh).empty();
2518 other.reconstruct();
2519 result.get_impl()->project(*other.get_impl(),fence);
2520 return result;
2521 }
2522
2523
2524 /// Computes the scalar/inner product between two functions
2525
2526 /// In Maple this would be \c int(conjugate(f(x))*g(x),x=-infinity..infinity)
2527 template <typename T, typename R, std::size_t NDIM>
2530 return f.inner(g);
2531 }
2532
2533
2534 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2535
2536 /// syntax similar to the inner product in tensor.h
2537 /// e.g result=inner<3>(f,g),{0},{1}) : r(x,y) = int f(x1,x) g(y,x1) dx1
2538 /// @param[in] task 0: everything, 1; prepare only (fence), 2: work only (no fence), 3: finalize only (fence)
2539 template<std::size_t NDIM, typename T, std::size_t LDIM, typename R, std::size_t KDIM,
2540 std::size_t CDIM = (KDIM + LDIM - NDIM) / 2>
2541 std::vector<Function<TENSOR_RESULT_TYPE(T, R), NDIM>>
2542 innerXX(const Function<T, LDIM>& f, const std::vector<Function<R, KDIM>>& vg, const std::array<int, CDIM> v1,
2543 const std::array<int, CDIM> v2, int task=0) {
2544 bool prepare = ((task==0) or (task==1));
2545 bool work = ((task==0) or (task==2));
2546 bool finish = ((task==0) or (task==3));
2547
2548 static_assert((KDIM + LDIM - NDIM) % 2 == 0, "faulty dimensions in inner (partial version)");
2549 static_assert(KDIM + LDIM - 2 * CDIM == NDIM, "faulty dimensions in inner (partial version)");
2550
2551 // contraction indices must be contiguous and either in the beginning or at the end
2552 for (int i=0; i<CDIM-1; ++i) MADNESS_CHECK((v1[i]+1)==v1[i+1]);
2553 MADNESS_CHECK((v1[0]==0) or (v1[CDIM-1]==LDIM-1));
2554
2555 for (int i=0; i<CDIM-1; ++i) MADNESS_CHECK((v2[i]+1)==v2[i+1]);
2556 MADNESS_CHECK((v2[0]==0) or (v2[CDIM-1]==KDIM-1));
2557
2558 MADNESS_CHECK(f.is_initialized());
2560 MADNESS_CHECK(f.world().id() == vg[0].world().id());
2561 // this needs to be run in a single world, so that all coefficients are local.
2562 // Use macrotasks if run on multiple processes.
2563 World& world=f.world();
2564 MADNESS_CHECK(world.size() == 1);
2565
2566 if (prepare) {
2567 f.change_tree_state(nonstandard);
2569 world.gop.fence();
2570 f.get_impl()->compute_snorm_and_dnorm(false);
2571 for (auto& g : vg) g.get_impl()->compute_snorm_and_dnorm(false);
2572 world.gop.fence();
2573 }
2574
2575 typedef TENSOR_RESULT_TYPE(T, R) resultT;
2576 std::vector<Function<resultT,NDIM>> result(vg.size());
2577 if (work) {
2579 for (int i=0; i<vg.size(); ++i) {
2581 .k(f.k()).thresh(f.thresh()).empty().nofence();
2582 result[i].get_impl()->partial_inner(*f.get_impl(),*(vg[i]).get_impl(),v1,v2);
2583 result[i].get_impl()->set_tree_state(nonstandard_after_apply);
2584 }
2585 world.gop.set_forbid_fence(false);
2586 }
2587
2588 if (finish) {
2589
2590 world.gop.fence();
2591// result.get_impl()->reconstruct(true);
2592
2594// result.reconstruct();
2595 // restore initial state of g and h
2596 auto erase_list = [] (const auto& funcimpl) {
2597 typedef typename std::decay_t<decltype(funcimpl)>::keyT keyTT;
2598 std::list<keyTT> to_be_erased;
2599 for (auto it=funcimpl.get_coeffs().begin(); it!=funcimpl.get_coeffs().end(); ++it) {
2600 const auto& key=it->first;
2601 const auto& node=it->second;
2602 if (not node.has_children()) to_be_erased.push_back(key);
2603 }
2604 return to_be_erased;
2605 };
2606
2607 FunctionImpl<T,LDIM>& f_nc=const_cast<FunctionImpl<T,LDIM>&>(*f.get_impl());
2608 for (auto& key : erase_list(f_nc)) f_nc.get_coeffs().erase(key);
2609 for (auto& g : vg) {
2610 FunctionImpl<R,KDIM>& g_nc=const_cast<FunctionImpl<R,KDIM>&>(*g.get_impl());
2611 for (auto& key : erase_list(g_nc)) g_nc.get_coeffs().erase(key);
2612 }
2613 world.gop.fence();
2615 f_nc.reconstruct(false);
2616 world.gop.fence();
2617
2618 }
2619
2620 return result;
2621 }
2622
2623
2624 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2625
2626 /// syntax similar to the inner product in tensor.h
2627 /// e.g result=inner<3>(f,g),{0},{1}) : r(x,y) = int f(x1,x) g(y,x1) dx1
2628 /// @param[in] task 0: everything, 1; prepare only (fence), 2: work only (no fence), 3: finalize only (fence)
2629 template<std::size_t NDIM, typename T, std::size_t LDIM, typename R, std::size_t KDIM,
2630 std::size_t CDIM = (KDIM + LDIM - NDIM) / 2>
2632 innerXX(const Function<T, LDIM>& f, const Function<R, KDIM>& g, const std::array<int, CDIM> v1,
2633 const std::array<int, CDIM> v2, int task=0) {
2634 return innerXX<NDIM,T,LDIM,R,KDIM>(f,std::vector<Function<R,KDIM>>({g}),v1,v2,task)[0];
2635 }
2636
2637 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2638
2639 /// syntax similar to the inner product in tensor.h
2640 /// e.g result=inner<3>(f,g),{0},{1}) : r(x,y) = int f(x1,x) g(y,x1) dx1
2641 template <typename T, std::size_t LDIM, typename R, std::size_t KDIM>
2643 inner(const Function<T,LDIM>& f, const Function<R,KDIM>& g, const std::tuple<int> v1, const std::tuple<int> v2) {
2644 return innerXX<KDIM+LDIM-2>(f,g,
2645 std::array<int,1>({std::get<0>(v1)}),
2646 std::array<int,1>({std::get<0>(v2)}));
2647 }
2648
2649 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2650
2651 /// syntax similar to the inner product in tensor.h
2652 /// e.g result=inner<3>(f,g),{0,1},{1,2}) : r(y) = int f(x1,x2) g(y,x1,x2) dx1 dx2
2653 template <typename T, std::size_t LDIM, typename R, std::size_t KDIM>
2655 inner(const Function<T,LDIM>& f, const Function<R,KDIM>& g, const std::tuple<int,int> v1, const std::tuple<int,int> v2) {
2656 return innerXX<KDIM+LDIM-4>(f,g,
2657 std::array<int,2>({std::get<0>(v1),std::get<1>(v1)}),
2658 std::array<int,2>({std::get<0>(v2),std::get<1>(v2)}));
2659 }
2660
2661 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2662
2663 /// syntax similar to the inner product in tensor.h
2664 /// e.g result=inner<3>(f,g),{1},{2}) : r(x,y,z) = int f(x,x1) g(y,z,x1) dx1
2665 template <typename T, std::size_t LDIM, typename R, std::size_t KDIM>
2667 inner(const Function<T,LDIM>& f, const Function<R,KDIM>& g, const std::tuple<int,int,int> v1, const std::tuple<int,int,int> v2) {
2668 return innerXX<KDIM+LDIM-6>(f,g,
2669 std::array<int,3>({std::get<0>(v1),std::get<1>(v1),std::get<2>(v1)}),
2670 std::array<int,3>({std::get<0>(v2),std::get<1>(v2),std::get<2>(v2)}));
2671 }
2672
2673
2674
2675 /// Computes the scalar/inner product between an MRA function and an external functor
2676
2677 /// Currently this defaults to inner_adaptive, which might be more expensive
2678 /// than inner_ext since it loops over all leaf nodes. If you feel inner_ext
2679 /// is more efficient you need to call it directly
2680 /// @param[in] f MRA function
2681 /// @param[in] g functor
2682 /// @result inner(f,g)
2683 template <typename T, typename opT, std::size_t NDIM>
2684 TENSOR_RESULT_TYPE(T,typename opT::value_type) inner(const Function<T,NDIM>& f, const opT& g) {
2686 std::shared_ptr< FunctionFunctorInterface<double,3> > func(new opT(g));
2687 return f.inner_adaptive(func);
2688 }
2689
2690 /// Computes the scalar/inner product between an MRA function and an external functor
2691
2692 /// Currently this defaults to inner_adaptive, which might be more expensive
2693 /// than inner_ext since it loops over all leaf nodes. If you feel inner_ext
2694 /// is more efficient you need to call it directly
2695 /// @param[in] g functor
2696 /// @param[in] f MRA function
2697 /// @result inner(f,g)
2698 template <typename T, typename opT, std::size_t NDIM>
2699 TENSOR_RESULT_TYPE(T,typename opT::value_type) inner(const opT& g, const Function<T,NDIM>& f) {
2700 return inner(f,g);
2701 }
2702
2703 template <typename T, typename R, std::size_t NDIM>
2706 return (f*R(1.0)).add_scalar(r);
2707 }
2708
2709 template <typename T, typename R, std::size_t NDIM>
2712 return (f*R(1.0)).add_scalar(r);
2713 }
2714
2715 template <typename T, typename R, std::size_t NDIM>
2718 return (f*R(1.0)).add_scalar(-r);
2719 }
2720
2721 template <typename T, typename R, std::size_t NDIM>
2724 return (f*R(-1.0)).add_scalar(r);
2725 }
2726
2727 namespace detail {
2728 template <std::size_t NDIM>
2729 struct realop {
2730 typedef double resultT;
2732 return real(t);
2733 }
2734
2735 template <typename Archive> void serialize (Archive& ar) {}
2736 };
2737
2738 template <std::size_t NDIM>
2739 struct imagop {
2740 typedef double resultT;
2742 return imag(t);
2743 }
2744
2745 template <typename Archive> void serialize (Archive& ar) {}
2746 };
2747
2748 template <std::size_t NDIM>
2749 struct abssqop {
2750 typedef double resultT;
2752 Tensor<double> r = abs(t);
2753 return r.emul(r);
2754 }
2755
2756 template <typename Archive> void serialize (Archive& ar) {}
2757 };
2758
2759 template <std::size_t NDIM>
2760 struct absop {
2761 typedef double resultT;
2763 Tensor<double> r = abs(t);
2764 return r;
2765 }
2766
2767 template <typename Archive> void serialize (Archive& ar) {}
2768 };
2769
2770 }
2771
2772 /// Returns a new function that is the real part of the input
2773 template <std::size_t NDIM>
2777
2778 /// Returns a new function that is the real part of the input
2779 template <std::size_t NDIM>
2781 return copy(z);
2782 }
2783
2784 /// Returns a new function that is the imaginary part of the input
2785 template <std::size_t NDIM>
2789
2790
2791 /// Create a new function that is the square of f - global comm only if not reconstructed
2792 template <typename T, std::size_t NDIM>
2795 Function<T,NDIM> result = copy(f,true); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2796 return result.square(true); //fence); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2797 }
2798
2799 /// Create a new function that is the abs of f - global comm only if not reconstructed
2800 template <typename T, std::size_t NDIM>
2801 Function<T,NDIM> abs(const Function<T,NDIM>& f, bool fence=true) {
2803 Function<T,NDIM> result = copy(f,true); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2804 return result.abs(true); //fence); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2805 }
2806
2807 /// Create a new function that is the abs_square of f - global comm only if not reconstructed
2808 template <typename T, std::size_t NDIM>
2809 typename std::enable_if<!TensorTypeData<T>::iscomplex, Function<T,NDIM> >::type
2810 abs_square(const Function<T,NDIM>& f, bool fence=true) {
2812 Function<T,NDIM> result = copy(f,true); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2813 return result.abs_square(true); //fence); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2814 }
2815
2816 /// Create a new function that is the abs_square of f - global comm only if not reconstructed
2817 template <typename T, std::size_t NDIM>
2818 typename std::enable_if<TensorTypeData<T>::iscomplex, Function<typename Tensor<T>::scalar_type,NDIM> >::type
2819 abs_square(const Function<T,NDIM>& f, bool fence=true) {
2821 }
2822
2823 /// Returns a new function that is the square of the absolute value of the input
2824 template <std::size_t NDIM>
2828
2829 /// Returns a new function that is the absolute value of the input
2830 template <std::size_t NDIM>
2834
2835 /// get tree state of a function
2836
2837 /// there is a corresponding function in vmra.h
2838 /// @param[in] f function
2839 /// @return TreeState::unknown if the function is not initialized
2840 template <typename T, std::size_t NDIM>
2842 if (f.is_initialized()) return f.get_impl()->get_tree_state();
2843 return TreeState::unknown;
2844 }
2845
2846 /// change tree state of a function
2847
2848 /// there is a corresponding function in vmra.h
2849 /// return this for chaining
2850 /// @param[in] f function
2851 /// @param[in] finalstate the new state
2852 /// @return this in the requested state
2853 template <typename T, std::size_t NDIM>
2858
2859 template <typename R, std::size_t MDIM>
2861 f1.impl.swap(f2.impl);
2862 }
2863
2864}
2865
2866#include <madness/mra/funcplot.h>
2867
2868namespace madness {
2869 namespace archive {
2870 template <class archiveT, class T, std::size_t NDIM>
2872 static inline void load(const ParallelInputArchive<archiveT>& ar, Function<T,NDIM>& f) {
2873 f.load(*ar.get_world(), ar);
2874 }
2875 };
2876
2877 template <class archiveT, class T, std::size_t NDIM>
2879 static inline void store(const ParallelOutputArchive<archiveT>& ar, const Function<T,NDIM>& f) {
2880 f.store(ar);
2881 }
2882 };
2883 }
2884
2885 template <class T, std::size_t NDIM>
2886 void save(const Function<T,NDIM>& f, const std::string name) {
2888 ar2 & f;
2889 }
2890
2891 template <class T, std::size_t NDIM>
2892 void load(Function<T,NDIM>& f, const std::string name) {
2894 ar2 & f;
2895 }
2896
2897}
2898
2899namespace madness {
2900 // type traits to check if a template parameter is a Function
2901 template<typename>
2902 struct is_madness_function : std::false_type {};
2903
2904 template<typename T, std::size_t NDIM>
2905 struct is_madness_function<madness::Function<T, NDIM>> : std::true_type {};
2906
2907}
2908
2909
2910/* @} */
2911
2912#include <madness/mra/derivative.h>
2913#include <madness/mra/operator.h>
2915#include <madness/mra/vmra.h>
2916// #include <madness/mra/mraimpl.h> !!!!!!!!!!!!! NOOOOOOOOOOOOOOOOOOOOOOOOOOOOOOOOOO !!!!!!!!!!!!!!!!!!
2917
2918#endif // MADNESS_MRA_MRA_H__INCLUDED
double q(double t)
Definition DKops.h:18
This header should include pretty much everything needed for the parallel runtime.
long dim(int i) const
Returns the size of dimension i.
Definition basetensor.h:147
This class is used to specify boundary conditions for all operators.
Definition bc.h:72
CompositeFunctorInterface implements a wrapper of holding several functions and functors.
Definition function_interface.h:165
FunctionDefaults holds default paramaters as static class members.
Definition funcdefaults.h:100
static const double & get_thresh()
Returns the default threshold.
Definition funcdefaults.h:177
static std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > & get_pmap()
Returns the default process map that was last initialized via set_default_pmap()
Definition funcdefaults.h:390
FunctionFactory implements the named-parameter idiom for Function.
Definition function_factory.h:86
virtual FunctionFactory & thresh(double thresh)
Definition function_factory.h:198
virtual FunctionFactory & k(int k)
Definition function_factory.h:193
Abstract base class interface required for functors used as input to Functions.
Definition function_interface.h:68
FunctionImpl holds all Function state to facilitate shallow copy semantics.
Definition funcimpl.h:945
const dcT & get_coeffs() const
Definition mraimpl.h:343
bool is_on_demand() const
Definition mraimpl.h:283
FunctionNode holds the coefficients, etc., at each node of the 2^NDIM-tree.
Definition funcimpl.h:127
A multiresolution adaptive numerical function.
Definition mra.h:139
void print_tree_json(std::ostream &os=std::cout) const
Definition mra.h:936
Tensor< T > eval_cube(const Tensor< double > &cell, const std::vector< long > &npt, bool eval_refine=false) const
Evaluates a cube/slice of points (probably for plotting) ... collective but no fence necessary.
Definition mra.h:355
TENSOR_RESULT_TYPE(T, R) inner(const Function< R
Returns the inner product.
Function< T, NDIM > & map_and_mirror(const Function< T, NDIM > &f, const std::vector< long > &map, const std::vector< long > &mirror, bool fence)
This is replaced with mirror(map(f)) ... private.
Definition mra.h:1736
void gaxpy_ext(const Function< L, NDIM > &left, T(*f)(const coordT &), T alpha, T beta, double tol, bool fence=true) const
Definition mra.h:1452
T inner_adaptive(const std::shared_ptr< FunctionFunctorInterface< T, NDIM > > f, const bool leaf_refine=true) const
Definition mra.h:1437
void unaryop_coeff(const opT &op, bool fence=true)
Unary operation applied inplace to the coefficients.
Definition mra.h:986
bool is_compressed() const
Returns true if compressed, false otherwise. No communication.
Definition mra.h:487
Function< T, NDIM/2 > dirac_convolution(const bool fence=true) const
Definition mra.h:1516
return impl inner_local * g())
void set_impl(const Function< R, NDIM > &f, bool zero=true)
Replace current FunctionImpl with a new one using the same parameters & map as f.
Definition mra.h:687
bool autorefine() const
Returns value of autorefine flag. No communication.
Definition mra.h:593
TreeState operating_state
Definition mra.h:1378
Function< T, NDIM > & add_scalar(T t, bool fence=true)
Inplace add scalar. No communication except for optional fence.
Definition mra.h:1035
void print_size(const std::string name) const
print some info about this
Definition mra.h:538
Function< T, NDIM > & fill_tree(const Function< R, NDIM > &g, bool fence=true)
With this being an on-demand function, fill the MRA tree according to different criteria.
Definition mra.h:1222
void broaden(const BoundaryConditions< NDIM > &bc=FunctionDefaults< NDIM >::get_bc(), bool fence=true) const
Inplace broadens support in scaling function basis.
Definition mra.h:908
void print_info() const
Print a summary of the load balancing info.
Definition mra.h:952
Function< T, NDIM > & scale(const Q q, bool fence=true)
Inplace, scale the function by a constant. No communication except for optional fence.
Definition mra.h:1025
void norm_tree(bool fence=true) const
Initializes information about the function norm at all length scales.
Definition mra.h:773
void load(World &world, Archive &ar)
Replaces this function with one loaded from an archive using the default processor map.
Definition mra.h:1533
Function< T, NDIM > & abs_square(bool fence=true)
Returns *this for chaining.
Definition mra.h:1162
void replicate(const DistributionType type, bool fence=true) const
Definition mra.h:711
double norm2sq_local() const
Returns the square of the norm of the local function ... no communication.
Definition mra.h:748
IsSupported< TensorTypeData< Q >, Function< T, NDIM > >::type & operator*=(const Q q)
Inplace scaling by a constant.
Definition mra.h:1134
void sum_down(bool fence=true) const
Sums scaling coeffs down tree restoring state with coeffs only at leaves. Optional fence....
Definition mra.h:870
Level depthpt(const coordT &xuser) const
Definition mra.h:431
bool is_redundant() const
Returns true if redundant, false otherwise. No communication.
Definition mra.h:517
Function< T, NDIM > & operator+=(const Function< Q, NDIM > &other)
Inplace addition of functions in the wavelet basis.
Definition mra.h:1094
Function< T, NDIM > & operator=(const Function< T, NDIM > &f)
Assignment is shallow. No communication, works in either basis.
Definition mra.h:191
void set_autorefine(bool value, bool fence=true)
Sets the value of the autorefine flag. Optional global fence.
Definition mra.h:603
World & world() const
Returns the world.
Definition mra.h:693
T trace() const
Returns global value of int(f(x),x) ... global comm required.
Definition mra.h:1184
T typeT
Definition mra.h:157
friend void swap(Function< R, MDIM > &f1, Function< R, MDIM > &f2)
implements swap algorithm
Definition mra.h:2860
Function< T, NDIM > & fill_tree(const opT &op, bool fence=true)
With this being an on-demand function, fill the MRA tree according to different criteria.
Definition mra.h:1240
const Function< T, NDIM > & change_tree_state(const TreeState finalstate, bool fence=true) const
changes tree state to given state
Definition mra.h:857
Function< typename opT::resultT, NDIM > & unary_op_coeffs(const Function< Q, NDIM > &func, const opT &op, bool fence)
This is replaced with left*right ... private.
Definition mra.h:1574
void print_tree_graphviz(std::ostream &os=std::cout) const
Process 0 prints a graphviz-formatted output of all nodes in the tree (collective)
Definition mra.h:942
double norm2() const
Returns the 2-norm of the function ... global sum ... works in either basis.
Definition mra.h:760
Function< T, NDIM > & fill_cuspy_tree(const opT &op, const bool fence=true)
Definition mra.h:1263
void change_tensor_type(const TensorArgs &targs, bool fence=true)
change the tensor type of the coefficients in the FunctionNode
Definition mra.h:1566
const Function< T, NDIM > & refine(bool fence=true) const
Inplace autorefines the function using same test as for squaring.
Definition mra.h:902
TENSOR_RESULT_TYPE(T, R) inner_on_demand(const Function< R
Returns the inner product for one on-demand function.
T operator()(const coordT &xuser) const
Evaluates the function at a point in user coordinates. Collective operation.
Definition mra.h:396
MADNESS_ASSERT(g.is_compressed())
Function< T, NDIM > & square(bool fence=true)
Inplace squaring of function ... global comm only if not reconstructed.
Definition mra.h:1144
Function< T, NDIM > & fill_tree(bool fence=true)
With this being an on-demand function, fill the MRA tree according to different criteria.
Definition mra.h:1251
void verify_tree() const
Verifies the tree data structure ... global sync implied.
Definition mra.h:478
Future< Level > evaldepthpt(const coordT &xuser) const
Definition mra.h:277
void do_hartree_product(const std::vector< std::shared_ptr< FunctionImpl< T, LDIM > > > left, const std::vector< std::shared_ptr< FunctionImpl< T, KDIM > > > right)
perform the hartree product of f*g, invoked by result
Definition mra.h:1332
Function< T, NDIM > & mapdim(const Function< T, NDIM > &f, const std::vector< long > &map, bool fence)
This is replaced with mapdim(f) ... private.
Definition mra.h:1702
impl world gop fence()
int k() const
Returns the number of multiwavelets (k). No communication.
Definition mra.h:631
const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > & get_pmap() const
Returns a shared pointer to the process map.
Definition mra.h:701
double thresh() const
Returns value of truncation threshold. No communication.
Definition mra.h:612
Function< T, NDIM > & gaxpy_oop(T alpha, const Function< L, NDIM > &left, T beta, const Function< R, NDIM > &right, bool fence)
This is replaced with alpha*left + beta*right ... private.
Definition mra.h:1688
void set_thresh(double value, bool fence=true)
Sets the value of the truncation threshold. Optional global fence.
Definition mra.h:622
Function< T, NDIM > & multiop_values(const opT &op, const std::vector< Function< T, NDIM > > &vf)
This is replaced with op(vector of functions) ... private.
Definition mra.h:1596
void distribute(std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > newmap) const
distribute this function according to newmap
Definition mra.h:739
void unaryop(const opT &op, bool fence=true)
Inplace unary operation on function values.
Definition mra.h:976
const std::shared_ptr< FunctionImpl< T, NDIM > > & get_impl() const
Returns a shared-pointer to the implementation.
Definition mra.h:659
void standard(bool fence=true)
Converts the function standard compressed form. Possible non-blocking comm.
Definition mra.h:820
void unaryop_node(const opT &op, bool fence=true)
Unary operation applied inplace to the nodes.
Definition mra.h:996
Function< T, NDIM > & truncate(double tol=0.0, bool fence=true)
Truncate the function with optional fence. Compresses with fence if not compressed.
Definition mra.h:647
Function< T, NDIM > & fill_cuspy_tree(const bool fence=true)
Special refinement on 6D boxes where the electrons come close (meet)
Definition mra.h:1276
std::size_t size() const
Returns the number of coefficients in the function ... collective global sum.
Definition mra.h:578
bool is_on_demand() const
Definition mra.h:681
Function< T, NDIM > & operator-=(const Function< Q, NDIM > &other)
Inplace subtraction of functions in the wavelet basis.
Definition mra.h:1113
bool compressed
Definition mra.h:1198
const Function< T, NDIM > & reconstruct(bool fence=true) const
Reconstructs the function, transforming into scaling function basis. Possible non-blocking comm.
Definition mra.h:847
Function< T, NDIM > & abs(bool fence=true)
Returns *this for chaining.
Definition mra.h:1153
Function< T, NDIM > & reduce_rank(const double thresh=0.0, const bool fence=true)
reduce the rank of the coefficient tensors
Definition mra.h:1766
Function< T, NDIM > & mirror(const Function< T, NDIM > &f, const std::vector< long > &mirrormap, bool fence)
This is replaced with mirror(f) ... private.
Definition mra.h:1717
std::size_t max_nodes() const
Returns the max number of nodes on a processor.
Definition mra.h:563
void do_hartree_product(const std::vector< std::shared_ptr< FunctionImpl< T, LDIM > > > left, const std::vector< std::shared_ptr< FunctionImpl< T, KDIM > > > right, const opT *op)
perform the hartree product of f*g, invoked by result
Definition mra.h:1318
std::shared_ptr< FunctionImpl< T, NDIM > > impl
Definition mra.h:146
void replicate_on_hosts(bool fence=true) const
replicate this function, one copy per host
Definition mra.h:732
void replicate(bool fence=true) const
replicate this function, generating a unique pmap
Definition mra.h:722
T trace_local() const
Returns local contribution to int(f(x),x) ... no communication.
Definition mra.h:1175
void mul_on_demand(const Function< L, NDIM > &f, const Function< R, NDIM > &g, bool fence=true)
Same as operator* but with optional fence and no automatic reconstruction.
Definition mra.h:1659
Function< T, NDIM > & gaxpy(const T &alpha, const Function< Q, NDIM > &other, const R &beta, bool fence=true)
Inplace, general bi-linear operation in wavelet basis. No communication except for optional fence.
Definition mra.h:1055
Vector< double, NDIM > coordT
Type of vector holding coordinates.
Definition mra.h:156
void store(Archive &ar) const
Stores the function to an archive.
Definition mra.h:1554
std::size_t max_local_depth() const
Returns the maximum local depth of the function tree ... no communications.
Definition mra.h:555
void vtransform(const std::vector< Function< R, NDIM > > &v, const Tensor< Q > &c, std::vector< Function< T, NDIM > > &vresult, double tol, bool fence=true)
sparse transformation of a vector of functions ... private
Definition mra.h:1677
TENSOR_RESULT_TYPE(T, R) local
std::size_t max_depth() const
Returns the maximum depth of the function tree ... collective global sum.
Definition mra.h:547
Function()
Default constructor makes uninitialized function. No communication.
Definition mra.h:174
std::size_t tree_size() const
Returns the number of nodes in the function tree ... collective global sum.
Definition mra.h:531
Function< TENSOR_RESULT_TYPE(T, R), NDIM-LDIM > project_out(const Function< R, LDIM > &g, const int dim) const
project this on the low-dim function g: h(x) = <f(x,y) | g(y)>
Definition mra.h:1492
static std::vector< std::shared_ptr< FunctionImpl< Q, D > > > vimpl(const std::vector< Function< Q, D > > &v)
Returns vector of FunctionImpl pointers corresponding to vector of functions.
Definition mra.h:1587
FunctionImpl< T, NDIM > implT
Definition mra.h:153
void clear(bool fence=true)
Clears the function as if constructed uninitialized. Optional fence.
Definition mra.h:919
void refine_general(const opT &op, bool fence=true) const
Inplace autorefines the function. Optional fence. Possible non-blocking comm.
Definition mra.h:883
static void doconj(const Key< NDIM >, Tensor< T > &t)
Definition mra.h:1006
std::pair< bool, T > eval_local_only(const Vector< double, NDIM > &xuser, Level maxlevel) const
Evaluate function only if point is local returning (true,value); otherwise return (false,...
Definition mra.h:245
void set_functor(const std::shared_ptr< FunctionFunctorInterface< T, NDIM > > functor)
Replace the current functor with the provided new one.
Definition mra.h:676
bool impl_initialized() const
Definition mra.h:149
Function< T, NDIM > & fill_nuclear_cuspy_tree(const opT &op, const size_t particle, const bool fence=true)
Definition mra.h:1291
bool is_redundant_after_merge() const
Returns true if redundant_after_merge, false otherwise. No communication.
Definition mra.h:525
return local
Definition mra.h:1393
auto func
Definition mra.h:1471
~Function()
Destruction of any underlying implementation is deferred to next global fence.
Definition mra.h:198
Function(const Function< T, NDIM > &f)
Copy constructor is shallow. No communication, works in either basis.
Definition mra.h:185
void set_impl(const std::shared_ptr< FunctionImpl< T, NDIM > > &impl)
Replace current FunctionImpl with provided new one.
Definition mra.h:666
T operator()(double x, double y=0, double z=0, double xx=0, double yy=0, double zz=0) const
Evaluates the function at a point in user coordinates. Collective operation.
Definition mra.h:410
T inner_ext_local(const std::shared_ptr< FunctionFunctorInterface< T, NDIM > > f, const bool leaf_refine=true, const bool keep_redundant=false) const
Definition mra.h:1404
std::size_t min_nodes() const
Returns the min number of nodes on a processor.
Definition mra.h:570
constexpr std::size_t LDIM
Definition mra.h:1470
static constexpr std::size_t dimT
Definition mra.h:158
change_tree_state(operating_state, false)
bool is_nonstandard() const
Returns true if nonstandard-compressed, false otherwise. No communication.
Definition mra.h:509
void verify() const
Asserts that the function is initialized.
Definition mra.h:162
double err(const funcT &func) const
Returns an estimate of the difference ||this-func|| ... global sum performed.
Definition mra.h:465
T inner_ext(const std::shared_ptr< FunctionFunctorInterface< T, NDIM > > f, const bool leaf_refine=true, const bool keep_redundant=false) const
Definition mra.h:1420
void make_redundant(bool fence=true)
Converts the function to redundant form, i.e. sum coefficients on all levels.
Definition mra.h:832
Function< T, NDIM > & fill_nuclear_cuspy_tree(const size_t particle, const bool fence=true)
Special refinement on 6D boxes for the nuclear potentials (regularized with cusp, non-regularized wit...
Definition mra.h:1304
double check_symmetry() const
check symmetry of a function by computing the 2nd derivative
Definition mra.h:1752
void multi_to_multi_op_values(const opT &op, const std::vector< Function< T, NDIM > > &vin, std::vector< Function< T, NDIM > > &vout, const bool fence=true)
apply op on the input vector yielding an output vector of functions
Definition mra.h:1615
FunctionFactory< T, NDIM > factoryT
Definition mra.h:155
std::size_t size_local() const
Return the number of coefficients in the function on this processor.
Definition mra.h:585
const Function< T, NDIM > & compress(bool fence=true) const
Compresses the function, transforming into wavelet basis. Possible non-blocking comm.
Definition mra.h:793
bool is_initialized() const
Returns true if the function is initialized.
Definition mra.h:167
bool is_reconstructed() const
Returns true if reconstructed, false otherwise. No communication.
Definition mra.h:498
Function< T, NDIM > & chop_at_level(const int n, const bool fence=true)
remove all nodes with level higher than n
Definition mra.h:1774
MADNESS_ASSERT(is_compressed())
void vmulXX(const Function< L, NDIM > &left, const std::vector< Function< R, NDIM > > &right, std::vector< Function< T, NDIM > > &result, double tol, bool fence)
Multiplication of function * vector of functions using recursive algorithm of mulxx.
Definition mra.h:1636
double errsq_local(const funcT &func) const
Returns an estimate of the difference ||this-func||^2 from local data.
Definition mra.h:450
void make_nonstandard(bool keepleaves, bool fence=true) const
Compresses the function retaining scaling function coeffs. Possible non-blocking comm.
Definition mra.h:806
if(VERIFY_TREE) verify_tree()
Future< T > eval(const coordT &xuser) const
Evaluates the function at a point in user coordinates. Possible non-blocking comm.
Definition mra.h:211
void print_tree(std::ostream &os=std::cout) const
Process 0 prints a summary of all nodes in the tree (collective)
Definition mra.h:929
FunctionNode< T, NDIM > nodeT
Definition mra.h:154
bool redundant
Definition mra.h:1199
Future< long > evalR(const coordT &xuser) const
Evaluates the function rank at a point in user coordinates. Possible non-blocking comm.
Definition mra.h:314
Function(const factoryT &factory)
Constructor from FunctionFactory provides named parameter idiom. Possible non-blocking communication.
Definition mra.h:178
NDIM &g const
Definition mra.h:1196
void unaryop(T(*f)(T))
Inplace unary operation on function values.
Definition mra.h:967
Function< T, NDIM > conj(bool fence=true)
Inplace complex conjugate. No communication except for optional fence.
Definition mra.h:1014
A future is a possibly yet unevaluated value.
Definition future.h:369
remote_refT remote_ref(World &world) const
Returns a structure used to pass references to another process.
Definition future.h:671
Key is the index for a node of the 2^NDIM-tree.
Definition key.h:69
Definition leafop.h:412
Definition leafop.h:281
Traits class to specify support of numeric types.
Definition type_data.h:56
A tensor is a multidimensional array.
Definition tensor.h:317
Tensor< T > & emul(const Tensor< T > &t)
Inplace multiply by corresponding elements of argument Tensor.
Definition tensor.h:1799
Tensor< T > & conj()
Inplace complex conjugate.
Definition tensor.h:716
A simple, fixed dimension vector.
Definition vector.h:64
void erase(const keyT &key)
Erases entry from container (non-blocking comm if remote)
Definition worlddc.h:1392
Interface to be provided by any process map.
Definition worlddc.h:122
void fence(bool debug=false)
Synchronizes all processes in communicator AND globally ensures no pending AM or tasks.
Definition worldgop.cc:161
bool set_forbid_fence(bool value)
Set forbid_fence flag to new value and return old value.
Definition worldgop.h:676
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
ProcessID size() const
Returns the number of processes in this World (same as MPI_Comm_size()).
Definition world.h:330
WorldGopInterface & gop
Global operations.
Definition world.h:207
World * get_world() const
Returns a pointer to the world.
Definition parallel_archive.h:130
An archive for storing local or parallel data, wrapping a BinaryFstreamInputArchive.
Definition parallel_archive.h:366
An archive for storing local or parallel data wrapping a BinaryFstreamOutputArchive.
Definition parallel_archive.h:321
Objects that implement their own parallel archive interface should derive from this class.
Definition parallel_archive.h:58
static const double R
Definition csqrt.cc:46
Declaration and initialization of tree traversal functions and generic derivative.
double(* f1)(const coord_3d &)
Definition derivatives.cc:55
double(* f2)(const coord_3d &)
Definition derivatives.cc:56
const double delta
Definition dielectric_external_field.cc:119
Provides FunctionDefaults and utilities for coordinate transformation.
Provides FunctionCommonData, FunctionImpl and FunctionFactory.
Defines/implements plotting interface for functions.
Provides typedefs to hide use of templates and to increase interoperability.
auto T(World &world, response_space &f) -> response_space
Definition global_functions.cc:28
const double beta
Definition gygi_soltion.cc:62
static const double v
Definition hatom_sf_dirac.cc:20
Provides IndexIterator.
Tensor< double > op(const Tensor< double > &x)
Definition kain.cc:508
Multidimension Key for MRA tree and associated iterators.
Implements (2nd generation) static load/data balancing for functions.
#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
Header to declare stuff which has not yet found a home.
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
double abs(double x)
Definition complexfun.h:48
Function< double, NDIM > abssq(const Function< double_complex, NDIM > &z, bool fence=true)
Returns a new function that is the square of the absolute value of the input.
Definition mra.h:2825
DistributionType
some introspection of how data is distributed
Definition worlddc.h:81
@ NodeReplicated
even if there are several ranks per node
Definition worlddc.h:84
@ RankReplicated
replicate the container over all world ranks
Definition worlddc.h:83
Function< TENSOR_RESULT_TYPE(typename opT::opT, R), NDIM > apply_1d_realspace_push(const opT &op, const Function< R, NDIM > &f, int axis, bool fence=true)
Definition mra.h:2343
Function< TENSOR_RESULT_TYPE(L, R), NDIM > sub(const Function< L, NDIM > &left, const Function< R, NDIM > &right, bool fence=true)
Same as operator- but with optional fence and no automatic compression.
Definition mra.h:2065
Function< TENSOR_RESULT_TYPE(L, R), NDIM > binary_op(const Function< L, NDIM > &left, const Function< R, NDIM > &right, const opT &op, bool fence=true)
Generate new function = op(left,right) where op acts on the function values.
Definition mra.h:1860
Function< Q, NDIM > convert(const Function< T, NDIM > &f, bool fence=true)
Type conversion implies a deep copy. No communication except for optional fence.
Definition mra.h:2143
Function< TENSOR_RESULT_TYPE(Q, T), NDIM > mul(const Q alpha, const Function< T, NDIM > &f, bool fence=true)
Returns new function equal to alpha*f(x) with optional fence.
Definition mra.h:1795
std::enable_if_t< NDIM%2==0, Function< T, NDIM > > swap_particles(const Function< T, NDIM > &f)
swap particles 1 and 2
Definition mra.h:2414
TreeState
Definition funcdefaults.h:59
@ nonstandard_after_apply
s and d coeffs, state after operator application
Definition funcdefaults.h:64
@ redundant_after_merge
s coeffs everywhere, must be summed up to yield the result
Definition funcdefaults.h:66
@ reconstructed
s coeffs at the leaves only
Definition funcdefaults.h:60
@ nonstandard
s and d coeffs in internal nodes
Definition funcdefaults.h:62
@ unknown
Definition funcdefaults.h:68
@ 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
std::vector< Function< TENSOR_RESULT_TYPE(T, R), NDIM > > innerXX(const Function< T, LDIM > &f, const std::vector< Function< R, KDIM > > &vg, const std::array< int, CDIM > v1, const std::array< int, CDIM > v2, int task=0)
Computes the partial scalar/inner product between two functions, returns a low-dim function.
Definition mra.h:2542
std::vector< CCPairFunction< T, NDIM > > operator*(const double fac, const std::vector< CCPairFunction< T, NDIM > > &arg)
Definition ccpairfunction.h:1089
int Level
Definition key.h:58
TreeState get_tree_state(const Function< T, NDIM > &f)
get tree state of a function
Definition mra.h:2841
std::vector< CCPairFunction< T, NDIM > > operator-(const std::vector< CCPairFunction< T, NDIM > > c1, const std::vector< CCPairFunction< T, NDIM > > &c2)
Definition ccpairfunction.h:1060
Function< TENSOR_RESULT_TYPE(L, R), NDIM > mul_sparse(const Function< L, NDIM > &left, const Function< R, NDIM > &right, double tol, bool fence=true)
Sparse multiplication — left and right must be reconstructed and if tol!=0 have tree of norms already...
Definition mra.h:1836
std::string get_mra_data_dir()
Definition startup.cc:209
Function< T, NDIM > gaxpy_oop_reconstructed(const double alpha, const Function< T, NDIM > &left, const double beta, const Function< T, NDIM > &right, const bool fence=true)
Returns new function alpha*left + beta*right optional fence, having both addends reconstructed.
Definition mra.h:2029
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
response_space apply(World &world, std::vector< std::vector< std::shared_ptr< real_convolution_3d > > > &op, response_space &f)
Definition basic_operators.cc:43
@ TT_2D
Definition gentensor.h:120
@ TT_FULL
Definition gentensor.h:120
NDIM & f
Definition mra.h:2528
Function< TENSOR_RESULT_TYPE(L, R), NDIM > add(const Function< L, NDIM > &left, const Function< R, NDIM > &right, bool fence=true)
Same as operator+ but with optional fence and no automatic compression.
Definition mra.h:2021
Function< T, NDIM > symmetrize(const Function< T, NDIM > &f, const std::string symmetry, bool fence=true)
symmetrize a function
Definition mra.h:2440
Function< TENSOR_RESULT_TYPE(typename opT::opT, R), NDIM > apply_only(const opT &op, const Function< R, NDIM > &f, bool fence=true)
Apply operator ONLY in non-standard form - required other steps missing !!
Definition mra.h:2231
double inner(response_space &a, response_space &b)
Definition response_functions.h:640
double imag(double x)
Definition complexfun.h:56
Function< typename opT::resultT, NDIM > unary_op(const Function< Q, NDIM > &func, const opT &op, bool fence=true)
Out of place application of unary operation to function values with optional fence.
Definition mra.h:1874
void startup(World &world, int argc, char **argv, bool doprint=false, bool make_stdcout_nice_to_reals=true)
initialize the internal state of the MADmra library
Definition startup.cc:64
std::string type(const PairType &n)
Definition PNOParameters.h:18
static bool print_timings
Definition SCF.cc:106
std::vector< CCPairFunction< T, NDIM > > operator+(const std::vector< CCPairFunction< T, NDIM > > c1, const std::vector< CCPairFunction< T, NDIM > > &c2)
Definition ccpairfunction.h:1052
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
void load(Function< T, NDIM > &f, const std::string name)
Definition mra.h:2892
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
double real(double x)
Definition complexfun.h:52
static XNonlinearSolver< std::vector< Function< T, NDIM > >, T, vector_function_allocator< T, NDIM > > nonlinear_vector_solver(World &world, const long nvec)
Definition nonlinsol.h:371
std::string name(const FuncType &type, const int ex=-1)
Definition ccpairfunction.h:28
void save(const Function< T, NDIM > &f, const std::string name)
Definition mra.h:2886
Function< T, KDIM+LDIM > hartree_product(const std::vector< Function< T, KDIM > > &left, const std::vector< Function< T, LDIM > > &right)
Performs a Hartree/outer product on the two given low-dimensional function vectors.
Definition mra.h:1926
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 double d
Definition nonlinschro.cc:121
Implements most functionality of separated operators.
Implements ParallelInputArchive and ParallelOutputArchive for parallel serialization of data.
double Q(double a)
Definition relops.cc:20
static const double c
Definition relops.cc:10
static const double L
Definition rk.cc:46
Definition test_ar.cc:204
Definition test_dc.cc:47
Definition leafop.h:142
void serialize(Archive &ar)
Definition mra.h:963
T(* f)(T)
Definition mra.h:958
void operator()(const Key< NDIM > &key, Tensor< T > &t) const
Definition mra.h:960
SimpleUnaryOpWrapper(T(*f)(T))
Definition mra.h:959
void serialize(Archive &ar)
Definition mra.h:896
bool operator()(implT *impl, const Key< NDIM > &key, const nodeT &t) const
Definition mra.h:892
Definition type_data.h:146
Definition leafop.h:194
Definition leafop.h:62
TensorArgs holds the arguments for creating a LowRankTensor.
Definition gentensor.h:134
static void load(const ParallelInputArchive< archiveT > &ar, Function< T, NDIM > &f)
Definition mra.h:2872
Default load of an object via serialize(ar, t).
Definition archive.h:667
static void store(const ParallelOutputArchive< archiveT > &ar, const Function< T, NDIM > &f)
Definition mra.h:2879
Default store of an object via serialize(ar, t).
Definition archive.h:612
Definition mra.h:2760
Tensor< double > operator()(const Key< NDIM > &key, const Tensor< double_complex > &t) const
Definition mra.h:2762
double resultT
Definition mra.h:2761
void serialize(Archive &ar)
Definition mra.h:2767
Definition mra.h:2749
void serialize(Archive &ar)
Definition mra.h:2756
Tensor< double > operator()(const Key< NDIM > &key, const Tensor< double_complex > &t) const
Definition mra.h:2751
double resultT
Definition mra.h:2750
Definition mra.h:2739
Tensor< double > operator()(const Key< NDIM > &key, const Tensor< double_complex > &t) const
Definition mra.h:2741
void serialize(Archive &ar)
Definition mra.h:2745
double resultT
Definition mra.h:2740
Definition mra.h:2729
Tensor< double > operator()(const Key< NDIM > &key, const Tensor< double_complex > &t) const
Definition mra.h:2731
double resultT
Definition mra.h:2730
void serialize(Archive &ar)
Definition mra.h:2735
Definition mra.h:127
Definition funcimpl.h:610
returns true if the result of a hartree_product is a leaf node (compute norm & error)
Definition funcimpl.h:500
Definition mra.h:2902
Definition mra.h:112
Definition mra.h:115
Definition funcimpl.h:564
Definition lowrankfunction.h:336
double real(double a)
Definition tdse4.cc:78
Defines and implements most of Tensor.
#define UNARY_OPTIMIZED_ITERATOR(X, x, exp)
Definition tensor_macros.h:658
template Tensor< Tensor< double_complex >::scalar_type > imag(const Tensor< double_complex > &t)
AtomicInt sum
Definition test_atomicint.cc:46
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
static const double alpha
Definition testcosine.cc:10
constexpr std::size_t NDIM
Definition testgconv.cc:54
std::size_t axis
Definition testpdiff.cc:59
#define TENSOR_RESULT_TYPE(L, R)
This macro simplifies access to TensorResultType.
Definition type_data.h:205
Defines operations on vectors of Functions.
Implements WorldContainer.
#define PROFILE_FUNC
Definition worldprofile.h:209
#define PROFILE_MEMBER_FUNC(classname)
Definition worldprofile.h:210