LinearOp< T > Class Template Reference

#include <LinearSolver.H>

Inheritance diagram for LinearOp< T >:

Inheritance graph
[legend]

List of all members.


Detailed Description

template<class T>
class LinearOp< T >

Operator class for Linear solver for solving L(phi) = rhs.

Public Member Functions

virtual ~LinearOp ()
virtual void residual (T &a_lhs, const T &a_phi, const T &a_rhs, bool a_homogeneous=false)=0
virtual void preCond (T &a_cor, const T &a_residual)=0
virtual void applyOp (T &a_lhs, const T &a_phi, bool a_homogeneous=false)=0
virtual void create (T &a_lhs, const T &a_rhs)=0
virtual void clear (T &a_lhs)
virtual void assign (T &a_lhs, const T &a_rhs)=0
virtual void assignLocal (T &a_lhs, const T &a_rhs)
virtual Real dotProduct (const T &a_1, const T &a_2)=0
virtual void mDotProduct (const T &a_1, const int a_sz, const T a_2[], Real a_mdots[])
virtual void incr (T &a_lhs, const T &a_x, Real a_scale)=0
virtual void axby (T &a_lhs, const T &a_x, const T &a_y, Real a_a, Real a_b)=0
virtual void scale (T &a_lhs, const Real &a_scale)=0
virtual Real norm (const T &a_rhs, int a_ord)=0
virtual Real dx () const
virtual void setToZero (T &a_lhs)=0
virtual void write (const T *a, const char *filename)

Constructor & Destructor Documentation

template<class T>
virtual LinearOp< T >::~LinearOp (  )  [inline, virtual]


Member Function Documentation

template<class T>
virtual void LinearOp< T >::residual ( T &  a_lhs,
const T &  a_phi,
const T &  a_rhs,
bool  a_homogeneous = false 
) [pure virtual]

Say you are solving L(phi) = rhs. Make a_lhs = L(a_phi) - a_rhs. If a_homogeneous is true, evaluate the operator using homogeneous boundary conditions.

Implemented in AMRNodeOp, AMRPoissonOp, MultilevelLinearOp< T >, NewPoissonOp, NewPoissonOp4, ViscousTensorOp, NWOViscousTensorOp, PoissonOp4, ResistivityOp, ViscousTensorOp, EBAMRPoissonOp, EBConductivityOp, EBPoissonOp, EBViscousTensorOp, NWOEBConductivityOp, NWOEBViscousTensorOp, and slowEBCO.

Referenced by MGLevelOp< LevelData< NodeFArrayBox > >::residualNF().

template<class T>
virtual void LinearOp< T >::preCond ( T &  a_cor,
const T &  a_residual 
) [pure virtual]

template<class T>
virtual void LinearOp< T >::applyOp ( T &  a_lhs,
const T &  a_phi,
bool  a_homogeneous = false 
) [pure virtual]

In the context of solving L(phi) = rhs, set a_lhs = L(a_phi). If a_homogeneous is true, evaluate the operator using homogeneous boundary conditions.

Implemented in AMRNodeOp, AMRPoissonOp, MultilevelLinearOp< T >, NewPoissonOp, NewPoissonOp4, ViscousTensorOp, NWOViscousTensorOp, PoissonOp4, ResistivityOp, ViscousTensorOp, EBAMRPoissonOp, EBConductivityOp, EBPoissonOp, EBViscousTensorOp, NWOEBConductivityOp, NWOEBViscousTensorOp, and slowEBCO.

template<class T>
virtual void LinearOp< T >::create ( T &  a_lhs,
const T &  a_rhs 
) [pure virtual]

template<class T>
virtual void LinearOp< T >::clear ( T &  a_lhs  )  [inline, virtual]

Opposite of create -- perform any operations required before lhs goes out of scope. In general, the only time this needs to be defined in a derived class is if the create() function called new. Otherwise, the default no-op function is sufficient.

Reimplemented in MultilevelLinearOp< T >.

Referenced by PetscSolver< T >::solve_mfree_private().

template<class T>
virtual void LinearOp< T >::assign ( T &  a_lhs,
const T &  a_rhs 
) [pure virtual]

template<class T>
virtual void LinearOp< T >::assignLocal ( T &  a_lhs,
const T &  a_rhs 
) [inline, virtual]

Reimplemented in AMRPoissonOp, and EBAMRPoissonOp.

template<class T>
virtual Real LinearOp< T >::dotProduct ( const T &  a_1,
const T &  a_2 
) [pure virtual]

template<class T>
virtual void LinearOp< T >::mDotProduct ( const T &  a_1,
const int  a_sz,
const T  a_2[],
Real  a_mdots[] 
) [inline, virtual]

Reimplemented in AMRPoissonOp.

template<class T>
virtual void LinearOp< T >::incr ( T &  a_lhs,
const T &  a_x,
Real  a_scale 
) [pure virtual]

template<class T>
virtual void LinearOp< T >::axby ( T &  a_lhs,
const T &  a_x,
const T &  a_y,
Real  a_a,
Real  a_b 
) [pure virtual]

template<class T>
virtual void LinearOp< T >::scale ( T &  a_lhs,
const Real a_scale 
) [pure virtual]

template<class T>
virtual Real LinearOp< T >::norm ( const T &  a_rhs,
int  a_ord 
) [pure virtual]

template<class T>
virtual Real LinearOp< T >::dx (  )  const [inline, virtual]

Return dx at this level of refinement

Reimplemented in AMRPoissonOp, ViscousTensorOp, NWOViscousTensorOp, ViscousTensorOp, EBAMRPoissonOp, EBConductivityOp, and NWOEBConductivityOp.

Referenced by PetscSolver< T >::define().

template<class T>
virtual void LinearOp< T >::setToZero ( T &  a_lhs  )  [pure virtual]

template<class T>
virtual void LinearOp< T >::write ( const T *  a,
const char *  filename 
) [inline, virtual]

Debugging aid for solvers. Print out a "T" to a file named "filename" default implementation is to print out a message saying "LinearOp::write not implemented"

Reimplemented in AMRPoissonOp, and MultilevelLinearOp< T >.


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

Generated on Fri Apr 5 04:25:09 2019 for Chombo + EB by  doxygen 1.5.5