Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions

Albany::KLResponseFunction Class Reference

A response function given by the KL decomposition of another response function. More...

#include <Albany_KLResponseFunction.hpp>

Inheritance diagram for Albany::KLResponseFunction:
Inheritance graph
[legend]
Collaboration diagram for Albany::KLResponseFunction:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 KLResponseFunction (const Teuchos::RCP< AbstractResponseFunction > &response, Teuchos::ParameterList &responseParams)
 Default constructor.
virtual ~KLResponseFunction ()
 Destructor.
virtual void setup ()
 Setup response function.
virtual Teuchos::RCP< const
Epetra_Map > 
responseMap () const
 Get the map associate with this response.
virtual bool isScalarResponse () const
 Is this response function "scalar" valued, i.e., has a replicated local response map.
virtual Teuchos::RCP
< Epetra_Operator > 
createGradientOp () const
 Create operator for gradient (e.g., dg/dx).
Deterministic evaluation functions

virtual void evaluateResponse (const double current_time, const Epetra_Vector *xdot, const Epetra_Vector *xdotdot, const Epetra_Vector &x, const Teuchos::Array< ParamVec > &p, Epetra_Vector &g)
 Evaluate responses.
virtual void evaluateTangent (const double alpha, const double beta, const double omega, const double current_time, bool sum_derivs, const Epetra_Vector *xdot, const Epetra_Vector *xdotdot, const Epetra_Vector &x, const Teuchos::Array< ParamVec > &p, ParamVec *deriv_p, const Epetra_MultiVector *Vxdot, const Epetra_MultiVector *Vxdotdot, const Epetra_MultiVector *Vx, const Epetra_MultiVector *Vp, Epetra_Vector *g, Epetra_MultiVector *gx, Epetra_MultiVector *gp)
 Evaluate tangent = dg/dx*dx/dp + dg/dxdot*dxdot/dp + dg/dp.
virtual void evaluateDerivative (const double current_time, const Epetra_Vector *xdot, const Epetra_Vector *xdotdot, const Epetra_Vector &x, const Teuchos::Array< ParamVec > &p, ParamVec *deriv_p, Epetra_Vector *g, const EpetraExt::ModelEvaluator::Derivative &dg_dx, const EpetraExt::ModelEvaluator::Derivative &dg_dxdot, const EpetraExt::ModelEvaluator::Derivative &dg_dxdotdot, const EpetraExt::ModelEvaluator::Derivative &dg_dp)
 Evaluate gradient = dg/dx, dg/dxdot, dg/dp.

Protected Member Functions

bool computeKL (const Stokhos::EpetraVectorOrthogPoly &sg_u, const int NumKL, Teuchos::Array< double > &evals, Teuchos::RCP< Epetra_MultiVector > &evecs)

Protected Attributes

Teuchos::RCP
< AbstractResponseFunction
response
 Response function we work with.
Teuchos::ParameterList responseParams
 Response parameters.
Teuchos::RCP
< Teuchos::FancyOStream > 
out
 Output stream;.
int num_kl
 Number of KL terms.

Private Member Functions

 KLResponseFunction (const KLResponseFunction &)
 Private to prohibit copying.
KLResponseFunctionoperator= (const KLResponseFunction &)
 Private to prohibit copying.

Detailed Description

A response function given by the KL decomposition of another response function.

It only defines the SG methods.

Definition at line 28 of file Albany_KLResponseFunction.hpp.


Constructor & Destructor Documentation

Albany::KLResponseFunction::KLResponseFunction ( const Teuchos::RCP< AbstractResponseFunction > &  response,
Teuchos::ParameterList &  responseParams 
)

Default constructor.

Albany::KLResponseFunction::~KLResponseFunction (  )  [virtual]

Destructor.

Definition at line 24 of file Albany_KLResponseFunction.cpp.

Albany::KLResponseFunction::KLResponseFunction ( const KLResponseFunction  )  [private]

Private to prohibit copying.


Member Function Documentation

virtual void Albany::KLResponseFunction::setup (  )  [inline, virtual]

Setup response function.

Implements Albany::AbstractResponseFunction.

Definition at line 40 of file Albany_KLResponseFunction.hpp.

Teuchos::RCP< const Epetra_Map > Albany::KLResponseFunction::responseMap (  )  const [virtual]

Get the map associate with this response.

Implements Albany::AbstractResponseFunction.

Definition at line 30 of file Albany_KLResponseFunction.cpp.

bool Albany::KLResponseFunction::isScalarResponse (  )  const [virtual]

Is this response function "scalar" valued, i.e., has a replicated local response map.

Implements Albany::AbstractResponseFunction.

Definition at line 44 of file Albany_KLResponseFunction.cpp.

Teuchos::RCP< Epetra_Operator > Albany::KLResponseFunction::createGradientOp (  )  const [virtual]

Create operator for gradient (e.g., dg/dx).

Implements Albany::AbstractResponseFunction.

Definition at line 37 of file Albany_KLResponseFunction.cpp.

void Albany::KLResponseFunction::evaluateResponse ( const double  current_time,
const Epetra_Vector *  xdot,
const Epetra_Vector *  xdotdot,
const Epetra_Vector &  x,
const Teuchos::Array< ParamVec > &  p,
Epetra_Vector &  g 
) [virtual]

Evaluate responses.

Implements Albany::AbstractResponseFunction.

Definition at line 51 of file Albany_KLResponseFunction.cpp.

void Albany::KLResponseFunction::evaluateTangent ( const double  alpha,
const double  beta,
const double  omega,
const double  current_time,
bool  sum_derivs,
const Epetra_Vector *  xdot,
const Epetra_Vector *  xdotdot,
const Epetra_Vector &  x,
const Teuchos::Array< ParamVec > &  p,
ParamVec deriv_p,
const Epetra_MultiVector *  Vxdot,
const Epetra_MultiVector *  Vxdotdot,
const Epetra_MultiVector *  Vx,
const Epetra_MultiVector *  Vp,
Epetra_Vector *  g,
Epetra_MultiVector *  gx,
Epetra_MultiVector *  gp 
) [virtual]

Evaluate tangent = dg/dx*dx/dp + dg/dxdot*dxdot/dp + dg/dp.

Implements Albany::AbstractResponseFunction.

Definition at line 63 of file Albany_KLResponseFunction.cpp.

void Albany::KLResponseFunction::evaluateDerivative ( const double  current_time,
const Epetra_Vector *  xdot,
const Epetra_Vector *  xdotdot,
const Epetra_Vector &  x,
const Teuchos::Array< ParamVec > &  p,
ParamVec deriv_p,
Epetra_Vector *  g,
const EpetraExt::ModelEvaluator::Derivative &  dg_dx,
const EpetraExt::ModelEvaluator::Derivative &  dg_dxdot,
const EpetraExt::ModelEvaluator::Derivative &  dg_dxdotdot,
const EpetraExt::ModelEvaluator::Derivative &  dg_dp 
) [virtual]

Evaluate gradient = dg/dx, dg/dxdot, dg/dp.

Implements Albany::AbstractResponseFunction.

Definition at line 88 of file Albany_KLResponseFunction.cpp.

KLResponseFunction& Albany::KLResponseFunction::operator= ( const KLResponseFunction  )  [private]

Private to prohibit copying.

bool Albany::KLResponseFunction::computeKL ( const Stokhos::EpetraVectorOrthogPoly &  sg_u,
const int  NumKL,
Teuchos::Array< double > &  evals,
Teuchos::RCP< Epetra_MultiVector > &  evecs 
) [protected]

Definition at line 271 of file Albany_KLResponseFunction.cpp.


Member Data Documentation

Response function we work with.

Definition at line 239 of file Albany_KLResponseFunction.hpp.

Teuchos::ParameterList Albany::KLResponseFunction::responseParams [protected]

Response parameters.

Definition at line 242 of file Albany_KLResponseFunction.hpp.

Teuchos::RCP<Teuchos::FancyOStream> Albany::KLResponseFunction::out [protected]

Output stream;.

Definition at line 245 of file Albany_KLResponseFunction.hpp.

Number of KL terms.

Definition at line 248 of file Albany_KLResponseFunction.hpp.


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