10#ifndef THYRA_MODEL_EVALUATOR_HELPERS_HPP 
   11#define THYRA_MODEL_EVALUATOR_HELPERS_HPP 
   14#include "Thyra_ModelEvaluator.hpp" 
   28RCP<ModelEvaluatorBase::InArgs<Scalar> >
 
   86  ,
const std::string &derivName
 
   95  ,
const std::string &derivName
 
  105template<
class Scalar>
 
  107  const std::string &modelEvalDescription,
 
  109  const std::string &deriv_name,
 
  111  const std::string &fnc_space_name,
 
  113  const std::string &var_space_name
 
  121template<
class Scalar>
 
  123  const std::string &modelEvalDescription,
 
  133template<
class Scalar>
 
  144template<
class Scalar>
 
  153template<
class Scalar>
 
  162template<
class Scalar>
 
  172template<
class Scalar>
 
  182template<
class Scalar>
 
  193template<
class Scalar>
 
  205template<
class Scalar>
 
  217template<
class Scalar>
 
  229template<
class Scalar>
 
  242#ifdef HAVE_THYRA_ME_POLYNOMIAL 
  246template<
class Scalar>
 
  256template<
class Scalar>
 
  277#include "Thyra_AssertOp.hpp" 
  278#include "Teuchos_Utils.hpp" 
  281template<
class Scalar>
 
  292template<
