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 /// Returns true if redundant, false otherwise. No communication.
510
511 /// If the function is not initialized, returns false.
512 bool is_redundant() const {
514 return impl ? impl->is_redundant() : false;
515 }
516
517 /// Returns true if redundant_after_merge, false otherwise. No communication.
518
519 /// If the function is not initialized, returns false.
522 return impl ? impl->is_redundant_after_merge() : false;
523 }
524
525 /// Returns the number of nodes in the function tree ... collective global sum
526 std::size_t tree_size() const {
528 if (!impl) return 0;
529 return impl->tree_size();
530 }
531
532 /// print some info about this
533 void print_size(const std::string name) const {
534 if (!impl) {
535 print("function",name,"not assigned yet");
536 } else {
537 impl->print_size(name);
538 }
539 }
540
541 /// Returns the maximum depth of the function tree ... collective global sum
542 std::size_t max_depth() const {
544 if (!impl) return 0;
545 return impl->max_depth();
546 }
547
548
549 /// Returns the maximum local depth of the function tree ... no communications
550 std::size_t max_local_depth() const {
552 if (!impl) return 0;
553 return impl->max_local_depth();
554 }
555
556
557 /// Returns the max number of nodes on a processor
558 std::size_t max_nodes() const {
560 if (!impl) return 0;
561 return impl->max_nodes();
562 }
563
564 /// Returns the min number of nodes on a processor
565 std::size_t min_nodes() const {
567 if (!impl) return 0;
568 return impl->min_nodes();
569 }
570
571
572 /// Returns the number of coefficients in the function ... collective global sum
573 std::size_t size() const {
575 if (!impl) return 0;
576 return impl->size();
577 }
578
579 /// Return the number of coefficients in the function on this processor
580 std::size_t size_local() const {
582 if (!impl) return 0;
583 return impl->size_local();
584 }
585
586
587 /// Returns value of autorefine flag. No communication.
588 bool autorefine() const {
590 if (!impl) return true;
591 return impl->get_autorefine();
592 }
593
594
595 /// Sets the value of the autorefine flag. Optional global fence.
596
597 /// A fence is required to ensure consistent global state.
598 void set_autorefine(bool value, bool fence = true) {
600 verify();
601 impl->set_autorefine(value);
602 if (fence) impl->world.gop.fence();
603 }
604
605
606 /// Returns value of truncation threshold. No communication.
607 double thresh() const {
609 if (!impl) return 0.0;
610 return impl->get_thresh();
611 }
612
613
614 /// Sets the value of the truncation threshold. Optional global fence.
615
616 /// A fence is required to ensure consistent global state.
617 void set_thresh(double value, bool fence = true) {
619 verify();
620 impl->set_thresh(value);
621 if (fence) impl->world.gop.fence();
622 }
623
624
625 /// Returns the number of multiwavelets (k). No communication.
626 int k() const {
628 verify();
629 return impl->get_k();
630 }
631
632
633 /// Truncate the function with optional fence. Compresses with fence if not compressed.
634
635 /// If the truncation threshold is less than or equal to zero the default value
636 /// specified when the function was created is used.
637 /// If the function is not initialized, it just returns.
638 ///
639 /// Returns this for chaining.
640 /// @param[in] tol Tolerance for truncating the coefficients. Default 0.0 means use the implementation's member value \c thresh instead.
641 /// @param[in] fence Do fence
642 Function<T,NDIM>& truncate(double tol = 0.0, bool fence = true) {
644 if (!impl) return *this;
645 verify();
646// if (!is_compressed()) compress();
647 impl->truncate(tol,fence);
649 return *this;
650 }
651
652
653 /// Returns a shared-pointer to the implementation
654 const std::shared_ptr< FunctionImpl<T,NDIM> >& get_impl() const {
656 verify();
657 return impl;
658 }
659
660 /// Replace current FunctionImpl with provided new one
661 void set_impl(const std::shared_ptr< FunctionImpl<T,NDIM> >& impl) {
663 this->impl = impl;
664 }
665
666
667 /// Replace the current functor with the provided new one
668
669 /// presumably the new functor will be a CompositeFunctor, which will
670 /// change the behavior of the function: multiply the functor with the function
671 void set_functor(const std::shared_ptr<FunctionFunctorInterface<T, NDIM> > functor) {
672 this->impl->set_functor(functor);
673 print("set functor in mra.h");
674 }
675
676 bool is_on_demand() const {return this->impl->is_on_demand();}
677
678 /// Replace current FunctionImpl with a new one using the same parameters & map as f
679
680 /// If zero is true the function is initialized to zero, otherwise it is empty
681 template <typename R>
682 void set_impl(const Function<R,NDIM>& f, bool zero = true) {
683 impl = std::shared_ptr<implT>(new implT(*f.get_impl(), f.get_pmap(), zero));
684 if (zero) world().gop.fence();
685 }
686
687 /// Returns the world
688 World& world() const {
690 verify();
691 return impl->world;
692 }
693
694
695 /// Returns a shared pointer to the process map
696 const std::shared_ptr< WorldDCPmapInterface< Key<NDIM> > >& get_pmap() const {
698 verify();
699 return impl->get_pmap();
700 }
701
702 /// replicate this function, generating a unique pmap
703 void replicate(bool fence=true) const {
704 verify();
705 impl->replicate(fence);
706 }
707
708 /// distribute this function according to newmap
709 void distribute(std::shared_ptr< WorldDCPmapInterface< Key<NDIM> > > newmap) const {
710 verify();
711 impl->distribute(newmap);
712 }
713
714
715 /// Returns the square of the norm of the local function ... no communication
716
717 /// Works in either basis
718 double norm2sq_local() const {
720 verify();
722 "function must be compressed or reconstructed for norm2sq_local");
723 return impl->norm2sq_local();
724 }
725
726
727 /// Returns the 2-norm of the function ... global sum ... works in either basis
728
729 /// See comments for err() w.r.t. applying to many functions.
730 double norm2() const {
732 verify();
734 double local = impl->norm2sq_local();
735
736 impl->world.gop.sum(local);
737 impl->world.gop.fence();
738 return sqrt(local);
739 }
740
741
742 /// Initializes information about the function norm at all length scales
743 void norm_tree(bool fence = true) const {
745 verify();
748 const_cast<Function<T,NDIM>*>(this)->impl->norm_tree(fence);
749 }
750
751
752 /// Compresses the function, transforming into wavelet basis. Possible non-blocking comm.
753
754 /// By default fence=true meaning that this operation completes before returning,
755 /// otherwise if fence=false it returns without fencing and the user must invoke
756 /// world.gop.fence() to assure global completion before using the function
757 /// for other purposes.
758 ///
759 /// Noop if already compressed or if not initialized.
760 ///
761 /// Since reconstruction/compression do not discard information we define them
762 /// as const ... "logical constness" not "bitwise constness".
763 const Function<T,NDIM>& compress(bool fence = true) const {
765 }
766
767
768 /// Compresses the function retaining scaling function coeffs. Possible non-blocking comm.
769
770 /// By default fence=true meaning that this operation completes before returning,
771 /// otherwise if fence=false it returns without fencing and the user must invoke
772 /// world.gop.fence() to assure global completion before using the function
773 /// for other purposes.
774 ///
775 /// Noop if already compressed or if not initialized.
776 void make_nonstandard(bool keepleaves, bool fence=true) const {
778 if (keepleaves) newstate=nonstandard_with_leaves;
779 change_tree_state(newstate,fence);
780 }
781
782 /// Converts the function standard compressed form. Possible non-blocking comm.
783
784 /// By default fence=true meaning that this operation completes before returning,
785 /// otherwise if fence=false it returns without fencing and the user must invoke
786 /// world.gop.fence() to assure global completion before using the function
787 /// for other purposes.
788 ///
789 /// Must be already compressed.
790 void standard(bool fence = true) {
792 }
793
794 /// Converts the function to redundant form, i.e. sum coefficients on all levels
795
796 /// By default fence=true meaning that this operation completes before returning,
797 /// otherwise if fence=false it returns without fencing and the user must invoke
798 /// world.gop.fence() to assure global completion before using the function
799 /// for other purposes.
800 ///
801 /// Must be already compressed.
802 void make_redundant(bool fence = true) {
804 }
805
806 /// Reconstructs the function, transforming into scaling function basis. Possible non-blocking comm.
807
808 /// By default fence=true meaning that this operation completes before returning,
809 /// otherwise if fence=false it returns without fencing and the user must invoke
810 /// world.gop.fence() to assure global completion before using the function
811 /// for other purposes.
812 ///
813 /// Noop if already reconstructed or if not initialized.
814 ///
815 /// Since reconstruction/compression do not discard information we define them
816 /// as const ... "logical constness" not "bitwise constness".
817 const Function<T,NDIM>& reconstruct(bool fence = true) const {
819 }
820
821 /// changes tree state to given state
822
823 /// Since reconstruction/compression do not discard information we define them
824 /// as const ... "logical constness" not "bitwise constness".
825 /// @param[in] finalstate The final state of the tree
826 /// @param[in] fence Fence after the operation (might not be respected!!!)
827 const Function<T,NDIM>& change_tree_state(const TreeState finalstate, bool fence = true) const {
829 if (not impl) return *this;
830 TreeState current_state = impl->get_tree_state();
831 if (finalstate == current_state) return *this;
832 MADNESS_CHECK_THROW(current_state != TreeState::unknown, "unknown tree state");
833
834 impl->change_tree_state(finalstate, fence);
835 if (fence && VERIFY_TREE) verify_tree();
836 return *this;
837 }
838
839 /// Sums scaling coeffs down tree restoring state with coeffs only at leaves. Optional fence. Possible non-blocking comm.
840 void sum_down(bool fence = true) const {
842 verify();
843 MADNESS_CHECK_THROW(impl->get_tree_state()==redundant_after_merge, "sum_down requires a redundant_after_merge state");
844 const_cast<Function<T,NDIM>*>(this)->impl->sum_down(fence);
845 const_cast<Function<T,NDIM>*>(this)->impl->set_tree_state(reconstructed);
846
847 if (fence && VERIFY_TREE) verify_tree(); // Must be after in case nonstandard
848 }
849
850
851 /// Inplace autorefines the function. Optional fence. Possible non-blocking comm.
852 template <typename opT>
853 void refine_general(const opT& op, bool fence = true) const {
855 verify();
857 impl->refine(op, fence);
858 }
859
860
862 bool operator()(implT* impl, const Key<NDIM>& key, const nodeT& t) const {
863 return impl->autorefine_square_test(key, t);
864 }
865
866 template <typename Archive> void serialize (Archive& ar) {}
867 };
868
869 /// Inplace autorefines the function using same test as for squaring.
870
871 /// return this for chaining
872 const Function<T,NDIM>& refine(bool fence = true) const {
874 return *this;
875 }
876
877 /// Inplace broadens support in scaling function basis
879 bool fence = true) const {
880 verify();
881 reconstruct();
882 impl->broaden(bc.is_periodic(), fence);
883 }
884
885
886 /// Clears the function as if constructed uninitialized. Optional fence.
887
888 /// Any underlying data will not be freed until the next global fence.
889 void clear(bool fence = true) {
891 if (impl) {
892 World& world = impl->world;
893 impl.reset();
894 if (fence) world.gop.fence();
895 }
896 }
897
898 /// Process 0 prints a summary of all nodes in the tree (collective)
899 void print_tree(std::ostream& os = std::cout) const {
901 if (impl) impl->print_tree(os);
902 }
903
904 /// same as print_tree() but produces JSON-formatted string
905 /// @warning enclose the result in braces to make it a valid JSON object
906 void print_tree_json(std::ostream& os = std::cout) const {
908 if (impl) impl->print_tree_json(os);
909 }
910
911 /// Process 0 prints a graphviz-formatted output of all nodes in the tree (collective)
912 void print_tree_graphviz(std::ostream& os = std::cout) const {
914 os << "digraph G {" << std::endl;
915 if (impl) impl->print_tree_graphviz(os);
916 os << "}" << std::endl;
917 }
918
919 /// Print a summary of the load balancing info
920
921 /// This is serial and VERY expensive
922 void print_info() const {
924 if (impl) impl->print_info();
925 }
926
928 T (*f)(T);
930 void operator()(const Key<NDIM>& key, Tensor<T>& t) const {
931 UNARY_OPTIMIZED_ITERATOR(T, t, *_p0 = f(*_p0));
932 }
933 template <typename Archive> void serialize(Archive& ar) {}
934 };
935
936 /// Inplace unary operation on function values
937 void unaryop(T (*f)(T)) {
938 // Must fence here due to temporary object on stack
939 // stopping us returning before complete
941 }
942
943
944 /// Inplace unary operation on function values
945 template <typename opT>
946 void unaryop(const opT& op, bool fence=true) {
948 verify();
949 reconstruct();
950 impl->unary_op_value_inplace(op, fence);
951 }
952
953
954 /// Unary operation applied inplace to the coefficients
955 template <typename opT>
956 void unaryop_coeff(const opT& op,
957 bool fence = true) {
959 verify();
960 impl->unary_op_coeff_inplace(op, fence);
961 }
962
963
964 /// Unary operation applied inplace to the nodes
965 template <typename opT>
966 void unaryop_node(const opT& op,
967 bool fence = true) {
969 verify();
970 impl->unary_op_node_inplace(op, fence);
971 }
972
973
974
975
976 static void doconj(const Key<NDIM>, Tensor<T>& t) {
978 t.conj();
979 }
980
981 /// Inplace complex conjugate. No communication except for optional fence.
982
983 /// Returns this for chaining. Works in either basis.
989
990
991 /// Inplace, scale the function by a constant. No communication except for optional fence.
992
993 /// Works in either basis. Returns reference to this for chaining.
994 template <typename Q>
995 Function<T,NDIM>& scale(const Q q, bool fence=true) {
997 verify();
999 impl->scale_inplace(q,fence);
1000 return *this;
1001 }
1002
1003
1004 /// Inplace add scalar. No communication except for optional fence.
1007 verify();
1008 if (VERIFY_TREE) verify_tree();
1009 impl->add_scalar_inplace(t,fence);
1010 return *this;
1011 }
1012
1013
1014 /// Inplace, general bi-linear operation in wavelet basis. No communication except for optional fence.
1015
1016 /// If the functions are not in the wavelet basis an exception is thrown since this routine
1017 /// is intended to be fast and unexpected compression is assumed to be a performance bug.
1018 ///
1019 /// Returns this for chaining, can be in states compressed of redundant_after_merge.
1020 ///
1021 /// this and other may have different distributions and may even live in different worlds
1022 ///
1023 /// this <-- this*alpha + other*beta
1024 template <typename Q, typename R>
1026 const Function<Q,NDIM>& other, const R& beta, bool fence=true) {
1028 verify();
1029 other.verify();
1030
1031 // operation is done either in compressed or reconstructed state
1032 TreeState operating_state=this->get_impl()->get_tensor_type()==TT_FULL ? compressed : reconstructed;
1033
1034 TreeState thisstate=impl->get_tree_state();
1035 TreeState otherstate=other.get_impl()->get_tree_state();
1036
1038 MADNESS_CHECK_THROW(thisstate==compressed, "gaxpy: this must be compressed");
1039 MADNESS_CHECK_THROW(otherstate==compressed, "gaxpy: other must be compressed");
1040 impl->gaxpy_inplace(alpha, *other.get_impl(), beta, fence);
1041
1042 } else if (operating_state==reconstructed) {
1043 // this works both in reconstructed and redundant_after_merge states
1045 "gaxpy: this must be reconstructed or redundant_after_merge");
1046 MADNESS_CHECK_THROW(otherstate==reconstructed or otherstate==redundant_after_merge,
1047 "gaxpy: other must be reconstructed or redundant_after_merge");
1048
1049 impl->gaxpy_inplace_reconstructed(alpha,*other.get_impl(),beta,fence);
1050 } else {
1051 MADNESS_EXCEPTION("unknown tree state",1);
1052 }
1053 return *this;
1054 }
1055
1056
1057 /// Inplace addition of functions in the wavelet basis
1058
1059 /// Using operator notation forces a global fence after every operation.
1060 /// Functions don't need to be compressed, it's the caller's responsibility
1061 /// to choose an appropriate state with performance, usually compressed for 3d,
1062 /// reconstructed for 6d)
1063 template <typename Q>
1066
1067 // do this in reconstructed or compressed form
1069 this->change_tree_state(operating_state);
1071
1072 MADNESS_ASSERT(impl->get_tree_state() == other.get_impl()->get_tree_state());
1073 if (VERIFY_TREE) verify_tree();
1074 if (VERIFY_TREE) other.verify_tree();
1075 return gaxpy(T(1.0), other, Q(1.0), true);
1076 }
1077
1078
1079 /// Inplace subtraction of functions in the wavelet basis
1080
1081 /// Using operator notation forces a global fence after every operation
1082 template <typename Q>
1085 if (NDIM<=3) {
1086 compress();
1087 other.compress();
1088 } else {
1089 reconstruct();
1090 other.reconstruct();
1091 }
1092 MADNESS_ASSERT(impl->get_tree_state() == other.get_impl()->get_tree_state());
1093 if (VERIFY_TREE) verify_tree();
1094 if (VERIFY_TREE) other.verify_tree();
1095 return gaxpy(T(1.0), other, Q(-1.0), true);
1096 }
1097
1098
1099 /// Inplace scaling by a constant
1100
1101 /// Using operator notation forces a global fence after every operation
1102 template <typename Q>
1104 operator*=(const Q q) {
1106 scale(q,true);
1107 return *this;
1108 }
1109
1110
1111 /// Inplace squaring of function ... global comm only if not reconstructed
1112
1113 /// Returns *this for chaining.
1116 if (!is_reconstructed()) reconstruct();
1117 if (VERIFY_TREE) verify_tree();
1118 impl->square_inplace(fence);
1119 return *this;
1120 }
1121
1122 /// Returns *this for chaining.
1125 if (!is_reconstructed()) reconstruct();
1126 if (VERIFY_TREE) verify_tree();
1127 impl->abs_inplace(fence);
1128 return *this;
1129 }
1130
1131 /// Returns *this for chaining.
1134 if (!is_reconstructed()) reconstruct();
1135 if (VERIFY_TREE) verify_tree();
1136 impl->abs_square_inplace(fence);
1137 return *this;
1138 }
1139
1140 /// Returns local contribution to \c int(f(x),x) ... no communication
1141
1142 /// In the wavelet basis this is just the coefficient of the first scaling
1143 /// function which is a constant. In the scaling function basis we
1144 /// must add up contributions from each box.
1145 T trace_local() const {
1147 if (!impl) return 0.0;
1148 if (VERIFY_TREE) verify_tree();
1149 return impl->trace_local();
1150 }
1151
1152
1153 /// Returns global value of \c int(f(x),x) ... global comm required
1154 T trace() const {
1156 if (!impl) return 0.0;
1157 T sum = impl->trace_local();
1158 impl->world.gop.sum(sum);
1159 impl->world.gop.fence();
1160 return sum;
1161 }
1162
1163
1164 /// Returns local part of inner product ... throws if both not compressed
1165 template <typename R>
1166 TENSOR_RESULT_TYPE(T,R) inner_local(const Function<R,NDIM>& g) const {
1173 return impl->inner_local(*(g.get_impl()));
1174 }
1175
1176 /// Returns local part of dot product ... throws if both not compressed
1177 template <typename R>
1178 TENSOR_RESULT_TYPE(T,R) dot_local(const Function<R,NDIM>& g) const {
1181 MADNESS_ASSERT(g.is_compressed());
1183 if (VERIFY_TREE) g.verify_tree();
1184 return impl->dot_local(*(g.get_impl()));
1185 }
1186
1187
1188 /// With this being an on-demand function, fill the MRA tree according to different criteria
1189
1190 /// @param[in] g the function after which the MRA structure is modeled (any basis works)
1191 template<typename R>
1193 MADNESS_ASSERT(g.is_initialized());
1195
1196 // clear what we have
1197 impl->get_coeffs().clear();
1198
1199 //leaf_op<T,NDIM> gnode_is_leaf(g.get_impl().get());
1200 Leaf_op_other<T,NDIM> gnode_is_leaf(g.get_impl().get());
1201 impl->make_Vphi(gnode_is_leaf,fence);
1202 return *this;
1203
1204 }
1205
1206 /// With this being an on-demand function, fill the MRA tree according to different criteria
1207
1208 /// @param[in] op the convolution operator for screening
1209 template<typename opT>
1210 Function<T,NDIM>& fill_tree(const opT& op, bool fence=true) {
1212 // clear what we have
1213 impl->get_coeffs().clear();
1216 impl ->make_Vphi(leaf_op,fence);
1217 return *this;
1218 }
1219
1220 /// With this being an on-demand function, fill the MRA tree according to different criteria
1223 // clear what we have
1224 impl->get_coeffs().clear();
1226 impl->make_Vphi(leaf_op,fence);
1227 return *this;
1228 }
1229
1230 /// Special refinement on 6D boxes where the electrons come close (meet)
1231 /// @param[in] op the convolution operator for screening
1232 template<typename opT>
1233 Function<T,NDIM>& fill_cuspy_tree(const opT& op,const bool fence=true){
1235 // clear what we have
1236 impl->get_coeffs().clear();
1238
1240 impl ->make_Vphi(leaf_op,fence);
1241
1242 return *this;
1243 }
1244
1245 /// Special refinement on 6D boxes where the electrons come close (meet)
1248 // clear what we have
1249 impl->get_coeffs().clear();
1251
1253 impl ->make_Vphi(leaf_op,fence);
1254
1255 return *this;
1256 }
1257
1258 /// Special refinement on 6D boxes for the nuclear potentials (regularized with cusp, non-regularized with singularity)
1259 /// @param[in] op the convolution operator for screening
1260 template<typename opT>
1261 Function<T,NDIM>& fill_nuclear_cuspy_tree(const opT& op,const size_t particle,const bool fence=true){
1263 // clear what we have
1264 impl->get_coeffs().clear();
1266
1268 impl ->make_Vphi(leaf_op,fence);
1269
1270 return *this;
1271 }
1272
1273 /// Special refinement on 6D boxes for the nuclear potentials (regularized with cusp, non-regularized with singularity)
1276 // clear what we have
1277 impl->get_coeffs().clear();
1279
1281 impl ->make_Vphi(leaf_op,fence);
1282
1283 return *this;
1284 }
1285
1286 /// perform the hartree product of f*g, invoked by result
1287 template<size_t LDIM, size_t KDIM, typename opT>
1288 void do_hartree_product(const std::vector<std::shared_ptr<FunctionImpl<T,LDIM>>> left,
1289 const std::vector<std::shared_ptr<FunctionImpl<T,KDIM>>> right,
1290 const opT* op) {
1291
1292 // get the right leaf operator
1294 impl->hartree_product(left,right,leaf_op,true);
1295 impl->finalize_sum();
1296// this->truncate();
1297
1298 }
1299
1300 /// perform the hartree product of f*g, invoked by result
1301 template<size_t LDIM, size_t KDIM>
1302 void do_hartree_product(const std::vector<std::shared_ptr<FunctionImpl<T,LDIM>>> left,
1303 const std::vector<std::shared_ptr<FunctionImpl<T,KDIM>>> right) {
1304
1305// hartree_leaf_op<T,KDIM+LDIM> leaf_op(impl.get(),cdata.s0);
1307 impl->hartree_product(left,right,leaf_op,true);
1308 impl->finalize_sum();
1309// this->truncate();
1310
1311 }
1312
1313 /// Returns the inner product
1314
1315 /// Not efficient for computing multiple inner products
1316 /// @param[in] g Function, optionally on-demand
1317 template <typename R>
1320
1321 // fast return if possible
1322 if (not this->is_initialized()) return 0.0;
1323 if (not g.is_initialized()) return 0.0;
1324
1325 // if this and g are the same, use norm2()
1326 if constexpr (std::is_same_v<T,R>) {
1327 if (this->get_impl() == g.get_impl()) {
1328 TreeState state = this->get_impl()->get_tree_state();
1329 if (not(state == reconstructed or state == compressed))
1331 double norm = this->norm2();
1332 return norm * norm;
1333 }
1334 }
1335
1336 // do it case-by-case
1337 if constexpr (std::is_same_v<R,T>) {
1338 if (this->is_on_demand())
1339 return g.inner_on_demand(*this);
1340 if (g.is_on_demand())
1341 return this->inner_on_demand(g);
1342 }
1343
1345 if (VERIFY_TREE) g.verify_tree();
1346
1347 // compute in compressed form if compression is fast, otherwise in redundant form
1349
1351 g.change_tree_state(operating_state,false);
1352 impl->world.gop.fence();
1353
1354 TENSOR_RESULT_TYPE(T,R) local = impl->inner_local(*g.get_impl());
1355 impl->world.gop.sum(local);
1356 impl->world.gop.fence();
1357
1358 // restore state -- no need for this
1359 // change_tree_state(state,false);
1360 // g.change_tree_state(gstate,false);
1361 // impl->world.gop.fence();
1362
1363 return local;
1364 }
1365
1366 /// Return the local part of inner product with external function ... no communication.
1367 /// If you are going to be doing a bunch of inner_ext calls, set
1368 /// keep_redundant to true and then manually undo_redundant when you
1369 /// are finished.
1370 /// @param[in] f Pointer to function of type T that take coordT arguments. This is the externally provided function
1371 /// @param[in] leaf_refine boolean switch to turn on/off refinement past leaf nodes
1372 /// @param[in] keep_redundant boolean switch to turn on/off undo_redundant
1373 /// @return Returns local part of the inner product, i.e. over the domain of all function nodes on this compute node.
1374 T inner_ext_local(const std::shared_ptr< FunctionFunctorInterface<T,NDIM> > f, const bool leaf_refine=true, const bool keep_redundant=false) const {
1377 T local = impl->inner_ext_local(f, leaf_refine);
1378 if (not keep_redundant) change_tree_state(reconstructed);
1379 return local;
1380 }
1381
1382 /// Return the inner product with external function ... requires communication.
1383 /// If you are going to be doing a bunch of inner_ext calls, set
1384 /// keep_redundant to true and then manually undo_redundant when you
1385 /// are finished.
1386 /// @param[in] f Reference to FunctionFunctorInterface. This is the externally provided function
1387 /// @param[in] leaf_refine boolean switch to turn on/off refinement past leaf nodes
1388 /// @param[in] keep_redundant boolean switch to turn on/off undo_redundant
1389 /// @return Returns the inner product
1390 T inner_ext(const std::shared_ptr< FunctionFunctorInterface<T,NDIM> > f, const bool leaf_refine=true, const bool keep_redundant=false) const {
1393 T local = impl->inner_ext_local(f, leaf_refine);
1394 impl->world.gop.sum(local);
1395 impl->world.gop.fence();
1396 if (not keep_redundant) change_tree_state(reconstructed);
1397 return local;
1398 }
1399
1400 /// Return the inner product with external function ... requires communication.
1401 /// If you are going to be doing a bunch of inner_ext calls, set
1402 /// keep_redundant to true and then manually undo_redundant when you
1403 /// are finished.
1404 /// @param[in] f Reference to FunctionFunctorInterface. This is the externally provided function
1405 /// @param[in] leaf_refine boolean switch to turn on/off refinement past leaf nodes
1406 /// @return Returns the inner product
1408 const bool leaf_refine=true) const {
1410 reconstruct();
1411 T local = impl->inner_adaptive_local(f, leaf_refine);
1412 impl->world.gop.sum(local);
1413 impl->world.gop.fence();
1414 return local;
1415 }
1416
1417 /// Return the local part of gaxpy with external function, this*alpha + f*beta ... no communication.
1418 /// @param[in] alpha prefactor for this Function
1419 /// @param[in] f Pointer to function of type T that take coordT arguments. This is the externally provided function
1420 /// @param[in] beta prefactor for f
1421 template <typename L>
1422 void gaxpy_ext(const Function<L,NDIM>& left, T (*f)(const coordT&), T alpha, T beta, double tol, bool fence=true) const {
1424 if (!left.is_reconstructed()) left.reconstruct();
1425 impl->gaxpy_ext(left.get_impl().get(), f, alpha, beta, tol, fence);
1426 }
1427
1428 /// Returns the inner product for one on-demand function
1429
1430 /// It does work, but it might not give you the precision you expect.
1431 /// The assumption is that the function g returns proper sum
1432 /// coefficients on the MRA tree of this. This might not be the case if
1433 /// g is constructed with an implicit multiplication, e.g.
1434 /// result = <this|g>, with g = 1/r12 | gg>
1435 /// @param[in] g on-demand function
1436 template<typename R>
1437 TENSOR_RESULT_TYPE(T, R) inner_on_demand(const Function<R, NDIM>& g) const {
1438 MADNESS_ASSERT(g.is_on_demand() and (not this->is_on_demand()));
1439
1440 constexpr std::size_t LDIM=std::max(NDIM/2,std::size_t(1));
1441 auto func=dynamic_cast<CompositeFunctorInterface<T,NDIM,LDIM>* >(g.get_impl()->get_functor().get());
1443 func->make_redundant(true);
1444 func->replicate_low_dim_functions(true);
1445 this->reconstruct(); // if this == &g we don't need g to be redundant
1446
1448
1449 TENSOR_RESULT_TYPE(T, R) local = impl->inner_local_on_demand(*g.get_impl());
1450 impl->world.gop.sum(local);
1451 impl->world.gop.fence();
1452
1453 return local;
1454 }
1455
1456 /// project this on the low-dim function g: h(x) = <f(x,y) | g(y)>
1457
1458 /// @param[in] g low-dim function
1459 /// @param[in] dim over which dimensions to be integrated: 0..LDIM-1 or LDIM..NDIM-1
1460 /// @return new function of dimension NDIM-LDIM
1461 template <typename R, size_t LDIM>
1463 if (NDIM<=LDIM) MADNESS_EXCEPTION("confused dimensions in project_out?",1);
1464 MADNESS_CHECK_THROW(dim==0 or dim==1,"dim must be 0 or 1 in project_out");
1465 verify();
1466 typedef TENSOR_RESULT_TYPE(T,R) resultT;
1467 static const size_t KDIM=NDIM-LDIM;
1468
1470 .k(g.k()).thresh(g.thresh());
1471 Function<resultT,KDIM> result=factory; // no empty() here!
1472
1474 g.change_tree_state(redundant,false);
1475 world().gop.fence();
1476 this->get_impl()->project_out(result.get_impl().get(),g.get_impl().get(),dim,true);
1477// result.get_impl()->project_out2(this->get_impl().get(),gimpl,dim);
1478 result.world().gop.fence();
1479 g.change_tree_state(reconstructed,false);
1480 result.get_impl()->trickle_down(false);
1481 result.get_impl()->set_tree_state(reconstructed);
1482 result.world().gop.fence();
1483 return result;
1484 }
1485
1486 Function<T,NDIM/2> dirac_convolution(const bool fence=true) const {
1487 constexpr std::size_t LDIM=NDIM/2;
1488 MADNESS_CHECK_THROW(NDIM==2*LDIM,"NDIM must be even");
1489// // this will be the result function
1491 Function<T,LDIM> f = factory;
1492 if(!is_reconstructed()) this->reconstruct();
1493 this->get_impl()->do_dirac_convolution(f.get_impl().get(),fence);
1494 return f;
1495 }
1496
1497 /// Replaces this function with one loaded from an archive using the default processor map
1498
1499 /// Archive can be sequential or parallel.
1500 ///
1501 /// The & operator for serializing will only work with parallel archives.
1502 template <typename Archive>
1503 void load(World& world, Archive& ar) {
1505 // Type checking since we are probably circumventing the archive's own type checking
1506 long magic = 0l, id = 0l, ndim = 0l, k = 0l;
1507 ar & magic & id & ndim & k;
1508 MADNESS_ASSERT(magic == 7776768); // Mellow Mushroom Pizza tel.# in Knoxville
1510 MADNESS_ASSERT(ndim == NDIM);
1511
1512 impl.reset(new implT(FunctionFactory<T,NDIM>(world).k(k).empty()));
1513
1514 impl->load(ar);
1515 }
1516
1517
1518 /// Stores the function to an archive
1519
1520 /// Archive can be sequential or parallel.
1521 ///
1522 /// The & operator for serializing will only work with parallel archives.
1523 template <typename Archive>
1524 void store(Archive& ar) const {
1526 verify();
1527 // For type checking, etc.
1528 ar & long(7776768) & long(TensorTypeData<T>::id) & long(NDIM) & long(k());
1529
1530 impl->store(ar);
1531 }
1532
1533 /// change the tensor type of the coefficients in the FunctionNode
1534
1535 /// @param[in] targs target tensor arguments (threshold and full/low rank)
1536 void change_tensor_type(const TensorArgs& targs, bool fence=true) {
1537 if (not impl) return;
1538 impl->change_tensor_type1(targs,fence);
1539 }
1540
1541
1542 /// This is replaced with left*right ... private
1543 template <typename Q, typename opT>
1545 const opT& op, bool fence) {
1547 func.verify();
1548 MADNESS_ASSERT(func.is_reconstructed());
1549 if (VERIFY_TREE) func.verify_tree();
1550 impl.reset(new implT(*func.get_impl(), func.get_pmap(), false));
1551 impl->unaryXX(func.get_impl().get(), op, fence);
1552 return *this;
1553 }
1554
1555 /// Returns vector of FunctionImpl pointers corresponding to vector of functions
1556 template <typename Q, std::size_t D>
1557 static std::vector< std::shared_ptr< FunctionImpl<Q,D> > > vimpl(const std::vector< Function<Q,D> >& v) {
1559 std::vector< std::shared_ptr< FunctionImpl<Q,D> > > r(v.size());
1560 for (unsigned int i=0; i<v.size(); ++i) r[i] = v[i].get_impl();
1561 return r;
1562 }
1563
1564 /// This is replaced with op(vector of functions) ... private
1565 template <typename opT>
1566 Function<T,NDIM>& multiop_values(const opT& op, const std::vector< Function<T,NDIM> >& vf) {
1567 std::vector<implT*> v(vf.size(),NULL);
1568 for (unsigned int i=0; i<v.size(); ++i) {
1569 if (vf[i].is_initialized()) v[i] = vf[i].get_impl().get();
1570 }
1571 impl->multiop_values(op, v);
1572 world().gop.fence();
1573 if (VERIFY_TREE) verify_tree();
1574
1575 return *this;
1576 }
1577
1578 /// apply op on the input vector yielding an output vector of functions
1579
1580 /// (*this) is just a dummy Function to be able to call internal methods in FuncImpl
1581 /// @param[in] op the operator working on vin
1582 /// @param[in] vin vector of input Functions
1583 /// @param[out] vout vector of output Functions vout = op(vin)
1584 template <typename opT>
1586 const std::vector< Function<T,NDIM> >& vin,
1587 std::vector< Function<T,NDIM> >& vout,
1588 const bool fence=true) {
1589 std::vector<implT*> vimplin(vin.size(),NULL);
1590 for (unsigned int i=0; i<vin.size(); ++i) {
1591 if (vin[i].is_initialized()) vimplin[i] = vin[i].get_impl().get();
1592 }
1593 std::vector<implT*> vimplout(vout.size(),NULL);
1594 for (unsigned int i=0; i<vout.size(); ++i) {
1595 if (vout[i].is_initialized()) vimplout[i] = vout[i].get_impl().get();
1596 }
1597
1598 impl->multi_to_multi_op_values(op, vimplin, vimplout, fence);
1599 if (VERIFY_TREE) verify_tree();
1600
1601 }
1602
1603
1604 /// Multiplication of function * vector of functions using recursive algorithm of mulxx
1605 template <typename L, typename R>
1606 void vmulXX(const Function<L,NDIM>& left,
1607 const std::vector< Function<R,NDIM> >& right,
1608 std::vector< Function<T,NDIM> >& result,
1609 double tol,
1610 bool fence) {
1612
1613 std::vector<FunctionImpl<T,NDIM>*> vresult(right.size());
1614 std::vector<const FunctionImpl<R,NDIM>*> vright(right.size());
1615 for (unsigned int i=0; i<right.size(); ++i) {
1616 result[i].set_impl(left,false);
1617 vresult[i] = result[i].impl.get();
1618 vright[i] = right[i].get_impl().get();
1619 }
1620
1621 left.world().gop.fence(); // Is this still essential? Yes.
1622 vresult[0]->mulXXvec(left.get_impl().get(), vright, vresult, tol, fence);
1623 }
1624
1625 /// Same as \c operator* but with optional fence and no automatic reconstruction
1626
1627 /// f or g are on-demand functions
1628 template<typename L, typename R>
1629 void mul_on_demand(const Function<L,NDIM>& f, const Function<R,NDIM>& g, bool fence=true) {
1630 const FunctionImpl<L,NDIM>* fimpl=f.get_impl().get();
1631 const FunctionImpl<R,NDIM>* gimpl=g.get_impl().get();
1632 if (fimpl->is_on_demand() and gimpl->is_on_demand()) {
1633 MADNESS_EXCEPTION("can't multiply two on-demand functions",1);
1634 }
1635
1636 if (fimpl->is_on_demand()) {
1637 leaf_op<R,NDIM> leaf_op1(gimpl);
1638 impl->multiply(leaf_op1,gimpl,fimpl,fence);
1639 } else {
1640 leaf_op<L,NDIM> leaf_op1(fimpl);
1641 impl->multiply(leaf_op1,fimpl,gimpl,fence);
1642 }
1643 }
1644
1645 /// sparse transformation of a vector of functions ... private
1646 template <typename R, typename Q>
1647 void vtransform(const std::vector< Function<R,NDIM> >& v,
1648 const Tensor<Q>& c,
1649 std::vector< Function<T,NDIM> >& vresult,
1650 double tol,
1651 bool fence=true) {
1653 vresult[0].impl->vtransform(vimpl(v), c, vimpl(vresult), tol, fence);
1654 }
1655
1656 /// This is replaced with alpha*left + beta*right ... private
1657 template <typename L, typename R>
1659 T beta, const Function<R,NDIM>& right, bool fence) {
1661 left.verify();
1662 right.verify();
1663 MADNESS_ASSERT(left.is_compressed() && right.is_compressed());
1664 if (VERIFY_TREE) left.verify_tree();
1665 if (VERIFY_TREE) right.verify_tree();
1666 impl.reset(new implT(*left.get_impl(), left.get_pmap(), false));
1667 impl->gaxpy(alpha,*left.get_impl(),beta,*right.get_impl(),fence);
1668 return *this;
1669 }
1670
1671 /// This is replaced with mapdim(f) ... private
1672 Function<T,NDIM>& mapdim(const Function<T,NDIM>& f, const std::vector<long>& map, bool fence) {
1674 f.verify();
1675 if (VERIFY_TREE) f.verify_tree();
1676 for (std::size_t i=0; i<NDIM; ++i) MADNESS_ASSERT(map[i]>=0 && static_cast<std::size_t>(map[i])<NDIM);
1677 impl.reset(new implT(*f.impl, f.get_pmap(), false));
1678 impl->mapdim(*f.impl,map,fence);
1679 return *this;
1680 }
1681
1682 /// This is replaced with mirror(f) ... private
1683
1684 /// similar to mapdim, but maps from x to -x, y to -y, and so on
1685 /// Example: mirror a 3d function on the xy plane: mirror={1,1,-1}
1686 /// @param[in] mirror array of -1 and 1, corresponding to mirror or not
1687 Function<T,NDIM>& mirror(const Function<T,NDIM>& f, const std::vector<long>& mirrormap, bool fence) {
1689 f.verify();
1690 if (VERIFY_TREE) f.verify_tree();
1691 for (std::size_t i=0; i<NDIM; ++i) MADNESS_ASSERT((mirrormap[i]==1) or (mirrormap[i]==-1));
1692 impl.reset(new implT(*f.impl, f.get_pmap(), false));
1693 impl->mirror(*f.impl,mirrormap,fence);
1694 return *this;
1695 }
1696
1697 /// This is replaced with mirror(map(f)) ... private
1698
1699 /// first map then mirror!
1700 /// mirror is similar to mapdim, but maps from x to -x, y to -y, and so on
1701 /// Example: mirror a 3d function on the xy plane: mirror={1,1,-1}
1702 /// Example: c4 rotation of a 3d function around the z axis:
1703 /// x->y, y->-x, z->z: map(1,0,2); mirror(-1,1,1)
1704 /// @param[in] map array holding dimensions
1705 /// @param[in] mirror array of -1 and 1, corresponding to mirror or not
1707 const std::vector<long>& map, const std::vector<long>& mirror,
1708 bool fence) {
1710 f.verify();
1711 if (VERIFY_TREE) f.verify_tree();
1712 for (std::size_t i=0; i<mirror.size(); ++i) MADNESS_ASSERT((mirror[i]==1) or (mirror[i]==-1));
1713 for (std::size_t i=0; i<map.size(); ++i) MADNESS_ASSERT(map[i]>=0 && static_cast<std::size_t>(map[i])<NDIM);
1714
1715 impl.reset(new implT(*f.impl, f.get_pmap(), false));
1716 impl->map_and_mirror(*f.impl,map,mirror,fence);
1717 return *this;
1718 }
1719
1720
1721 /// check symmetry of a function by computing the 2nd derivative
1722 double check_symmetry() const {
1723
1725 if (VERIFY_TREE) verify_tree();
1726 double local = impl->check_symmetry_local();
1727 impl->world.gop.sum(local);
1728 impl->world.gop.fence();
1729 double asy=sqrt(local);
1730 if (this->world().rank()==0) print("asymmetry wrt particle",asy);
1732 return asy;
1733 }
1734
1735 /// reduce the rank of the coefficient tensors
1736 Function<T,NDIM>& reduce_rank(const double thresh=0.0, const bool fence=true) {
1737 verify();
1738 double thresh1= (thresh==0.0) ? impl->get_tensor_args().thresh : thresh;
1739 impl->reduce_rank(thresh1,fence);
1740 return *this;
1741 }
1742
1743 /// remove all nodes with level higher than n
1744 Function<T,NDIM>& chop_at_level(const int n, const bool fence=true) {
1745 verify();
1747 impl->chop_at_level(n,true);
1749 return *this;
1750 }
1751 };
1752
1753// template <typename T, typename opT, std::size_t NDIM>
1754 template <typename T, typename opT, std::size_t NDIM>
1755 Function<T,NDIM> multiop_values(const opT& op, const std::vector< Function<T,NDIM> >& vf) {
1757 r.set_impl(vf[0], false);
1758 r.multiop_values(op, vf);
1759 return r;
1760 }
1761
1762 /// Returns new function equal to alpha*f(x) with optional fence
1763 template <typename Q, typename T, std::size_t NDIM>
1765 mul(const Q alpha, const Function<T,NDIM>& f, bool fence=true) {
1767 f.verify();
1768 if (VERIFY_TREE) f.verify_tree();
1770 result.set_impl(f, false);
1771 result.get_impl()->scale_oop(alpha,*f.get_impl(),fence);
1772 return result;
1773 }
1774
1775
1776 /// Returns new function equal to f(x)*alpha with optional fence
1777 template <typename Q, typename T, std::size_t NDIM>
1779 mul(const Function<T,NDIM>& f, const Q alpha, bool fence=true) {
1781 return mul(alpha,f,fence);
1782 }
1783
1784
1785 /// Returns new function equal to f(x)*alpha
1786
1787 /// Using operator notation forces a global fence after each operation
1788 template <typename Q, typename T, std::size_t NDIM>
1791 return mul(alpha, f, true);
1792 }
1793
1794 /// Returns new function equal to alpha*f(x)
1795
1796 /// Using operator notation forces a global fence after each operation
1797 template <typename Q, typename T, std::size_t NDIM>
1800 return mul(alpha, f, true);
1801 }
1802
1803 /// Sparse multiplication --- left and right must be reconstructed and if tol!=0 have tree of norms already created
1804 template <typename L, typename R,std::size_t NDIM>
1806 mul_sparse(const Function<L,NDIM>& left, const Function<R,NDIM>& right, double tol, bool fence=true) {
1808 left.verify();
1809 right.verify();
1811 if (VERIFY_TREE) left.verify_tree();
1812 if (VERIFY_TREE) right.verify_tree();
1813
1815 result.set_impl(left, false);
1816 result.get_impl()->mulXX(left.get_impl().get(), right.get_impl().get(), tol, fence);
1817 return result;
1818 }
1819
1820 /// Same as \c operator* but with optional fence and no automatic reconstruction
1821 template <typename L, typename R,std::size_t NDIM>
1823 mul(const Function<L,NDIM>& left, const Function<R,NDIM>& right, bool fence=true) {
1824 return mul_sparse(left,right,0.0,fence);
1825 }
1826
1827 /// Generate new function = op(left,right) where op acts on the function values
1828 template <typename L, typename R, typename opT, std::size_t NDIM>
1830 binary_op(const Function<L,NDIM>& left, const Function<R,NDIM>& right, const opT& op, bool fence=true) {
1832 if (!left.is_reconstructed()) left.reconstruct();
1833 if (!right.is_reconstructed()) right.reconstruct();
1834
1836 result.set_impl(left, false);
1837 result.get_impl()->binaryXX(left.get_impl().get(), right.get_impl().get(), op, fence);
1838 return result;
1839 }
1840
1841 /// Out of place application of unary operation to function values with optional fence
1842 template <typename Q, typename opT, std::size_t NDIM>
1843 Function<typename opT::resultT, NDIM>
1844 unary_op(const Function<Q,NDIM>& func, const opT& op, bool fence=true) {
1845 if (!func.is_reconstructed()) func.reconstruct();
1848 result.set_impl(func, false);
1849 result.get_impl()->unaryXXvalues(func.get_impl().get(), op, fence);
1850 return result;
1851 }
1852
1853
1854 /// Out of place application of unary operation to scaling function coefficients with optional fence
1855 template <typename Q, typename opT, std::size_t NDIM>
1856 Function<typename opT::resultT, NDIM>
1857 unary_op_coeffs(const Function<Q,NDIM>& func, const opT& op, bool fence=true) {
1858 if (!func.is_reconstructed()) func.reconstruct();
1860 return result.unary_op_coeffs(func,op,fence);
1861 }
1862
1863 /// Use the vmra/mul(...) interface instead
1864
1865 /// This so that we don't have to have friend functions in a different header.
1866 ///
1867 /// If using sparsity (tol != 0) you must have created the tree of norms
1868 /// already for both left and right.
1869 template <typename L, typename R, std::size_t D>
1870 std::vector< Function<TENSOR_RESULT_TYPE(L,R),D> >
1871 vmulXX(const Function<L,D>& left, const std::vector< Function<R,D> >& vright, double tol, bool fence=true) {
1872 if (vright.size() == 0) return std::vector< Function<TENSOR_RESULT_TYPE(L,R),D> >();
1873 std::vector< Function<TENSOR_RESULT_TYPE(L,R),D> > vresult(vright.size());
1874 vresult[0].vmulXX(left, vright, vresult, tol, fence);
1875 return vresult;
1876 }
1877
1878 /// Multiplies two functions with the new result being of type TensorResultType<L,R>
1879
1880 /// Using operator notation forces a global fence after each operation but also
1881 /// enables us to automatically reconstruct the input functions as required.
1882 template <typename L, typename R, std::size_t NDIM>
1884 operator*(const Function<L,NDIM>& left, const Function<R,NDIM>& right) {
1885 if (!left.is_reconstructed()) left.reconstruct();
1886 if (!right.is_reconstructed()) right.reconstruct();
1887 MADNESS_ASSERT(not (left.is_on_demand() or right.is_on_demand()));
1888 return mul(left,right,true);
1889 }
1890
1891 /// Performs a Hartree/outer product on the two given low-dimensional function vectors
1892
1893 /// @return result(x,y) = \sum_i f_i(x) g_i(y)
1894 template<typename T, std::size_t KDIM, std::size_t LDIM>
1895 Function<T,KDIM+LDIM>
1896 hartree_product(const std::vector<Function<T,KDIM>>& left, const std::vector<Function<T,LDIM>>& right) {
1897
1898 MADNESS_CHECK_THROW(left.size()==right.size(), "hartree_product: left and right must have same size");
1899 if (left.size()==0) return Function<T,KDIM+LDIM>();
1900
1902
1904 .k(left.front().k()).thresh(thresh);
1905 Function<T,KDIM+LDIM> result=factory.empty();
1906
1907 // some prep work
1910 std::vector<std::shared_ptr<FunctionImpl<T,KDIM>>> vleft=get_impl(left);
1911 std::vector<std::shared_ptr<FunctionImpl<T,LDIM>>> vright=get_impl(right);
1912
1913 result.do_hartree_product(vleft,vright);
1914
1915 return result;
1916
1917 }
1918
1919 /// Performs a Hartree product on the two given low-dimensional functions
1920 template<typename T, std::size_t KDIM, std::size_t LDIM>
1921 Function<T,KDIM+LDIM>
1923 typedef std::vector<Function<T,KDIM>> vector;
1924 return hartree_product(vector({left2}),vector({right2}));
1925 }
1926
1927 /// Performs a Hartree product on the two given low-dimensional functions
1928 template<typename T, std::size_t KDIM, std::size_t LDIM, typename opT>
1929 Function<T,KDIM+LDIM>
1931 const opT& op) {
1932
1933 // we need both sum and difference coeffs for error estimation
1934 Function<T,KDIM>& left = const_cast< Function<T,KDIM>& >(left2);
1935 Function<T,LDIM>& right = const_cast< Function<T,LDIM>& >(right2);
1936
1938
1940 .k(left.k()).thresh(thresh);
1941 Function<T,KDIM+LDIM> result=factory.empty();
1942
1943 if (result.world().rank()==0) {
1944 print("incomplete FunctionFactory in Function::hartree_product");
1945 print("thresh: ", thresh);
1946 }
1947 bool same=(left2.get_impl()==right2.get_impl());
1948
1949 // some prep work
1950 left.make_nonstandard(true, true);
1951 right.make_nonstandard(true, true);
1952
1953 std::vector<std::shared_ptr<FunctionImpl<T,KDIM>>> vleft;
1954 std::vector<std::shared_ptr<FunctionImpl<T,LDIM>>> vright;
1955 vleft.push_back(left.get_impl());
1956 vright.push_back(right.get_impl());
1957 result.do_hartree_product(vleft,right,&op);
1958
1959 left.standard(false);
1960 if (not same) right.standard(false);
1961 left2.world().gop.fence();
1962
1963 return result;
1964 }
1965
1966 /// adds beta*right only left: alpha*left + beta*right optional fence and no automatic compression
1967
1968 /// left and right might live in different worlds, the accumulation is non-blocking
1969 template <typename L, typename R,std::size_t NDIM>
1970 void
1972 TENSOR_RESULT_TYPE(L,R) beta, const Function<R,NDIM>& right, bool fence=true) {
1975 left.gaxpy(alpha, right, beta, fence);
1976 }
1977
1978 /// Returns new function alpha*left + beta*right optional fence and no automatic compression
1979 template <typename L, typename R,std::size_t NDIM>
1982 TENSOR_RESULT_TYPE(L,R) beta, const Function<R,NDIM>& right, bool fence=true) {
1985 return result.gaxpy_oop(alpha, left, beta, right, fence);
1986 }
1987
1988 /// Same as \c operator+ but with optional fence and no automatic compression
1989 template <typename L, typename R,std::size_t NDIM>
1991 add(const Function<L,NDIM>& left, const Function<R,NDIM>& right, bool fence=true) {
1992 return gaxpy_oop(TENSOR_RESULT_TYPE(L,R)(1.0), left,
1993 TENSOR_RESULT_TYPE(L,R)(1.0), right, fence);
1994 }
1995
1996
1997 /// Returns new function alpha*left + beta*right optional fence, having both addends reconstructed
1998 template<typename T, std::size_t NDIM>
2000 const double beta, const Function<T,NDIM>& right, const bool fence=true) {
2001 Function<T,NDIM> result;
2002 result.set_impl(right,false);
2003
2006 result.get_impl()->gaxpy_oop_reconstructed(alpha,*left.get_impl(),beta,*right.get_impl(),fence);
2007 return result;
2008
2009 }
2010
2011 /// Adds two functions with the new result being of type TensorResultType<L,R>
2012
2013 /// Using operator notation forces a global fence after each operation
2014 template <typename L, typename R, std::size_t NDIM>
2016 operator+(const Function<L,NDIM>& left, const Function<R,NDIM>& right) {
2017 if (VERIFY_TREE) left.verify_tree();
2018 if (VERIFY_TREE) right.verify_tree();
2019
2020 // no compression for high-dimensional functions
2021 if (NDIM==6) {
2022 left.reconstruct();
2023 right.reconstruct();
2024 return gaxpy_oop_reconstructed(1.0,left,1.0,right,true);
2025 } else {
2026 if (!left.is_compressed()) left.compress();
2027 if (!right.is_compressed()) right.compress();
2028 return add(left,right,true);
2029 }
2030 }
2031
2032 /// Same as \c operator- but with optional fence and no automatic compression
2033 template <typename L, typename R,std::size_t NDIM>
2035 sub(const Function<L,NDIM>& left, const Function<R,NDIM>& right, bool fence=true) {
2036 return gaxpy_oop(TENSOR_RESULT_TYPE(L,R)(1.0), left,
2037 TENSOR_RESULT_TYPE(L,R)(-1.0), right, fence);
2038 }
2039
2040
2041 /// Subtracts two functions with the new result being of type TensorResultType<L,R>
2042
2043 /// Using operator notation forces a global fence after each operation
2044 template <typename L, typename R, std::size_t NDIM>
2046 operator-(const Function<L,NDIM>& left, const Function<R,NDIM>& right) {
2048 // no compression for high-dimensional functions
2049 if (NDIM==6) {
2050 left.reconstruct();
2051 right.reconstruct();
2052 return gaxpy_oop_reconstructed(1.0,left,-1.0,right,true);
2053 } else {
2054 if (!left.is_compressed()) left.compress();
2055 if (!right.is_compressed()) right.compress();
2056 return sub(left,right,true);
2057 }
2058 }
2059
2060 /// Create a new copy of the function with different distribution and optional fence
2061
2062 /// Works in either basis. Different distributions imply
2063 /// asynchronous communication and the optional fence is
2064 /// collective.
2065 template <typename T, std::size_t NDIM>
2067 const std::shared_ptr< WorldDCPmapInterface< Key<NDIM> > >& pmap,
2068 bool fence = true) {
2070 f.verify();
2071 Function<T,NDIM> result;
2073 result.set_impl(std::shared_ptr<implT>(new implT(*f.get_impl(), pmap, false)));
2074 result.get_impl()->copy_coeffs(*f.get_impl(), fence);
2075 if (VERIFY_TREE) result.verify_tree();
2076 return result;
2077 }
2078
2079 /// Create a new copy of the function with the same distribution and optional fence
2080 template <typename T, std::size_t NDIM>
2083 return copy(f, f.get_pmap(), fence);
2084 }
2085
2086 /// Type conversion implies a deep copy. No communication except for optional fence.
2087
2088 /// Works in either basis but any loss of precision may result in different errors
2089 /// in applied in a different basis.
2090 ///
2091 /// The new function is formed with the options from the default constructor.
2092 ///
2093 /// There is no automatic type conversion since this is generally a rather dangerous
2094 /// thing and because there would be no way to make the fence optional.
2095 template <typename T, typename Q, std::size_t NDIM>
2098 f.verify();
2099 Function<Q,NDIM> result;
2100 result.set_impl(f, false);
2101 result.get_impl()->copy_coeffs(*f.get_impl(), fence);
2102 return result;
2103 }
2104
2105
2106 /// Return the complex conjugate of the input function with the same distribution and optional fence
2107
2108 /// !!! The fence is actually not optional in the current implementation !!!
2109 template <typename T, std::size_t NDIM>
2112 Function<T,NDIM> result = copy(f,true);
2113 return result.conj(fence);
2114 }
2115
2116 /// Apply operator on a hartree product of two low-dimensional functions
2117
2118 /// Supposed to be something like result= G( f(1)*f(2))
2119 /// the hartree product is never constructed explicitly, but its coeffs are
2120 /// constructed on the fly and processed immediately.
2121 /// @param[in] op the operator
2122 /// @param[in] f1 function of particle 1
2123 /// @param[in] f2 function of particle 2
2124 /// @param[in] fence if we shall fence
2125 /// @return a function of dimension NDIM=LDIM+LDIM
2126 template <typename opT, typename T, std::size_t LDIM>
2127 Function<TENSOR_RESULT_TYPE(typename opT::opT,T), LDIM+LDIM>
2128 apply(const opT& op, const std::vector<Function<T,LDIM>>& f1, const std::vector<Function<T,LDIM>>& f2, bool fence=true) {
2129
2130 World& world=f1.front().world();
2131
2132 typedef TENSOR_RESULT_TYPE(T,typename opT::opT) resultT;
2133 typedef std::vector<Function<T,LDIM>> vecfuncL;
2134
2135 vecfuncL& ff1 = const_cast< vecfuncL& >(f1);
2136 vecfuncL& ff2 = const_cast< vecfuncL& >(f2);
2137
2138 bool same=(ff1[0].get_impl()==ff2[0].get_impl());
2139
2140 reconstruct(world,f1,false);
2141 reconstruct(world,f2,false);
2142 world.gop.fence();
2143 // keep the leaves! They are assumed to be there later
2144 // even for modified op we need NS form for the hartree_leaf_op
2145 for (auto& f : f1) f.make_nonstandard(true,false);
2146 for (auto& f : f2) f.make_nonstandard(true,false);
2147 world.gop.fence();
2148
2149
2152 Function<resultT,LDIM+LDIM> result=factory.empty().fence();
2153
2154 result.get_impl()->reset_timer();
2155 op.reset_timer();
2156
2157 // will fence here
2158 for (size_t i=0; i<f1.size(); ++i)
2159 result.get_impl()->recursive_apply(op, f1[i].get_impl().get(),f2[i].get_impl().get(),false);
2160 world.gop.fence();
2161
2162 if (op.print_timings) {
2163 result.get_impl()->print_timer();
2164 op.print_timer();
2165 }
2166
2167 result.get_impl()->finalize_apply(); // need fence before reconstruct
2168
2169 if (op.modified()) {
2170 result.get_impl()->trickle_down(true);
2171 } else {
2172 result.get_impl()->reconstruct(true);
2173 }
2174 standard(world,ff1,false);
2175 if (not same) standard(world,ff2,false);
2176
2177 return result;
2178 }
2179
2180
2181 /// Apply operator ONLY in non-standard form - required other steps missing !!
2182 template <typename opT, typename R, std::size_t NDIM>
2183 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM>
2184 apply_only(const opT& op, const Function<R,NDIM>& f, bool fence=true) {
2185 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> result;
2186
2187 constexpr std::size_t OPDIM=opT::opdim;
2188 constexpr bool low_dim=(OPDIM*2==NDIM); // apply on some dimensions only
2189
2190 // specialized version for 3D
2191 if (NDIM <= 3 and (not low_dim)) {
2192 result.set_impl(f, false);
2193 result.get_impl()->apply(op, *f.get_impl(), fence);
2194
2195 } else { // general version for higher dimension
2196 //bool print_timings=false;
2197 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> r1;
2198
2199 result.set_impl(f, false);
2200 r1.set_impl(f, false);
2201
2202 result.get_impl()->reset_timer();
2203 op.reset_timer();
2204
2205 result.get_impl()->apply_source_driven(op, *f.get_impl(), fence);
2206
2207 // recursive_apply is about 20% faster than apply_source_driven
2208 //result.get_impl()->recursive_apply(op, f.get_impl().get(),
2209 // r1.get_impl().get(),true); // will fence here
2210
2211 }
2212
2213 return result;
2214 }
2215
2216 /// Apply operator in non-standard form
2217
2218 /// Returns a new function with the same distribution
2219 ///
2220 /// !!! For the moment does NOT respect fence option ... always fences
2221 /// if the operator acts on one particle only the result will be sorted as
2222 /// g.particle=1: g(f) = \int g(x,x') f(x',y) dx' = result(x,y)
2223 /// g.particle=2: g(f) = \int g(y,y') f(x,y') dy' = result(x,y)
2224 /// for the second case it will notably *not* be as it is implemented in the partial inner product!
2225 /// g.particle=2 g(f) = result(x,y)
2226 /// inner(g(y,y'),f(x,y'),1,1) = result(y,x)
2227 /// also note the confusion with the counting of the particles/integration variables
2228 template <typename opT, typename R, std::size_t NDIM>
2229 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM>
2230 apply(const opT& op, const Function<R,NDIM>& f, bool fence=true) {
2231
2232 typedef TENSOR_RESULT_TYPE(typename opT::opT,R) resultT;
2233 Function<R,NDIM>& ff = const_cast< Function<R,NDIM>& >(f);
2235
2236 MADNESS_ASSERT(not f.is_on_demand());
2237 bool print_timings=op.print_timings;
2238
2239 if (VERIFY_TREE) ff.verify_tree();
2240 ff.reconstruct();
2241 if (print_timings) ff.print_size("ff in apply after reconstruct");
2242
2243 if (op.modified()) {
2244
2246// ff.get_impl()->make_redundant(true);
2247 result = apply_only(op, ff, fence);
2248 ff.get_impl()->undo_redundant(false);
2249 result.get_impl()->trickle_down(true);
2250
2251 } else {
2252
2253 // saves the standard() step, which is very expensive in 6D
2254// Function<R,NDIM> fff=copy(ff);
2255 Function<R,NDIM> fff=(ff);
2256 fff.make_nonstandard(op.doleaves, true);
2257 if (print_timings) fff.print_size("ff in apply after make_nonstandard");
2258 if ((print_timings) and (f.world().rank()==0)) {
2259 fff.get_impl()->timer_filter.print("filter");
2260 fff.get_impl()->timer_compress_svd.print("compress_svd");
2261 }
2262 result = apply_only(op, fff, fence);
2263 result.get_impl()->set_tree_state(nonstandard_after_apply);
2264 ff.world().gop.fence();
2265 if (print_timings) result.print_size("result after apply_only");
2266
2267 // svd-tensors need some post-processing
2268 if (result.get_impl()->get_tensor_type()==TT_2D) {
2269 double elapsed=result.get_impl()->finalize_apply();
2270 if (print_timings) printf("time in finalize_apply %8.2f\n",elapsed);
2271 }
2272 if (print_timings) {
2273 result.get_impl()->print_timer();
2274 op.print_timer();
2275 }
2276
2277 result.get_impl()->reconstruct(true);
2278
2279// fff.clear();
2280 if (op.destructive()) {
2281 ff.world().gop.fence();
2282 ff.clear();
2283 } else {
2284 // ff.standard();
2285 ff.reconstruct();
2286 }
2287
2288 }
2289 if (print_timings) result.print_size("result after reconstruction");
2290 return result;
2291 }
2292
2293
2294 template <typename opT, typename R, std::size_t NDIM>
2295 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM>
2296 apply_1d_realspace_push(const opT& op, const Function<R,NDIM>& f, int axis, bool fence=true) {
2298 Function<R,NDIM>& ff = const_cast< Function<R,NDIM>& >(f);
2299 if (VERIFY_TREE) ff.verify_tree();
2300 ff.reconstruct();
2301
2302 Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> result;
2303
2304 result.set_impl(ff, false);
2305 result.get_impl()->apply_1d_realspace_push(op, ff.get_impl().get(), axis, fence);
2306 result.get_impl()->set_tree_state(redundant_after_merge);
2307 return result;
2308 }
2309
2310
2311 /// Generate a new function by reordering dimensions ... optional fence
2312
2313 /// You provide an array of dimension NDIM that maps old to new dimensions
2314 /// according to
2315 /// \code
2316 /// newdim = mapdim[olddim]
2317 /// \endcode
2318 /// Works in either scaling function or wavelet basis.
2319 ///
2320 /// Would be easy to modify this to also change the procmap here
2321 /// if desired but presently it uses the same procmap as f.
2322 template <typename T, std::size_t NDIM>
2323 Function<T,NDIM>
2324 mapdim(const Function<T,NDIM>& f, const std::vector<long>& map, bool fence=true) {
2326 Function<T,NDIM> result;
2327 return result.mapdim(f,map,fence);
2328 }
2329
2330 /// Generate a new function by mirroring within the dimensions .. optional fence
2331
2332 /// similar to mapdim
2333 /// @param[in] mirror array with -1 and 1, corresponding to mirror this dimension or not
2334 template <typename T, std::size_t NDIM>
2335 Function<T,NDIM>
2336 mirror(const Function<T,NDIM>& f, const std::vector<long>& mirrormap, bool fence=true) {
2338 Function<T,NDIM> result;
2339 return result.mirror(f,mirrormap,fence);
2340 }
2341
2342 /// This is replaced with mirror(map(f)), optional fence
2343
2344 /// first map then mirror!
2345 /// mirror is similar to mapdim, but maps from x to -x, y to -y, and so on
2346 /// Example: mirror a 3d function on the xy plane: mirror={1,1,-1}
2347 /// Example: c4 rotation of a 3d function around the z axis:
2348 /// x->y, y->-x, z->z: map(1,0,2); mirror(-1,1,1)
2349 /// @param[in] map array holding dimensions
2350 /// @param[in] mirror array of -1 and 1, corresponding to mirror or not
2351 template <typename T, std::size_t NDIM>
2352 Function<T,NDIM>
2353 map_and_mirror(const Function<T,NDIM>& f, const std::vector<long>& map,
2354 const std::vector<long>& mirror, bool fence=true) {
2356 Function<T,NDIM> result;
2357 return result.map_and_mirror(f,map,mirror,fence);
2358 }
2359
2360
2361 /// swap particles 1 and 2
2362
2363 /// param[in] f a function of 2 particles f(1,2)
2364 /// return the input function with particles swapped g(1,2) = f(2,1)
2365 template <typename T, std::size_t NDIM>
2366 typename std::enable_if_t<NDIM%2==0, Function<T,NDIM>>
2368 // this could be done more efficiently for SVD, but it works decently
2369 std::vector<long> map(NDIM);
2370 constexpr std::size_t LDIM=NDIM/2;
2371 static_assert(LDIM*2==NDIM);
2372 for (auto d=0; d<LDIM; ++d) {
2373 map[d]=d+LDIM;
2374 map[d+LDIM]=d;
2375 }
2376// map[0]=3;
2377// map[1]=4;
2378// map[2]=5; // 2 -> 1
2379// map[3]=0;
2380// map[4]=1;
2381// map[5]=2; // 1 -> 2
2382 return mapdim(f,map);
2383 }
2384
2385 /// symmetrize a function
2386
2387 /// @param[in] symmetry possibilities are:
2388 /// (anti-) symmetric particle permutation ("sy_particle", "antisy_particle")
2389 /// symmetric mirror plane ("xy", "xz", "yz")
2390 /// @return a new function symmetrized according to the input parameter
2391 template <typename T, std::size_t NDIM>
2392 Function<T,NDIM>
2393 symmetrize(const Function<T,NDIM>& f, const std::string symmetry, bool fence=true) {
2394 Function<T,NDIM> result;
2395
2396 MADNESS_ASSERT(NDIM==6); // works only for pair functions
2397 std::vector<long> map(NDIM);
2398
2399 // symmetric particle permutation
2400 if (symmetry=="sy_particle") {
2401 map[0]=3; map[1]=4; map[2]=5;
2402 map[3]=0; map[4]=1; map[5]=2;
2403 } else if (symmetry=="cx") {
2404 map[0]=0; map[1]=2; map[2]=1;
2405 map[3]=3; map[4]=5; map[5]=4;
2406
2407 } else if (symmetry=="cy") {
2408 map[0]=2; map[1]=1; map[2]=0;
2409 map[3]=5; map[4]=4; map[5]=3;
2410
2411 } else if (symmetry=="cz") {
2412 map[0]=1; map[1]=0; map[2]=2;
2413 map[3]=4; map[4]=3; map[5]=5;
2414
2415 } else {
2416 if (f.world().rank()==0) {
2417 print("unknown parameter in symmetrize:",symmetry);
2418 }
2419 MADNESS_EXCEPTION("unknown parameter in symmetrize",1);
2420 }
2421
2422 result.mapdim(f,map,true); // need to fence here
2423 result.get_impl()->average(*f.get_impl());
2424
2425 return result;
2426 }
2427
2428
2429
2430 /// multiply a high-dimensional function with a low-dimensional function
2431
2432 /// @param[in] f NDIM function of 2 particles: f=f(1,2)
2433 /// @param[in] g LDIM function of 1 particle: g=g(1) or g=g(2)
2434 /// @param[in] particle if g=g(1) or g=g(2)
2435 /// @return h(1,2) = f(1,2) * g(p)
2436 template<typename T, std::size_t NDIM, std::size_t LDIM>
2437 Function<T,NDIM> multiply(const Function<T,NDIM> f, const Function<T,LDIM> g, const int particle, const bool fence=true) {
2438
2439 static_assert(LDIM+LDIM==NDIM);
2441
2442 Function<T,NDIM> result;
2443 result.set_impl(f, false);
2444
2445// Function<T,NDIM>& ff = const_cast< Function<T,NDIM>& >(f);
2446// Function<T,LDIM>& gg = const_cast< Function<T,LDIM>& >(g);
2447
2448 f.change_tree_state(redundant,false);
2449 g.change_tree_state(redundant);
2450 FunctionImpl<T,NDIM>* fimpl=f.get_impl().get();
2451 FunctionImpl<T,LDIM>* gimpl=g.get_impl().get();
2452
2453 result.get_impl()->multiply(fimpl,gimpl,particle);
2454 result.world().gop.fence();
2455
2456 f.change_tree_state(reconstructed,false);
2457 g.change_tree_state(reconstructed);
2458 return result;
2459 }
2460
2461
2462 template <typename T, std::size_t NDIM>
2463 Function<T,NDIM>
2467 bool fence=true)
2468 {
2471 other.reconstruct();
2472 result.get_impl()->project(*other.get_impl(),fence);
2473 return result;
2474 }
2475
2476
2477 /// Computes the scalar/inner product between two functions
2478
2479 /// In Maple this would be \c int(conjugate(f(x))*g(x),x=-infinity..infinity)
2480 template <typename T, typename R, std::size_t NDIM>
2483 return f.inner(g);
2484 }
2485
2486
2487 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2488
2489 /// syntax similar to the inner product in tensor.h
2490 /// e.g result=inner<3>(f,g),{0},{1}) : r(x,y) = int f(x1,x) g(y,x1) dx1
2491 /// @param[in] task 0: everything, 1; prepare only (fence), 2: work only (no fence), 3: finalize only (fence)
2492 template<std::size_t NDIM, typename T, std::size_t LDIM, typename R, std::size_t KDIM,
2493 std::size_t CDIM = (KDIM + LDIM - NDIM) / 2>
2494 std::vector<Function<TENSOR_RESULT_TYPE(T, R), NDIM>>
2495 innerXX(const Function<T, LDIM>& f, const std::vector<Function<R, KDIM>>& vg, const std::array<int, CDIM> v1,
2496 const std::array<int, CDIM> v2, int task=0) {
2497 bool prepare = ((task==0) or (task==1));
2498 bool work = ((task==0) or (task==2));
2499 bool finish = ((task==0) or (task==3));
2500
2501 static_assert((KDIM + LDIM - NDIM) % 2 == 0, "faulty dimensions in inner (partial version)");
2502 static_assert(KDIM + LDIM - 2 * CDIM == NDIM, "faulty dimensions in inner (partial version)");
2503
2504 // contraction indices must be contiguous and either in the beginning or at the end
2505 for (int i=0; i<CDIM-1; ++i) MADNESS_CHECK((v1[i]+1)==v1[i+1]);
2506 MADNESS_CHECK((v1[0]==0) or (v1[CDIM-1]==LDIM-1));
2507
2508 for (int i=0; i<CDIM-1; ++i) MADNESS_CHECK((v2[i]+1)==v2[i+1]);
2509 MADNESS_CHECK((v2[0]==0) or (v2[CDIM-1]==KDIM-1));
2510
2511 MADNESS_CHECK(f.is_initialized());
2513 MADNESS_CHECK(f.world().id() == vg[0].world().id());
2514 // this needs to be run in a single world, so that all coefficients are local.
2515 // Use macrotasks if run on multiple processes.
2516 World& world=f.world();
2517 MADNESS_CHECK(world.size() == 1);
2518
2519 if (prepare) {
2520 f.change_tree_state(nonstandard);
2522 world.gop.fence();
2523 f.get_impl()->compute_snorm_and_dnorm(false);
2524 for (auto& g : vg) g.get_impl()->compute_snorm_and_dnorm(false);
2525 world.gop.fence();
2526 }
2527
2528 typedef TENSOR_RESULT_TYPE(T, R) resultT;
2529 std::vector<Function<resultT,NDIM>> result(vg.size());
2530 if (work) {
2532 for (int i=0; i<vg.size(); ++i) {
2534 .k(f.k()).thresh(f.thresh()).empty().nofence();
2535 result[i].get_impl()->partial_inner(*f.get_impl(),*(vg[i]).get_impl(),v1,v2);
2536 result[i].get_impl()->set_tree_state(nonstandard_after_apply);
2537 }
2538 world.gop.set_forbid_fence(false);
2539 }
2540
2541 if (finish) {
2542
2543 world.gop.fence();
2544// result.get_impl()->reconstruct(true);
2545
2547// result.reconstruct();
2548 // restore initial state of g and h
2549 auto erase_list = [] (const auto& funcimpl) {
2550 typedef typename std::decay_t<decltype(funcimpl)>::keyT keyTT;
2551 std::list<keyTT> to_be_erased;
2552 for (auto it=funcimpl.get_coeffs().begin(); it!=funcimpl.get_coeffs().end(); ++it) {
2553 const auto& key=it->first;
2554 const auto& node=it->second;
2555 if (not node.has_children()) to_be_erased.push_back(key);
2556 }
2557 return to_be_erased;
2558 };
2559
2560 FunctionImpl<T,LDIM>& f_nc=const_cast<FunctionImpl<T,LDIM>&>(*f.get_impl());
2561 for (auto& key : erase_list(f_nc)) f_nc.get_coeffs().erase(key);
2562 for (auto& g : vg) {
2563 FunctionImpl<R,KDIM>& g_nc=const_cast<FunctionImpl<R,KDIM>&>(*g.get_impl());
2564 for (auto& key : erase_list(g_nc)) g_nc.get_coeffs().erase(key);
2565 }
2566 world.gop.fence();
2568 f_nc.reconstruct(false);
2569 world.gop.fence();
2570
2571 }
2572
2573 return result;
2574 }
2575
2576
2577 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2578
2579 /// syntax similar to the inner product in tensor.h
2580 /// e.g result=inner<3>(f,g),{0},{1}) : r(x,y) = int f(x1,x) g(y,x1) dx1
2581 /// @param[in] task 0: everything, 1; prepare only (fence), 2: work only (no fence), 3: finalize only (fence)
2582 template<std::size_t NDIM, typename T, std::size_t LDIM, typename R, std::size_t KDIM,
2583 std::size_t CDIM = (KDIM + LDIM - NDIM) / 2>
2585 innerXX(const Function<T, LDIM>& f, const Function<R, KDIM>& g, const std::array<int, CDIM> v1,
2586 const std::array<int, CDIM> v2, int task=0) {
2587 return innerXX<NDIM,T,LDIM,R,KDIM>(f,std::vector<Function<R,KDIM>>({g}),v1,v2,task)[0];
2588 }
2589
2590 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2591
2592 /// syntax similar to the inner product in tensor.h
2593 /// e.g result=inner<3>(f,g),{0},{1}) : r(x,y) = int f(x1,x) g(y,x1) dx1
2594 template <typename T, std::size_t LDIM, typename R, std::size_t KDIM>
2596 inner(const Function<T,LDIM>& f, const Function<R,KDIM>& g, const std::tuple<int> v1, const std::tuple<int> v2) {
2597 return innerXX<KDIM+LDIM-2>(f,g,
2598 std::array<int,1>({std::get<0>(v1)}),
2599 std::array<int,1>({std::get<0>(v2)}));
2600 }
2601
2602 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2603
2604 /// syntax similar to the inner product in tensor.h
2605 /// e.g result=inner<3>(f,g),{0,1},{1,2}) : r(y) = int f(x1,x2) g(y,x1,x2) dx1 dx2
2606 template <typename T, std::size_t LDIM, typename R, std::size_t KDIM>
2608 inner(const Function<T,LDIM>& f, const Function<R,KDIM>& g, const std::tuple<int,int> v1, const std::tuple<int,int> v2) {
2609 return innerXX<KDIM+LDIM-4>(f,g,
2610 std::array<int,2>({std::get<0>(v1),std::get<1>(v1)}),
2611 std::array<int,2>({std::get<0>(v2),std::get<1>(v2)}));
2612 }
2613
2614 /// Computes the partial scalar/inner product between two functions, returns a low-dim function
2615
2616 /// syntax similar to the inner product in tensor.h
2617 /// e.g result=inner<3>(f,g),{1},{2}) : r(x,y,z) = int f(x,x1) g(y,z,x1) dx1
2618 template <typename T, std::size_t LDIM, typename R, std::size_t KDIM>
2620 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) {
2621 return innerXX<KDIM+LDIM-6>(f,g,
2622 std::array<int,3>({std::get<0>(v1),std::get<1>(v1),std::get<2>(v1)}),
2623 std::array<int,3>({std::get<0>(v2),std::get<1>(v2),std::get<2>(v2)}));
2624 }
2625
2626
2627
2628 /// Computes the scalar/inner product between an MRA function and an external functor
2629
2630 /// Currently this defaults to inner_adaptive, which might be more expensive
2631 /// than inner_ext since it loops over all leaf nodes. If you feel inner_ext
2632 /// is more efficient you need to call it directly
2633 /// @param[in] f MRA function
2634 /// @param[in] g functor
2635 /// @result inner(f,g)
2636 template <typename T, typename opT, std::size_t NDIM>
2637 TENSOR_RESULT_TYPE(T,typename opT::value_type) inner(const Function<T,NDIM>& f, const opT& g) {
2639 std::shared_ptr< FunctionFunctorInterface<double,3> > func(new opT(g));
2640 return f.inner_adaptive(func);
2641 }
2642
2643 /// Computes the scalar/inner product between an MRA function and an external functor
2644
2645 /// Currently this defaults to inner_adaptive, which might be more expensive
2646 /// than inner_ext since it loops over all leaf nodes. If you feel inner_ext
2647 /// is more efficient you need to call it directly
2648 /// @param[in] g functor
2649 /// @param[in] f MRA function
2650 /// @result inner(f,g)
2651 template <typename T, typename opT, std::size_t NDIM>
2652 TENSOR_RESULT_TYPE(T,typename opT::value_type) inner(const opT& g, const Function<T,NDIM>& f) {
2653 return inner(f,g);
2654 }
2655
2656 template <typename T, typename R, std::size_t NDIM>
2657 typename IsSupported<TensorTypeData<R>, Function<TENSOR_RESULT_TYPE(T,R),NDIM> >::type
2659 return (f*R(1.0)).add_scalar(r);
2660 }
2661
2662 template <typename T, typename R, std::size_t NDIM>
2663 typename IsSupported<TensorTypeData<R>, Function<TENSOR_RESULT_TYPE(T,R),NDIM> >::type
2665 return (f*R(1.0)).add_scalar(r);
2666 }
2667
2668 template <typename T, typename R, std::size_t NDIM>
2669 typename IsSupported<TensorTypeData<R>, Function<TENSOR_RESULT_TYPE(T,R),NDIM> >::type
2671 return (f*R(1.0)).add_scalar(-r);
2672 }
2673
2674 template <typename T, typename R, std::size_t NDIM>
2675 typename IsSupported<TensorTypeData<R>, Function<TENSOR_RESULT_TYPE(T,R),NDIM> >::type
2677 return (f*R(-1.0)).add_scalar(r);
2678 }
2679
2680 namespace detail {
2681 template <std::size_t NDIM>
2682 struct realop {
2683 typedef double resultT;
2685 return real(t);
2686 }
2687
2688 template <typename Archive> void serialize (Archive& ar) {}
2689 };
2690
2691 template <std::size_t NDIM>
2692 struct imagop {
2693 typedef double resultT;
2695 return imag(t);
2696 }
2697
2698 template <typename Archive> void serialize (Archive& ar) {}
2699 };
2700
2701 template <std::size_t NDIM>
2702 struct abssqop {
2703 typedef double resultT;
2705 Tensor<double> r = abs(t);
2706 return r.emul(r);
2707 }
2708
2709 template <typename Archive> void serialize (Archive& ar) {}
2710 };
2711
2712 template <std::size_t NDIM>
2713 struct absop {
2714 typedef double resultT;
2716 Tensor<double> r = abs(t);
2717 return r;
2718 }
2719
2720 template <typename Archive> void serialize (Archive& ar) {}
2721 };
2722
2723 }
2724
2725 /// Returns a new function that is the real part of the input
2726 template <std::size_t NDIM>
2730
2731 /// Returns a new function that is the real part of the input
2732 template <std::size_t NDIM>
2734 return copy(z);
2735 }
2736
2737 /// Returns a new function that is the imaginary part of the input
2738 template <std::size_t NDIM>
2742
2743
2744 /// Create a new function that is the square of f - global comm only if not reconstructed
2745 template <typename T, std::size_t NDIM>
2748 Function<T,NDIM> result = copy(f,true); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2749 return result.square(true); //fence); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2750 }
2751
2752 /// Create a new function that is the abs of f - global comm only if not reconstructed
2753 template <typename T, std::size_t NDIM>
2754 Function<T,NDIM> abs(const Function<T,NDIM>& f, bool fence=true) {
2756 Function<T,NDIM> result = copy(f,true); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2757 return result.abs(true); //fence); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2758 }
2759
2760 /// Create a new function that is the abs_square of f - global comm only if not reconstructed
2761 template <typename T, std::size_t NDIM>
2762 typename std::enable_if<!TensorTypeData<T>::iscomplex, Function<T,NDIM> >::type
2763 abs_square(const Function<T,NDIM>& f, bool fence=true) {
2765 Function<T,NDIM> result = copy(f,true); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2766 return result.abs_square(true); //fence); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2767 }
2768
2769 /// Create a new function that is the abs_square of f - global comm only if not reconstructed
2770 template <typename T, std::size_t NDIM>
2771 typename std::enable_if<TensorTypeData<T>::iscomplex, Function<typename Tensor<T>::scalar_type,NDIM> >::type
2772 abs_square(const Function<T,NDIM>& f, bool fence=true) {
2774 }
2775
2776 /// Returns a new function that is the square of the absolute value of the input
2777 template <std::size_t NDIM>
2781
2782 /// Returns a new function that is the absolute value of the input
2783 template <std::size_t NDIM>
2787
2788 /// get tree state of a function
2789
2790 /// there is a corresponding function in vmra.h
2791 /// @param[in] f function
2792 /// @return TreeState::unknown if the function is not initialized
2793 template <typename T, std::size_t NDIM>
2795 if (f.is_initialized()) return f.get_impl()->get_tree_state();
2796 return TreeState::unknown;
2797 }
2798
2799 /// change tree state of a function
2800
2801 /// there is a corresponding function in vmra.h
2802 /// return this for chaining
2803 /// @param[in] f function
2804 /// @param[in] finalstate the new state
2805 /// @return this in the requested state
2806 template <typename T, std::size_t NDIM>
2808 const TreeState finalstate, bool fence=true) {
2809 return f.change_tree_state(finalstate,fence);
2810 }
2811
2812
2813}
2814
2815#include <madness/mra/funcplot.h>
2816
2817namespace madness {
2818 namespace archive {
2819 template <class archiveT, class T, std::size_t NDIM>
2821 static inline void load(const ParallelInputArchive<archiveT>& ar, Function<T,NDIM>& f) {
2822 f.load(*ar.get_world(), ar);
2823 }
2824 };
2825
2826 template <class archiveT, class T, std::size_t NDIM>
2828 static inline void store(const ParallelOutputArchive<archiveT>& ar, const Function<T,NDIM>& f) {
2829 f.store(ar);
2830 }
2831 };
2832 }
2833
2834 template <class T, std::size_t NDIM>
2835 void save(const Function<T,NDIM>& f, const std::string name) {
2837 ar2 & f;
2838 }
2839
2840 template <class T, std::size_t NDIM>
2841 void load(Function<T,NDIM>& f, const std::string name) {
2843 ar2 & f;
2844 }
2845
2846}
2847
2848namespace madness {
2849 // type traits to check if a template parameter is a Function
2850 template<typename>
2851 struct is_madness_function : std::false_type {};
2852
2853 template<typename T, std::size_t NDIM>
2854 struct is_madness_function<madness::Function<T, NDIM>> : std::true_type {};
2855
2856 template<typename>
2857 struct is_madness_function_vector : std::false_type {
2858 };
2859
2860 template<typename T, std::size_t NDIM>
2861 struct is_madness_function_vector<std::vector<typename madness::Function<T, NDIM>>> : std::true_type {
2862};
2863
2864}
2865
2866
2867/* @} */
2868
2869#include <madness/mra/derivative.h>
2870#include <madness/mra/operator.h>
2872#include <madness/mra/vmra.h>
2873// #include <madness/mra/mraimpl.h> !!!!!!!!!!!!! NOOOOOOOOOOOOOOOOOOOOOOOOOOOOOOOOOO !!!!!!!!!!!!!!!!!!
2874
2875#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:282
virtual FunctionFactory & thresh(double thresh)
Definition function_factory.h:198
FunctionFactory & fence(bool fence=true)
Definition function_factory.h:276
virtual FunctionFactory & k(int k)
Definition function_factory.h:193
FunctionFactory & empty()
Definition function_factory.h:246
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:339
void reconstruct(bool fence)
reconstruct this tree – respects fence
Definition mraimpl.h:1464
bool is_on_demand() const
Definition mraimpl.h:279
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:906
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:1706
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:1422
T inner_adaptive(const std::shared_ptr< FunctionFunctorInterface< T, NDIM > > f, const bool leaf_refine=true) const
Definition mra.h:1407
void unaryop_coeff(const opT &op, bool fence=true)
Unary operation applied inplace to the coefficients.
Definition mra.h:956
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:1486
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:682
bool autorefine() const
Returns value of autorefine flag. No communication.
Definition mra.h:588
TreeState operating_state
Definition mra.h:1348
Function< T, NDIM > & add_scalar(T t, bool fence=true)
Inplace add scalar. No communication except for optional fence.
Definition mra.h:1005
void print_size(const std::string name) const
print some info about this
Definition mra.h:533
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:1192
void broaden(const BoundaryConditions< NDIM > &bc=FunctionDefaults< NDIM >::get_bc(), bool fence=true) const
Inplace broadens support in scaling function basis.
Definition mra.h:878
void print_info() const
Print a summary of the load balancing info.
Definition mra.h:922
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:995
void norm_tree(bool fence=true) const
Initializes information about the function norm at all length scales.
Definition mra.h:743
void load(World &world, Archive &ar)
Replaces this function with one loaded from an archive using the default processor map.
Definition mra.h:1503
Function< T, NDIM > & abs_square(bool fence=true)
Returns *this for chaining.
Definition mra.h:1132
double norm2sq_local() const
Returns the square of the norm of the local function ... no communication.
Definition mra.h:718
IsSupported< TensorTypeData< Q >, Function< T, NDIM > >::type & operator*=(const Q q)
Inplace scaling by a constant.
Definition mra.h:1104
void sum_down(bool fence=true) const
Sums scaling coeffs down tree restoring state with coeffs only at leaves. Optional fence....
Definition mra.h:840
Level depthpt(const coordT &xuser) const
Definition mra.h:426
bool is_redundant() const
Returns true if redundant, false otherwise. No communication.
Definition mra.h:512
Function< T, NDIM > & operator+=(const Function< Q, NDIM > &other)
Inplace addition of functions in the wavelet basis.
Definition mra.h:1064
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:598
World & world() const
Returns the world.
Definition mra.h:688
T trace() const
Returns global value of int(f(x),x) ... global comm required.
Definition mra.h:1154
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:1210
const Function< T, NDIM > & change_tree_state(const TreeState finalstate, bool fence=true) const
changes tree state to given state
Definition mra.h:827
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:1544
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:912
double norm2() const
Returns the 2-norm of the function ... global sum ... works in either basis.
Definition mra.h:730
Function< T, NDIM > & fill_cuspy_tree(const opT &op, const bool fence=true)
Definition mra.h:1233
void change_tensor_type(const TensorArgs &targs, bool fence=true)
change the tensor type of the coefficients in the FunctionNode
Definition mra.h:1536
const Function< T, NDIM > & refine(bool fence=true) const
Inplace autorefines the function using same test as for squaring.
Definition mra.h:872
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:1114
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:1221
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:1302
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:1672
impl world gop fence()
int k() const
Returns the number of multiwavelets (k). No communication.
Definition mra.h:626
const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > & get_pmap() const
Returns a shared pointer to the process map.
Definition mra.h:696
double thresh() const
Returns value of truncation threshold. No communication.
Definition mra.h:607
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:1658
void set_thresh(double value, bool fence=true)
Sets the value of the truncation threshold. Optional global fence.
Definition mra.h:617
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:1566
void distribute(std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > newmap) const
distribute this function according to newmap
Definition mra.h:709
void unaryop(const opT &op, bool fence=true)
Inplace unary operation on function values.
Definition mra.h:946
const std::shared_ptr< FunctionImpl< T, NDIM > > & get_impl() const
Returns a shared-pointer to the implementation.
Definition mra.h:654
void standard(bool fence=true)
Converts the function standard compressed form. Possible non-blocking comm.
Definition mra.h:790
void unaryop_node(const opT &op, bool fence=true)
Unary operation applied inplace to the nodes.
Definition mra.h:966
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:642
Function< T, NDIM > & fill_cuspy_tree(const bool fence=true)
Special refinement on 6D boxes where the electrons come close (meet)
Definition mra.h:1246
std::size_t size() const
Returns the number of coefficients in the function ... collective global sum.
Definition mra.h:573
bool is_on_demand() const
Definition mra.h:676
Function< T, NDIM > & operator-=(const Function< Q, NDIM > &other)
Inplace subtraction of functions in the wavelet basis.
Definition mra.h:1083
bool compressed
Definition mra.h:1168
const Function< T, NDIM > & reconstruct(bool fence=true) const
Reconstructs the function, transforming into scaling function basis. Possible non-blocking comm.
Definition mra.h:817
Function< T, NDIM > & abs(bool fence=true)
Returns *this for chaining.
Definition mra.h:1123
Function< T, NDIM > & reduce_rank(const double thresh=0.0, const bool fence=true)
reduce the rank of the coefficient tensors
Definition mra.h:1736
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:1687
std::size_t max_nodes() const
Returns the max number of nodes on a processor.
Definition mra.h:558
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:1288
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:703
T trace_local() const
Returns local contribution to int(f(x),x) ... no communication.
Definition mra.h:1145
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:1629
Function< T, NDIM > & gaxpy(const T &alpha, const Function< Q, NDIM > &other, const R &beta, bool fence=true)
Inplace, general bi-linear operation in wavelet basis. No communication except for optional fence.
Definition mra.h:1025
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:1524
std::size_t max_local_depth() const
Returns the maximum local depth of the function tree ... no communications.
Definition mra.h:550
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:1647
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:542
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:526
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:1462
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:1557
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:889
void refine_general(const opT &op, bool fence=true) const
Inplace autorefines the function. Optional fence. Possible non-blocking comm.
Definition mra.h:853
static void doconj(const Key< NDIM >, Tensor< T > &t)
Definition mra.h:976
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:671
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:1261
bool is_redundant_after_merge() const
Returns true if redundant_after_merge, false otherwise. No communication.
Definition mra.h:520
return local
Definition mra.h:1363
auto func
Definition mra.h:1441
~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:661
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:1374
std::size_t min_nodes() const
Returns the min number of nodes on a processor.
Definition mra.h:565
constexpr std::size_t LDIM
Definition mra.h:1440
static constexpr std::size_t dimT
Definition mra.h:158
change_tree_state(operating_state, false)
bool is_nonstandard() const
Returns true if nonstandard-compressed, false otherwise. No communication.
Definition mra.h: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:1390
void make_redundant(bool fence=true)
Converts the function to redundant form, i.e. sum coefficients on all levels.
Definition mra.h:802
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:1274
double check_symmetry() const
check symmetry of a function by computing the 2nd derivative
Definition mra.h:1722
void multi_to_multi_op_values(const opT &op, const std::vector< Function< T, NDIM > > &vin, std::vector< Function< T, NDIM > > &vout, const bool fence=true)
apply op on the input vector yielding an output vector of functions
Definition mra.h:1585
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:580
const Function< T, NDIM > & compress(bool fence=true) const
Compresses the function, transforming into wavelet basis. Possible non-blocking comm.
Definition mra.h:763
bool is_initialized() const
Returns true if the function is initialized.
Definition mra.h:167
bool is_reconstructed() const
Returns true if reconstructed, false otherwise. No communication.
Definition mra.h: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:1744
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:1606
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:776
if(VERIFY_TREE) verify_tree()
Future< T > eval(const coordT &xuser) const
Evaluates the function at a point in user coordinates. Possible non-blocking comm.
Definition mra.h: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:899
FunctionNode< T, NDIM > nodeT
Definition mra.h:154
bool redundant
Definition mra.h:1169
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
NDIM &g const
Definition mra.h:1166
void unaryop(T(*f)(T))
Inplace unary operation on function values.
Definition mra.h:937
Function< T, NDIM > conj(bool fence=true)
Inplace complex conjugate. No communication except for optional fence.
Definition mra.h:984
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:69
Definition leafop.h:412
Definition leafop.h:281
Traits class to specify support of numeric types.
Definition type_data.h:56
A tensor is a multidimensional array.
Definition tensor.h:317
Tensor< T > & emul(const Tensor< T > &t)
Inplace multiply by corresponding elements of argument Tensor.
Definition tensor.h:1799
Tensor< T > & conj()
Inplace complex conjugate.
Definition tensor.h:716
A simple, fixed dimension vector.
Definition vector.h:64
void erase(const keyT &key)
Erases entry from container (non-blocking comm if remote)
Definition worlddc.h: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
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
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:2778
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:2296
Function< TENSOR_RESULT_TYPE(L, R), NDIM > sub(const Function< L, NDIM > &left, const Function< R, NDIM > &right, bool fence=true)
Same as operator- but with optional fence and no automatic compression.
Definition mra.h:2035
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:1830
Function< Q, NDIM > convert(const Function< T, NDIM > &f, bool fence=true)
Type conversion implies a deep copy. No communication except for optional fence.
Definition mra.h:2096
Function< TENSOR_RESULT_TYPE(Q, T), NDIM > mul(const Q alpha, const Function< T, NDIM > &f, bool fence=true)
Returns new function equal to alpha*f(x) with optional fence.
Definition mra.h:1765
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:2367
TreeState
Definition funcdefaults.h:59
@ nonstandard_after_apply
s and d coeffs, state after operator application
Definition funcdefaults.h:64
@ redundant_after_merge
s coeffs everywhere, must be summed up to yield the result
Definition funcdefaults.h:66
@ reconstructed
s coeffs at the leaves only
Definition funcdefaults.h:60
@ nonstandard
s and d coeffs in internal nodes
Definition funcdefaults.h:62
@ unknown
Definition funcdefaults.h:68
@ nonstandard_with_leaves
like nonstandard, with s coeffs at the leaves
Definition funcdefaults.h:63
static void user_to_sim(const Vector< double, NDIM > &xuser, Vector< double, NDIM > &xsim)
Convert user coords (cell[][]) to simulation coords ([0,1]^ndim)
Definition funcdefaults.h:425
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:2495
std::vector< CCPairFunction< T, NDIM > > operator*(const double fac, const std::vector< CCPairFunction< T, NDIM > > &arg)
Definition ccpairfunction.h:1089
int Level
Definition key.h:58
TreeState get_tree_state(const Function< T, NDIM > &f)
get tree state of a function
Definition mra.h:2794
std::vector< CCPairFunction< T, NDIM > > operator-(const std::vector< CCPairFunction< T, NDIM > > c1, const std::vector< CCPairFunction< T, NDIM > > &c2)
Definition ccpairfunction.h:1060
Function< TENSOR_RESULT_TYPE(L, R), NDIM > mul_sparse(const Function< L, NDIM > &left, const Function< R, NDIM > &right, double tol, bool fence=true)
Sparse multiplication — left and right must be reconstructed and if tol!=0 have tree of norms already...
Definition mra.h:1806
std::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:1999
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
@ TT_FULL
Definition gentensor.h:120
NDIM & f
Definition mra.h:2481
Function< TENSOR_RESULT_TYPE(L, R), NDIM > add(const Function< L, NDIM > &left, const Function< R, NDIM > &right, bool fence=true)
Same as operator+ but with optional fence and no automatic compression.
Definition mra.h:1991
Function< T, NDIM > symmetrize(const Function< T, NDIM > &f, const std::string symmetry, bool fence=true)
symmetrize a function
Definition mra.h:2393
Function< TENSOR_RESULT_TYPE(typename opT::opT, R), NDIM > apply_only(const opT &op, const Function< R, NDIM > &f, bool fence=true)
Apply operator ONLY in non-standard form - required other steps missing !!
Definition mra.h:2184
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:1844
void startup(World &world, int argc, char **argv, bool doprint=false, bool make_stdcout_nice_to_reals=true)
initialize the internal state of the MADmra library
Definition startup.cc:64
std::string type(const PairType &n)
Definition PNOParameters.h:18
static bool print_timings
Definition SCF.cc: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:2437
void load(Function< T, NDIM > &f, const std::string name)
Definition mra.h:2841
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:2464
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:2835
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:1896
Function< T, NDIM > copy(const Function< T, NDIM > &f, const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > &pmap, bool fence=true)
Create a new copy of the function with different distribution and optional fence.
Definition mra.h:2066
Definition mraimpl.h:50
static const double d
Definition nonlinschro.cc:121
Implements most functionality of separated operators.
Implements ParallelInputArchive and ParallelOutputArchive for parallel serialization of data.
double Q(double a)
Definition relops.cc:20
static const double c
Definition relops.cc:10
static const double L
Definition rk.cc:46
Definition test_ar.cc:204
Definition leafop.h:142
void serialize(Archive &ar)
Definition mra.h:933
T(* f)(T)
Definition mra.h:928
void operator()(const Key< NDIM > &key, Tensor< T > &t) const
Definition mra.h:930
SimpleUnaryOpWrapper(T(*f)(T))
Definition mra.h:929
void serialize(Archive &ar)
Definition mra.h:866
bool operator()(implT *impl, const Key< NDIM > &key, const nodeT &t) const
Definition mra.h:862
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:2821
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:2828
Default store of an object via serialize(ar, t).
Definition archive.h:611
Definition mra.h:2713
Tensor< double > operator()(const Key< NDIM > &key, const Tensor< double_complex > &t) const
Definition mra.h:2715
double resultT
Definition mra.h:2714
void serialize(Archive &ar)
Definition mra.h:2720
Definition mra.h:2702
void serialize(Archive &ar)
Definition mra.h:2709
Tensor< double > operator()(const Key< NDIM > &key, const Tensor< double_complex > &t) const
Definition mra.h:2704
double resultT
Definition mra.h:2703
Definition mra.h:2692
Tensor< double > operator()(const Key< NDIM > &key, const Tensor< double_complex > &t) const
Definition mra.h:2694
void serialize(Archive &ar)
Definition mra.h:2698
double resultT
Definition mra.h:2693
Definition mra.h:2682
Tensor< double > operator()(const Key< NDIM > &key, const Tensor< double_complex > &t) const
Definition mra.h:2684
double resultT
Definition mra.h:2683
void serialize(Archive &ar)
Definition mra.h:2688
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:2851
Definition mra.h:112
Definition mra.h:115
Definition funcimpl.h:564
Definition lowrankfunction.h:332
double real(double a)
Definition tdse4.cc:78
Defines and implements most of Tensor.
#define UNARY_OPTIMIZED_ITERATOR(X, x, exp)
Definition tensor_macros.h:658
template Tensor< Tensor< double_complex >::scalar_type > imag(const Tensor< double_complex > &t)
AtomicInt sum
Definition test_atomicint.cc:46
double norm(const T i1)
Definition test_cloud.cc: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