00001
00002
00003
00004
00005
00006
00007 #ifndef ALBANY_FIELD_MANAGER_SCALAR_RESPONSE_FUNCTION_HPP
00008 #define ALBANY_FIELD_MANAGER_SCALAR_RESPONSE_FUNCTION_HPP
00009
00010 #include "Albany_ScalarResponseFunction.hpp"
00011 #include "Albany_Application.hpp"
00012 #include "Albany_AbstractProblem.hpp"
00013 #include "Albany_StateManager.hpp"
00014 #include "Albany_StateInfoStruct.hpp"
00015 #include "PHAL_AlbanyTraits.hpp"
00016 #include "Phalanx.hpp"
00017
00018 namespace Albany {
00019
00023 class FieldManagerScalarResponseFunction :
00024 public ScalarResponseFunction {
00025 public:
00026
00028 FieldManagerScalarResponseFunction(
00029 const Teuchos::RCP<Albany::Application>& application,
00030 const Teuchos::RCP<Albany::AbstractProblem>& problem,
00031 const Teuchos::RCP<Albany::MeshSpecsStruct>& ms,
00032 const Teuchos::RCP<Albany::StateManager>& stateMgr,
00033 Teuchos::ParameterList& responseParams);
00034
00036 virtual ~FieldManagerScalarResponseFunction();
00037
00039 virtual unsigned int numResponses() const;
00040
00042 virtual void
00043 evaluateResponse(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 Epetra_Vector& g);
00049
00051 virtual void
00052 evaluateTangent(const double alpha,
00053 const double beta,
00054 const double omega,
00055 const double current_time,
00056 bool sum_derivs,
00057 const Epetra_Vector* xdot,
00058 const Epetra_Vector* xdotdot,
00059 const Epetra_Vector& x,
00060 const Teuchos::Array<ParamVec>& p,
00061 ParamVec* deriv_p,
00062 const Epetra_MultiVector* Vxdot,
00063 const Epetra_MultiVector* Vxdotdot,
00064 const Epetra_MultiVector* Vx,
00065 const Epetra_MultiVector* Vp,
00066 Epetra_Vector* g,
00067 Epetra_MultiVector* gx,
00068 Epetra_MultiVector* gp);
00069
00071 virtual void
00072 evaluateGradient(const double current_time,
00073 const Epetra_Vector* xdot,
00074 const Epetra_Vector* xdotdot,
00075 const Epetra_Vector& x,
00076 const Teuchos::Array<ParamVec>& p,
00077 ParamVec* deriv_p,
00078 Epetra_Vector* g,
00079 Epetra_MultiVector* dg_dx,
00080 Epetra_MultiVector* dg_dxdot,
00081 Epetra_MultiVector* dg_dxdotdot,
00082 Epetra_MultiVector* dg_dp);
00083
00085
00086
00087 #ifdef ALBANY_SG_MP
00088
00089 virtual void init_sg(
00090 const Teuchos::RCP<const Stokhos::OrthogPolyBasis<int,double> >& basis,
00091 const Teuchos::RCP<const Stokhos::Quadrature<int,double> >& quad,
00092 const Teuchos::RCP<Stokhos::OrthogPolyExpansion<int,double> >& expansion,
00093 const Teuchos::RCP<const EpetraExt::MultiComm>& multiComm) {}
00094
00096 virtual void evaluateSGResponse(
00097 const double curr_time,
00098 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00099 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00100 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00101 const Teuchos::Array<ParamVec>& p,
00102 const Teuchos::Array<int>& sg_p_index,
00103 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00104 Stokhos::EpetraVectorOrthogPoly& sg_g);
00105
00107 virtual void evaluateSGTangent(
00108 const double alpha,
00109 const double beta,
00110 const double omega,
00111 const double current_time,
00112 bool sum_derivs,
00113 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00114 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00115 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00116 const Teuchos::Array<ParamVec>& p,
00117 const Teuchos::Array<int>& sg_p_index,
00118 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00119 ParamVec* deriv_p,
00120 const Epetra_MultiVector* Vx,
00121 const Epetra_MultiVector* Vxdot,
00122 const Epetra_MultiVector* Vxdotdot,
00123 const Epetra_MultiVector* Vp,
00124 Stokhos::EpetraVectorOrthogPoly* sg_g,
00125 Stokhos::EpetraMultiVectorOrthogPoly* sg_JV,
00126 Stokhos::EpetraMultiVectorOrthogPoly* sg_gp);
00127
00129 virtual void evaluateSGGradient(
00130 const double current_time,
00131 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00132 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00133 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00134 const Teuchos::Array<ParamVec>& p,
00135 const Teuchos::Array<int>& sg_p_index,
00136 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00137 ParamVec* deriv_p,
00138 Stokhos::EpetraVectorOrthogPoly* sg_g,
00139 Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dx,
00140 Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dxdot,
00141 Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dxdotdot,
00142 Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dp);
00143 #endif //ALBANY_SG_MP
00144
00146
00148
00149
00150 #ifdef ALBANY_SG_MP
00151
00152 virtual void evaluateMPResponse(
00153 const double curr_time,
00154 const Stokhos::ProductEpetraVector* mp_xdot,
00155 const Stokhos::ProductEpetraVector* mp_xdotdot,
00156 const Stokhos::ProductEpetraVector& mp_x,
00157 const Teuchos::Array<ParamVec>& p,
00158 const Teuchos::Array<int>& mp_p_index,
00159 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00160 Stokhos::ProductEpetraVector& mp_g);
00161
00163 virtual void evaluateMPTangent(
00164 const double alpha,
00165 const double beta,
00166 const double omega,
00167 const double current_time,
00168 bool sum_derivs,
00169 const Stokhos::ProductEpetraVector* mp_xdot,
00170 const Stokhos::ProductEpetraVector* mp_xdotdot,
00171 const Stokhos::ProductEpetraVector& mp_x,
00172 const Teuchos::Array<ParamVec>& p,
00173 const Teuchos::Array<int>& mp_p_index,
00174 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00175 ParamVec* deriv_p,
00176 const Epetra_MultiVector* Vx,
00177 const Epetra_MultiVector* Vxdot,
00178 const Epetra_MultiVector* Vxdotdot,
00179 const Epetra_MultiVector* Vp,
00180 Stokhos::ProductEpetraVector* mp_g,
00181 Stokhos::ProductEpetraMultiVector* mp_JV,
00182 Stokhos::ProductEpetraMultiVector* mp_gp);
00183
00185 virtual void evaluateMPGradient(
00186 const double current_time,
00187 const Stokhos::ProductEpetraVector* mp_xdot,
00188 const Stokhos::ProductEpetraVector* mp_xdotdot,
00189 const Stokhos::ProductEpetraVector& mp_x,
00190 const Teuchos::Array<ParamVec>& p,
00191 const Teuchos::Array<int>& mp_p_index,
00192 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00193 ParamVec* deriv_p,
00194 Stokhos::ProductEpetraVector* mp_g,
00195 Stokhos::ProductEpetraMultiVector* mp_dg_dx,
00196 Stokhos::ProductEpetraMultiVector* mp_dg_dxdot,
00197 Stokhos::ProductEpetraMultiVector* mp_dg_dxdotdot,
00198 Stokhos::ProductEpetraMultiVector* mp_dg_dp);
00199 #endif //ALBANY_SG_MP
00200
00202
00203 private:
00204
00206 FieldManagerScalarResponseFunction(const FieldManagerScalarResponseFunction&);
00207
00209 FieldManagerScalarResponseFunction& operator=(const FieldManagerScalarResponseFunction&);
00210
00211 protected:
00212
00214
00217 FieldManagerScalarResponseFunction(
00218 const Teuchos::RCP<Albany::Application>& application,
00219 const Teuchos::RCP<Albany::AbstractProblem>& problem,
00220 const Teuchos::RCP<Albany::MeshSpecsStruct>& ms,
00221 const Teuchos::RCP<Albany::StateManager>& stateMgr);
00222
00224 void setup(Teuchos::ParameterList& responseParams);
00225
00227 template <typename EvalT>
00228 void visResponseGraph(const std::string& res_type);
00229
00230 protected:
00231
00233 Teuchos::RCP<Albany::Application> application;
00234
00236 Teuchos::RCP<Albany::AbstractProblem> problem;
00237
00239 Teuchos::RCP<Albany::MeshSpecsStruct> meshSpecs;
00240
00242 Teuchos::RCP<Albany::StateManager> stateMgr;
00243
00245 Teuchos::RCP<PHX::FieldManager<PHAL::AlbanyTraits> > rfm;
00246
00248 unsigned int num_responses;
00249
00251 int vis_response_graph;
00252
00254 std::string vis_response_name;
00255
00256 };
00257
00258 template <typename EvalT>
00259 void
00260 Albany::FieldManagerScalarResponseFunction::
00261 visResponseGraph(const std::string& res_type) {
00262
00263 static bool first = true;
00264 if (first && vis_response_graph > 0) {
00265 bool detail = false; if (vis_response_graph > 1) detail=true;
00266 Teuchos::RCP<Teuchos::FancyOStream> out =
00267 Teuchos::VerboseObjectBase::getDefaultOStream();
00268 *out << "Phalanx writing graphviz file for graph of Response fill "
00269 << "(detail = "<< vis_response_graph << ")"<< std::endl;
00270 std::string detail_name =
00271 "responses_graph_" + vis_response_name + res_type;
00272 *out << "Process using 'dot -Tpng -O ' " << detail_name << "\n" << std::endl;
00273 rfm->writeGraphvizFile<EvalT>(detail_name,detail,detail);
00274 first = false;
00275 }
00276 }
00277
00278 }
00279
00280 #endif // ALBANY_FIELD_MANAGER_SCALAR_RESPONSE_FUNCTION_HPP