00001
00002
00003
00004
00005
00006
00007 #ifndef ALBANY_SAMPLING_BASED_SCALAR_RESPONSE_FUNCTION_HPP
00008 #define ALBANY_SAMPLING_BASED_SCALAR_RESPONSE_FUNCTION_HPP
00009
00010 #include "Albany_ScalarResponseFunction.hpp"
00011 #include "Stokhos_Quadrature.hpp"
00012
00013 namespace Albany {
00014
00019 class SamplingBasedScalarResponseFunction :
00020 public ScalarResponseFunction {
00021 public:
00022
00024 SamplingBasedScalarResponseFunction(
00025 const Teuchos::RCP<const Epetra_Comm>& comm) :
00026 ScalarResponseFunction(comm) {};
00027
00029 virtual ~SamplingBasedScalarResponseFunction() {};
00030
00032
00033
00034 #ifdef ALBANY_SG_MP
00035
00036 virtual void init_sg(
00037 const Teuchos::RCP<const Stokhos::OrthogPolyBasis<int,double> >& basis,
00038 const Teuchos::RCP<const Stokhos::Quadrature<int,double> >& quad,
00039 const Teuchos::RCP<Stokhos::OrthogPolyExpansion<int,double> >& expansion,
00040 const Teuchos::RCP<const EpetraExt::MultiComm>& multiComm);
00041
00043 virtual void evaluateSGResponse(
00044 const double curr_time,
00045 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00046 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00047 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00048 const Teuchos::Array<ParamVec>& p,
00049 const Teuchos::Array<int>& sg_p_index,
00050 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00051 Stokhos::EpetraVectorOrthogPoly& sg_g);
00052
00054 virtual void evaluateSGTangent(
00055 const double alpha,
00056 const double beta,
00057 const double omega,
00058 const double current_time,
00059 bool sum_derivs,
00060 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00061 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00062 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00063 const Teuchos::Array<ParamVec>& p,
00064 const Teuchos::Array<int>& sg_p_index,
00065 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00066 ParamVec* deriv_p,
00067 const Epetra_MultiVector* Vx,
00068 const Epetra_MultiVector* Vxdot,
00069 const Epetra_MultiVector* Vxdotdot,
00070 const Epetra_MultiVector* Vp,
00071 Stokhos::EpetraVectorOrthogPoly* sg_g,
00072 Stokhos::EpetraMultiVectorOrthogPoly* sg_JV,
00073 Stokhos::EpetraMultiVectorOrthogPoly* sg_gp);
00074
00076 virtual void evaluateSGGradient(
00077 const double current_time,
00078 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00079 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00080 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00081 const Teuchos::Array<ParamVec>& p,
00082 const Teuchos::Array<int>& sg_p_index,
00083 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00084 ParamVec* deriv_p,
00085 Stokhos::EpetraVectorOrthogPoly* sg_g,
00086 Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dx,
00087 Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dxdot,
00088 Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dxdotdot,
00089 Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dp);
00090 #endif //ALBANY_SG_MP
00091
00093
00095
00096
00097 #ifdef ALBANY_SG_MP
00098
00099 virtual void evaluateMPResponse(
00100 const double curr_time,
00101 const Stokhos::ProductEpetraVector* mp_xdot,
00102 const Stokhos::ProductEpetraVector* mp_xdotdot,
00103 const Stokhos::ProductEpetraVector& mp_x,
00104 const Teuchos::Array<ParamVec>& p,
00105 const Teuchos::Array<int>& mp_p_index,
00106 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00107 Stokhos::ProductEpetraVector& mp_g);
00108
00110 virtual void evaluateMPTangent(
00111 const double alpha,
00112 const double beta,
00113 const double omega,
00114 const double current_time,
00115 bool sum_derivs,
00116 const Stokhos::ProductEpetraVector* mp_xdot,
00117 const Stokhos::ProductEpetraVector* mp_xdotdot,
00118 const Stokhos::ProductEpetraVector& mp_x,
00119 const Teuchos::Array<ParamVec>& p,
00120 const Teuchos::Array<int>& mp_p_index,
00121 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00122 ParamVec* deriv_p,
00123 const Epetra_MultiVector* Vx,
00124 const Epetra_MultiVector* Vxdot,
00125 const Epetra_MultiVector* Vxdotdot,
00126 const Epetra_MultiVector* Vp,
00127 Stokhos::ProductEpetraVector* mp_g,
00128 Stokhos::ProductEpetraMultiVector* mp_JV,
00129 Stokhos::ProductEpetraMultiVector* mp_gp);
00130
00132 virtual void evaluateMPGradient(
00133 const double current_time,
00134 const Stokhos::ProductEpetraVector* mp_xdot,
00135 const Stokhos::ProductEpetraVector* mp_xdotdot,
00136 const Stokhos::ProductEpetraVector& mp_x,
00137 const Teuchos::Array<ParamVec>& p,
00138 const Teuchos::Array<int>& mp_p_index,
00139 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00140 ParamVec* deriv_p,
00141 Stokhos::ProductEpetraVector* mp_g,
00142 Stokhos::ProductEpetraMultiVector* mp_dg_dx,
00143 Stokhos::ProductEpetraMultiVector* mp_dg_dxdot,
00144 Stokhos::ProductEpetraMultiVector* mp_dg_dxdotdot,
00145 Stokhos::ProductEpetraMultiVector* mp_dg_dp);
00146 #endif //ALBANY_SG_MP
00147
00149
00150 private:
00151
00153 SamplingBasedScalarResponseFunction(
00154 const SamplingBasedScalarResponseFunction&);
00155
00157 SamplingBasedScalarResponseFunction& operator=(
00158 const SamplingBasedScalarResponseFunction&);
00159
00160 protected:
00161
00162 #ifdef ALBANY_SG_MP
00163
00164 Teuchos::RCP<const Stokhos::OrthogPolyBasis<int,double> > basis;
00165
00167 Teuchos::RCP<const Stokhos::Quadrature<int,double> > quad;
00168 #endif //ALBANY_SG_MP
00169
00170 };
00171
00172 }
00173
00174 #endif // ALBANY_SAMPLING_BASED_SCALAR_RESPONSE_FUNCTION_HPP