MADNESS  0.10.1
Classes | Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Private Attributes | List of all members
madness::Function< T, NDIM > Class Template Reference

A multiresolution adaptive numerical function. More...

#include <mra.h>

Inheritance diagram for madness::Function< T, NDIM >:
Inheritance graph
[legend]
Collaboration diagram for madness::Function< T, NDIM >:
Collaboration graph
[legend]

Classes

struct  autorefine_square_op
 
struct  SimpleUnaryOpWrapper
 

Public Types

typedef Vector< double, NDIMcoordT
 Type of vector holding coordinates. More...
 
typedef FunctionFactory< T, NDIMfactoryT
 
typedef FunctionImpl< T, NDIMimplT
 
typedef FunctionNode< T, NDIMnodeT
 

Public Member Functions

 Function ()
 Default constructor makes uninitialized function. No communication. More...
 
 Function (const factoryT &factory)
 Constructor from FunctionFactory provides named parameter idiom. Possible non-blocking communication. More...
 
 Function (const Function< T, NDIM > &f)
 Copy constructor is shallow. No communication, works in either basis. More...
 
 ~Function ()
 Destruction of any underlying implementation is deferred to next global fence. More...
 
Function< T, NDIM > & abs (bool fence=true)
 Returns *this for chaining. More...
 
Function< T, NDIM > & abs_square (bool fence=true)
 Returns *this for chaining. More...
 
Function< T, NDIM > & add_scalar (T t, bool fence=true)
 Inplace add scalar. No communication except for optional fence. More...
 
bool autorefine () const
 Retunrs. More...
 
void broaden (const BoundaryConditions< NDIM > &bc=FunctionDefaults< NDIM >::get_bc(), bool fence=true) const
 Inplace broadens support in scaling function basis. More...
 
void change_tensor_type (const TensorArgs &targs, bool fence=true)
 change the tensor type of the coefficients in the FunctionNode More...
 
const Function< T, NDIM > & change_tree_state (const TreeState finalstate, bool fence=true) const
 changes tree state to given state More...
 
g change_tree_state (gstate, false)
 
g change_tree_state (redundant, false)
 
 change_tree_state (state, false)
 
double check_symmetry () const
 check symmetry of a function by computing the 2nd derivative More...
 
Function< T, NDIM > & chop_at_level (const int n, const bool fence=true)
 remove all nodes with level higher than n More...
 
void clear (bool fence=true)
 Clears the function as if constructed uninitialized. Optional fence. More...
 
const Function< T, NDIM > & compress (bool fence=true) const
 Compresses the function, transforming into wavelet basis. Possible non-blocking comm. More...
 
Function< T, NDIMconj (bool fence=true)
 Inplace complex conjugate. No communication except for optional fence. More...
 
Level depthpt (const coordT &xuser) const
 
Function< T, NDIM/2 > dirac_convolution (const bool fence=true) const
 
void distribute (std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > newmap) const
 distribute this function according to newmap More...
 
template<size_t LDIM, size_t KDIM>
void do_hartree_product (const std::vector< std::shared_ptr< FunctionImpl< T, LDIM >>> left, const std::vector< std::shared_ptr< FunctionImpl< T, KDIM >>> right)
 perform the hartree product of f*g, invoked by result More...
 
template<size_t LDIM, size_t KDIM, typename opT >
void do_hartree_product (const std::vector< std::shared_ptr< FunctionImpl< T, LDIM >>> left, const std::vector< std::shared_ptr< FunctionImpl< T, KDIM >>> right, const opT *op)
 perform the hartree product of f*g, invoked by result More...
 
template<typename funcT >
double err (const funcT &func) const
 Returns an estimate of the difference ||this-func|| ... global sum performed. More...
 
template<typename funcT >
double errsq_local (const funcT &func) const
 Returns an estimate of the difference ||this-func||^2 from local data. More...
 
Future< Teval (const coordT &xuser) const
 Evaluates the function at a point in user coordinates. Possible non-blocking comm. More...
 
Tensor< Teval_cube (const Tensor< double > &cell, const std::vector< long > &npt, bool eval_refine=false) const
 Evaluates a cube/slice of points (probably for plotting) ... collective but no fence necessary. More...
 
std::pair< bool, Teval_local_only (const Vector< double, NDIM > &xuser, Level maxlevel) const
 Evaluate function only if point is local returning (true,value); otherwise return (false,0.0) More...
 
Future< Levelevaldepthpt (const coordT &xuser) const
 
Future< long > evalR (const coordT &xuser) const
 Evaluates the function rank at a point in user coordinates. Possible non-blocking comm. More...
 
impl world gop fence ()
 
impl world gop fence ()
 
impl world gop fence ()
 
impl world gop fence ()
 
Function< T, NDIM > & fill_cuspy_tree (const bool fence=true)
 Special refinement on 6D boxes where the electrons come close (meet) More...
 
template<typename opT >
Function< T, NDIM > & fill_cuspy_tree (const opT &op, const bool fence=true)
 
template<typename opT >
Function< T, NDIM > & fill_nuclear_cuspy_tree (const opT &op, const size_t particle, const bool fence=true)
 
Function< T, NDIM > & fill_nuclear_cuspy_tree (const size_t particle, const bool fence=true)
 Special refinement on 6D boxes for the nuclear potentials (regularized with cusp, non-regularized with singularity) More...
 
Function< T, NDIM > & fill_tree (bool fence=true)
 With this being an on-demand function, fill the MRA tree according to different criteria. More...
 
template<typename R >
Function< T, NDIM > & fill_tree (const Function< R, NDIM > &g, bool fence=true)
 With this being an on-demand function, fill the MRA tree according to different criteria. More...
 
template<typename opT >
Function< T, NDIM > & fill_tree (const opT &op, bool fence=true)
 With this being an on-demand function, fill the MRA tree according to different criteria. More...
 
return impl inner_local * g ())
 
template<typename Q , typename R >
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. More...
 
template<typename L >
void gaxpy_ext (const Function< L, NDIM > &left, T(*f)(const coordT &), T alpha, T beta, double tol, bool fence=true) const
 
template<typename L , typename R >
Function< T, NDIM > & gaxpy_oop (T alpha, const Function< L, NDIM > &left, T beta, const Function< R, NDIM > &right, bool fence)
 This is replaced with alpha*left + beta*right ... private. More...
 
const std::shared_ptr< FunctionImpl< T, NDIM > > & get_impl () const
 Returns a shared-pointer to the implementation. More...
 
const std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > > & get_pmap () const
 Returns a shared pointer to the process map. More...
 
 if (g.is_on_demand()) return this -> inner_on_demand(g)
 
 if (NDIM<=3)
 
 if (not g.is_initialized()) return 0.0
 
 if (not this->is_initialized()) return 0.0
 
 if (this->get_impl()==g.get_impl())
 
 if (this->is_compressed() and g.is_compressed())
 
 if (this->is_on_demand()) return g.inner_on_demand(*this)
 
 if (VERIFY_TREE) g.verify_tree()
 
 if (VERIFY_TREE) g.verify_tree()
 
 if (VERIFY_TREE) verify_tree()
 
 if (VERIFY_TREE) verify_tree()
 
 if (VERIFY_TREE) verify_tree()
 
bool impl_initialized () const
 
T inner_adaptive (const std::shared_ptr< FunctionFunctorInterface< T, NDIM > > f, const bool leaf_refine=true) const
 
T inner_ext (const std::shared_ptr< FunctionFunctorInterface< T, NDIM > > f, const bool leaf_refine=true, const bool keep_redundant=false) const
 
T inner_ext_local (const std::shared_ptr< FunctionFunctorInterface< T, NDIM > > f, const bool leaf_refine=true, const bool keep_redundant=false) const
 
bool is_compressed () const
 Returns true if compressed, false otherwise. No communication. More...
 
bool is_initialized () const
 Returns true if the function is initialized. More...
 
bool is_nonstandard () const
 Returns true if nonstandard-compressed, false otherwise. No communication. More...
 
bool is_on_demand () const
 
bool is_reconstructed () const
 Returns true if reconstructed, false otherwise. No communication. More...
 
int k () const
 Returns the number of multiwavelets (k). No communication. More...
 
template<typename Archive >
void load (World &world, Archive &ar)
 Replaces this function with one loaded from an archive using the default processor map. More...
 
 MADNESS_ASSERT (func)
 
 MADNESS_ASSERT (g.is_compressed())
 
 MADNESS_ASSERT (is_compressed())
 
void make_nonstandard (bool keepleaves, bool fence=true) const
 Compresses the function retaining scaling function coeffs. Possible non-blocking comm. More...
 
void make_redundant (bool fence=true)
 Converts the function to redundant form, i.e. sum coefficients on all levels. More...
 
func make_redundant (true)
 
Function< T, NDIM > & map_and_mirror (const Function< T, NDIM > &f, const std::vector< long > &map, const std::vector< long > &mirror, bool fence)
 This is replaced with mirror(map(f)) ... private. More...
 
Function< T, NDIM > & mapdim (const Function< T, NDIM > &f, const std::vector< long > &map, bool fence)
 This is replaced with mapdim(f) ... private. More...
 
std::size_t max_depth () const
 Returns the maximum depth of the function tree ... collective global sum. More...
 
std::size_t max_local_depth () const
 Returns the maximum local depth of the function tree ... no communications. More...
 
std::size_t max_nodes () const
 Returns the max number of nodes on a processor. More...
 
std::size_t min_nodes () const
 Returns the min number of nodes on a processor. More...
 
Function< T, NDIM > & mirror (const Function< T, NDIM > &f, const std::vector< long > &mirrormap, bool fence)
 This is replaced with mirror(f) ... private. More...
 
template<typename L , typename R >
void mul_on_demand (const Function< L, NDIM > &f, const Function< R, NDIM > &g, bool fence=true)
 Same as operator* but with optional fence and no automatic reconstruction. More...
 
template<typename opT >
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 More...
 
template<typename opT >
Function< T, NDIM > & multiop_values (const opT &op, const std::vector< Function< T, NDIM > > &vf)
 This is replaced with op(vector of functions) ... private. More...
 
