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 if(v.empty()) return v;
478
479 Tensor<T> U;
481 syev(ovlp,U,s);
482
483 // transform s to s^{-1}
484 for(size_t i=0;i<v.size();++i) s(i)=1.0/(sqrt(s(i)));
485
486 // save Ut before U gets modified with s^{-1}
487 const Tensor<T> Ut=transpose(U);
488 for(size_t i=0;i<v.size();++i){
489 for(size_t j=0;j<v.size();++j){
490 U(i,j)=U(i,j)*s(j);
491 }
492 }
493
494 Tensor<T> X=inner(U,Ut,1,0);
495
496 World& world=v.front().world();
497 return transform(world,v,X);
498
499 }
500 /// convenience routine for symmetric orthonormalization (see e.g. Szabo/Ostlund)
501 /// overlap matrix is calculated
502 /// @param[in] the vector to orthonormalize
503 template <typename T, std::size_t NDIM>
504 std::vector<Function<T,NDIM> > orthonormalize_symmetric(const std::vector<Function<T,NDIM> >& v){
505 if(v.empty()) return v;
506
507
508 World& world=v.front().world();
509 Tensor<T> ovlp = matrix_inner(world, v, v);
510
511 return orthonormalize_symmetric(v,ovlp);
512 }
513
514 /// canonical orthonormalization (see e.g. Szabo/Ostlund)
515 /// @param[in] the vector to orthonormalize
516 /// @param[in] overlap matrix
517 /// @param[in] lindep linear dependency threshold relative to largest eigenvalue
518 template <typename T, std::size_t NDIM>
519 std::vector<Function<T,NDIM> > orthonormalize_canonical(
520 const std::vector<Function<T,NDIM> >& v,
521 const Tensor<T>& ovlp,
522 double lindep) {
523
524 if(v.empty()) return v;
525
526 Tensor<T> U;
528 syev(ovlp,U,s);
529 lindep*=s(s.size()-1); // eigenvalues are in ascending order
530
531 // transform s to s^{-1}
532 int rank=0,lo=0;
534 for(size_t i=0;i<v.size();++i) {
535 if (s(i)>lindep) {
536 sqrts(i)=1.0/(sqrt(s(i)));
537 rank++;
538 } else {
539 sqrts(i)=0.0;
540 lo++;
541 }
542 }
543 MADNESS_ASSERT(size_t(lo+rank)==v.size());
544
545 for(size_t i=0;i<v.size();++i){
546 for(size_t j=0;j<v.size();++j){
547 U(i,j)=U(i,j)*(sqrts(j));
548 }
549 }
550 Tensor<T> X=U(_,Slice(lo,-1));
551
552 World& world=v.front().world();
553 return transform(world,v,X);
554
555 }
556
557 /// convenience routine for canonical routine for symmetric orthonormalization (see e.g. Szabo/Ostlund)
558 /// overlap matrix is calculated
559 /// @param[in] the vector to orthonormalize
560 template <typename T, std::size_t NDIM>
561 std::vector<Function<T,NDIM> > orthonormalize_canonical(const std::vector<Function<T,NDIM> >& v,
562 const double lindep){
563 if(v.empty()) return v;
564
565 World& world=v.front().world();
566 Tensor<T> ovlp = matrix_inner(world, v, v);
567
568 return orthonormalize_canonical(v,ovlp,lindep);
569 }
570
571 /// cholesky orthonormalization without pivoting
572 /// @param[in] the vector to orthonormalize
573 /// @param[in] overlap matrix, destroyed on return!
574 template <typename T, std::size_t NDIM>
575 std::vector<Function<T,NDIM> > orthonormalize_cd(
576 const std::vector<Function<T,NDIM> >& v,
577 Tensor<T>& ovlp) {
578
579 if (v.empty()) return v;
580
581 cholesky(ovlp); // destroys ovlp and gives back Upper ∆ Matrix from CD
582
583 Tensor<T> L = transpose(ovlp);
584 Tensor<T> Linv = inverse(L);
585 Tensor<T> U = transpose(Linv);
586
587 World& world=v.front().world();
588 return transform(world, v, U);
589
590 }
591
592 /// convenience routine for cholesky orthonormalization without pivoting
593 /// @param[in] the vector to orthonormalize
594 /// @param[in] overlap matrix
595 template <typename T, std::size_t NDIM>
596 std::vector<Function<T,NDIM> > orthonormalize_cd(const std::vector<Function<T,NDIM> >& v){
597 if(v.empty()) return v;
598
599 World& world=v.front().world();
600 Tensor<T> ovlp = matrix_inner(world, v, v);
601
602 return orthonormalize_cd(v,ovlp);
603 }
604
605 /// @param[in] the vector to orthonormalize
606 /// @param[in] overlap matrix, will be destroyed on return!
607 /// @param[in] tolerance for numerical rank reduction
608 /// @param[out] pivoting vector, no allocation on input needed
609 /// @param[out] rank
610 /// @return orthonrormalized vector (may or may not be truncated)
611 template <typename T, std::size_t NDIM>
612 std::vector<Function<T,NDIM> > orthonormalize_rrcd(
613 const std::vector<Function<T,NDIM> >& v,
614 Tensor<T>& ovlp,
615 const double tol,
616 Tensor<integer>& piv,
617 int& rank) {
618
619 if (v.empty()) {
620 return v;
621 }
622
623 rr_cholesky(ovlp,tol,piv,rank); // destroys ovlp and gives back Upper ∆ Matrix from CCD
624
625 // rearrange and truncate the functions according to the pivoting of the rr_cholesky
626 std::vector<Function<T,NDIM> > pv(rank);
627 for(integer i=0;i<rank;++i){
628 pv[i]=v[piv[i]];
629 }
630 ovlp=ovlp(Slice(0,rank-1),Slice(0,rank-1));
631
632 Tensor<T> L = transpose(ovlp);
633 Tensor<T> Linv = inverse(L);
634 Tensor<T> U = transpose(Linv);
635
636 World& world=v.front().world();
637 return transform(world, pv, U);
638 }
639
640 /// convenience routine for orthonromalize_cholesky: orthonromalize_cholesky without information on pivoting and rank
641 /// @param[in] the vector to orthonormalize
642 /// @param[in] overlap matrix
643 /// @param[in] tolerance for numerical rank reduction
644 template <typename T, std::size_t NDIM>
645 std::vector<Function<T,NDIM> > orthonormalize_rrcd(const std::vector<Function<T,NDIM> >& v, Tensor<T> ovlp , const double tol) {
646 Tensor<integer> piv;
647 int rank;
648 return orthonormalize_rrcd(v,ovlp,tol,piv,rank);
649 }
650
651 /// convenience routine for orthonromalize_cholesky: computes the overlap matrix and then calls orthonromalize_cholesky
652 /// @param[in] the vector to orthonormalize
653 /// @param[in] tolerance for numerical rank reduction
654 template <typename T, std::size_t NDIM>
655 std::vector<Function<T,NDIM> > orthonormalize_rrcd(const std::vector<Function<T,NDIM> >& v, const double tol) {
656 if (v.empty()) {
657 return v;
658 }
659 // compute overlap
660 World& world=v.front().world();
661 Tensor<T> ovlp = matrix_inner(world, v, v);
662 return orthonormalize_rrcd(v,ovlp,tol);
663 }
664
665 /// combine two vectors
666 template <typename T, std::size_t NDIM>
667 std::vector<Function<T,NDIM> > append(const std::vector<Function<T,NDIM> > & lhs, const std::vector<Function<T,NDIM> > & rhs){
668 std::vector<Function<T,NDIM> > v=lhs;
669 for (std::size_t i = 0; i < rhs.size(); ++i) v.push_back(rhs[i]);
670 return v;
671 }
672
673 template <typename T, std::size_t NDIM>
674 std::vector<Function<T,NDIM> > flatten(const std::vector< std::vector<Function<T,NDIM> > >& vv){
675 std::vector<Function<T,NDIM> >result;
676 for(const auto& x:vv) result=append(result,x);
677 return result;
678 }
679
680 template<typename T, std::size_t NDIM>
681 std::vector<std::shared_ptr<FunctionImpl<T,NDIM>>> get_impl(const std::vector<Function<T,NDIM>>& v) {
682 std::vector<std::shared_ptr<FunctionImpl<T,NDIM>>> result;
683 for (auto& f : v) result.push_back(f.get_impl());
684 return result;
685 }
686
687 template<typename T, std::size_t NDIM>
688 void set_impl(std::vector<Function<T,NDIM>>& v, const std::vector<std::shared_ptr<FunctionImpl<T,NDIM>>> vimpl) {
689 MADNESS_CHECK(vimpl.size()==v.size());
690 for (std::size_t i=0; i<vimpl.size(); ++i) v[i].set_impl(vimpl[i]);
691 }
692
693 template<typename T, std::size_t NDIM>
694 std::vector<Function<T,NDIM>> impl2function(const std::vector<std::shared_ptr<FunctionImpl<T,NDIM>>> vimpl) {
695 std::vector<Function<T,NDIM>> v(vimpl.size());
696 for (std::size_t i=0; i<vimpl.size(); ++i) v[i].set_impl(vimpl[i]);
697 return v;
698 }
699
700
701 /// Transforms a vector of functions according to new[i] = sum[j] old[j]*c[j,i]
702
703 /// Uses sparsity in the transformation matrix --- set small elements to
704 /// zero to take advantage of this.
705 template <typename T, typename R, std::size_t NDIM>
706 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> >
708 const std::vector< Function<T,NDIM> >& v,
709 const Tensor<R>& c,
710 bool fence=true) {
711
712 PROFILE_BLOCK(Vtransformsp);
713 typedef TENSOR_RESULT_TYPE(T,R) resultT;
714 int n = v.size(); // n is the old dimension
715 int m = c.dim(1); // m is the new dimension
716 MADNESS_CHECK(n==c.dim(0));
717
718 std::vector< Function<resultT,NDIM> > vc = zero_functions_compressed<resultT,NDIM>(world, m);
719 compress(world, v);
720
721 for (int i=0; i<m; ++i) {
722 for (int j=0; j<n; ++j) {
723 if (c(j,i) != R(0.0)) vc[i].gaxpy(resultT(1.0),v[j],resultT(c(j,i)),false);
724 }
725 }
726
727 if (fence) world.gop.fence();
728 return vc;
729 }
730
731 /// Transforms a vector of functions according to new[i] = sum[j] old[j]*c[j,i]
732
733 /// all trees are in reconstructed state, final trees have to be summed down if no fence is present
734 template <typename T, typename R, std::size_t NDIM>
735 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> >
737 const std::vector< Function<T,NDIM> >& v,
738 const Tensor<R>& c,
739 bool fence=true) {
740
741 PROFILE_BLOCK(Vtransformsp);
742 typedef TENSOR_RESULT_TYPE(T,R) resultT;
743 int n = v.size(); // n is the old dimension
744 int m = c.dim(1); // m is the new dimension
745 MADNESS_CHECK(n==c.dim(0));
746
747 // if we fence set the right tree state here, otherwise it has to be correct from the start.
749 for (const auto& vv : v) MADNESS_CHECK_THROW(
750 vv.get_impl()->get_tree_state()==reconstructed,"trees have to be reconstructed in transform_reconstructed");
751
752 std::vector< Function<resultT,NDIM> > result = zero_functions<resultT,NDIM>(world, m);
753
754 for (int i=0; i<m; ++i) {
755 result[i].get_impl()->set_tree_state(redundant_after_merge);
756 for (int j=0; j<n; ++j) {
757 if (c(j,i) != R(0.0)) v[j].get_impl()->accumulate_trees(*(result[i].get_impl()),resultT(c(j,i)),true);
758 }
759 }
760
761 // if we fence we can as well finish the job here. Otherwise no harm done, as the tree state is well-defined.
762 if (fence) {
763 world.gop.fence();
764 // for (auto& r : vc) r.sum_down(false);
765 for (auto& r : result) r.get_impl()->finalize_sum();
766 world.gop.fence();
767 }
768 return result;
769 }
770
771 /// this version of transform uses Function::vtransform and screens
772 /// using both elements of `c` and `v`
773 template <typename L, typename R, std::size_t NDIM>
774 std::vector< Function<TENSOR_RESULT_TYPE(L,R),NDIM> >
775 transform(World& world, const std::vector< Function<L,NDIM> >& v,
776 const Tensor<R>& c, double tol, bool fence) {
777 PROFILE_BLOCK(Vtransform);
778 MADNESS_ASSERT(v.size() == (unsigned int)(c.dim(0)));
779
780 std::vector< Function<TENSOR_RESULT_TYPE(L,R),NDIM> > vresult
781 = zero_functions_compressed<TENSOR_RESULT_TYPE(L,R),NDIM>(world, c.dim(1));
782
783 compress(world, v, true);
784 vresult[0].vtransform(v, c, vresult, tol, fence);
785 return vresult;
786 }
787
788 template <typename T, typename R, std::size_t NDIM>
789 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> >
791 const std::vector< Function<T,NDIM> >& v,
792 const DistributedMatrix<R>& c,
793 bool fence=true) {
795
796 typedef TENSOR_RESULT_TYPE(T,R) resultT;
797 long n = v.size(); // n is the old dimension
798 long m = c.rowdim(); // m is the new dimension
799 MADNESS_ASSERT(n==c.coldim());
800
801 // new(i) = sum(j) old(j) c(j,i)
802
803 Tensor<T> tmp(n,m);
804 c.copy_to_replicated(tmp); // for debugging
805 tmp = transpose(tmp);
806
807 std::vector< Function<resultT,NDIM> > vc = zero_functions_compressed<resultT,NDIM>(world, m);
808 compress(world, v);
809
810 for (int i=0; i<m; ++i) {
811 for (int j=0; j<n; ++j) {
812 if (tmp(j,i) != R(0.0)) vc[i].gaxpy(1.0,v[j],tmp(j,i),false);
813 }
814 }
815
816 if (fence) world.gop.fence();
817 return vc;
818 }
819
820
821 /// Scales inplace a vector of functions by distinct values
822 template <typename T, typename Q, std::size_t NDIM>
823 void scale(World& world,
824 std::vector< Function<T,NDIM> >& v,
825 const std::vector<Q>& factors,
826 bool fence=true) {
827 PROFILE_BLOCK(Vscale);
828 for (unsigned int i=0; i<v.size(); ++i) v[i].scale(factors[i],false);
829 if (fence) world.gop.fence();
830 }
831
832 /// Scales inplace a vector of functions by the same
833 template <typename T, typename Q, std::size_t NDIM>
834 void scale(World& world,
835 std::vector< Function<T,NDIM> >& v,
836 const Q factor,
837 bool fence=true) {
838 PROFILE_BLOCK(Vscale);
839 for (unsigned int i=0; i<v.size(); ++i) v[i].scale(factor,false);
840 if (fence) world.gop.fence();
841 }
842
843 /// Computes the 2-norms of a vector of functions
844 template <typename T, std::size_t NDIM>
845 std::vector<double> norm2s(World& world,
846 const std::vector< Function<T,NDIM> >& v) {
847 PROFILE_BLOCK(Vnorm2);
848 std::vector<double> norms(v.size());
850 for (unsigned int i=0; i<v.size(); ++i) norms[i] = v[i].norm2sq_local();
851 world.gop.sum(&norms[0], norms.size());
852 for (unsigned int i=0; i<v.size(); ++i) norms[i] = sqrt(norms[i]);
853 world.gop.fence();
854 return norms;
855 }
856 /// Computes the 2-norms of a vector of functions
857 template <typename T, std::size_t NDIM>
858 Tensor<double> norm2s_T(World& world, const std::vector<Function<T, NDIM>>& v) {
859 PROFILE_BLOCK(Vnorm2);
860 Tensor<double> norms(v.size());
862 for (unsigned int i = 0; i < v.size(); ++i) norms[i] = v[i].norm2sq_local();
863 world.gop.sum(&norms[0], norms.size());
864 for (unsigned int i = 0; i < v.size(); ++i) norms[i] = sqrt(norms[i]);
865 world.gop.fence();
866 return norms;
867 }
868
869 /// Computes the 2-norm of a vector of functions
870 template <typename T, std::size_t NDIM>
871 double norm2(World& world,const std::vector< Function<T,NDIM> >& v) {
872 PROFILE_BLOCK(Vnorm2);
873 if (v.size()==0) return 0.0;
875 std::vector<double> norms(v.size());
876 for (unsigned int i=0; i<v.size(); ++i) norms[i] = v[i].norm2sq_local();
877 world.gop.sum(&norms[0], norms.size());
878 for (unsigned int i=1; i<v.size(); ++i) norms[0] += norms[i];
879 world.gop.fence();
880 return sqrt(norms[0]);
881 }
882
883 inline double conj(double x) {
884 return x;
885 }
886
887 inline double conj(float x) {
888 return x;
889 }
890
891// !!! FIXME: this task is broken because FunctionImpl::inner_local forces a
892// future on return from WorldTaskQueue::reduce, which will causes a deadlock if
893// run inside a task. This behavior must be changed before this task can be used
894// again.
895//
896// template <typename T, typename R, std::size_t NDIM>
897// struct MatrixInnerTask : public TaskInterface {
898// Tensor<TENSOR_RESULT_TYPE(T,R)> result; // Must be a copy
899// const Function<T,NDIM>& f;
900// const std::vector< Function<R,NDIM> >& g;
901// long jtop;
902//
903// MatrixInnerTask(const Tensor<TENSOR_RESULT_TYPE(T,R)>& result,
904// const Function<T,NDIM>& f,
905// const std::vector< Function<R,NDIM> >& g,
906// long jtop)
907// : result(result), f(f), g(g), jtop(jtop) {}
908//
909// void run(World& world) {
910// for (long j=0; j<jtop; ++j) {
911// result(j) = f.inner_local(g[j]);
912// }
913// }
914//
915// private:
916// /// Get the task id
917//
918// /// \param id The id to set for this task
919// virtual void get_id(std::pair<void*,unsigned short>& id) const {
920// PoolTaskInterface::make_id(id, *this);
921// }
922// }; // struct MatrixInnerTask
923
924
925
926 template <typename T, std::size_t NDIM>
928 const std::vector< Function<T,NDIM> >& f,
929 const std::vector< Function<T,NDIM> >& g,
930 bool sym=false)
931 {
934 const int64_t n = A.coldim();
935 const int64_t m = A.rowdim();
936 MADNESS_ASSERT(int64_t(f.size()) == n && int64_t(g.size()) == m);
937
938 // Assume we can always create an ichunk*jchunk matrix locally
939 const int ichunk = 1000;
940 const int jchunk = 1000; // 1000*1000*8 = 8 MBytes
941 for (int64_t ilo=0; ilo<n; ilo+=ichunk) {
942 int64_t ihi = std::min(ilo + ichunk, n);
943 std::vector< Function<T,NDIM> > ivec(f.begin()+ilo, f.begin()+ihi);
944 for (int64_t jlo=0; jlo<m; jlo+=jchunk) {
945 int64_t jhi = std::min(jlo + jchunk, m);
946 std::vector< Function<T,NDIM> > jvec(g.begin()+jlo, g.begin()+jhi);
947
948 Tensor<T> P = matrix_inner(A.get_world(), ivec, jvec);
949 A.copy_from_replicated_patch(ilo, ihi - 1, jlo, jhi - 1, P);
950 }
951 }
952 return A;
953 }
954
955 /// Computes the matrix inner product of two function vectors - q(i,j) = inner(f[i],g[j])
956
957 /// For complex types symmetric is interpreted as Hermitian.
958
959 /// The current parallel loop is non-optimal but functional.
960 template <typename T, typename R, std::size_t NDIM>
962 const std::vector< Function<T,NDIM> >& f,
963 const std::vector< Function<R,NDIM> >& g,
964 bool sym=false)
965 {
966 world.gop.fence();
967 auto tensor_type = [](const std::vector<Function<T,NDIM>>& v) {
968 return v.front().get_impl()->get_tensor_type();
969 };
970 TreeState operating_state=tensor_type(f)==TT_FULL ? compressed : redundant;
971 ensure_tree_state_respecting_fence(f,operating_state,true);
972 ensure_tree_state_respecting_fence(g,operating_state,true);
973
974 std::vector<const FunctionImpl<T,NDIM>*> left(f.size());
975 std::vector<const FunctionImpl<R,NDIM>*> right(g.size());
976 for (unsigned int i=0; i<f.size(); i++) left[i] = f[i].get_impl().get();
977 for (unsigned int i=0; i<g.size(); i++) right[i]= g[i].get_impl().get();
978
980
981 world.gop.fence();
982 world.gop.sum(r.ptr(),f.size()*g.size());
983
984 return r;
985 }
986
987 /// Computes the matrix inner product of two function vectors - q(i,j) = inner(f[i],g[j])
988
989 /// For complex types symmetric is interpreted as Hermitian.
990 ///
991 /// The current parallel loop is non-optimal but functional.
992 template <typename T, typename R, std::size_t NDIM>
994 const std::vector< Function<T,NDIM> >& f,
995 const std::vector< Function<R,NDIM> >& g,
996 bool sym=false) {
997 PROFILE_BLOCK(Vmatrix_inner);
998 long n=f.size(), m=g.size();
999 Tensor< TENSOR_RESULT_TYPE(T,R) > r(n,m);
1000 if (sym) MADNESS_ASSERT(n==m);
1001
1002 world.gop.fence();
1003 compress(world, f);
1004 if ((void*)(&f) != (void*)(&g)) compress(world, g);
1005
1006 for (long i=0; i<n; ++i) {
1007 long jtop = m;
1008 if (sym) jtop = i+1;
1009 for (long j=0; j<jtop; ++j) {
1010 r(i,j) = f[i].inner_local(g[j]);
1011 if (sym) r(j,i) = conj(r(i,j));
1012 }
1013 }
1014
1015// for (long i=n-1; i>=0; --i) {
1016// long jtop = m;
1017// if (sym) jtop = i+1;
1018// world.taskq.add(new MatrixInnerTask<T,R,NDIM>(r(i,_), f[i], g, jtop));
1019// }
1020 world.gop.fence();
1021 world.gop.sum(r.ptr(),n*m);
1022
1023// if (sym) {
1024// for (int i=0; i<n; ++i) {
1025// for (int j=0; j<i; ++j) {
1026// r(j,i) = conj(r(i,j));
1027// }
1028// }
1029// }
1030 return r;
1031 }
1032
1033 /// Computes the element-wise inner product of two function vectors - q(i) = inner(f[i],g[i])
1034
1035 /// works in reconstructed or compressed state, state is chosen based on TensorType
1036 template <typename T, typename R, std::size_t NDIM>
1038 const std::vector< Function<T,NDIM> >& f,
1039 const std::vector< Function<R,NDIM> >& g) {
1040 PROFILE_BLOCK(Vinnervv);
1041 long n=f.size(), m=g.size();
1042 MADNESS_CHECK(n==m);
1043 Tensor< TENSOR_RESULT_TYPE(T,R) > r(n);
1044
1045 auto tensor_type = [](const std::vector<Function<T,NDIM>>& v) {
1046 return v.front().get_impl()->get_tensor_type();
1047 };
1048 TreeState operating_state=tensor_type(f)==TT_FULL ? compressed : redundant;
1049 ensure_tree_state_respecting_fence(f,operating_state,true);
1050 ensure_tree_state_respecting_fence(g,operating_state,true);
1051
1052 for (long i=0; i<n; ++i) r(i) = f[i].inner_local(g[i]);
1053
1054 world.taskq.fence();
1055 world.gop.sum(r.ptr(),n);
1056 world.gop.fence();
1057 return r;
1058 }
1059
1060
1061 /// Computes the inner product of a function with a function vector - q(i) = inner(f,g[i])
1062
1063 /// works in reconstructed or compressed state, state is chosen based on TensorType
1064 template <typename T, typename R, std::size_t NDIM>
1066 const Function<T,NDIM>& f,
1067 const std::vector< Function<R,NDIM> >& g) {
1068 PROFILE_BLOCK(Vinner);
1069 long n=g.size();
1070 Tensor< TENSOR_RESULT_TYPE(T,R) > r(n);
1071
1072 auto tensor_type = [](const std::vector<Function<T,NDIM>>& v) {
1073 return v.front().get_impl()->get_tensor_type();
1074 };
1075 TreeState operating_state=tensor_type(g)==TT_FULL ? compressed : redundant;
1076 f.change_tree_state(operating_state,false);
1077 ensure_tree_state_respecting_fence(g,operating_state,true);
1078 world.gop.fence();
1079
1080 for (long i=0; i<n; ++i) {
1081 r(i) = f.inner_local(g[i]);
1082 }
1083
1084 world.taskq.fence();
1085 world.gop.sum(r.ptr(),n);
1086 world.gop.fence();
1087 return r;
1088 }
1089
1090 /// inner function with right signature for the nonlinear solver
1091 /// this is needed for the KAIN solvers and other functions
1092 template <typename T, typename R, std::size_t NDIM>
1094 const std::vector< Function<R,NDIM> >& g){
1095 MADNESS_ASSERT(f.size()==g.size());
1096 if(f.empty()) return 0.0;
1097 else return inner(f[0].world(),f,g).sum();
1098 }
1099
1100
1101 /// Multiplies a function against a vector of functions --- q[i] = a * v[i]
1102 template <typename T, typename R, std::size_t NDIM>
1103 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1104 mul(World& world,
1105 const Function<T,NDIM>& a,
1106 const std::vector< Function<R,NDIM> >& v,
1107 bool fence=true) {
1108 PROFILE_BLOCK(Vmul);
1109 a.reconstruct(false);
1110 reconstruct(world, v, false);
1111 world.gop.fence();
1112 return vmulXX(a, v, 0.0, fence);
1113 }
1114
1115 /// Multiplies a function against a vector of functions using sparsity of a and v[i] --- q[i] = a * v[i]
1116 template <typename T, typename R, std::size_t NDIM>
1117 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1119 const Function<T,NDIM>& a,
1120 const std::vector< Function<R,NDIM> >& v,
1121 double tol,
1122 bool fence=true) {
1123 PROFILE_BLOCK(Vmulsp);
1124 a.reconstruct(false);
1125 reconstruct(world, v, false);
1126 world.gop.fence();
1127 for (unsigned int i=0; i<v.size(); ++i) {
1128 v[i].norm_tree(false);
1129 }
1130 a.norm_tree();
1131 return vmulXX(a, v, tol, fence);
1132 }
1133
1134
1135 /// Outer product of a vector of functions with a vector of functions using sparsity
1136
1137 /// \tparam T type parameter for first factor
1138 /// \tparam R type parameter for second factor
1139 /// \tparam NDIM dimension of first and second factors
1140 /// \param world the world
1141 /// \param f first vector of functions
1142 /// \param g second vector of functions
1143 /// \param tol threshold for multiplication
1144 /// \param fence force fence (will always fence if necessary)
1145 /// \param symm if true, only compute f(i) * g(j) for j<=i
1146 /// \return fg(i,j) = f(i) * g(j), as a vector of vectors
1147 template <typename T, typename R, std::size_t NDIM>
1148 std::vector<std::vector<Function<TENSOR_RESULT_TYPE(T, R), NDIM> > >
1150 const std::vector<Function<R, NDIM> > &f,
1151 const std::vector<Function<R, NDIM> > &g,
1152 double tol,
1153 bool fence = true,
1154 bool symm = false) {
1155 PROFILE_BLOCK(Vmulsp);
1156 bool same=(&f == &g);
1157 reconstruct(world, f, false);
1158 if (not same) reconstruct(world, g, false);
1159 world.gop.fence();
1160 for (auto& ff : f) ff.norm_tree(false);
1161 if (not same) for (auto& gg : g) gg.norm_tree(false);
1162 world.gop.fence();
1163
1164 std::vector<std::vector<Function<R,NDIM> > >result(f.size());
1165 std::vector<Function<R,NDIM>> g_i;
1166 for (int64_t i=f.size()-1; i>=0; --i) {
1167 if (!symm)
1168 result[i]= vmulXX(f[i], g, tol, false);
1169 else {
1170 if (g_i.empty()) g_i = g;
1171 g_i.resize(i+1); // this shrinks g_i down to single function for i=0
1172 result[i]= vmulXX(f[i], g_i, tol, false);
1173 }
1174 }
1175 if (fence) world.gop.fence();
1176 return result;
1177 }
1178
1179 /// Makes the norm tree for all functions in a vector
1180 template <typename T, std::size_t NDIM>
1181 void norm_tree(World& world,
1182 const std::vector< Function<T,NDIM> >& v,
1183 bool fence=true)
1184 {
1185 PROFILE_BLOCK(Vnorm_tree);
1186 for (unsigned int i=0; i<v.size(); ++i) {
1187 v[i].norm_tree(false);
1188 }
1189 if (fence) world.gop.fence();
1190 }
1191
1192 /// Multiplies two vectors of functions q[i] = a[i] * b[i]
1193 template <typename T, typename R, std::size_t NDIM>
1194 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1195 mul(World& world,
1196 const std::vector< Function<T,NDIM> >& a,
1197 const std::vector< Function<R,NDIM> >& b,
1198 bool fence=true) {
1199 PROFILE_BLOCK(Vmulvv);
1200 reconstruct(world, a, true);
1201 reconstruct(world, b, true);
1202// if (&a != &b) reconstruct(world, b, true); // fails if type(a) != type(b)
1203
1204 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> > q(a.size());
1205 for (unsigned int i=0; i<a.size(); ++i) {
1206 q[i] = mul(a[i], b[i], false);
1207 }
1208 if (fence) world.gop.fence();
1209 return q;
1210 }
1211
1212
1213 /// multiply a high-dimensional function with a low-dimensional function
1214
1215 /// @param[in] f NDIM function of NDIM dimensions
1216 /// @param[in] g LDIM function of LDIM
1217 /// @param[in] v dimension indices of f to multiply
1218 /// @return h[i](0,1,2,3) = f(0,1,2,3) * g[i](1,2,3) for v={1,2,3}
1219 template<typename T, std::size_t NDIM, std::size_t LDIM>
1220 std::vector<Function<T,NDIM> > partial_mul(const Function<T,NDIM> f, const std::vector<Function<T,LDIM> > g,
1221 const int particle) {
1222
1223 World& world=f.world();
1224 std::vector<Function<T,NDIM> > result(g.size());
1225 for (auto& r : result) r.set_impl(f, false);
1226
1227 FunctionImpl<T,NDIM>* fimpl=f.get_impl().get();
1228// fimpl->make_redundant(false);
1229 fimpl->change_tree_state(redundant,false);
1230 make_redundant(world,g,false);
1231 world.gop.fence();
1232
1233 for (std::size_t i=0; i<result.size(); ++i) {
1234 FunctionImpl<T,LDIM>* gimpl=g[i].get_impl().get();
1235 result[i].get_impl()->multiply(fimpl,gimpl,particle); // stupid naming inconsistency
1236 }
1237 world.gop.fence();
1238
1239 fimpl->undo_redundant(false);
1240 for (auto& ig : g) ig.get_impl()->undo_redundant(false);
1241 world.gop.fence();
1242 return result;
1243 }
1244
1245 template<typename T, std::size_t NDIM, std::size_t LDIM>
1246 std::vector<Function<T,NDIM> > multiply(const Function<T,NDIM> f, const std::vector<Function<T,LDIM> > g,
1247 const std::tuple<int,int,int> v) {
1248 return partial_mul<T,NDIM,LDIM>(f,g,std::array<int,3>({std::get<0>(v),std::get<1>(v),std::get<2>(v)}));
1249 }
1250
1251
1252/// Computes the square of a vector of functions --- q[i] = v[i]**2
1253 template <typename T, std::size_t NDIM>
1254 std::vector< Function<T,NDIM> >
1256 const std::vector< Function<T,NDIM> >& v,
1257 bool fence=true) {
1258 return mul<T,T,NDIM>(world, v, v, fence);
1259// std::vector< Function<T,NDIM> > vsq(v.size());
1260// for (unsigned int i=0; i<v.size(); ++i) {
1261// vsq[i] = square(v[i], false);
1262// }
1263// if (fence) world.gop.fence();
1264// return vsq;
1265 }
1266
1267
1268 /// Computes the square of a vector of functions --- q[i] = abs(v[i])**2
1269 template <typename T, std::size_t NDIM>
1270 std::vector< Function<typename Tensor<T>::scalar_type,NDIM> >
1271 abssq(World& world,
1272 const std::vector< Function<T,NDIM> >& v,
1273 bool fence=true) {
1274 typedef typename Tensor<T>::scalar_type scalartype;
1275 reconstruct(world,v);
1276 std::vector<Function<scalartype,NDIM> > result(v.size());
1277 for (size_t i=0; i<v.size(); ++i) result[i]=abs_square(v[i],false);
1278 if (fence) world.gop.fence();
1279 return result;
1280 }
1281
1282
1283 /// Sets the threshold in a vector of functions
1284 template <typename T, std::size_t NDIM>
1285 void set_thresh(World& world, std::vector< Function<T,NDIM> >& v, double thresh, bool fence=true) {
1286 for (unsigned int j=0; j<v.size(); ++j) {
1287 v[j].set_thresh(thresh,false);
1288 }
1289 if (fence) world.gop.fence();
1290 }
1291
1292 /// Returns the complex conjugate of the vector of functions
1293 template <typename T, std::size_t NDIM>
1294 std::vector< Function<T,NDIM> >
1295 conj(World& world,
1296 const std::vector< Function<T,NDIM> >& v,
1297 bool fence=true) {
1298 PROFILE_BLOCK(Vconj);
1299 std::vector< Function<T,NDIM> > r = copy(world, v); // Currently don't have oop conj
1300 for (unsigned int i=0; i<v.size(); ++i) {
1301 r[i].conj(false);
1302 }
1303 if (fence) world.gop.fence();
1304 return r;
1305 }
1306
1307 /// Returns a deep copy of a vector of functions
1308 template <typename T, typename R, std::size_t NDIM>
1309 std::vector< Function<R,NDIM> > convert(World& world,
1310 const std::vector< Function<T,NDIM> >& v, bool fence=true) {
1311 PROFILE_BLOCK(Vcopy);
1312 std::vector< Function<R,NDIM> > r(v.size());
1313 for (unsigned int i=0; i<v.size(); ++i) {
1314 r[i] = convert<T,R,NDIM>(v[i], false);
1315 }
1316 if (fence) world.gop.fence();
1317 return r;
1318 }
1319
1320
1321 /// Returns a deep copy of a vector of functions
1322 template <typename T, std::size_t NDIM>
1323 std::vector< Function<T,NDIM> >
1324 copy(World& world,
1325 const std::vector< Function<T,NDIM> >& v,
1326 bool fence=true) {
1327 PROFILE_BLOCK(Vcopy);
1328 std::vector< Function<T,NDIM> > r(v.size());
1329 for (unsigned int i=0; i<v.size(); ++i) {
1330 r[i] = copy(v[i], false);
1331 }
1332 if (fence) world.gop.fence();
1333 return r;
1334 }
1335
1336
1337 /// Returns a deep copy of a vector of functions
1338 template <typename T, std::size_t NDIM>
1339 std::vector< Function<T,NDIM> >
1340 copy(const std::vector< Function<T,NDIM> >& v, bool fence=true) {
1341 PROFILE_BLOCK(Vcopy);
1342 std::vector< Function<T,NDIM> > r(v.size());
1343 if (v.size()>0) r=copy(v.front().world(),v,fence);
1344 return r;
1345 }
1346
1347 /// Returns a vector of `n` deep copies of a function
1348 template <typename T, std::size_t NDIM>
1349 std::vector< Function<T,NDIM> >
1350 copy(World& world,
1351 const Function<T,NDIM>& v,
1352 const unsigned int n,
1353 bool fence=true) {
1354 PROFILE_BLOCK(Vcopy1);
1355 std::vector< Function<T,NDIM> > r(n);
1356 for (unsigned int i=0; i<n; ++i) {
1357 r[i] = copy(v, false);
1358 }
1359 if (fence) world.gop.fence();
1360 return r;
1361 }
1362
1363 /// Create a new copy of the function with different distribution and optional
1364 /// fence
1365
1366 /// Works in either basis. Different distributions imply
1367 /// asynchronous communication and the optional fence is
1368 /// collective.
1369 //
1370 /// Returns a deep copy of a vector of functions
1371
1372 template <typename T, std::size_t NDIM>
1373 std::vector<Function<T, NDIM>> copy(World& world,
1374 const std::vector<Function<T, NDIM>>& v,
1375 const std::shared_ptr<WorldDCPmapInterface<Key<NDIM>>>& pmap,
1376 bool fence = true) {
1377 PROFILE_BLOCK(Vcopy);
1378 std::vector<Function<T, NDIM>> r(v.size());
1379 for (unsigned int i = 0; i < v.size(); ++i) {
1380 r[i] = copy(v[i], pmap, false);
1381 }
1382 if (fence) world.gop.fence();
1383 return r;
1384 }
1385
1386 /// Returns new vector of functions --- q[i] = a[i] + b[i]
1387 template <typename T, typename R, std::size_t NDIM>
1388 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1389 add(World& world,
1390 const std::vector< Function<T,NDIM> >& a,
1391 const std::vector< Function<R,NDIM> >& b,
1392 bool fence=true) {
1393 PROFILE_BLOCK(Vadd);
1394 MADNESS_ASSERT(a.size() == b.size());
1395 compress(world, a);
1396 compress(world, b);
1397
1398 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> > r(a.size());
1399 for (unsigned int i=0; i<a.size(); ++i) {
1400 r[i] = add(a[i], b[i], false);
1401 }
1402 if (fence) world.gop.fence();
1403 return r;
1404 }
1405
1406 /// Returns new vector of functions --- q[i] = a + b[i]
1407 template <typename T, typename R, std::size_t NDIM>
1408 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1409 add(World& world,
1410 const Function<T,NDIM> & a,
1411 const std::vector< Function<R,NDIM> >& b,
1412 bool fence=true) {
1413 PROFILE_BLOCK(Vadd1);
1414 a.compress();
1415 compress(world, b);
1416
1417 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> > r(b.size());
1418 for (unsigned int i=0; i<b.size(); ++i) {
1419 r[i] = add(a, b[i], false);
1420 }
1421 if (fence) world.gop.fence();
1422 return r;
1423 }
1424 template <typename T, typename R, std::size_t NDIM>
1425 inline std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1426 add(World& world,
1427 const std::vector< Function<R,NDIM> >& b,
1428 const Function<T,NDIM> & a,
1429 bool fence=true) {
1430 return add(world, a, b, fence);
1431 }
1432
1433 /// Returns new vector of functions --- q[i] = a[i] - b[i]
1434 template <typename T, typename R, std::size_t NDIM>
1435 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1436 sub(World& world,
1437 const std::vector< Function<T,NDIM> >& a,
1438 const std::vector< Function<R,NDIM> >& b,
1439 bool fence=true) {
1440 PROFILE_BLOCK(Vsub);
1441 MADNESS_ASSERT(a.size() == b.size());
1442 compress(world, a);
1443 compress(world, b);
1444
1445 std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> > r(a.size());
1446 for (unsigned int i=0; i<a.size(); ++i) {
1447 r[i] = sub(a[i], b[i], false);
1448 }
1449 if (fence) world.gop.fence();
1450 return r;
1451 }
1452
1453 /// Returns new function --- q = sum_i f[i]
1454 template <typename T, std::size_t NDIM>
1455 Function<T, NDIM> sum(World& world, const std::vector<Function<T,NDIM> >& f,
1456 bool fence=true) {
1457
1458 compress(world, f);
1460
1461 for (unsigned int i=0; i<f.size(); ++i) r.gaxpy(1.0,f[i],1.0,false);
1462 if (fence) world.gop.fence();
1463 return r;
1464 }
1465
1466 template <typename T, std::size_t NDIM>
1468 const std::vector<Function<T, NDIM>>& f,
1469 const std::vector<Function<T, NDIM>>& g,
1470 bool sym=false)
1471 {
1474 const int64_t n = A.coldim();
1475 const int64_t m = A.rowdim();
1476 MADNESS_ASSERT(int64_t(f.size()) == n && int64_t(g.size()) == m);
1477
1478 // Assume we can always create an ichunk*jchunk matrix locally
1479 const int ichunk = 1000;
1480 const int jchunk = 1000; // 1000*1000*8 = 8 MBytes
1481 for (int64_t ilo = 0; ilo < n; ilo += ichunk) {
1482 int64_t ihi = std::min(ilo + ichunk, n);
1483 std::vector<Function<T, NDIM>> ivec(f.begin() + ilo, f.begin() + ihi);
1484 for (int64_t jlo = 0; jlo < m; jlo += jchunk) {
1485 int64_t jhi = std::min(jlo + jchunk, m);
1486 std::vector<Function<T, NDIM>> jvec(g.begin() + jlo, g.begin() + jhi);
1487
1488 Tensor<T> P = matrix_dot(A.get_world(), ivec, jvec, sym);
1489 A.copy_from_replicated_patch(ilo, ihi - 1, jlo, jhi - 1, P);
1490 }
1491 }
1492 return A;
1493 }
1494
1495 /// Computes the matrix dot product of two function vectors - q(i,j) = dot(f[i],g[j])
1496
1497 /// For complex types symmetric is interpreted as Hermitian.
1498 ///
1499 /// The current parallel loop is non-optimal but functional.
1500 template <typename T, typename R, std::size_t NDIM>
1502 const std::vector<Function<T, NDIM>>& f,
1503 const std::vector<Function<R, NDIM>>& g,
1504 bool sym=false)
1505 {
1506 world.gop.fence();
1507 compress(world, f);
1508 // if ((void*)(&f) != (void*)(&g)) compress(world, g);
1509 compress(world, g);
1510
1511 std::vector<const FunctionImpl<T, NDIM>*> left(f.size());
1512 std::vector<const FunctionImpl<R, NDIM>*> right(g.size());
1513 for (unsigned int i = 0; i < f.size(); i++) left[i] = f[i].get_impl().get();
1514 for (unsigned int i = 0; i < g.size(); i++) right[i] = g[i].get_impl().get();
1515
1517
1518 world.gop.fence();
1519 world.gop.sum(r.ptr(), f.size() * g.size());
1520
1521 return r;
1522 }
1523
1524 /// Computes the matrix dot product of two function vectors - q(i,j) = dot(f[i],g[j])
1525
1526 /// For complex types symmetric is interpreted as Hermitian.
1527 ///
1528 /// The current parallel loop is non-optimal but functional.
1529 template <typename T, typename R, std::size_t NDIM>
1531 const std::vector< Function<T,NDIM> >& f,
1532 const std::vector< Function<R,NDIM> >& g,
1533 bool sym=false) {
1534 PROFILE_BLOCK(Vmatrix_dot);
1535 long n=f.size(), m=g.size();
1536 Tensor< TENSOR_RESULT_TYPE(T,R) > r(n,m);
1537 if (sym) MADNESS_ASSERT(n==m);
1538
1539 world.gop.fence();
1540 compress(world, f);
1541 if ((void*)(&f) != (void*)(&g)) compress(world, g);
1542
1543 for (long i=0; i<n; ++i) {
1544 long jtop = m;
1545 if (sym) jtop = i+1;
1546 for (long j=0; j<jtop; ++j) {
1547 if (sym) {
1548 r(j,i) = f[i].dot_local(g[j]);
1549 if (i != j)
1550 r(i,j) = conj(r(j,i));
1551 } else
1552 r(i,j) = f[i].dot_local(g[j]);
1553 }
1554 }
1555
1556 world.gop.fence();
1557 world.gop.sum(r.ptr(),n*m);
1558
1559 return r;
1560 }
1561
1562 /// Multiplies and sums two vectors of functions r = \sum_i a[i] * b[i]
1563 template <typename T, typename R, std::size_t NDIM>
1564 Function<TENSOR_RESULT_TYPE(T,R), NDIM>
1565 dot(World& world,
1566 const std::vector< Function<T,NDIM> >& a,
1567 const std::vector< Function<R,NDIM> >& b,
1568 bool fence=true) {
1569 MADNESS_CHECK(a.size()==b.size());
1570 return sum(world,mul(world,a,b,true),fence);
1571 }
1572
1573
1574
1575 /// out-of-place gaxpy for two vectors: result[i] = alpha * a[i] + beta * b[i]
1576 template <typename T, typename Q, typename R, std::size_t NDIM>
1577 std::vector<Function<TENSOR_RESULT_TYPE(Q,TENSOR_RESULT_TYPE(T,R)),NDIM> >
1579 const std::vector< Function<T,NDIM> >& a,
1580 Q beta,
1581 const std::vector< Function<R,NDIM> >& b,
1582 bool fence=true) {
1583
1584 MADNESS_ASSERT(a.size() == b.size());
1585 typedef TENSOR_RESULT_TYPE(Q,TENSOR_RESULT_TYPE(T,R)) resultT;
1586 if (a.size()==0) return std::vector<Function<resultT,NDIM> >();
1587
1588 auto tensor_type = [](const std::vector<Function<T,NDIM>>& v) {
1589 return v.front().get_impl()->get_tensor_type();
1590 };
1591
1592 // gaxpy can be done either in reconstructed or in compressed state
1593 World& world=a[0].world();
1594 std::vector<Function<resultT,NDIM> > result(a.size());
1595
1596 TreeState operating_state=tensor_type(a)==TT_FULL ? compressed : reconstructed;
1597 try {
1598 ensure_tree_state_respecting_fence(a,operating_state,fence);
1599 ensure_tree_state_respecting_fence(b,operating_state,fence);
1600 } catch (...) {
1601 print("could not respect fence in gaxpy");
1602 change_tree_state(a,operating_state,true);
1603 change_tree_state(b,operating_state,true);
1604 }
1605
1606 if (operating_state==compressed) {
1607 for (unsigned int i=0; i<a.size(); ++i) result[i]=gaxpy_oop(alpha, a[i], beta, b[i], false);
1608 } else {
1609 for (unsigned int i=0; i<a.size(); ++i) result[i]=gaxpy_oop_reconstructed(alpha, a[i], beta, b[i], false);
1610 }
1611
1612 if (fence) world.gop.fence();
1613 return result;
1614 }
1615
1616
1617 /// out-of-place gaxpy for a vectors and a function: result[i] = alpha * a[i] + beta * b
1618 template <typename T, typename Q, typename R, std::size_t NDIM>
1619 std::vector<Function<TENSOR_RESULT_TYPE(Q,TENSOR_RESULT_TYPE(T,R)),NDIM> >
1621 const std::vector< Function<T,NDIM> >& a,
1622 Q beta,
1623 const Function<R,NDIM>& b,
1624 bool fence=true) {
1625
1626 typedef TENSOR_RESULT_TYPE(Q,TENSOR_RESULT_TYPE(T,R)) resultT;
1627 if (a.size()==0) return std::vector<Function<resultT,NDIM> >();
1628
1629 World& world=a[0].world();
1630 try {
1632 // ensure_tree_state_respecting_fence({b},compressed,fence);
1634 } catch (...) {
1635 print("could not respect fence in gaxpy_oop");
1636 compress(world,a);
1637 b.compress();
1638 }
1639 std::vector<Function<resultT,NDIM> > result(a.size());
1640 for (unsigned int i=0; i<a.size(); ++i) {
1641 result[i]=gaxpy_oop(alpha, a[i], beta, b, false);
1642 }
1643 if (fence) world.gop.fence();
1644 return result;
1645 }
1646
1647
1648 /// Generalized A*X+Y for vectors of functions ---- a[i] = alpha*a[i] + beta*b[i]
1649 template <typename T, typename Q, typename R, std::size_t NDIM>
1650 void gaxpy(Q alpha, std::vector<Function<T,NDIM>>& a, Q beta, const std::vector<Function<R,NDIM>>& b, const bool fence) {
1651 if (a.size() == 0) return;
1652 World& world=a.front().world();
1653 gaxpy(world,alpha,a,beta,b,fence);
1654 }
1655
1656 /// Generalized A*X+Y for vectors of functions ---- a[i] = alpha*a[i] + beta*b[i]
1657 template <typename T, typename Q, typename R, std::size_t NDIM>
1658 void gaxpy(World& world,
1659 Q alpha,
1660 std::vector< Function<T,NDIM> >& a,
1661 Q beta,
1662 const std::vector< Function<R,NDIM> >& b,
1663 bool fence=true) {
1664 PROFILE_BLOCK(Vgaxpy);
1665 MADNESS_ASSERT(a.size() == b.size());
1666 if (a.empty()) return;
1667
1668 auto tensor_type = [](const std::vector<Function<T,NDIM>>& v) {
1669 return v.front().get_impl()->get_tensor_type();
1670 };
1671
1672 // gaxpy can be done either in reconstructed or in compressed state
1673 bool do_in_reconstructed_state=tensor_type(a)!=TT_FULL;
1674 TreeState operating_state=do_in_reconstructed_state ? reconstructed : compressed;
1675
1676 if (operating_state==compressed) {
1677 // this is strict: both vectors have to be compressed
1678 try {
1679 ensure_tree_state_respecting_fence(a,operating_state,fence);
1680 ensure_tree_state_respecting_fence(b,operating_state,fence);
1681 } catch (...) {
1682 print("could not respect fence in gaxpy");
1683 change_tree_state(a,operating_state,true);
1684 change_tree_state(b,operating_state,true);
1685 }
1686 MADNESS_CHECK_THROW(get_tree_state(a)==get_tree_state(b),"gaxpy requires same tree state for all functions");
1687 MADNESS_CHECK_THROW(get_tree_state(a)==operating_state,"gaxpy requires reconstructed/compressed tree state for all functions");
1688 } else {
1689 // both vectors can be reconstructed or redundant_after_merge, and they don't have to be the same
1690 TreeState astate=get_tree_state(a);
1691 TreeState bstate=get_tree_state(b);
1692 if (not (astate==reconstructed or astate==redundant_after_merge)) {
1693 try {
1694 ensure_tree_state_respecting_fence(a,operating_state,fence);
1695 } catch (...) {
1696 print("could not respect fence in gaxpy for a");
1698 }
1699 }
1700 if (not (bstate==reconstructed or bstate==redundant_after_merge)) {
1701 try {
1702 ensure_tree_state_respecting_fence(b,operating_state,fence);
1703 } catch (...) {
1704 print("could not respect fence in gaxpy for b");
1706 }
1707 }
1708 }
1709
1710 // finally do the work
1711 for (unsigned int i=0; i<a.size(); ++i) {
1712 a[i].gaxpy(alpha, b[i], beta, false);
1713 }
1714 if (fence and (get_tree_state(a)==redundant_after_merge)) {
1715 for (unsigned int i=0; i<a.size(); ++i) a[i].get_impl()->finalize_sum();
1716 }
1717
1718 if (fence) world.gop.fence();
1719 }
1720
1721
1722 /// Applies a vector of operators to a vector of functions --- q[i] = apply(op[i],f[i])
1723 template <typename opT, typename R, std::size_t NDIM>
1724 std::vector< Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> >
1725 apply(World& world,
1726 const std::vector< std::shared_ptr<opT> >& op,
1727 const std::vector< Function<R,NDIM> > f) {
1728
1729 PROFILE_BLOCK(Vapplyv);
1730 MADNESS_ASSERT(f.size()==op.size());
1731
1732 std::vector< Function<R,NDIM> >& ncf = *const_cast< std::vector< Function<R,NDIM> >* >(&f);
1733
1734// reconstruct(world, f);
1735 make_nonstandard(world, ncf);
1736
1737 std::vector< Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> > result(f.size());
1738 for (unsigned int i=0; i<f.size(); ++i) {
1739 result[i] = apply_only(*op[i], f[i], false);
1740 result[i].get_impl()->set_tree_state(nonstandard_after_apply);
1741 }
1742
1743 world.gop.fence();
1744
1745 standard(world, ncf, false); // restores promise of logical constness
1746 reconstruct(result);
1747 world.gop.fence();
1748
1749 return result;
1750 }
1751
1752
1753 /// Applies an operator to a vector of functions --- q[i] = apply(op,f[i])
1754 template <typename T, typename R, std::size_t NDIM, std::size_t KDIM>
1755 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1757 const std::vector< Function<R,NDIM> > f) {
1758 return apply(op.get_world(),op,f);
1759 }
1760
1761
1762 /// Applies an operator to a vector of functions --- q[i] = apply(op,f[i])
1763 template <typename T, typename R, std::size_t NDIM, std::size_t KDIM>
1764 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
1765 apply(World& world,
1767 const std::vector< Function<R,NDIM> > f) {
1768 PROFILE_BLOCK(Vapply);
1769
1770 std::vector< Function<R,NDIM> >& ncf = *const_cast< std::vector< Function<R,NDIM> >* >(&f);
1771 bool print_timings=(NDIM==6) and (world.rank()==0) and op.print_timings;
1772
1773 double wall0=wall_time();
1774// reconstruct(world, f);
1775 make_nonstandard(world, ncf);
1776 double wall1=wall_time();
1777 if (print_timings) printf("timer: %20.20s %8.2fs\n", "make_nonstandard", wall1-wall0);
1778
1779 std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> > result(f.size());
1780 for (unsigned int i=0; i<f.size(); ++i) {
1781 result[i] = apply_only(op, f[i], false);
1782 }
1783
1784 world.gop.fence();
1785
1786 // restores promise of logical constness
1787 if (op.destructive()) {
1788 for (auto& ff : ncf) ff.clear(false);
1789 world.gop.fence();
1790 } else {
1791 reconstruct(world,f);
1792 }
1793
1794 // svd-tensor requires some cleanup after apply
1795 if (result[0].get_impl()->get_tensor_type()==TT_2D) {
1796 for (auto& r : result) r.get_impl()->finalize_apply();
1797 }
1798
1799 if (print_timings) {
1800 for (auto& r : result) r.get_impl()->print_timer();
1801 op.print_timer();
1802 }
1803 reconstruct(world, result);
1804
1805 return result;
1806 }
1807
1808 /// Normalizes a vector of functions --- v[i] = v[i].scale(1.0/v[i].norm2())
1809 template <typename T, std::size_t NDIM>
1810 void normalize(World& world, std::vector< Function<T,NDIM> >& v, bool fence=true) {
1811 PROFILE_BLOCK(Vnormalize);
1812 std::vector<double> nn = norm2s(world, v);
1813 for (unsigned int i=0; i<v.size(); ++i) v[i].scale(1.0/nn[i],false);
1814 if (fence) world.gop.fence();
1815 }
1816
1817 template <typename T, std::size_t NDIM>
1818 void print_size(World &world, const std::vector<Function<T,NDIM> > &v, const std::string &msg = "vectorfunction" ){
1819 if(v.empty()){
1820 if(world.rank()==0) std::cout << "print_size: " << msg << " is empty" << std::endl;
1821 }else if(v.size()==1){
1822 v.front().print_size(msg);
1823 }else{
1824 for(auto x:v){
1825 x.print_size(msg);
1826 }
1827 }
1828 }
1829
1830 /// return the size of a vector of functions for each rank
1831 template <typename T, std::size_t NDIM>
1832 double get_size_local(World& world, const std::vector< Function<T,NDIM> >& v){
1833 double size=0.0;
1834 for(auto x:v){
1835 if (x.is_initialized()) size+=x.size_local();
1836 }
1837 const double d=sizeof(T);
1838 const double fac=1024*1024*1024;
1839 return size/fac*d;
1840 }
1841
1842 /// return the size of a function for each rank
1843 template <typename T, std::size_t NDIM>
1845 return get_size_local(f.world(),std::vector<Function<T,NDIM> >(1,f));
1846 }
1847
1848
1849 // gives back the size in GB
1850 template <typename T, std::size_t NDIM>
1851 double get_size(World& world, const std::vector< Function<T,NDIM> >& v){
1852
1853 if (v.empty()) return 0.0;
1854
1855 const double d=sizeof(T);
1856 const double fac=1024*1024*1024;
1857
1858 double size=0.0;
1859 for(unsigned int i=0;i<v.size();i++){
1860 if (v[i].is_initialized()) size+=v[i].size();
1861 }
1862
1863 return size/fac*d;
1864
1865 }
1866
1867 // gives back the size in GB
1868 template <typename T, std::size_t NDIM>
1869 double get_size(const Function<T,NDIM> & f){
1870 const double d=sizeof(T);
1871 const double fac=1024*1024*1024;
1872 double size=f.size();
1873 return size/fac*d;
1874 }
1875
1876 /// apply op on the input vector yielding an output vector of functions
1877
1878 /// @param[in] op the operator working on vin
1879 /// @param[in] vin vector of input Functions; needs to be refined to common level!
1880 /// @return vector of output Functions vout = op(vin)
1881 template <typename T, typename opT, std::size_t NDIM>
1882 std::vector<Function<T,NDIM> > multi_to_multi_op_values(const opT& op,
1883 const std::vector< Function<T,NDIM> >& vin,
1884 const bool fence=true) {
1885 MADNESS_ASSERT(vin.size()>0);
1886 MADNESS_ASSERT(vin[0].is_initialized()); // might be changed
1887 World& world=vin[0].world();
1888 Function<T,NDIM> dummy;
1889 dummy.set_impl(vin[0], false);
1890 std::vector<Function<T,NDIM> > vout=zero_functions<T,NDIM>(world, op.get_result_size());
1891 for (auto& out : vout) out.set_impl(vin[0],false);
1892 dummy.multi_to_multi_op_values(op, vin, vout, fence);
1893 return vout;
1894 }
1895
1896
1897
1898
1899 // convenience operators
1900
1901 /// result[i] = a[i] + b[i]
1902 template <typename T, std::size_t NDIM>
1903 std::vector<Function<T,NDIM> > operator+(const std::vector<Function<T,NDIM> >& lhs,
1904 const std::vector<Function<T,NDIM>>& rhs) {
1905 MADNESS_CHECK(lhs.size() == rhs.size());
1906 return gaxpy_oop(1.0,lhs,1.0,rhs);
1907 }
1908
1909 /// result[i] = a[i] - b[i]
1910 template <typename T, std::size_t NDIM>
1911 std::vector<Function<T,NDIM> > operator-(const std::vector<Function<T,NDIM> >& lhs,
1912 const std::vector<Function<T,NDIM> >& rhs) {
1913 MADNESS_CHECK(lhs.size() == rhs.size());
1914 return gaxpy_oop(1.0,lhs,-1.0,rhs);
1915 }
1916
1917 /// result[i] = a[i] + b
1918 template <typename T, std::size_t NDIM>
1919 std::vector<Function<T,NDIM> > operator+(const std::vector<Function<T,NDIM> >& lhs,
1920 const Function<T,NDIM>& rhs) {
1921 // MADNESS_CHECK(lhs.size() == rhs.size()); // no!!
1922 return gaxpy_oop(1.0,lhs,1.0,rhs);
1923 }
1924
1925 /// result[i] = a[i] - b
1926 template <typename T, std::size_t NDIM>
1927 std::vector<Function<T,NDIM> > operator-(const std::vector<Function<T,NDIM> >& lhs,
1928 const Function<T,NDIM>& rhs) {
1929 // MADNESS_CHECK(lhs.size() == rhs.size()); // no
1930 return gaxpy_oop(1.0,lhs,-1.0,rhs);
1931 }
1932
1933 /// result[i] = a + b[i]
1934 template <typename T, std::size_t NDIM>
1935 std::vector<Function<T,NDIM> > operator+(const Function<T,NDIM>& lhs,
1936 const std::vector<Function<T,NDIM> >& rhs) {
1937 // MADNESS_CHECK(lhs.size() == rhs.size()); // no
1938 return gaxpy_oop(1.0,rhs,1.0,lhs);
1939 }
1940
1941 /// result[i] = a - b[i]
1942 template <typename T, std::size_t NDIM>
1943 std::vector<Function<T,NDIM> > operator-(const Function<T,NDIM>& lhs,
1944 const std::vector<Function<T,NDIM> >& rhs) {
1945// MADNESS_CHECK(lhs.size() == rhs.size()); // no
1946 return gaxpy_oop(-1.0,rhs,1.0,lhs);
1947 }
1948
1949
1950 template <typename T, typename R, std::size_t NDIM>
1951 std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> > operator*(const R fac,
1952 const std::vector<Function<T,NDIM> >& rhs) {
1953 if (rhs.size()>0) {
1954 std::vector<Function<T,NDIM> > tmp=copy(rhs[0].world(),rhs);
1955 scale(tmp[0].world(),tmp,TENSOR_RESULT_TYPE(T,R)(fac));
1956 return tmp;
1957 }
1958 return std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> >();
1959 }
1960
1961 template <typename T, typename R, std::size_t NDIM>
1962 std::vector<Function<T,NDIM> > operator*(const std::vector<Function<T,NDIM> >& rhs,
1963 const R fac) {
1964 if (rhs.size()>0) {
1965 std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> > tmp=copy(rhs[0].world(),rhs);
1966 scale(tmp[0].world(),tmp,TENSOR_RESULT_TYPE(T,R)(fac));
1967 return tmp;
1968 }
1969 return std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> >();
1970 }
1971
1972 /// multiply a vector of functions with a function: r[i] = v[i] * a
1973 template <typename T, typename R, std::size_t NDIM>
1975 const std::vector<Function<R,NDIM> >& v) {
1976 if (v.size()>0) return mul(v[0].world(),a,v,true);
1977 return std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> >();
1978 }
1979
1980
1981 /// multiply a vector of functions with a function: r[i] = a * v[i]
1982 template <typename T, typename R, std::size_t NDIM>
1983 std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> > operator*(const std::vector<Function<T,NDIM> >& v,
1984 const Function<R,NDIM>& a) {
1985 if (v.size()>0) return mul(v[0].world(),a,v,true);
1986 return std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> >();
1987 }
1988
1989
1990 template <typename T, std::size_t NDIM>
1991 std::vector<Function<T,NDIM> > operator+=(std::vector<Function<T,NDIM> >& lhs, const std::vector<Function<T,NDIM> >& rhs) {
1992 MADNESS_CHECK(lhs.size() == rhs.size());
1993 if (lhs.size() > 0) gaxpy(lhs.front().world(), 1.0, lhs, 1.0, rhs);
1994 return lhs;
1995 }
1996
1997 template <typename T, std::size_t NDIM>
1998 std::vector<Function<T,NDIM> > operator-=(std::vector<Function<T,NDIM> >& lhs,
1999 const std::vector<Function<T,NDIM> >& rhs) {
2000 MADNESS_CHECK(lhs.size() == rhs.size());
2001 if (lhs.size() > 0) gaxpy(lhs.front().world(), 1.0, lhs, -1.0, rhs);
2002 return lhs;
2003 }
2004
2005 /// return the real parts of the vector's function (if complex)
2006 template <typename T, std::size_t NDIM>
2007 std::vector<Function<typename Tensor<T>::scalar_type,NDIM> >
2008 real(const std::vector<Function<T,NDIM> >& v, bool fence=true) {
2009 std::vector<Function<typename Tensor<T>::scalar_type,NDIM> > result(v.size());
2010 for (std::size_t i=0; i<v.size(); ++i) result[i]=real(v[i],false);
2011 if (fence and result.size()>0) result[0].world().gop.fence();
2012 return result;
2013 }
2014
2015 /// return the imaginary parts of the vector's function (if complex)
2016 template <typename T, std::size_t NDIM>
2017 std::vector<Function<typename Tensor<T>::scalar_type,NDIM> >
2018 imag(const std::vector<Function<T,NDIM> >& v, bool fence=true) {
2019 std::vector<Function<typename Tensor<T>::scalar_type,NDIM> > result(v.size());
2020 for (std::size_t i=0; i<v.size(); ++i) result[i]=imag(v[i],false);
2021 if (fence and result.size()>0) result[0].world().gop.fence();
2022 return result;
2023 }
2024
2025 /// shorthand gradient operator
2026
2027 /// returns the differentiated function f in all NDIM directions
2028 /// @param[in] f the function on which the grad operator works on
2029 /// @param[in] refine refinement before diff'ing makes the result more accurate
2030 /// @param[in] fence fence after completion; if reconstruction is needed always fence
2031 /// @return the vector \frac{\partial}{\partial x_i} f
2032 template <typename T, std::size_t NDIM>
2033 std::vector<Function<T,NDIM> > grad(const Function<T,NDIM>& f,
2034 bool refine=false, bool fence=true) {
2035
2036 World& world=f.world();
2037 f.reconstruct();
2038 if (refine) f.refine(); // refine to make result more precise
2039
2040 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2041 gradient_operator<T,NDIM>(world);
2042
2043 std::vector<Function<T,NDIM> > result(NDIM);
2044 for (size_t i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),f,false);
2045 if (fence) world.gop.fence();
2046 return result;
2047 }
2048
2049 // BLM first derivative
2050 template <typename T, std::size_t NDIM>
2051 std::vector<Function<T,NDIM> > grad_ble_one(const Function<T,NDIM>& f,
2052 bool refine=false, bool fence=true) {
2053
2054 World& world=f.world();
2055 f.reconstruct();
2056 if (refine) f.refine(); // refine to make result more precise
2057
2058 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2059 gradient_operator<T,NDIM>(world);
2060
2061 // Read in new coeff for each operator
2062 for (unsigned int i=0; i<NDIM; ++i) (*grad[i]).set_ble1();
2063
2064 std::vector<Function<T,NDIM> > result(NDIM);
2065 for (unsigned int i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),f,false);
2066 if (fence) world.gop.fence();
2067 return result;
2068 }
2069
2070 // BLM second derivative
2071 template <typename T, std::size_t NDIM>
2072 std::vector<Function<T,NDIM> > grad_ble_two(const Function<T,NDIM>& f,
2073 bool refine=false, bool fence=true) {
2074
2075 World& world=f.world();
2076 f.reconstruct();
2077 if (refine) f.refine(); // refine to make result more precise
2078
2079 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2080 gradient_operator<T,NDIM>(world);
2081
2082 // Read in new coeff for each operator
2083 for (unsigned int i=0; i<NDIM; ++i) (*grad[i]).set_ble2();
2084
2085 std::vector<Function<T,NDIM> > result(NDIM);
2086 for (unsigned int i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),f,false);
2087 if (fence) world.gop.fence();
2088 return result;
2089 }
2090
2091 // Bspline first derivative
2092 template <typename T, std::size_t NDIM>
2093 std::vector<Function<T,NDIM> > grad_bspline_one(const Function<T,NDIM>& f,
2094 bool refine=false, bool fence=true) {
2095
2096 World& world=f.world();
2097 f.reconstruct();
2098 if (refine) f.refine(); // refine to make result more precise
2099
2100 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2101 gradient_operator<T,NDIM>(world);
2102
2103 // Read in new coeff for each operator
2104 for (unsigned int i=0; i<NDIM; ++i) (*grad[i]).set_bspline1();
2105
2106 std::vector<Function<T,NDIM> > result(NDIM);
2107 for (unsigned int i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),f,false);
2108 if (fence) world.gop.fence();
2109 return result;
2110 }
2111
2112 // Bpsline second derivative
2113 template <typename T, std::size_t NDIM>
2114 std::vector<Function<T,NDIM> > grad_bpsline_two(const Function<T,NDIM>& f,
2115 bool refine=false, bool fence=true) {
2116
2117 World& world=f.world();
2118 f.reconstruct();
2119 if (refine) f.refine(); // refine to make result more precise
2120
2121 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2122 gradient_operator<T,NDIM>(world);
2123
2124 // Read in new coeff for each operator
2125 for (unsigned int i=0; i<NDIM; ++i) (*grad[i]).set_bspline2();
2126
2127 std::vector<Function<T,NDIM> > result(NDIM);
2128 for (unsigned int i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),f,false);
2129 if (fence) world.gop.fence();
2130 return result;
2131 }
2132
2133 // Bspline third derivative
2134 template <typename T, std::size_t NDIM>
2135 std::vector<Function<T,NDIM> > grad_bspline_three(const Function<T,NDIM>& f,
2136 bool refine=false, bool fence=true) {
2137
2138 World& world=f.world();
2139 f.reconstruct();
2140 if (refine) f.refine(); // refine to make result more precise
2141
2142 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2143 gradient_operator<T,NDIM>(world);
2144
2145 // Read in new coeff for each operator
2146 for (unsigned int i=0; i<NDIM; ++i) (*grad[i]).set_bspline3();
2147
2148 std::vector<Function<T,NDIM> > result(NDIM);
2149 for (unsigned int i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),f,false);
2150 if (fence) world.gop.fence();
2151 return result;
2152 }
2153
2154
2155
2156 /// shorthand div operator
2157
2158 /// returns the dot product of nabla with a vector f
2159 /// @param[in] f the vector of functions on which the div operator works on
2160 /// @param[in] refine refinement before diff'ing makes the result more accurate
2161 /// @param[in] fence fence after completion; currently always fences
2162 /// @return the vector \frac{\partial}{\partial x_i} f
2163 /// TODO: add this to operator fusion
2164 template <typename T, std::size_t NDIM>
2166 bool do_refine=false, bool fence=true) {
2167
2168 MADNESS_ASSERT(v.size()>0);
2169 World& world=v[0].world();
2170 reconstruct(world,v);
2171 if (do_refine) refine(world,v); // refine to make result more precise
2172
2173 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2174 gradient_operator<T,NDIM>(world);
2175
2176 std::vector<Function<T,NDIM> > result(NDIM);
2177 for (size_t i=0; i<NDIM; ++i) result[i]=apply(*(grad[i]),v[i],false);
2178 world.gop.fence();
2179 return sum(world,result,fence);
2180 }
2181
2182 /// shorthand rot operator
2183
2184 /// returns the cross product of nabla with a vector f
2185 /// @param[in] f the vector of functions on which the rot operator works on
2186 /// @param[in] refine refinement before diff'ing makes the result more accurate
2187 /// @param[in] fence fence after completion; currently always fences
2188 /// @return the vector \frac{\partial}{\partial x_i} f
2189 /// TODO: add this to operator fusion
2190 template <typename T, std::size_t NDIM>
2191 std::vector<Function<T,NDIM> > rot(const std::vector<Function<T,NDIM> >& v,
2192 bool do_refine=false, bool fence=true) {
2193
2194 MADNESS_ASSERT(v.size()==3);
2195 World& world=v[0].world();
2196 reconstruct(world,v);
2197 if (do_refine) refine(world,v); // refine to make result more precise
2198
2199 std::vector< std::shared_ptr< Derivative<T,NDIM> > > grad=
2200 gradient_operator<T,NDIM>(world);
2201
2202 std::vector<Function<T,NDIM> > d(NDIM),dd(NDIM);
2203 d[0]=apply(*(grad[1]),v[2],false); // Dy z
2204 d[1]=apply(*(grad[2]),v[0],false); // Dz x
2205 d[2]=apply(*(grad[0]),v[1],false); // Dx y
2206 dd[0]=apply(*(grad[2]),v[1],false); // Dz y
2207 dd[1]=apply(*(grad[0]),v[2],false); // Dx z
2208 dd[2]=apply(*(grad[1]),v[0],false); // Dy x
2209 world.gop.fence();
2210
2211 compress(world,d,false);
2212 compress(world,dd,false);
2213 world.gop.fence();
2214 d[0].gaxpy(1.0,dd[0],-1.0,false);
2215 d[1].gaxpy(1.0,dd[1],-1.0,false);
2216 d[2].gaxpy(1.0,dd[2],-1.0,false);
2217
2218 world.gop.fence();
2219 reconstruct(d);
2220 return d;
2221 }
2222
2223 /// shorthand cross operator
2224
2225 /// returns the cross product of vectors f and g
2226 /// @param[in] f the vector of functions on which the rot operator works on
2227 /// @param[in] g the vector of functions on which the rot operator works on
2228 /// @param[in] fence fence after completion; currently always fences
2229 /// @return the vector \frac{\partial}{\partial x_i} f
2230 /// TODO: add this to operator fusion
2231 template <typename T, typename R, std::size_t NDIM>
2232 std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> > cross(const std::vector<Function<T,NDIM> >& f,
2233 const std::vector<Function<R,NDIM> >& g,
2234 bool do_refine=false, bool fence=true) {
2235
2236 MADNESS_ASSERT(f.size()==3);
2237 MADNESS_ASSERT(g.size()==3);
2238 World& world=f[0].world();
2239 reconstruct(world,f,false);
2240 reconstruct(world,g);
2241
2242 std::vector<Function<TENSOR_RESULT_TYPE(T,R),NDIM> > d(f.size()),dd(f.size());
2243
2244 d[0]=mul(f[1],g[2],false);
2245 d[1]=mul(f[2],g[0],false);
2246 d[2]=mul(f[0],g[1],false);
2247
2248 dd[0]=mul(f[2],g[1],false);
2249 dd[1]=mul(f[0],g[2],false);
2250 dd[2]=mul(f[1],g[0],false);
2251 world.gop.fence();
2252
2253 compress(world,d,false);
2254 compress(world,dd);
2255
2256 d[0].gaxpy(1.0,dd[0],-1.0,false);
2257 d[1].gaxpy(1.0,dd[1],-1.0,false);
2258 d[2].gaxpy(1.0,dd[2],-1.0,false);
2259
2260
2261 world.gop.fence();
2262 return d;
2263 }
2264
2265 template<typename T, std::size_t NDIM>
2266 void load_balance(World& world, std::vector<Function<T,NDIM> >& vf) {
2267
2268 struct LBCost {
2269 LBCost() = default;
2270 double operator()(const Key<NDIM>& key, const FunctionNode<T,NDIM>& node) const {
2271 return node.coeff().size();
2272 }
2273 };
2274
2275 LoadBalanceDeux<6> lb(world);
2276 for (const auto& f : vf) lb.add_tree(f, LBCost());
2278
2279 }
2280
2281 /// load a vector of functions
2282 template<typename T, size_t NDIM>
2283 void load_function(World& world, std::vector<Function<T,NDIM> >& f,
2284 const std::string name) {
2285 if (world.rank()==0) print("loading vector of functions",name);
2287 std::size_t fsize=0;
2288 ar & fsize;
2289 f.resize(fsize);
2290 for (std::size_t i=0; i<fsize; ++i) ar & f[i];
2291 }
2292
2293 /// save a vector of functions
2294 template<typename T, size_t NDIM>
2295 void save_function(const std::vector<Function<T,NDIM> >& f, const std::string name) {
2296 if (f.size()>0) {
2297 World& world=f.front().world();
2298 if (world.rank()==0) print("saving vector of functions",name);
2300 std::size_t fsize=f.size();
2301 ar & fsize;
2302 for (std::size_t i=0; i<fsize; ++i) ar & f[i];
2303 }
2304 }
2305
2306
2307}
2308#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:414
static TensorType get_tensor_type()
Returns the default tensor type.
Definition funcdefaults.h:324
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:6000
void undo_redundant(const bool fence)
convert this from redundant to standard reconstructed form
Definition mraimpl.h:1534
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:3580
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:1403
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:6052
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:688
Function< T, NDIM > & gaxpy(const T &alpha, const Function< Q, NDIM > &other, const R &beta, bool fence=true)
Inplace, general bi-linear operation in wavelet basis. No communication except for optional fence.
Definition mra.h:1025
void set_impl(const std::shared_ptr< FunctionImpl< T, NDIM > > &impl)
Replace current FunctionImpl with provided new one.
Definition mra.h:661
void multi_to_multi_op_values(const opT &op, const std::vector< Function< T, NDIM > > &vin, std::vector< Function< T, NDIM > > &vout, const bool fence=true)
apply op on the input vector yielding an output vector of functions
Definition mra.h:1585
bool is_initialized() const
Returns true if the function is initialized.
Definition mra.h:167
long size() const
Definition lowranktensor.h:482
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:1825
A simple, fixed dimension vector.
Definition vector.h:64
Interface to be provided by any process map.
Definition worlddc.h:82
void fence(bool debug=false)
Synchronizes all processes in communicator AND globally ensures no pending AM or tasks.
Definition worldgop.cc:161
void sum(T *buf, size_t nelem)
Inplace global sum while still processing AM & tasks.
Definition worldgop.h:870
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
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.
static double lo
Definition dirac-hatom.cc:23
auto T(World &world, response_space &f) -> response_space
Definition global_functions.cc:34
const double beta
Definition gygi_soltion.cc:62
static const double v
Definition hatom_sf_dirac.cc:20
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:2295
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:612
Function< double, NDIM > abssq(const Function< double_complex, NDIM > &z, bool fence=true)
Returns a new function that is the square of the absolute value of the input.
Definition mra.h:2778
Function< TENSOR_RESULT_TYPE(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:1981
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:2746
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:2035
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:858
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:845
std::vector< Function< T, NDIM > > grad_bspline_one(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
Definition vmra.h:2093
void set_impl(std::vector< Function< T, NDIM > > &v, const std::vector< std::shared_ptr< FunctionImpl< T, NDIM > > > vimpl)
Definition vmra.h:688
Function< Q, NDIM > convert(const Function< T, NDIM > &f, bool fence=true)
Type conversion implies a deep copy. No communication except for optional fence.
Definition mra.h:2096
Function< TENSOR_RESULT_TYPE(Q, T), NDIM > mul(const Q alpha, const Function< T, NDIM > &f, bool fence=true)
Returns new function equal to alpha*f(x) with optional fence.
Definition mra.h:1765
void truncate(World &world, response_space &v, double tol, bool fence)
Definition basic_operators.cc:30
std::vector< std::shared_ptr< FunctionImpl< T, NDIM > > > get_impl(const std::vector< Function< T, NDIM > > &v)
Definition vmra.h:681
Function< T, NDIM > div(const std::vector< Function< T, NDIM > > &v, bool do_refine=false, bool fence=true)
shorthand div operator
Definition vmra.h:2165
std::vector< Function< T, NDIM > > orthonormalize_cd(const std::vector< Function< T, NDIM > > &v, Tensor< T > &ovlp)
Definition vmra.h:575
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:1181
tensorT Q2(const tensorT &s)
Given overlap matrix, return rotation with 2nd order error to orthonormalize the vectors.
Definition SCF.cc:135
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:707
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:2232
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:2110
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:1149
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:694
std::vector< Function< T, NDIM > > grad_bpsline_two(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
Definition vmra.h:2114
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:1285
double norm2(World &world, const std::vector< Function< T, NDIM > > &v)
Computes the 2-norm of a vector of functions.
Definition vmra.h:871
std::vector< Function< T, NDIM > > flatten(const std::vector< std::vector< Function< T, NDIM > > > &vv)
Definition vmra.h:674
std::vector< Function< T, NDIM > > orthonormalize_canonical(const std::vector< Function< T, NDIM > > &v, const Tensor< T > &ovlp, double lindep)
Definition vmra.h:519
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:72
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:1882
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:2794
std::vector< CCPairFunction< T, NDIM > > operator-(const std::vector< CCPairFunction< T, NDIM > > c1, const std::vector< CCPairFunction< T, NDIM > > &c2)
Definition ccpairfunction.h:1060
Function< TENSOR_RESULT_TYPE(L, R), NDIM > mul_sparse(const Function< L, NDIM > &left, const Function< R, NDIM > &right, double tol, bool fence=true)
Sparse multiplication — left and right must be reconstructed and if tol!=0 have tree of norms already...
Definition mra.h:1806
std::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:1220
Function< T, NDIM > gaxpy_oop_reconstructed(const double alpha, const Function< T, NDIM > &left, const double beta, const Function< T, NDIM > &right, const bool fence=true)
Returns new function alpha*left + beta*right optional fence, having both addends reconstructed.
Definition mra.h:1999
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:736
void print(const T &t, const Ts &... ts)
Print items to std::cout (items separated by spaces) and terminate with a new line.
Definition print.h:225
response_space apply(World &world, std::vector< std::vector< std::shared_ptr< real_convolution_3d > > > &op, response_space &f)
Definition basic_operators.cc:39
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:667
@ 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:1818
NDIM & f
Definition mra.h:2481
Function< TENSOR_RESULT_TYPE(L, R), NDIM > add(const Function< L, NDIM > &left, const Function< R, NDIM > &right, bool fence=true)
Same as operator+ but with optional fence and no automatic compression.
Definition mra.h:1991
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:2807
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:1565
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:2481
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:2051
Function< TENSOR_RESULT_TYPE(typename opT::opT, R), NDIM > apply_only(const opT &op, const Function< R, NDIM > &f, bool fence=true)
Apply operator ONLY in non-standard form - required other steps missing !!
Definition mra.h:2184
std::vector< Function< T, NDIM > > grad_bspline_three(const Function< T, NDIM > &f, bool refine=false, bool fence=true)
Definition vmra.h:2135
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:442
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:2283
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:1871
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:2072
static bool print_timings
Definition SCF.cc:104
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:1810
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:2033
Function< T, NDIM > multiply(const Function< T, NDIM > f, const Function< T, LDIM > g, const int particle, const bool fence=true)
multiply a high-dimensional function with a low-dimensional function
Definition mra.h:2437
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:1467
std::vector< Function< T, NDIM > > orthonormalize_symmetric(const std::vector< Function< T, NDIM > > &v, const Tensor< T > &ovlp)
symmetric orthonormalization (see e.g. Szabo/Ostlund)
Definition vmra.h:474
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
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:1530
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:1851
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:1832
Function< T, NDIM > copy(const Function< T, NDIM > &f, const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > &pmap, bool fence=true)
Create a new copy of the function with different distribution and optional fence.
Definition mra.h:2066
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:993
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:245
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:332
Definition dirac-hatom.cc:108
std::string ok(const bool b)
Definition test6.cc:43
AtomicInt sum
Definition test_atomicint.cc:46
int P
Definition test_binsorter.cc:9
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