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

Albany_SolutionAverageResponseFunction.cpp

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 #include "Albany_SolutionAverageResponseFunction.hpp"
00008 
00009 Albany::SolutionAverageResponseFunction::
00010 SolutionAverageResponseFunction(const Teuchos::RCP<const Epetra_Comm>& comm) :
00011   ScalarResponseFunction(comm)
00012 {
00013 }
00014 
00015 Albany::SolutionAverageResponseFunction::
00016 ~SolutionAverageResponseFunction()
00017 {
00018 }
00019 
00020 unsigned int
00021 Albany::SolutionAverageResponseFunction::
00022 numResponses() const 
00023 {
00024   return 1;
00025 }
00026 
00027 void
00028 Albany::SolutionAverageResponseFunction::
00029 evaluateResponse(const double current_time,
00030      const Epetra_Vector* xdot,
00031      const Epetra_Vector* xdotdot,
00032      const Epetra_Vector& x,
00033      const Teuchos::Array<ParamVec>& p,
00034      Epetra_Vector& g)
00035 {
00036   x.MeanValue(&g[0]);
00037 }
00038 
00039 void
00040 Albany::SolutionAverageResponseFunction::
00041 evaluateTangent(const double alpha, 
00042     const double beta,
00043     const double omega,
00044     const double current_time,
00045     bool sum_derivs,
00046     const Epetra_Vector* xdot,
00047     const Epetra_Vector* xdotdot,
00048     const Epetra_Vector& x,
00049     const Teuchos::Array<ParamVec>& p,
00050     ParamVec* deriv_p,
00051     const Epetra_MultiVector* Vxdot,
00052     const Epetra_MultiVector* Vxdotdot,
00053     const Epetra_MultiVector* Vx,
00054     const Epetra_MultiVector* Vp,
00055     Epetra_Vector* g,
00056     Epetra_MultiVector* gx,
00057     Epetra_MultiVector* gp)
00058 {
00059   // Evaluate response g
00060   if (g != NULL)
00061     x.MeanValue(&(*g)[0]);
00062 
00063   // Evaluate tangent of g = dg/dx*Vx + dg/dxdot*Vxdot + dg/dp*Vp
00064   // If Vx == NULL, Vx is the identity
00065   if (gx != NULL) {
00066     if (Vx != NULL)
00067       for (int j=0; j<Vx->NumVectors(); j++)
00068   (*Vx)(j)->MeanValue(&(*gx)[j][0]);
00069     else
00070       gx->PutScalar(1.0/x.GlobalLength());
00071     gx->Scale(alpha);
00072   }
00073   if (gp != NULL)
00074     gp->PutScalar(0.0);
00075 }
00076 
00077 void
00078 Albany::SolutionAverageResponseFunction::
00079 evaluateGradient(const double current_time,
00080      const Epetra_Vector* xdot,
00081      const Epetra_Vector* xdotdot,
00082      const Epetra_Vector& x,
00083      const Teuchos::Array<ParamVec>& p,
00084      ParamVec* deriv_p,
00085      Epetra_Vector* g,
00086      Epetra_MultiVector* dg_dx,
00087      Epetra_MultiVector* dg_dxdot,
00088      Epetra_MultiVector* dg_dxdotdot,
00089      Epetra_MultiVector* dg_dp)
00090 {
00091 
00092   // Evaluate response g
00093   if (g != NULL)
00094     x.MeanValue(&(*g)[0]);
00095 
00096   // Evaluate dg/dx
00097   if (dg_dx != NULL)
00098     dg_dx->PutScalar(1.0 / x.GlobalLength());
00099 
00100   // Evaluate dg/dxdot
00101   if (dg_dxdot != NULL)
00102     dg_dxdot->PutScalar(0.0);
00103   if (dg_dxdotdot != NULL)
00104     dg_dxdotdot->PutScalar(0.0);
00105 
00106   // Evaluate dg/dp
00107   if (dg_dp != NULL)
00108     dg_dp->PutScalar(0.0);
00109 }
00110 
00111 #ifdef ALBANY_SG_MP
00112 void
00113 Albany::SolutionAverageResponseFunction::
00114 evaluateSGResponse(
00115   const double current_time,
00116   const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00117   const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00118   const Stokhos::EpetraVectorOrthogPoly& sg_x,
00119   const Teuchos::Array<ParamVec>& p,
00120   const Teuchos::Array<int>& sg_p_index,
00121   const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00122   Stokhos::EpetraVectorOrthogPoly& sg_g)
00123 {
00124   for (int i=0; i<sg_x.size(); i++)
00125     sg_x[i].MeanValue(&sg_g[i][0]);
00126 }
00127 
00128 void
00129 Albany::SolutionAverageResponseFunction::
00130 evaluateSGTangent(
00131   const double alpha, 
00132   const double beta, 
00133   const double omega, 
00134   const double current_time,
00135   bool sum_derivs,
00136   const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00137   const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00138   const Stokhos::EpetraVectorOrthogPoly& sg_x,
00139   const Teuchos::Array<ParamVec>& p,
00140   const Teuchos::Array<int>& sg_p_index,
00141   const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00142   ParamVec* deriv_p,
00143   const Epetra_MultiVector* Vx,
00144   const Epetra_MultiVector* Vxdot,
00145   const Epetra_MultiVector* Vxdotdot,
00146   const Epetra_MultiVector* Vp,
00147   Stokhos::EpetraVectorOrthogPoly* sg_g,
00148   Stokhos::EpetraMultiVectorOrthogPoly* sg_JV,
00149   Stokhos::EpetraMultiVectorOrthogPoly* sg_gp)
00150 {
00151   // Evaluate response g
00152   if (sg_g != NULL)
00153     for (int i=0; i<sg_x.size(); i++)
00154       sg_x[i].MeanValue(&(*sg_g)[i][0]);
00155 
00156   // Evaluate tangent of g = dg/dx*Vx + dg/dxdot*Vxdot + dg/dp*Vp
00157   // If Vx == NULL, Vx is the identity
00158   if (sg_JV != NULL) {
00159     sg_JV->init(0.0);
00160     if (Vx != NULL)
00161       for (int j=0; j<Vx->NumVectors(); j++)
00162   (*Vx)(j)->MeanValue(&(*sg_JV)[0][j][0]);
00163     else
00164       (*sg_JV)[0].PutScalar(alpha/sg_x[0].GlobalLength());
00165   }
00166   if (sg_gp != NULL)
00167     sg_gp->init(0.0);
00168 }
00169 
00170 void
00171 Albany::SolutionAverageResponseFunction::
00172 evaluateSGGradient(
00173   const double current_time,
00174   const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00175   const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00176   const Stokhos::EpetraVectorOrthogPoly& sg_x,
00177   const Teuchos::Array<ParamVec>& p,
00178   const Teuchos::Array<int>& sg_p_index,
00179   const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00180   ParamVec* deriv_p,
00181   Stokhos::EpetraVectorOrthogPoly* sg_g,
00182   Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dx,
00183   Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dxdot,
00184   Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dxdotdot,
00185   Stokhos::EpetraMultiVectorOrthogPoly* sg_dg_dp)
00186 {
00187   // Evaluate response g
00188   if (sg_g != NULL)
00189     for (int i=0; i<sg_x.size(); i++)
00190       sg_x[i].MeanValue(&(*sg_g)[i][0]);
00191 
00192   // Evaluate dg/dx
00193   if (sg_dg_dx != NULL)
00194     (*sg_dg_dx)[0].PutScalar(1.0 / sg_x[0].GlobalLength());
00195 
00196   // Evaluate dg/dxdot
00197   if (sg_dg_dxdot != NULL)
00198     sg_dg_dxdot->init(0.0);
00199   if (sg_dg_dxdotdot != NULL)
00200     sg_dg_dxdotdot->init(0.0);
00201 
00202   // Evaluate dg/dp
00203   if (sg_dg_dp != NULL)
00204     sg_dg_dp->init(0.0);
00205 }
00206 
00207 void
00208 Albany::SolutionAverageResponseFunction::
00209 evaluateMPResponse(
00210   const double current_time,
00211   const Stokhos::ProductEpetraVector* mp_xdot,
00212   const Stokhos::ProductEpetraVector* mp_xdotdot,
00213   const Stokhos::ProductEpetraVector& mp_x,
00214   const Teuchos::Array<ParamVec>& p,
00215   const Teuchos::Array<int>& mp_p_index,
00216   const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00217   Stokhos::ProductEpetraVector& mp_g)
00218 {
00219   for (int i=0; i<mp_x.size(); i++)
00220     mp_x[i].MeanValue(&mp_g[i][0]);
00221 }
00222 
00223 void
00224 Albany::SolutionAverageResponseFunction::
00225 evaluateMPTangent(
00226   const double alpha, 
00227   const double beta, 
00228   const double omega, 
00229   const double current_time,
00230   bool sum_derivs,
00231   const Stokhos::ProductEpetraVector* mp_xdot,
00232   const Stokhos::ProductEpetraVector* mp_xdotdot,
00233   const Stokhos::ProductEpetraVector& mp_x,
00234   const Teuchos::Array<ParamVec>& p,
00235   const Teuchos::Array<int>& mp_p_index,
00236   const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00237   ParamVec* deriv_p,
00238   const Epetra_MultiVector* Vx,
00239   const Epetra_MultiVector* Vxdot,
00240   const Epetra_MultiVector* Vxdotdot,
00241   const Epetra_MultiVector* Vp,
00242   Stokhos::ProductEpetraVector* mp_g,
00243   Stokhos::ProductEpetraMultiVector* mp_JV,
00244   Stokhos::ProductEpetraMultiVector* mp_gp)
00245 {
00246   // Evaluate response g
00247   if (mp_g != NULL)
00248     for (int i=0; i<mp_x.size(); i++)
00249       mp_x[i].MeanValue(&(*mp_g)[i][0]);
00250 
00251   // Evaluate tangent of g = dg/dx*Vx + dg/dxdot*Vxdot + dg/dp*Vp
00252   // If Vx == NULL, Vx is the identity
00253   if (mp_JV != NULL) {
00254     if (Vx != NULL)
00255       for (int i=0; i<mp_x.size(); i++)
00256   for (int j=0; j<Vx->NumVectors(); j++)
00257     (*Vx)(j)->MeanValue(&(*mp_JV)[i][j][0]);
00258     else
00259       for (int i=0; i<mp_x.size(); i++)
00260   (*mp_JV)[i].PutScalar(alpha/mp_x[0].GlobalLength());
00261   }
00262   if (mp_gp != NULL)
00263     mp_gp->init(0.0);
00264 }
00265 
00266 void
00267 Albany::SolutionAverageResponseFunction::
00268 evaluateMPGradient(
00269   const double current_time,
00270   const Stokhos::ProductEpetraVector* mp_xdot,
00271   const Stokhos::ProductEpetraVector* mp_xdotdot,
00272   const Stokhos::ProductEpetraVector& mp_x,
00273   const Teuchos::Array<ParamVec>& p,
00274   const Teuchos::Array<int>& mp_p_index,
00275   const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00276   ParamVec* deriv_p,
00277   Stokhos::ProductEpetraVector* mp_g,
00278   Stokhos::ProductEpetraMultiVector* mp_dg_dx,
00279   Stokhos::ProductEpetraMultiVector* mp_dg_dxdot,
00280   Stokhos::ProductEpetraMultiVector* mp_dg_dxdotdot,
00281   Stokhos::ProductEpetraMultiVector* mp_dg_dp)
00282 {
00283   // Evaluate response g
00284   if (mp_g != NULL)
00285     for (int i=0; i<mp_x.size(); i++)
00286       mp_x[i].MeanValue(&(*mp_g)[i][0]);
00287 
00288   // Evaluate dg/dx
00289   if (mp_dg_dx != NULL)
00290     for (int i=0; i<mp_x.size(); i++)
00291       (*mp_dg_dx)[i].PutScalar(1.0 / mp_x[0].GlobalLength());
00292 
00293   // Evaluate dg/dxdot
00294   if (mp_dg_dxdot != NULL)
00295     mp_dg_dxdot->init(0.0);
00296   if (mp_dg_dxdotdot != NULL)
00297     mp_dg_dxdotdot->init(0.0);
00298 
00299   // Evaluate dg/dp
00300   if (mp_dg_dp != NULL)
00301     mp_dg_dp->init(0.0);
00302 }
00303 #endif //ALBANY_SG_MP

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