double norm2 () const
 Returns the 2-norm of the function ... global sum ... works in either basis. More...
 
double norm2sq_local () const
 Returns the square of the norm of the local function ... no communication. More...
 
void norm_tree (bool fence=true) const
 Initializes information about the function norm at all length scales. More...
 
T operator() (const coordT &xuser) const
 Evaluates the function at a point in user coordinates. Collective operation. More...
 
T operator() (double x, double y=0, double z=0, double xx=0, double yy=0, double zz=0) const
 Evaluates the function at a point in user coordinates. Collective operation. More...
 
template<typename Q >
IsSupported< TensorTypeData< Q >, Function< T, NDIM > >::typeoperator*= (const Q q)
 Inplace scaling by a constant. More...
 
template<typename Q >
Function< T, NDIM > & operator+= (const Function< Q, NDIM > &other)
 Inplace addition of functions in the wavelet basis. More...
 
template<typename Q >
Function< T, NDIM > & operator-= (const Function< Q, NDIM > &other)
 Inplace subtraction of functions in the wavelet basis. More...
 
Function< T, NDIM > & operator= (const Function< T, NDIM > &f)
 Assignment is shallow. No communication, works in either basis. More...
 
void print_info () const
 Print a summary of the load balancing info. More...
 
void print_size (const std::string name) const
 print some info about this More...
 
void print_tree (std::ostream &os=std::cout) const
 Process 0 prints a summary of all nodes in the tree (collective) More...
 
void print_tree_graphviz (std::ostream &os=std::cout) const
 Process 0 prints a graphviz-formatted output of all nodes in the tree (collective) More...
 
void print_tree_json (std::ostream &os=std::cout) const
 
template<typename R , size_t LDIM>
Function< TENSOR_RESULT_TYPE(T, R), NDIM-LDIMproject_out (const Function< R, LDIM > &g, const int dim) const
 project this on the low-dim function g: h(x) = <f(x,y) | g(y)> More...
 
this reconstruct ()
 
const Function< T, NDIM > & reconstruct (bool fence=true) const
 Reconstructs the function, transforming into scaling function basis. Possible non-blocking comm. More...
 
Function< T, NDIM > & reduce_rank (const double thresh=0.0, const bool fence=true)
 reduce the rank of the coefficient tensors More...
 
const Function< T, NDIM > & refine (bool fence=true) const
 Inplace autorefines the function using same test as for squaring. More...
 
template<typename opT >
void refine_general (const opT &op, bool fence=true) const
 Inplace autorefines the function. Optional fence. Possible non-blocking comm. More...
 
void replicate (bool fence=true) const
 replicate this function, generating a unique pmap More...
 
func replicate_low_dim_functions (true)
 
template<typename Q >
Function< T, NDIM > & scale (const Q q, bool fence=true)
 Inplace, scale the function by a constant. No communication except for optional fence. More...
 
void set_autorefine (bool value, bool fence=true)
 Sets the value of the autorefine flag. Optional global fence. More...
 
void set_functor (const std::shared_ptr< FunctionFunctorInterface< T, NDIM > > functor)
 Replace the current functor with the provided new one. More...
 
template<typename R >
void set_impl (const Function< R, NDIM > &f, bool zero=true)
 Replace current FunctionImpl with a new one using the same parameters & map as f. More...
 
void set_impl (const std::shared_ptr< FunctionImpl< T, NDIM > > &impl)
 Replace current FunctionImpl with provided new one. More...
 
void set_thresh (double value, bool fence=true)
 Sets the value of the truncation threshold. Optional global fence. More...
 
std::size_t size () const
 Returns the number of coefficients in the function ... collective global sum. More...
 
Function< T, NDIM > & square (bool fence=true)
 Inplace squaring of function ... global comm only if not reconstructed. More...
 
void standard (bool fence=true)
 Converts the function standard compressed form. Possible non-blocking comm. More...
 
template<typename Archive >
void store (Archive &ar) const
 Stores the function to an archive. More...
 
impl world gop sum (local)
 
impl world gop sum (local)
 
void sum_down (bool fence=true) const
 Sums scaling coeffs down tree restoring state with coeffs only at leaves. Optional fence. Possible non-blocking comm. More...
 
