00001
00002
00003
00004
00005
00006
00007 #ifndef ALBANY_APPLICATION_HPP
00008 #define ALBANY_APPLICATION_HPP
00009
00010 #include <vector>
00011
00012 #include "Teuchos_RCP.hpp"
00013 #include "Teuchos_ArrayRCP.hpp"
00014 #include "Teuchos_ParameterList.hpp"
00015 #include "Teuchos_SerialDenseMatrix.hpp"
00016 #include "Teuchos_VerboseObject.hpp"
00017 #include "Teuchos_TimeMonitor.hpp"
00018
00019 #include "Epetra_Map.h"
00020 #include "Epetra_Vector.h"
00021 #include "Epetra_CrsMatrix.h"
00022 #include "Epetra_Import.h"
00023 #include "Epetra_Export.h"
00024
00025 #include "Albany_AbstractDiscretization.hpp"
00026 #include "Albany_AbstractProblem.hpp"
00027 #include "Albany_AbstractResponseFunction.hpp"
00028 #include "Albany_StateManager.hpp"
00029 #include "AAdapt_AdaptiveSolutionManager.hpp"
00030
00031 #ifdef ALBANY_CUTR
00032 #include "CUTR_CubitMeshMover.hpp"
00033 #include "STKMeshData.hpp"
00034 #endif
00035
00036 #include "Sacado_ScalarParameterLibrary.hpp"
00037 #include "Sacado_ScalarParameterVector.hpp"
00038 #include "Sacado_ParameterAccessor.hpp"
00039 #include "Sacado_ParameterRegistration.hpp"
00040
00041 #include "PHAL_AlbanyTraits.hpp"
00042 #include "Phalanx.hpp"
00043
00044 #include "Stokhos_OrthogPolyExpansion.hpp"
00045 #include "Stokhos_Quadrature.hpp"
00046 #include "Stokhos_EpetraVectorOrthogPoly.hpp"
00047 #include "Stokhos_EpetraMultiVectorOrthogPoly.hpp"
00048 #include "EpetraExt_MultiComm.h"
00049
00050 #include "LOCA_Epetra_Group.H"
00051
00052 #include "Teko_InverseLibrary.hpp"
00053
00054 #ifdef ALBANY_MOR
00055 #include "MOR/Albany_MORFacade.hpp"
00056 #endif
00057
00058 namespace Albany {
00059
00060 class Application :
00061 public Sacado::ParameterAccessor<PHAL::AlbanyTraits::Residual, SPL_Traits> {
00062 public:
00063
00064 enum SolutionMethod {Steady, Transient, Continuation, Eigensolve};
00065
00067 Application(const Teuchos::RCP<const Epetra_Comm>& comm,
00068 const Teuchos::RCP<Teuchos::ParameterList>& params,
00069 const Teuchos::RCP<const Epetra_Vector>& initial_guess =
00070 Teuchos::null);
00071
00073 ~Application();
00074
00076 Teuchos::RCP<Albany::AbstractDiscretization> getDiscretization() const;
00077
00079 Teuchos::RCP<Albany::AbstractProblem> getProblem() const;
00080
00082 Teuchos::RCP<const Epetra_Comm> getComm() const;
00083
00085 Teuchos::RCP<const Epetra_Map> getMap() const;
00086
00088 Teuchos::RCP<const Epetra_CrsGraph> getJacobianGraph() const;
00089
00091 Teuchos::RCP<Epetra_Operator> getPreconditioner();
00092
00094 Teuchos::RCP<AAdapt::AdaptiveSolutionManager> getAdaptSolMgr(){ return solMgr;}
00095
00097 Teuchos::RCP<ParamLib> getParamLib();
00098
00100 SolutionMethod getSolutionMethod() const {return solMethod; }
00101
00103 int getNumResponses() const;
00104
00106 Teuchos::RCP<AbstractResponseFunction> getResponse(int i) const;
00107
00109 bool suppliesPreconditioner() const;
00110
00112 Teuchos::RCP<Stokhos::OrthogPolyExpansion<int,double> >
00113 getStochasticExpansion();
00114
00116 #ifdef ALBANY_SG_MP
00117 void init_sg(
00118 const Teuchos::RCP<const Stokhos::OrthogPolyBasis<int,double> >& basis,
00119 const Teuchos::RCP<const Stokhos::Quadrature<int,double> >& quad,
00120 const Teuchos::RCP<Stokhos::OrthogPolyExpansion<int,double> >& expansion,
00121 const Teuchos::RCP<const EpetraExt::MultiComm>& multiComm);
00122 #endif //ALBANY_SG_MP
00123
00125
00128 void computeGlobalResidual(const double current_time,
00129 const Epetra_Vector* xdot,
00130 const Epetra_Vector* xdotdot,
00131 const Epetra_Vector& x,
00132 const Teuchos::Array<ParamVec>& p,
00133 Epetra_Vector& f);
00134
00136
00139 void computeGlobalJacobian(const double alpha,
00140 const double beta,
00141 const double omega,
00142 const double current_time,
00143 const Epetra_Vector* xdot,
00144 const Epetra_Vector* xdotdot,
00145 const Epetra_Vector& x,
00146 const Teuchos::Array<ParamVec>& p,
00147 Epetra_Vector* f,
00148 Epetra_CrsMatrix& jac);
00149
00151
00154 void computeGlobalPreconditioner(const Teuchos::RCP<Epetra_CrsMatrix>& jac,
00155 const Teuchos::RCP<Epetra_Operator>& prec);
00156
00158
00161 void computeGlobalTangent(const double alpha,
00162 const double beta,
00163 const double omega,
00164 const double current_time,
00165 bool sum_derivs,
00166 const Epetra_Vector* xdot,
00167 const Epetra_Vector* xdotdot,
00168 const Epetra_Vector& x,
00169 const Teuchos::Array<ParamVec>& p,
00170 ParamVec* deriv_p,
00171 const Epetra_MultiVector* Vx,
00172 const Epetra_MultiVector* Vxdot,
00173 const Epetra_MultiVector* Vxdotdot,
00174 const Epetra_MultiVector* Vp,
00175 Epetra_Vector* f,
00176 Epetra_MultiVector* JV,
00177 Epetra_MultiVector* fp);
00178
00180
00183 void evaluateResponse(
00184 int response_index,
00185 const double current_time,
00186 const Epetra_Vector* xdot,
00187 const Epetra_Vector* xdotdot,
00188 const Epetra_Vector& x,
00189 const Teuchos::Array<ParamVec>& p,
00190 Epetra_Vector& g);
00191
00193
00196 void evaluateResponseTangent(
00197 int response_index,
00198 const double alpha,
00199 const double beta,
00200 const double omega,
00201 const double current_time,
00202 bool sum_derivs,
00203 const Epetra_Vector* xdot,
00204 const Epetra_Vector* xdotdot,
00205 const Epetra_Vector& x,
00206 const Teuchos::Array<ParamVec>& p,
00207 ParamVec* deriv_p,
00208 const Epetra_MultiVector* Vxdot,
00209 const Epetra_MultiVector* Vxdotdot,
00210 const Epetra_MultiVector* Vx,
00211 const Epetra_MultiVector* Vp,
00212 Epetra_Vector* g,
00213 Epetra_MultiVector* gx,
00214 Epetra_MultiVector* gp);
00215
00217
00220 void evaluateResponseDerivative(
00221 int response_index,
00222 const double current_time,
00223 const Epetra_Vector* xdot,
00224 const Epetra_Vector* xdotdot,
00225 const Epetra_Vector& x,
00226 const Teuchos::Array<ParamVec>& p,
00227 ParamVec* deriv_p,
00228 Epetra_Vector* g,
00229 const EpetraExt::ModelEvaluator::Derivative& dg_dx,
00230 const EpetraExt::ModelEvaluator::Derivative& dg_dxdot,
00231 const EpetraExt::ModelEvaluator::Derivative& dg_dxdotdot,
00232 const EpetraExt::ModelEvaluator::Derivative& dg_dp);
00233
00234 #ifdef ALBANY_SG_MP
00235
00236
00239 void computeGlobalSGResidual(
00240 const double current_time,
00241 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00242 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00243 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00244 const Teuchos::Array<ParamVec>& p,
00245 const Teuchos::Array<int>& sg_p_index,
00246 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00247 Stokhos::EpetraVectorOrthogPoly& sg_f);
00248
00250
00253 void computeGlobalSGJacobian(
00254 const double alpha,
00255 const double beta,
00256 const double omega,
00257 const double current_time,
00258 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00259 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00260 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00261 const Teuchos::Array<ParamVec>& p,
00262 const Teuchos::Array<int>& sg_p_index,
00263 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00264 Stokhos::EpetraVectorOrthogPoly* sg_f,
00265 Stokhos::VectorOrthogPoly<Epetra_CrsMatrix>& sg_jac);
00266
00268
00271 void computeGlobalSGTangent(
00272 const double alpha,
00273 const double beta,
00274 const double omega,
00275 const double current_time,
00276 bool sum_derivs,
00277 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00278 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00279 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00280 const Teuchos::Array<ParamVec>& p,
00281 const Teuchos::Array<int>& sg_p_index,
00282 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00283 ParamVec* deriv_p,
00284 const Epetra_MultiVector* Vx,
00285 const Epetra_MultiVector* Vxdot,
00286 const Epetra_MultiVector* Vxdotdot,
00287 const Epetra_MultiVector* Vp,
00288 Stokhos::EpetraVectorOrthogPoly* sg_f,
00289 Stokhos::EpetraMultiVectorOrthogPoly* sg_JVx,
00290 Stokhos::EpetraMultiVectorOrthogPoly* sg_fVp);
00291
00293
00296 void evaluateSGResponse(
00297 int response_index,
00298 const double curr_time,
00299 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00300 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00301 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00302 const Teuchos::Array<ParamVec>& p,
00303 const Teuchos::Array<int>& sg_p_index,
00304 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00305 Stokhos::EpetraVectorOrthogPoly& sg_g);
00306
00308
00311 void
00312 evaluateSGResponseTangent(
00313 int response_index,
00314 const double alpha,
00315 const double beta,
00316 const double omega,
00317 const double current_time,
00318 bool sum_derivs,
00319 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00320 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00321 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00322 const Teuchos::Array<ParamVec>& p,
00323 const Teuchos::Array<int>& sg_p_index,
00324 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00325 ParamVec* deriv_p,
00326 const Epetra_MultiVector* Vx,
00327 const Epetra_MultiVector* Vxdot,
00328 const Epetra_MultiVector* Vxdotdot,
00329 const Epetra_MultiVector* Vp,
00330 Stokhos::EpetraVectorOrthogPoly* sg_g,
00331 Stokhos::EpetraMultiVectorOrthogPoly* sg_JV,
00332 Stokhos::EpetraMultiVectorOrthogPoly* sg_gp);
00333
00335
00338 void
00339 evaluateSGResponseDerivative(
00340 int response_index,
00341 const double current_time,
00342 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00343 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00344 const Stokhos::EpetraVectorOrthogPoly& sg_x,
00345 const Teuchos::Array<ParamVec>& p,
00346 const Teuchos::Array<int>& sg_p_index,
00347 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00348 ParamVec* deriv_p,
00349 Stokhos::EpetraVectorOrthogPoly* sg_g,
00350 const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dx,
00351 const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dxdot,
00352 const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dxdotdot,
00353 const EpetraExt::ModelEvaluator::SGDerivative& sg_dg_dp);
00354
00356
00359 void computeGlobalMPResidual(
00360 const double current_time,
00361 const Stokhos::ProductEpetraVector* mp_xdot,
00362 const Stokhos::ProductEpetraVector* mp_xdotdot,
00363 const Stokhos::ProductEpetraVector& mp_x,
00364 const Teuchos::Array<ParamVec>& p,
00365 const Teuchos::Array<int>& mp_p_index,
00366 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00367 Stokhos::ProductEpetraVector& mp_f);
00368
00370
00373 void computeGlobalMPJacobian(
00374 const double alpha,
00375 const double beta,
00376 const double omega,
00377 const double current_time,
00378 const Stokhos::ProductEpetraVector* mp_xdot,
00379 const Stokhos::ProductEpetraVector* mp_xdotdot,
00380 const Stokhos::ProductEpetraVector& mp_x,
00381 const Teuchos::Array<ParamVec>& p,
00382 const Teuchos::Array<int>& mp_p_index,
00383 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00384 Stokhos::ProductEpetraVector* mp_f,
00385 Stokhos::ProductContainer<Epetra_CrsMatrix>& mp_jac);
00386
00388
00391 void computeGlobalMPTangent(
00392 const double alpha,
00393 const double beta,
00394 const double omega,
00395 const double current_time,
00396 bool sum_derivs,
00397 const Stokhos::ProductEpetraVector* mp_xdot,
00398 const Stokhos::ProductEpetraVector* mp_xdotdot,
00399 const Stokhos::ProductEpetraVector& mp_x,
00400 const Teuchos::Array<ParamVec>& p,
00401 const Teuchos::Array<int>& mp_p_index,
00402 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00403 ParamVec* deriv_p,
00404 const Epetra_MultiVector* Vx,
00405 const Epetra_MultiVector* Vxdot,
00406 const Epetra_MultiVector* Vxdotdot,
00407 const Epetra_MultiVector* Vp,
00408 Stokhos::ProductEpetraVector* mp_f,
00409 Stokhos::ProductEpetraMultiVector* mp_JVx,
00410 Stokhos::ProductEpetraMultiVector* mp_fVp);
00411
00413
00416 void evaluateMPResponse(
00417 int response_index,
00418 const double curr_time,
00419 const Stokhos::ProductEpetraVector* mp_xdot,
00420 const Stokhos::ProductEpetraVector* mp_xdotdot,
00421 const Stokhos::ProductEpetraVector& mp_x,
00422 const Teuchos::Array<ParamVec>& p,
00423 const Teuchos::Array<int>& mp_p_index,
00424 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00425 Stokhos::ProductEpetraVector& mp_g);
00426
00428
00431 void
00432 evaluateMPResponseTangent(
00433 int response_index,
00434 const double alpha,
00435 const double beta,
00436 const double omega,
00437 const double current_time,
00438 bool sum_derivs,
00439 const Stokhos::ProductEpetraVector* mp_xdot,
00440 const Stokhos::ProductEpetraVector* mp_xdotdot,
00441 const Stokhos::ProductEpetraVector& mp_x,
00442 const Teuchos::Array<ParamVec>& p,
00443 const Teuchos::Array<int>& mp_p_index,
00444 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00445 ParamVec* deriv_p,
00446 const Epetra_MultiVector* Vx,
00447 const Epetra_MultiVector* Vxdot,
00448 const Epetra_MultiVector* Vxdotdot,
00449 const Epetra_MultiVector* Vp,
00450 Stokhos::ProductEpetraVector* mp_g,
00451 Stokhos::ProductEpetraMultiVector* mp_JV,
00452 Stokhos::ProductEpetraMultiVector* mp_gp);
00453
00455
00458 void
00459 evaluateMPResponseDerivative(
00460 int response_index,
00461 const double current_time,
00462 const Stokhos::ProductEpetraVector* mp_xdot,
00463 const Stokhos::ProductEpetraVector* mp_xdotdot,
00464 const Stokhos::ProductEpetraVector& mp_x,
00465 const Teuchos::Array<ParamVec>& p,
00466 const Teuchos::Array<int>& mp_p_index,
00467 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00468 ParamVec* deriv_p,
00469 Stokhos::ProductEpetraVector* mp_g,
00470 const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dx,
00471 const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dxdot,
00472 const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dxdotdot,
00473 const EpetraExt::ModelEvaluator::MPDerivative& mp_dg_dp);
00474 #endif //ALBANY_SG_MP
00475
00477 PHAL::AlbanyTraits::Residual::ScalarT& getValue(const std::string &n);
00478
00480 StateManager& getStateMgr() {return stateMgr; }
00481
00483 void evaluateStateFieldManager(const double current_time,
00484 const Epetra_Vector* xdot,
00485 const Epetra_Vector* xdotdot,
00486 const Epetra_Vector& x);
00487
00489 int getNumWorksets() {
00490 return disc->getWsElNodeEqID().size();
00491 }
00492
00493 bool is_adjoint;
00494
00495 private:
00496
00498 Application(const Application&);
00499
00501 Application& operator=(const Application&);
00502
00504 Teuchos::RCP<Epetra_Operator> buildWrappedOperator(
00505 const Teuchos::RCP<Epetra_Operator>& Jac,
00506 const Teuchos::RCP<Epetra_Operator>& wrapInput,
00507 bool reorder=false) const;
00508
00510 void registerShapeParameters();
00511
00512 void defineTimers();
00513
00514 public:
00515
00517 template <typename EvalT>
00518 void loadWorksetBucketInfo(PHAL::Workset& workset, const int& ws);
00519
00521 void loadBasicWorksetInfo(
00522 PHAL::Workset& workset,
00523 double current_time);
00524
00525 void loadWorksetJacobianInfo(PHAL::Workset& workset,
00526 const double& alpha, const double& beta, const double& omega);
00527
00529 void loadWorksetNodesetInfo(PHAL::Workset& workset);
00530
00532 void loadWorksetSidesetInfo(PHAL::Workset& workset, const int ws);
00533
00534 void setupBasicWorksetInfo(
00535 PHAL::Workset& workset,
00536 double current_time,
00537 const Epetra_Vector* xdot,
00538 const Epetra_Vector* xdotdot,
00539 const Epetra_Vector* x,
00540 const Teuchos::Array<ParamVec>& p);
00541
00542 #ifdef ALBANY_SG_MP
00543 void setupBasicWorksetInfo(
00544 PHAL::Workset& workset,
00545 double current_time,
00546 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00547 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00548 const Stokhos::EpetraVectorOrthogPoly* sg_x,
00549 const Teuchos::Array<ParamVec>& p,
00550 const Teuchos::Array<int>& sg_p_index,
00551 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals);
00552
00553 void setupBasicWorksetInfo(
00554 PHAL::Workset& workset,
00555 double current_time,
00556 const Stokhos::ProductEpetraVector* mp_xdot,
00557 const Stokhos::ProductEpetraVector* mp_xdotdot,
00558 const Stokhos::ProductEpetraVector* mp_x,
00559 const Teuchos::Array<ParamVec>& p,
00560 const Teuchos::Array<int>& mp_p_index,
00561 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals);
00562 #endif //ALBANY_SG_MP
00563
00564 void setupTangentWorksetInfo(
00565 PHAL::Workset& workset,
00566 double current_time,
00567 bool sum_derivs,
00568 const Epetra_Vector* xdot,
00569 const Epetra_Vector* xdotdot,
00570 const Epetra_Vector* x,
00571 const Teuchos::Array<ParamVec>& p,
00572 ParamVec* deriv_p,
00573 const Epetra_MultiVector* Vxdot,
00574 const Epetra_MultiVector* Vxdotdot,
00575 const Epetra_MultiVector* Vx,
00576 const Epetra_MultiVector* Vp);
00577
00578 #ifdef ALBANY_SG_MP
00579 void setupTangentWorksetInfo(
00580 PHAL::Workset& workset,
00581 double current_time,
00582 bool sum_derivs,
00583 const Stokhos::EpetraVectorOrthogPoly* sg_xdot,
00584 const Stokhos::EpetraVectorOrthogPoly* sg_xdotdot,
00585 const Stokhos::EpetraVectorOrthogPoly* sg_x,
00586 const Teuchos::Array<ParamVec>& p,
00587 ParamVec* deriv_p,
00588 const Teuchos::Array<int>& sg_p_index,
00589 const Teuchos::Array< Teuchos::Array<SGType> >& sg_p_vals,
00590 const Epetra_MultiVector* Vxdot,
00591 const Epetra_MultiVector* Vxdotdot,
00592 const Epetra_MultiVector* Vx,
00593 const Epetra_MultiVector* Vp);
00594
00595 void setupTangentWorksetInfo(
00596 PHAL::Workset& workset,
00597 double current_time,
00598 bool sum_derivs,
00599 const Stokhos::ProductEpetraVector* mp_xdot,
00600 const Stokhos::ProductEpetraVector* mp_xdotdot,
00601 const Stokhos::ProductEpetraVector* mp_x,
00602 const Teuchos::Array<ParamVec>& p,
00603 ParamVec* deriv_p,
00604 const Teuchos::Array<int>& mp_p_index,
00605 const Teuchos::Array< Teuchos::Array<MPType> >& mp_p_vals,
00606 const Epetra_MultiVector* Vxdot,
00607 const Epetra_MultiVector* Vxdotdot,
00608 const Epetra_MultiVector* Vx,
00609 const Epetra_MultiVector* Vp);
00610 #endif //ALBANY_SG_MP
00611
00612 void postRegSetup(std::string eval);
00613
00614 #ifdef ALBANY_MOR
00615 Teuchos::RCP<MORFacade> getMorFacade();
00616 #endif
00617
00618 protected:
00619
00621 Teuchos::RCP<const Epetra_Comm> comm;
00622
00624 Teuchos::RCP<Teuchos::FancyOStream> out;
00625
00627 Teuchos::RCP<Albany::AbstractDiscretization> disc;
00628
00630 Teuchos::RCP<Albany::AbstractProblem> problem;
00631
00633 Teuchos::RCP<ParamLib> paramLib;
00634
00636 Teuchos::RCP<AAdapt::AdaptiveSolutionManager> solMgr;
00637
00639 Teuchos::Array< Teuchos::RCP<Albany::AbstractResponseFunction> > responses;
00640
00642 Teuchos::ArrayRCP<Teuchos::RCP<PHX::FieldManager<PHAL::AlbanyTraits> > > fm;
00643
00645 Teuchos::RCP<PHX::FieldManager<PHAL::AlbanyTraits> > dfm;
00646
00648 Teuchos::ArrayRCP<Teuchos::RCP<PHX::FieldManager<PHAL::AlbanyTraits> > > nfm;
00649
00651 Teuchos::Array< Teuchos::RCP<PHX::FieldManager<PHAL::AlbanyTraits> > > sfm;
00652
00654 Teuchos::RCP<const Stokhos::OrthogPolyBasis<int,double> > sg_basis;
00655
00657 Teuchos::RCP<const Stokhos::Quadrature<int,double> > sg_quad;
00658
00660 Teuchos::RCP<Stokhos::OrthogPolyExpansion<int,double> > sg_expansion;
00661
00663 Teuchos::RCP<const EpetraExt::MultiComm> product_comm;
00664
00666 Teuchos::RCP<const Epetra_BlockMap> sg_overlap_map;
00667
00669 Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > sg_overlapped_x;
00670
00672 Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > sg_overlapped_xdot;
00673 Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > sg_overlapped_xdotdot;
00674
00676 Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > sg_overlapped_f;
00677
00679 Teuchos::RCP< Stokhos::VectorOrthogPoly<Epetra_CrsMatrix> > sg_overlapped_jac;
00680
00682 Teuchos::RCP< Stokhos::ProductEpetraVector > mp_overlapped_x;
00683
00685 Teuchos::RCP< Stokhos::ProductEpetraVector > mp_overlapped_xdot;
00686 Teuchos::RCP< Stokhos::ProductEpetraVector > mp_overlapped_xdotdot;
00687
00689 Teuchos::RCP< Stokhos::ProductEpetraVector > mp_overlapped_f;
00690
00692 Teuchos::RCP< Stokhos::ProductContainer<Epetra_CrsMatrix> > mp_overlapped_jac;
00693
00695 bool physicsBasedPreconditioner;
00696 Teuchos::RCP<Teuchos::ParameterList> tekoParams;
00697
00699 SolutionMethod solMethod;
00700
00702
00703
00704
00705
00706
00707 int writeToMatrixMarketJac;
00708 int writeToMatrixMarketRes;
00710 int writeToCoutJac;
00711 int writeToCoutRes;
00712
00714 bool shapeParamsHaveBeenReset;
00715 std::vector<RealType> shapeParams;
00716 std::vector<std::string> shapeParamNames;
00717 #ifdef ALBANY_CUTR
00718 Teuchos::RCP<CUTR::CubitMeshMover> meshMover;
00719 #endif
00720
00721 unsigned int neq;
00722
00724 Teuchos::RCP<Teko::InverseLibrary> inverseLib;
00725 Teuchos::RCP<Teko::InverseFactory> inverseFac;
00726 Teuchos::RCP<Epetra_Operator> wrappedJac;
00727 std::vector<int> blockDecomp;
00728
00729 std::set<std::string> setupSet;
00730 mutable int phxGraphVisDetail;
00731 mutable int stateGraphVisDetail;
00732
00733 StateManager stateMgr;
00734
00735 bool morphFromInit;
00736 bool ignore_residual_in_jacobian;
00737
00739
00740 double perturbBetaForDirichlets;
00741
00742 void determinePiroSolver(const Teuchos::RCP<Teuchos::ParameterList>& topLevelParams);
00743
00744 #ifdef ALBANY_MOR
00745 Teuchos::RCP<MORFacade> morFacade;
00746 #endif
00747 };
00748 }
00749
00750 template <typename EvalT>
00751 void Albany::Application::loadWorksetBucketInfo(PHAL::Workset& workset,
00752 const int& ws)
00753 {
00754
00755 const WorksetArray<Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::ArrayRCP<int> > > >::type&
00756 wsElNodeEqID = disc->getWsElNodeEqID();
00757 const WorksetArray<Teuchos::ArrayRCP<Teuchos::ArrayRCP<int> > >::type&
00758 wsElNodeID = disc->getWsElNodeID();
00759 const WorksetArray<Teuchos::ArrayRCP<Teuchos::ArrayRCP<double*> > >::type&
00760 coords = disc->getCoords();
00761 const WorksetArray<Teuchos::ArrayRCP<Teuchos::ArrayRCP<double> > >::type&
00762 sHeight = disc->getSurfaceHeight();
00763 const WorksetArray<Teuchos::ArrayRCP<double> >::type&
00764 temperature = disc->getTemperature();
00765 const WorksetArray<Teuchos::ArrayRCP<Teuchos::ArrayRCP<double> > >::type&
00766 basalFriction = disc->getBasalFriction();
00767 const WorksetArray<Teuchos::ArrayRCP<Teuchos::ArrayRCP<double> > >::type&
00768 thickness = disc->getThickness();
00769 const WorksetArray<Teuchos::ArrayRCP<double> >::type&
00770 flowFactor = disc->getFlowFactor();
00771 const WorksetArray<Teuchos::ArrayRCP<Teuchos::ArrayRCP<double*> > >::type&
00772 surfaceVelocity = disc->getSurfaceVelocity();
00773 const WorksetArray<Teuchos::ArrayRCP<Teuchos::ArrayRCP<double*> > >::type&
00774 velocityRMS = disc->getVelocityRMS();
00775 const WorksetArray<std::string>::type& wsEBNames = disc->getWsEBNames();
00776
00777 workset.numCells = wsElNodeEqID[ws].size();
00778 workset.wsElNodeEqID = wsElNodeEqID[ws];
00779 workset.wsElNodeID = wsElNodeID[ws];
00780 workset.wsCoords = coords[ws];
00781 workset.wsSHeight = sHeight[ws];
00782 workset.wsTemperature = temperature[ws];
00783 workset.wsBasalFriction = basalFriction[ws];
00784 workset.wsThickness = thickness[ws];
00785 workset.wsFlowFactor = flowFactor[ws];
00786 workset.wsSurfaceVelocity = surfaceVelocity[ws];
00787 workset.wsVelocityRMS = velocityRMS[ws];
00788 workset.EBName = wsEBNames[ws];
00789 workset.wsIndex = ws;
00790
00791
00792
00793
00794 loadWorksetSidesetInfo(workset, ws);
00795
00796 workset.stateArrayPtr = &stateMgr.getStateArray(Albany::StateManager::ELEM, ws);
00797 workset.eigenDataPtr = stateMgr.getEigenData();
00798 workset.auxDataPtr = stateMgr.getAuxData();
00799
00800 PHAL::BuildSerializer<EvalT> bs(workset);
00801 }
00802
00803 #endif // ALBANY_APPLICATION_HPP