MADNESS  version 0.9
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 
122 #include <madness/mra/mra.h>
123 #include <madness/mra/derivative.h>
124 #include <cstdio>
125 
126 namespace madness {
127 
129  template <typename T, std::size_t NDIM>
130  void compress(World& world,
131  const std::vector< Function<T,NDIM> >& v,
132  bool fence=true) {
133 
134  PROFILE_BLOCK(Vcompress);
135  bool must_fence = false;
136  for (unsigned int i=0; i<v.size(); ++i) {
137  if (!v[i].is_compressed()) {
138  v[i].compress(false);
139  must_fence = true;
140  }
141  }
142 
143  if (fence && must_fence) world.gop.fence();
144  }
145 
146 
148  template <typename T, std::size_t NDIM>
149  void reconstruct(World& world,
150  const std::vector< Function<T,NDIM> >& v,
151  bool fence=true) {
152  PROFILE_BLOCK(Vreconstruct);
153  bool must_fence = false;
154  for (unsigned int i=0; i<v.size(); ++i) {
155  if (v[i].is_compressed()) {
156  v[i].reconstruct(false);
157  must_fence = true;
158  }
159  }
160 
161  if (fence && must_fence) world.gop.fence();
162  }
163 
164 
166  template <typename T, std::size_t NDIM>
167  void nonstandard(World& world,
168  std::vector< Function<T,NDIM> >& v,
169  bool fence=true) {
170  PROFILE_BLOCK(Vnonstandard);
171  reconstruct(world, v);
172  for (unsigned int i=0; i<v.size(); ++i) {
173  v[i].nonstandard(false,false);
174  }
175  if (fence) world.gop.fence();
176  }
177 
178 
180  template <typename T, std::size_t NDIM>
181  void standard(World& world,
182  std::vector< Function<T,NDIM> >& v,
183  bool fence=true) {
184  PROFILE_BLOCK(Vstandard);
185  for (unsigned int i=0; i<v.size(); ++i) {
186  v[i].standard(false);
187  }
188  if (fence) world.gop.fence();
189  }
190 
191 
193  template <typename T, std::size_t NDIM>
194  void truncate(World& world,
195  std::vector< Function<T,NDIM> >& v,
196  double tol=0.0,
197  bool fence=true) {
198  PROFILE_BLOCK(Vtruncate);
199 
200  compress(world, v);
201 
202  for (unsigned int i=0; i<v.size(); ++i) {
203  v[i].truncate(tol, false);
204  }
205 
206  if (fence) world.gop.fence();
207  }
208 
210  template <typename T, std::size_t NDIM>
211  std::vector< Function<T,NDIM> >
212  apply(World& world,
213  const Derivative<T,NDIM>& D,
214  const std::vector< Function<T,NDIM> >& v,
215  bool fence=true)
216  {
217  reconstruct(world, v);
218  std::vector< Function<T,NDIM> > df(v.size());
219  for (unsigned int i=0; i<v.size(); ++i) {
220  df[i] = D(v[i],false);
221  }
222  if (fence) world.gop.fence();
223  return df;
224  }
225 
227  template <typename T, std::size_t NDIM>
228  std::vector< Function<T,NDIM> >
229  zero_functions(World& world, int n) {
230  PROFILE_BLOCK(Vzero_functions);
231  std::vector< Function<T,NDIM> > r(n);
232  for (int i=0; i<n; ++i)
234 
235  return r;
236  }
237 
238 
240 
243  template <typename T, typename R, std::size_t NDIM>
244  std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> >
245  transform(World& world,
246  const std::vector< Function<T,NDIM> >& v,
247  const Tensor<R>& c,
248  bool fence=true) {
249 
250  PROFILE_BLOCK(Vtransformsp);
251  typedef TENSOR_RESULT_TYPE(T,R) resultT;
252  int n = v.size(); // n is the old dimension
253  int m = c.dim(1); // m is the new dimension
254  MADNESS_ASSERT(n==c.dim(0));
255 
256  std::vector< Function<resultT,NDIM> > vc = zero_functions<resultT,NDIM>(world, m);
257  compress(world, v);
258  compress(world, vc);
259 
260  for (int i=0; i<m; ++i) {
261  for (int j=0; j<n; ++j) {
262  if (c(j,i) != R(0.0)) vc[i].gaxpy(1.0,v[j],c(j,i),false);
263  }
264  }
265 
266  if (fence) world.gop.fence();
267  return vc;
268  }
269 
270 
271  template <typename L, typename R, std::size_t NDIM>
272  std::vector< Function<TENSOR_RESULT_TYPE(L,R),NDIM> >
273  transform(World& world, const std::vector< Function<L,NDIM> >& v, const Tensor<R>& c, double tol, bool fence) {
274  PROFILE_BLOCK(Vtransform);
275  MADNESS_ASSERT(v.size() == (unsigned int)(c.dim(0)));
276 
277  std::vector< Function<TENSOR_RESULT_TYPE(L,R),NDIM> > vresult(c.dim(1));
278  for (int i=0; i<c.dim(1); ++i) {
280  }
281  compress(world, v, false);
282  compress(world, vresult, false);
283  world.gop.fence();
284  vresult[0].vtransform(v, c, vresult, tol, fence);
285  return vresult;
286  }
287 
289  template <typename T, typename Q, std::size_t NDIM>
290  void scale(World& world,
291  std::vector< Function<T,NDIM> >& v,
292  const std::vector<Q>& factors,
293  bool fence=true) {
294  PROFILE_BLOCK(Vscale);
295  for (unsigned int i=0; i<v.size(); ++i) v[i].scale(factors[i],false);
296  if (fence) world.gop.fence();
297  }
298 
300  template <typename T, typename Q, std::size_t NDIM>
301  void scale(World& world,
302  std::vector< Function<T,NDIM> >& v,
303  const Q factor,
304  bool fence=true) {
305  PROFILE_BLOCK(Vscale);
306  for (unsigned int i=0; i<v.size(); ++i) v[i].scale(factor,false);
307  if (fence) world.gop.fence();
308  }
309 
311  template <typename T, std::size_t NDIM>
312  std::vector<double> norm2s(World& world,
313  const std::vector< Function<T,NDIM> >& v) {
314  PROFILE_BLOCK(Vnorm2);
315  std::vector<double> norms(v.size());
316  for (unsigned int i=0; i<v.size(); ++i) norms[i] = v[i].norm2sq_local();
317  world.gop.sum(&norms[0], norms.size());
318  for (unsigned int i=0; i<v.size(); ++i) norms[i] = sqrt(norms[i]);
319  world.gop.fence();
320  return norms;
321  }
322 
324  template <typename T, std::size_t NDIM>
325  double norm2(World& world,
326  const std::vector< Function<T,NDIM> >& v) {
327  PROFILE_BLOCK(Vnorm2);
328  std::vector<double> norms(v.size());
329  for (unsigned int i=0; i<v.size(); ++i) norms[i] = v[i].norm2sq_local();
330  world.gop.sum(&norms[0], norms.size());
331  for (unsigned int i=1; i<v.size(); ++i) norms[0] += norms[i];
332  world.gop.fence();
333  return sqrt(norms[0]);
334  }
335 
336  inline double conj(double x) {
337  return x;
338  }
339 
340  inline double conj(float x) {
341  return x;
342  }
343 
344 // !!! FIXME: this task is broken because FunctionImpl::inner_local forces a
345 // future on return from WorldTaskQueue::reduce, which will causes a deadlock if
346 // run inside a task. This behavior must be changed before this task can be used
347 // again.
348 //
349 // template <typename T, typename R, std::size_t NDIM>
350 // struct MatrixInnerTask : public TaskInterface {
351 // Tensor<TENSOR_RESULT_TYPE(T,R)> result; // Must be a copy
352 // const Function<T,NDIM>& f;
353 // const std::vector< Function<R,NDIM> >& g;
354 // long jtop;
355 //
356 // MatrixInnerTask(const Tensor<TENSOR_RESULT_TYPE(T,R)>& result,
357 // const Function<T,NDIM>& f,
358 // const std::vector< Function<R,NDIM> >& g,
359 // long jtop)
360 // : result(result), f(f), g(g), jtop(jtop) {}
361 //
362 // void run(World& world) {
363 // for (long j=0; j<jtop; ++j) {
364 // result(j) = f.inner_local(g[j]);
365 // }
366 // }
367 //
368 // private:
369 // /// Get the task id
370 //
371 // /// \param id The id to set for this task
372 // virtual void get_id(std::pair<void*,unsigned short>& id) const {
373 // PoolTaskInterface::make_id(id, *this);
374 // }
375 // }; // struct MatrixInnerTask
376 
377 
379 
383  template <typename T, typename R, std::size_t NDIM>
384  Tensor< TENSOR_RESULT_TYPE(T,R) > matrix_inner(World& world,
385  const std::vector< Function<T,NDIM> >& f,
386  const std::vector< Function<R,NDIM> >& g,
387  bool sym=false) {
388  PROFILE_BLOCK(Vmatrix_inner);
389  long n=f.size(), m=g.size();
390  Tensor< TENSOR_RESULT_TYPE(T,R) > r(n,m);
391  if (sym) MADNESS_ASSERT(n==m);
392 
393  world.gop.fence();
394  compress(world, f);
395  if (&f != &g) compress(world, g);
396 
397  for (long i=0; i<n; ++i) {
398  long jtop = m;
399  if (sym) jtop = i+1;
400  for (long j=0; j<jtop; ++j) {
401  r(i,j) = f[i].inner_local(g[j]);
402  if (sym) r(j,i) = conj(r(i,j));
403  }
404  }
405 
406 // for (long i=n-1; i>=0; --i) {
407 // long jtop = m;
408 // if (sym) jtop = i+1;
409 // world.taskq.add(new MatrixInnerTask<T,R,NDIM>(r(i,_), f[i], g, jtop));
410 // }
411  world.gop.fence();
412  world.gop.sum(r.ptr(),n*m);
413 
414 // if (sym) {
415 // for (int i=0; i<n; ++i) {
416 // for (int j=0; j<i; ++j) {
417 // r(j,i) = conj(r(i,j));
418 // }
419 // }
420 // }
421  return r;
422  }
423 
425  template <typename T, typename R, std::size_t NDIM>
426  Tensor< TENSOR_RESULT_TYPE(T,R) > inner(World& world,
427  const std::vector< Function<T,NDIM> >& f,
428  const std::vector< Function<R,NDIM> >& g) {
429  PROFILE_BLOCK(Vinnervv);
430  long n=f.size(), m=g.size();
431  MADNESS_ASSERT(n==m);
432  Tensor< TENSOR_RESULT_TYPE(T,R) > r(n);
433 
434  compress(world, f);
435  compress(world, g);
436 
437  for (long i=0; i<n; ++i) {
438  r(i) = f[i].inner_local(g[i]);
439  }
440 
441  world.taskq.fence();
442  world.gop.sum(r.ptr(),n);
443  world.gop.fence();
444  return r;
445  }
446 
447 
449  template <typename T, typename R, std::size_t NDIM>
450  Tensor< TENSOR_RESULT_TYPE(T,R) > inner(World& world,
451  const Function<T,NDIM>& f,
452  const std::vector< Function<R,NDIM> >& g) {
453  PROFILE_BLOCK(Vinner);
454  long n=g.size();
455  Tensor< TENSOR_RESULT_TYPE(T,R) > r(n);
456 
457  f.compress();
458  compress(world, g);
459 
460  for (long i=0; i<n; ++i) {
461  r(i) = f.inner_local(g[i]);
462  }
463 
464  world.taskq.fence();
465  world.gop.sum(r.ptr(),n);
466  world.gop.fence();
467  return r;
468  }
469 
470 
472  template <typename T, typename R, std::size_t NDIM>
473  std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
474  mul(World& world,
475  const Function<T,NDIM>& a,
476  const std::vector< Function<R,NDIM> >& v,
477  bool fence=true) {
478  PROFILE_BLOCK(Vmul);
479  a.reconstruct(false);
480  reconstruct(world, v, false);
481  world.gop.fence();
482  return vmulXX(a, v, 0.0, fence);
483  }
484 
486  template <typename T, typename R, std::size_t NDIM>
487  std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
489  const Function<T,NDIM>& a,
490  const std::vector< Function<R,NDIM> >& v,
491  double tol,
492  bool fence=true) {
493  PROFILE_BLOCK(Vmulsp);
494  a.reconstruct(false);
495  reconstruct(world, v, false);
496  world.gop.fence();
497  for (unsigned int i=0; i<v.size(); ++i) {
498  v[i].norm_tree(false);
499  }
500  a.norm_tree();
501  return vmulXX(a, v, tol, fence);
502  }
503 
505  template <typename T, std::size_t NDIM>
506  void norm_tree(World& world,
507  const std::vector< Function<T,NDIM> >& v,
508  bool fence=true)
509  {
510  PROFILE_BLOCK(Vnorm_tree);
511  for (unsigned int i=0; i<v.size(); ++i) {
512  v[i].norm_tree(false);
513  }
514  if (fence) world.gop.fence();
515  }
516 
518  template <typename T, typename R, std::size_t NDIM>
519  std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
520  mul(World& world,
521  const std::vector< Function<T,NDIM> >& a,
522  const std::vector< Function<R,NDIM> >& b,
523  bool fence=true) {
524  PROFILE_BLOCK(Vmulvv);
525  reconstruct(world, a, false);
526  if (&a != &b) reconstruct(world, b, false);
527  world.gop.fence();
528 
529  std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> > q(a.size());
530  for (unsigned int i=0; i<a.size(); ++i) {
531  q[i] = mul(a[i], b[i], false);
532  }
533  if (fence) world.gop.fence();
534  return q;
535  }
536 
537 
539  template <typename T, std::size_t NDIM>
540  std::vector< Function<T,NDIM> >
541  square(World& world,
542  const std::vector< Function<T,NDIM> >& v,
543  bool fence=true) {
544  return mul<T,T,NDIM>(world, v, v, fence);
545 // std::vector< Function<T,NDIM> > vsq(v.size());
546 // for (unsigned int i=0; i<v.size(); ++i) {
547 // vsq[i] = square(v[i], false);
548 // }
549 // if (fence) world.gop.fence();
550 // return vsq;
551  }
552 
553 
555  template <typename T, std::size_t NDIM>
556  void set_thresh(World& world, std::vector< Function<T,NDIM> >& v, double thresh, bool fence=true) {
557  for (unsigned int j=0; j<v.size(); ++j) {
558  v[j].set_thresh(thresh,false);
559  }
560  if (fence) world.gop.fence();
561  }
562 
564  template <typename T, std::size_t NDIM>
565  std::vector< Function<T,NDIM> >
566  conj(World& world,
567  const std::vector< Function<T,NDIM> >& v,
568  bool fence=true) {
569  PROFILE_BLOCK(Vconj);
570  std::vector< Function<T,NDIM> > r = copy(world, v); // Currently don't have oop conj
571  for (unsigned int i=0; i<v.size(); ++i) {
572  r[i].conj(false);
573  }
574  if (fence) world.gop.fence();
575  return r;
576  }
577 
578 
580  template <typename T, std::size_t NDIM>
581  std::vector< Function<T,NDIM> >
582  copy(World& world,
583  const std::vector< Function<T,NDIM> >& v,
584  bool fence=true) {
585  PROFILE_BLOCK(Vcopy);
586  std::vector< Function<T,NDIM> > r(v.size());
587  for (unsigned int i=0; i<v.size(); ++i) {
588  r[i] = copy(v[i], false);
589  }
590  if (fence) world.gop.fence();
591  return r;
592  }
593 
595  template <typename T, std::size_t NDIM>
596  std::vector< Function<T,NDIM> >
597  copy(World& world,
598  const Function<T,NDIM>& v,
599  const unsigned int n,
600  bool fence=true) {
601  PROFILE_BLOCK(Vcopy1);
602  std::vector< Function<T,NDIM> > r(n);
603  for (unsigned int i=0; i<n; ++i) {
604  r[i] = copy(v, false);
605  }
606  if (fence) world.gop.fence();
607  return r;
608  }
609 
611  template <typename T, typename R, std::size_t NDIM>
612  std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
613  add(World& world,
614  const std::vector< Function<T,NDIM> >& a,
615  const std::vector< Function<R,NDIM> >& b,
616  bool fence=true) {
617  PROFILE_BLOCK(Vadd);
618  MADNESS_ASSERT(a.size() == b.size());
619  compress(world, a);
620  compress(world, b);
621 
622  std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> > r(a.size());
623  for (unsigned int i=0; i<a.size(); ++i) {
624  r[i] = add(a[i], b[i], false);
625  }
626  if (fence) world.gop.fence();
627  return r;
628  }
629 
631  template <typename T, typename R, std::size_t NDIM>
632  std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
633  add(World& world,
634  const Function<T,NDIM> & a,
635  const std::vector< Function<R,NDIM> >& b,
636  bool fence=true) {
637  PROFILE_BLOCK(Vadd1);
638  a.compress();
639  compress(world, b);
640 
641  std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> > r(b.size());
642  for (unsigned int i=0; i<b.size(); ++i) {
643  r[i] = add(a, b[i], false);
644  }
645  if (fence) world.gop.fence();
646  return r;
647  }
648  template <typename T, typename R, std::size_t NDIM>
649  inline std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
650  add(World& world,
651  const std::vector< Function<R,NDIM> >& b,
652  const Function<T,NDIM> & a,
653  bool fence=true) {
654  return add(world, a, b, fence);
655  }
656 
658  template <typename T, typename R, std::size_t NDIM>
659  std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
660  sub(World& world,
661  const std::vector< Function<T,NDIM> >& a,
662  const std::vector< Function<R,NDIM> >& b,
663  bool fence=true) {
664  PROFILE_BLOCK(Vsub);
665  MADNESS_ASSERT(a.size() == b.size());
666  compress(world, a);
667  compress(world, b);
668 
669  std::vector< Function<TENSOR_RESULT_TYPE(T,R),NDIM> > r(a.size());
670  for (unsigned int i=0; i<a.size(); ++i) {
671  r[i] = sub(a[i], b[i], false);
672  }
673  if (fence) world.gop.fence();
674  return r;
675  }
676 
677 
679  template <typename T, typename Q, typename R, std::size_t NDIM>
680  void gaxpy(World& world,
681  Q alpha,
682  std::vector< Function<T,NDIM> >& a,
683  Q beta,
684  const std::vector< Function<R,NDIM> >& b,
685  bool fence=true) {
686  PROFILE_BLOCK(Vgaxpy);
687  MADNESS_ASSERT(a.size() == b.size());
688  compress(world, a);
689  compress(world, b);
690 
691  for (unsigned int i=0; i<a.size(); ++i) {
692  a[i].gaxpy(alpha, b[i], beta, false);
693  }
694  if (fence) world.gop.fence();
695  }
696 
697 
699  template <typename opT, typename R, std::size_t NDIM>
700  std::vector< Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> >
701  apply(World& world,
702  const std::vector< std::shared_ptr<opT> >& op,
703  const std::vector< Function<R,NDIM> > f) {
704 
705  PROFILE_BLOCK(Vapplyv);
706  MADNESS_ASSERT(f.size()==op.size());
707 
708  std::vector< Function<R,NDIM> >& ncf = *const_cast< std::vector< Function<R,NDIM> >* >(&f);
709 
710  reconstruct(world, f);
711  nonstandard(world, ncf);
712 
713  std::vector< Function<TENSOR_RESULT_TYPE(typename opT::opT,R), NDIM> > result(f.size());
714  for (unsigned int i=0; i<f.size(); ++i) {
715  MADNESS_ASSERT(not op[i]->is_slaterf12);
716  result[i] = apply_only(*op[i], f[i], false);
717  }
718 
719  world.gop.fence();
720 
721  standard(world, ncf, false); // restores promise of logical constness
722  world.gop.fence();
723  reconstruct(world, result);
724 
725  return result;
726  }
727 
728 
730  template <typename T, typename R, std::size_t NDIM>
731  std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> >
732  apply(World& world,
734  const std::vector< Function<R,NDIM> > f) {
735  PROFILE_BLOCK(Vapply);
736 
737  std::vector< Function<R,NDIM> >& ncf = *const_cast< std::vector< Function<R,NDIM> >* >(&f);
738 
739  reconstruct(world, f);
740  nonstandard(world, ncf);
741 
742  std::vector< Function<TENSOR_RESULT_TYPE(T,R), NDIM> > result(f.size());
743  for (unsigned int i=0; i<f.size(); ++i) {
744  result[i] = apply_only(op, f[i], false);
745  }
746 
747  world.gop.fence();
748 
749  standard(world, ncf, false); // restores promise of logical constness
750  reconstruct(world, result);
751 
752  if (op.is_slaterf12) {
753  MADNESS_ASSERT(not op.destructive());
754  for (unsigned int i=0; i<f.size(); ++i) {
755  double trace=f[i].trace();
756  result[i]=(result[i]-trace).scale(-0.5/op.mu());
757  }
758  }
759 
760  return result;
761  }
762 
764  template <typename T, std::size_t NDIM>
765  void normalize(World& world, std::vector< Function<T,NDIM> >& v, bool fence=true) {
766  PROFILE_BLOCK(Vnormalize);
767  std::vector<double> nn = norm2s(world, v);
768  for (unsigned int i=0; i<v.size(); ++i) v[i].scale(1.0/nn[i],false);
769  if (fence) world.gop.fence();
770  }
771 
772  // gives back the size in GB
773  template <typename T, std::size_t NDIM>
774  double get_size(World& world, const std::vector< Function<T,NDIM> >& v){
775 
776  if (v.empty()) return 0.0;
777 
778  const double d=sizeof(T);
779  const double fac=1024*1024*1024;
780 
781  double size=0.0;
782  for(unsigned int i=0;i<v.size();i++){
783  size+=v[i].size();
784  }
785 
786  return size/fac*d;
787 
788  }
789 
790 }
791 #endif // MADNESS_MRA_VMRA_H__INCLUDED
void matrix_inner(DistributedMatrix< T > &A, const std::vector< Function< T, NDIM > > &f, const std::vector< Function< T, NDIM > > &g, bool sym=false)
Definition: chem/distpm.cc:38
WorldGopInterface & gop
Global operations.
Definition: worldfwd.h:462
Definition: shared_ptr_bits.h:38
void reconstruct(World &world, const std::vector< Function< T, NDIM > > &v, bool fence=true)
Reconstruct a vector of functions.
Definition: vmra.h:149
NDIM const Function< R, NDIM > & g
Definition: mra.h:2179
const double R
Definition: dielectric.cc:191
Main include file for MADNESS and defines Function interface.
void gaxpy(World &world, Q alpha, std::vector< Function< T, NDIM > > &a, Q beta, const std::vector< Function< R, NDIM > > &b, bool fence=true)
Generalized A*X+Y for vectors of functions -— a[i] = alpha*a[i] + beta*b[i].
Definition: vmra.h:680
const Function< T, NDIM > & compress(bool fence=true) const
Compresses the function, transforming into wavelet basis. Possible non-blocking comm.
Definition: mra.h:683
Implements derivatives operators with variety of boundary conditions on simulation domain...
Definition: derivative.h:272
void truncate(World &world, std::vector< Function< T, NDIM > > &v, double tol=0.0, bool fence=true)
Truncates a vector of functions.
Definition: vmra.h:194
void nonstandard(World &world, std::vector< Function< T, NDIM > > &v, bool fence=true)
Generates non-standard form of a vector of functions.
Definition: vmra.h:167
NDIM & f
Definition: mra.h:2179
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:1942
Declaration and initialization of tree traversal functions and generic derivative.
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:1734
TENSOR_RESULT_TYPE(T, R) inner(const Function<T
Computes the scalar/inner product between two functions.
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:506
const double beta
Definition: gygi_soltion.cc:63
T inner(const vecfunc< T, NDIM > &a, const vecfunc< T, NDIM > &b)
the non-linear solver requires an inner product
Definition: nemo.h:112
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:312
bool is_slaterf12
Definition: operator.h:134
#define PROFILE_BLOCK(name)
Definition: worldprofile.h:197
const T1 &f1 return GTEST_2_TUPLE_() T(f0, f1)
void scale(World &world, std::vector< Function< T, NDIM > > &v, const std::vector< Q > &factors, bool fence=true)
Scales inplace a vector of functions by distinct values.
Definition: vmra.h:290
void reconstruct(bool fence=true) const
Reconstructs the function, transforming into scaling function basis. Possible non-blocking comm...
Definition: mra.h:737
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:1835
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:1778
FLOAT a(int j, FLOAT z)
Definition: y1.cc:86
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:1634
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:556
void sum(T *buf, size_t nelem)
Inplace global sum while still processing AM & tasks.
Definition: worldgop.h:767
void standard(World &world, std::vector< Function< T, NDIM > > &v, bool fence=true)
Generates standard form of a vector of functions.
Definition: vmra.h:181
tensorT sqrt(const tensorT &s, double tol=1e-8)
Computes matrix square root (not used any more?)
Definition: DFcode/moldft.cc:446
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:1806
A parallel world with full functionality wrapping an MPI communicator.
Definition: worldfwd.h:416
void fence()
Synchronizes all processes in communicator AND globally ensures no pending AM or tasks.
Definition: worldgop.cc:52
double norm2(World &world, const std::vector< Function< T, NDIM > > &v)
Computes the 2-norm of a vector of functions.
Definition: vmra.h:325
A multiresolution adaptive numerical function.
Definition: derivative.h:61
void fence()
Returns after all local tasks have completed.
Definition: worldtask.h:1086
const double m
Definition: gfit.cc:199
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 alread...
Definition: mra.h:1569
WorldTaskQueue & taskq
Task queue.
Definition: worldfwd.h:461
bool & destructive()
Definition: operator.h:161
double get_size(World &world, const std::vector< Function< T, NDIM > > &v)
Definition: vmra.h:774
Tensor< double > op(const Tensor< double > &x)
Definition: kain.cc:508
void norm_tree(bool fence=true) const
Initializes information about the function norm at all length scales.
Definition: mra.h:663
const double & mu() const
Definition: operator.h:165
FunctionFactory implements the named-parameter idiom for Function.
Definition: funcimpl.h:70
Holds machinery to set up Functions/FuncImpls using various Factories and Interfaces.
Definition: chem/atomutil.cc:45
Function< T, NDIM > apply(const Derivative< T, NDIM > &D, const Function< T, NDIM > &f, bool fence=true)
Applies derivative operator to function (for syntactic equivalence to integral operator apply) ...
Definition: derivative.h:613
const double c
Definition: gfit.cc:200
std::vector< Function< TENSOR_RESULT_TYPE(T, R), NDIM > > transform(World &world, const std::vector< Function< T, NDIM > > &v, const DistributedMatrix< R > &c, bool fence=true)
Definition: chem/SCF.cc:86
FLOAT b(int j, FLOAT z)
Definition: y1.cc:79
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:765
std::vector< Function< T, NDIM > > zero_functions(World &world, int n)
Generates a vector of zero functions.
Definition: vmra.h:229
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:1528
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:1879
void compress(World &world, const std::vector< Function< T, NDIM > > &v, bool fence=true)
Compress a vector of functions.
Definition: vmra.h:130