template<typename R >
 TENSOR_RESULT_TYPE (T, R) inner(const Function< R
 Returns the inner product. More...
 
template<typename R >
 TENSOR_RESULT_TYPE (T, R) inner_local(const Function< R
 Returns local part of inner product ... throws if both not compressed. More...
 
template<typename R >
 TENSOR_RESULT_TYPE (T, R) inner_on_demand(const Function< R
 Returns the inner product for one on-demand function. More...
 
 TENSOR_RESULT_TYPE (T, R) local
 
 TENSOR_RESULT_TYPE (T, R) local
 
double thresh () const
 Returns value of truncation threshold. No communication. More...
 
T trace () const
 Returns global value of int(f(x),x) ... global comm required. More...
 
T trace_local () const
 Returns local contribution to int(f(x),x) ... no communication. More...
 
std::size_t tree_size () const
 Returns the number of nodes in the function tree ... collective global sum. More...
 
Function< T, NDIM > & truncate (double tol=0.0, bool fence=true)
 Truncate the function with optional fence. Compresses with fence if not compressed. More...
 
template<typename Q , typename opT >
Function< typename opT::resultT, NDIM > & unary_op_coeffs (const Function< Q, NDIM > &func, const opT &op, bool fence)
 This is replaced with left*right ... private. More...
 
template<typename opT >
void unaryop (const opT &op, bool fence=true)
 Inplace unary operation on function values. More...
 
void unaryop (T(*f)(T))
 Inplace unary operation on function values. More...
 
template<typename opT >
void unaryop_coeff (const opT &op, bool fence=true)
 Unary operation applied inplace to the coefficients. More...
 
template<typename opT >
void unaryop_node (const opT &op, bool fence=true)
 Unary operation applied inplace to the nodes. More...
 
void verify () const
 Asserts that the function is initialized. More...
 
void verify_tree () const
 Verifies the tree data structure ... global sync implied. More...
 
template<typename L , typename R >
void vmulXX (const Function< L, NDIM > &left, const std::vector< Function< R, NDIM > > &right, std::vector< Function< T, NDIM > > &result, double tol, bool fence)
 Multiplication of function * vector of functions using recursive algorithm of mulxx. More...
 
template<typename R , typename Q >
void vtransform (const std::vector< Function< R, NDIM > > &v, const Tensor< Q > &c, std::vector< Function< T, NDIM > > &vresult, double tol, bool fence=true)
 sparse transformation of a vector of functions ... private More...
 
Worldworld () const
 Returns the world. More...
 

Static Public Member Functions

static void doconj (const Key< NDIM >, Tensor< T > &t)
 
template<typename Q , std::size_t D>
static std::vector< std::shared_ptr< FunctionImpl< Q, D > > > vimpl (const std::vector< Function< Q, D > > &v)
 Returns vector of FunctionImpl pointers corresponding to vector of functions. More...
 

Public Attributes

NDIM &g const
 
 else
 
auto func =dynamic_cast<CompositeFunctorInterface<T,NDIM,LDIM>* >(g.get_impl()->get_functor().get())
 
TreeState gstate =g.get_impl()->get_tree_state()
 
constexpr std::size_t LDIM =std::max(NDIM/2,std::size_t(1))
 
return local
 
TreeState state =this->get_impl()->get_tree_state()
 

Private Attributes

std::shared_ptr< FunctionImpl< T, NDIM > > impl
 

Detailed Description

template<typename T, std::size_t NDIM>
class madness::Function< T, NDIM >

A multiresolution adaptive numerical function.

Member Typedef Documentation

◆ coordT

template<typename T , std::size_t NDIM>
typedef Vector<double,NDIM> madness::Function< T, NDIM >::coordT

Type of vector holding coordinates.

◆ factoryT

template<typename T , std::size_t NDIM>
typedef FunctionFactory<T,NDIM> madness::Function< T, NDIM >::factoryT

◆ implT

template<typename T , std::size_t NDIM>
typedef FunctionImpl<T,NDIM> madness::Function< T, NDIM >::implT

◆ nodeT

template<typename T , std::size_t NDIM>
typedef FunctionNode<T,NDIM> madness::Function< T, NDIM >::nodeT

Constructor & Destructor Documentation

◆ Function() [1/3]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::Function ( )
inline

Default constructor makes uninitialized function. No communication.

An uninitialized function can only be assigned to. Any other operation will throw.

◆ Function() [2/3]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::Function ( const factoryT factory)
inline

Constructor from FunctionFactory provides named parameter idiom. Possible non-blocking communication.

References PROFILE_MEMBER_FUNC.

◆ Function() [3/3]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::Function ( const Function< T, NDIM > &  f)
inline

Copy constructor is shallow. No communication, works in either basis.

◆ ~Function()

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::~Function ( )
inline

Destruction of any underlying implementation is deferred to next global fence.

Member Function Documentation

◆ abs()

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::abs ( bool  fence = true)
inline

◆ abs_square()

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::abs_square ( bool  fence = true)
inline

◆ add_scalar()

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::add_scalar ( T  t,
bool  fence = true 
)
inline

◆ autorefine()

template<typename T , std::size_t NDIM>
bool madness::Function< T, NDIM >::autorefine ( ) const
inline

Retunrs.

Returns value of autorefine flag. No communication.

References madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

Referenced by madness::smooth< T, NDIM >::make_sigma().

◆ broaden()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::broaden ( const BoundaryConditions< NDIM > &  bc = FunctionDefaults<NDIM>::get_bc(),
bool  fence = true 
) const
inline

◆ change_tensor_type()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::change_tensor_type ( const TensorArgs targs,
bool  fence = true 
)
inline

change the tensor type of the coefficients in the FunctionNode

Parameters
[in]targstarget tensor arguments (threshold and full/low rank)

References madness::Function< T, NDIM >::fence(), and madness::Function< T, NDIM >::impl.

◆ change_tree_state() [1/4]

template<typename T , std::size_t NDIM>
const Function<T,NDIM>& madness::Function< T, NDIM >::change_tree_state ( const TreeState  finalstate,
bool  fence = true 
) const
inline

◆ change_tree_state() [2/4]

template<typename T , std::size_t NDIM>
g madness::Function< T, NDIM >::change_tree_state ( gstate  ,
false   
)

◆ change_tree_state() [3/4]

template<typename T , std::size_t NDIM>
g madness::Function< T, NDIM >::change_tree_state ( redundant  ,
false   
)

◆ change_tree_state() [4/4]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::change_tree_state ( state  ,
false   
)

◆ check_symmetry()

template<typename T , std::size_t NDIM>
double madness::Function< T, NDIM >::check_symmetry ( ) const
inline

◆ chop_at_level()

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::chop_at_level ( const int  n,
const bool  fence = true 
)
inline

◆ clear()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::clear ( bool  fence = true)
inline

◆ compress()

template<typename T , std::size_t NDIM>
const Function<T,NDIM>& madness::Function< T, NDIM >::compress ( bool  fence = true) const
inline

Compresses the function, transforming into wavelet basis. Possible non-blocking comm.

By default fence=true meaning that this operation completes before returning, otherwise if fence=false it returns without fencing and the user must invoke world.gop.fence() to assure global completion before using the function for other purposes.

Noop if already compressed or if not initialized.

Since reconstruction/compression do not discard information we define them as const ... "logical constness" not "bitwise constness".

References madness::Function< T, NDIM >::change_tree_state(), madness::compressed, and madness::Function< T, NDIM >::fence().

Referenced by RealFuncLinearOp::action(), madness::NuclearCorrelationFactor::apply_U(), madness::Solver< T, NDIM >::compute_rho(), energy(), main(), madness::SCF::make_density(), MiniDFT::make_density(), madness::GTHPseudopotential< Q >::make_pseudo_potential(), Calculation::make_reference(), makeao(), moments(), madness::operator+(), madness::Function< T, NDIM >::operator+=(), madness::operator-(), madness::Function< T, NDIM >::operator-=(), propagate(), q_c(), realfunc2(), test_add(), and madness::NonlinearSolverND< NDIM >::update().

◆ conj()

template<typename T , std::size_t NDIM>
Function<T,NDIM> madness::Function< T, NDIM >::conj ( bool  fence = true)
inline

Inplace complex conjugate. No communication except for optional fence.

Returns this for chaining. Works in either basis.

References madness::Function< T, NDIM >::fence(), PROFILE_MEMBER_FUNC, and madness::Function< T, NDIM >::unaryop_coeff().

Referenced by madness::conj().

◆ depthpt()

template<typename T , std::size_t NDIM>
Level madness::Function< T, NDIM >::depthpt ( const coordT xuser) const
inline

Throws if function is not initialized.

This function mimics operator() by going through the tree looking for the depth of the tree at the point. It blocks until the result is available and then broadcasts the result to everyone. Therefore, if you are evaluating many points in parallel it is vastly less efficient than calling evaldepthpt directly, saving the futures, and then forcing all of the results.

References madness::Function< T, NDIM >::evaldepthpt(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_reconstructed(), PROFILE_MEMBER_FUNC, madness::Function< T, NDIM >::reconstruct(), and madness::Function< T, NDIM >::verify().

◆ dirac_convolution()

template<typename T , std::size_t NDIM>
Function<T,NDIM/2> madness::Function< T, NDIM >::dirac_convolution ( const bool  fence = true) const
inline

◆ distribute()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::distribute ( std::shared_ptr< WorldDCPmapInterface< Key< NDIM > > >  newmap) const
inline

distribute this function according to newmap

References madness::Function< T, NDIM >::impl, and madness::Function< T, NDIM >::verify().

Referenced by test_replicate().

◆ do_hartree_product() [1/2]

template<typename T , std::size_t NDIM>
template<size_t LDIM, size_t KDIM>
void madness::Function< T, NDIM >::do_hartree_product ( const std::vector< std::shared_ptr< FunctionImpl< T, LDIM >>>  left,
const std::vector< std::shared_ptr< FunctionImpl< T, KDIM >>>  right 
)
inline

perform the hartree product of f*g, invoked by result

References madness::Function< T, NDIM >::impl, and madness::Function< T, NDIM >::k().

◆ do_hartree_product() [2/2]

template<typename T , std::size_t NDIM>
template<size_t LDIM, size_t KDIM, typename opT >
void madness::Function< T, NDIM >::do_hartree_product ( const std::vector< std::shared_ptr< FunctionImpl< T, LDIM >>>  left,
const std::vector< std::shared_ptr< FunctionImpl< T, KDIM >>>  right,
const opT *  op 
)
inline

perform the hartree product of f*g, invoked by result

References madness::Function< T, NDIM >::impl, and op().

Referenced by madness::hartree_product().

◆ doconj()

template<typename T , std::size_t NDIM>
static void madness::Function< T, NDIM >::doconj ( const Key< NDIM ,
Tensor< T > &  t 
)
inlinestatic

◆ err()

template<typename T , std::size_t NDIM>
template<typename funcT >
double madness::Function< T, NDIM >::err ( const funcT &  func) const
inline

Returns an estimate of the difference ||this-func|| ... global sum performed.

If the function is compressed, it is reconstructed first. For efficient use especially with many functions, reconstruct them all first, and use errsq_local instead so you can perform a global sum on all at the same time.

Parameters
[in]funcTemplated interface to the a user specified function

References madness::Function< T, NDIM >::func, madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_reconstructed(), madness::Function< T, NDIM >::local, PROFILE_MEMBER_FUNC, madness::Function< T, NDIM >::reconstruct(), madness::Function< T, NDIM >::verify(), VERIFY_TREE, and madness::Function< T, NDIM >::verify_tree().

Referenced by main(), test_add(), test_bsh(), test_coulomb(), test_diff(), test_math(), test_op(), test_opdir(), and testgradG().

◆ errsq_local()

template<typename T , std::size_t NDIM>
template<typename funcT >
double madness::Function< T, NDIM >::errsq_local ( const funcT &  func) const
inline

Returns an estimate of the difference ||this-func||^2 from local data.

No communication is performed. If the function is not reconstructed, it throws an exception. To get the global value either do a global sum of the local values or call errsq

Parameters
[in]funcTemplated interface to the a user specified function

References madness::Function< T, NDIM >::func, madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_reconstructed(), MADNESS_EXCEPTION, PROFILE_MEMBER_FUNC, and madness::Function< T, NDIM >::verify().

◆ eval()

template<typename T , std::size_t NDIM>
Future<T> madness::Function< T, NDIM >::eval ( const coordT xuser) const
inline

Evaluates the function at a point in user coordinates. Possible non-blocking comm.

Only the invoking process will receive the result via the future though other processes may be involved in the evaluation.

Throws if function is not initialized.

References d(), e(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_reconstructed(), madness::Function< T, NDIM >::MADNESS_ASSERT(), MADNESS_EXCEPTION, NDIM, PROFILE_MEMBER_FUNC, madness::Future< T >::remote_ref(), madness::user_to_sim(), and madness::Function< T, NDIM >::verify().

Referenced by madness::Function< T, NDIM >::operator()(), and test_apply_push_1d().

◆ eval_cube()

template<typename T , std::size_t NDIM>
Tensor<T> madness::Function< T, NDIM >::eval_cube ( const Tensor< double > &  cell,
const std::vector< long > &  npt,
bool  eval_refine = false 
) const
inline

Evaluates a cube/slice of points (probably for plotting) ... collective but no fence necessary.

All processes recieve the entire result (which is a rather severe limit on the size of the cube that is possible). Set eval_refine=true to return the refinment levels of the given function.

Parameters
[in]cellA Tensor describe the cube where the function to be evaluated in
[in]nptHow many points to evaluate in each dimension
[in]eval_refineWether to return the refinment levels of the given function

References d(), delta, madness::BaseTensor::dim(), e(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::MADNESS_ASSERT(), NDIM, PROFILE_MEMBER_FUNC, madness::Function< T, NDIM >::reconstruct(), madness::user_to_sim(), and madness::Function< T, NDIM >::verify().

◆ eval_local_only()

template<typename T , std::size_t NDIM>
std::pair<bool,T> madness::Function< T, NDIM >::eval_local_only ( const Vector< double, NDIM > &  xuser,
Level  maxlevel 
) const
inline

Evaluate function only if point is local returning (true,value); otherwise return (false,0.0)

maxlevel is the maximum depth to search down to — the max local depth can be computed with max_local_depth();

References d(), e(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_reconstructed(), madness::Function< T, NDIM >::MADNESS_ASSERT(), MADNESS_EXCEPTION, NDIM, madness::user_to_sim(), and madness::Function< T, NDIM >::verify().

◆ evaldepthpt()

template<typename T , std::size_t NDIM>
Future<Level> madness::Function< T, NDIM >::evaldepthpt ( const coordT xuser) const
inline

Only the invoking process will receive the result via the future though other processes may be involved in the evaluation.

Throws if function is not initialized.

This function is a minimally-modified version of eval()

References d(), e(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_reconstructed(), madness::Function< T, NDIM >::MADNESS_ASSERT(), MADNESS_EXCEPTION, NDIM, PROFILE_MEMBER_FUNC, madness::Future< T >::remote_ref(), madness::user_to_sim(), and madness::Function< T, NDIM >::verify().

Referenced by madness::Function< T, NDIM >::depthpt().

◆ evalR()

template<typename T , std::size_t NDIM>
Future<long> madness::Function< T, NDIM >::evalR ( const coordT xuser) const
inline

Evaluates the function rank at a point in user coordinates. Possible non-blocking comm.

Only the invoking process will receive the result via the future though other processes may be involved in the evaluation.

Throws if function is not initialized.

References d(), e(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_reconstructed(), madness::Function< T, NDIM >::MADNESS_ASSERT(), MADNESS_EXCEPTION, NDIM, PROFILE_MEMBER_FUNC, madness::Future< T >::remote_ref(), madness::user_to_sim(), and madness::Function< T, NDIM >::verify().

◆ fence() [1/4]

template<typename T , std::size_t NDIM>
impl world gop madness::Function< T, NDIM >::fence ( )

Referenced by madness::Function< T, NDIM >::abs(), madness::Function< T, NDIM >::abs_square(), madness::Function< T, NDIM >::add_scalar(), madness::Function< T, NDIM >::broaden(), madness::Function< T, NDIM >::change_tensor_type(), madness::Function< T, NDIM >::change_tree_state(), madness::Function< T, NDIM >::clear(), madness::Function< T, NDIM >::compress(), madness::Function< T, NDIM >::conj(), madness::Function< T, NDIM >::dirac_convolution(), madness::Function< T, NDIM >::fill_cuspy_tree(), madness::Function< T, NDIM >::fill_nuclear_cuspy_tree(), madness::Function< T, NDIM >::fill_tree(), madness::Function< T, NDIM >::gaxpy(), madness::Function< T, NDIM >::gaxpy_ext(), madness::Function< T, NDIM >::gaxpy_oop(), madness::Function< T, NDIM >::make_nonstandard(), madness::Function< T, NDIM >::make_redundant(), madness::Function< T, NDIM >::map_and_mirror(), madness::Function< T, NDIM >::mapdim(), madness::Function< T, NDIM >::mirror(), madness::Function< T, NDIM >::mul_on_demand(), madness::Function< T, NDIM >::multi_to_multi_op_values(), madness::Function< T, NDIM >::norm_tree(), madness::Function< T, NDIM >::reconstruct(), madness::Function< T, NDIM >::reduce_rank(), madness::Function< T, NDIM >::refine(), madness::Function< T, NDIM >::refine_general(), madness::Function< T, NDIM >::replicate(), madness::Function< T, NDIM >::scale(), madness::Function< T, NDIM >::set_autorefine(), madness::Function< T, NDIM >::set_thresh(), madness::Function< T, NDIM >::square(), madness::Function< T, NDIM >::standard(), madness::Function< T, NDIM >::sum_down(), madness::Function< T, NDIM >::truncate(), madness::Function< T, NDIM >::unary_op_coeffs(), madness::Function< T, NDIM >::unaryop(), madness::Function< T, NDIM >::unaryop_coeff(), madness::Function< T, NDIM >::unaryop_node(), madness::Function< T, NDIM >::vmulXX(), and madness::Function< T, NDIM >::vtransform().

◆ fence() [2/4]

template<typename T , std::size_t NDIM>
impl world gop madness::Function< T, NDIM >::fence ( )

◆ fence() [3/4]

template<typename T , std::size_t NDIM>
impl world gop madness::Function< T, NDIM >::fence ( )

◆ fence() [4/4]

template<typename T , std::size_t NDIM>
impl world gop madness::Function< T, NDIM >::fence ( )

◆ fill_cuspy_tree() [1/2]

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::fill_cuspy_tree ( const bool  fence = true)
inline

◆ fill_cuspy_tree() [2/2]

template<typename T , std::size_t NDIM>
template<typename opT >
Function<T,NDIM>& madness::Function< T, NDIM >::fill_cuspy_tree ( const opT &  op,
const bool  fence = true 
)
inline

Special refinement on 6D boxes where the electrons come close (meet)

Parameters
[in]opthe convolution operator for screening

References madness::Function< T, NDIM >::fence(), madness::Function< T, NDIM >::get_impl(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_on_demand(), madness::Function< T, NDIM >::MADNESS_ASSERT(), and op().

Referenced by madness::MP2::make_Uphi0(), and test_vector_composite().

◆ fill_nuclear_cuspy_tree() [1/2]

template<typename T , std::size_t NDIM>
template<typename opT >
Function<T,NDIM>& madness::Function< T, NDIM >::fill_nuclear_cuspy_tree ( const opT &  op,
const size_t  particle,
const bool  fence = true 
)
inline

Special refinement on 6D boxes for the nuclear potentials (regularized with cusp, non-regularized with singularity)

Parameters
[in]opthe convolution operator for screening

References madness::Function< T, NDIM >::fence(), madness::Function< T, NDIM >::get_impl(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_on_demand(), madness::Function< T, NDIM >::MADNESS_ASSERT(), and op().

◆ fill_nuclear_cuspy_tree() [2/2]

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::fill_nuclear_cuspy_tree ( const size_t  particle,
const bool  fence = true 
)
inline

Special refinement on 6D boxes for the nuclear potentials (regularized with cusp, non-regularized with singularity)

References madness::Function< T, NDIM >::fence(), madness::Function< T, NDIM >::get_impl(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_on_demand(), and madness::Function< T, NDIM >::MADNESS_ASSERT().

◆ fill_tree() [1/3]

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::fill_tree ( bool  fence = true)
inline

◆ fill_tree() [2/3]

template<typename T , std::size_t NDIM>
template<typename R >
Function<T,NDIM>& madness::Function< T, NDIM >::fill_tree ( const Function< R, NDIM > &  g,
bool  fence = true 
)
inline

◆ fill_tree() [3/3]

template<typename T , std::size_t NDIM>
template<typename opT >
Function<T,NDIM>& madness::Function< T, NDIM >::fill_tree ( const opT &  op,
bool  fence = true 
)
inline

With this being an on-demand function, fill the MRA tree according to different criteria.

Parameters
[in]opthe convolution operator for screening

References madness::Function< T, NDIM >::fence(), madness::Function< T, NDIM >::get_impl(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_on_demand(), madness::Function< T, NDIM >::MADNESS_ASSERT(), and op().

◆ g()

template<typename T , std::size_t NDIM>
return impl inner_local* madness::Function< T, NDIM >::g ( )

◆ gaxpy()

template<typename T , std::size_t NDIM>
template<typename Q , typename R >
Function<T,NDIM>& madness::Function< T, NDIM >::gaxpy ( const T alpha,
const Function< Q, NDIM > &  other,
const R beta,
bool  fence = true 
)
inline

◆ gaxpy_ext()

template<typename T , std::size_t NDIM>
template<typename L >
void madness::Function< T, NDIM >::gaxpy_ext ( const Function< L, NDIM > &  left,
T(*)(const coordT &)  f,
T  alpha,
T  beta,
double  tol,
bool  fence = true 
) const
inline

Return the local part of gaxpy with external function, this*alpha + f*beta ... no communication.

Parameters
[in]alphaprefactor for this Function
[in]fPointer to function of type T that take coordT arguments. This is the externally provided function
[in]betaprefactor for f

References alpha, beta, madness::f, madness::Function< T, NDIM >::fence(), madness::Function< T, NDIM >::get_impl(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_reconstructed(), PROFILE_MEMBER_FUNC, and madness::Function< T, NDIM >::reconstruct().

Referenced by main().

◆ gaxpy_oop()

template<typename T , std::size_t NDIM>
template<typename L , typename R >
Function<T,NDIM>& madness::Function< T, NDIM >::gaxpy_oop ( T  alpha,
const Function< L, NDIM > &  left,
T  beta,
const Function< R, NDIM > &  right,
bool  fence 
)
inline

◆ get_impl()

template<typename T , std::size_t NDIM>
const std::shared_ptr< FunctionImpl<T,NDIM> >& madness::Function< T, NDIM >::get_impl ( ) const
inline

◆ get_pmap()

template<typename T , std::size_t NDIM>
const std::shared_ptr< WorldDCPmapInterface< Key<NDIM> > >& madness::Function< T, NDIM >::get_pmap ( ) const
inline

◆ if() [1/12]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::if ( g.  is_on_demand()) -> inner_on_demand(g)

◆ if() [2/12]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::if ( NDIM<=  3)
inline

◆ if() [3/12]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::if ( not g.  is_initialized())

◆ if() [4/12]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::if ( not this->  is_initialized())

◆ if() [5/12]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::if ( this->  get_impl() = =g.get_impl())
inline

◆ if() [6/12]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::if ( this->  is_compressed) and g.is_compressed()
inline

◆ if() [7/12]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::if ( this->  is_on_demand())

◆ if() [8/12]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::if ( VERIFY_TREE  )

◆ if() [9/12]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::if ( VERIFY_TREE  )

◆ if() [10/12]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::if ( VERIFY_TREE  )

◆ if() [11/12]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::if ( VERIFY_TREE  )

◆ if() [12/12]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::if ( VERIFY_TREE  )

◆ impl_initialized()

template<typename T , std::size_t NDIM>
bool madness::Function< T, NDIM >::impl_initialized ( ) const
inline

◆ inner_adaptive()

template<typename T , std::size_t NDIM>
T madness::Function< T, NDIM >::inner_adaptive ( const std::shared_ptr< FunctionFunctorInterface< T, NDIM > >  f,
const bool  leaf_refine = true 
) const
inline

Return the inner product with external function ... requires communication. If you are going to be doing a bunch of inner_ext calls, set keep_redundant to true and then manually undo_redundant when you are finished.

Parameters
[in]fReference to FunctionFunctorInterface. This is the externally provided function
[in]leaf_refineboolean switch to turn on/off refinement past leaf nodes
Returns
Returns the inner product

References madness::f, madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::local, PROFILE_MEMBER_FUNC, madness::Function< T, NDIM >::reconstruct(), and T().

◆ inner_ext()

template<typename T , std::size_t NDIM>
T madness::Function< T, NDIM >::inner_ext ( const std::shared_ptr< FunctionFunctorInterface< T, NDIM > >  f,
const bool  leaf_refine = true,
const bool  keep_redundant = false 
) const
inline

Return the inner product with external function ... requires communication. If you are going to be doing a bunch of inner_ext calls, set keep_redundant to true and then manually undo_redundant when you are finished.

Parameters
[in]fReference to FunctionFunctorInterface. This is the externally provided function
[in]leaf_refineboolean switch to turn on/off refinement past leaf nodes
[in]keep_redundantboolean switch to turn on/off undo_redundant
Returns
Returns the inner product

References madness::Function< T, NDIM >::change_tree_state(), madness::f, madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::local, PROFILE_MEMBER_FUNC, madness::reconstructed, madness::redundant, and T().

Referenced by main().

◆ inner_ext_local()

template<typename T , std::size_t NDIM>
T madness::Function< T, NDIM >::inner_ext_local ( const std::shared_ptr< FunctionFunctorInterface< T, NDIM > >  f,
const bool  leaf_refine = true,
const bool  keep_redundant = false 
) const
inline

Return the local part of inner product with external function ... no communication. If you are going to be doing a bunch of inner_ext calls, set keep_redundant to true and then manually undo_redundant when you are finished.

Parameters
[in]fPointer to function of type T that take coordT arguments. This is the externally provided function
[in]leaf_refineboolean switch to turn on/off refinement past leaf nodes
[in]keep_redundantboolean switch to turn on/off undo_redundant
Returns
Returns local part of the inner product, i.e. over the domain of all function nodes on this compute node.

References madness::Function< T, NDIM >::change_tree_state(), madness::f, madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::local, PROFILE_MEMBER_FUNC, madness::reconstructed, madness::redundant, and T().

◆ is_compressed()

template<typename T , std::size_t NDIM>
bool madness::Function< T, NDIM >::is_compressed ( ) const
inline

Returns true if compressed, false otherwise. No communication.

If the function is not initialized, returns false.

References madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

Referenced by madness::Function< T, NDIM >::gaxpy(), madness::Function< T, NDIM >::gaxpy_oop(), madness::operator+(), madness::operator-(), and test_add().

◆ is_initialized()

template<typename T , std::size_t NDIM>
bool madness::Function< T, NDIM >::is_initialized ( ) const
inline

◆ is_nonstandard()

template<typename T , std::size_t NDIM>
bool madness::Function< T, NDIM >::is_nonstandard ( ) const
inline

Returns true if nonstandard-compressed, false otherwise. No communication.

If the function is not initialized, returns false.

References madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

◆ is_on_demand()

template<typename T , std::size_t NDIM>
bool madness::Function< T, NDIM >::is_on_demand ( ) const
inline

◆ is_reconstructed()

template<typename T , std::size_t NDIM>
bool madness::Function< T, NDIM >::is_reconstructed ( ) const
inline

◆ k()

template<typename T , std::size_t NDIM>
int madness::Function< T, NDIM >::k ( ) const
inline

◆ load()

template<typename T , std::size_t NDIM>
template<typename Archive >
void madness::Function< T, NDIM >::load ( World world,
Archive &  ar 
)
inline

Replaces this function with one loaded from an archive using the default processor map.

Archive can be sequential or parallel.

The & operator for serializing will only work with parallel archives.

References madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::k(), madness::Function< T, NDIM >::MADNESS_ASSERT(), NDIM, PROFILE_MEMBER_FUNC, and madness::Function< T, NDIM >::world().

◆ MADNESS_ASSERT() [1/3]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::MADNESS_ASSERT ( func  )

◆ MADNESS_ASSERT() [2/3]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::MADNESS_ASSERT ( g.  is_compressed())

◆ MADNESS_ASSERT() [3/3]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::MADNESS_ASSERT ( is_compressed()  )

◆ make_nonstandard()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::make_nonstandard ( bool  keepleaves,
bool  fence = true 
) const
inline

Compresses the function retaining scaling function coeffs. Possible non-blocking comm.

By default fence=true meaning that this operation completes before returning, otherwise if fence=false it returns without fencing and the user must invoke world.gop.fence() to assure global completion before using the function for other purposes.

Noop if already compressed or if not initialized.

References madness::Function< T, NDIM >::change_tree_state(), madness::Function< T, NDIM >::fence(), madness::nonstandard, and madness::nonstandard_with_leaves.

Referenced by madness::apply(), and madness::hartree_product().

◆ make_redundant() [1/2]

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::make_redundant ( bool  fence = true)
inline

Converts the function to redundant form, i.e. sum coefficients on all levels.

By default fence=true meaning that this operation completes before returning, otherwise if fence=false it returns without fencing and the user must invoke world.gop.fence() to assure global completion before using the function for other purposes.

Must be already compressed.

References madness::Function< T, NDIM >::change_tree_state(), madness::Function< T, NDIM >::fence(), and madness::redundant.

◆ make_redundant() [2/2]

template<typename T , std::size_t NDIM>
func madness::Function< T, NDIM >::make_redundant ( true  )

◆ map_and_mirror()

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::map_and_mirror ( const Function< T, NDIM > &  f,
const std::vector< long > &  map,
const std::vector< long > &  mirror,
bool  fence 
)
inline

This is replaced with mirror(map(f)) ... private.

first map then mirror! mirror is similar to mapdim, but maps from x to -x, y to -y, and so on Example: mirror a 3d function on the xy plane: mirror={1,1,-1} Example: c4 rotation of a 3d function around the z axis: x->y, y->-x, z->z: map(1,0,2); mirror(-1,1,1)

Parameters
[in]maparray holding dimensions
[in]mirrorarray of -1 and 1, corresponding to mirror or not

References madness::f, madness::Function< T, NDIM >::fence(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::MADNESS_ASSERT(), madness::Function< T, NDIM >::mirror(), NDIM, PROFILE_MEMBER_FUNC, and VERIFY_TREE.

Referenced by madness::map_and_mirror().

◆ mapdim()

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::mapdim ( const Function< T, NDIM > &  f,
const std::vector< long > &  map,
bool  fence 
)
inline

◆ max_depth()

template<typename T , std::size_t NDIM>
std::size_t madness::Function< T, NDIM >::max_depth ( ) const
inline

Returns the maximum depth of the function tree ... collective global sum.

References madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

Referenced by testNavierStokes().

◆ max_local_depth()

template<typename T , std::size_t NDIM>
std::size_t madness::Function< T, NDIM >::max_local_depth ( ) const
inline

Returns the maximum local depth of the function tree ... no communications.

References madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

◆ max_nodes()

template<typename T , std::size_t NDIM>
std::size_t madness::Function< T, NDIM >::max_nodes ( ) const
inline

Returns the max number of nodes on a processor.

References madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

Referenced by test_gaussian_num_coeffs().

◆ min_nodes()

template<typename T , std::size_t NDIM>
std::size_t madness::Function< T, NDIM >::min_nodes ( ) const
inline

Returns the min number of nodes on a processor.

References madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

◆ mirror()

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::mirror ( const Function< T, NDIM > &  f,
const std::vector< long > &  mirrormap,
bool  fence 
)
inline

This is replaced with mirror(f) ... private.

similar to mapdim, but maps from x to -x, y to -y, and so on Example: mirror a 3d function on the xy plane: mirror={1,1,-1}

Parameters
[in]mirrorarray of -1 and 1, corresponding to mirror or not

References madness::f, madness::Function< T, NDIM >::fence(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::MADNESS_ASSERT(), NDIM, PROFILE_MEMBER_FUNC, and VERIFY_TREE.

Referenced by madness::Function< T, NDIM >::map_and_mirror(), and madness::mirror().

◆ mul_on_demand()

template<typename T , std::size_t NDIM>
template<typename L , typename R >
void madness::Function< T, NDIM >::mul_on_demand ( const Function< L, NDIM > &  f,
const Function< R, NDIM > &  g,
bool  fence = true 
)
inline

Same as operator* but with optional fence and no automatic reconstruction.

f or g are on-demand functions

References madness::f, madness::Function< T, NDIM >::fence(), madness::Function< T, NDIM >::g(), madness::Function< T, NDIM >::impl, madness::FunctionImpl< T, NDIM >::is_on_demand(), and MADNESS_EXCEPTION.

◆ multi_to_multi_op_values()

template<typename T , std::size_t NDIM>
template<typename opT >
void madness::Function< T, NDIM >::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 
)
inline

apply op on the input vector yielding an output vector of functions

(*this) is just a dummy Function to be able to call internal methods in FuncImpl

Parameters
[in]opthe operator working on vin
[in]vinvector of input Functions
[out]voutvector of output Functions vout = op(vin)

References madness::Function< T, NDIM >::fence(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_initialized(), op(), VERIFY_TREE, and madness::Function< T, NDIM >::verify_tree().

Referenced by madness::multi_to_multi_op_values().

◆ multiop_values()

template<typename T , std::size_t NDIM>
template<typename opT >
Function<T,NDIM>& madness::Function< T, NDIM >::multiop_values ( const opT &  op,
const std::vector< Function< T, NDIM > > &  vf 
)
inline

◆ norm2()

template<typename T , std::size_t NDIM>
double madness::Function< T, NDIM >::norm2 ( ) const
inline

Returns the 2-norm of the function ... global sum ... works in either basis.

See comments for err() w.r.t. applying to many functions.

References madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::local, PROFILE_MEMBER_FUNC, madness::Function< T, NDIM >::verify(), VERIFY_TREE, and madness::Function< T, NDIM >::verify_tree().

Referenced by madness::XCOperator< T, NDIM >::XCOperator(), madness::Solver< T, NDIM >::apply_hf_exchange(), check(), check_operator_multiplications_3d(), compute_energy(), converge(), converge2s(), doit(), GygiPot::ESP(), sgl_guess::get_wf(), madness::Function< T, NDIM >::if(), madness::MP2::increment(), iterate(), iterate_excite(), iterate_ground(), madness::CC2::iterate_pair(), main(), madness::EigSolver< T, NDIM >::multi_solve(), MiniDFT::orthogonalize(), orthogonalize(), SurfaceMoleculeInteraction::perturbed_molecular_pot(), Calculation::project_potential_basis(), propagate(), run(), madness::OEP::selftest(), simple_example(), madness::EigSolver< T, NDIM >::solve(), madness::MP2::solve_residual_equations(), test(), test_adaptive_tree(), test_add(), test_combined_operators(), madness::Znemo::test_compute_current_density(), test_converged_function(), test_convolution(), test_coulomb(), test_exchange(), test_hf_be(), test_integration(), test_Kcommutator(), test_kinetic(), test_lowrank_function(), madness::Diamagnetic_potential_factor::test_lz_commutator(), test_op(), test_opdir(), test_recursive_application(), test_rot(), test_U_el(), test_vector_composite(), and testNavierStokes().

◆ norm2sq_local()

template<typename T , std::size_t NDIM>
double madness::Function< T, NDIM >::norm2sq_local ( ) const
inline

Returns the square of the norm of the local function ... no communication.

Works in either basis

References madness::Function< T, NDIM >::impl, PROFILE_MEMBER_FUNC, and madness::Function< T, NDIM >::verify().

◆ norm_tree()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::norm_tree ( bool  fence = true) const
inline

◆ operator()() [1/2]

template<typename T , std::size_t NDIM>
T madness::Function< T, NDIM >::operator() ( const coordT xuser) const
inline

Evaluates the function at a point in user coordinates. Collective operation.

Throws if function is not initialized.

This function calls eval, blocks until the result is available and then broadcasts the result to everyone. Therefore, if you are evaluating many points in parallel it is vastly less efficient than calling eval directly, saving the futures, and then forcing all of the results.

References madness::Function< T, NDIM >::eval(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::is_reconstructed(), PROFILE_MEMBER_FUNC, madness::Function< T, NDIM >::reconstruct(), T(), and madness::Function< T, NDIM >::verify().

◆ operator()() [2/2]

template<typename T , std::size_t NDIM>
T madness::Function< T, NDIM >::operator() ( double  x,
double  y = 0,
double  z = 0,
double  xx = 0,
double  yy = 0,
double  zz = 0 
) const
inline

Evaluates the function at a point in user coordinates. Collective operation.

See "operator()(const coordT& xuser)" for more info

References NDIM.

◆ operator*=()

template<typename T , std::size_t NDIM>
template<typename Q >
IsSupported<TensorTypeData<Q>, Function<T,NDIM> >::type& madness::Function< T, NDIM >::operator*= ( const Q  q)
inline

Inplace scaling by a constant.

Using operator notation forces a global fence after every operation

References PROFILE_MEMBER_FUNC, q(), and madness::Function< T, NDIM >::scale().

◆ operator+=()

template<typename T , std::size_t NDIM>
template<typename Q >
Function<T,NDIM>& madness::Function< T, NDIM >::operator+= ( const Function< Q, NDIM > &  other)
inline

Inplace addition of functions in the wavelet basis.

Using operator notation forces a global fence after every operation. Functions don't need to be compressed, it's the caller's responsibility to choose an appropriate state with performance, usually compressed for 3d, reconstructed for 6d)

References madness::Function< T, NDIM >::compress(), madness::Function< T, NDIM >::gaxpy(), madness::Function< T, NDIM >::get_impl(), madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::MADNESS_ASSERT(), NDIM, PROFILE_MEMBER_FUNC, Q(), madness::Function< T, NDIM >::reconstruct(), T(), VERIFY_TREE, and madness::Function< T, NDIM >::verify_tree().

◆ operator-=()

template<typename T , std::size_t NDIM>
template<typename Q >
Function<T,NDIM>& madness::Function< T, NDIM >::operator-= ( const Function< Q, NDIM > &  other)
inline

◆ operator=()

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::operator= ( const Function< T, NDIM > &  f)
inline

Assignment is shallow. No communication, works in either basis.

References madness::f, madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

◆ print_info()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::print_info ( ) const
inline

Print a summary of the load balancing info.

This is serial and VERY expensive

References madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

◆ print_size()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::print_size ( const std::string  name) const
inline

◆ print_tree()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::print_tree ( std::ostream &  os = std::cout) const
inline

Process 0 prints a summary of all nodes in the tree (collective)

References madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

Referenced by main().

◆ print_tree_graphviz()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::print_tree_graphviz ( std::ostream &  os = std::cout) const
inline

Process 0 prints a graphviz-formatted output of all nodes in the tree (collective)

References madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

◆ print_tree_json()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::print_tree_json ( std::ostream &  os = std::cout) const
inline

same as print_tree() but produces JSON-formatted string

Warning
enclose the result in braces to make it a valid JSON object

References madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

◆ project_out()

template<typename T , std::size_t NDIM>
template<typename R , size_t LDIM>
Function<TENSOR_RESULT_TYPE(T,R),NDIM-LDIM> madness::Function< T, NDIM >::project_out ( const Function< R, LDIM > &  g,
const int  dim 
) const
inline

◆ reconstruct() [1/2]

template<typename T , std::size_t NDIM>
this madness::Function< T, NDIM >::reconstruct ( )

◆ reconstruct() [2/2]

template<typename T , std::size_t NDIM>
const Function<T,NDIM>& madness::Function< T, NDIM >::reconstruct ( bool  fence = true) const
inline

Reconstructs the function, transforming into scaling function basis. Possible non-blocking comm.

By default fence=true meaning that this operation completes before returning, otherwise if fence=false it returns without fencing and the user must invoke world.gop.fence() to assure global completion before using the function for other purposes.

Noop if already reconstructed or if not initialized.

Since reconstruction/compression do not discard information we define them as const ... "logical constness" not "bitwise constness".

References madness::Function< T, NDIM >::change_tree_state(), madness::Function< T, NDIM >::fence(), and madness::reconstructed.

Referenced by madness::Derivative< T, NDIM >::Derivative(), APPLY(), madness::APPLY(), madness::SCF::APPLY(), madness::apply(), madness::apply_1d_realspace_push(), apply_potential(), madness::binary_op(), Calculation::calc_optimal_coeffs(), chin_chen(), madness::SCF::do_plots(), MiniDFT::doit(), dostuff(), madness::Function< T, NDIM >::gaxpy_ext(), madness::Solver< T, NDIM >::initial_guess(), madness::SCF::initial_guess(), main(), madness::smooth< T, NDIM >::make_density_mask(), madness::SCF::make_lda_potential(), madness::smooth< T, NDIM >::make_mask(), madness::PotentialManager::make_nuclear_potential(), Calculation::make_reference(), Calculation::make_Upsi(), madness::operator*(), madness::operator+(), madness::Function< T, NDIM >::operator+=(), madness::operator-(), madness::Function< T, NDIM >::operator-=(), madness::XCOperator< T, NDIM >::prep_xc_args(), madness::project(), propagate(), run(), solve(), test_add(), test_coulomb(), test_he_potential(), test_hf_be(), test_hf_he(), test_periodic(), test_periodic1(), test_periodic2(), test_periodic_bsh(), test_xc2(), testbsh(), testPeriodicCoulomb3d(), and madness::CCConvolutionOperator< T, NDIM >::update_elements().

◆ reduce_rank()

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::reduce_rank ( const double  thresh = 0.0,
const bool  fence = true 
)
inline

◆ refine()

template<typename T , std::size_t NDIM>
const Function<T,NDIM>& madness::Function< T, NDIM >::refine ( bool  fence = true) const
inline

Inplace autorefines the function using same test as for squaring.

return this for chaining

References madness::Function< T, NDIM >::fence(), and madness::Function< T, NDIM >::refine_general().

Referenced by dostuff(), madness::Nemo::gradient(), and madness::smooth< T, NDIM >::smooth_density_from_orbitals().

◆ refine_general()

template<typename T , std::size_t NDIM>
template<typename opT >
void madness::Function< T, NDIM >::refine_general ( const opT &  op,
bool  fence = true 
) const
inline

◆ replicate()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::replicate ( bool  fence = true) const
inline

replicate this function, generating a unique pmap

References madness::Function< T, NDIM >::fence(), madness::Function< T, NDIM >::impl, and madness::Function< T, NDIM >::verify().

Referenced by test_replicate().

◆ replicate_low_dim_functions()

template<typename T , std::size_t NDIM>
func madness::Function< T, NDIM >::replicate_low_dim_functions ( true  )

◆ scale()

template<typename T , std::size_t NDIM>
template<typename Q >
Function<T,NDIM>& madness::Function< T, NDIM >::scale ( const Q  q,
bool  fence = true 
)
inline

Inplace, scale the function by a constant. No communication except for optional fence.

Works in either basis. Returns reference to this for chaining.

References madness::Function< T, NDIM >::fence(), madness::Function< T, NDIM >::impl, PROFILE_MEMBER_FUNC, q(), madness::Function< T, NDIM >::verify(), VERIFY_TREE, and madness::Function< T, NDIM >::verify_tree().

Referenced by DirichletCondIntOp::action(), WSTFunctional::apply_xc(), madness::DFT< T, NDIM >::calculate_tot_xc_energy(), madness::Coulomb< T, NDIM >::compute_density(), madness::Znemo::compute_nemo_density(), madness::Coulomb< T, NDIM >::compute_potential(), madness::Zcis::compute_potentials(), madness::Diamagnetic_potential_factor::compute_R_times_T_commutator_scalar_term_numerically(), madness::Solver< T, NDIM >::compute_rho_slow(), converge(), converge2s(), madness::PotentialManager::core_projection(), madness::SCF::do_plots(), MiniDFT::doit(), doit(), madness::Nemo::get_coulomb_potential(), sgl_guess::get_wf(), madness::Nemo::gradient(), madness::MP2::increment(), madness::Solver< T, NDIM >::initial_guess(), madness::SCF::initial_guess(), madness::Zcis::iterate(), iterate(), iterate_excite(), iterate_ground(), iterate_xy(), main(), make_density(), DF::make_fermi_potential(), madness::Zcis::make_guess(), SurfaceMoleculeInteraction::make_surfcharge(), madness::XCFunctionalLDA< T, NDIM >::op_r(), madness::Function< T, NDIM >::operator*=(), MiniDFT::orthogonalize(), orthogonalize(), run(), madness::EigSolver< T, NDIM >::solve(), madness::SCF::solve(), madness::MP2::solve_coupled_equations(), madness::MP2::solve_residual_equations(), squaremod(), test_adaptive_tree(), test_add(), madness::Znemo::test_compute_current_density(), test_he_potential(), test_hf_be(), test_Kcommutator(), test_kinetic(), test_modified(), and test_recursive_application().

◆ set_autorefine()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::set_autorefine ( bool  value,
bool  fence = true 
)
inline

Sets the value of the autorefine flag. Optional global fence.

A fence is required to ensure consistent global state.

References madness::Function< T, NDIM >::fence(), madness::Function< T, NDIM >::impl, PROFILE_MEMBER_FUNC, and madness::Function< T, NDIM >::verify().

◆ set_functor()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::set_functor ( const std::shared_ptr< FunctionFunctorInterface< T, NDIM > >  functor)
inline

Replace the current functor with the provided new one.

presumably the new functor will be a CompositeFunctor, which will change the behavior of the function: multiply the functor with the function

References madness::print().

◆ set_impl() [1/2]

template<typename T , std::size_t NDIM>
template<typename R >
void madness::Function< T, NDIM >::set_impl ( const Function< R, NDIM > &  f,
bool  zero = true 
)
inline

Replace current FunctionImpl with a new one using the same parameters & map as f.

If zero is true the function is initialized to zero, otherwise it is empty

References madness::f, madness::WorldGopInterface::fence(), madness::World::gop, madness::Function< T, NDIM >::impl, and madness::Function< T, NDIM >::world().

◆ set_impl() [2/2]

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::set_impl ( const std::shared_ptr< FunctionImpl< T, NDIM > > &  impl)
inline

◆ set_thresh()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::set_thresh ( double  value,
bool  fence = true 
)
inline

◆ size()

template<typename T , std::size_t NDIM>
std::size_t madness::Function< T, NDIM >::size ( ) const
inline

Returns the number of coefficients in the function ... collective global sum.

References madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

Referenced by madness::hartree_product(), main(), run(), solve(), and trotter().

◆ square()

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::square ( bool  fence = true)
inline

◆ standard()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::standard ( bool  fence = true)
inline

Converts the function standard compressed form. Possible non-blocking comm.

By default fence=true meaning that this operation completes before returning, otherwise if fence=false it returns without fencing and the user must invoke world.gop.fence() to assure global completion before using the function for other purposes.

Must be already compressed.

References madness::Function< T, NDIM >::change_tree_state(), madness::compressed, and madness::Function< T, NDIM >::fence().

Referenced by madness::apply(), and madness::hartree_product().

◆ store()

template<typename T , std::size_t NDIM>
template<typename Archive >
void madness::Function< T, NDIM >::store ( Archive &  ar) const
inline

Stores the function to an archive.

Archive can be sequential or parallel.

The & operator for serializing will only work with parallel archives.

References madness::Function< T, NDIM >::impl, madness::Function< T, NDIM >::k(), NDIM, PROFILE_MEMBER_FUNC, and madness::Function< T, NDIM >::verify().

◆ sum() [1/2]

template<typename T , std::size_t NDIM>
impl world gop madness::Function< T, NDIM >::sum ( local  )

◆ sum() [2/2]

template<typename T , std::size_t NDIM>
impl world gop madness::Function< T, NDIM >::sum ( local  )

◆ sum_down()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::sum_down ( bool  fence = true) const
inline

◆ TENSOR_RESULT_TYPE() [1/5]

template<typename T , std::size_t NDIM>
template<typename R >
madness::Function< T, NDIM >::TENSOR_RESULT_TYPE ( T  ,
R   
) const

Returns the inner product.

Not efficient for computing multiple inner products

Parameters
[in]gFunction, optionally on-demand

References PROFILE_MEMBER_FUNC.

◆ TENSOR_RESULT_TYPE() [2/5]

template<typename T , std::size_t NDIM>
template<typename R >
madness::Function< T, NDIM >::TENSOR_RESULT_TYPE ( T  ,
R   
) const

Returns local part of inner product ... throws if both not compressed.

Referenced by madness::Function< T, NDIM >::project_out().

◆ TENSOR_RESULT_TYPE() [3/5]

template<typename T , std::size_t NDIM>
template<typename R >
madness::Function< T, NDIM >::TENSOR_RESULT_TYPE ( T  ,
R   
) const

Returns the inner product for one on-demand function.

It does work, but it might not give you the precision you expect. The assumption is that the function g returns proper sum coefficients on the MRA tree of this. This might not be the case if g is constructed with an implicit multiplication, e.g. result = <this|g>, with g = 1/r12 | gg>

Parameters
[in]gon-demand function

References madness::Function< T, NDIM >::g(), madness::Function< T, NDIM >::is_on_demand(), and madness::Function< T, NDIM >::MADNESS_ASSERT().

◆ TENSOR_RESULT_TYPE() [4/5]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::TENSOR_RESULT_TYPE ( T  ,
R   
)

◆ TENSOR_RESULT_TYPE() [5/5]

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::TENSOR_RESULT_TYPE ( T  ,
R   
)

◆ thresh()

template<typename T , std::size_t NDIM>
double madness::Function< T, NDIM >::thresh ( ) const
inline

◆ trace()

template<typename T , std::size_t NDIM>
T madness::Function< T, NDIM >::trace ( ) const
inline

◆ trace_local()

template<typename T , std::size_t NDIM>
T madness::Function< T, NDIM >::trace_local ( ) const
inline

Returns local contribution to int(f(x),x) ... no communication.

In the wavelet basis this is just the coefficient of the first scaling function which is a constant. In the scaling function basis we must add up contributions from each box.

References madness::Function< T, NDIM >::impl, PROFILE_MEMBER_FUNC, VERIFY_TREE, and madness::Function< T, NDIM >::verify_tree().

◆ tree_size()

template<typename T , std::size_t NDIM>
std::size_t madness::Function< T, NDIM >::tree_size ( ) const
inline

Returns the number of nodes in the function tree ... collective global sum.

References madness::Function< T, NDIM >::impl, and PROFILE_MEMBER_FUNC.

Referenced by main(), solve(), and test_conversion().

◆ truncate()

template<typename T , std::size_t NDIM>
Function<T,NDIM>& madness::Function< T, NDIM >::truncate ( double  tol = 0.0,
bool  fence = true 
)
inline

Truncate the function with optional fence. Compresses with fence if not compressed.

If the truncation threshold is less than or equal to zero the default value specified when the function was created is used. If the function is not initialized, it just returns.

Returns this for chaining.

Parameters
[in]tolTolerance for truncating the coefficients. Default 0.0 means use the implementation's member value thresh instead.
[in]fenceDo fence

References madness::Function< T, NDIM >::fence(), madness::Function< T, NDIM >::impl, PROFILE_MEMBER_FUNC, madness::Function< T, NDIM >::verify(), VERIFY_TREE, and madness::Function< T, NDIM >::verify_tree().

Referenced by SVPESolver::SVPESolver(), DirichletCondIntOp::action(), RealFuncLinearOp::action(), madness::Solver< T, NDIM >::apply_hf_exchange3(), madness::Solver< T, NDIM >::apply_hf_exchange4(), apply_potential(), madness::SCF::apply_potential(), madness::NuclearCorrelationFactor::apply_U(), madness::PseudoNuclearCorrelationFactor::apply_U(), apply_U_ncf(), madness::XCOperator< T, NDIM >::apply_xc_kernel(), applyexpLt(), applyN(), Calculation::calc_optimal_coeffs(), madness::HartreeFock< T, NDIM >::calculate_coulomb_energy(), chin_chen(), madness::Coulomb< T, NDIM >::compute_density(), madness::Znemo::compute_density(), madness::Zcis::compute_potentials(), madness::Diamagnetic_potential_factor::compute_R_times_T_commutator_scalar_term_numerically(), madness::Solver< T, NDIM >::compute_rho(), madness::EigSolver< T, NDIM >::compute_rho(), madness::Solver< T, NDIM >::compute_rho_slow(), madness::Znemo::compute_spin_density(), madness::Diamagnetic_potential_factor::compute_U2(), converge(), converge2s(), madness::CCPairFunction< T, NDIM >::convert_to_pure_no_op_inplace(), madness::F12Potentials::convolve_with_fg(), cplxfunc1(), cplxfunc2(), madness::SCF::do_plots(), doit(), expV(), madness::MP2::increment(), madness::SCF::initial_guess(), madness::NuclearCorrelationFactor::initialize(), DF::iterate(), iterate(), iterate_excite(), iterate_ground(), madness::CC2::iterate_pair(), iterate_xy(), main(), madness::SCF::make_density(), make_density(), make_exp(), madness::Solver< T, NDIM >::make_nuclear_charge_density_impl(), madness::PotentialManager::make_nuclear_potential(), Calculation::make_reference(), Calculation::make_Upsi(), madness::HartreeFockExchangeOp< T, NDIM >::op_o(), madness::CCPairFunction< T, NDIM >::op_pure_to_pure(), madness::CCConvolutionOperator< T, NDIM >::operator()(), madness::StrongOrthogonalityProjector< T, NDIM >::operator()(), Calculation::project_potential_basis(), propagate(), q_c(), realfunc1(), realfunc2(), run(), SVPESolver::solve(), madness::EigSolver< T, NDIM >::solve(), madness::SCF::solve(), madness::MP2::solve_coupled_equations(), DF::solve_occupied(), madness::MP2::solve_residual_equations(), test_he_potential(), test_hf_he(), test_hydro(), test_modified(), test_multiply(), testbsh(), testPeriodicCoulomb3d(), trotter(), madness::NonlinearSolverND< NDIM >::update(), and madness::CCConvolutionOperator< T, NDIM >::update_elements().

◆ unary_op_coeffs()

template<typename T , std::size_t NDIM>
template<typename Q , typename opT >
Function<typename opT::resultT,NDIM>& madness::Function< T, NDIM >::unary_op_coeffs ( const Function< Q, NDIM > &  func,
const opT &  op,
bool  fence 
)
inline

◆ unaryop() [1/2]

template<typename T , std::size_t NDIM>
template<typename opT >
void madness::Function< T, NDIM >::unaryop ( const opT &  op,
bool  fence = true 
)
inline

◆ unaryop() [2/2]

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::unaryop ( T(*)(T f)
inline

◆ unaryop_coeff()

template<typename T , std::size_t NDIM>
template<typename opT >
void madness::Function< T, NDIM >::unaryop_coeff ( const opT &  op,
bool  fence = true 
)
inline

◆ unaryop_node()

template<typename T , std::size_t NDIM>
template<typename opT >
void madness::Function< T, NDIM >::unaryop_node ( const opT &  op,
bool  fence = true 
)
inline

◆ verify()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::verify ( ) const
inline

Asserts that the function is initialized.

References madness::Function< T, NDIM >::impl, and madness::Function< T, NDIM >::MADNESS_ASSERT().

Referenced by madness::Function< T, NDIM >::add_scalar(), madness::Function< T, NDIM >::broaden(), madness::Function< T, NDIM >::chop_at_level(), madness::Function< T, NDIM >::depthpt(), madness::Function< T, NDIM >::distribute(), madness::Function< T, NDIM >::err(), madness::Function< T, NDIM >::errsq_local(), madness::Function< T, NDIM >::eval(), madness::Function< T, NDIM >::eval_cube(), madness::Function< T, NDIM >::eval_local_only(), madness::Function< T, NDIM >::evaldepthpt(), madness::Function< T, NDIM >::evalR(), madness::Function< T, NDIM >::gaxpy(), madness::Function< T, NDIM >::gaxpy_oop(), madness::Function< T, NDIM >::get_impl(), madness::Function< T, NDIM >::get_pmap(), madness::Function< T, NDIM >::k(), madness::mul_sparse(), madness::Function< T, NDIM >::norm2(), madness::Function< T, NDIM >::norm2sq_local(), madness::Function< T, NDIM >::norm_tree(), madness::Function< T, NDIM >::operator()(), madness::Function< T, NDIM >::project_out(), madness::Function< T, NDIM >::reduce_rank(), madness::Function< T, NDIM >::refine_general(), madness::Function< T, NDIM >::replicate(), madness::Function< T, NDIM >::scale(), madness::Function< T, NDIM >::set_autorefine(), madness::Function< T, NDIM >::set_thresh(), madness::Function< T, NDIM >::store(), madness::Function< T, NDIM >::sum_down(), madness::Function< T, NDIM >::truncate(), madness::Function< T, NDIM >::unaryop(), madness::Function< T, NDIM >::unaryop_coeff(), madness::Function< T, NDIM >::unaryop_node(), and madness::Function< T, NDIM >::world().

◆ verify_tree()

template<typename T , std::size_t NDIM>
void madness::Function< T, NDIM >::verify_tree ( ) const
inline

◆ vimpl()

template<typename T , std::size_t NDIM>
template<typename Q , std::size_t D>
static std::vector< std::shared_ptr< FunctionImpl<Q,D> > > madness::Function< T, NDIM >::vimpl ( const std::vector< Function< Q, D > > &  v)
inlinestatic

Returns vector of FunctionImpl pointers corresponding to vector of functions.

References madness::Function< T, NDIM >::get_impl(), PROFILE_MEMBER_FUNC, and v.

Referenced by madness::Function< T, NDIM >::vtransform().

◆ vmulXX()

template<typename T , std::size_t NDIM>
template<typename L , typename R >
void madness::Function< T, NDIM >::vmulXX ( const Function< L, NDIM > &  left,
const std::vector< Function< R, NDIM > > &  right,
std::vector< Function< T, NDIM > > &  result,
double  tol,
bool  fence 
)
inline

◆ vtransform()

template<typename T , std::size_t NDIM>
template<typename R , typename Q >
void madness::Function< T, NDIM >::vtransform ( const std::vector< Function< R, NDIM > > &  v,
const Tensor< Q > &  c,
std::vector< Function< T, NDIM > > &  vresult,
double  tol,
bool  fence = true 
)
inline

sparse transformation of a vector of functions ... private

References c, madness::Function< T, NDIM >::fence(), PROFILE_MEMBER_FUNC, v, and madness::Function< T, NDIM >::vimpl().

◆ world()

template<typename T , std::size_t NDIM>
World& madness::Function< T, NDIM >::world ( ) const
inline

Member Data Documentation

◆ const

template<typename T , std::size_t NDIM>
NDIM &g madness::Function< T, NDIM >::const
Initial value:
{
Function()
Default constructor makes uninitialized function. No communication.
Definition: mra.h:154
#define PROFILE_MEMBER_FUNC(classname)
Definition: worldprofile.h:210

◆ else

template<typename T , std::size_t NDIM>
madness::Function< T, NDIM >::else
Initial value:
{
const Function< T, NDIM > & change_tree_state(const TreeState finalstate, bool fence=true) const
changes tree state to given state
Definition: mra.h:785
@ redundant
s coeffs everywhere
Definition: funcdefaults.h:64

◆ func

template<typename T , std::size_t NDIM>
auto madness::Function< T, NDIM >::func =dynamic_cast<CompositeFunctorInterface<T,NDIM,LDIM>* >(g.get_impl()->get_functor().get())

◆ gstate

template<typename T , std::size_t NDIM>
TreeState madness::Function< T, NDIM >::gstate =g.get_impl()->get_tree_state()

◆ impl

template<typename T , std::size_t NDIM>
std::shared_ptr< FunctionImpl<T,NDIM> > madness::Function< T, NDIM >::impl
private

Referenced by madness::Function< T, NDIM >::abs(), madness::Function< T, NDIM >::abs_square(), madness::Function< T, NDIM >::add_scalar(), madness::Function< T, NDIM >::autorefine(), madness::Function< T, NDIM >::broaden(), madness::Function< T, NDIM >::change_tensor_type(), madness::Function< T, NDIM >::change_tree_state(), madness::Function< T, NDIM >::check_symmetry(), madness::Function< T, NDIM >::chop_at_level(), madness::Function< T, NDIM >::clear(), madness::Function< T, NDIM >::depthpt(), madness::Function< T, NDIM >::distribute(), madness::Function< T, NDIM >::do_hartree_product(), madness::Function< T, NDIM >::err(), madness::Function< T, NDIM >::errsq_local(), madness::Function< T, NDIM >::eval(), madness::Function< T, NDIM >::eval_cube(), madness::Function< T, NDIM >::eval_local_only(), madness::Function< T, NDIM >::evaldepthpt(), madness::Function< T, NDIM >::evalR(), madness::Function< T, NDIM >::fill_cuspy_tree(), madness::Function< T, NDIM >::fill_nuclear_cuspy_tree(), madness::Function< T, NDIM >::fill_tree(), madness::Function< T, NDIM >::gaxpy(), madness::Function< T, NDIM >::gaxpy_ext(), madness::Function< T, NDIM >::gaxpy_oop(), madness::Function< T, NDIM >::get_impl(), madness::Function< T, NDIM >::get_pmap(), madness::Function< T, NDIM >::if(), madness::Function< T, NDIM >::impl_initialized(), madness::Function< T, NDIM >::inner_adaptive(), madness::Function< T, NDIM >::inner_ext(), madness::Function< T, NDIM >::inner_ext_local(), madness::Function< T, NDIM >::is_compressed(), madness::Function< T, NDIM >::is_initialized(), madness::Function< T, NDIM >::is_nonstandard(), madness::Function< T, NDIM >::is_reconstructed(), madness::Function< T, NDIM >::k(), madness::Function< T, NDIM >::load(), madness::Function< T, NDIM >::map_and_mirror(), madness::Function< T, NDIM >::mapdim(), madness::Function< T, NDIM >::max_depth(), madness::Function< T, NDIM >::max_local_depth(), madness::Function< T, NDIM >::max_nodes(), madness::Function< T, NDIM >::min_nodes(), madness::Function< T, NDIM >::mirror(), madness::Function< T, NDIM >::mul_on_demand(), madness::Function< T, NDIM >::multi_to_multi_op_values(), madness::Function< T, NDIM >::multiop_values(), madness::Function< T, NDIM >::norm2(), madness::Function< T, NDIM >::norm2sq_local(), madness::Function< T, NDIM >::norm_tree(), madness::Function< T, NDIM >::operator()(), madness::Function< T, NDIM >::autorefine_square_op::operator()(), madness::Function< T, NDIM >::operator+=(), madness::Function< T, NDIM >::operator-=(), madness::Function< T, NDIM >::operator=(), madness::Function< T, NDIM >::print_info(), madness::Function< T, NDIM >::print_size(), madness::Function< T, NDIM >::print_tree(), madness::Function< T, NDIM >::print_tree_graphviz(), madness::Function< T, NDIM >::print_tree_json(), madness::Function< T, NDIM >::reduce_rank(), madness::Function< T, NDIM >::refine_general(), madness::Function< T, NDIM >::replicate(), madness::Function< T, NDIM >::scale(), madness::Function< T, NDIM >::set_autorefine(), madness::Function< T, NDIM >::set_impl(), madness::Function< T, NDIM >::set_thresh(), madness::Function< T, NDIM >::size(), madness::Function< T, NDIM >::square(), madness::Function< T, NDIM >::store(), madness::Function< T, NDIM >::sum_down(), madness::Function< T, NDIM >::thresh(), madness::Function< T, NDIM >::trace(), madness::Function< T, NDIM >::trace_local(), madness::Function< T, NDIM >::tree_size(), madness::Function< T, NDIM >::truncate(), madness::Function< T, NDIM >::unary_op_coeffs(), madness::Function< T, NDIM >::unaryop(), madness::Function< T, NDIM >::unaryop_coeff(), madness::Function< T, NDIM >::unaryop_node(), madness::Function< T, NDIM >::verify(), madness::Function< T, NDIM >::verify_tree(), and madness::Function< T, NDIM >::world().

◆ LDIM

template<typename T , std::size_t NDIM>
constexpr std::size_t madness::Function< T, NDIM >::LDIM =std::max(NDIM/2,std::size_t(1))
constexpr

◆ local

template<typename T , std::size_t NDIM>
return madness::Function< T, NDIM >::local

◆ state

template<typename T , std::size_t NDIM>
TreeState madness::Function< T, NDIM >::state =this->get_impl()->get_tree_state()

The documentation for this class was generated from the following files: