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

Albany_AbstractResponseFunction.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_ABSTRACTRESPONSEFUNCTION_HPP
00008 #define ALBANY_ABSTRACTRESPONSEFUNCTION_HPP
00009 
00010 #include "Teuchos_Array.hpp"
00011 #include "Teuchos_RCP.hpp"
00012 #include "Epetra_Map.h"
00013 #include "Epetra_Vector.h"
00014 #include "Epetra_MultiVector.h"
00015 #include "EpetraExt_ModelEvaluator.h"
00016 #include "Stokhos_OrthogPolyBasis.hpp"
00017 #include "Stokhos_OrthogPolyExpansion.hpp"
00018 #include "Stokhos_Quadrature.hpp"
00019 #include "EpetraExt_MultiComm.h"
00020 #include "Stokhos_EpetraVectorOrthogPoly.hpp"
00021 #include "Stokhos_EpetraMultiVectorOrthogPoly.hpp"
00022 #include "Stokhos_ProductEpetraVector.hpp"
00023 #include "Stokhos_ProductEpetraMultiVector.hpp"
00024 #include "PHAL_AlbanyTraits.hpp"
00025 
00026 namespace Albany {
00027 
00031   class AbstractResponseFunction {
00032   public:
00033   
00035     AbstractResponseFunction() {};
00036 
00038     virtual ~AbstractResponseFunction() {};
00039 
00041     virtual void setup() = 0;
00042 
00044     virtual Teuchos::RCP<const Epetra_Map> responseMap() const = 0;
00045 
00050     virtual bool isScalarResponse() const = 0;
00051 
00053     virtual Teuchos::RCP<Epetra_Operator> createGradientOp() const = 0;
00054 
00056 
00057 
00059     virtual void evaluateResponse(
00060       const double current_time,
00061       const Epetra_Vector* xdot,
00062       const Epetra_Vector* xdotdot,
00063       const Epetra_Vector& x,
00064       const Teuchos::Array<ParamVec>& p,
00065       Epetra_Vector& g) = 0;
00066 
00068     virtual void evaluateTangent(
00069       const double alpha, 
00070       const double beta,
00071       const double omega,
00072       const double current_time,
00073       bool sum_derivs,
00074       const Epetra_Vector* xdot,
00075       const Epetra_Vector* xdotdot,
00076       const Epetra_Vector& x,
00077       const Teuchos::Array<ParamVec>& p,
00078       ParamVec* deriv_p,
00079       const Epetra_MultiVector* Vxdot,
00080       const Epetra_MultiVector* Vxdotdot,
00081       const Epetra_MultiVector* Vx,
00082       const Epetra_MultiVector* Vp,
00083       Epetra_Vector* g,
00084       Epetra_MultiVector* gx,
00085       Epetra_MultiVector* gp) = 0;
00086 
00088     virtual void evaluateDerivative(
00089       const double current_time,
00090       const Epetra_Vector* xdot,
00091       const Epetra_Vector* xdotdot,
00092       const Epetra_Vector& x,
00093       const Teuchos::Array<ParamVec>& p,
00094       ParamVec* deriv_p,
00095       Epetra_Vector* g,
00096       const EpetraExt::ModelEvaluator::Derivative& dg_dx,
00097       const EpetraExt::ModelEvaluator::Derivative& dg_dxdot,
00098       const EpetraExt::ModelEvaluator::Derivative& dg_dxdotdot,
00099       const EpetraExt::ModelEvaluator::Derivative& dg_dp) = 0;
00100 
00102 
00104 
00105 
00106 #ifdef ALBANY_SG_MP
00107 
00108     virtual void init_sg(
00109       const Teuchos::RCP<const Stokhos::OrthogPolyBasis<int,double> >& basis,
00110       const Teuchos::RCP<const Stokhos::Quadrature<int,double> >& quad,
00111       const Teuchos::RCP<Stokhos::OrthogPolyExpansion<int,double> >& expansion,
00112       const Teuchos::RCP<const EpetraExt::MultiComm>& multiComm) = 0;
00113 
00115     virtual void evaluateSGResponse(
00116       const double curr_time,
00117       const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00118       const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00119       const Stokhos::EpetraVectorOrthogPoly& sg_x,
00120       const Teuchos::Array<ParamVec>& p,
00121       const Teuchos::Array<int>& sg_p_index,
00122       const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00123       Stokhos::EpetraVectorOrthogPoly& sg_g) = 0;
00124 
00126     virtual void evaluateSGTangent(
00127       const double alpha, 
00128       const double beta, 
00129       const double omega, 
00130       const double current_time,
00131       bool sum_derivs,
00132       const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00133       const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00134       const Stokhos::EpetraVectorOrthogPoly& sg_x,
00135       const Teuchos::Array<ParamVec>& p,
00136       const Teuchos::Array<int>& sg_p_index,
00137       const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00138       ParamVec* deriv_p,
00139       const Epetra_MultiVector* Vx,
00140       const Epetra_MultiVector* Vxdot,
00141       const Epetra_MultiVector* Vxdotdot,
00142       const Epetra_MultiVector* Vp,
00143       Stokhos::EpetraVectorOrthogPoly* sg_g,
00144       Stokhos::EpetraMultiVectorOrthogPoly* sg_JV,
00145       Stokhos::EpetraMultiVectorOrthogPoly* sg_gp) = 0;
00146 
00148     virtual void evaluateSGDerivative(
00149       const double current_time,
00150       const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00151       const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00152       const Stokhos::EpetraVectorOrthogPoly& sg_x,
00153       const Teuchos::Array<ParamVec>& p,
00154       const Teuchos::Array<int>& sg_p_index,
00155       const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00156       ParamVec* deriv_p,
00157       Stokhos::EpetraVectorOrthogPoly* sg_g,
00158       const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dx,
00159       const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dxdot,
00160       const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dxdotdot,
00161       const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dp) = 0;
00162 #endif //ALBANY_SG_MP
00163 
00165 
00167 
00168 
00169 #ifdef ALBANY_SG_MP
00170 
00171     virtual void evaluateMPResponse(
00172       const double curr_time,
00173       const Stokhos::ProductEpetraVector* mp_xdot,
00174       const Stokhos::ProductEpetraVector* mp_xdotdot,
00175       const Stokhos::ProductEpetraVector& mp_x,
00176       const Teuchos::Array<ParamVec>& p,
00177       const Teuchos::Array<int>& mp_p_index,
00178       const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00179       Stokhos::ProductEpetraVector& mp_g) = 0;
00180 
00182     virtual void evaluateMPTangent(
00183       const double alpha, 
00184       const double beta, 
00185       const double omega, 
00186       const double current_time,
00187       bool sum_derivs,
00188       const Stokhos::ProductEpetraVector* mp_xdot,
00189       const Stokhos::ProductEpetraVector* mp_xdotdot,
00190       const Stokhos::ProductEpetraVector& mp_x,
00191       const Teuchos::Array<ParamVec>& p,
00192       const Teuchos::Array<int>& mp_p_index,
00193       const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00194       ParamVec* deriv_p,
00195       const Epetra_MultiVector* Vx,
00196       const Epetra_MultiVector* Vxdot,
00197       const Epetra_MultiVector* Vxdotdot,
00198       const Epetra_MultiVector* Vp,
00199       Stokhos::ProductEpetraVector* mp_g,
00200       Stokhos::ProductEpetraMultiVector* mp_JV,
00201       Stokhos::ProductEpetraMultiVector* mp_gp) = 0;
00202 
00204     virtual void evaluateMPDerivative(
00205       const double current_time,
00206       const Stokhos::ProductEpetraVector* mp_xdot,
00207       const Stokhos::ProductEpetraVector* mp_xdotdot,
00208       const Stokhos::ProductEpetraVector& mp_x,
00209       const Teuchos::Array<ParamVec>& p,
00210       const Teuchos::Array<int>& mp_p_index,
00211       const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00212       ParamVec* deriv_p,
00213       Stokhos::ProductEpetraVector* mp_g,
00214       const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dx,
00215       const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dxdot,
00216       const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dxdotdot,
00217       const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dp) = 0;
00218 #endif //ALBANY_SG_MP
00219 
00221 
00222   private:
00223 
00225     AbstractResponseFunction(const AbstractResponseFunction&);
00226     
00228     AbstractResponseFunction& operator=(const AbstractResponseFunction&);
00229 
00230   };
00231 
00232 }
00233 
00234 #endif // ALBANY_ABSTRACTRESPONSEFUNCTION_HPP

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