• Main Page
  • Related Pages
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

Albany_DistributedResponseFunction.hpp

Go to the documentation of this file.
00001 //*****************************************************************//
00002 //    Albany 2.0:  Copyright 2012 Sandia Corporation               //
00003 //    This Software is released under the BSD license detailed     //
00004 //    in the file "license.txt" in the top-level Albany directory  //
00005 //*****************************************************************//
00006 
00007 #ifndef ALBANY_DISTRIBUTED_RESPONSE_FUNCTION_HPP
00008 #define ALBANY_DISTRIBUTED_RESPONSE_FUNCTION_HPP
00009 
00010 #include "Albany_AbstractResponseFunction.hpp"
00011 #include "Stokhos_ProductEpetraOperator.hpp"
00012 #include "Stokhos_EpetraOperatorOrthogPoly.hpp"
00013 
00014 namespace Albany {
00015 
00022   class DistributedResponseFunction : 
00023     public AbstractResponseFunction {
00024   public:
00025   
00027     DistributedResponseFunction() {};
00028 
00030     virtual ~DistributedResponseFunction() {};
00031 
00033     virtual void evaluateGradient(
00034       const double current_time,
00035       const Epetra_Vector* xdot,
00036       const Epetra_Vector* xdotdot,
00037       const Epetra_Vector& x,
00038       const Teuchos::Array<ParamVec>& p,
00039       ParamVec* deriv_p,
00040       Epetra_Vector* g,
00041       Epetra_Operator* dg_dx,
00042       Epetra_Operator* dg_dxdot,
00043       Epetra_Operator* dg_dxdotdot,
00044       Epetra_MultiVector* dg_dp) = 0;
00045 
00046 #ifdef ALBANY_SG_MP
00047 
00048     virtual void evaluateSGGradient(
00049       const double current_time,
00050       const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00051       const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00052       const Stokhos::EpetraVectorOrthogPoly& sg_x,
00053       const Teuchos::Array<ParamVec>& p,
00054       const Teuchos::Array<int>& sg_p_index,
00055       const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00056       ParamVec* deriv_p,
00057       Stokhos::EpetraVectorOrthogPoly* sg_g,
00058       Stokhos::EpetraOperatorOrthogPoly* sg_dg_dx,
00059       Stokhos::EpetraOperatorOrthogPoly* sg_dg_dxdot,
00060       Stokhos::EpetraOperatorOrthogPoly* sg_dg_dxdotdot,
00061       Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dp) = 0;
00062 
00064     virtual void evaluateMPGradient(
00065       const double current_time,
00066       const Stokhos::ProductEpetraVector* mp_xdot,
00067       const Stokhos::ProductEpetraVector* mp_xdotdot,
00068       const Stokhos::ProductEpetraVector& mp_x,
00069       const Teuchos::Array<ParamVec>& p,
00070       const Teuchos::Array<int>& mp_p_index,
00071       const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00072       ParamVec* deriv_p,
00073       Stokhos::ProductEpetraVector* mp_g,
00074       Stokhos::ProductEpetraOperator* mp_dg_dx,
00075       Stokhos::ProductEpetraOperator* mp_dg_dxdot,
00076       Stokhos::ProductEpetraOperator* mp_dg_dxdotdot,
00077       Stokhos::ProductEpetraMultiVector* mp_dg_dp) = 0;
00078 #endif //ALBANY_SG_MP
00079 
00081 
00082 
00087     virtual bool isScalarResponse() const { return false; }
00088 
00090     virtual void evaluateDerivative(
00091       const double current_time,
00092       const Epetra_Vector* xdot,
00093       const Epetra_Vector* xdotdot,
00094       const Epetra_Vector& x,
00095       const Teuchos::Array<ParamVec>& p,
00096       ParamVec* deriv_p,
00097       Epetra_Vector* g,
00098       const EpetraExt::ModelEvaluator::Derivative& dg_dx,
00099       const EpetraExt::ModelEvaluator::Derivative& dg_dxdot,
00100       const EpetraExt::ModelEvaluator::Derivative& dg_dxdotdot,
00101       const EpetraExt::ModelEvaluator::Derivative& dg_dp);
00102 
00103 #ifdef ALBANY_SG_MP
00104 
00105     virtual void evaluateSGDerivative(
00106       const double current_time,
00107       const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00108       const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00109       const Stokhos::EpetraVectorOrthogPoly& sg_x,
00110       const Teuchos::Array<ParamVec>& p,
00111       const Teuchos::Array<int>& sg_p_index,
00112       const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00113       ParamVec* deriv_p,
00114       Stokhos::EpetraVectorOrthogPoly* sg_g,
00115       const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dx,
00116       const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dxdot,
00117       const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dxdotdot,
00118       const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dp);
00119 
00121     virtual void evaluateMPDerivative(
00122       const double current_time,
00123       const Stokhos::ProductEpetraVector* mp_xdot,
00124       const Stokhos::ProductEpetraVector* mp_xdotdot,
00125       const Stokhos::ProductEpetraVector& mp_x,
00126       const Teuchos::Array<ParamVec>& p,
00127       const Teuchos::Array<int>& mp_p_index,
00128       const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00129       ParamVec* deriv_p,
00130       Stokhos::ProductEpetraVector* mp_g,
00131       const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dx,
00132       const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dxdot,
00133       const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dxdotdot,
00134       const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dp);
00135 #endif //ALBANY_SG_MP
00136 
00138 
00139 
00140   private:
00141 
00143     DistributedResponseFunction(const DistributedResponseFunction&);
00144     
00146     DistributedResponseFunction& operator=(const DistributedResponseFunction&);
00147 
00148   protected:
00149 
00151     Teuchos::RCP<const Epetra_Comm> comm;
00152 
00153   };
00154 
00155 }
00156 
00157 #endif // ALBANY_SCALAR_RESPONSE_FUNCTION_HPP

Generated on Wed Mar 26 2014 18:36:36 for Albany: a Trilinos-based PDE code by  doxygen 1.7.1