00001
00002
00003
00004
00005
00006
00007 #ifndef MULTISCALESTRESS_HPP
00008 #define MULTISCALESTRESS_HPP
00009
00010 #include "Phalanx_ConfigDefs.hpp"
00011 #include "Phalanx_Evaluator_WithBaseImpl.hpp"
00012 #include "Phalanx_Evaluator_Derived.hpp"
00013 #include "Phalanx_MDField.hpp"
00014
00015 namespace LCM {
00022
00023 enum MessageType {STRESS_TENSOR, STRAIN_TENSOR, TANGENT, DIE};
00024
00025
00026 template<typename EvalT, typename Traits>
00027 class MultiScaleStressBase : public PHX::EvaluatorWithBaseImpl<Traits>,
00028 public PHX::EvaluatorDerived<EvalT, Traits> {
00029
00030 public:
00031
00032 MultiScaleStressBase(const Teuchos::ParameterList& p);
00033
00034 void postRegistrationSetup(typename Traits::SetupData d,
00035 PHX::FieldManager<Traits>& vm);
00036
00037 void evaluateFields(typename Traits::EvalData d);
00038
00039 protected:
00040
00041 void calcStress(typename Traits::EvalData workset);
00042
00043
00044 void mesoBridgeStressRealType(PHX::MDField<RealType, Cell, QuadPoint, Dim, Dim>& stressFieldOut,
00045 PHX::MDField<RealType, Cell, QuadPoint, Dim, Dim>& stressFieldIn,
00046 typename Traits::EvalData workset);
00047
00048 PHX::MDField<RealType, Cell, QuadPoint, Dim, Dim> stressFieldRealType;
00049
00050
00051 struct MesoPt {
00052
00053 int cell;
00054 int qp;
00055
00056 };
00057
00058 typedef typename EvalT::ScalarT ScalarT;
00059 typedef typename EvalT::MeshScalarT MeshScalarT;
00060
00061
00062 PHX::MDField<ScalarT, Cell, QuadPoint, Dim, Dim> strain;
00063 PHX::MDField<ScalarT, Cell, QuadPoint> elasticModulus;
00064 PHX::MDField<ScalarT, Cell, QuadPoint> poissonsRatio;
00065
00066 unsigned int numQPs;
00067 unsigned int numDims;
00068
00069 int numMesoPEs;
00070 std::vector<double> exchanged_stresses;
00071 std::vector<MesoPt> loc_data;
00072 Teuchos::RCP<MPI_Comm> interCommunicator;
00073
00074
00075
00076 PHX::MDField<ScalarT, Cell, QuadPoint, Dim, Dim> stress;
00077
00078
00079
00080 void sendCellQPData(int cell, int qp, int toProc, MessageType type,
00081 PHX::MDField<RealType, Cell, QuadPoint, Dim, Dim>& stressFieldIn);
00082
00083 void rcvCellQPData(int procIDReached, MessageType type,
00084 PHX::MDField<RealType, Cell, QuadPoint, Dim, Dim>& stressFieldOut);
00085
00086 };
00087
00088
00089 template<typename EvalT, typename Traits> class MultiScaleStress;
00090
00091
00092
00093 template<typename EvalT, typename Traits>
00094 class MultiScaleStress : public MultiScaleStressBase<EvalT, Traits> {
00095 public:
00096 MultiScaleStress(Teuchos::ParameterList& p) : MultiScaleStressBase<EvalT, Traits>(p) {};
00097 };
00098
00099
00100
00101 template<typename Traits>
00102 class MultiScaleStress<PHAL::AlbanyTraits::Residual, Traits>
00103 : public MultiScaleStressBase<PHAL::AlbanyTraits::Residual, Traits> {
00104 public:
00105 MultiScaleStress(Teuchos::ParameterList& p) : MultiScaleStressBase<PHAL::AlbanyTraits::Residual, Traits>(p) {};
00106 void evaluateFields(typename Traits::EvalData d);
00107 };
00108
00109
00110 template<typename Traits>
00111 class MultiScaleStress<PHAL::AlbanyTraits::Jacobian, Traits>
00112 : public MultiScaleStressBase<PHAL::AlbanyTraits::Jacobian, Traits> {
00113 public:
00114 MultiScaleStress(Teuchos::ParameterList& p) : MultiScaleStressBase<PHAL::AlbanyTraits::Jacobian, Traits>(p) {};
00115 void evaluateFields(typename Traits::EvalData d);
00116 };
00117
00118
00119 template<typename Traits>
00120 class MultiScaleStress<PHAL::AlbanyTraits::Tangent, Traits>
00121 : public MultiScaleStressBase<PHAL::AlbanyTraits::Tangent, Traits> {
00122 public:
00123 MultiScaleStress(Teuchos::ParameterList& p) : MultiScaleStressBase<PHAL::AlbanyTraits::Tangent, Traits>(p) {};
00124 void evaluateFields(typename Traits::EvalData d);
00125 };
00126
00127 }
00128
00129 #endif