00001
00002
00003
00004
00005
00006
00007 #ifndef ALBANY_SCALAR_RESPONSE_FUNCTION_HPP
00008 #define ALBANY_SCALAR_RESPONSE_FUNCTION_HPP
00009
00010 #include "Albany_AbstractResponseFunction.hpp"
00011
00012 namespace Albany {
00013
00022 class ScalarResponseFunction :
00023 public AbstractResponseFunction {
00024 public:
00025
00027 ScalarResponseFunction(const Teuchos::RCP<const Epetra_Comm>& comm_) :
00028 comm(comm_) {};
00029
00031 virtual ~ScalarResponseFunction() {};
00032
00034 virtual unsigned int numResponses() const = 0;
00035
00037 virtual Teuchos::RCP<const Epetra_Comm> getComm() const {
00038 return comm;
00039 }
00040
00042 virtual void evaluateGradient(
00043 const double current_time,
00044 const Epetra_Vector* xdot,
00045 const Epetra_Vector* xdotdot,
00046 const Epetra_Vector& x,
00047 const Teuchos::Array<ParamVec>& p,
00048 ParamVec* deriv_p,
00049 Epetra_Vector* g,
00050 Epetra_MultiVector* dg_dx,
00051 Epetra_MultiVector* dg_dxdot,
00052 Epetra_MultiVector* dg_dxdotdot,
00053 Epetra_MultiVector* dg_dp) = 0;
00054
00055 #ifdef ALBANY_SG_MP
00056
00057 virtual void evaluateSGGradient(
00058 const double current_time,
00059 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00060 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00061 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00062 const Teuchos::Array<ParamVec>& p,
00063 const Teuchos::Array<int>& sg_p_index,
00064 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00065 ParamVec* deriv_p,
00066 Stokhos::EpetraVectorOrthogPoly* sg_g,
00067 Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dx,
00068 Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dxdot,
00069 Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dxdotdot,
00070 Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dp) = 0;
00071
00073 virtual void evaluateMPGradient(
00074 const double current_time,
00075 const Stokhos::ProductEpetraVector* mp_xdot,
00076 const Stokhos::ProductEpetraVector* mp_xdotdot,
00077 const Stokhos::ProductEpetraVector& mp_x,
00078 const Teuchos::Array<ParamVec>& p,
00079 const Teuchos::Array<int>& mp_p_index,
00080 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00081 ParamVec* deriv_p,
00082 Stokhos::ProductEpetraVector* mp_g,
00083 Stokhos::ProductEpetraMultiVector* mp_dg_dx,
00084 Stokhos::ProductEpetraMultiVector* mp_dg_dxdot,
00085 Stokhos::ProductEpetraMultiVector* mp_dg_dxdotdot,
00086 Stokhos::ProductEpetraMultiVector* mp_dg_dp) = 0;
00087 #endif //ALBANY_SG_MP
00088
00090
00091
00093 virtual void setup() {}
00094
00099 virtual bool isScalarResponse() const { return true; }
00100
00102
00106 virtual Teuchos::RCP<Epetra_Operator> createGradientOp() const;
00107
00109 virtual Teuchos::RCP<const Epetra_Map> responseMap() const;
00110
00112 virtual void evaluateDerivative(
00113 const double current_time,
00114 const Epetra_Vector* xdot,
00115 const Epetra_Vector* xdotdot,
00116 const Epetra_Vector& x,
00117 const Teuchos::Array<ParamVec>& p,
00118 ParamVec* deriv_p,
00119 Epetra_Vector* g,
00120 const EpetraExt::ModelEvaluator::Derivative& dg_dx,
00121 const EpetraExt::ModelEvaluator::Derivative& dg_dxdot,
00122 const EpetraExt::ModelEvaluator::Derivative& dg_dxdotdot,
00123 const EpetraExt::ModelEvaluator::Derivative& dg_dp);
00124
00125 #ifdef ALBANY_SG_MP
00126
00127 virtual void evaluateSGDerivative(
00128 const double current_time,
00129 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00130 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00131 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00132 const Teuchos::Array<ParamVec>& p,
00133 const Teuchos::Array<int>& sg_p_index,
00134 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00135 ParamVec* deriv_p,
00136 Stokhos::EpetraVectorOrthogPoly* sg_g,
00137 const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dx,
00138 const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dxdot,
00139 const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dxdotdot,
00140 const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dp);
00141
00143 virtual void evaluateMPDerivative(
00144 const double current_time,
00145 const Stokhos::ProductEpetraVector* mp_xdot,
00146 const Stokhos::ProductEpetraVector* mp_xdotdot,
00147 const Stokhos::ProductEpetraVector& mp_x,
00148 const Teuchos::Array<ParamVec>& p,
00149 const Teuchos::Array<int>& mp_p_index,
00150 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00151 ParamVec* deriv_p,
00152 Stokhos::ProductEpetraVector* mp_g,
00153 const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dx,
00154 const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dxdot,
00155 const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dxdotdot,
00156 const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dp);
00157 #endif //ALBANY_SG_MP
00158
00160
00161
00162 private:
00163
00165 ScalarResponseFunction(const ScalarResponseFunction&);
00166
00168 ScalarResponseFunction& operator=(const ScalarResponseFunction&);
00169
00170 protected:
00171
00173 Teuchos::RCP<const Epetra_Comm> comm;
00174
00175 };
00176
00177 }
00178
00179 #endif // ALBANY_SCALAR_RESPONSE_FUNCTION_HPP