class Scalar>
 
  294Thyra::derivativeGradient(
 
  295  const RCP<MultiVectorBase<Scalar> > &grad
 
  298  return ModelEvaluatorBase::Derivative<Scalar>(
 
  300    ModelEvaluatorBase::DERIV_MV_GRADIENT_FORM
 
  305template<
class Scalar>
 
  307Thyra::create_DfDp_mv(
 
  308  const ModelEvaluator<Scalar>& model,
 
  310  ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
 
  314  return createMembers( model.get_f_space(), model.get_p_space(l)->dim() );
 
  318template<
class Scalar>
 
  320Thyra::create_DgDx_dot_mv(
 
  321  const ModelEvaluator<Scalar>& model,
 
  323  ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
 
  326  typedef ModelEvaluatorBase MEB;
 
  327  switch(orientation) {
 
  328    case MEB::DERIV_MV_BY_COL:
 
  330        MEB::DerivativeMultiVector<Scalar>(
 
  331          createMembers( model.get_g_space(j), model.get_x_space()->dim() )
 
  332          ,MEB::DERIV_MV_BY_COL
 
  334    case MEB::DERIV_TRANS_MV_BY_ROW:
 
  336        MEB::DerivativeMultiVector<Scalar>(
 
  337          createMembers( model.get_x_space(), model.get_g_space(j)->dim() )
 
  338          ,MEB::DERIV_TRANS_MV_BY_ROW
 
  347template<
class Scalar>
 
  349Thyra::create_DgDx_mv(
 
  350  const ModelEvaluator<Scalar>& model,
 
  352  ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
 
  355  return create_DgDx_dot_mv(model,j,orientation);
 
  359template<
class Scalar>
 
  361Thyra::create_DgDp_mv(
 
  362  const ModelEvaluator<Scalar>& model,
 
  365  ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
 
  368  typedef ModelEvaluatorBase MEB;
 
  369  switch(orientation) {
 
  370    case MEB::DERIV_MV_BY_COL:
 
  372        MEB::DerivativeMultiVector<Scalar>(
 
  373          createMembers( model.get_g_space(j), model.get_p_space(l)->dim() )
 
  374          ,MEB::DERIV_MV_BY_COL
 
  376    case MEB::DERIV_TRANS_MV_BY_ROW:
 
  378        MEB::DerivativeMultiVector<Scalar>(
 
  379          createMembers( model.get_p_space(l), model.get_g_space(j)->dim() )
 
  380          ,MEB::DERIV_TRANS_MV_BY_ROW
 
  389template<
class Scalar>
 
  392  const ModelEvaluatorBase::Derivative<Scalar> &deriv
 
  393  ,
const std::string &derivName
 
  397    deriv.getLinearOp().get()!=NULL, std::logic_error
 
  398    ,
"Error, LinearOpBase type not expected for " << derivName <<
"!" 
  400  return deriv.getDerivativeMultiVector();
 
  404template<
class Scalar>
 
  407  const ModelEvaluatorBase::Derivative<Scalar> &deriv
 
  408  ,
const std::string &derivName
 
  409  ,ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
 
  412  typedef ModelEvaluatorBase MEB;
 
  414    deriv.getLinearOp().get()!=NULL, std::logic_error
 
  415    ,
"Error, LinearOpBase type not expected for " << derivName <<
"!" 
  417  MEB::DerivativeMultiVector<Scalar>
 
  418    dmv = deriv.getDerivativeMultiVector();
 
  419  RCP<MultiVectorBase<Scalar> >
 
  420    mv = dmv.getMultiVector();
 
  423      dmv.getOrientation() != orientation, std::logic_error
 
  424      ,
"Error, the orientation " << 
toString(dmv.getOrientation()) << 
" is not the" 
  425      " expected orientation of " << 
toString(orientation)
 
  426      << 
" for " << derivName << 
"!" 
  433template<
class Scalar>
 
  434void Thyra::assertDerivSpaces(
 
  435  const std::string &modelEvalDescription,
 
  436  const ModelEvaluatorBase::Derivative<Scalar> &deriv,
 
  437  const std::string &deriv_name,
 
  438  const VectorSpaceBase<Scalar> &fnc_space,
 
  439  const std::string &fnc_space_name,
 
  440  const VectorSpaceBase<Scalar> &var_space,
 
  441  const std::string &var_space_name
 
  444  typedef ModelEvaluatorBase MEB;
 
  445  if (!
is_null(deriv.getLinearOp())) {
 
  446    const RCP<const LinearOpBase<Scalar> > lo = deriv.getLinearOp();
 
  449        modelEvalDescription,
 
  450        *lo->range(), deriv_name + 
".range()",
 
  451        fnc_space, fnc_space_name
 
  454        modelEvalDescription,
 
  455        *lo->domain(), deriv_name + 
".domain()",
 
  456        var_space, var_space_name
 
  460  else if(!
is_null(deriv.getMultiVector())) {
 
  461    const RCP<const LinearOpBase<Scalar> > mv = deriv.getMultiVector();
 
  462    switch(deriv.getMultiVectorOrientation()) {
 
  463      case MEB::DERIV_MV_BY_COL: {
 
  465          modelEvalDescription,
 
  466          *mv->range(), deriv_name + 
".range()",
 
  467          fnc_space, fnc_space_name
 
  470          modelEvalDescription,
 
  471          *mv->domain(), deriv_name + 
".domain()",
 
  472          var_space, var_space_name
 
  476      case MEB::DERIV_TRANS_MV_BY_ROW: {
 
  478          modelEvalDescription,
 
  479          *mv->range(), deriv_name + 
"^T.range()",
 
  480          var_space, var_space_name
 
  483          modelEvalDescription,
 
  484          *mv->domain(), deriv_name + 
"^T.domain()",
 
  485          fnc_space, fnc_space_name
 
  498template<
class Scalar>
 
  499void Thyra::assertInArgsOutArgsSetup(
 
  500  const std::string &modelEvalDescription,
 
  501  const ModelEvaluatorBase::InArgs<Scalar> &inArgs,
 
  502  const ModelEvaluatorBase::OutArgs<Scalar> &outArgs
 
  506  typedef ModelEvaluatorBase MEB;
 
  508  const int Ng = outArgs.Ng();
 
  509  const int Np = outArgs.Np();
 
  517    inArgs.Np() != outArgs.Np(), std::logic_error,
 
  518    "Error: The underlying model " << modelEvalDescription << 
" incorrectly\n" 
  519    "set inArgs.Np() = "<<inArgs.Np()<<
" != outArgs.Np() = " 
  525    inArgs.supports(MEB::IN_ARG_x_dot) && !inArgs.supports(MEB::IN_ARG_x),
 
  527    "Error: The underlying model " << modelEvalDescription << 
" supports\n" 
  528    "x_dot but does not support x!" 
  533    inArgs.supports(MEB::IN_ARG_x_dot) && !inArgs.supports(MEB::IN_ARG_t),
 
  535    "Error: The underlying model " << modelEvalDescription << 
" supports\n" 
  536    "x_dot but does not support t!" 
  542      ( outArgs.supports(MEB::OUT_ARG_W) || outArgs.supports(MEB::OUT_ARG_W_op) )
 
  544      !inArgs.supports(MEB::IN_ARG_x)
 
  547    "Error: The underlying model " << modelEvalDescription << 
" says that\n" 
  548    "it supports W and/or W_op but it does not support x!" 
  552      ( outArgs.supports(MEB::OUT_ARG_W) || outArgs.supports(MEB::OUT_ARG_W_op) )
 
  554      inArgs.supports(MEB::IN_ARG_x_dot)
 
  556      !( inArgs.supports(MEB::IN_ARG_alpha) && inArgs.supports(MEB::IN_ARG_beta) )
 
  559    "Error: The underlying model " << modelEvalDescription << 
" supports W and/or W_op\n" 
  560    "and x_dot but it does not support alpha and beta as InArgs! \n" 
  561    "If the model supports W and x_dot, then it can be interpreted as an implicit \n" 
  562    "ODE/DAE and therefore D(f)/D(x_dot) must be non-zero and therefore alpha must \n" 
  563    "be supported.  If, however, the model can be interpreted as an explicit ODE \n" 
  564    "then x_dot should not be supported at all." 
  567  for ( 
int l = 0; l < Np; ++l ) {
 
  571    for ( 
int j = 0; j < Ng; ++j ) {
 
  575        ( !outArgs.supports(MEB::OUT_ARG_DgDx_dot,j).none()
 
  576          && !inArgs.supports(MEB::IN_ARG_x_dot) ),
 
  578        "Error: The underlying model " << modelEvalDescription << 
" says that\n" 
  579        "it supports DgDx_dot("<<j<<
") but it does not support x_dot!" 
  584        ( !outArgs.supports(MEB::OUT_ARG_DgDx,j).none()
 
  585          && !inArgs.supports(MEB::IN_ARG_x) ),
 
  587        "Error: The underlying model " << modelEvalDescription << 
" says that\n" 
  588        "it supports DgDx("<<j<<
") but it does not support x!" 
  600template<
class Scalar>
 
  601void Thyra::assertInArgsEvalObjects(
 
  602  const ModelEvaluator<Scalar> &model,
 
  603  const ModelEvaluatorBase::InArgs<Scalar> &inArgs
 
  607  typedef ModelEvaluatorBase MEB;
 
  609  const std::string description = model.description();
 
  610  const int Np = inArgs.Np();
 
  612  model.createInArgs().assertSameSupport(inArgs);
 
  615  if ( inArgs.supports(MEB::IN_ARG_x_dot) && !
is_null(inArgs.get_x_dot()) ) {
 
  617      description, *inArgs.get_x_dot()->space(), *model.get_x_space() );
 
  621  if ( inArgs.supports(MEB::IN_ARG_x) && !
is_null(inArgs.get_x()) ) {
 
  623      description, *inArgs.get_x()->space(), *model.get_x_space() );
 
  627  for ( 
int l = 0; l < Np; ++l ) {
 
  628    if (!
is_null(inArgs.get_p(l))) {
 
  630        description, *inArgs.get_p(l)->space(), *model.get_p_space(l) );
 
  637template<
class Scalar>
 
  638void Thyra::assertOutArgsEvalObjects(
 
  639  const ModelEvaluator<Scalar> &model,
 
  640  const ModelEvaluatorBase::OutArgs<Scalar> &outArgs,
 
  641  const ModelEvaluatorBase::InArgs<Scalar> *inArgs
 
  645  typedef ScalarTraits<Scalar> ST;
 
  647  typedef ModelEvaluatorBase MEB;
 
  649  const std::string description = model.description();
 
  650  const int Ng = outArgs.Ng();
 
  651  const int Np = outArgs.Np();
 
  657  model.createOutArgs().assertSameSupport(outArgs);
 
  660  if ( outArgs.supports(MEB::OUT_ARG_f) && !
is_null(outArgs.get_f()) ) {
 
  662      description, *outArgs.get_f()->space(), *model.get_f_space() );
 
  666  if ( outArgs.supports(MEB::OUT_ARG_W) && !
is_null(outArgs.get_W()) ) {
 
  667    if (!
is_null(outArgs.get_W()->range())) {
 
  669        description, *outArgs.get_W()->range(), *model.get_f_space() );
 
  671        description, *outArgs.get_W()->domain(), *model.get_x_space() );
 
  676  if ( outArgs.supports(MEB::OUT_ARG_W_op) && !
is_null(outArgs.get_W_op()) ) {
 
  677    if (!
is_null(outArgs.get_W_op()->range())) {
 
  679        description, *outArgs.get_W_op()->range(), *model.get_f_space() );
 
  681        description, *outArgs.get_W_op()->domain(), *model.get_x_space() );
 
  691      ( outArgs.supports(MEB::OUT_ARG_W) && !
is_null(outArgs.get_W()) )
 
  693      ( outArgs.supports(MEB::OUT_ARG_W_op) && !
is_null(outArgs.get_W_op()) )
 
  697    if ( inArgs->supports(MEB::IN_ARG_alpha) && inArgs->supports(MEB::IN_ARG_beta) ) {
 
  705    else if ( inArgs->supports(MEB::IN_ARG_beta) ) {
 
  711  if (outArgs.supports(MEB::OUT_ARG_f)) {
 
  712    for ( 
int l = 0; l < Np; ++l ) {
 
  713      if (!outArgs.supports(MEB::OUT_ARG_DfDp,l).none()) {
 
  716          outArgs.get_DfDp(l), 
"DfDp("+TU::toString(l)+
")",
 
  717          *model.get_f_space(), 
"f_space",
 
  718          *model.get_p_space(l), 
"p_space("+TU::toString(l)+
")" 
  725  for ( 
int j = 0; j < Ng; ++j ) {
 
  726    if (!
is_null(outArgs.get_g(j))) {
 
  728        description, *outArgs.get_g(j)->space(), *model.get_g_space(j) );
 
  733  for ( 
int j = 0; j < Ng; ++j ) {
 
  734    if (!outArgs.supports(MEB::OUT_ARG_DgDx_dot,j).none()) {
 
  737        outArgs.get_DgDx_dot(j), 
"DgDx_dot("+TU::toString(j)+
")",
 
  738        *model.get_g_space(j), 
"g_space("+TU::toString(j)+
")",
 
  739        *model.get_x_space(), 
"x_space" 
  745  for ( 
int j = 0; j < Ng; ++j ) {
 
  746    if (!outArgs.supports(MEB::OUT_ARG_DgDx,j).none()) {
 
  749        outArgs.get_DgDx(j), 
"DgDx("+TU::toString(j)+
")",
 
  750        *model.get_g_space(j), 
"g_space("+TU::toString(j)+
")",
 
  751        *model.get_x_space(), 
"x_space" 
  757  for ( 
int j = 0; j < Ng; ++j ) {
 
  758    for ( 
int l = 0; l < Np; ++l ) {
 
  759      if (!outArgs.supports(MEB::OUT_ARG_DgDp,j,l).none()) {
 
  760        const std::string j_str = TU::toString(j);
 
  761        const std::string l_str = TU::toString(l);
 
  764          outArgs.get_DgDp(j,l), 
"DgDp("+j_str+
","+l_str+
")",
 
  765          *model.get_g_space(j), 
"g_space("+j_str+
")",
 
  766          *model.get_p_space(l), 
"p_space("+l_str+
")" 
  775template<
class Scalar>
 
  777  const ModelEvaluator<Scalar> &model
 
  778  ,
const VectorBase<Scalar> &x
 
  779  ,VectorBase<Scalar> *f
 
  782  typedef ModelEvaluatorBase MEB;
 
  783  MEB::InArgs<Scalar> inArgs = model.createInArgs();
 
  784  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
 
  787  model.evalModel(inArgs,outArgs);
 
  791template<
class Scalar>
 
  793  const ModelEvaluator<Scalar> &model
 
  794  ,
const VectorBase<Scalar> &x
 
  795  ,VectorBase<Scalar> *f
 
  796  ,LinearOpWithSolveBase<Scalar> *W
 
  800  typedef ModelEvaluatorBase MEB;
 
  802  MEB::InArgs<Scalar> inArgs = model.createInArgs();
 
  803  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
 
  810  model.evalModel(inArgs,outArgs);
 
  815template<
class Scalar>
 
  817  const ModelEvaluator<Scalar> &model
 
  818  ,
const VectorBase<Scalar> &x
 
  820  ,VectorBase<Scalar> *f
 
  823  typedef ModelEvaluatorBase MEB;
 
  824  MEB::InArgs<Scalar> inArgs = model.createInArgs();
 
  825  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
 
  827  if(inArgs.supports(MEB::IN_ARG_t)) inArgs.set_t(t);
 
  829  model.evalModel(inArgs,outArgs);
 
  833template<
class Scalar>
 
  835  const ModelEvaluator<Scalar> &model,
 
  837  const VectorBase<Scalar> &p_l,
 
  839  const Ptr<VectorBase<Scalar> > &g_j
 
  842  typedef ModelEvaluatorBase MEB;
 
  843  MEB::InArgs<Scalar> inArgs = model.createInArgs();
 
  844  MEB::OutArgs<Scalar> outArgs= model.createOutArgs();
 
  845  inArgs.set_p(l, Teuchos::rcpFromRef(p_l));
 
  846  outArgs.set_g(j, Teuchos::rcpFromRef(*g_j));
 
  847  model.evalModel(inArgs,outArgs);
 
  851template<
class Scalar>
 
  853  const ModelEvaluator<Scalar> &model,
 
  855  const VectorBase<Scalar> &p_l,
 
  858  VectorBase<Scalar> *g_j
 
  861  typedef ModelEvaluatorBase MEB;
 
  862  MEB::InArgs<Scalar> inArgs = model.createInArgs();
 
  863  MEB::OutArgs<Scalar> outArgs= model.createOutArgs();
 
  867  model.evalModel(inArgs,outArgs);
 
  871template<
class Scalar>
 
  872void Thyra::eval_g_DgDp(
 
  873  const ModelEvaluator<Scalar> &model,
 
  875  const VectorBase<Scalar> &p_l,
 
  877  const Ptr<VectorBase<Scalar> > &g_j,
 
  878  const ModelEvaluatorBase::Derivative<Scalar> &DgDp_j_l
 
  881  typedef ModelEvaluatorBase MEB;
 
  882  MEB::InArgs<Scalar> inArgs = model.createInArgs();
 
  883  MEB::OutArgs<Scalar> outArgs= model.createOutArgs();
 
  884  inArgs.set_p(l, Teuchos::rcpFromRef(p_l));
 
  886    outArgs.set_g(j, Teuchos::rcpFromPtr(g_j));
 
  888  if (!DgDp_j_l.isEmpty()) {
 
  889    outArgs.set_DgDp(j, l, DgDp_j_l);
 
  891  model.evalModel(inArgs,outArgs);
 
  895template<
class Scalar>
 
  897  const ModelEvaluator<Scalar> &model
 
  898  ,
const VectorBase<Scalar> &x_dot
 
  899  ,
const VectorBase<Scalar> &x
 
  900  ,
const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
 
  901  ,VectorBase<Scalar> *f
 
  905  typedef ModelEvaluatorBase MEB;
 
  907  MEB::InArgs<Scalar> inArgs = model.createInArgs();
 
  908  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
 
  912  if(inArgs.supports(MEB::IN_ARG_t))
 
  917  model.evalModel(inArgs,outArgs);
 
  922template<
class Scalar>
 
  924  const ModelEvaluator<Scalar> &model
 
  925  ,
const VectorBase<Scalar> &x_dot
 
  926  ,
const VectorBase<Scalar> &x
 
  927  ,
const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
 
  930  ,VectorBase<Scalar> *f
 
  931  ,LinearOpWithSolveBase<Scalar> *W
 
  935  typedef ModelEvaluatorBase MEB;
 
  937  MEB::InArgs<Scalar> inArgs = model.createInArgs();
 
  938  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
 
  942  if(inArgs.supports(MEB::IN_ARG_t))
 
  944  inArgs.set_alpha(alpha);
 
  945  inArgs.set_beta(beta);
 
  950  model.evalModel(inArgs,outArgs);
 
  955#ifdef HAVE_THYRA_ME_POLYNOMIAL 
  958template<
class Scalar>
 
  959void Thyra::eval_f_poly(
 
  960  const ModelEvaluator<Scalar> &model
 
  962  ,
const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
 
  967  typedef ModelEvaluatorBase MEB;
 
  969  MEB::InArgs<Scalar> inArgs = model.createInArgs();
 
  970  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
 
  973  if(inArgs.supports(MEB::IN_ARG_t))
 
  978  model.evalModel(inArgs,outArgs);
 
  983template<
class Scalar>
 
  984void Thyra::eval_f_poly(
 
  985  const ModelEvaluator<Scalar> &model
 
  987  ,
const VectorBase<Scalar> &x_poly
 
  988  ,
const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
 
  993  typedef ModelEvaluatorBase MEB;
 
  995  MEB::InArgs<Scalar> inArgs = model.createInArgs();
 
  996  MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
 
 1000  if(inArgs.supports(MEB::IN_ARG_t))
 
 1005  model.evalModel(inArgs,outArgs);
 
Base class for all linear operators that can support a high-level solve operation.
 
Simple aggregate class for a derivative object represented as a column-wise multi-vector or its trans...
 
Simple aggregate class that stores a derivative object as a general linear operator or as a multi-vec...
 
Concrete aggregate class for all input arguments computable by a ModelEvaluator subclass object.
 
Teuchos::ScalarTraits< Scalar >::magnitudeType ScalarMag
.
 
Concrete aggregate class for all output arguments computable by a ModelEvaluator subclass object.
 
EDerivativeMultiVectorOrientation
 
void assertOutArgsEvalObjects(const ModelEvaluator< Scalar > &model, const ModelEvaluatorBase::OutArgs< Scalar > &outArgs, const ModelEvaluatorBase::InArgs< Scalar > *inArgs=0)
Assert that the objects in an OutArgs object match a given model.
 
RCP< ModelEvaluatorBase::InArgs< Scalar > > clone(const ModelEvaluatorBase::InArgs< Scalar > &inArgs)
Create a clone of an InArgs object.
 
void assertInArgsOutArgsSetup(const std::string &modelEvalDescription, const ModelEvaluatorBase::InArgs< Scalar > &inArgs, const ModelEvaluatorBase::OutArgs< Scalar > &outArgs)
Assert that an InArgs and OutArgs object are setup consistently.
 
void assertInArgsEvalObjects(const ModelEvaluator< Scalar > &model, const ModelEvaluatorBase::InArgs< Scalar > &inArgs)
Assert that the objects in an InArgs object match a given model.
 
void assertDerivSpaces(const std::string &modelEvalDescription, const ModelEvaluatorBase::Derivative< Scalar > &deriv, const std::string &deriv_name, const VectorSpaceBase< Scalar > &fnc_space, const std::string &fnc_space_name, const VectorSpaceBase< Scalar > &var_space, const std::string &var_space_name)
Assert that that Thyra objects imbedded in a Derivative object matches its function and variable spac...
 
Pure abstract base interface for evaluating a stateless "model" that can be mapped into a number of d...
 
Interface for a collection of column vectors called a multi-vector.
 
Abstract interface for finite-dimensional dense vectors.
 
Abstract interface for objects that represent a space for vectors.
 
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
 
bool is_null(const std::shared_ptr< T > &p)
 
#define THYRA_ASSERT_VEC_SPACES(FUNC_NAME, VS1, VS2)
This is a very useful macro that should be used to validate that two vector spaces are compatible.
 
#define THYRA_ASSERT_VEC_SPACES_NAMES(FUNC_NAME, VS1, VS1_NAME, VS2, VS2_NAME)
Helper assertion macro.
 
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
 
std::string toString(const T &t)
 
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)