MADNESS 0.10.1
vmra.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 $Id$
32*/
33#ifndef MADNESS_MRA_VMRA_H__INCLUDED
34#define MADNESS_MRA_VMRA_H__INCLUDED
35
36/*!
37 \file vmra.h
38 \brief Defines operations on vectors of Functions
39 \ingroup mra
40
41 This file defines a number of operations on vectors of functions.
42 Assume v is a vector of NDIM-D functions of a certain type.
43
44
45 Operations on array of functions
46
47 *) copying: deep copying of vectors of functions to vector of functions
48 \code
49 vector2 = copy(world, vector1,fence);
50 \endcode
51
52 *) compress: convert multiwavelet representation to legendre representation
53 \code
54 compress(world, vector, fence);
55 \endcode
56
57 *) reconstruct: convert representation to multiwavelets
58 \code
59 reconstruct(world, vector, fence);
60 \endcode
61
62 *) make_nonstandard: convert to non-standard form
63 \code
64 make_nonstandard(world, v, fence);
65 \endcode
66
67 *) standard: convert to standard form
68 \code
69 standard(world, v, fence);
70 \endcode
71
72 *) truncate: truncating vectors of functions to desired precision
73 \code
74 truncate(world, v, tolerance, fence);
75 \endcode
76
77
78 *) zero function: create a vector of zero functions of length n
79 \code
80 v=zero(world, n);
81 \endcode
82
83 *) transform: transform a representation from one basis to another
84 \code
85 transform(world, vector, tensor, tolerance, fence )
86 \endcode
87
88 Setting thresh-hold for precision
89
90 *) set_thresh: setting a finite thresh-hold for a vector of functions
91 \code
92 void set_thresh(World& world, std::vector< Function<T,NDIM> >& v, double thresh, bool fence=true);
93 \endcode
94
95 Arithmetic Operations on arrays of functions
96
97 *) conjugation: conjugate a vector of complex functions
98
99 *) add
100 *) sub
101 *) mul
102 - mul_sparse
103 *) square
104 *) gaxpy
105 *) apply
106
107 Norms, inner-products, blas-1 like operations on vectors of functions
108
109 *) inner
110 *) matrix_inner
111 *) norm_tree
112 *) normalize
113 *) norm2
114 - norm2s
115 *) scale(world, v, alpha);
116
117
118
119
120*/
121
122#include <madness/mra/mra.h>
125#include <cstdio>
126
127namespace madness {
128
129
130 /// get tree state of a vector of functions
131
132 /// @return TreeState::unknown if the vector is empty or if the functions have different tree states
133 template <typename T, std::size_t NDIM>
135 if (v.size()==0) return TreeState::unknown;
136 // return unknown if any function is not initialized
137 if (std::any_of(v.begin(), v.end(), [](const Function<T,NDIM>& f) {return not f.is_initialized();})) {
138 return TreeState::unknown;
139 }
140 TreeState state=v[0].get_impl()->get_tree_state();
141 for (const auto& f : v) {
142 if (f.get_impl()->get_tree_state()!=state) state=TreeState::unknown;
143 }
144 return state;
145 }
146
147 /// Compress a vector of functions
148 template <typename T, std::size_t NDIM>
149 void compress(World& world,
150 const std::vector< Function<T,NDIM> >& v,
151 bool fence=true) {
152 PROFILE_BLOCK(Vcompress);
154 }
155
156
157 /// reconstruct a vector of functions
158
159 /// implies fence
160 /// return v for chaining
161 template <typename T, std::size_t NDIM>
162 const std::vector< Function<T,NDIM> >& reconstruct(const std::vector< Function<T,NDIM> >& v) {
164 }
165
166 /// compress a vector of functions
167
168 /// implies fence
169 /// return v for chaining
170 template <typename T, std::size_t NDIM>
171 const std::vector< Function<T,NDIM> >& compress(const std::vector< Function<T,NDIM> >& v) {
173 }
174
175 /// Reconstruct a vector of functions
176 template <typename T, std::size_t NDIM>
177 void reconstruct(World& world,
178 const std::vector< Function<T,NDIM> >& v,
179 bool fence=true) {
180 PROFILE_BLOCK(Vreconstruct);
182 }
183
184 /// change tree_state of a vector of functions to redundant
185 template <typename T, std::size_t NDIM>
187 const std::vector< Function<T,NDIM> >& v,
188 bool fence=true) {
189
190 PROFILE_BLOCK(Vcompress);
192 }
193
194 /// refine the functions according to the autorefine criteria
195 template <typename T, std::size_t NDIM>
196 void refine(World& world, const std::vector<Function<T,NDIM> >& vf,
197 bool fence=true) {
198 for (const auto& f : vf) f.refine(false);
199 if (fence) world.gop.fence();
200 }
201
202 /// refine all functions to a common (finest) level
203
204 /// if functions are not initialized (impl==NULL) they are ignored
205 template <typename T, std::size_t NDIM>
206 void refine_to_common_level(World& world, std::vector<Function<T,NDIM> >& vf,
207 bool fence=true) {
208
209 reconstruct(world,vf);
211 std::vector<FunctionImpl<T,NDIM>*> v_ptr;
212
213 // push initialized function pointers into the vector v_ptr
214 for (unsigned int i=0; i<vf.size(); ++i) {
215 if (vf[i].is_initialized()) v_ptr.push_back(vf[i].get_impl().get());
216 }
217
218 // sort and remove duplicates to not confuse the refining function
219 std::sort(v_ptr.begin(),v_ptr.end());
220 typename std::vector<FunctionImpl<T, NDIM>*>::iterator it;
221 it = std::unique(v_ptr.begin(), v_ptr.end());
222 v_ptr.resize( std::distance(v_ptr.begin(),it) );
223
224 std::vector< Tensor<T> > c(v_ptr.size());
225 v_ptr[0]->refine_to_common_level(v_ptr, c, key0);
226 if (fence) v_ptr[0]->world.gop.fence();
227 if (VERIFY_TREE)
228 for (unsigned int i=0; i<vf.size(); i++) vf[i].verify_tree();
229 }
230
231 /// Generates non-standard form of a vector of functions
232 template <typename T, std::size_t NDIM>
234 std::vector< Function<T,NDIM> >& v,
235 bool fence= true) {
236 PROFILE_BLOCK(Vnonstandard);
238 }
239
240
241 /// Generates standard form of a vector of functions
242 template <typename T, std::size_t NDIM>
243 void standard(World& world,
244 std::vector< Function<T,NDIM> >& v,
245 bool fence=true) {
246 PROFILE_BLOCK(Vstandard);
248 }
249
250
251 /// change tree state of the functions
252
253 /// might not respect fence
254 /// @return v for chaining
255 template <typename T, std::size_t NDIM>
256 const std::vector<Function<T,NDIM>>& change_tree_state(const std::vector<Function<T,NDIM>>& v,
257 const TreeState finalstate,
258 const bool fence=true) {
259 // fast return
260 if (v.size()==0) return v;
261 if (get_tree_state(v)==finalstate) return v;
262
263 // find initialized function with world
264 Function<T,NDIM> dummy;
265 for (const auto& f : v)
266 if (f.is_initialized()) {
267 dummy=f;
268 break;
269 }
270 if (not dummy.is_initialized()) return v;
271 World& world=dummy.world();
272
273
274 // if a tree state cannot directly be changed to finalstate, we need to go via intermediate
275 auto change_initial_to_intermediate =[](const std::vector<Function<T,NDIM>>& v,
276 const TreeState initialstate,
277 const TreeState intermediatestate) {
278 int must_fence=0;
279 for (auto& f : v) {
280 if (f.is_initialized() and f.get_impl()->get_tree_state()==initialstate) {
281 f.change_tree_state(intermediatestate,false);
282 must_fence=1;
283 }
284 }
285 return must_fence;
286 };
287
288 int do_fence=0;
289 if (finalstate==compressed) {
290 do_fence+=change_initial_to_intermediate(v,redundant,TreeState::reconstructed);
291 }
292 if (finalstate==nonstandard) {
293 do_fence+=change_initial_to_intermediate(v,compressed,TreeState::reconstructed);
294 do_fence+=change_initial_to_intermediate(v,redundant,TreeState::reconstructed);
295 }
296 if (finalstate==nonstandard_with_leaves) {
297 do_fence+=change_initial_to_intermediate(v,compressed,TreeState::reconstructed);
298 do_fence+=change_initial_to_intermediate(v,nonstandard,TreeState::reconstructed);
299 do_fence+=change_initial_to_intermediate(v,redundant,TreeState::reconstructed);
300 }
301 if (finalstate==redundant) {
302 do_fence+=change_initial_to_intermediate(v,compressed,TreeState::reconstructed);
303 do_fence+=change_initial_to_intermediate(v,nonstandard,TreeState::reconstructed);
304 do_fence+=change_initial_to_intermediate(v,nonstandard_with_leaves,TreeState::reconstructed);
305 }
306 if (do_fence>0) world.gop.fence();
307
308 for (unsigned int i=0; i<v.size(); ++i) v[i].change_tree_state(finalstate,fence);
309 if (fence) world.gop.fence();
310
311 return v;
312 }
313
314 /// ensure v has the requested tree state, change the tree state of v if necessary and no fence is given
315 template<typename T, std::size_t NDIM>
317 const TreeState state, bool fence) {
318 // fast return
319 if (get_tree_state(v)==state) return true;;
320
321 // if there is a fence we can simply change the tree state, might be a no-op
322 if (fence) change_tree_state(v,state,true);
323
324 // check success, throw if not
325 bool ok=get_tree_state(v)==state;
326 if (not ok) {
327 print("ensure_tree_state_respecting_fence failed");
328 throw std::runtime_error("ensure_tree_state_respecting_fence failed");
329 }
330 return ok;
331 }
332
333 /// Truncates a vector of functions
334 template <typename T, std::size_t NDIM>
335 void truncate(World& world,
336 std::vector< Function<T,NDIM> >& v,
337 double tol=0.0,
338 bool fence=true) {
339 PROFILE_BLOCK(Vtruncate);
340
341 // truncate in compressed form only for low-dimensional functions
342 // compression is very expensive if low-rank tensor approximations are used
343 if (NDIM<4) compress(world, v);
344
345 for (auto& vv: v) {
346 vv.truncate(tol, false);
347 }
348
349 if (fence) world.gop.fence();
350 }
351
352 /// Truncates a vector of functions
353
354 /// @return the truncated vector for chaining
355 template <typename T, std::size_t NDIM>
356 std::vector< Function<T,NDIM> > truncate(std::vector< Function<T,NDIM> > v,
357 double tol=0.0, bool fence=true) {
358 if (v.size()>0) truncate(v[0].world(),v,tol,fence);
359 return v;
360 }
361
362 /// reduces the tensor rank of the coefficient tensor (if applicable)
363
364 /// @return the vector for chaining
365 template <typename T, std::size_t NDIM>
366 std::vector< Function<T,NDIM> > reduce_rank(std::vector< Function<T,NDIM> > v,
367 double thresh=0.0, bool fence=true) {
368 if (v.size()==0) return v;
369 for (auto& vv : v) vv.reduce_rank(thresh,false);
370 if (fence) v[0].world().gop.fence();
371 return v;
372 }
373
374
375 /// Applies a derivative operator to a vector of functions
376 template <typename T, std::size_t NDIM>
377 std::vector< Function<T,NDIM> >
378 apply(World& world,
379 const Derivative<T,NDIM>& D,
380 const std::vector< Function<T,NDIM> >& v,
381 bool fence=true)
382 {
383 reconstruct(world, v);
384 std::vector< Function<T,NDIM> > df(v.size());
385 for (unsigned int i=0; i<v.size(); ++i) {
386 df[i] = D(v[i],false);
387 }
388 if (fence) world.gop.fence();
389 return df;
390 }
391
392 /// Generates a vector of zero functions with a given tree state
393 template <typename T, std::size_t NDIM>
394 std::vector< Function<T,NDIM> >
395 zero_functions_tree_state(World& world, int n, const TreeState state, bool fence=true) {
396 std::vector< Function<T,NDIM> > r(n);
397 for (int i=0; i<n; ++i) {
398 if (state==compressed)
399 r[i] = Function<T,NDIM>(FunctionFactory<T,NDIM>(world).fence(false).compressed(true).initial_level(1));
400 else if (state==reconstructed)
401 r[i] = Function<T,NDIM>(FunctionFactory<T,NDIM>(world).fence(false));
402 else {
403 print("zero_functions_tree_state: unknown tree state");
404 throw std::runtime_error("zero_functions_tree_state: unknown tree state");
405 }
406 }
407
408 if (n && fence) world.gop.fence();
409 return r;
410
411 }
412
413 /// Generates a vector of zero functions (reconstructed)
414 template <typename T, std::size_t NDIM>
415 std::vector< Function<T,NDIM> >
416 zero_functions(World& world, int n, bool fence=true) {
417 return zero_functions_tree_state<T,NDIM>(world,n,reconstructed,fence);
418 }
419
420 /// Generates a vector of zero functions (compressed)
421 template <typename T, std::size_t NDIM>
422 std::vector< Function<T,NDIM> >
423 zero_functions_compressed(World& world, int n, bool fence=true) {
424 return zero_functions_tree_state<T,NDIM>(world,n,compressed,fence);
425 }
426
427 /// Generates a vector of zero functions, either compressed or reconstructed, depending on tensor type
428 template <typename T, std::size_t NDIM>
429 std::vector< Function<T,NDIM> >
430 zero_functions_auto_tree_state(World& world, int n, bool fence=true) {
432 return zero_functions_tree_state<T,NDIM>(world,n,state,fence);
433 }
434
435
436
437 /// orthonormalize the vectors
438 template<typename T, std::size_t NDIM>
439 std::vector<Function<T,NDIM>> orthonormalize(const std::vector<Function<T,NDIM> >& vf_in) {
440 if (vf_in.size()==0) return std::vector<Function<T,NDIM>>();
441 World& world=vf_in.front().world();
442 auto vf=copy(world,vf_in);
443 normalize(world,vf);
444 if (vf.size()==1) return copy(world,vf_in);
445 double maxq;
446 double trantol=0.0;
447 auto Q2=[](const Tensor<T>& s) {
448 Tensor<T> Q = -0.5*s;
449 for (int i=0; i<s.dim(0); ++i) Q(i,i) += 1.5;
450 return Q;
451 };
452
453 do {
454 Tensor<T> Q = Q2(matrix_inner(world, vf, vf));
455 maxq=0.0;
456 for (int i=0; i<Q.dim(0); ++i)
457 for (int j=0; j<i; ++j)
458 maxq = std::max(maxq,std::abs(Q(i,j)));
459
460 vf = transform(world, vf, Q, trantol, true);
461 truncate(world, vf);
462
463 } while (maxq>0.01);
464 normalize(world,vf);
465 return vf;
466 }
467
468
469 /// symmetric orthonormalization (see e.g. Szabo/Ostlund)
470
471 /// @param[in] the vector to orthonormalize
472 /// @param[in] overlap matrix
473 template <typename T, std::size_t NDIM>
474 std::vector<Function<T,NDIM> > orthonormalize_symmetric(
475 const std::vector<Function<T,NDIM> >& v,
476 const Tensor<T>& ovlp,
477 double lindep = 1e-12) {
478 if(v.empty()) return v;
479
480 World& world = v.front().world();
481 const size_t n = v.size();
482
483 Tensor<T> U;
485 syev(ovlp, U, s);
486 lindep *= s(s.size() - 1); // eigenvalues are in ascending order
487
488 // transform s to s^{-1/2} in-place
489 int rank = 0, nlindep = 0;
490 for(size_t i = 0; i < n; ++i) {
491 const auto s_i = s(i);
492 s(i) = 1.0 / sqrt(s_i);
493 (s_i > lindep) ? rank++ : nlindep++;
494 }
495 MADNESS_ASSERT(size_t(nlindep + rank) == n);
496
497 // warn of linearly dependent vectors and values
498 if (nlindep > 0) {
499 if (world.rank() == 0)
500 print("WARNING: linear dependencies detected in ", nlindep,
501 " functions, rank = ", rank);
502 }
503
504 // save Ut before U gets modified with s^{-1/2}
505 const Tensor<T> Ut = conj_transpose(U);
506
507 for(size_t i = 0; i < n; ++i){
508 for(size_t j = 0; j < n; ++j){
509 U(i, j) = U(i, j) * s(j);
510 }
511 }
512
513 Tensor<T> X = inner(U, Ut, 1, 0);
514
515 return transform(world, v, X);
516 }
517
518 /// convenience routine for symmetric orthonormalization (see e.g. Szabo/Ostlund)
519 /// overlap matrix is calculated
520 /// @param[in] the vector to orthonormalize
521 template <typename T, std::size_t NDIM>
522 std::vector<Function<T,NDIM> > orthonormalize_symmetric(const std::vector<Function<T,NDIM> >& v,
523 double lindep = 1e-12){
524 if(v.empty()) return v;
525
526 Tensor<T> ovlp = matrix_inner(v.front().world(), v, v, /* sym= */ true);
527
528 return orthonormalize_symmetric(v, ovlp, lindep);
529 }
530
531 /// canonical orthonormalization (see e.g. Szabo/Ostlund)
532 /// @param[in] the vector to orthonormalize
533 /// @param[in] overlap matrix
534 /// @param[in] lindep linear dependency threshold relative to largest eigenvalue
535 template <typename T, std::size_t NDIM>
536 std::vector<Function<T,NDIM> > orthonormalize_canonical(
537 const std::vector<Function<T,NDIM> >& v,
538 const Tensor<T>& ovlp,
539 double lindep = 1e-12) {
540 if(v.empty()) return v;
541
542 World& world = v.front().world();
543 const size_t n = v.size();
544
545 Tensor<T> U;
547 syev(ovlp, U, s);
548 lindep *= s(s.size() - 1); // eigenvalues are in ascending order
549
550 // transform s to s^{-1/2} in-place
551 int rank = 0, nlindep = 0;
552 for(size_t i = 0; i < n; ++i) {
553 const auto s_i = s(i);
554 if (s_i > lindep) {
555 s(i) = 1.0 / sqrt(s_i);
556 rank++;
557 } else {
558 nlindep++;
559 }
560 }
561 MADNESS_ASSERT(size_t(nlindep + rank) == n);
562
563 // remove linearly dependent vectors and values
564 if (nlindep > 0) {
565 if (world.rank() == 0)
566 print("Linear dependencies detected: removed ", nlindep,
567 " functions, rank = ", rank);
568 U = U(_, Slice(nlindep, -1));
569 s = s(Slice(nlindep, -1));
570 }
571
572 // modify U in-place, U is now transformation matrix (U * s^{-1/2})
573 for(size_t i = 0; i < n; ++i){
574 for(size_t j = 0; j < rank; ++j){
575 U(i, j) = U(i, j) * s(j);
576 }
577 }
578
579 return transform(world, v, U);
580 }
581
582 /// convenience routine for canonical routine for symmetric orthonormalization (see e.g. Szabo/Ostlund)
583 /// overlap matrix is calculated
584 /// @param[in] the vector to orthonormalize
585 template <typename T, std::size_t NDIM>
586 std::vector<Function<T,NDIM> > orthonormalize_canonical(const std::vector<Function<T,NDIM> >& v,
587 double lindep = 1e-12){
588 if(v.empty()) return v;
589
590 Tensor<T> ovlp = matrix_inner(v.front().world(), v, v, /* sym= */ true);
591
592 return orthonormalize_canonical(v, ovlp, lindep);
593 }
594
595 /// cholesky orthonormalization without pivoting
596 /// @param[in] the vector to orthonormalize
597 /// @param[in] overlap matrix, destroyed on return!
598 template <typename T, std::size_t NDIM>
599 std::vector<Function<T,NDIM> > orthonormalize_cd(
600 const std::vector<Function<T,NDIM> >& v,
601 Tensor<T>& ovlp) {
602
603 if (v.empty()) return v;
604
605 cholesky(ovlp); // destroys ovlp and gives back Upper ∆ Matrix from CD
606
607 Tensor<T> L = transpose(ovlp);
608 Tensor<T> Linv = inverse(L);
609 Tensor<T> U = transpose(Linv);
610
611 World& world=v.front().world();
612 return transform(world, v, U);
613
614 }
615
616 /// convenience routine for cholesky orthonormalization without pivoting
617 /// @param[in] the vector to orthonormalize
618 /// @param[in] overlap matrix
619 template <typename T, std::size_t NDIM>
620 std::vector<Function<T,NDIM> > orthonormalize_cd(const std::vector<Function<T,NDIM> >& v){
621 if(v.empty()) return v;
622
623 World& world=v.front().world();
624 Tensor<T> ovlp = matrix_inner(world, v, v, /* sym= */ true);
625
626 return orthonormalize_cd(v,ovlp);
627 }
628
629 /// @param[in] the vector to orthonormalize
630 /// @param[in] overlap matrix, will be destroyed on return!
631 /// @param[in] tolerance for numerical rank reduction
632 /// @param[out] pivoting vector, no allocation on input needed
633 /// @param[out] rank
634 /// @return orthonormalized vector (may or may not be truncated)
635 template <typename T, std::size_t NDIM>
636 std::vector<Function<T,NDIM> > orthonormalize_rrcd(
637 const std::vector<Function<T,NDIM> >& v,
638 Tensor<T>& ovlp,
639 const double tol,
640 Tensor<integer>& piv,
641 int& rank) {
642
643 if (v.empty()) {
644 return v;
645 }
646
647 rr_cholesky(ovlp,tol,piv,rank); // destroys ovlp and gives back Upper ∆ Matrix from CCD
648
649 // rearrange and truncate the functions according to the pivoting of the rr_cholesky
650 std::vector<Function<T,NDIM> > pv(rank);
651 for(integer i=0;i<rank;++i){
652 pv[i]=v[piv[i]];
653 }
654 ovlp=ovlp(Slice(0,rank-1),Slice(0,rank-1));
655
656 Tensor<T> L = transpose(ovlp);
657 Tensor<T> Linv = inverse(L);
658 Tensor<T> U = transpose(Linv);
659
660 World& world=v.front().world();
661 return transform(world, pv, U);
662 }
663
664 /// convenience routine for orthonormalize_cholesky: orthonormalize_cholesky without information on pivoting and rank
665 /// @param[in] the vector to orthonormalize
666 /// @param[in] overlap matrix
667 /// @param[in] tolerance for numerical rank reduction
668 template <typename T, std::size_t NDIM>
669 std::vector<Function<T,NDIM> > orthonormalize_rrcd(const std::vector<Function<T,NDIM> >& v, Tensor<T> ovlp , const double tol) {
670 Tensor<integer> piv;
671 int rank;
672 return orthonormalize_rrcd(v,ovlp,tol,piv,rank);
673 }
674
675 /// convenience routine for orthonormalize_cholesky: computes the overlap matrix and then calls orthonormalize_cholesky
676 /// @param[in] the vector to orthonormalize
677 /// @param[in] tolerance for numerical rank reduction
678 template <typename T, std::size_t NDIM>
679 std::vector<Function<T,NDIM> > orthonormalize_rrcd(const std::vector<Function<T,NDIM> >& v, const double tol) {
680 if (v.empty()) {
681 return v;
682 }
683 // compute overlap
684 World& world=v.front().world();
685 Tensor<T> ovlp = matrix_inner(world, v, v, /* sym= */ true);
686 return orthonormalize_rrcd(v,ovlp,tol);
687 }
688
689 /// combine two vectors
690 template <typename T, std::size_t NDIM>
691 std::vector<Function<T,NDIM> > append(const std::vector<Function<T,NDIM> > & lhs, const std::vector<Function<T,NDIM> > & rhs){
692 std::vector<Function<T,NDIM> > v=lhs;
693 for (std::size_t i = 0; i < rhs.size(); ++i) v.push_back(rhs[i]);
694 return v;
695 }
696
697 template <typename T, std::size_t NDIM>
698 std::vector<Function<T,NDIM> > flatten(const std::vector< std::vector<Function<T,NDIM> > >& vv){
699 std::vector<Function<T,NDIM> >result;
700 for(const auto& x:vv) result=append(result,x);
701 return result;
702 }
703
704 template<typename T, std::size_t NDIM>
705 std::vector<std::shared_ptr<FunctionImpl<T,NDIM>>> get_impl(const std::vector<Function<T,NDIM>>& v) {
706 std::vector<std::shared_ptr<FunctionImpl<T,NDIM>>> result;
707 for (auto& f : v) result.push_back(f.get_impl());
708 return result;
709 }
710
711 template<typename T, std::size_t NDIM>
712 void set_impl(std::vector<Function<T,NDIM>>& v, const std::vector<std::shared_ptr<FunctionImpl<T,NDIM>>> vimpl) {
713 MADNESS_CHECK(vimpl.size()==v.size());
714 for (std::size_t i=0; i<vimpl.size(); ++i) v[i].set_impl(vimpl[i]);
715 }
716
717 template<typename T, std::size_t NDIM>
718 std::vector<Function<T,NDIM>> impl2function(const std::vector<std::shared_ptr<FunctionImpl<T,NDIM>>> vimpl) {
719 std::vector<Function<T,NDIM>> v(vimpl.size());
720 for (std::size_t i=0; i<vimpl.size(); ++i) v[i].set_impl(vimpl[i]);
721 return v;
722 }
723
724
725 /// Transforms a vector of functions according to new[i] = sum[j] old[j]*c[j,i]
726
727 /// Uses sparsity in the transformation matrix --- set small elements to
728 /// zero to take advantage of this.
729 template <typename T, typename R, std::size_t NDIM>
730 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> >
732 const std::vector< Function<T,NDIM> >& v,
733 const Tensor<R>& c,
734 bool fence=true) {
735
736 PROFILE_BLOCK(Vtransformsp);
737 typedef TENSOR_RESULT_TYPE(T,R) resultT;
738 int n = v.size(); // n is the old dimension
739 int m = c.dim(1); // m is the new dimension
740 MADNESS_CHECK(n==c.dim(0));
741
742 std::vector< Function<resultT,NDIM> > vc = zero_functions_compressed<resultT,NDIM>(world, m);
743 compress(world, v);
744
745 for (int i=0; i<m; ++i) {
746 for (int j=0; j<n; ++j) {
747 if (c(j,i) != R(0.0)) vc[i].gaxpy(resultT(1.0),v[j],resultT(c(j,i)),false);
748 }
749 }
750
751 if (fence) world.gop.fence();
752 return vc;
753 }
754
755 /// Transforms a vector of functions according to new[i] = sum[j] old[j]*c[j,i]
756
757 /// all trees are in reconstructed state, final trees have to be summed down if no fence is present
758 template <typename T, typename R, std::size_t NDIM>
759 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> >
761 const std::vector< Function<T,NDIM> >& v,
762 const Tensor<R>& c,
763 bool fence=true) {
764
765 PROFILE_BLOCK(Vtransformsp);
766 typedef TENSOR_RESULT_TYPE(T,R) resultT;
767 int n = v.size(); // n is the old dimension
768 int m = c.dim(1); // m is the new dimension
769 MADNESS_CHECK(n==c.dim(0));
770
771 // if we fence set the right tree state here, otherwise it has to be correct from the start.
773 for (const auto& vv : v) MADNESS_CHECK_THROW(
774 vv.get_impl()->get_tree_state()==reconstructed,"trees have to be reconstructed in transform_reconstructed");
775
776 std::vector< Function<resultT,NDIM> > result = zero_functions<resultT,NDIM>(world, m);
777
778 for (int i=0; i<m; ++i) {
779 result[i].get_impl()->set_tree_state(redundant_after_merge);
780 for (int j=0; j<n; ++j) {
781 if (c(j,i) != R(0.0)) v[j].get_impl()->accumulate_trees(*(result[i].get_impl()),resultT(c(j,i)),true);
782 }
783 }
784
785 // if we fence we can as well finish the job here. Otherwise no harm done, as the tree state is well-defined.
786 if (fence) {
787 world.gop.fence();
788 // for (auto& r : vc) r.sum_down(false);
789 for (auto& r : result) r.get_impl()->finalize_sum();
790 world.gop.fence();
791 }
792 return result;
793 }
794
795 /// this version of transform uses Function::vtransform and screens
796 /// using both elements of `c` and `v`
797 template <typename L, typename R, std::size_t NDIM>
798 std::vector< Function<TENSOR_RESULT_TYPE(L,R),NDIM> >
799 transform(World& world, const std::vector< Function<L,NDIM> >& v,
800 const Tensor<R>& c, double tol, bool fence) {
801 PROFILE_BLOCK(Vtransform);
802 MADNESS_ASSERT(v.size() == (unsigned int)(c.dim(0)));
803
804 std::vector< Function<TENSOR_RESULT_TYPE(L,R),NDIM> > vresult
805 = zero_functions_compressed<TENSOR_RESULT_TYPE(L,R),NDIM>(world, c.dim(1));
806
807 compress(world, v, true);
808 vresult[0].vtransform(v, c, vresult, tol, fence);
809 return vresult;
810 }
811
812 template <typename T, typename R, std::size_t NDIM>
813 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> >
815 const std::vector< Function<T,NDIM> >& v,
816 const DistributedMatrix<R>& c,
817 bool fence=true) {
819
820 typedef TENSOR_RESULT_TYPE(T,R) resultT;
821 long n = v.size(); // n is the old dimension
822 long m = c.rowdim(); // m is the new dimension
823 MADNESS_ASSERT(n==c.coldim());
824
825 // new(i) = sum(j) old(j) c(j,i)
826
827 Tensor<T> tmp(n,m);
828 c.copy_to_replicated(tmp); // for debugging
829 tmp = transpose(tmp);
830
831 std::vector< Function<resultT,NDIM> > vc = zero_functions_compressed<resultT,NDIM>(world, m);
832 compress(world, v);
833
834 for (int i=0; i<m; ++i) {
835 for (int j=0; j<n; ++j) {
836 if (tmp(j,i) != R(0.0)) vc[i].gaxpy(1.0,v[j],tmp(j,i),false);
837 }
838 }
839
840 if (fence) world.gop.fence();
841 return vc;
842 }
843
844
845 /// Scales inplace a vector of functions by distinct values
846 template <typename T, typename Q, std::size_t NDIM>
847 void scale(World& world,
848 std::vector< Function<T,NDIM> >& v,
849 const std::vector<Q>& factors,
850 bool fence=true) {
851 PROFILE_BLOCK(Vscale);
852 for (unsigned int i=0; i<v.size(); ++i) v[i].scale(factors[i],false);
853 if (fence) world.gop.fence();
854 }
855
856 /// Scales inplace a vector of functions by the same
857 template <typename T, typename Q, std::size_t NDIM>
858 void scale(World& world,
859 std::vector< Function<T,NDIM> >& v,
860 const Q factor,
861 bool fence=true) {
862 PROFILE_BLOCK(Vscale);
863 for (unsigned int i=0; i<v.size(); ++i) v[i].scale(factor,false);
864 if (fence) world.gop.fence();
865 }
866
867 /// Computes the 2-norms of a vector of functions
868 template <typename T, std::size_t NDIM>
869 std::vector<double> norm2s(World& world,
870 const std::vector< Function<T,NDIM> >& v) {
871 PROFILE_BLOCK(Vnorm2);
872 std::vector<double> norms(v.size());
874 for (unsigned int i=0; i<v.size(); ++i) norms[i] = v[i].norm2sq_local();
875 world.gop.sum(&norms[0], norms.size());
876 for (unsigned int i=0; i<v.size(); ++i) norms[i] = sqrt(norms[i]);
877 world.gop.fence();
878 return norms;
879 }
880 /// Computes the 2-norms of a vector of functions
881 template <typename T, std::size_t NDIM>
882 Tensor<double> norm2s_T(World& world, const std::vector<Function<T, NDIM>>& v) {
883 PROFILE_BLOCK(Vnorm2);
884 Tensor<double> norms(v.size());
886 for (unsigned int i = 0; i < v.size(); ++i) norms[i] = v[i].norm2sq_local();
887 world.gop.sum(&norms[0], norms.size());
888 for (unsigned int i = 0; i < v.size(); ++i) norms[i] = sqrt(norms[i]);
889 world.gop.fence();
890 return norms;
891 }
892
893 /// Computes the 2-norm of a vector of functions
894 template <typename T, std::size_t NDIM>
895 double norm2(World& world,const std::vector< Function<T,NDIM> >& v) {
896 PROFILE_BLOCK(Vnorm2);
897 if (v.size()==0) return 0.0;
899 std::vector<double> norms(v.size());
900 for (unsigned int i=0; i<v.size(); ++i) norms[i] = v[i].norm2sq_local();
901 world.gop.sum(&norms[0], norms.size());
902 for (unsigned int i=1; i<v.size(); ++i) norms[0] += norms[i];
903 world.gop.fence();
904 return sqrt(norms[0]);
905 }
906
907 inline double conj(double x) {
908 return x;
909 }
910
911 inline double conj(float x) {
912 return x;
913 }
914
915// !!! FIXME: this task is broken because FunctionImpl::inner_local forces a
916// future on return from WorldTaskQueue::reduce, which will causes a deadlock if
917// run inside a task. This behavior must be changed before this task can be used
918// again.
919//
920// template <typename T, typename R, std::size_t NDIM>
921// struct MatrixInnerTask : public TaskInterface {
922// Tensor<TENSOR_RESULT_TYPE(T,R)> result; // Must be a copy
923// const Function<T,NDIM>& f;
924// const std::vector< Function<R,NDIM> >& g;
925// long jtop;
926//
927// MatrixInnerTask(const Tensor<TENSOR_RESULT_TYPE(T,R)>& result,
928// const Function<T,NDIM>& f,
929// const std::vector< Function<R,NDIM> >& g,
930// long jtop)
931// : result(result), f(f), g(g), jtop(jtop) {}
932//
933// void run(World& world) {
934// for (long j=0; j<jtop; ++j) {
935// result(j) = f.inner_local(g[j]);
936// }
937// }
938//
939// private:
940// /// Get the task id
941//
942// /// \param id The id to set for this task
943// virtual void get_id(std::pair<void*,unsigned short>& id) const {
944// PoolTaskInterface::make_id(id, *this);
945// }
946// }; // struct MatrixInnerTask
947
948
949
950 template <typename T, std::size_t NDIM>
952 const std::vector< Function<T,NDIM> >& f,
953 const std::vector< Function<T,NDIM> >& g,
954 bool sym=false)
955 {
958 const int64_t n = A.coldim();
959 const int64_t m = A.rowdim();
960 MADNESS_ASSERT(int64_t(f.size()) == n && int64_t(g.size()) == m);
961
962 // Assume we can always create an ichunk*jchunk matrix locally
963 const int ichunk = 1000;
964 const int jchunk = 1000; // 1000*1000*8 = 8 MBytes
965 for (int64_t ilo=0; ilo<n; ilo+=ichunk) {
966 int64_t ihi = std::min(ilo + ichunk, n);
967 std::vector< Function<T,NDIM> > ivec(f.begin()+ilo, f.begin()+ihi);
968 for (int64_t jlo=0; jlo<m; jlo+=jchunk) {
969 int64_t jhi = std::min(jlo + jchunk, m);
970 std::vector< Function<T,NDIM> > jvec(g.begin()+jlo, g.begin()+jhi);
971
972 Tensor<T> P = matrix_inner(A.get_world(), ivec, jvec);
973 A.copy_from_replicated_patch(ilo, ihi - 1, jlo, jhi - 1, P);
974 }
975 }
976 return A;
977 }
978
979 /// Computes the matrix inner product of two function vectors - q(i,j) = inner(f[i],g[j])
980
981 /// For complex types symmetric is interpreted as Hermitian.
982
983 /// The current parallel loop is non-optimal but functional.
984 template <typename T, typename R, std::size_t NDIM>
986 const std::vector< Function<T,NDIM> >& f,
987 const std::vector< Function<R,NDIM> >& g,
988 bool sym=false)
989 {
990 world.gop.fence();
991 auto tensor_type = [](const std::vector<Function<T,NDIM>>& v) {
992 return v.front().get_impl()->get_tensor_type();
993 };
994 TreeState operating_state=tensor_type(f)==TT_FULL ? compressed : redundant;
995 ensure_tree_state_respecting_fence(f,operating_state,true);
996 ensure_tree_state_respecting_fence(g,operating_state,true);
997
998 std::vector<const FunctionImpl<T,NDIM>*> left(f.size());
999 std::vector<const FunctionImpl<R,NDIM>*> right(g.size());
1000 for (unsigned int i=0; i<f.size(); i++) left[i] = f[i].get_impl().get();
1001 for (unsigned int i=0; i<g.size(); i++) right[i]= g[i].get_impl().get();
1002
1004
1005 world.gop.fence();
1006 world.gop.sum(r.ptr(),f.size()*g.size());
1007
1008 return r;
1009 }
1010
1011 /// Computes the matrix inner product of two function vectors - q(i,j) = inner(f[i],g[j])
1012
1013 /// For complex types symmetric is interpreted as Hermitian.
1014 ///
1015 /// The current parallel loop is non-optimal but functional.
1016 template <typename T, typename R, std::size_t NDIM>
1018 const std::vector< Function<T,NDIM> >& f,
1019 const std::vector< Function<R,NDIM> >& g,
1020 bool sym=false) {
1021 PROFILE_BLOCK(Vmatrix_inner);
1022 long n=f.size(), m=g.size();
1023 Tensor< TENSOR_RESULT_TYPE(T,R) > r(n,m);
1024 if (sym) MADNESS_ASSERT(n==m);
1025
1026 world.gop.fence();
1027 compress(world, f);
1028 if ((void*)(&f) != (void*)(&g)) compress(world, g);
1029
1030 for (long i=0; i<n; ++i) {
1031 long jtop = m;
1032 if (sym) jtop = i+1;
1033 for (long j=0; j<jtop; ++j) {
1034 r(i,j) = f[i].inner_local(g[j]);
1035 if (sym) r(j,i) = conj(r(i,j));
1036 }
1037 }
1038
1039// for (long i=n-1; i>=0; --i) {
1040// long jtop = m;
1041// if (sym) jtop = i+1;
1042// world.taskq.add(new MatrixInnerTask<T,R,NDIM>(r(i,_), f[i], g, jtop));
1043// }
1044 world.gop.fence();
1045 world.gop.sum(r.ptr(),n*m);
1046
1047// if (sym) {
1048// for (int i=0; i<n; ++i) {
1049// for (int j=0; j<i; ++j) {
1050// r(j,i) = conj(r(i,j));
1051// }
1052// }
1053// }
1054 return r;
1055 }
1056
1057 /// Computes the element-wise inner product of two function vectors - q(i) = inner(f[i],g[i])
1058
1059 /// works in reconstructed or compressed state, state is chosen based on TensorType
1060 template <typename T, typename R, std::size_t NDIM>
1062 const std::vector< Function<T,NDIM> >& f,
1063 const std::vector< Function<R,NDIM> >& g) {
1064 PROFILE_BLOCK(Vinnervv);
1065 long n=f.size(), m=g.size();
1066 MADNESS_CHECK(n==m);
1067 Tensor< TENSOR_RESULT_TYPE(T,R) > r(n);
1068
1069 auto tensor_type = [](const std::vector<Function<T,NDIM>>& v) {
1070 return v.front().get_impl()->get_tensor_type();
1071 };
1072 TreeState operating_state=tensor_type(f)==TT_FULL ? compressed : redundant;
1073 ensure_tree_state_respecting_fence(f,operating_state,true);
1074 ensure_tree_state_respecting_fence(g,operating_state,true);
1075
1076 for (long i=0; i<n; ++i) r(i) = f[i].inner_local(g[i]);
1077
1078 world.taskq.fence();
1079 world.gop.sum(r.ptr(),n);
1080 world.gop.fence();
1081 return r;
1082 }
1083
1084
1085 /// Computes the inner product of a function with a function vector - q(i) = inner(f,g[i])
1086
1087 /// works in reconstructed or compressed state, state is chosen based on TensorType
1088 template <typename T, typename R, std::size_t NDIM>
1090 const Function<T,NDIM>& f,
1091 const std::vector< Function<R,NDIM> >& g) {
1092 PROFILE_BLOCK(Vinner);
1093 long n=g.size();
1094 Tensor< TENSOR_RESULT_TYPE(T,R) > r(n);
1095
1096 auto tensor_type = [](const std::vector<Function<T,NDIM>>& v) {
1097 return v.front().get_impl()->get_tensor_type();
1098 };
1099 TreeState operating_state=tensor_type(g)==TT_FULL ? compressed : redundant;
1100 f.change_tree_state(operating_state,false);
1101 ensure_tree_state_respecting_fence(g,operating_state,true);
1102 world.gop.fence();
1103
1104 for (long i=0; i<n; ++i) {
1105 r(i) = f.inner_local(g[i]);
1106 }
1107
1108 world.taskq.fence();
1109 world.gop.sum(r.ptr(),n);
1110 world.gop.fence();
1111 return r;
1112 }
1113
1114 /// inner function with right signature for the nonlinear solver
1115 /// this is needed for the KAIN solvers and other functions
1116 template <typename T, typename R, std::size_t NDIM>
1118 const std::vector< Function<R,NDIM> >& g){
1119 MADNESS_ASSERT(f.size()==g.size());
1120 if(f.empty()) return 0.0;
1121 else return inner(f[0].world(),f,g).sum();
1122 }
1123
1124
1125 /// Multiplies a function against a vector of functions --- q[i] = a * v[i]
1126 template <typename T, typename R, std::size_t NDIM>
1127 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1128 mul(World& world,
1129 const Function<T,NDIM>& a,
1130 const std::vector< Function<R,NDIM> >& v,
1131 bool fence=true) {
1132 PROFILE_BLOCK(Vmul);
1133 a.reconstruct(false);
1134 reconstruct(world, v, false);
1135 world.gop.fence();
1136 return vmulXX(a, v, 0.0, fence);
1137 }
1138
1139 /// Multiplies a function against a vector of functions using sparsity of a and v[i] --- q[i] = a * v[i]
1140 template <typename T, typename R, std::size_t NDIM>
1141 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1143 const Function<T,NDIM>& a,
1144 const std::vector< Function<R,NDIM> >& v,
1145 double tol,
1146 bool fence=true) {
1147 PROFILE_BLOCK(Vmulsp);
1148 a.reconstruct(false);
1149 reconstruct(world, v, false);
1150 world.gop.fence();
1151 for (unsigned int i=0; i<v.size(); ++i) {
1152 v[i].norm_tree(false);
1153 }
1154 a.norm_tree();
1155 return vmulXX(a, v, tol, fence);
1156 }
1157
1158
1159 /// Outer product of a vector of functions with a vector of functions using sparsity
1160
1161 /// \tparam T type parameter for first factor
1162 /// \tparam R type parameter for second factor
1163 /// \tparam NDIM dimension of first and second factors
1164 /// \param world the world
1165 /// \param f first vector of functions
1166 /// \param g second vector of functions
1167 /// \param tol threshold for multiplication
1168 /// \param fence force fence (will always fence if necessary)
1169 /// \param symm if true, only compute f(i) * g(j) for j<=i
1170 /// \return fg(i,j) = f(i) * g(j), as a vector of vectors
1171 template <typename T, typename R, std::size_t NDIM>
1172 std::vector<std::vector<Function<TENSOR_RESULT_TYPE(T, R), NDIM> > >
1174 const std::vector<Function<R, NDIM> > &f,
1175 const std::vector<Function<R, NDIM> > &g,
1176 double tol,
1177 bool fence = true,
1178 bool symm = false) {
1179 PROFILE_BLOCK(Vmulsp);
1180 bool same=(&f == &g);
1181 reconstruct(world, f, false);
1182 if (not same) reconstruct(world, g, false);
1183 world.gop.fence();
1184 for (auto& ff : f) ff.norm_tree(false);
1185 if (not same) for (auto& gg : g) gg.norm_tree(false);
1186 world.gop.fence();
1187
1188 std::vector<std::vector<Function<R,NDIM> > >result(f.size());
1189 std::vector<Function<R,NDIM>> g_i;
1190 for (int64_t i=f.size()-1; i>=0; --i) {
1191 if (!symm)
1192 result[i]= vmulXX(f[i], g, tol, false);
1193 else {
1194 if (g_i.empty()) g_i = g;
1195 g_i.resize(i+1); // this shrinks g_i down to single function for i=0
1196 result[i]= vmulXX(f[i], g_i, tol, false);
1197 }
1198 }
1199 if (fence) world.gop.fence();
1200 return result;
1201 }
1202
1203 /// Makes the norm tree for all functions in a vector
1204 template <typename T, std::size_t NDIM>
1205 void norm_tree(World& world,
1206 const std::vector< Function<T,NDIM> >& v,
1207 bool fence=true)
1208 {
1209 PROFILE_BLOCK(Vnorm_tree);
1210 for (unsigned int i=0; i<v.size(); ++i) {
1211 v[i].norm_tree(false);
1212 }
1213 if (fence) world.gop.fence();
1214 }
1215
1216 /// Multiplies two vectors of functions q[i] = a[i] * b[i]
1217 template <typename T, typename R, std::size_t NDIM>
1218 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1219 mul(World& world,
1220 const std::vector< Function<T,NDIM> >& a,
1221 const std::vector< Function<R,NDIM> >& b,
1222 bool fence=true) {
1223 PROFILE_BLOCK(Vmulvv);
1224 reconstruct(world, a, true);
1225 reconstruct(world, b, true);
1226// if (&a != &b) reconstruct(world, b, true); // fails if type(a) != type(b)
1227
1228 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> > q(a.size());
1229 for (unsigned int i=0; i<a.size(); ++i) {
1230 q[i] = mul(a[i], b[i], false);
1231 }
1232 if (fence) world.gop.fence();
1233 return q;
1234 }
1235
1236
1237 /// multiply a high-dimensional function with a low-dimensional function
1238
1239 /// @param[in] f NDIM function of NDIM dimensions
1240 /// @param[in] g LDIM function of LDIM
1241 /// @param[in] v dimension indices of f to multiply
1242 /// @return h[i](0,1,2,3) = f(0,1,2,3) * g[i](1,2,3) for v={1,2,3}
1243 template<typename T, std::size_t NDIM, std::size_t LDIM>
1244 std::vector<Function<T,NDIM> > partial_mul(const Function<T,NDIM> f, const std::vector<Function<T,LDIM> > g,
1245 const int particle) {
1246
1247 World& world=f.world();
1248 std::vector<Function<T,NDIM> > result(g.size());
1249 for (auto& r : result) r.set_impl(f, false);
1250
1251 FunctionImpl<T,NDIM>* fimpl=f.get_impl().get();
1252// fimpl->make_redundant(false);
1253 fimpl->change_tree_state(redundant,false);
1254 make_redundant(world,g,false);
1255 world.gop.fence();
1256
1257 for (std::size_t i=0; i<result.size(); ++i) {
1258 FunctionImpl<T,LDIM>* gimpl=g[i].get_impl().get();
1259 result[i].get_impl()->multiply(fimpl,gimpl,particle); // stupid naming inconsistency
1260 }
1261 world.gop.fence();
1262
1263 fimpl->undo_redundant(false);
1264 for (auto& ig : g) ig.get_impl()->undo_redundant(false);
1265 world.gop.fence();
1266 return result;
1267 }
1268
1269 template<typename T, std::size_t NDIM, std::size_t LDIM>
1270 std::vector<Function<T,NDIM> > multiply(const Function<T,NDIM> f, const std::vector<Function<T,LDIM> > g,
1271 const std::tuple<int,int,int> v) {
1272 return partial_mul<T,NDIM,LDIM>(f,g,std::array<int,3>({std::get<0>(v),std::get<1>(v),std::get<2>(v)}));
1273 }
1274
1275
1276/// Computes the square of a vector of functions --- q[i] = v[i]**2
1277 template <typename T, std::size_t NDIM>
1278 std::vector< Function<T,NDIM> >
1280 const std::vector< Function<T,NDIM> >& v,
1281 bool fence=true) {
1282 return mul<T,T,NDIM>(world, v, v, fence);
1283// std::vector< Function<T,NDIM> > vsq(v.size());
1284// for (unsigned int i=0; i<v.size(); ++i) {
1285// vsq[i] = square(v[i], false);
1286// }
1287// if (fence) world.gop.fence();
1288// return vsq;
1289 }
1290
1291
1292 /// Computes the square of a vector of functions --- q[i] = abs(v[i])**2
1293 template <typename T, std::size_t NDIM>
1294 std::vector< Function<typename Tensor<T>::scalar_type,NDIM> >
1295 abssq(World& world,
1296 const std::vector< Function<T,NDIM> >& v,
1297 bool fence=true) {
1298 typedef typename Tensor<T>::scalar_type scalartype;
1299 reconstruct(world,v);
1300 std::vector<Function<scalartype,NDIM> > result(v.size());
1301 for (size_t i=0; i<v.size(); ++i) result[i]=abs_square(v[i],false);
1302 if (fence) world.gop.fence();
1303 return result;
1304 }
1305
1306
1307 /// Sets the threshold in a vector of functions
1308 template <typename T, std::size_t NDIM>
1309 void set_thresh(World& world, std::vector< Function<T,NDIM> >& v, double thresh, bool fence=true) {
1310 for (unsigned int j=0; j<v.size(); ++j) {
1311 v[j].set_thresh(thresh,false);
1312 }
1313 if (fence) world.gop.fence();
1314 }
1315
1316 /// Returns the complex conjugate of the vector of functions
1317 template <typename T, std::size_t NDIM>
1318 std::vector< Function<T,NDIM> >
1319 conj(World& world,
1320 const std::vector< Function<T,NDIM> >& v,
1321 bool fence=true) {
1322 PROFILE_BLOCK(Vconj);
1323 std::vector< Function<T,NDIM> > r = copy(world, v); // Currently don't have oop conj
1324 for (unsigned int i=0; i<v.size(); ++i) {
1325 r[i].conj(false);
1326 }
1327 if (fence) world.gop.fence();
1328 return r;
1329 }
1330
1331 /// Returns a deep copy of a vector of functions
1332 template <typename T, typename R, std::size_t NDIM>
1333 std::vector< Function<R,NDIM> > convert(World& world,
1334 const std::vector< Function<T,NDIM> >& v, bool fence=true) {
1335 PROFILE_BLOCK(Vcopy);
1336 std::vector< Function<R,NDIM> > r(v.size());
1337 for (unsigned int i=0; i<v.size(); ++i) {
1338 r[i] = convert<T,R,NDIM>(v[i], false);
1339 }
1340 if (fence) world.gop.fence();
1341 return r;
1342 }
1343
1344
1345 /// Returns a deep copy of a vector of functions
1346 template <typename T, std::size_t NDIM>
1347 std::vector< Function<T,NDIM> >
1348 copy(World& world,
1349 const std::vector< Function<T,NDIM> >& v,
1350 bool fence=true) {
1351 PROFILE_BLOCK(Vcopy);
1352 std::vector< Function<T,NDIM> > r(v.size());
1353 for (unsigned int i=0; i<v.size(); ++i) {
1354 r[i] = copy(v[i], false);
1355 }
1356 if (fence) world.gop.fence();
1357 return r;
1358 }
1359
1360
1361 /// Returns a deep copy of a vector of functions
1362 template <typename T, std::size_t NDIM>
1363 std::vector< Function<T,NDIM> >
1364 copy(const std::vector< Function<T,NDIM> >& v, bool fence=true) {
1365 PROFILE_BLOCK(Vcopy);
1366 std::vector< Function<T,NDIM> > r(v.size());
1367 if (v.size()>0) r=copy(v.front().world(),v,fence);
1368 return r;
1369 }
1370
1371 /// Returns a vector of `n` deep copies of a function
1372 template <typename T, std::size_t NDIM>
1373 std::vector< Function<T,NDIM> >
1375 const Function<T,NDIM>& v,
1376 const unsigned int n,
1377 bool fence=true) {
1378 PROFILE_BLOCK(Vcopy1);
1379 std::vector< Function<T,NDIM> > r(n);
1380 for (unsigned int i=0; i<n; ++i) {
1381 r[i] = copy(v, false);
1382 }
1383 if (fence) world.gop.fence();
1384 return r;
1385 }
1386
1387 /// Create a new copy of the function with different distribution and optional
1388 /// fence
1389
1390 /// Works in either basis. Different distributions imply
1391 /// asynchronous communication and the optional fence is
1392 /// collective.
1393 //
1394 /// Returns a deep copy of a vector of functions
1395
1396 template <typename T, std::size_t NDIM>
1397 std::vector<Function<T, NDIM>> copy(World& world,
1398 const std::vector<Function<T, NDIM>>& v,
1399 const std::shared_ptr<WorldDCPmapInterface<Key<NDIM>>>& pmap,
1400 bool fence = true) {
1401 PROFILE_BLOCK(Vcopy);
1402 std::vector<Function<T, NDIM>> r(v.size());
1403 for (unsigned int i = 0; i < v.size(); ++i) {
1404 r[i] = copy(v[i], pmap, false);
1405 }
1406 if (fence) world.gop.fence();
1407 return r;
1408 }
1409
1410 /// Returns new vector of functions --- q[i] = a[i] + b[i]
1411 template <typename T, typename R, std::size_t NDIM>
1412 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1413 add(World& world,
1414 const std::vector< Function<T,NDIM> >& a,
1415 const std::vector< Function<R,NDIM> >& b,
1416 bool fence=true) {
1417 PROFILE_BLOCK(Vadd);
1418 MADNESS_ASSERT(a.size() == b.size());
1419 compress(world, a);
1420 compress(world, b);
1421
1422 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> > r(a.size());
1423 for (unsigned int i=0; i<a.size(); ++i) {
1424 r[i] = add(a[i], b[i], false);
1425 }
1426 if (fence) world.gop.fence();
1427 return r;
1428 }
1429
1430 /// Returns new vector of functions --- q[i] = a + b[i]
1431 template <typename T, typename R, std::size_t NDIM>
1432 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1433 add(World& world,
1434 const Function<T,NDIM> & a,
1435 const std::vector< Function<R,NDIM> >& b,
1436 bool fence=true) {
1437 PROFILE_BLOCK(Vadd1);
1438 a.compress();
1439 compress(world, b);
1440
1441 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> > r(b.size());
1442 for (unsigned int i=0; i<b.size(); ++i) {
1443 r[i] = add(a, b[i], false);
1444 }
1445 if (fence) world.gop.fence();
1446 return r;
1447 }
1448 template <typename T, typename R, std::size_t NDIM>
1449 inline std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1450 add(World& world,
1451 const std::vector< Function<R,NDIM> >& b,
1452 const Function<T,NDIM> & a,
1453 bool fence=true) {
1454 return add(world, a, b, fence);
1455 }
1456
1457 /// Returns new vector of functions --- q[i] = a[i] - b[i]
1458 template <typename T, typename R, std::size_t NDIM>
1459 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1460 sub(World& world,
1461 const std::vector< Function<T,NDIM> >& a,
1462 const std::vector< Function<R,NDIM> >& b,
1463 bool fence=true) {
1464 PROFILE_BLOCK(Vsub);
1465 MADNESS_ASSERT(a.size() == b.size());
1466 compress(world, a);
1467 compress(world, b);
1468
1469 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> > r(a.size());
1470 for (unsigned int i=0; i<a.size(); ++i) {
1471 r[i] = sub(a[i], b[i], false);
1472 }
1473 if (fence) world.gop.fence();
1474 return r;
1475 }
1476
1477 /// Returns new function --- q = sum_i f[i]
1478 template <typename T, std::size_t NDIM>
1479 Function<T, NDIM> sum(World& world, const std::vector<Function<T,NDIM> >& f,
1480 bool fence=true) {
1481
1482 compress(world, f);
1484
1485 for (unsigned int i=0; i<f.size(); ++i) r.gaxpy(1.0,f[i],1.0,false);
1486 if (fence) world.gop.fence();
1487 return r;
1488 }
1489
1490 template <typename T, std::size_t NDIM>
1492 const std::vector<Function<T, NDIM>>& f,
1493 const std::vector<Function<T, NDIM>>& g,
1494 bool sym=false)
1495 {
1498 const int64_t n = A.coldim();
1499 const int64_t m = A.rowdim();
1500 MADNESS_ASSERT(int64_t(f.size()) == n && int64_t(g.size()) == m);
1501
1502 // Assume we can always create an ichunk*jchunk matrix locally
1503 const int ichunk = 1000;
1504 const int jchunk = 1000; // 1000*1000*8 = 8 MBytes
1505 for (int64_t ilo = 0; ilo < n; ilo += ichunk) {
1506 int64_t ihi = std::min(ilo + ichunk, n);
1507 std::vector<Function<T, NDIM>> ivec(f.begin() + ilo, f.begin() + ihi);
1508 for (int64_t jlo = 0; jlo < m; jlo += jchunk) {
1509 int64_t jhi = std::min(jlo + jchunk, m);
1510 std::vector<Function<T, NDIM>> jvec(g.begin() + jlo, g.begin() + jhi);
1511
1512 Tensor<T> P = matrix_dot(A.get_world(), ivec, jvec, sym);
1513 A.copy_from_replicated_patch(ilo, ihi - 1, jlo, jhi - 1, P);
1514 }
1515 }
1516 return A;
1517 }
1518
1519 /// Computes the matrix dot product of two function vectors - q(i,j) = dot(f[i],g[j])
1520
1521 /// For complex types symmetric is interpreted as Hermitian.
1522 ///
1523 /// The current parallel loop is non-optimal but functional.
1524 template <typename T, typename R, std::size_t NDIM>
1526 const std::vector<Function<T, NDIM>>& f,
1527 const std::vector<Function<R, NDIM>>& g,
1528 bool sym=false)
1529 {
1530 world.gop.fence();
1531 compress(world, f);
1532 // if ((void*)(&f) != (void*)(&g)) compress(world, g);
1533 compress(world, g);
1534
1535 std::vector<const FunctionImpl<T, NDIM>*> left(f.size());
1536 std::vector<const FunctionImpl<R, NDIM>*> right(g.size());
1537 for (unsigned int i = 0; i < f.size(); i++) left[i] = f[i].get_impl().get();
1538 for (unsigned int i = 0; i < g.size(); i++) right[i] = g[i].get_impl().get();
1539
1541
1542 world.gop.fence();
1543 world.gop.sum(r.ptr(), f.size() * g.size());
1544
1545 return r;
1546 }
1547
1548 /// Computes the matrix dot product of two function vectors - q(i,j) = dot(f[i],g[j])
1549
1550 /// For complex types symmetric is interpreted as Hermitian.
1551 ///
1552 /// The current parallel loop is non-optimal but functional.
1553 template <typename T, typename R, std::size_t NDIM>
1555 const std::vector< Function<T,NDIM> >& f,
1556 const std::vector< Function<R,NDIM> >& g,
1557 bool sym=false) {
1558 PROFILE_BLOCK(Vmatrix_dot);
1559 long n=f.size(), m=g.size();
1560 Tensor< TENSOR_RESULT_TYPE(T,R) > r(n,m);
1561 if (sym) MADNESS_ASSERT(n==m);
1562
1563 world.gop.fence();
1564 compress(world, f);
1565 if ((void*)(&f) != (void*)(&g)) compress(world, g);
1566
1567 for (long i=0; i<n; ++i) {
1568 long jtop = m;
1569 if (sym) jtop = i+1;
1570 for (long j=0; j<jtop; ++j) {
1571 if (sym) {
1572 r(j,i) = f[i].dot_local(g[j]);
1573 if (i != j)
1574 r(i,j) = conj(r(j,i));
1575 } else
1576 r(i,j) = f[i].dot_local(g[j]);
1577 }
1578 }
1579
1580 world.gop.fence();
1581 world.gop.sum(r.ptr(),n*m);
1582
1583 return r;
1584 }
1585
1586 /// Multiplies and sums two vectors of functions r = \sum_i a[i] * b[i]
1587 template <typename T, typename R, std::size_t NDIM>
1588 Function<TENSOR_RESULT_TYPE(T,R), NDIM>
1589 dot(World& world,
1590 const std::vector< Function<T,NDIM> >& a,
1591 const std::vector< Function<R,NDIM> >& b,
1592 bool fence=true) {
1593 MADNESS_CHECK(a.size()==b.size());
1594 return sum(world,mul(world,a,b,true),fence);
1595 }
1596
1597
1598
1599 /// out-of-place gaxpy for two vectors: result[i] = alpha * a[i] + beta * b[i]
1600 template <typename T, typename Q, typename R, std::size_t NDIM>
1601 std::vector<Function<TENSOR_RESULT_TYPE(Q,TENSOR_RESULT_TYPE(T,R)),NDIM> >
1603 const std::vector< Function<T,NDIM> >& a,
1604 Q beta,
1605 const std::vector< Function<R,NDIM> >& b,
1606 bool fence=true) {
1607
1608 MADNESS_ASSERT(a.size() == b.size());
1609 typedef TENSOR_RESULT_TYPE(Q,TENSOR_RESULT_TYPE(T,R)) resultT;
1610 if (a.size()==0) return std::vector<Function<resultT,NDIM> >();
1611
1612 auto tensor_type = [](const std::vector<Function<T,NDIM>>& v) {
1613 return v.front().get_impl()->get_tensor_type();
1614 };
1615
1616 // gaxpy can be done either in reconstructed or in compressed state
1617 World& world=a[0].world();
1618 std::vector<Function<resultT,NDIM> > result(a.size());
1619
1620 TreeState operating_state=tensor_type(a)==TT_FULL ? compressed : reconstructed;
1621 try {
1622 ensure_tree_state_respecting_fence(a,operating_state,fence);
1623 ensure_tree_state_respecting_fence(b,operating_state,fence);
1624 } catch (...) {
1625 print("could not respect fence in gaxpy");
1626 change_tree_state(a,operating_state,true);
1627 change_tree_state(b,operating_state,true);
1628 }
1629
1630 if (operating_state==compressed) {
1631 for (unsigned int i=0; i<a.size(); ++i) result[i]=gaxpy_oop(alpha, a[i], beta, b[i], false);
1632 } else {
1633 for (unsigned int i=0; i<a.size(); ++i) result[i]=gaxpy_oop_reconstructed(alpha, a[i], beta, b[i], false);
1634 }
1635
1636 if (fence) world.gop.fence();
1637 return result;
1638 }
1639
1640
1641 /// out-of-place gaxpy for a vectors and a function: result[i] = alpha * a[i] + beta * b
1642 template <typename T, typename Q, typename R, std::size_t NDIM>
1643 std::vector<Function<TENSOR_RESULT_TYPE(Q,TENSOR_RESULT_TYPE(T,R)),NDIM> >
1645 const std::vector< Function<T,NDIM> >& a,
1646 Q beta,
1647 const Function<R,NDIM>& b,
1648 bool fence=true) {
1649
1650 typedef TENSOR_RESULT_TYPE(Q,TENSOR_RESULT_TYPE(T,R)) resultT;
1651 if (a.size()==0) return std::vector<Function<resultT,NDIM> >();
1652
1653 World& world=a[0].world();
1654 try {
1656 // ensure_tree_state_respecting_fence({b},compressed,fence);
1658 } catch (...) {
1659 print("could not respect fence in gaxpy_oop");
1660 compress(world,a);
1661 b.compress();
1662 }
1663 std::vector<Function<resultT,NDIM> > result(a.size());
1664 for (unsigned int i=0; i<a.size(); ++i) {
1665 result[i]=gaxpy_oop(alpha, a[i], beta, b, false);
1666 }
1667 if (fence) world.gop.fence();
1668 return result;
1669 }
1670
1671
1672 /// Generalized A*X+Y for vectors of functions ---- a[i] = alpha*a[i] + beta*b[i]
1673 template <typename T, typename Q, typename R, std::size_t NDIM>
1674 void gaxpy(Q alpha, std::vector<Function<T,NDIM>>& a, Q beta, const std::vector<Function<R,NDIM>>& b, const bool fence) {
1675 if (a.size() == 0) return;
1676 World& world=a.front().world();
1677 gaxpy(world,alpha,a,beta,b,fence);
1678 }
1679
1680 /// Generalized A*X+Y for vectors of functions ---- a[i] = alpha*a[i] + beta*b[i]
1681 template <typename T, typename Q, typename R, std::size_t NDIM>
1682 void gaxpy(World& world,
1683 Q alpha,
1684 std::vector< Function<T,NDIM> >& a,
1685 Q beta,
1686 const std::vector< Function<R,NDIM> >& b,
1687 bool fence=true) {
1688 PROFILE_BLOCK(Vgaxpy);
1689 MADNESS_ASSERT(a.size() == b.size());
1690 if (a.empty()) return;
1691
1692 auto tensor_type = [](const std::vector<Function<T,NDIM>>& v) {
1693 return v.front().get_impl()->get_tensor_type();
1694 };
1695
1696 // gaxpy can be done either in reconstructed or in compressed state
1697 bool do_in_reconstructed_state=tensor_type(a)!=TT_FULL;
1698 TreeState operating_state=do_in_reconstructed_state ? reconstructed : compressed;
1699
1700 if (operating_state==compressed) {
1701 // this is strict: both vectors have to be compressed
1702 try {
1703 ensure_tree_state_respecting_fence(a,operating_state,fence);
1704 ensure_tree_state_respecting_fence(b,operating_state,fence);
1705 } catch (...) {
1706 print("could not respect fence in gaxpy");
1707 change_tree_state(a,operating_state,true);
1708 change_tree_state(b,operating_state,true);
1709 }
1710 MADNESS_CHECK_THROW(get_tree_state(a)==get_tree_state(b),"gaxpy requires same tree state for all functions");
1711 MADNESS_CHECK_THROW(get_tree_state(a)==operating_state,"gaxpy requires reconstructed/compressed tree state for all functions");
1712 } else {
1713 // both vectors can be reconstructed or redundant_after_merge, and they don't have to be the same
1714 TreeState astate=get_tree_state(a);
1715 TreeState bstate=get_tree_state(b);
1716 if (not (astate==reconstructed or astate==redundant_after_merge)) {
1717 try {
1718 ensure_tree_state_respecting_fence(a,operating_state,fence);
1719 } catch (...) {
1720 print("could not respect fence in gaxpy for a");
1722 }
1723 }
1724 if (not (bstate==reconstructed or bstate==redundant_after_merge)) {
1725 try {
1726 ensure_tree_state_respecting_fence(b,operating_state,fence);
1727 } catch (...) {
1728 print("could not respect fence in gaxpy for b");
1730 }
1731 }
1732 }
1733
1734 // finally do the work
1735 for (unsigned int i=0; i<a.size(); ++i) {
1736 a[i].gaxpy(alpha, b[i], beta, false);
1737 }
1738 if (fence and (get_tree_state(a)==redundant_after_merge)) {
1739 for (unsigned int i=0; i<a.size(); ++i) a[i].get_impl()->finalize_sum();
1740 }
1741
1742 if (fence) world.gop.fence();
1743 }
1744
1745
1746 /// Applies a vector of operators to a vector of functions --- q[i] = apply(op[i],f[i])
1747 template <typename opT, typename R, std::size_t NDIM>
1748 std::vector< Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> >
1749 apply(World& world,
1750 const std::vector< std::shared_ptr<opT> >& op,
1751 const std::vector< Function<R,NDIM> > f) {
1752
1753 PROFILE_BLOCK(Vapplyv);
1754 MADNESS_ASSERT(f.size()==op.size());
1755
1756 std::vector< Function<R,NDIM> >& ncf = *const_cast< std::vector< Function<R,NDIM> >* >(&f);
1757
1758// reconstruct(world, f);
1759 make_nonstandard(world, ncf);
1760
1761 std::vector< Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> > result(f.size());
1762 for (unsigned int i=0; i<f.size(); ++i) {
1763 result[i] = apply_only(*op[i], f[i], false);
1764 result[i].get_impl()->set_tree_state(nonstandard_after_apply);
1765 }
1766
1767 world.gop.fence();
1768
1769 standard(world, ncf, false); // restores promise of logical constness
1770 reconstruct(result);
1771 world.gop.fence();
1772
1773 return result;
1774 }
1775
1776
1777 /// Applies an operator to a vector of functions --- q[i] = apply(op,f[i])
1778 template <typename T, typename R, std::size_t NDIM, std::size_t KDIM>
1779 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1781 const std::vector< Function<R,NDIM> > f) {
1782 return apply(op.get_world(),op,f);
1783 }
1784
1785
1786 /// Applies an operator to a vector of functions --- q[i] = apply(op,f[i])
1787 template <typename T, typename R, std::size_t NDIM, std::size_t KDIM>
1788 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1789 apply(World& world,
1791 const std::vector< Function<R,NDIM> > f) {
1792 PROFILE_BLOCK(Vapply);
1793
1794 std::vector< Function<R,NDIM> >& ncf = *const_cast< std::vector< Function<R,NDIM> >* >(&f);
1795 bool print_timings=(NDIM==6) and (world.rank()==0) and op.print_timings;
1796
1797 double wall0=wall_time();
1798// reconstruct(world, f);
1799 make_nonstandard(world, ncf);
1800 double wall1=wall_time();
1801 if (print_timings) printf("timer: %20.20s %8.2fs\n", "make_nonstandard", wall1-wall0);
1802
1803 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> > result(f.size());
1804 for (unsigned int i=0; i<f.size(); ++i) {
1805 result[i] = apply_only(op, f[i], false);
1806 }
1807
1808 world.gop.fence();
1809
1810 // restores promise of logical constness
1811 if (op.destructive()) {
1812 for (auto& ff : ncf) ff.clear(false);
1813 world.gop.fence();
1814 } else {
1815 reconstruct(world,f);
1816 }
1817
1818 // svd-tensor requires some cleanup after apply
1819 if (result[0].get_impl()->get_tensor_type()==TT_2D) {
1820 for (auto& r : result) r.get_impl()->finalize_apply();
1821 }
1822
1823 if (print_timings) {
1824 for (auto& r : result) r.get_impl()->print_timer();
1825 op.print_timer();
1826 }
1827 reconstruct(world, result);
1828
1829 return result;
1830 }
1831
1832 /// Normalizes a vector of functions --- v[i] = v[i].scale(1.0/v[i].norm2())
1833 template <typename T, std::size_t NDIM>
1834 void normalize(World& world, std::vector< Function<T,NDIM> >& v, bool fence=true) {
1835 PROFILE_BLOCK(Vnormalize);
1836 std::vector<double> nn = norm2s(world, v);
1837 for (unsigned int i=0; i<v.size(); ++i) v[i].scale(1.0/nn[i],false);
1838 if (fence) world.gop.fence();
1839 }
1840
1841 template <typename T, std::size_t NDIM>
1842 void print_size(World &world, const std::vector<Function<T,NDIM> > &v, const std::string &msg = "vectorfunction" ){
1843 if(v.empty()){
1844 if(world.rank()==0) std::cout << "print_size: " << msg << " is empty" << std::endl;
1845 }else if(v.size()==1){
1846 v.front().print_size(msg);
1847 }else{
1848 for(auto x:v){
1849 // print("impl",x.get_impl().get());
1850 x.print_size(msg);
1851 }
1852 }
1853 }
1854
1855 /// return the size of a vector of functions for each rank
1856 template <typename T, std::size_t NDIM>
1857 double get_size_local(World& world, const std::vector< Function<T,NDIM> >& v){
1858 double size=0.0;
1859 for(auto x:v){
1860 if (x.is_initialized()) size+=x.size_local();
1861 }
1862 const double d=sizeof(T);
1863 const double fac=1024*1024*1024;
1864 return size/fac*d;
1865 }
1866
1867 /// return the size of a function for each rank
1868 template <typename T, std::size_t NDIM>
1870 return get_size_local(f.world(),std::vector<Function<T,NDIM> >(1,f));
1871 }
1872
1873
1874 // gives back the size in GB
1875 template <typename T, std::size_t NDIM>
1876 double get_size(World& world, const std::vector< Function<T,NDIM> >& v){
1877
1878 if (v.empty()) return 0.0;
1879
1880 const double d=sizeof(T);
1881 const double fac=1024*1024*1024;
1882
1883 double size=0.0;
1884 for(unsigned int i=0;i<v.size();i++){
1885 if (v[i].is_initialized()) size+=v[i].size();
1886 }
1887
1888 return size/fac*d;
1889
1890 }
1891
1892 // gives back the size in GB
1893 template <typename T, std::size_t NDIM>
1894 double get_size(const Function<T,NDIM> & f){
1895 const double d=sizeof(T);
1896 const double fac=1024*1024*1024;
1897 double size=f.size();
1898 return size/fac*d;
1899 }
1900
1901 /// apply op on the input vector yielding an output vector of functions
1902
1903 /// @param[in] op the operator working on vin
1904 /// @param[in] vin vector of input Functions; needs to be refined to common level!
1905 /// @return vector of output Functions vout = op(vin)
1906 template <typename T, typename opT, std::size_t NDIM>
1907 std::vector<Function<T,NDIM> > multi_to_multi_op_values(const opT& op,
1908 const std::vector< Function<T,NDIM> >& vin,
1909 const bool fence=true) {
1910 MADNESS_ASSERT(vin.size()>0);
1911 MADNESS_ASSERT(vin[0].is_initialized()); // might be changed
1912 World& world=vin[0].world();
1913 Function<T,NDIM> dummy;
1914 dummy.set_impl(vin[0], false);
1915 std::vector<Function<T,NDIM> > vout=zero_functions<T,NDIM>(world, op.get_result_size());
1916 for (auto& out : vout) out.set_impl(vin[0],false);
1917 dummy.multi_to_multi_op_values(op, vin, vout, fence);
1918 return vout;
1919 }
1920
1921
1922
1923
1924 // convenience operators
1925
1926 /// result[i] = a[i] + b[i]
1927 template <typename T, std::size_t NDIM>
1928 std::vector<Function<T,NDIM> > operator+(const std::vector<Function<T,NDIM> >& lhs,
1929 const std::vector<Function<T,NDIM>>& rhs) {
1930 MADNESS_CHECK(lhs.size() == rhs.size());
1931 return gaxpy_oop(1.0,lhs,1.0,rhs);
1932 }
1933
1934 /// result[i] = a[i] - b[i]
1935 template <typename T, std::size_t NDIM>
1936 std::vector<Function<T,NDIM> > operator-(const std::vector<Function<T,NDIM> >& lhs,
1937 const std::vector<Function<T,NDIM> >& rhs) {
1938 MADNESS_CHECK(lhs.size() == rhs.size());
1939 return gaxpy_oop(1.0,lhs,-1.0,rhs);
1940 }
1941
1942 /// result[i] = a[i] + b
1943 template <typename T, std::size_t NDIM>
1944 std::vector<Function<T,NDIM> > operator+(const std::vector<Function<T,NDIM> >& lhs,
1945 const Function<T,NDIM>& rhs) {
1946 // MADNESS_CHECK(lhs.size() == rhs.size()); // no!!
1947 return gaxpy_oop(1.0,lhs,1.0,rhs);
1948 }
1949
1950 /// result[i] = a[i] - b
1951 template <typename T, std::size_t NDIM>
1952 std::vector<Function<T,NDIM> > operator-(const std::vector<Function<T,NDIM> >& lhs,
1953 const Function<T,NDIM>& rhs) {
1954 // MADNESS_CHECK(lhs.size() == rhs.size()); // no
1955 return gaxpy_oop(1.0,lhs,-1.0,rhs);
1956 }
1957
1958 /// result[i] = a + b[i]
1959 template <typename T, std::size_t NDIM>
1960 std::vector<Function<T,NDIM> > operator+(const Function<T,NDIM>& lhs,
1961 const std::vector<Function<T,NDIM> >& rhs) {
1962 // MADNESS_CHECK(lhs.size() == rhs.size()); // no
1963 return gaxpy_oop(1.0,rhs,1.0,lhs);
1964 }
1965
1966 /// result[i] = a - b[i]
1967 template <typename T, std::size_t NDIM>
1968 std::vector<Function<T,NDIM> > operator-(const Function<T,NDIM>& lhs,
1969 const std::vector<Function<T,NDIM> >& rhs) {
1970// MADNESS_CHECK(lhs.size() == rhs.size()); // no
1971 return gaxpy_oop(-1.0,rhs,1.0,lhs);
1972 }
1973
1974
1975 template <typename T, typename R, std::size_t NDIM>
1976 std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> > operator*(const R fac,
1977 const std::vector<Function<T,NDIM> >& rhs) {
1978 if (rhs.size()>0) {
1979 std::vector<Function<T,NDIM> > tmp=copy(rhs[0].world(),rhs);
1980 scale(tmp[0].world(),tmp,TENSOR_RESULT_TYPE(T,R)(fac));
1981 return tmp;
1982 }
1983 return std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> >();
1984 }
1985
1986 template <typename T, typename R, std::size_t NDIM>
1987 std::vector<Function<T,NDIM> > operator*(const std::vector<Function<T,NDIM> >& rhs,
1988 const R fac) {
1989 if (rhs.size()>0) {
1990 std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> > tmp=copy(rhs[0].world(),rhs);
1991 scale(tmp[0].world(),tmp,TENSOR_RESULT_TYPE(T,R)(fac));
1992 return tmp;
1993 }
1994 return std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> >();
1995 }
1996
1997 /// multiply a vector of functions with a function: r[i] = v[i] * a
1998 template <typename T, typename R, std::size_t NDIM>
2000 const std::vector<Function<R,NDIM> >& v) {
2001 if (v.size()>0) return mul(v[0].world(),a,v,true);
2002 return std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> >();
2003 }
2004
2005
2006 /// multiply a vector of functions with a function: r[i] = a * v[i]
2007 template <typename T, typename R, std::size_t NDIM>
2008 std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> > operator*(const std::vector<Function<T,NDIM> >& v,
2009 const Function<R,NDIM>& a) {
2010 if (v.size()>0) return mul(v[0].world(),a,v,true);
2011 return std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> >();
2012 }
2013
2014
2015 template <typename T, std::size_t NDIM>
2016 std::vector<Function<T,NDIM> > operator+=(std::vector<Function<T,NDIM> >& lhs, const std::vector<Function<T,NDIM> >& rhs) {
2017 MADNESS_CHECK(lhs.size() == rhs.size());
2018 if (lhs.size() > 0) gaxpy(lhs.front().world(), 1.0, lhs, 1.0, rhs);
2019 return lhs;
2020 }
2021
2022 template <typename T, std::size_t NDIM>
2023 std::vector<Function<T,NDIM> > operator-=(std::vector<Function<T,NDIM> >& lhs,
2024 const std::vector<Function<T,NDIM> >& rhs) {
2025 MADNESS_CHECK(lhs.size() == rhs.size());
2026 if (lhs.size() > 0) gaxpy(lhs.front().world(), 1.0, lhs, -1.0, rhs);
2027 return lhs;
2028 }
2029
2030 /// return the real parts of the vector's function (if complex)
2031 template <typename T, std::size_t NDIM>
2032 std::vector<Function<typename Tensor<T>::scalar_type,NDIM> >
2033 real(const std::vector<Function<T,NDIM> >& v, bool fence=true) {
2034 std::vector<Function<typename Tensor<T>::scalar_type,NDIM> > result(v.size());
2035 for (std::size_t i=0; i<v.size(); ++i) result[i]=real(v[i],false);
2036 if (fence and result.size()>0) result[0].world().gop.fence();
2037 return result;
2038 }
2039
2040 /// return the imaginary parts of the vector's function (if complex)
2041 template <typename T, std::size_t NDIM>
2042 std::vector<Function<typename Tensor<T>::scalar_type,NDIM> >
2043 imag(const std::vector<Function<T,NDIM> >& v, bool fence=true) {
2044 std::vector<Function<typename Tensor<T>::scalar_type,NDIM> > result(v.size());
2045 for (std::size_t i=0; i<v.size(); ++i) result[i]=imag(v[i],false);
2046 if (fence and result.size()>0) result[0].world().gop.fence();
2047 return result;
2048 }
2049
2050 /// shorthand gradient operator
2051
2052 /// returns the differentiated function f in all NDIM directions
2053 /// @param[in] f the function on which the grad operator works on
2054 /// @param[in] refine refinement before diff'ing makes the result more accurate
2055 /// @param[in] fence fence after completion; if reconstruction is needed always fence
2056 /// @return the vector \frac{\partial}{\partial x_i} f
2057 template <typename T, std::size_t NDIM>
2058 std::vector<Function<T,NDIM> > grad(const Function<T,NDIM>& f,
2059 bool refine=false, bool fence=true) {
2060
2061 World& world=f.world();
2062 f.reconstruct();
2063 if (refine) f.refine(); // refine to make result more precise
2064
2065 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2066 gradient_operator<T,NDIM>(world);
2067
2068 std::vector<Function<T,NDIM> > result(NDIM);
2069 for (size_t i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),f,false);
2070 if (fence) world.gop.fence();
2071 return result;
2072 }
2073
2074 // BLM first derivative
2075 template <typename T, std::size_t NDIM>
2076 std::vector<Function<T,NDIM> > grad_ble_one(const Function<T,NDIM>& f,
2077 bool refine=false, bool fence=true) {
2078
2079 World& world=f.world();
2080 f.reconstruct();
2081 if (refine) f.refine(); // refine to make result more precise
2082
2083 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2084 gradient_operator<T,NDIM>(world);
2085
2086 // Read in new coeff for each operator
2087 for (unsigned int i=0; i<NDIM; ++i) (*grad[i]).set_ble1();
2088
2089 std::vector<Function<T,NDIM> > result(NDIM);
2090 for (unsigned int i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),f,false);
2091 if (fence) world.gop.fence();
2092 return result;
2093 }
2094
2095 // BLM second derivative
2096 template <typename T, std::size_t NDIM>
2097 std::vector<Function<T,NDIM> > grad_ble_two(const Function<T,NDIM>& f,
2098 bool refine=false, bool fence=true) {
2099
2100 World& world=f.world();
2101 f.reconstruct();
2102 if (refine) f.refine(); // refine to make result more precise
2103
2104 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2105 gradient_operator<T,NDIM>(world);
2106
2107 // Read in new coeff for each operator
2108 for (unsigned int i=0; i<NDIM; ++i) (*grad[i]).set_ble2();
2109
2110 std::vector<Function<T,NDIM> > result(NDIM);
2111 for (unsigned int i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),f,false);
2112 if (fence) world.gop.fence();
2113 return result;
2114 }
2115
2116 // Bspline first derivative
2117 template <typename T, std::size_t NDIM>
2118 std::vector<Function<T,NDIM> > grad_bspline_one(const Function<T,NDIM>& f,
2119 bool refine=false, bool fence=true) {
2120
2121 World& world=f.world();
2122 f.reconstruct();
2123 if (refine) f.refine(); // refine to make result more precise
2124
2125 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2126 gradient_operator<T,NDIM>(world);
2127
2128 // Read in new coeff for each operator
2129 for (unsigned int i=0; i<NDIM; ++i) (*grad[i]).set_bspline1();
2130
2131 std::vector<Function<T,NDIM> > result(NDIM);
2132 for (unsigned int i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),f,false);
2133 if (fence) world.gop.fence();
2134 return result;
2135 }
2136
2137 // Bpsline second derivative
2138 template <typename T, std::size_t NDIM>
2139 std::vector<Function<T,NDIM> > grad_bpsline_two(const Function<T,NDIM>& f,
2140 bool refine=false, bool fence=true) {
2141
2142 World& world=f.world();
2143 f.reconstruct();
2144 if (refine) f.refine(); // refine to make result more precise
2145
2146 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2147 gradient_operator<T,NDIM>(world);
2148
2149 // Read in new coeff for each operator
2150 for (unsigned int i=0; i<NDIM; ++i) (*grad[i]).set_bspline2();
2151
2152 std::vector<Function<T,NDIM> > result(NDIM);
2153 for (unsigned int i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),f,false);
2154 if (fence) world.gop.fence();
2155 return result;
2156 }
2157
2158 // Bspline third derivative
2159 template <typename T, std::size_t NDIM>
2160 std::vector<Function<T,NDIM> > grad_bspline_three(const Function<T,NDIM>& f,
2161 bool refine=false, bool fence=true) {
2162
2163 World& world=f.world();
2164 f.reconstruct();
2165 if (refine) f.refine(); // refine to make result more precise
2166
2167 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2168 gradient_operator<T,NDIM>(world);
2169
2170 // Read in new coeff for each operator
2171 for (unsigned int i=0; i<NDIM; ++i) (*grad[i]).set_bspline3();
2172
2173 std::vector<Function<T,NDIM> > result(NDIM);
2174 for (unsigned int i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),f,false);
2175 if (fence) world.gop.fence();
2176 return result;
2177 }
2178
2179
2180
2181 /// shorthand div operator
2182
2183 /// returns the dot product of nabla with a vector f
2184 /// @param[in] f the vector of functions on which the div operator works on
2185 /// @param[in] refine refinement before diff'ing makes the result more accurate
2186 /// @param[in] fence fence after completion; currently always fences
2187 /// @return the vector \frac{\partial}{\partial x_i} f
2188 /// TODO: add this to operator fusion
2189 template <typename T, std::size_t NDIM>
2191 bool do_refine=false, bool fence=true) {
2192
2193 MADNESS_ASSERT(v.size()>0);
2194 World& world=v[0].world();
2195 reconstruct(world,v);
2196 if (do_refine) refine(world,v); // refine to make result more precise
2197
2198 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2199 gradient_operator<T,NDIM>(world);
2200
2201 std::vector<Function<T,NDIM> > result(NDIM);
2202 for (size_t i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),v[i],false);
2203 world.gop.fence();
2204 return sum(world,result,fence);
2205 }
2206
2207 /// shorthand rot operator
2208
2209 /// returns the cross product of nabla with a vector f
2210 /// @param[in] f the vector of functions on which the rot operator works on
2211 /// @param[in] refine refinement before diff'ing makes the result more accurate
2212 /// @param[in] fence fence after completion; currently always fences
2213 /// @return the vector \frac{\partial}{\partial x_i} f
2214 /// TODO: add this to operator fusion
2215 template <typename T, std::size_t NDIM>
2216 std::vector<Function<T,NDIM> > rot(const std::vector<Function<T,NDIM> >& v,
2217 bool do_refine=false, bool fence=true) {
2218
2219 MADNESS_ASSERT(v.size()==3);
2220 World& world=v[0].world();
2221 reconstruct(world,v);
2222 if (do_refine) refine(world,v); // refine to make result more precise
2223
2224 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2225 gradient_operator<T,NDIM>(world);
2226
2227 std::vector<Function<T,NDIM> > d(NDIM),dd(NDIM);
2228 d[0]=apply(*(grad[1]),v[2],false); // Dy z
2229 d[1]=apply(*(grad[2]),v[0],false); // Dz x
2230 d[2]=apply(*(grad[0]),v[1],false); // Dx y
2231 dd[0]=apply(*(grad[2]),v[1],false); // Dz y
2232 dd[1]=apply(*(grad[0]),v[2],false); // Dx z
2233 dd[2]=apply(*(grad[1]),v[0],false); // Dy x
2234 world.gop.fence();
2235
2236 compress(world,d,false);
2237 compress(world,dd,false);
2238 world.gop.fence();
2239 d[0].gaxpy(1.0,dd[0],-1.0,false);
2240 d[1].gaxpy(1.0,dd[1],-1.0,false);
2241 d[2].gaxpy(1.0,dd[2],-1.0,false);
2242
2243 world.gop.fence();
2244 reconstruct(d);
2245 return d;
2246 }
2247
2248 /// shorthand cross operator
2249
2250 /// returns the cross product of vectors f and g
2251 /// @param[in] f the vector of functions on which the rot operator works on
2252 /// @param[in] g the vector of functions on which the rot operator works on
2253 /// @param[in] fence fence after completion; currently always fences
2254 /// @return the vector \frac{\partial}{\partial x_i} f
2255 /// TODO: add this to operator fusion
2256 template <typename T, typename R, std::size_t NDIM>
2257 std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> > cross(const std::vector<Function<T,NDIM> >& f,
2258 const std::vector<Function<R,NDIM> >& g,
2259 bool do_refine=false, bool fence=true) {
2260
2261 MADNESS_ASSERT(f.size()==3);
2262 MADNESS_ASSERT(g.size()==3);
2263 World& world=f[0].world();
2264 reconstruct(world,f,false);
2265 reconstruct(world,g);
2266
2267 std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> > d(f.size()),dd(f.size());
2268
2269 d[0]=mul(f[1],g[2],false);
2270 d[1]=mul(f[2],g[0],false);
2271 d[2]=mul(f[0],g[1],false);
2272
2273 dd[0]=mul(f[2],g[1],false);
2274 dd[1]=mul(f[0],g[2],false);
2275 dd[2]=mul(f[1],g[0],false);
2276 world.gop.fence();
2277
2278 compress(world,d,false);
2279 compress(world,dd);
2280
2281 d[0].gaxpy(1.0,dd[0],-1.0,false);
2282 d[1].gaxpy(1.0,dd[1],-1.0,false);
2283 d[2].gaxpy(1.0,dd[2],-1.0,false);
2284
2285
2286 world.gop.fence();
2287 return d;
2288 }
2289
2290 template<typename T, std::size_t NDIM>
2291 void load_balance(World& world, std::vector<Function<T,NDIM> >& vf) {
2292
2293 struct LBCost {
2294 LBCost() = default;
2295 double operator()(const Key<NDIM>& key, const FunctionNode<T,NDIM>& node) const {
2296 return node.coeff().size();
2297 }
2298 };
2299
2300 LoadBalanceDeux<6> lb(world);
2301 for (const auto& f : vf) lb.add_tree(f, LBCost());
2303
2304 }
2305
2306 /// load a vector of functions
2307 template<typename T, size_t NDIM>
2308 void load_function(World& world, std::vector<Function<T,NDIM> >& f,
2309 const std::string name) {
2310 if (world.rank()==0) print("loading vector of functions",name);
2312 std::size_t fsize=0;
2313 ar & fsize;
2314 f.resize(fsize);
2315 for (std::size_t i=0; i<fsize; ++i) ar & f[i];
2316 }
2317
2318 /// save a vector of functions
2319 template<typename T, size_t NDIM>
2320 void save_function(const std::vector<Function<T,NDIM> >& f, const std::string name) {
2321 if (f.size()>0) {
2322 World& world=f.front().world();
2323 if (world.rank()==0) print("saving vector of functions",name);
2325 std::size_t fsize=f.size();
2326 ar & fsize;
2327 for (std::size_t i=0; i<fsize; ++i) ar & f[i];
2328 }
2329 }
2330
2331
2332}
2333#endif // MADNESS_MRA_VMRA_H__INCLUDED
double q(double t)
Definition DKops.h:18
Definition test_ar.cc:118
long size() const
Returns the number of elements in the tensor.
Definition basetensor.h:138
Implements derivatives operators with variety of boundary conditions on simulation domain.
Definition derivative.h:266
Definition distributed_matrix.h:68
Manages data associated with a row/column/block distributed array.
Definition distributed_matrix.h:388
static void redistribute(World &world, const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > &newpmap)
Sets the default process map and redistributes all functions using the old map.
Definition funcdefaults.h:433
static TensorType get_tensor_type()
Returns the default tensor type.
Definition funcdefaults.h:325
FunctionFactory implements the named-parameter idiom for Function.
Definition function_factory.h:86
FunctionFactory & compressed(bool value=true)
Definition function_factory.h:168
FunctionImpl holds all Function state to facilitate shallow copy semantics.
Definition funcimpl.h:945
static Tensor< TENSOR_RESULT_TYPE(T, R) > inner_local(const std::vector< const FunctionImpl< T, NDIM > * > &left, const std::vector< const FunctionImpl< R, NDIM > * > &right, bool sym)
Definition funcimpl.h:6081
void undo_redundant(const bool fence)
convert this from redundant to standard reconstructed form
Definition mraimpl.h:1538
void multiply(const implT *f, const FunctionImpl< T, LDIM > *g, const int particle)
multiply f (a pair function of NDIM) with an orbital g (LDIM=NDIM/2)
Definition funcimpl.h:3659
void change_tree_state(const TreeState finalstate, bool fence=true)
change the tree state of this function, might or might not respect fence!
Definition mraimpl.h:1407
static Tensor< TENSOR_RESULT_TYPE(T, R)> dot_local(const std::vector< const FunctionImpl< T, NDIM > * > &left, const std::vector< const FunctionImpl< R, NDIM > * > &right, bool sym)
Definition funcimpl.h:6133
FunctionNode holds the coefficients, etc., at each node of the 2^NDIM-tree.
Definition funcimpl.h:127
coeffT & coeff()
Returns a non-const reference to the tensor containing the coeffs.
Definition funcimpl.h:227
A multiresolution adaptive numerical function.
Definition mra.h:139
World & world() const
Returns the world.
Definition mra.h:693
Function< T, NDIM > & gaxpy(const T &alpha, const Function< Q, NDIM > &other, const R &beta, bool fence=true)
Inplace, general bi-linear operation in wavelet basis. No communication except for optional fence.
Definition mra.h:1055
void set_impl(const std::shared_ptr< FunctionImpl< T, NDIM > > &impl)
Replace current FunctionImpl with provided new one.
Definition mra.h:666
void multi_to_multi_op_values(const opT &op, const std::vector< Function< T, NDIM > > &vin, std::vector< Function< T, NDIM > > &vout, const bool fence=true)
apply op on the input vector yielding an output vector of functions
Definition mra.h:1615
bool is_initialized() const
Returns true if the function is initialized.
Definition mra.h:167
long size() const
Definition lowranktensor.h:488
Key is the index for a node of the 2^NDIM-tree.
Definition key.h:69
Definition lbdeux.h:228
std::shared_ptr< WorldDCPmapInterface< keyT > > load_balance(double fac=1.0, bool printstuff=false)
Actually does the partitioning of the tree.
Definition lbdeux.h:320
void add_tree(const Function< T, NDIM > &f, const costT &costfn, bool fence=false)
Accumulates cost from a function.
Definition lbdeux.h:289
Convolutions in separated form (including Gaussian)
Definition operator.h:139
A slice defines a sub-range or patch of a dimension.
Definition slice.h:103
A tensor is a multidimensional array.
Definition tensor.h:317
TensorTypeData< T >::scalar_type scalar_type
C++ typename of the real type associated with a complex type.
Definition tensor.h:409
T * ptr()
Returns a pointer to the internal data.
Definition tensor.h:1840
A simple, fixed dimension vector.
Definition vector.h:64
Interface to be provided by any process map.
Definition worlddc.h:122
void fence(bool debug=false)
Synchronizes all processes in communicator AND globally ensures no pending AM or tasks.
Definition worldgop.cc:161
void sum(T *buf, size_t nelem)
Inplace global sum while still processing AM & tasks.
Definition worldgop.h:872
void fence()
Returns after all local tasks have completed.
Definition world_task_queue.h:1384
A parallel world class.
Definition world.h:132
WorldTaskQueue & taskq
Task queue.
Definition world.h:206
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
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
int integer
Definition crayio.c:25
static const double R
Definition csqrt.cc:46
Declaration and initialization of tree traversal functions and generic derivative.
auto T(World &world, response_space &f) -> response_space
Definition global_functions.cc:28
Tensor< T > conj_transpose(const Tensor< T > &t)
Returns a new deep copy of the complex conjugate transpose of the input tensor.
Definition tensor.h:2043
const double beta
Definition gygi_soltion.cc:62
static const double v
Definition hatom_sf_dirac.cc:20
Tensor< double > op(const Tensor< double > &x)
Definition kain.cc:508
#define rot(x, k)
Definition lookup3.c:72
#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_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
Main include file for MADNESS and defines Function interface.
static const bool VERIFY_TREE
Definition mra.h:57
Namespace for all elements and tools of MADNESS.
Definition DFParameters.h:10
void save_function(const std::vector< Function< T, NDIM > > &f, const std::string name)
save a vector of functions
Definition vmra.h:2320
bool ensure_tree_state_respecting_fence(const std::vector< Function< T, NDIM > > &v, const TreeState state, bool fence)
ensure v has the requested tree state, change the tree state of v if necessary and no fence is given
Definition vmra.h:316
void rr_cholesky(Tensor< T > &A, typename Tensor< T >::scalar_type tol, Tensor< integer > &piv, int &rank)
Compute the rank-revealing Cholesky factorization.
Definition lapack.cc:1203
void make_redundant(World &world, const std::vector< Function< T, NDIM > > &v, bool fence=true)
change tree_state of a vector of functions to redundant
Definition vmra.h:186
std::vector< Function< T, NDIM > > orthonormalize_rrcd(const std::vector< Function< T, NDIM > > &v, Tensor< T > &ovlp, const double tol, Tensor< integer > &piv, int &rank)
Definition vmra.h:636
Function< double, NDIM > abssq(const Function< double_complex, NDIM > &z, bool fence=true)
Returns a new function that is the square of the absolute value of the input.
Definition mra.h:2825
Function< TENSOR_RESULT_TYPE(L, R), NDIM > gaxpy_oop(TENSOR_RESULT_TYPE(L, R) alpha, const Function< L, NDIM > &left, TENSOR_RESULT_TYPE(L, R) beta, const Function< R, NDIM > &right, bool fence=true)
Returns new function alpha*left + beta*right optional fence and no automatic compression.
Definition mra.h:2011
Function< typename TensorTypeData< Q >::scalar_type, NDIM > abs_square(const Function< Q, NDIM > &func)
Definition complexfun.h:121
Function< T, NDIM > square(const Function< T, NDIM > &f, bool fence=true)
Create a new function that is the square of f - global comm only if not reconstructed.
Definition mra.h:2793
response_space scale(response_space a, double b)
Function< TENSOR_RESULT_TYPE(L, R), NDIM > sub(const Function< L, NDIM > &left, const Function< R, NDIM > &right, bool fence=true)
Same as operator- but with optional fence and no automatic compression.
Definition mra.h:2065
std::vector< Function< T, NDIM > > reduce_rank(std::vector< Function< T, NDIM > > v, double thresh=0.0, bool fence=true)
reduces the tensor rank of the coefficient tensor (if applicable)
Definition vmra.h:366
Tensor< double > norm2s_T(World &world, const std::vector< Function< T, NDIM > > &v)
Computes the 2-norms of a vector of functions.
Definition vmra.h:882
std::vector< double > norm2s(World &world, const std::vector< Function< T, NDIM > > &v)
Computes the 2-norms of a vector of functions.
Definition vmra.h:869
std::vector< Function< T, NDIM > > grad_bspline_one(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
Definition vmra.h:2118
void set_impl(std::vector< Function< T, NDIM > > &v, const std::vector< std::shared_ptr< FunctionImpl< T, NDIM > > > vimpl)
Definition vmra.h:712
Function< Q, NDIM > convert(const Function< T, NDIM > &f, bool fence=true)
Type conversion implies a deep copy. No communication except for optional fence.
Definition mra.h:2143
Function< TENSOR_RESULT_TYPE(Q, T), NDIM > mul(const Q alpha, const Function< T, NDIM > &f, bool fence=true)
Returns new function equal to alpha*f(x) with optional fence.
Definition mra.h:1795
void truncate(World &world, response_space &v, double tol, bool fence)
Definition basic_operators.cc:31
std::vector< Function< T, NDIM > > orthonormalize_symmetric(const std::vector< Function< T, NDIM > > &v, const Tensor< T > &ovlp, double lindep=1e-12)
symmetric orthonormalization (see e.g. Szabo/Ostlund)
Definition vmra.h:474
std::vector< std::shared_ptr< FunctionImpl< T, NDIM > > > get_impl(const std::vector< Function< T, NDIM > > &v)
Definition vmra.h:705
Function< T, NDIM > div(const std::vector< Function< T, NDIM > > &v, bool do_refine=false, bool fence=true)
shorthand div operator
Definition vmra.h:2190
std::vector< Function< T, NDIM > > orthonormalize_cd(const std::vector< Function< T, NDIM > > &v, Tensor< T > &ovlp)
Definition vmra.h:599
std::vector< Function< T, NDIM > > copy_n(World &world, const Function< T, NDIM > &v, const unsigned int n, bool fence=true)
Returns a vector of n deep copies of a function.
Definition vmra.h:1374
void norm_tree(World &world, const std::vector< Function< T, NDIM > > &v, bool fence=true)
Makes the norm tree for all functions in a vector.
Definition vmra.h:1205
tensorT Q2(const tensorT &s)
Given overlap matrix, return rotation with 2nd order error to orthonormalize the vectors.
Definition SCF.cc:137
std::vector< Function< TENSOR_RESULT_TYPE(T, R), NDIM > > transform(World &world, const std::vector< Function< T, NDIM > > &v, const Tensor< R > &c, bool fence=true)
Transforms a vector of functions according to new[i] = sum[j] old[j]*c[j,i].
Definition vmra.h:731
std::vector< Function< TENSOR_RESULT_TYPE(T, R), NDIM > > cross(const std::vector< Function< T, NDIM > > &f, const std::vector< Function< R, NDIM > > &g, bool do_refine=false, bool fence=true)
shorthand cross operator
Definition vmra.h:2257
TreeState
Definition funcdefaults.h:59
@ nonstandard_after_apply
s and d coeffs, state after operator application
Definition funcdefaults.h:64
@ redundant_after_merge
s coeffs everywhere, must be summed up to yield the result
Definition funcdefaults.h:66
@ reconstructed
s coeffs at the leaves only
Definition funcdefaults.h:60
@ nonstandard
s and d coeffs in internal nodes
Definition funcdefaults.h:62
@ unknown
Definition funcdefaults.h:68
@ compressed
d coeffs in internal nodes, s and d coeffs at the root
Definition funcdefaults.h:61
@ redundant
s coeffs everywhere
Definition funcdefaults.h:65
@ nonstandard_with_leaves
like nonstandard, with s coeffs at the leaves
Definition funcdefaults.h:63
Function< T, NDIM > conj(const Function< T, NDIM > &f, bool fence=true)
Return the complex conjugate of the input function with the same distribution and optional fence.
Definition mra.h:2157
void cholesky(Tensor< T > &A)
Compute the Cholesky factorization.
Definition lapack.cc:1174
std::vector< std::vector< Function< TENSOR_RESULT_TYPE(T, R), NDIM > > > matrix_mul_sparse(World &world, const std::vector< Function< R, NDIM > > &f, const std::vector< Function< R, NDIM > > &g, double tol, bool fence=true, bool symm=false)
Outer product of a vector of functions with a vector of functions using sparsity.
Definition vmra.h:1173
void standard(World &world, std::vector< Function< T, NDIM > > &v, bool fence=true)
Generates standard form of a vector of functions.
Definition vmra.h:243
void compress(World &world, const std::vector< Function< T, NDIM > > &v, bool fence=true)
Compress a vector of functions.
Definition vmra.h:149
const std::vector< Function< T, NDIM > > & reconstruct(const std::vector< Function< T, NDIM > > &v)
reconstruct a vector of functions
Definition vmra.h:162
response_space transpose(response_space &f)
Definition basic_operators.cc:10
std::vector< Function< T, NDIM > > impl2function(const std::vector< std::shared_ptr< FunctionImpl< T, NDIM > > > vimpl)
Definition vmra.h:718
std::vector< Function< T, NDIM > > grad_bpsline_two(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
Definition vmra.h:2139
void set_thresh(World &world, std::vector< Function< T, NDIM > > &v, double thresh, bool fence=true)
Sets the threshold in a vector of functions.
Definition vmra.h:1309
double norm2(World &world, const std::vector< Function< T, NDIM > > &v)
Computes the 2-norm of a vector of functions.
Definition vmra.h:895
std::vector< Function< T, NDIM > > flatten(const std::vector< std::vector< Function< T, NDIM > > > &vv)
Definition vmra.h:698
std::vector< CCPairFunction< T, NDIM > > operator*(const double fac, const std::vector< CCPairFunction< T, NDIM > > &arg)
Definition ccpairfunction.h:1089
static void verify_tree(World &world, const std::vector< Function< T, NDIM > > &v)
Definition SCF.cc:74
std::vector< Function< T, NDIM > > multi_to_multi_op_values(const opT &op, const std::vector< Function< T, NDIM > > &vin, const bool fence=true)
apply op on the input vector yielding an output vector of functions
Definition vmra.h:1907
static const Slice _(0,-1, 1)
Tensor< T > inverse(const Tensor< T > &a_in)
invert general square matrix A
Definition lapack.cc:832
void load_balance(const real_function_6d &f, const bool leaf)
do some load-balancing
Definition madness/chem/mp2.cc:70
TreeState get_tree_state(const Function< T, NDIM > &f)
get tree state of a function
Definition mra.h:2841
std::vector< CCPairFunction< T, NDIM > > operator-(const std::vector< CCPairFunction< T, NDIM > > c1, const std::vector< CCPairFunction< T, NDIM > > &c2)
Definition ccpairfunction.h:1060
Function< TENSOR_RESULT_TYPE(L, R), NDIM > mul_sparse(const Function< L, NDIM > &left, const Function< R, NDIM > &right, double tol, bool fence=true)
Sparse multiplication — left and right must be reconstructed and if tol!=0 have tree of norms already...
Definition mra.h:1836
std::vector< Function< T, NDIM > > partial_mul(const Function< T, NDIM > f, const std::vector< Function< T, LDIM > > g, const int particle)
multiply a high-dimensional function with a low-dimensional function
Definition vmra.h:1244
Function< T, NDIM > gaxpy_oop_reconstructed(const double alpha, const Function< T, NDIM > &left, const double beta, const Function< T, NDIM > &right, const bool fence=true)
Returns new function alpha*left + beta*right optional fence, having both addends reconstructed.
Definition mra.h:2029
std::vector< Function< TENSOR_RESULT_TYPE(T, R), NDIM > > transform_reconstructed(World &world, const std::vector< Function< T, NDIM > > &v, const Tensor< R > &c, bool fence=true)
Transforms a vector of functions according to new[i] = sum[j] old[j]*c[j,i].
Definition vmra.h:760
void print(const T &t, const Ts &... ts)
Print items to std::cout (items separated by spaces) and terminate with a new line.
Definition print.h:226
response_space apply(World &world, std::vector< std::vector< std::shared_ptr< real_convolution_3d > > > &op, response_space &f)
Definition basic_operators.cc:43
std::vector< Function< T, NDIM > > append(const std::vector< Function< T, NDIM > > &lhs, const std::vector< Function< T, NDIM > > &rhs)
combine two vectors
Definition vmra.h:691
@ TT_2D
Definition gentensor.h:120
@ TT_FULL
Definition gentensor.h:120
void refine(World &world, const std::vector< Function< T, NDIM > > &vf, bool fence=true)
refine the functions according to the autorefine criteria
Definition vmra.h:196
void print_size(World &world, const std::vector< Function< T, NDIM > > &v, const std::string &msg="vectorfunction")
Definition vmra.h:1842
NDIM & f
Definition mra.h:2528
Function< TENSOR_RESULT_TYPE(L, R), NDIM > add(const Function< L, NDIM > &left, const Function< R, NDIM > &right, bool fence=true)
Same as operator+ but with optional fence and no automatic compression.
Definition mra.h:2021
const Function< T, NDIM > & change_tree_state(const Function< T, NDIM > &f, const TreeState finalstate, bool fence=true)
change tree state of a function
Definition mra.h:2854
Function< TENSOR_RESULT_TYPE(T, R), NDIM > dot(World &world, const std::vector< Function< T, NDIM > > &a, const std::vector< Function< R, NDIM > > &b, bool fence=true)
Multiplies and sums two vectors of functions r = \sum_i a[i] * b[i].
Definition vmra.h:1589
std::vector< CCPairFunction< T, NDIM > > & operator-=(std::vector< CCPairFunction< T, NDIM > > &rhs, const std::vector< CCPairFunction< T, NDIM > > &lhs)
Definition ccpairfunction.h:1082
std::vector< Function< T, NDIM > > orthonormalize(const std::vector< Function< T, NDIM > > &vf_in)
orthonormalize the vectors
Definition vmra.h:439
NDIM const Function< R, NDIM > & g
Definition mra.h:2528
double wall_time()
Returns the wall time in seconds relative to an arbitrary origin.
Definition timers.cc:48
std::vector< Function< T, NDIM > > grad_ble_one(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
Definition vmra.h:2076
Function< TENSOR_RESULT_TYPE(typename opT::opT, R), NDIM > apply_only(const opT &op, const Function< R, NDIM > &f, bool fence=true)
Apply operator ONLY in non-standard form - required other steps missing !!
Definition mra.h:2231
std::vector< Function< T, NDIM > > grad_bspline_three(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
Definition vmra.h:2160
std::vector< Function< T, NDIM > > zero_functions_compressed(World &world, int n, bool fence=true)
Generates a vector of zero functions (compressed)
Definition vmra.h:423
double inner(response_space &a, response_space &b)
Definition response_functions.h:639
double imag(double x)
Definition complexfun.h:56
void load_function(World &world, std::vector< Function< T, NDIM > > &f, const std::string name)
load a vector of functions
Definition vmra.h:2308
std::vector< Function< TENSOR_RESULT_TYPE(L, R), D > > vmulXX(const Function< L, D > &left, const std::vector< Function< R, D > > &vright, double tol, bool fence=true)
Use the vmra/mul(...) interface instead.
Definition mra.h:1901
void refine_to_common_level(World &world, std::vector< Function< T, NDIM > > &vf, bool fence=true)
refine all functions to a common (finest) level
Definition vmra.h:206
std::vector< Function< T, NDIM > > grad_ble_two(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
Definition vmra.h:2097
static bool print_timings
Definition SCF.cc:106
void normalize(World &world, std::vector< Function< T, NDIM > > &v, bool fence=true)
Normalizes a vector of functions — v[i] = v[i].scale(1.0/v[i].norm2())
Definition vmra.h:1834
std::vector< Function< T, NDIM > > zero_functions(World &world, int n, bool fence=true)
Generates a vector of zero functions (reconstructed)
Definition vmra.h:416
std::vector< CCPairFunction< T, NDIM > > operator+(const std::vector< CCPairFunction< T, NDIM > > c1, const std::vector< CCPairFunction< T, NDIM > > &c2)
Definition ccpairfunction.h:1052
std::vector< Function< T, NDIM > > grad(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
shorthand gradient operator
Definition vmra.h:2058
Function< T, NDIM > multiply(const Function< T, NDIM > f, const Function< T, LDIM > g, const int particle, const bool fence=true)
multiply a high-dimensional function with a low-dimensional function
Definition mra.h:2484
std::vector< Function< T, NDIM > > zero_functions_auto_tree_state(World &world, int n, bool fence=true)
Generates a vector of zero functions, either compressed or reconstructed, depending on tensor type.
Definition vmra.h:430
DistributedMatrix< T > matrix_dot(const DistributedMatrixDistribution &d, const std::vector< Function< T, NDIM > > &f, const std::vector< Function< T, NDIM > > &g, bool sym=false)
Definition vmra.h:1491
std::vector< Function< T, NDIM > > zero_functions_tree_state(World &world, int n, const TreeState state, bool fence=true)
Generates a vector of zero functions with a given tree state.
Definition vmra.h:395
double real(double x)
Definition complexfun.h:52
std::vector< CCPairFunction< T, NDIM > > & operator+=(std::vector< CCPairFunction< T, NDIM > > &lhs, const CCPairFunction< T, NDIM > &rhs)
Definition ccpairfunction.h:1068
std::vector< Function< T, NDIM > > orthonormalize_canonical(const std::vector< Function< T, NDIM > > &v, const Tensor< T > &ovlp, double lindep=1e-12)
Definition vmra.h:536
Tensor< TENSOR_RESULT_TYPE(T, R) > matrix_dot_old(World &world, const std::vector< Function< T, NDIM > > &f, const std::vector< Function< R, NDIM > > &g, bool sym=false)
Computes the matrix dot product of two function vectors - q(i,j) = dot(f[i],g[j])
Definition vmra.h:1554
std::string name(const FuncType &type, const int ex=-1)
Definition ccpairfunction.h:28
void matrix_inner(DistributedMatrix< T > &A, const std::vector< Function< T, NDIM > > &f, const std::vector< Function< T, NDIM > > &g, bool sym=false)
Definition distpm.cc:46
double get_size(World &world, const std::vector< Function< T, NDIM > > &v)
Definition vmra.h:1876
double get_size_local(World &world, const std::vector< Function< T, NDIM > > &v)
return the size of a vector of functions for each rank
Definition vmra.h:1857
Function< T, NDIM > copy(const Function< T, NDIM > &f, const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > &pmap, bool fence=true)
Create a new copy of the function with different distribution and optional fence.
Definition mra.h:2096
void syev(const Tensor< T > &A, Tensor< T > &V, Tensor< typename Tensor< T >::scalar_type > &e)
Real-symmetric or complex-Hermitian eigenproblem.
Definition lapack.cc:969
Tensor< TENSOR_RESULT_TYPE(T, R) > matrix_inner_old(World &world, const std::vector< Function< T, NDIM > > &f, const std::vector< Function< R, NDIM > > &g, bool sym=false)
Computes the matrix inner product of two function vectors - q(i,j) = inner(f[i],g[j])
Definition vmra.h:1017
void make_nonstandard(World &world, std::vector< Function< T, NDIM > > &v, bool fence=true)
Generates non-standard form of a vector of functions.
Definition vmra.h:233
void gaxpy(const double a, ScalarResult< T > &left, const double b, const T &right, const bool fence=true)
the result type of a macrotask must implement gaxpy
Definition macrotaskq.h:244
int distance(const madness::Hash_private::HashIterator< hashT > &it, const madness::Hash_private::HashIterator< hashT > &jt)
Definition worldhashmap.h:616
static long abs(long a)
Definition tensor.h:218
static const double b
Definition nonlinschro.cc:119
static const double d
Definition nonlinschro.cc:121
static const double a
Definition nonlinschro.cc:118
double Q(double a)
Definition relops.cc:20
static const double c
Definition relops.cc:10
static const double m
Definition relops.cc:9
static const double L
Definition rk.cc:46
static const double thresh
Definition rk.cc:45
Definition test_ar.cc:204
Definition mp2.h:63
double operator()(const Key< 6 > &key, const FunctionNode< double, 6 > &node) const
Definition mp2.h:70
Definition lowrankfunction.h:336
Definition dirac-hatom.cc:112
std::string ok(const bool b)
Definition test6.cc:43
AtomicInt sum
Definition test_atomicint.cc:46
int P
Definition test_binsorter.cc:9
void e()
Definition test_sig.cc:75
static const double alpha
Definition testcosine.cc:10
constexpr std::size_t NDIM
Definition testgconv.cc:54
#define TENSOR_RESULT_TYPE(L, R)
This macro simplifies access to TensorResultType.
Definition type_data.h:205
#define PROFILE_FUNC
Definition worldprofile.h:209
#define PROFILE_BLOCK(name)
Definition worldprofile.h:208