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

Albany_SolutionTwoNormResponseFunction.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_SolutionTwoNormResponseFunction.hpp"
00008 #include "Epetra_Comm.h"
00009 
00010 Albany::SolutionTwoNormResponseFunction::
00011 SolutionTwoNormResponseFunction(const Teuchos::RCP<const Epetra_Comm>& comm) :
00012   SamplingBasedScalarResponseFunction(comm)
00013 {
00014 }
00015 
00016 Albany::SolutionTwoNormResponseFunction::
00017 ~SolutionTwoNormResponseFunction()
00018 {
00019 }
00020 
00021 unsigned int
00022 Albany::SolutionTwoNormResponseFunction::
00023 numResponses() const 
00024 {
00025   return 1;
00026 }
00027 
00028 void
00029 Albany::SolutionTwoNormResponseFunction::
00030 evaluateResponse(const double current_time,
00031      const Epetra_Vector* xdot,
00032      const Epetra_Vector* xdotdot,
00033      const Epetra_Vector& x,
00034      const Teuchos::Array<ParamVec>& p,
00035      Epetra_Vector& g)
00036 {
00037   x.Norm2(&g[0]);
00038 }
00039 
00040 void
00041 Albany::SolutionTwoNormResponseFunction::
00042 evaluateTangent(const double alpha, 
00043     const double beta,
00044     const double omega,
00045     const double current_time,
00046     bool sum_derivs,
00047     const Epetra_Vector* xdot,
00048     const Epetra_Vector* xdotdot,
00049     const Epetra_Vector& x,
00050     const Teuchos::Array<ParamVec>& p,
00051     ParamVec* deriv_p,
00052     const Epetra_MultiVector* Vxdot,
00053     const Epetra_MultiVector* Vxdotdot,
00054     const Epetra_MultiVector* Vx,
00055     const Epetra_MultiVector* Vp,
00056     Epetra_Vector* g,
00057     Epetra_MultiVector* gx,
00058     Epetra_MultiVector* gp)
00059 {
00060   double nrm;
00061   x.Norm2(&nrm);
00062 
00063   // Evaluate response g
00064   if (g != NULL)
00065     (*g)[0] = nrm;
00066 
00067   // Evaluate tangent of g = dg/dx*dx/dp + dg/dxdot*dxdot/dp + dg/dp
00068   // dg/dx = 1/||x|| * x^T
00069   if (gx != NULL) {
00070     if (Vx != NULL)
00071       gx->Multiply('T','N',alpha/nrm,x,*Vx,0.0);
00072     else
00073       gx->Update(alpha/nrm, x, 0.0);
00074   }
00075 
00076   if (gp != NULL)
00077     gp->PutScalar(0.0);
00078 }
00079 
00080 void
00081 Albany::SolutionTwoNormResponseFunction::
00082 evaluateGradient(const double current_time,
00083      const Epetra_Vector* xdot,
00084      const Epetra_Vector* xdotdot,
00085      const Epetra_Vector& x,
00086      const Teuchos::Array<ParamVec>& p,
00087      ParamVec* deriv_p,
00088      Epetra_Vector* g,
00089      Epetra_MultiVector* dg_dx,
00090      Epetra_MultiVector* dg_dxdot,
00091      Epetra_MultiVector* dg_dxdotdot,
00092      Epetra_MultiVector* dg_dp)
00093 {
00094 
00095   // Evaluate response g
00096   if (g != NULL)
00097     x.Norm2(&(*g)[0]);
00098 
00099   // Evaluate dg/dx
00100   if (dg_dx != NULL) {
00101     double nrm;
00102     if (g != NULL)
00103       nrm = (*g)[0];
00104     else
00105       x.Norm2(&nrm);
00106     dg_dx->Scale(1.0/nrm,x);
00107   }
00108 
00109   // Evaluate dg/dxdot
00110   if (dg_dxdot != NULL)
00111     dg_dxdot->PutScalar(0.0);
00112   if (dg_dxdotdot != NULL)
00113     dg_dxdotdot->PutScalar(0.0);
00114 
00115   // Evaluate dg/dp
00116   if (dg_dp != NULL)
00117     dg_dp->PutScalar(0.0);
00118 }

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