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.
178 Function(const factoryT& factory)
179 : impl(new FunctionImpl<T,NDIM>(factory)) {
181 }
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 /// Evaluates the function at a point in user coordinates. Possible non-blocking comm.
201
202 /// Only the invoking process will receive the result via the future
203 /// though other processes may be involved in the evaluation.
204 ///
205 /// Throws if function is not initialized.
206 Future<T> eval(const coordT& xuser) const {
208 const double eps=1e-15;
209 verify();
211 coordT xsim;
212 user_to_sim(xuser,xsim);
213 // If on the boundary, move the point just inside the
214 // volume so that the evaluation logic does not fail
215 for (std::size_t d=0; d<NDIM; ++d) {
216 if (xsim[d] < -eps) {
217 MADNESS_EXCEPTION("eval: coordinate lower-bound error in dimension", d);
218 }
219 else if (xsim[d] < eps) {
220 xsim[d] = eps;
221 }
222
223 if (xsim[d] > 1.0+eps) {
224 MADNESS_EXCEPTION("eval: coordinate upper-bound error in dimension", d);
225 }
226 else if (xsim[d] > 1.0-eps) {
227 xsim[d] = 1.0-eps;
228 }
229 }
230
231 Future<T> result;
232 impl->eval(xsim, impl->key0(), result.remote_ref(impl->world));
233 return result;
234 }
235
236 /// Evaluate function only if point is local returning (true,value); otherwise return (false,0.0)
237
238 /// maxlevel is the maximum depth to search down to --- the max local depth can be
239 /// computed with max_local_depth();
240 std::pair<bool,T> eval_local_only(const Vector<double,NDIM>& xuser, Level maxlevel) const {
241 const double eps=1e-15;
242 verify();
244 coordT xsim;
245 user_to_sim(xuser,xsim);
246 // If on the boundary, move the point just inside the
247 // volume so that the evaluation logic does not fail
248 for (std::size_t d=0; d<NDIM; ++d) {
249 if (xsim[d] < -eps) {
250 MADNESS_EXCEPTION("eval: coordinate lower-bound error in dimension", d);
251 }
252 else if (xsim[d] < eps) {
253 xsim[d] = eps;
254 }
255
256 if (xsim[d] > 1.0+eps) {
257 MADNESS_EXCEPTION("eval: coordinate upper-bound error in dimension", d);
258 }
259 else if (xsim[d] > 1.0-eps) {
260 xsim[d] = 1.0-eps;
261 }
262 }
263 return impl->eval_local_only(xsim,maxlevel);
264 }
265
266 /// Only the invoking process will receive the result via the future
267 /// though other processes may be involved in the evaluation.
268 ///
269 /// Throws if function is not initialized.
270 ///
271 /// This function is a minimally-modified version of eval()
272 Future<Level> evaldepthpt(const coordT& xuser) const {
274 const double eps=1e-15;
275 verify();
277 coordT xsim;
278 user_to_sim(xuser,xsim);
279 // If on the boundary, move the point just inside the
280 // volume so that the evaluation logic does not fail
281 for (std::size_t d=0; d<NDIM; ++d) {
282 if (xsim[d] < -eps) {
283 MADNESS_EXCEPTION("eval: coordinate lower-bound error in dimension", d);
284 }
285 else if (xsim[d] < eps) {
286 xsim[d] = eps;
287 }
288
289 if (xsim[d] > 1.0+eps) {
290 MADNESS_EXCEPTION("eval: coordinate upper-bound error in dimension", d);
291 }
292 else if (xsim[d] > 1.0-eps) {
293 xsim[d] = 1.0-eps;
294 }
295 }
296
297 Future<Level> result;
298 impl->evaldepthpt(xsim, impl->key0(), result.remote_ref(impl->world));
299 return result;
300 }
301
302
303 /// Evaluates the function rank at a point in user coordinates. Possible non-blocking comm.
304
305 /// Only the invoking process will receive the result via the future
306 /// though other processes may be involved in the evaluation.
307 ///
308 /// Throws if function is not initialized.
309 Future<long> evalR(const coordT& xuser) const {
311 const double eps=1e-15;
312 verify();
314 coordT xsim;
315 user_to_sim(xuser,xsim);
316 // If on the boundary, move the point just inside the
317 // volume so that the evaluation logic does not fail
318 for (std::size_t d=0; d<NDIM; ++d) {
319 if (xsim[d] < -eps) {
320 MADNESS_EXCEPTION("eval: coordinate lower-bound error in dimension", d);
321 }
322 else if (xsim[d] < eps) {
323 xsim[d] = eps;
324 }
325
326 if (xsim[d] > 1.0+eps) {
327 MADNESS_EXCEPTION("eval: coordinate upper-bound error in dimension", d);
328 }
329 else if (xsim[d] > 1.0-eps) {
330 xsim[d] = 1.0-eps;
331 }
332 }
333
334 Future<long> result;
335 impl->evalR(xsim, impl->key0(), result.remote_ref(impl->world));
336 return result;
337 }
338
339 /// Evaluates a cube/slice of points (probably for plotting) ... collective but no fence necessary
340
341 /// All processes recieve the entire result (which is a rather severe limit
342 /// on the size of the cube that is possible).
343
344 /// Set eval_refine=true to return the refinment levels of
345 /// the given function.
346
347 /// @param[in] cell A Tensor describe the cube where the function to be evaluated in
348 /// @param[in] npt How many points to evaluate in each dimension
349 /// @param[in] eval_refine Wether to return the refinment levels of the given function
351 const std::vector<long>& npt,
352 bool eval_refine = false) const {
353 MADNESS_ASSERT(static_cast<std::size_t>(cell.dim(0))>=NDIM && cell.dim(1)==2 && npt.size()>=NDIM);
355 const double eps=1e-14;
356 verify();
357 reconstruct();
358 coordT simlo, simhi;
359 for (std::size_t d=0; d<NDIM; ++d) {
360 simlo[d] = cell(d,0);
361 simhi[d] = cell(d,1);
362 }
363 user_to_sim(simlo, simlo);
364 user_to_sim(simhi, simhi);
365
366 // Move the bounding box infintesimally inside dyadic
367 // points so that the evaluation logic does not fail
368 for (std::size_t d=0; d<NDIM; ++d) {
369 MADNESS_ASSERT(simhi[d] >= simlo[d]);
370 MADNESS_ASSERT(simlo[d] >= 0.0);
371 MADNESS_ASSERT(simhi[d] <= 1.0);
372
373 double delta = eps*(simhi[d]-simlo[d]);
374 simlo[d] += delta;
375 simhi[d] -= 2*delta; // deliberate asymmetry
376 }
377 return impl->eval_plot_cube(simlo, simhi, npt, eval_refine);
378 }
379
380
381 /// Evaluates the function at a point in user coordinates. Collective operation.
382
383 /// Throws if function is not initialized.
384 ///
385 /// This function calls eval, blocks until the result is
386 /// available and then broadcasts the result to everyone.
387 /// Therefore, if you are evaluating many points in parallel
388 /// it is \em vastly less efficient than calling eval
389 /// directly, saving the futures, and then forcing all of the
390 /// results.
391 T operator()(const coordT& xuser) const {
393 verify();
395 T result;
396 if (impl->world.rank() == 0) result = eval(xuser).get();
397 impl->world.gop.broadcast(result);
398 //impl->world.gop.fence();
399 return result;
400 }
401
402 /// Evaluates the function at a point in user coordinates. Collective operation.
403
404 /// See "operator()(const coordT& xuser)" for more info
405 T operator()(double x, double y=0, double z=0, double xx=0, double yy=0, double zz=0) const {
406 coordT r;
407 r[0] = x;
408 if (NDIM>=2) r[1] = y;
409 if (NDIM>=3) r[2] = z;
410 if (NDIM>=4) r[3] = xx;
411 if (NDIM>=5) r[4] = yy;
412 if (NDIM>=6) r[5] = zz;
413 return (*this)(r);
414 }
415
416 /// Throws if function is not initialized.
417 ///
418 /// This function mimics operator() by going through the
419 /// tree looking for the depth of the tree at the point.
420 /// It blocks until the result is
421 /// available and then broadcasts the result to everyone.
422 /// Therefore, if you are evaluating many points in parallel
423 /// it is \em vastly less efficient than calling evaldepthpt
424 /// directly, saving the futures, and then forcing all of the
425 /// results.
426 Level depthpt(const coordT& xuser) const {
428 verify();
430 Level result;
431 if (impl->world.rank() == 0) result = evaldepthpt(xuser).get();
432 impl->world.gop.broadcast(result);
433 //impl->world.gop.fence();
434 return result;
435 }
436
437 /// Returns an estimate of the difference ||this-func||^2 from local data
438
439 /// No communication is performed. If the function is not
440 /// reconstructed, it throws an exception. To get the global
441 /// value either do a global sum of the local values or call
442 /// errsq
443 /// @param[in] func Templated interface to the a user specified function
444 template <typename funcT>
445 double errsq_local(const funcT& func) const {
447 verify();
448 if (!is_reconstructed()) MADNESS_EXCEPTION("Function:errsq_local:not reconstructed",0);
449 return impl->errsq_local(func);
450 }
451
452
453 /// Returns an estimate of the difference ||this-func|| ... global sum performed
454
455 /// If the function is compressed, it is reconstructed first. For efficient use
456 /// especially with many functions, reconstruct them all first, and use errsq_local
457 /// instead so you can perform a global sum on all at the same time.
458 /// @param[in] func Templated interface to the a user specified function
459 template <typename funcT>
460 double err(const funcT& func) const {
462 verify();
466 double local = impl->errsq_local(func);
467 impl->world.gop.sum(local);
468 impl->world.gop.fence();
469 return sqrt(local);
470 }
471
472 /// Verifies the tree data structure ... global sync implied
473 void verify_tree() const {
475 if (impl) impl->verify_tree();
476 }
477
478
479 /// Returns true if compressed, false otherwise. No communication.
480
481 /// If the function is not initialized, returns false.
482 bool is_compressed() const {
484 if (impl)
485 return impl->is_compressed();
486 else
487 return false;
488 }
489
490 /// Returns true if reconstructed, false otherwise. No communication.
491
492 /// If the function is not initialized, returns false.
493 bool is_reconstructed() const {
495 if (impl)
496 return impl->is_reconstructed();
497 else
498 return false;
499 }
500
501 /// Returns true if nonstandard-compressed, false otherwise. No communication.
502
503 /// If the function is not initialized, returns false.
504 bool is_nonstandard() const {
506 return impl ? impl->is_nonstandard() : false;
507 }
508
509
510 /// Returns the number of nodes in the function tree ... collective global sum
511 std::size_t tree_size() const {
513 if (!impl) return 0;
514 return impl->tree_size();
515 }
516
517 /// print some info about this
518 void print_size(const std::string name) const {
519 if (!impl) {
520 print("function",name,"not assigned yet");
521 } else {
522 impl->print_size(name);
523 }
524 }
525
526 /// Returns the maximum depth of the function tree ... collective global sum
527 std::size_t max_depth() const {
529 if (!impl) return 0;
530 return impl->max_depth();
531 }
532
533
534 /// Returns the maximum local depth of the function tree ... no communications
535 std::size_t max_local_depth() const {
537 if (!impl) return 0;
538 return impl->max_local_depth();
539 }
540
541
542 /// Returns the max number of nodes on a processor
543 std::size_t max_nodes() const {
545 if (!impl) return 0;
546 return impl->max_nodes();
547 }
548
549 /// Returns the min number of nodes on a processor
550 std::size_t min_nodes() const {
552 if (!impl) return 0;
553 return impl->min_nodes();
554 }
555
556
557 /// Returns the number of coefficients in the function ... collective global sum
558 std::size_t size() const {
560 if (!impl) return 0;
561 return impl->size();
562 }
563
564 /// Return the number of coefficients in the function on this processor
565 std::size_t size_local() const {
567 if (!impl) return 0;
568 return impl->size_local();
569 }
570
571
572 /// Returns value of autorefine flag. No communication.
573 bool autorefine() const {
575 if (!impl) return true;
576 return impl->get_autorefine();
577 }
578
579
580 /// Sets the value of the autorefine flag. Optional global fence.
581
582 /// A fence is required to ensure consistent global state.
583 void set_autorefine(bool value, bool fence = true) {
585 verify();
586 impl->set_autorefine(value);
587 if (fence) impl->world.gop.fence();
588 }
589
590
591 /// Returns value of truncation threshold. No communication.
592 double thresh() const {
594 if (!impl) return 0.0;
595 return impl->get_thresh();
596 }
597
598
599 /// Sets the value of the truncation threshold. Optional global fence.
600
601 /// A fence is required to ensure consistent global state.
602 void set_thresh(double value, bool fence = true) {
604 verify();
605 impl->set_thresh(value);
606 if (fence) impl->world.gop.fence();
607 }
608
609
610 /// Returns the number of multiwavelets (k). No communication.
611 int k() const {
613 verify();
614 return impl->get_k();
615 }
616
617
618 /// Truncate the function with optional fence. Compresses with fence if not compressed.
619
620 /// If the truncation threshold is less than or equal to zero the default value
621 /// specified when the function was created is used.
622 /// If the function is not initialized, it just returns.
623 ///
624 /// Returns this for chaining.
625 /// @param[in] tol Tolerance for truncating the coefficients. Default 0.0 means use the implementation's member value \c thresh instead.
626 /// @param[in] fence Do fence
627 Function<T,NDIM>& truncate(double tol = 0.0, bool fence = true) {
629 if (!impl) return *this;
630 verify();
631// if (!is_compressed()) compress();
632 impl->truncate(tol,fence);
634 return *this;
635 }
636
637
638 /// Returns a shared-pointer to the implementation
639 const std::shared_ptr< FunctionImpl<T,NDIM> >& get_impl() const {
641 verify();
642 return impl;
643 }
644
645 /// Replace current FunctionImpl with provided new one
646 void set_impl(const std::shared_ptr< FunctionImpl<T,NDIM> >& impl) {
648 this->impl = impl;
649 }
650
651
652 /// Replace the current functor with the provided new one
653
654 /// presumably the new functor will be a CompositeFunctor, which will
655 /// change the behavior of the function: multiply the functor with the function
656 void set_functor(const std::shared_ptr<FunctionFunctorInterface<T, NDIM> > functor) {
657 this->impl->set_functor(functor);
658 print("set functor in mra.h");
659 }
660
661 bool is_on_demand() const {return this->impl->is_on_demand();}
662
663 /// Replace current FunctionImpl with a new one using the same parameters & map as f
664
665 /// If zero is true the function is initialized to zero, otherwise it is empty
666 template <typename R>
667 void set_impl(const Function<R,NDIM>& f, bool zero = true) {
668 impl = std::shared_ptr<implT>(new implT(*f.get_impl(), f.get_pmap(), zero));
669 if (zero) world().gop.fence();
670 }
671
672 /// Returns the world
673 World& world() const {
675 verify();
676 return impl->world;
677 }
678
679
680 /// Returns a shared pointer to the process map
681 const std::shared_ptr< WorldDCPmapInterface< Key<NDIM> > >& get_pmap() const {
683 verify();
684 return impl->get_pmap();
685 }
686
687 /// replicate this function, generating a unique pmap
688 void replicate(bool fence=true) const {
689 verify();
690 impl->replicate(fence);
691 }
692
693 /// distribute this function according to newmap
694 void distribute(std::shared_ptr< WorldDCPmapInterface< Key<NDIM> > > newmap) const {
695 verify();
696 impl->distribute(newmap);
697 }
698
699
700 /// Returns the square of the norm of the local function ... no communication
701
702 /// Works in either basis
703 double norm2sq_local() const {
705 verify();
706 return impl->norm2sq_local();
707 }
708
709
710 /// Returns the 2-norm of the function ... global sum ... works in either basis
711
712 /// See comments for err() w.r.t. applying to many functions.
713 double norm2() const {
715 verify();
717 double local = impl->norm2sq_local();
718
719 impl->world.gop.sum(local);
720 impl->world.gop.fence();
721 return sqrt(local);
722 }
723
724
725 /// Initializes information about the function norm at all length scales
726 void norm_tree(bool fence = true) const {
728 verify();
731 const_cast<Function<T,NDIM>*>(this)->impl->norm_tree(fence);
732 }
733
734
735 /// Compresses the function, transforming into wavelet basis. Possible non-blocking comm.
736
737 /// By default fence=true meaning that this operation completes before returning,
738 /// otherwise if fence=false it returns without fencing and the user must invoke
739 /// world.gop.fence() to assure global completion before using the function
740 /// for other purposes.
741 ///
742 /// Noop if already compressed or if not initialized.
743 ///
744 /// Since reconstruction/compression do not discard information we define them
745 /// as const ... "logical constness" not "bitwise constness".
746 const Function<T,NDIM>& compress(bool fence = true) const {
748 }
749
750
751 /// Compresses the function retaining scaling function coeffs. Possible non-blocking comm.
752
753 /// By default fence=true meaning that this operation completes before returning,
754 /// otherwise if fence=false it returns without fencing and the user must invoke
755 /// world.gop.fence() to assure global completion before using the function
756 /// for other purposes.
757 ///
758 /// Noop if already compressed or if not initialized.
759 void make_nonstandard(bool keepleaves, bool fence=true) const {
761 if (keepleaves) newstate=nonstandard_with_leaves;
762 change_tree_state(newstate,fence);
763 }
764
765 /// Converts the function standard compressed form. Possible non-blocking comm.
766
767 /// By default fence=true meaning that this operation completes before returning,
768 /// otherwise if fence=false it returns without fencing and the user must invoke
769 /// world.gop.fence() to assure global completion before using the function
770 /// for other purposes.
771 ///
772 /// Must be already compressed.
773 void standard(bool fence = true) {
775 }
776
777 /// Converts the function to redundant form, i.e. sum coefficients on all levels
778
779 /// By default fence=true meaning that this operation completes before returning,
780 /// otherwise if fence=false it returns without fencing and the user must invoke
781 /// world.gop.fence() to assure global completion before using the function
782 /// for other purposes.
783 ///
784 /// Must be already compressed.
785 void make_redundant(bool fence = true) {
787 }
788
789 /// Reconstructs the function, transforming into scaling function basis. Possible non-blocking comm.
790
791 /// By default fence=true meaning that this operation completes before returning,
792 /// otherwise if fence=false it returns without fencing and the user must invoke
793 /// world.gop.fence() to assure global completion before using the function
794 /// for other purposes.
795 ///
796 /// Noop if already reconstructed or if not initialized.
797 ///
798 /// Since reconstruction/compression do not discard information we define them
799 /// as const ... "logical constness" not "bitwise constness".
800 const Function<T,NDIM>& reconstruct(bool fence = true) const {
802 }
803
804 /// changes tree state to given state
805
806 /// Since reconstruction/compression do not discard information we define them
807 /// as const ... "logical constness" not "bitwise constness".
808 /// @param[in] finalstate The final state of the tree
809 /// @param[in] fence Fence after the operation (might not be respected!!!)
810 const Function<T,NDIM>& change_tree_state(const TreeState finalstate, bool fence = true) const {
812 if (not impl) return *this;
813 TreeState current_state = impl->get_tree_state();
814 if (finalstate == current_state) return *this;
815 MADNESS_CHECK_THROW(current_state != TreeState::unknown, "unknown tree state");
816
817 impl->change_tree_state(finalstate, fence);
818 if (fence && VERIFY_TREE) verify_tree();
819 return *this;
820 }
821
822 /// Sums scaling coeffs down tree restoring state with coeffs only at leaves. Optional fence. Possible non-blocking comm.
823 void sum_down(bool fence = true) const {
825 verify();
826 MADNESS_CHECK_THROW(impl->get_tree_state()==redundant_after_merge, "sum_down requires a redundant_after_merge state");
827 const_cast<Function<T,NDIM>*>(this)->impl->sum_down(fence);
828 const_cast<Function<T,NDIM>*>(this)->impl->set_tree_state(reconstructed);
829
830 if (fence && VERIFY_TREE) verify_tree(); // Must be after in case nonstandard
831 }
832
833
834 /// Inplace autorefines the function. Optional fence. Possible non-blocking comm.
835 template <typename opT>
836 void refine_general(const opT& op, bool fence = true) const {
838 verify();
840 impl->refine(op, fence);
841 }
842
843
845 bool operator()(implT* impl, const Key<NDIM>& key, const nodeT& t) const {
846 return impl->autorefine_square_test(key, t);
847 }
848
849 template <typename Archive> void serialize (Archive& ar) {}
850 };
851
852 /// Inplace autorefines the function using same test as for squaring.
853
854 /// return this for chaining
855 const Function<T,NDIM>& refine(bool fence = true) const {
857 return *this;
858 }
859
860 /// Inplace broadens support in scaling function basis
862 bool fence = true) const {
863 verify();
864 reconstruct();
865 impl->broaden(bc.is_periodic(), fence);
866 }
867
868
869 /// Clears the function as if constructed uninitialized. Optional fence.
870
871 /// Any underlying data will not be freed until the next global fence.
872 void clear(bool fence = true) {
874 if (impl) {
875 World& world = impl->world;
876 impl.reset();
877 if (fence) world.gop.fence();
878 }
879 }
880
881 /// Process 0 prints a summary of all nodes in the tree (collective)
882 void print_tree(std::ostream& os = std::cout) const {
884 if (impl) impl->print_tree(os);
885 }
886
887 /// same as print_tree() but produces JSON-formatted string
888 /// @warning enclose the result in braces to make it a valid JSON object
889 void print_tree_json(std::ostream& os = std::cout) const {
891 if (impl) impl->print_tree_json(os);
892 }
893
894 /// Process 0 prints a graphviz-formatted output of all nodes in the tree (collective)
895 void print_tree_graphviz(std::ostream& os = std::cout) const {
897 os << "digraph G {" << std::endl;
898 if (impl) impl->print_tree_graphviz(os);
899 os << "}" << std::endl;
900 }
901
902 /// Print a summary of the load balancing info
903
904 /// This is serial and VERY expensive
905 void print_info() const {
907 if (impl) impl->print_info();
908 }
909
911 T (*f)(T);
913 void operator()(const Key<NDIM>& key, Tensor<T>& t) const {
914 UNARY_OPTIMIZED_ITERATOR(T, t, *_p0 = f(*_p0));
915 }
916 template <typename Archive> void serialize(Archive& ar) {}
917 };
918
919 /// Inplace unary operation on function values
920 void unaryop(T (*f)(T)) {
921 // Must fence here due to temporary object on stack
922 // stopping us returning before complete
924 }
925
926
927 /// Inplace unary operation on function values
928 template <typename opT>
929 void unaryop(const opT& op, bool fence=true) {
931 verify();
932 reconstruct();
933 impl->unary_op_value_inplace(op, fence);
934 }
935
936
937 /// Unary operation applied inplace to the coefficients
938 template <typename opT>
939 void unaryop_coeff(const opT& op,
940 bool fence = true) {
942 verify();
943 impl->unary_op_coeff_inplace(op, fence);
944 }
945
946
947 /// Unary operation applied inplace to the nodes
948 template <typename opT>
949 void unaryop_node(const opT& op,
950 bool fence = true) {
952 verify();
953 impl->unary_op_node_inplace(op, fence);
954 }
955
956
957
958
959 static void doconj(const Key<NDIM>, Tensor<T>& t) {
961 t.conj();
962 }
963
964 /// Inplace complex conjugate. No communication except for optional fence.
965
966 /// Returns this for chaining. Works in either basis.
972
973
974 /// Inplace, scale the function by a constant. No communication except for optional fence.
975
976 /// Works in either basis. Returns reference to this for chaining.
977 template <typename Q>
978 Function<T,NDIM>& scale(const Q q, bool fence=true) {
980 verify();
982 impl->scale_inplace(q,fence);
983 return *this;
984 }
985
986
987 /// Inplace add scalar. No communication except for optional fence.
990 verify();
992 impl->add_scalar_inplace(t,fence);
993 return *this;
994 }
995
996
997 /// Inplace, general bi-linear operation in wavelet basis. No communication except for optional fence.
998
999 /// If the functions are not in the wavelet basis an exception is thrown since this routine
1000 /// is intended to be fast and unexpected compression is assumed to be a performance bug.
1001 ///
1002 /// Returns this for chaining, can be in states compressed of redundant_after_merge.
1003 ///
1004 /// this <-- this*alpha + other*beta
1005 template <typename Q, typename R>
1007 const Function<Q,NDIM>& other, const R& beta, bool fence=true) {
1009 verify();
1010 other.verify();
1011 MADNESS_CHECK_THROW(impl->get_tree_state() == other.get_impl()->get_tree_state(),
1012 "gaxpy requires both functions to be in the same tree state");
1013 MADNESS_CHECK_THROW(impl->get_tree_state()==reconstructed or impl->get_tree_state()==compressed,
1014 "gaxpy requires the tree state to be reconstructed or compressed");
1015
1016 bool same_world=this->world().id()==other.world().id();
1017 MADNESS_CHECK(same_world or is_compressed());
1018
1019 if (not same_world) {
1020 impl->gaxpy_inplace(alpha,*other.get_impl(),beta,fence);
1021 } else {
1022 if (is_compressed()) impl->gaxpy_inplace(alpha, *other.get_impl(), beta, fence);
1023 if (is_reconstructed()) impl->gaxpy_inplace_reconstructed(alpha,*other.get_impl(),beta,fence);
1024 }
1025 return *this;
1026 }
1027
1028
1029 /// Inplace addition of functions in the wavelet basis
1030
1031 /// Using operator notation forces a global fence after every operation.
1032 /// Functions don't need to be compressed, it's the caller's responsibility
1033 /// to choose an appropriate state with performance, usually compressed for 3d,
1034 /// reconstructed for 6d)
1035 template <typename Q>
1038 if (NDIM<=3) {
1039 compress();
1040 other.compress();
1041 } else {
1042 reconstruct();
1043 other.reconstruct();
1044 }
1045 MADNESS_ASSERT(impl->get_tree_state() == other.get_impl()->get_tree_state());
1046 if (VERIFY_TREE) verify_tree();
1047 if (VERIFY_TREE) other.verify_tree();
1048 return gaxpy(T(1.0), other, Q(1.0), true);
1049 }
1050
1051
1052 /// Inplace subtraction of functions in the wavelet basis
1053
1054 /// Using operator notation forces a global fence after every operation
1055 template <typename Q>
1058 if (NDIM<=3) {
1059 compress();
1060 other.compress();
1061 } else {
1062 reconstruct();
1063 other.reconstruct();
1064 }
1065 MADNESS_ASSERT(impl->get_tree_state() == other.get_impl()->get_tree_state());
1066 if (VERIFY_TREE) verify_tree();
1067 if (VERIFY_TREE) other.verify_tree();
1068 return gaxpy(T(1.0), other, Q(-1.0), true);
1069 }
1070
1071
1072 /// Inplace scaling by a constant
1073
1074 /// Using operator notation forces a global fence after every operation
1075 template <typename Q>
1077 operator*=(const Q q) {
1079 scale(q,true);
1080 return *this;
1081 }
1082
1083
1084 /// Inplace squaring of function ... global comm only if not reconstructed
1085
1086 /// Returns *this for chaining.
1089 if (!is_reconstructed()) reconstruct();
1090 if (VERIFY_TREE) verify_tree();
1091 impl->square_inplace(fence);
1092 return *this;
1093 }
1094
1095 /// Returns *this for chaining.
1098 if (!is_reconstructed()) reconstruct();
1099 if (VERIFY_TREE) verify_tree();
1100 impl->abs_inplace(fence);
1101 return *this;
1102 }
1103
1104 /// Returns *this for chaining.
1107 if (!is_reconstructed()) reconstruct();
1108 if (VERIFY_TREE) verify_tree();
1109 impl->abs_square_inplace(fence);
1110 return *this;
1111 }
1112
1113 /// Returns local contribution to \c int(f(x),x) ... no communication
1114
1115 /// In the wavelet basis this is just the coefficient of the first scaling
1116 /// function which is a constant. In the scaling function basis we
1117 /// must add up contributions from each box.
1118 T trace_local() const {
1120 if (!impl) return 0.0;
1121 if (VERIFY_TREE) verify_tree();
1122 return impl->trace_local();
1123 }
1124
1125
1126 /// Returns global value of \c int(f(x),x) ... global comm required
1127 T trace() const {
1129 if (!impl) return 0.0;
1130 T sum = impl->trace_local();
1131 impl->world.gop.sum(sum);
1132 impl->world.gop.fence();
1133 return sum;
1134 }
1135
1136
1137 /// Returns local part of inner product ... throws if both not compressed
1138 template <typename R>
1139 TENSOR_RESULT_TYPE(T,R) inner_local(const Function<R,NDIM>& g) const {
1142 MADNESS_ASSERT(g.is_compressed());
1144 if (VERIFY_TREE) g.verify_tree();
1145 return impl->inner_local(*(g.get_impl()));
1146 }
1147
1148
1149 /// With this being an on-demand function, fill the MRA tree according to different criteria
1150
1151 /// @param[in] g the function after which the MRA structure is modeled (any basis works)
1152 template<typename R>
1154 MADNESS_ASSERT(g.is_initialized());
1156
1157 // clear what we have
1158 impl->get_coeffs().clear();
1159
1160 //leaf_op<T,NDIM> gnode_is_leaf(g.get_impl().get());
1161 Leaf_op_other<T,NDIM> gnode_is_leaf(g.get_impl().get());
1162 impl->make_Vphi(gnode_is_leaf,fence);
1163 return *this;
1164
1165 }
1166
1167 /// With this being an on-demand function, fill the MRA tree according to different criteria
1168
1169 /// @param[in] op the convolution operator for screening
1170 template<typename opT>
1171 Function<T,NDIM>& fill_tree(const opT& op, bool fence=true) {
1173 // clear what we have
1174 impl->get_coeffs().clear();
1177 impl ->make_Vphi(leaf_op,fence);
1178 return *this;
1179 }
1180
1181 /// With this being an on-demand function, fill the MRA tree according to different criteria
1184 // clear what we have
1185 impl->get_coeffs().clear();
1187 impl->make_Vphi(leaf_op,fence);
1188 return *this;
1189 }
1190
1191 /// Special refinement on 6D boxes where the electrons come close (meet)
1192 /// @param[in] op the convolution operator for screening
1193 template<typename opT>
1194 Function<T,NDIM>& fill_cuspy_tree(const opT& op,const bool fence=true){
1196 // clear what we have
1197 impl->get_coeffs().clear();
1199
1201 impl ->make_Vphi(leaf_op,fence);
1202
1203 return *this;
1204 }
1205
1206 /// Special refinement on 6D boxes where the electrons come close (meet)
1209 // clear what we have
1210 impl->get_coeffs().clear();
1212
1214 impl ->make_Vphi(leaf_op,fence);
1215
1216 return *this;
1217 }
1218
1219 /// Special refinement on 6D boxes for the nuclear potentials (regularized with cusp, non-regularized with singularity)
1220 /// @param[in] op the convolution operator for screening
1221 template<typename opT>
1222 Function<T,NDIM>& fill_nuclear_cuspy_tree(const opT& op,const size_t particle,const bool fence=true){
1224 // clear what we have
1225 impl->get_coeffs().clear();
1227
1229 impl ->make_Vphi(leaf_op,fence);
1230
1231 return *this;
1232 }
1233
1234 /// Special refinement on 6D boxes for the nuclear potentials (regularized with cusp, non-regularized with singularity)
1237 // clear what we have
1238 impl->get_coeffs().clear();
1240
1242 impl ->make_Vphi(leaf_op,fence);
1243
1244 return *this;
1245 }
1246
1247 /// perform the hartree product of f*g, invoked by result
1248 template<size_t LDIM, size_t KDIM, typename opT>
1249 void do_hartree_product(const std::vector<std::shared_ptr<FunctionImpl<T,LDIM>>> left,
1250 const std::vector<std::shared_ptr<FunctionImpl<T,KDIM>>> right,
1251 const opT* op) {
1252
1253 // get the right leaf operator
1255 impl->hartree_product(left,right,leaf_op,true);
1256 impl->finalize_sum();
1257// this->truncate();
1258
1259 }
1260
1261 /// perform the hartree product of f*g, invoked by result
1262 template<size_t LDIM, size_t KDIM>
1263 void do_hartree_product(const std::vector<std::shared_ptr<FunctionImpl<T,LDIM>>> left,
1264 const std::vector<std::shared_ptr<FunctionImpl<T,KDIM>>> right) {
1265
1266// hartree_leaf_op<T,KDIM+LDIM> leaf_op(impl.get(),cdata.s0);
1268 impl->hartree_product(left,right,leaf_op,true);
1269 impl->finalize_sum();
1270// this->truncate();
1271
1272 }
1273
1274 /// Returns the inner product
1275
1276 /// Not efficient for computing multiple inner products
1277 /// @param[in] g Function, optionally on-demand
1278 template <typename R>
1281
1282 // fast return if possible
1283 if (not this->is_initialized()) return 0.0;
1284 if (not g.is_initialized()) return 0.0;
1285
1286 // if this and g are the same, use norm2()
1287 if constexpr (std::is_same_v<T,R>) {
1288 if (this->get_impl() == g.get_impl()) {
1289 TreeState state = this->get_impl()->get_tree_state();
1290 if (not(state == reconstructed or state == compressed))
1292 double norm = this->norm2();
1293 return norm * norm;
1294 }
1295 }
1296
1297 // do it case-by-case
1298 if constexpr (std::is_same_v<R,T>) {
1299 if (this->is_on_demand())
1300 return g.inner_on_demand(*this);
1301 if (g.is_on_demand())
1302 return this->inner_on_demand(g);
1303 }
1304
1306 if (VERIFY_TREE) g.verify_tree();
1307
1308 // compression is more efficient for 3D
1309 TreeState state=this->get_impl()->get_tree_state();
1310 TreeState gstate=g.get_impl()->get_tree_state();
1311 if (NDIM<=3) {
1313 g.change_tree_state(compressed,false);
1314 impl->world.gop.fence();
1315 }
1316
1317 if (this->is_compressed() and g.is_compressed()) {
1318 } else {
1320 g.change_tree_state(redundant,false);
1321 impl->world.gop.fence();
1322 }
1323
1324
1325 TENSOR_RESULT_TYPE(T,R) local = impl->inner_local(*g.get_impl());
1326 impl->world.gop.sum(local);
1327 impl->world.gop.fence();
1328
1329 // restore state
1331 g.change_tree_state(gstate,false);
1332 impl->world.gop.fence();
1333
1334 return local;
1335 }
1336
1337 /// Return the local part of inner product with external function ... no communication.
1338 /// If you are going to be doing a bunch of inner_ext calls, set
1339 /// keep_redundant to true and then manually undo_redundant when you
1340 /// are finished.
1341 /// @param[in] f Pointer to function of type T that take coordT arguments. This is the externally provided function
1342 /// @param[in] leaf_refine boolean switch to turn on/off refinement past leaf nodes
1343 /// @param[in] keep_redundant boolean switch to turn on/off undo_redundant
1344 /// @return Returns local part of the inner product, i.e. over the domain of all function nodes on this compute node.
1345 T inner_ext_local(const std::shared_ptr< FunctionFunctorInterface<T,NDIM> > f, const bool leaf_refine=true, const bool keep_redundant=false) const {
1348 T local = impl->inner_ext_local(f, leaf_refine);
1349 if (not keep_redundant) change_tree_state(reconstructed);
1350 return local;
1351 }
1352
1353 /// Return the inner product with external function ... requires communication.
1354 /// If you are going to be doing a bunch of inner_ext calls, set
1355 /// keep_redundant to true and then manually undo_redundant when you
1356 /// are finished.
1357 /// @param[in] f Reference to FunctionFunctorInterface. This is the externally provided function
1358 /// @param[in] leaf_refine boolean switch to turn on/off refinement past leaf nodes
1359 /// @param[in] keep_redundant boolean switch to turn on/off undo_redundant
1360 /// @return Returns the inner product
1361 T inner_ext(const std::shared_ptr< FunctionFunctorInterface<T,NDIM> > f, const bool leaf_refine=true, const bool keep_redundant=false) const {
1364 T local = impl->inner_ext_local(f, leaf_refine);
1365 impl->world.gop.sum(local);
1366 impl->world.gop.fence();
1367 if (not keep_redundant) change_tree_state(reconstructed);
1368 return local;
1369 }
1370
1371 /// Return the inner product with external function ... requires communication.
1372 /// If you are going to be doing a bunch of inner_ext calls, set
1373 /// keep_redundant to true and then manually undo_redundant when you
1374 /// are finished.
1375 /// @param[in] f Reference to FunctionFunctorInterface. This is the externally provided function
1376 /// @param[in] leaf_refine boolean switch to turn on/off refinement past leaf nodes
1377 /// @return Returns the inner product
1379 const bool leaf_refine=true) const {
1381 reconstruct();
1382 T local = impl->inner_adaptive_local(f, leaf_refine);
1383 impl->world.gop.sum(local);
1384 impl->world.gop.fence();
1385 return local;
1386 }
1387
1388 /// Return the local part of gaxpy with external function, this*alpha + f*beta ... no communication.
1389 /// @param[in] alpha prefactor for this Function
1390 /// @param[in] f Pointer to function of type T that take coordT arguments. This is the externally provided function
1391 /// @param[in] beta prefactor for f
1392 template <typename L>
1393 void gaxpy_ext(const Function<L,NDIM>& left, T (*f)(const coordT&), T alpha, T beta, double tol, bool fence=true) const {
1395 if (!left.is_reconstructed()) left.reconstruct();
1396 impl->gaxpy_ext(left.get_impl().get(), f, alpha, beta, tol, fence);
1397 }
1398
1399 /// Returns the inner product for one on-demand function
1400
1401 /// It does work, but it might not give you the precision you expect.
1402 /// The assumption is that the function g returns proper sum
1403 /// coefficients on the MRA tree of this. This might not be the case if
1404 /// g is constructed with an implicit multiplication, e.g.
1405 /// result = <this|g>, with g = 1/r12 | gg>
1406 /// @param[in] g on-demand function
1407 template<typename R>
1408 TENSOR_RESULT_TYPE(T, R) inner_on_demand(const Function<R, NDIM>& g) const {
1409 MADNESS_ASSERT(g.is_on_demand() and (not this->is_on_demand()));
1410
1411 constexpr std::size_t LDIM=std::max(NDIM/2,std::size_t(1));
1412 auto func=dynamic_cast<CompositeFunctorInterface<T,NDIM,LDIM>* >(g.get_impl()->get_functor().get());
1414 func->make_redundant(true);
1415 func->replicate_low_dim_functions(true);
1416 this->reconstruct(); // if this == &g we don't need g to be redundant
1417
1419
1420 TENSOR_RESULT_TYPE(T, R) local = impl->inner_local_on_demand(*g.get_impl());
1421 impl->world.gop.sum(local);
1422 impl->world.gop.fence();
1423
1424 return local;
1425 }
1426
1427 /// project this on the low-dim function g: h(x) = <f(x,y) | g(y)>
1428
1429 /// @param[in] g low-dim function
1430 /// @param[in] dim over which dimensions to be integrated: 0..LDIM-1 or LDIM..NDIM-1
1431 /// @return new function of dimension NDIM-LDIM
1432 template <typename R, size_t LDIM>
1434 if (NDIM<=LDIM) MADNESS_EXCEPTION("confused dimensions in project_out?",1);
1435 MADNESS_CHECK_THROW(dim==0 or dim==1,"dim must be 0 or 1 in project_out");
1436 verify();
1437 typedef TENSOR_RESULT_TYPE(T,R) resultT;
1438 static const size_t KDIM=NDIM-LDIM;
1439
1441 .k(g.k()).thresh(g.thresh());
1442 Function<resultT,KDIM> result=factory; // no empty() here!
1443
1445 g.change_tree_state(redundant,false);
1446 world().gop.fence();
1447 this->get_impl()->project_out(result.get_impl().get(),g.get_impl().get(),dim,true);
1448// result.get_impl()->project_out2(this->get_impl().get(),gimpl,dim);
1449 result.world().gop.fence();
1450 g.change_tree_state(reconstructed,false);
1451 result.get_impl()->trickle_down(false);
1452 result.get_impl()->set_tree_state(reconstructed);
1453 result.world().gop.fence();
1454 return result;
1455 }
1456
1457 Function<T,NDIM/2> dirac_convolution(const bool fence=true) const {
1458 constexpr std::size_t LDIM=NDIM/2;
1459 MADNESS_CHECK_THROW(NDIM==2*LDIM,"NDIM must be even");
1460// // this will be the result function
1462 Function<T,LDIM> f = factory;
1463 if(!is_reconstructed()) this->reconstruct();
1464 this->get_impl()->do_dirac_convolution(f.get_impl().get(),fence);
1465 return f;
1466 }
1467
1468 /// Replaces this function with one loaded from an archive using the default processor map
1469
1470 /// Archive can be sequential or parallel.
1471 ///
1472 /// The & operator for serializing will only work with parallel archives.
1473 template <typename Archive>
1474 void load(World& world, Archive& ar) {
1476 // Type checking since we are probably circumventing the archive's own type checking
1477 long magic = 0l, id = 0l, ndim = 0l, k = 0l;
1478 ar & magic & id & ndim & k;
1479 MADNESS_ASSERT(magic == 7776768); // Mellow Mushroom Pizza tel.# in Knoxville
1481 MADNESS_ASSERT(ndim == NDIM);
1482
1483 impl.reset(new implT(FunctionFactory<T,NDIM>(world).k(k).empty()));
1484
1485 impl->load(ar);
1486 }
1487
1488
1489 /// Stores the function to an archive
1490
1491 /// Archive can be sequential or parallel.
1492 ///
1493 /// The & operator for serializing will only work with parallel archives.
1494 template <typename Archive>
1495 void store(Archive& ar) const {
1497 verify();
1498 // For type checking, etc.
1499 ar & long(7776768) & long(TensorTypeData<T>::id) & long(NDIM) & long(k());
1500
1501 impl->store(ar);
1502 }
1503
1504 /// change the tensor type of the coefficients in the FunctionNode
1505
1506 /// @param[in] targs target tensor arguments (threshold and full/low rank)
1507 void change_tensor_type(const TensorArgs& targs, bool fence=true) {
1508 if (not impl) return;
1509 impl->change_tensor_type1(targs,fence);
1510 }
1511
1512
1513 /// This is replaced with left*right ... private
1514 template <typename Q, typename opT>
1516 const opT& op, bool fence) {
1518 func.verify();
1519 MADNESS_ASSERT(func.is_reconstructed());
1520 if (VERIFY_TREE) func.verify_tree();
1521 impl.reset(new implT(*func.get_impl(), func.get_pmap(), false));
1522 impl->unaryXX(func.get_impl().get(), op, fence);
1523 return *this;
1524 }
1525
1526 /// Returns vector of FunctionImpl pointers corresponding to vector of functions
1527 template <typename Q, std::size_t D>
1528 static std::vector< std::shared_ptr< FunctionImpl<Q,D> > > vimpl(const std::vector< Function<Q,D> >& v) {
1530 std::vector< std::shared_ptr< FunctionImpl<Q,D> > > r(v.size());
1531 for (unsigned int i=0; i<v.size(); ++i) r[i] = v[i].get_impl();
1532 return r;
1533 }
1534
1535 /// This is replaced with op(vector of functions) ... private
1536 template <typename opT>
1537 Function<T,NDIM>& multiop_values(const opT& op, const std::vector< Function<T,NDIM> >& vf) {
1538 std::vector<implT*> v(vf.size(),NULL);
1539 for (unsigned int i=0; i<v.size(); ++i) {
1540 if (vf[i].is_initialized()) v[i] = vf[i].get_impl().get();
1541 }
1542 impl->multiop_values(op, v);
1543 world().gop.fence();
1544 if (VERIFY_TREE) verify_tree();
1545
1546 return *this;
1547 }
1548
1549 /// apply op on the input vector yielding an output vector of functions
1550
1551 /// (*this) is just a dummy Function to be able to call internal methods in FuncImpl
1552 /// @param[in] op the operator working on vin
1553 /// @param[in] vin vector of input Functions
1554 /// @param[out] vout vector of output Functions vout = op(vin)
1555 template <typename opT>
1557 const std::vector< Function<T,NDIM> >& vin,
1558 std::vector< Function<T,NDIM> >& vout,
1559 const bool fence=true) {
1560 std::vector<implT*> vimplin(vin.size(),NULL);
1561 for (unsigned int i=0; i<vin.size(); ++i) {
1562 if (vin[i].is_initialized()) vimplin[i] = vin[i].get_impl().get();
1563 }
1564 std::vector<implT*> vimplout(vout.size(),NULL);
1565 for (unsigned int i=0; i<vout.size(); ++i) {
1566 if (vout[i].is_initialized()) vimplout[i] = vout[i].get_impl().get();
1567 }
1568
1569 impl->multi_to_multi_op_values(op, vimplin, vimplout, fence);
1570 if (VERIFY_TREE) verify_tree();
1571
1572 }
1573
1574
1575 /// Multiplication of function * vector of functions using recursive algorithm of mulxx
1576 template <typename L, typename R>
1577 void vmulXX(const Function<L,NDIM>& left,
1578 const std::vector< Function<R,NDIM> >& right,
1579 std::vector< Function<T,NDIM> >& result,
1580 double tol,
1581 bool fence) {
1583
1584 std::vector<FunctionImpl<T,NDIM>*> vresult(right.size());
1585 std::vector<const FunctionImpl<R,NDIM>*> vright(right.size());
1586 for (unsigned int i=0; i<right.size(); ++i) {
1587 result[i].set_impl(left,false);
1588 vresult[i] = result[i].impl.get();
1589 vright[i] = right[i].get_impl().get();
1590 }
1591
1592 left.world().gop.fence(); // Is this still essential? Yes.
1593 vresult[0]->mulXXvec(left.get_impl().get(), vright, vresult, tol, fence);
1594 }
1595
1596 /// Same as \c operator* but with optional fence and no automatic reconstruction
1597
1598 /// f or g are on-demand functions
1599 template<typename L, typename R>
1600 void mul_on_demand(const Function<L,NDIM>& f, const Function<R,NDIM>& g, bool fence=true) {
1601 const FunctionImpl<L,NDIM>* fimpl=f.get_impl().get();
1602 const FunctionImpl<R,NDIM>* gimpl=g.get_impl().get();
1603 if (fimpl->is_on_demand() and gimpl->is_on_demand()) {
1604 MADNESS_EXCEPTION("can't multiply two on-demand functions",1);
1605 }
1606
1607 if (fimpl->is_on_demand()) {
1608 leaf_op<R,NDIM> leaf_op1(gimpl);
1609 impl->multiply(leaf_op1,gimpl,fimpl,fence);
1610 } else {
1611 leaf_op<L,NDIM> leaf_op1(fimpl);
1612 impl->multiply(leaf_op1,fimpl,gimpl,fence);
1613 }
1614 }
1615
1616 /// sparse transformation of a vector of functions ... private
1617 template <typename R, typename Q>
1618 void vtransform(const std::vector< Function<R,NDIM> >& v,
1619 const Tensor<Q>& c,
1620 std::vector< Function<T,NDIM> >& vresult,
1621 double tol,
1622 bool fence=true) {
1624 vresult[0].impl->vtransform(vimpl(v), c, vimpl(vresult), tol, fence);
1625 }
1626
1627 /// This is replaced with alpha*left + beta*right ... private
1628 template <typename L, typename R>
1630 T beta, const Function<R,NDIM>& right, bool fence) {
1632 left.verify();
1633 right.verify();
1634 MADNESS_ASSERT(left.is_compressed() && right.is_compressed());
1635 if (VERIFY_TREE) left.verify_tree();
1636 if (VERIFY_TREE) right.verify_tree();
1637 impl.reset(new implT(*left.get_impl(), left.get_pmap(), false));
1638 impl->gaxpy(alpha,*left.get_impl(),beta,*right.get_impl(),fence);
1639 return *this;
1640 }
1641
1642 /// This is replaced with mapdim(f) ... private
1643 Function<T,NDIM>& mapdim(const Function<T,NDIM>& f, const std::vector<long>& map, bool fence) {
1645 f.verify();
1646 if (VERIFY_TREE) f.verify_tree();
1647 for (std::size_t i=0; i<NDIM; ++i) MADNESS_ASSERT(map[i]>=0 && static_cast<std::size_t>(map[i])<NDIM);
1648 impl.reset(new implT(*f.impl, f.get_pmap(), false));
1649 impl->mapdim(*f.impl,map,fence);
1650 return *this;
1651 }
1652
1653 /// This is replaced with mirror(f) ... private
1654
1655 /// similar to mapdim, but maps from x to -x, y to -y, and so on
1656 /// Example: mirror a 3d function on the xy plane: mirror={1,1,-1}
1657 /// @param[in] mirror array of -1 and 1, corresponding to mirror or not
1658 Function<T,NDIM>& mirror(const Function<T,NDIM>& f, const std::vector<long>& mirrormap, bool fence) {
1660 f.verify();
1661 if (VERIFY_TREE) f.verify_tree();
1662 for (std::size_t i=0; i<NDIM; ++i) MADNESS_ASSERT((mirrormap[i]==1) or (mirrormap[i]==-1));
1663 impl.reset(new implT(*f.impl, f.get_pmap(), false));
1664 impl->mirror(*f.impl,mirrormap,fence);
1665 return *this;
1666 }
1667
1668 /// This is replaced with mirror(map(f)) ... private
1669
1670 /// first map then mirror!
1671 /// mirror is similar to mapdim, but maps from x to -x, y to -y, and so on
1672 /// Example: mirror a 3d function on the xy plane: mirror={1,1,-1}
1673 /// Example: c4 rotation of a 3d function around the z axis:
1674 /// x->y, y->-x, z->z: map(1,0,2); mirror(-1,1,1)
1675 /// @param[in] map array holding dimensions
1676 /// @param[in] mirror array of -1 and 1, corresponding to mirror or not
1678 const std::vector<long>& map, const std::vector<long>& mirror,
1679 bool fence) {
1681 f.verify();
1682 if (VERIFY_TREE) f.verify_tree();
1683 for (std::size_t i=0; i<mirror.size(); ++i) MADNESS_ASSERT((mirror[i]==1) or (mirror[i]==-1));
1684 for (std::size_t i=0; i<map.size(); ++i) MADNESS_ASSERT(map[i]>=0 && static_cast<std::size_t>(map[i])<NDIM);
1685
1686 impl.reset(new implT(*f.impl, f.get_pmap(), false));
1687 impl->map_and_mirror(*f.impl,map,mirror,fence);
1688 return *this;
1689 }
1690
1691
1692 /// check symmetry of a function by computing the 2nd derivative
1693 double check_symmetry() const {
1694
1696 if (VERIFY_TREE) verify_tree();
1697 double local = impl->check_symmetry_local();
1698 impl->world.gop.sum(local);
1699 impl->world.gop.fence();
1700 double asy=sqrt(local);
1701 if (this->world().rank()==0) print("asymmetry wrt particle",asy);
1703 return asy;
1704 }
1705
1706 /// reduce the rank of the coefficient tensors
1707 Function<T,NDIM>& reduce_rank(const double thresh=0.0, const bool fence=true) {
1708 verify();
1709 double thresh1= (thresh==0.0) ? impl->get_tensor_args().thresh : thresh;
1710 impl->reduce_rank(thresh1,fence);
1711 return *this;
1712 }
1713
1714 /// remove all nodes with level higher than n
1715 Function<T,NDIM>& chop_at_level(const int n, const bool fence=true) {
1716 verify();
1718 impl->chop_at_level(n,true);
1720 return *this;
1721 }
1722 };
1723
1724// template <typename T, typename opT, std::size_t NDIM>
1725 template <typename T, typename opT, std::size_t NDIM>
1726 Function<T,NDIM> multiop_values(const opT& op, const std::vector< Function<T,NDIM> >& vf) {
1728 r.set_impl(vf[0], false);
1729 r.multiop_values(op, vf);
1730 return r;
1731 }
1732
1733 /// Returns new function equal to alpha*f(x) with optional fence
1734 template <typename Q, typename T, std::size_t NDIM>
1735 Function<TENSOR_RESULT_TYPE(Q,T),NDIM>
1736 mul(const Q alpha, const Function<T,NDIM>& f, bool fence=true) {
1738 f.verify();
1739 if (VERIFY_TREE) f.verify_tree();
1741 result.set_impl(f, false);
1742 result.get_impl()->scale_oop(alpha,*f.get_impl(),fence);
1743 return result;
1744 }
1745
1746
1747 /// Returns new function equal to f(x)*alpha with optional fence
1748 template <typename Q, typename T, std::size_t NDIM>
1749 Function<TENSOR_RESULT_TYPE(Q,T),NDIM>
1750 mul(const Function<T,NDIM>& f, const Q alpha, bool fence=true) {
1752 return mul(alpha,f,fence);
1753 }
1754
1755
1756 /// Returns new function equal to f(x)*alpha
1757
1758 /// Using operator notation forces a global fence after each operation
1759 template <typename Q, typename T, std::size_t NDIM>
1760 Function<TENSOR_RESULT_TYPE(Q,T),NDIM>
1762 return mul(alpha, f, true);
1763 }
1764
1765 /// Returns new function equal to alpha*f(x)
1766
1767 /// Using operator notation forces a global fence after each operation
1768 template <typename Q, typename T, std::size_t NDIM>
1769 Function<TENSOR_RESULT_TYPE(Q,T),NDIM>
1771 return mul(alpha, f, true);
1772 }
1773
1774 /// Sparse multiplication --- left and right must be reconstructed and if tol!=0 have tree of norms already created
1775 template <typename L, typename R,std::size_t NDIM>
1776 Function<TENSOR_RESULT_TYPE(L,R),NDIM>
1777 mul_sparse(const Function<L,NDIM>& left, const Function<R,NDIM>& right, double tol, bool fence=true) {
1779 left.verify();
1780 right.verify();
1782 if (VERIFY_TREE) left.verify_tree();
1783 if (VERIFY_TREE) right.verify_tree();
1784
1786 result.set_impl(left, false);
1787 result.get_impl()->mulXX(left.get_impl().get(), right.get_impl().get(), tol, fence);
1788 return result;
1789 }
1790
1791 /// Same as \c operator* but with optional fence and no automatic reconstruction
1792 template <typename L, typename R,std::size_t NDIM>
1793 Function<TENSOR_RESULT_TYPE(L,R),NDIM>
1794 mul(const Function<L,NDIM>& left, const Function<R,NDIM>& right, bool fence=true) {
1795 return mul_sparse(left,right,0.0,fence);
1796 }
1797
1798 /// Generate new function = op(left,right) where op acts on the function values
1799 template <typename L, typename R, typename opT, std::size_t NDIM>
1800 Function<TENSOR_RESULT_TYPE(L,R),NDIM>
1801 binary_op(const Function<L,NDIM>& left, const Function<R,NDIM>& right, const opT& op, bool fence=true) {
1803 if (!left.is_reconstructed()) left.reconstruct();
1804 if (!right.is_reconstructed()) right.reconstruct();
1805
1807 result.set_impl(left, false);
1808 result.get_impl()->binaryXX(left.get_impl().get(), right.get_impl().get(), op, fence);
1809 return result;
1810 }
1811
1812 /// Out of place application of unary operation to function values with optional fence
1813 template <typename Q, typename opT, std::size_t NDIM>
1814 Function<typename opT::resultT, NDIM>
1815 unary_op(const Function<Q,NDIM>& func, const opT& op, bool fence=true) {
1816 if (!func.is_reconstructed()) func.reconstruct();
1819 result.set_impl(func, false);
1820 result.get_impl()->unaryXXvalues(func.get_impl().get(), op, fence);
1821 return result;
1822 }
1823
1824
1825 /// Out of place application of unary operation to scaling function coefficients with optional fence
1826 template <typename Q, typename opT, std::size_t NDIM>
1827 Function<typename opT::resultT, NDIM>
1828 unary_op_coeffs(const Function<Q,NDIM>& func, const opT& op, bool fence=true) {
1829 if (!func.is_reconstructed()) func.reconstruct();
1831 return result.unary_op_coeffs(func,op,fence);
1832 }
1833
1834 /// Use the vmra/mul(...) interface instead
1835
1836 /// This so that we don't have to have friend functions in a different header.
1837 ///
1838 /// If using sparsity (tol != 0) you must have created the tree of norms
1839 /// already for both left and right.
1840 template <typename L, typename R, std::size_t D>
1841 std::vector< Function<TENSOR_RESULT_TYPE(L,R),D> >
1842 vmulXX(const Function<L,D>& left, const std::vector< Function<R,D> >& vright, double tol, bool fence=true) {
1843 if (vright.size() == 0) return std::vector< Function<TENSOR_RESULT_TYPE(L,R),D> >();
1844 std::vector< Function<TENSOR_RESULT_TYPE(L,R),D> > vresult(vright.size());
1845 vresult[0].vmulXX(left, vright, vresult, tol, fence);
1846 return vresult;
1847 }
1848
1849 /// Multiplies two functions with the new result being of type TensorResultType<L,R>
1850
1851 /// Using operator notation forces a global fence after each operation but also
1852 /// enables us to automatically reconstruct the input functions as required.
1853 template <typename L, typename R, std::size_t NDIM>
1854 Function<TENSOR_RESULT_TYPE(L,R), NDIM>
1855 operator*(const Function<L,NDIM>& left, const Function<R,NDIM>& right) {
1856 if (!left.is_reconstructed()) left.reconstruct();
1857 if (!right.is_reconstructed()) right.reconstruct();
1858 MADNESS_ASSERT(not (left.is_on_demand() or right.is_on_demand()));
1859 return mul(left,right,true);
1860 }
1861
1862 /// Performs a Hartree/outer product on the two given low-dimensional function vectors
1863
1864 /// @return result(x,y) = \sum_i f_i(x) g_i(y)
1865 template<typename T, std::size_t KDIM, std::size_t LDIM>
1866 Function<T,KDIM+LDIM>
1867 hartree_product(const std::vector<Function<T,KDIM>>& left, const std::vector<Function<T,LDIM>>& right) {
1868
1869 MADNESS_CHECK_THROW(left.size()==right.size(), "hartree_product: left and right must have same size");
1870 if (left.size()==0) return Function<T,KDIM+LDIM>();
1871
1873
1875 .k(left.front().k()).thresh(thresh);
1876 Function<T,KDIM+LDIM> result=factory.empty();
1877
1878 // some prep work
1881 std::vector<std::shared_ptr<FunctionImpl<T,KDIM>>> vleft=get_impl(left);
1882 std::vector<std::shared_ptr<FunctionImpl<T,LDIM>>> vright=get_impl(right);
1883
1884 result.do_hartree_product(vleft,vright);
1885
1886 return result;
1887
1888 }
1889
1890 /// Performs a Hartree product on the two given low-dimensional functions
1891 template<typename T, std::size_t KDIM, std::size_t LDIM>
1892 Function<T,KDIM+LDIM>
1894 typedef std::vector<Function<T,KDIM>> vector;
1895 return hartree_product(vector({left2}),vector({right2}));
1896 }
1897
1898 /// Performs a Hartree product on the two given low-dimensional functions
1899 template<typename T, std::size_t KDIM, std::size_t LDIM, typename opT>
1900 Function<T,KDIM+LDIM>
1902 const opT& op) {
1903
1904 // we need both sum and difference coeffs for error estimation
1905 Function<T,KDIM>& left = const_cast< Function<T,KDIM>& >(left2);
1906 Function<T,LDIM>& right = const_cast< Function<T,LDIM>& >(right2);
1907
1909
1911 .k(left.k()).thresh(thresh);
1912 Function<T,KDIM+LDIM> result=factory.empty();
1913
1914 if (result.world().rank()==0) {
1915 print("incomplete FunctionFactory in Function::hartree_product");
1916 print("thresh: ", thresh);
1917 }
1918 bool same=(left2.get_impl()==right2.get_impl());
1919
1920 // some prep work
1921 left.make_nonstandard(true, true);
1922 right.make_nonstandard(true, true);
1923
1924 std::vector<std::shared_ptr<FunctionImpl<T,KDIM>>> vleft;
1925 std::vector<std::shared_ptr<FunctionImpl<T,LDIM>>> vright;
1926 vleft.push_back(left.get_impl());
1927 vright.push_back(right.get_impl());
1928 result.do_hartree_product(vleft,right,&op);
1929
1930 left.standard(false);
1931 if (not same) right.standard(false);
1932 left2.world().gop.fence();
1933
1934 return result;
1935 }
1936
1937 /// adds beta*right only left: alpha*left + beta*right optional fence and no automatic compression
1938
1939 /// left and right might live in different worlds, the accumulation is non-blocking
1940 template <typename L, typename R,std::size_t NDIM>
1941 void
1943 TENSOR_RESULT_TYPE(L,R) beta, const Function<R,NDIM>& right, bool fence=true) {
1946 left.gaxpy(alpha, right, beta, fence);
1947 }
1948
1949 /// Returns new function alpha*left + beta*right optional fence and no automatic compression
1950 template <typename L, typename R,std::size_t NDIM>
1951 Function<TENSOR_RESULT_TYPE(L,R),NDIM>
1953 TENSOR_RESULT_TYPE(L,R) beta, const Function<R,NDIM>& right, bool fence=true) {
1956 return result.gaxpy_oop(alpha, left, beta, right, fence);
1957 }
1958
1959 /// Same as \c operator+ but with optional fence and no automatic compression
1960 template <typename L, typename R,std::size_t NDIM>
1961 Function<TENSOR_RESULT_TYPE(L,R),NDIM>
1962 add(const Function<L,NDIM>& left, const Function<R,NDIM>& right, bool fence=true) {
1963 return gaxpy_oop(TENSOR_RESULT_TYPE(L,R)(1.0), left,
1964 TENSOR_RESULT_TYPE(L,R)(1.0), right, fence);
1965 }
1966
1967
1968 /// Returns new function alpha*left + beta*right optional fence, having both addends reconstructed
1969 template<typename T, std::size_t NDIM>
1971 const double beta, const Function<T,NDIM>& right, const bool fence=true) {
1972 Function<T,NDIM> result;
1973 result.set_impl(right,false);
1974
1977 result.get_impl()->gaxpy_oop_reconstructed(alpha,*left.get_impl(),beta,*right.get_impl(),fence);
1978 return result;
1979
1980 }
1981
1982 /// Adds two functions with the new result being of type TensorResultType<L,R>
1983
1984 /// Using operator notation forces a global fence after each operation
1985 template <typename L, typename R, std::size_t NDIM>
1986 Function<TENSOR_RESULT_TYPE(L,R), NDIM>
1987 operator+(const Function<L,NDIM>& left, const Function<R,NDIM>& right) {
1988 if (VERIFY_TREE) left.verify_tree();
1989 if (VERIFY_TREE) right.verify_tree();
1990
1991 // no compression for high-dimensional functions
1992 if (NDIM==6) {
1993 left.reconstruct();
1994 right.reconstruct();
1995 return gaxpy_oop_reconstructed(1.0,left,1.0,right,true);
1996 } else {
1997 if (!left.is_compressed()) left.compress();
1998 if (!right.is_compressed()) right.compress();
1999 return add(left,right,true);
2000 }
2001 }
2002
2003 /// Same as \c operator- but with optional fence and no automatic compression
2004 template <typename L, typename R,std::size_t NDIM>
2005 Function<TENSOR_RESULT_TYPE(L,R),NDIM>
2006 sub(const Function<L,NDIM>& left, const Function<R,NDIM>& right, bool fence=true) {
2007 return gaxpy_oop(TENSOR_RESULT_TYPE(L,R)(1.0), left,
2008 TENSOR_RESULT_TYPE(L,R)(-1.0), right, fence);
2009 }
2010
2011
2012 /// Subtracts two functions with the new result being of type TensorResultType<L,R>
2013
2014 /// Using operator notation forces a global fence after each operation
2015 template <typename L, typename R, std::size_t NDIM>
2016 Function<TENSOR_RESULT_TYPE(L,R), NDIM>
2017 operator-(const Function<L,NDIM>& left, const Function<R,NDIM>& right) {
2019 // no compression for high-dimensional functions
2020 if (NDIM==6) {
2021 left.reconstruct();
2022 right.reconstruct();
2023 return gaxpy_oop_reconstructed(1.0,left,-1.0,right,true);
2024 } else {
2025 if (!left.is_compressed()) left.compress();
2026 if (!right.is_compressed()) right.compress();
2027 return sub(left,right,true);
2028 }
2029 }
2030
2031 /// Create a new copy of the function with different distribution and optional fence
2032
2033 /// Works in either basis. Different distributions imply
2034 /// asynchronous communication and the optional fence is
2035 /// collective.
2036 template <typename T, std::size_t NDIM>
2038 const std::shared_ptr< WorldDCPmapInterface< Key<NDIM> > >& pmap,
2039 bool fence = true) {
2041 f.verify();
2042 Function<T,NDIM> result;
2043 typedef FunctionImpl<T,NDIM> implT;
2044 result.set_impl(std::shared_ptr<implT>(new implT(*f.get_impl(), pmap, false)));
2045 result.get_impl()->copy_coeffs(*f.get_impl(), fence);
2046 if (VERIFY_TREE) result.verify_tree();
2047 return result;
2048 }
2049
2050 /// Create a new copy of the function with the same distribution and optional fence
2051 template <typename T, std::size_t NDIM>
2052 Function<T,NDIM> copy(const Function<T,NDIM>& f, bool fence = true) {
2054 return copy(f, f.get_pmap(), fence);
2055 }
2056
2057 /// Type conversion implies a deep copy. No communication except for optional fence.
2058
2059 /// Works in either basis but any loss of precision may result in different errors
2060 /// in applied in a different basis.
2061 ///
2062 /// The new function is formed with the options from the default constructor.
2063 ///
2064 /// There is no automatic type conversion since this is generally a rather dangerous
2065 /// thing and because there would be no way to make the fence optional.
2066 template <typename T, typename Q, std::size_t NDIM>
2067 Function<Q,NDIM> convert(const Function<T,NDIM>& f, bool fence = true) {
2069 f.verify();
2070 Function<Q,NDIM> result;
2071 result.set_impl(f, false);
2072 result.get_impl()->copy_coeffs(*f.get_impl(), fence);
2073 return result;
2074 }
2075
2076
2077 /// Return the complex conjugate of the input function with the same distribution and optional fence
2078
2079 /// !!! The fence is actually not optional in the current implementation !!!
2080 template <typename T, std::size_t NDIM>
2081 Function<T,NDIM> conj(const Function<T,NDIM>& f, bool fence = true) {
2083 Function<T,NDIM> result = copy(f,true);
2084 return result.conj(fence);
2085 }
2086
2087 /// Apply operator on a hartree product of two low-dimensional functions
2088
2089 /// Supposed to be something like result= G( f(1)*f(2))
2090 /// the hartree product is never constructed explicitly, but its coeffs are
2091 /// constructed on the fly and processed immediately.
2092 /// @param[in] op the operator
2093 /// @param[in] f1 function of particle 1
2094 /// @param[in] f2 function of particle 2
2095 /// @param[in] fence if we shall fence
2096 /// @return a function of dimension NDIM=LDIM+LDIM
2097 template <typename opT, typename T, std::size_t LDIM>
2098 Function<TENSOR_RESULT_TYPE(typename opT::opT,T), LDIM+LDIM>
2099 apply(const opT& op, const std::vector<Function<T,LDIM>>& f1, const std::vector<Function<T,LDIM>>& f2, bool fence=true) {
2100
2101 World& world=f1.front().world();
2102
2103 typedef TENSOR_RESULT_TYPE(T,typename opT::opT) resultT;
2104 typedef std::vector<Function<T,LDIM>> vecfuncL;
2105
2106 vecfuncL& ff1 = const_cast< vecfuncL& >(f1);
2107 vecfuncL& ff2 = const_cast< vecfuncL& >(f2);
2108
2109 bool same=(ff1[0].get_impl()==ff2[0].get_impl());
2110
2111 reconstruct(world,f1,false);
2112 reconstruct(world,f2,false);
2113 world.gop.fence();
2114 // keep the leaves! They are assumed to be there later
2115 // even for modified op we need NS form for the hartree_leaf_op
2116 for (auto& f : f1) f.make_nonstandard(true,false);
2117 for (auto& f : f2) f.make_nonstandard(true,false);
2118 world.gop.fence();
2119
2120
2123 Function<resultT,LDIM+LDIM> result=factory.empty().fence();
2124
2125 result.get_impl()->reset_timer();
2126 op.reset_timer();
2127
2128 // will fence here
2129 for (size_t i=0; i<f1.size(); ++i)
2130 result.get_impl()->recursive_apply(op, f1[i].get_impl().get(),f2[i].get_impl().get(),false);
2131 world.gop.fence();
2132
2133 if (op.print_timings) {
2134 result.get_impl()->print_timer();
2135 op.print_timer();
2136 }
2137
2138 result.get_impl()->finalize_apply(); // need fence before reconstruct
2139
2140 if (op.modified()) {
2141 result.get_impl()->trickle_down(true);
2142 } else {
2143 result.get_impl()->reconstruct(true);
2144 }
2145 standard(world,ff1,false);
2146 if (not same) standard(world,ff2,false);
2147
2148 return result;
2149 }
2150
2151
2152 /// Apply operator ONLY in non-standard form - required other steps missing !!
2153 template <typename opT, typename R, std::size_t NDIM>
2154 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM>
2155 apply_only(const opT& op, const Function<R,NDIM>& f, bool fence=true) {
2156 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> result;
2157
2158 constexpr std::size_t OPDIM=opT::opdim;
2159 constexpr bool low_dim=(OPDIM*2==NDIM); // apply on some dimensions only
2160
2161 // specialized version for 3D
2162 if (NDIM <= 3 and (not low_dim)) {
2163 result.set_impl(f, false);
2164 result.get_impl()->apply(op, *f.get_impl(), fence);
2165
2166 } else { // general version for higher dimension
2167 //bool print_timings=false;
2168 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> r1;
2169
2170 result.set_impl(f, false);
2171 r1.set_impl(f, false);
2172
2173 result.get_impl()->reset_timer();
2174 op.reset_timer();
2175
2176 result.get_impl()->apply_source_driven(op, *f.get_impl(), fence);
2177
2178 // recursive_apply is about 20% faster than apply_source_driven
2179 //result.get_impl()->recursive_apply(op, f.get_impl().get(),
2180 // r1.get_impl().get(),true); // will fence here
2181
2182 }
2183
2184 return result;
2185 }
2186
2187 /// Apply operator in non-standard form
2188
2189 /// Returns a new function with the same distribution
2190 ///
2191 /// !!! For the moment does NOT respect fence option ... always fences
2192 /// if the operator acts on one particle only the result will be sorted as
2193 /// g.particle=1: g(f) = \int g(x,x') f(x',y) dx' = result(x,y)
2194 /// g.particle=2: g(f) = \int g(y,y') f(x,y') dy' = result(x,y)
2195 /// for the second case it will notably *not* be as it is implemented in the partial inner product!
2196 /// g.particle=2 g(f) = result(x,y)
2197 /// inner(g(y,y'),f(x,y'),1,1) = result(y,x)
2198 /// also note the confusion with the counting of the particles/integration variables
2199 template <typename opT, typename R, std::size_t NDIM>
2200 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM>
2201 apply(const opT& op, const Function<R,NDIM>& f, bool fence=true) {
2202
2203 typedef TENSOR_RESULT_TYPE(typename opT::opT,R) resultT;
2204 Function<R,NDIM>& ff = const_cast< Function<R,NDIM>& >(f);
2206
2207 MADNESS_ASSERT(not f.is_on_demand());
2208 bool print_timings=op.print_timings;
2209
2210 if (VERIFY_TREE) ff.verify_tree();
2211 ff.reconstruct();
2212 if (print_timings) ff.print_size("ff in apply after reconstruct");
2213
2214 if (op.modified()) {
2215
2217// ff.get_impl()->make_redundant(true);
2218 result = apply_only(op, ff, fence);
2219 ff.get_impl()->undo_redundant(false);
2220 result.get_impl()->trickle_down(true);
2221
2222 } else {
2223
2224 // saves the standard() step, which is very expensive in 6D
2225// Function<R,NDIM> fff=copy(ff);
2226 Function<R,NDIM> fff=(ff);
2227 fff.make_nonstandard(op.doleaves, true);
2228 if (print_timings) fff.print_size("ff in apply after make_nonstandard");
2229 if ((print_timings) and (f.world().rank()==0)) {
2230 fff.get_impl()->timer_filter.print("filter");
2231 fff.get_impl()->timer_compress_svd.print("compress_svd");
2232 }
2233 result = apply_only(op, fff, fence);
2234 result.get_impl()->set_tree_state(nonstandard_after_apply);
2235 ff.world().gop.fence();
2236 if (print_timings) result.print_size("result after apply_only");
2237
2238 // svd-tensors need some post-processing
2239 if (result.get_impl()->get_tensor_type()==TT_2D) {
2240 double elapsed=result.get_impl()->finalize_apply();
2241 if (print_timings) printf("time in finalize_apply %8.2f\n",elapsed);
2242 }
2243 if (print_timings) {
2244 result.get_impl()->print_timer();
2245 op.print_timer();
2246 }
2247
2248 result.get_impl()->reconstruct(true);
2249
2250// fff.clear();
2251 if (op.destructive()) {
2252 ff.world().gop.fence();
2253 ff.clear();
2254 } else {
2255 ff.standard();
2256 }
2257
2258 }
2259 if (print_timings) result.print_size("result after reconstruction");
2260 return result;
2261 }
2262
2263
2264 template <typename opT, typename R, std::size_t NDIM>
2265 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM>
2266 apply_1d_realspace_push(const opT& op, const Function<R,NDIM>& f, int axis, bool fence=true) {
2268 Function<R,NDIM>& ff = const_cast< Function<R,NDIM>& >(f);
2269 if (VERIFY_TREE) ff.verify_tree();
2270 ff.reconstruct();
2271
2272 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> result;
2273
2274 result.set_impl(ff, false);
2275 result.get_impl()->apply_1d_realspace_push(op, ff.get_impl().get(), axis, fence);
2276 result.get_impl()->set_tree_state(redundant_after_merge);
2277 return result;
2278 }
2279
2280
2281 /// Generate a new function by reordering dimensions ... optional fence
2282
2283 /// You provide an array of dimension NDIM that maps old to new dimensions
2284 /// according to
2285 /// \code
2286 /// newdim = mapdim[olddim]
2287 /// \endcode
2288 /// Works in either scaling function or wavelet basis.
2289 ///
2290 /// Would be easy to modify this to also change the procmap here
2291 /// if desired but presently it uses the same procmap as f.
2292 template <typename T, std::size_t NDIM>
2293 Function<T,NDIM>
2294 mapdim(const Function<T,NDIM>& f, const std::vector<long>& map, bool fence=true) {
2296 Function<T,NDIM> result;
2297 return result.mapdim(f,map,fence);
2298 }
2299
2300 /// Generate a new function by mirroring within the dimensions .. optional fence
2301
2302 /// similar to mapdim
2303 /// @param[in] mirror array with -1 and 1, corresponding to mirror this dimension or not
2304 template <typename T, std::size_t NDIM>
2305 Function<T,NDIM>
2306 mirror(const Function<T,NDIM>& f, const std::vector<long>& mirrormap, bool fence=true) {
2308 Function<T,NDIM> result;
2309 return result.mirror(f,mirrormap,fence);
2310 }
2311
2312 /// This is replaced with mirror(map(f)), optional fence
2313
2314 /// first map then mirror!
2315 /// mirror is similar to mapdim, but maps from x to -x, y to -y, and so on
2316 /// Example: mirror a 3d function on the xy plane: mirror={1,1,-1}
2317 /// Example: c4 rotation of a 3d function around the z axis:
2318 /// x->y, y->-x, z->z: map(1,0,2); mirror(-1,1,1)
2319 /// @param[in] map array holding dimensions
2320 /// @param[in] mirror array of -1 and 1, corresponding to mirror or not
2321 template <typename T, std::size_t NDIM>
2322 Function<T,NDIM>
2323 map_and_mirror(const Function<T,NDIM>& f, const std::vector<long>& map,
2324 const std::vector<long>& mirror, bool fence=true) {
2326 Function<T,NDIM> result;
2327 return result.map_and_mirror(f,map,mirror,fence);
2328 }
2329
2330
2331 /// swap particles 1 and 2
2332
2333 /// param[in] f a function of 2 particles f(1,2)
2334 /// return the input function with particles swapped g(1,2) = f(2,1)
2335 template <typename T, std::size_t NDIM>
2336 typename std::enable_if_t<NDIM%2==0, Function<T,NDIM>>
2338 // this could be done more efficiently for SVD, but it works decently
2339 std::vector<long> map(NDIM);
2340 constexpr std::size_t LDIM=NDIM/2;
2341 static_assert(LDIM*2==NDIM);
2342 for (auto d=0; d<LDIM; ++d) {
2343 map[d]=d+LDIM;
2344 map[d+LDIM]=d;
2345 }
2346// map[0]=3;
2347// map[1]=4;
2348// map[2]=5; // 2 -> 1
2349// map[3]=0;
2350// map[4]=1;
2351// map[5]=2; // 1 -> 2
2352 return mapdim(f,map);
2353 }
2354
2355 /// symmetrize a function
2356
2357 /// @param[in] symmetry possibilities are:
2358 /// (anti-) symmetric particle permutation ("sy_particle", "antisy_particle")
2359 /// symmetric mirror plane ("xy", "xz", "yz")
2360 /// @return a new function symmetrized according to the input parameter
2361 template <typename T, std::size_t NDIM>
2362 Function<T,NDIM>
2363 symmetrize(const Function<T,NDIM>& f, const std::string symmetry, bool fence=true) {
2364 Function<T,NDIM> result;
2365
2366 MADNESS_ASSERT(NDIM==6); // works only for pair functions
2367 std::vector<long> map(NDIM);
2368
2369 // symmetric particle permutation
2370 if (symmetry=="sy_particle") {
2371 map[0]=3; map[1]=4; map[2]=5;
2372 map[3]=0; map[4]=1; map[5]=2;
2373 } else if (symmetry=="cx") {
2374 map[0]=0; map[1]=2; map[2]=1;
2375 map[3]=3; map[4]=5; map[5]=4;
2376
2377 } else if (symmetry=="cy") {
2378 map[0]=2; map[1]=1; map[2]=0;
2379 map[3]=5; map[4]=4; map[5]=3;
2380
2381 } else if (symmetry=="cz") {
2382 map[0]=1; map[1]=0; map[2]=2;
2383 map[3]=4; map[4]=3; map[5]=5;
2384
2385 } else {
2386 if (f.world().rank()==0) {
2387 print("unknown parameter in symmetrize:",symmetry);
2388 }
2389 MADNESS_EXCEPTION("unknown parameter in symmetrize",1);
2390 }
2391
2392 result.mapdim(f,map,true); // need to fence here
2393 result.get_impl()->average(*f.get_impl());
2394
2395 return result;
2396 }
2397
2398
2399
2400 /// multiply a high-dimensional function with a low-dimensional function
2401
2402 /// @param[in] f NDIM function of 2 particles: f=f(1,2)
2403 /// @param[in] g LDIM function of 1 particle: g=g(1) or g=g(2)
2404 /// @param[in] particle if g=g(1) or g=g(2)
2405 /// @return h(1,2) = f(1,2) * g(p)
2406 template<typename T, std::size_t NDIM, std::size_t LDIM>
2407 Function<T,NDIM> multiply(const Function<T,NDIM> f, const Function<T,LDIM> g, const int particle, const bool fence=true) {
2408
2409 static_assert(LDIM+LDIM==NDIM);
2411
2412 Function<T,NDIM> result;
2413 result.set_impl(f, false);
2414
2415// Function<T,NDIM>& ff = const_cast< Function<T,NDIM>& >(f);
2416// Function<T,LDIM>& gg = const_cast< Function<T,LDIM>& >(g);
2417
2418 f.change_tree_state(redundant,false);
2419 g.change_tree_state(redundant);
2420 FunctionImpl<T,NDIM>* fimpl=f.get_impl().get();
2421 FunctionImpl<T,LDIM>* gimpl=g.get_impl().get();
2422
2423 result.get_impl()->multiply(fimpl,gimpl,particle);
2424 result.world().gop.fence();
2425
2426 f.change_tree_state(reconstructed,false);
2427 g.change_tree_state(reconstructed);
2428 return result;
2429 }
2430
2431
2432 template <typename T, std::size_t NDIM>
2433 Function<T,NDIM>
2437 bool fence=true)
2438 {
2441 other.reconstruct();
2442 result.get_impl()->project(*other.get_impl(),fence);
2443 return result;
2444 }
2445
2446
2447 /// Computes the scalar/inner product between two functions
2448
2449 /// In Maple this would be \c int(conjugate(f(x))*g(x),x=-infinity..infinity)
2450 template <typename T, typename R, std::size_t NDIM>
2453 return f.inner(g);
2454 }
2455
2456
2457 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2458
2459 /// syntax similar to the inner product in tensor.h
2460 /// e.g result=inner<3>(f,g),{0},{1}) : r(x,y) = int f(x1,x) g(y,x1) dx1
2461 /// @param[in] task 0: everything, 1; prepare only (fence), 2: work only (no fence), 3: finalize only (fence)
2462 template<std::size_t NDIM, typename T, std::size_t LDIM, typename R, std::size_t KDIM,
2463 std::size_t CDIM = (KDIM + LDIM - NDIM) / 2>
2464 std::vector<Function<TENSOR_RESULT_TYPE(T, R), NDIM>>
2465 innerXX(const Function<T, LDIM>& f, const std::vector<Function<R, KDIM>>& vg, const std::array<int, CDIM> v1,
2466 const std::array<int, CDIM> v2, int task=0) {
2467 bool prepare = ((task==0) or (task==1));
2468 bool work = ((task==0) or (task==2));
2469 bool finish = ((task==0) or (task==3));
2470
2471 static_assert((KDIM + LDIM - NDIM) % 2 == 0, "faulty dimensions in inner (partial version)");
2472 static_assert(KDIM + LDIM - 2 * CDIM == NDIM, "faulty dimensions in inner (partial version)");
2473
2474 // contraction indices must be contiguous and either in the beginning or at the end
2475 for (int i=0; i<CDIM-1; ++i) MADNESS_CHECK((v1[i]+1)==v1[i+1]);
2476 MADNESS_CHECK((v1[0]==0) or (v1[CDIM-1]==LDIM-1));
2477
2478 for (int i=0; i<CDIM-1; ++i) MADNESS_CHECK((v2[i]+1)==v2[i+1]);
2479 MADNESS_CHECK((v2[0]==0) or (v2[CDIM-1]==KDIM-1));
2480
2481 MADNESS_CHECK(f.is_initialized());
2482 MADNESS_CHECK(vg[0].is_initialized());
2483 MADNESS_CHECK(f.world().id() == vg[0].world().id());
2484 // this needs to be run in a single world, so that all coefficients are local.
2485 // Use macrotasks if run on multiple processes.
2486 World& world=f.world();
2487 MADNESS_CHECK(world.size() == 1);
2488
2489 if (prepare) {
2490 f.change_tree_state(nonstandard);
2492 world.gop.fence();
2493 f.get_impl()->compute_snorm_and_dnorm(false);
2494 for (auto& g : vg) g.get_impl()->compute_snorm_and_dnorm(false);
2495 world.gop.fence();
2496 }
2497
2498 typedef TENSOR_RESULT_TYPE(T, R) resultT;
2499 std::vector<Function<resultT,NDIM>> result(vg.size());
2500 if (work) {
2501 world.gop.set_forbid_fence(true);
2502 for (int i=0; i<vg.size(); ++i) {
2503 result[i]=FunctionFactory<resultT,NDIM>(world)
2504 .k(f.k()).thresh(f.thresh()).empty().nofence();
2505 result[i].get_impl()->partial_inner(*f.get_impl(),*(vg[i]).get_impl(),v1,v2);
2506 result[i].get_impl()->set_tree_state(nonstandard_after_apply);
2507 }
2508 world.gop.set_forbid_fence(false);
2509 }
2510
2511 if (finish) {
2512
2513 world.gop.fence();
2514// result.get_impl()->reconstruct(true);
2515
2517// result.reconstruct();
2518 // restore initial state of g and h
2519 auto erase_list = [] (const auto& funcimpl) {
2520 typedef typename std::decay_t<decltype(funcimpl)>::keyT keyTT;
2521 std::list<keyTT> to_be_erased;
2522 for (auto it=funcimpl.get_coeffs().begin(); it!=funcimpl.get_coeffs().end(); ++it) {
2523 const auto& key=it->first;
2524 const auto& node=it->second;
2525 if (not node.has_children()) to_be_erased.push_back(key);
2526 }
2527 return to_be_erased;
2528 };
2529
2530 FunctionImpl<T,LDIM>& f_nc=const_cast<FunctionImpl<T,LDIM>&>(*f.get_impl());
2531 for (auto& key : erase_list(f_nc)) f_nc.get_coeffs().erase(key);
2532 for (auto& g : vg) {
2533 FunctionImpl<R,KDIM>& g_nc=const_cast<FunctionImpl<R,KDIM>&>(*g.get_impl());
2534 for (auto& key : erase_list(g_nc)) g_nc.get_coeffs().erase(key);
2535 }
2536 world.gop.fence();
2538 f_nc.reconstruct(false);
2539 world.gop.fence();
2540
2541 }
2542
2543 return result;
2544 }
2545
2546
2547 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2548
2549 /// syntax similar to the inner product in tensor.h
2550 /// e.g result=inner<3>(f,g),{0},{1}) : r(x,y) = int f(x1,x) g(y,x1) dx1
2551 /// @param[in] task 0: everything, 1; prepare only (fence), 2: work only (no fence), 3: finalize only (fence)
2552 template<std::size_t NDIM, typename T, std::size_t LDIM, typename R, std::size_t KDIM,
2553 std::size_t CDIM = (KDIM + LDIM - NDIM) / 2>
2554 Function<TENSOR_RESULT_TYPE(T, R), NDIM>
2555 innerXX(const Function<T, LDIM>& f, const Function<R, KDIM>& g, const std::array<int, CDIM> v1,
2556 const std::array<int, CDIM> v2, int task=0) {
2557 return innerXX<NDIM,T,LDIM,R,KDIM>(f,std::vector<Function<R,KDIM>>({g}),v1,v2,task)[0];
2558 }
2559
2560 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2561
2562 /// syntax similar to the inner product in tensor.h
2563 /// e.g result=inner<3>(f,g),{0},{1}) : r(x,y) = int f(x1,x) g(y,x1) dx1
2564 template <typename T, std::size_t LDIM, typename R, std::size_t KDIM>
2565 Function<TENSOR_RESULT_TYPE(T,R),KDIM+LDIM-2>
2566 inner(const Function<T,LDIM>& f, const Function<R,KDIM>& g, const std::tuple<int> v1, const std::tuple<int> v2) {
2567 return innerXX<KDIM+LDIM-2>(f,g,
2568 std::array<int,1>({std::get<0>(v1)}),
2569 std::array<int,1>({std::get<0>(v2)}));
2570 }
2571
2572 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2573
2574 /// syntax similar to the inner product in tensor.h
2575 /// e.g result=inner<3>(f,g),{0,1},{1,2}) : r(y) = int f(x1,x2) g(y,x1,x2) dx1 dx2
2576 template <typename T, std::size_t LDIM, typename R, std::size_t KDIM>
2577 Function<TENSOR_RESULT_TYPE(T,R),KDIM+LDIM-4>
2578 inner(const Function<T,LDIM>& f, const Function<R,KDIM>& g, const std::tuple<int,int> v1, const std::tuple<int,int> v2) {
2579 return innerXX<KDIM+LDIM-4>(f,g,
2580 std::array<int,2>({std::get<0>(v1),std::get<1>(v1)}),
2581 std::array<int,2>({std::get<0>(v2),std::get<1>(v2)}));
2582 }
2583
2584 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2585
2586 /// syntax similar to the inner product in tensor.h
2587 /// e.g result=inner<3>(f,g),{1},{2}) : r(x,y,z) = int f(x,x1) g(y,z,x1) dx1
2588 template <typename T, std::size_t LDIM, typename R, std::size_t KDIM>
2589 Function<TENSOR_RESULT_TYPE(T,R),KDIM+LDIM-6>
2590 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) {
2591 return innerXX<KDIM+LDIM-6>(f,g,
2592 std::array<int,3>({std::get<0>(v1),std::get<1>(v1),std::get<2>(v1)}),
2593 std::array<int,3>({std::get<0>(v2),std::get<1>(v2),std::get<2>(v2)}));
2594 }
2595
2596
2597
2598 /// Computes the scalar/inner product between an MRA function and an external functor
2599
2600 /// Currently this defaults to inner_adaptive, which might be more expensive
2601 /// than inner_ext since it loops over all leaf nodes. If you feel inner_ext
2602 /// is more efficient you need to call it directly
2603 /// @param[in] f MRA function
2604 /// @param[in] g functor
2605 /// @result inner(f,g)
2606 template <typename T, typename opT, std::size_t NDIM>
2607 TENSOR_RESULT_TYPE(T,typename opT::value_type) inner(const Function<T,NDIM>& f, const opT& g) {
2609 std::shared_ptr< FunctionFunctorInterface<double,3> > func(new opT(g));
2610 return f.inner_adaptive(func);
2611 }
2612
2613 /// Computes the scalar/inner product between an MRA function and an external functor
2614
2615 /// Currently this defaults to inner_adaptive, which might be more expensive
2616 /// than inner_ext since it loops over all leaf nodes. If you feel inner_ext
2617 /// is more efficient you need to call it directly
2618 /// @param[in] g functor
2619 /// @param[in] f MRA function
2620 /// @result inner(f,g)
2621 template <typename T, typename opT, std::size_t NDIM>
2622 TENSOR_RESULT_TYPE(T,typename opT::value_type) inner(const opT& g, const Function<T,NDIM>& f) {
2623 return inner(f,g);
2624 }
2625
2626 template <typename T, typename R, std::size_t NDIM>
2627 typename IsSupported<TensorTypeData<R>, Function<TENSOR_RESULT_TYPE(T,R),NDIM> >::type
2629 return (f*R(1.0)).add_scalar(r);
2630 }
2631
2632 template <typename T, typename R, std::size_t NDIM>
2633 typename IsSupported<TensorTypeData<R>, Function<TENSOR_RESULT_TYPE(T,R),NDIM> >::type
2635 return (f*R(1.0)).add_scalar(r);
2636 }
2637
2638 template <typename T, typename R, std::size_t NDIM>
2639 typename IsSupported<TensorTypeData<R>, Function<TENSOR_RESULT_TYPE(T,R),NDIM> >::type
2641 return (f*R(1.0)).add_scalar(-r);
2642 }
2643
2644 template <typename T, typename R, std::size_t NDIM>
2645 typename IsSupported<TensorTypeData<R>, Function<TENSOR_RESULT_TYPE(T,R),NDIM> >::type
2647 return (f*R(-1.0)).add_scalar(r);
2648 }
2649
2650 namespace detail {
2651 template <std::size_t NDIM>
2652 struct realop {
2653 typedef double resultT;
2655 return real(t);
2656 }
2657
2658 template <typename Archive> void serialize (Archive& ar) {}
2659 };
2660
2661 template <std::size_t NDIM>
2662 struct imagop {
2663 typedef double resultT;
2665 return imag(t);
2666 }
2667
2668 template <typename Archive> void serialize (Archive& ar) {}
2669 };
2670
2671 template <std::size_t NDIM>
2672 struct abssqop {
2673 typedef double resultT;
2675 Tensor<double> r = abs(t);
2676 return r.emul(r);
2677 }
2678
2679 template <typename Archive> void serialize (Archive& ar) {}
2680 };
2681
2682 template <std::size_t NDIM>
2683 struct absop {
2684 typedef double resultT;
2686 Tensor<double> r = abs(t);
2687 return r;
2688 }
2689
2690 template <typename Archive> void serialize (Archive& ar) {}
2691 };
2692
2693 }
2694
2695 /// Returns a new function that is the real part of the input
2696 template <std::size_t NDIM>
2698 return unary_op_coeffs(z, detail::realop<NDIM>(), fence);
2699 }
2700
2701 /// Returns a new function that is the real part of the input
2702 template <std::size_t NDIM>
2704 return copy(z);
2705 }
2706
2707 /// Returns a new function that is the imaginary part of the input
2708 template <std::size_t NDIM>
2710 return unary_op_coeffs(z, detail::imagop<NDIM>(), fence);
2711 }
2712
2713
2714 /// Create a new function that is the square of f - global comm only if not reconstructed
2715 template <typename T, std::size_t NDIM>
2716 Function<T,NDIM> square(const Function<T,NDIM>& f, bool fence=true) {
2718 Function<T,NDIM> result = copy(f,true); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2719 return result.square(true); //fence); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2720 }
2721
2722 /// Create a new function that is the abs of f - global comm only if not reconstructed
2723 template <typename T, std::size_t NDIM>
2724 Function<T,NDIM> abs(const Function<T,NDIM>& f, bool fence=true) {
2726 Function<T,NDIM> result = copy(f,true); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2727 return result.abs(true); //fence); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2728 }
2729
2730 /// Create a new function that is the abs_square of f - global comm only if not reconstructed
2731 template <typename T, std::size_t NDIM>
2732 typename std::enable_if<!TensorTypeData<T>::iscomplex, Function<T,NDIM> >::type
2733 abs_square(const Function<T,NDIM>& f, bool fence=true) {
2735 Function<T,NDIM> result = copy(f,true); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2736 return result.abs_square(true); //fence); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2737 }
2738
2739 /// Create a new function that is the abs_square of f - global comm only if not reconstructed
2740 template <typename T, std::size_t NDIM>
2741 typename std::enable_if<TensorTypeData<T>::iscomplex, Function<typename Tensor<T>::scalar_type,NDIM> >::type
2742 abs_square(const Function<T,NDIM>& f, bool fence=true) {
2743 return unary_op(f, detail::abssqop<NDIM>(), fence);
2744 }
2745
2746 /// Returns a new function that is the square of the absolute value of the input
2747 template <std::size_t NDIM>
2749 return unary_op(z, detail::abssqop<NDIM>(), fence);
2750 }
2751
2752 /// Returns a new function that is the absolute value of the input
2753 template <std::size_t NDIM>
2755 return unary_op(z, detail::absop<NDIM>(), fence);
2756 }
2757
2758}
2759
2760#include <madness/mra/funcplot.h>
2761
2762namespace madness {
2763 namespace archive {
2764 template <class archiveT, class T, std::size_t NDIM>
2766 static inline void load(const ParallelInputArchive<archiveT>& ar, Function<T,NDIM>& f) {
2767 f.load(*ar.get_world(), ar);
2768 }
2769 };
2770
2771 template <class archiveT, class T, std::size_t NDIM>
2773 static inline void store(const ParallelOutputArchive<archiveT>& ar, const Function<T,NDIM>& f) {
2774 f.store(ar);
2775 }
2776 };
2777 }
2778
2779 template <class T, std::size_t NDIM>
2780 void save(const Function<T,NDIM>& f, const std::string name) {
2782 ar2 & f;
2783 }
2784
2785 template <class T, std::size_t NDIM>
2786 void load(Function<T,NDIM>& f, const std::string name) {
2788 ar2 & f;
2789 }
2790
2791}
2792
2793namespace madness {
2794 // type traits to check if a template parameter is a Function
2795 template<typename>
2796 struct is_madness_function : std::false_type {};
2797
2798 template<typename T, std::size_t NDIM>
2799 struct is_madness_function<madness::Function<T, NDIM>> : std::true_type {};
2800}
2801
2802
2803/* @} */
2804
2805#include <madness/mra/derivative.h>
2806#include <madness/mra/operator.h>
2808#include <madness/mra/vmra.h>
2809// #include <madness/mra/mraimpl.h> !!!!!!!!!!!!! NOOOOOOOOOOOOOOOOOOOOOOOOOOOOOOOOOO !!!!!!!!!!!!!!!!!!
2810
2811#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:176
FunctionFactory implements the named-parameter idiom for Function.
Definition function_factory.h:86
FunctionFactory & nofence()
Definition function_factory.h:275
virtual FunctionFactory & thresh(double thresh)
Definition function_factory.h:191
FunctionFactory & fence(bool fence=true)
Definition function_factory.h:269
virtual FunctionFactory & k(int k)
Definition function_factory.h:186
FunctionFactory & empty()
Definition function_factory.h:239
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:322
void reconstruct(bool fence)
reconstruct this tree – respects fence
Definition mraimpl.h:1495
bool is_on_demand() const
Definition mraimpl.h:262
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:889
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:350
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:1677
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:1393
T inner_adaptive(const std::shared_ptr< FunctionFunctorInterface< T, NDIM > > f, const bool leaf_refine=true) const
Definition mra.h:1378
void unaryop_coeff(const opT &op, bool fence=true)
Unary operation applied inplace to the coefficients.
Definition mra.h:939
bool is_compressed() const
Returns true if compressed, false otherwise. No communication.
Definition mra.h:482
Function< T, NDIM/2 > dirac_convolution(const bool fence=true) const
Definition mra.h:1457
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:667
bool autorefine() const
Returns value of autorefine flag. No communication.
Definition mra.h:573
Function< T, NDIM > & add_scalar(T t, bool fence=true)
Inplace add scalar. No communication except for optional fence.
Definition mra.h:988
void print_size(const std::string name) const
print some info about this
Definition mra.h:518
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:1153
void broaden(const BoundaryConditions< NDIM > &bc=FunctionDefaults< NDIM >::get_bc(), bool fence=true) const
Inplace broadens support in scaling function basis.
Definition mra.h:861
void print_info() const
Print a summary of the load balancing info.
Definition mra.h:905
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:978
void norm_tree(bool fence=true) const
Initializes information about the function norm at all length scales.
Definition mra.h:726
void load(World &world, Archive &ar)
Replaces this function with one loaded from an archive using the default processor map.
Definition mra.h:1474
Function< T, NDIM > & abs_square(bool fence=true)
Returns *this for chaining.
Definition mra.h:1105
double norm2sq_local() const
Returns the square of the norm of the local function ... no communication.
Definition mra.h:703
IsSupported< TensorTypeData< Q >, Function< T, NDIM > >::type & operator*=(const Q q)
Inplace scaling by a constant.
Definition mra.h:1077
void sum_down(bool fence=true) const
Sums scaling coeffs down tree restoring state with coeffs only at leaves. Optional fence....
Definition mra.h:823
Level depthpt(const coordT &xuser) const
Definition mra.h:426
Function< T, NDIM > & operator+=(const Function< Q, NDIM > &other)
Inplace addition of functions in the wavelet basis.
Definition mra.h:1036
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:583
World & world() const
Returns the world.
Definition mra.h:673
T trace() const
Returns global value of int(f(x),x) ... global comm required.
Definition mra.h:1127
T typeT
Definition mra.h:157
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:1171
const Function< T, NDIM > & change_tree_state(const TreeState finalstate, bool fence=true) const
changes tree state to given state
Definition mra.h:810
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:1515
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:895
double norm2() const
Returns the 2-norm of the function ... global sum ... works in either basis.
Definition mra.h:713
Function< T, NDIM > & fill_cuspy_tree(const opT &op, const bool fence=true)
Definition mra.h:1194
void change_tensor_type(const TensorArgs &targs, bool fence=true)
change the tensor type of the coefficients in the FunctionNode
Definition mra.h:1507
const Function< T, NDIM > & refine(bool fence=true) const
Inplace autorefines the function using same test as for squaring.
Definition mra.h:855
TreeState gstate
Definition mra.h:1310
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:391
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:1087
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:1182
void verify_tree() const
Verifies the tree data structure ... global sync implied.
Definition mra.h:473
Future< Level > evaldepthpt(const coordT &xuser) const
Definition mra.h:272
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:1263
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:1643
impl world gop fence()
int k() const
Returns the number of multiwavelets (k). No communication.
Definition mra.h:611
const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > & get_pmap() const
Returns a shared pointer to the process map.
Definition mra.h:681
double thresh() const
Returns value of truncation threshold. No communication.
Definition mra.h:592
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:1629
void set_thresh(double value, bool fence=true)
Sets the value of the truncation threshold. Optional global fence.
Definition mra.h:602
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:1537
void distribute(std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > newmap) const
distribute this function according to newmap
Definition mra.h:694
void unaryop(const opT &op, bool fence=true)
Inplace unary operation on function values.
Definition mra.h:929
const std::shared_ptr< FunctionImpl< T, NDIM > > & get_impl() const
Returns a shared-pointer to the implementation.
Definition mra.h:639
void standard(bool fence=true)
Converts the function standard compressed form. Possible non-blocking comm.
Definition mra.h:773
void unaryop_node(const opT &op, bool fence=true)
Unary operation applied inplace to the nodes.
Definition mra.h:949
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:627
Function< T, NDIM > & fill_cuspy_tree(const bool fence=true)
Special refinement on 6D boxes where the electrons come close (meet)
Definition mra.h:1207
std::size_t size() const
Returns the number of coefficients in the function ... collective global sum.
Definition mra.h:558
bool is_on_demand() const
Definition mra.h:661
Function< T, NDIM > & operator-=(const Function< Q, NDIM > &other)
Inplace subtraction of functions in the wavelet basis.
Definition mra.h:1056
const Function< T, NDIM > & reconstruct(bool fence=true) const
Reconstructs the function, transforming into scaling function basis. Possible non-blocking comm.
Definition mra.h:800
Function< T, NDIM > & abs(bool fence=true)
Returns *this for chaining.
Definition mra.h:1096
Function< T, NDIM > & reduce_rank(const double thresh=0.0, const bool fence=true)
reduce the rank of the coefficient tensors
Definition mra.h:1707
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:1658
std::size_t max_nodes() const
Returns the max number of nodes on a processor.
Definition mra.h:543
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:1249
std::shared_ptr< FunctionImpl< T, NDIM > > impl
Definition mra.h:146
void replicate(bool fence=true) const
replicate this function, generating a unique pmap
Definition mra.h:688
T trace_local() const
Returns local contribution to int(f(x),x) ... no communication.
Definition mra.h:1118
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:1600
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:1006
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:1495
std::size_t max_local_depth() const
Returns the maximum local depth of the function tree ... no communications.
Definition mra.h:535
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:1618
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:527
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:511
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:1433
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:1528
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:872
void refine_general(const opT &op, bool fence=true) const
Inplace autorefines the function. Optional fence. Possible non-blocking comm.
Definition mra.h:836
static void doconj(const Key< NDIM >, Tensor< T > &t)
Definition mra.h:959
TreeState state
Definition mra.h:1309
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:240
void set_functor(const std::shared_ptr< FunctionFunctorInterface< T, NDIM > > functor)
Replace the current functor with the provided new one.
Definition mra.h:656
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:1222
return local
Definition mra.h:1334
auto func
Definition mra.h:1412
~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:646
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:405
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:1345
std::size_t min_nodes() const
Returns the min number of nodes on a processor.
Definition mra.h:550
constexpr std::size_t LDIM
Definition mra.h:1411
static constexpr std::size_t dimT
Definition mra.h:158
bool is_nonstandard() const
Returns true if nonstandard-compressed, false otherwise. No communication.
Definition mra.h:504
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:460
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:1361
void make_redundant(bool fence=true)
Converts the function to redundant form, i.e. sum coefficients on all levels.
Definition mra.h:785
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:1235
double check_symmetry() const
check symmetry of a function by computing the 2nd derivative
Definition mra.h:1693
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:1556
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:565
const Function< T, NDIM > & compress(bool fence=true) const
Compresses the function, transforming into wavelet basis. Possible non-blocking comm.
Definition mra.h:746
bool is_initialized() const
Returns true if the function is initialized.
Definition mra.h:167
change_tree_state(state, false)
bool is_reconstructed() const
Returns true if reconstructed, false otherwise. No communication.
Definition mra.h:493
Function< T, NDIM > & chop_at_level(const int n, const bool fence=true)
remove all nodes with level higher than n
Definition mra.h:1715
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:1577
double errsq_local(const funcT &func) const
Returns an estimate of the difference ||this-func||^2 from local data.
Definition mra.h:445
void make_nonstandard(bool keepleaves, bool fence=true) const
Compresses the function retaining scaling function coeffs. Possible non-blocking comm.
Definition mra.h:759
Future< T > eval(const coordT &xuser) const
Evaluates the function at a point in user coordinates. Possible non-blocking comm.
Definition mra.h:206
void print_tree(std::ostream &os=std::cout) const
Process 0 prints a summary of all nodes in the tree (collective)
Definition mra.h:882
FunctionNode< T, NDIM > nodeT
Definition mra.h:154
Future< long > evalR(const coordT &xuser) const
Evaluates the function rank at a point in user coordinates. Possible non-blocking comm.
Definition mra.h:309
Function(const factoryT &factory)
Constructor from FunctionFactory provides named parameter idiom. Possible non-blocking communication.
Definition mra.h:178
void unaryop(T(*f)(T))
Inplace unary operation on function values.
Definition mra.h:920
Function< T, NDIM > conj(bool fence=true)
Inplace complex conjugate. No communication except for optional fence.
Definition mra.h:967
A future is a possibly yet unevaluated value.
Definition future.h:373
remote_refT remote_ref(World &world) const
Returns a structure used to pass references to another process.
Definition future.h:675
Key is the index for a node of the 2^NDIM-tree.
Definition key.h:68
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:1798
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:1105
Interface to be provided by any process map.
Definition worlddc.h:82
void fence(bool debug=false)
Synchronizes all processes in communicator AND globally ensures no pending AM or tasks.
Definition worldgop.cc:161
bool set_forbid_fence(bool value)
Set forbid_fence flag to new value and return old value.
Definition worldgop.h:674
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
unsigned long id() const
Definition world.h:315
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:34
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
Namespace for all elements and tools of MADNESS.
Definition DFParameters.h:10
Function< T, NDIM > map_and_mirror(const Function< T, NDIM > &f, const std::vector< long > &map, const std::vector< long > &mirror, bool fence=true)
This is replaced with mirror(map(f)), optional fence.
Definition mra.h:2323
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:2748
Function< TENSOR_RESULT_TYPE(L, R), NDIM > gaxpy_oop(TENSOR_RESULT_TYPE(L, R) alpha, const Function< L, NDIM > &left, TENSOR_RESULT_TYPE(L, R) beta, const Function< R, NDIM > &right, bool fence=true)
Returns new function alpha*left + beta*right optional fence and no automatic compression.
Definition mra.h:1952
Function< typename TensorTypeData< Q >::scalar_type, NDIM > abs_square(const Function< Q, NDIM > &func)
Definition complexfun.h:121
Function< T, NDIM > mapdim(const Function< T, NDIM > &f, const std::vector< long > &map, bool fence=true)
Generate a new function by reordering dimensions ... optional fence.
Definition mra.h:2294
const std::vector< Function< T, NDIM > > & change_tree_state(const std::vector< Function< T, NDIM > > &v, const TreeState finalstate, const bool fence=true)
change tree state of the functions
Definition vmra.h:252
Function< T, NDIM > square(const Function< T, NDIM > &f, bool fence=true)
Create a new function that is the square of f - global comm only if not reconstructed.
Definition mra.h:2716
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:2266
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:2006
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:1801
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:2067
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:1736
std::vector< std::shared_ptr< FunctionImpl< T, NDIM > > > get_impl(const std::vector< Function< T, NDIM > > &v)
Definition vmra.h:638
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:2337
TreeState
Definition funcdefaults.h:59
@ nonstandard_after_apply
s and d coeffs, state after operator application
Definition funcdefaults.h:64
@ redundant_after_merge
s coeffs everywhere, must be summed up to yield the result
Definition funcdefaults.h:66
@ reconstructed
s coeffs at the leaves only
Definition funcdefaults.h:60
@ nonstandard
s and d coeffs in internal nodes
Definition funcdefaults.h:62
@ unknown
Definition funcdefaults.h:68
@ compressed
d coeffs in internal nodes, s and d coeffs at the root
Definition funcdefaults.h:61
@ redundant
s coeffs everywhere
Definition funcdefaults.h:65
@ nonstandard_with_leaves
like nonstandard, with s coeffs at the leaves
Definition funcdefaults.h:63
Function< T, NDIM > conj(const Function< T, NDIM > &f, bool fence=true)
Return the complex conjugate of the input function with the same distribution and optional fence.
Definition mra.h:2081
static void user_to_sim(const Vector< double, NDIM > &xuser, Vector< double, NDIM > &xsim)
Convert user coords (cell[][]) to simulation coords ([0,1]^ndim)
Definition funcdefaults.h:425
Function< T, NDIM > multiop_values(const opT &op, const std::vector< Function< T, NDIM > > &vf)
Definition mra.h:1726
void standard(World &world, std::vector< Function< T, NDIM > > &v, bool fence=true)
Generates standard form of a vector of functions.
Definition vmra.h:239
const std::vector< Function< T, NDIM > > & reconstruct(const std::vector< Function< T, NDIM > > &v)
reconstruct a vector of functions
Definition vmra.h:158
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:2465
Function< T, NDIM > mirror(const Function< T, NDIM > &f, const std::vector< long > &mirrormap, bool fence=true)
Generate a new function by mirroring within the dimensions .. optional fence.
Definition mra.h:2306
std::vector< CCPairFunction< T, NDIM > > operator*(const double fac, const std::vector< CCPairFunction< T, NDIM > > &arg)
Definition ccpairfunction.h:1089
std::shared_ptr< FunctionFunctorInterface< double, 3 > > func(new opT(g))
int Level
Definition key.h:57
Function< typename opT::resultT, NDIM > unary_op_coeffs(const Function< Q, NDIM > &func, const opT &op, bool fence=true)
Out of place application of unary operation to scaling function coefficients with optional fence.
Definition mra.h:1828
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:1777
std::string get_mra_data_dir()
Definition startup.cc:208
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:1970
void print(const T &t, const Ts &... ts)
Print items to std::cout (items separated by spaces) and terminate with a new line.
Definition print.h:225
response_space apply(World &world, std::vector< std::vector< std::shared_ptr< real_convolution_3d > > > &op, response_space &f)
Definition basic_operators.cc:39
@ TT_2D
Definition gentensor.h:120
NDIM & f
Definition mra.h:2451
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:1962
Function< T, NDIM > symmetrize(const Function< T, NDIM > &f, const std::string symmetry, bool fence=true)
symmetrize a function
Definition mra.h:2363
NDIM const Function< R, NDIM > & g
Definition mra.h:2451
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:2155
double inner(response_space &a, response_space &b)
Definition response_functions.h:442
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:1815
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
std::vector< Function< TENSOR_RESULT_TYPE(L, R), D > > vmulXX(const Function< L, D > &left, const std::vector< Function< R, D > > &vright, double tol, bool fence=true)
Use the vmra/mul(...) interface instead.
Definition mra.h:1842
static bool print_timings
Definition SCF.cc:104
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:2407
void load(Function< T, NDIM > &f, const std::string name)
Definition mra.h:2786
Function< T, NDIM > project(const Function< T, NDIM > &other, int k=FunctionDefaults< NDIM >::get_k(), double thresh=FunctionDefaults< NDIM >::get_thresh(), bool fence=true)
Definition mra.h:2434
double real(double x)
Definition complexfun.h:52
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:2780
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:1867
Function< T, NDIM > copy(const Function< T, NDIM > &f, const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > &pmap, bool fence=true)
Create a new copy of the function with different distribution and optional fence.
Definition mra.h:2037
void gaxpy(const double a, ScalarResult< T > &left, const double b, const T &right, const bool fence=true)
the result type of a macrotask must implement gaxpy
Definition macrotaskq.h:246
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
static const double thresh
Definition rk.cc:45
static const long k
Definition rk.cc:44
Definition test_ar.cc:204
Definition leafop.h:142
void serialize(Archive &ar)
Definition mra.h:916
T(* f)(T)
Definition mra.h:911
void operator()(const Key< NDIM > &key, Tensor< T > &t) const
Definition mra.h:913
SimpleUnaryOpWrapper(T(*f)(T))
Definition mra.h:912
void serialize(Archive &ar)
Definition mra.h:849
bool operator()(implT *impl, const Key< NDIM > &key, const nodeT &t) const
Definition mra.h:845
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:2766
Default load of an object via serialize(ar, t).
Definition archive.h:666
static void store(const ParallelOutputArchive< archiveT > &ar, const Function< T, NDIM > &f)
Definition mra.h:2773
Default store of an object via serialize(ar, t).
Definition archive.h:611
Definition mra.h:2683
Tensor< double > operator()(const Key< NDIM > &key, const Tensor< double_complex > &t) const
Definition mra.h:2685
double resultT
Definition mra.h:2684
void serialize(Archive &ar)
Definition mra.h:2690
Definition mra.h:2672
void serialize(Archive &ar)
Definition mra.h:2679
Tensor< double > operator()(const Key< NDIM > &key, const Tensor< double_complex > &t) const
Definition mra.h:2674
double resultT
Definition mra.h:2673
Definition mra.h:2662
Tensor< double > operator()(const Key< NDIM > &key, const Tensor< double_complex > &t) const
Definition mra.h:2664
void serialize(Archive &ar)
Definition mra.h:2668
double resultT
Definition mra.h:2663
Definition mra.h:2652
Tensor< double > operator()(const Key< NDIM > &key, const Tensor< double_complex > &t) const
Definition mra.h:2654
double resultT
Definition mra.h:2653
void serialize(Archive &ar)
Definition mra.h:2658
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:2796
Definition mra.h:112
Definition mra.h:115
Definition funcimpl.h:564
Definition lowrankfunction.h:332
Defines and implements most of Tensor.
#define UNARY_OPTIMIZED_ITERATOR(X, x, exp)
Definition tensor_macros.h:658
AtomicInt sum
Definition test_atomicint.cc:46
double norm(const T i1)
Definition test_cloud.cc:72
int task(int i)
Definition test_runtime.cpp:4
void e()
Definition test_sig.cc:75
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