00001
00002
00003
00004
00005
00006
00007 #include <Intrepid_MiniTensor.h>
00008 #include <Teuchos_TestForException.hpp>
00009 #include <Phalanx_DataLayout.hpp>
00010 #include <typeinfo>
00011
00012 namespace LCM
00013 {
00014
00015
00016 template<typename EvalT, typename Traits>
00017 AnisotropicHyperelasticDamageModel<EvalT, Traits>::
00018 AnisotropicHyperelasticDamageModel(Teuchos::ParameterList* p,
00019 const Teuchos::RCP<Albany::Layouts>& dl) :
00020 LCM::ConstitutiveModel<EvalT, Traits>(p, dl),
00021 k_f1_(p->get<RealType>("Fiber 1 k", 1.0)),
00022 q_f1_(p->get<RealType>("Fiber 1 q", 1.0)),
00023 volume_fraction_f1_(p->get<RealType>("Fiber 1 volume fraction", 0.0)),
00024 max_damage_f1_(p->get<RealType>("Fiber 1 maximum damage", 1.0)),
00025 saturation_f1_(p->get<RealType>("Fiber 1 damage saturation", 0.0)),
00026 k_f2_(p->get<RealType>("Fiber 2 k", 1.0)),
00027 q_f2_(p->get<RealType>("Fiber 2 q", 1.0)),
00028 volume_fraction_f2_(p->get<RealType>("Fiber 2 volume fraction", 0.0)),
00029 max_damage_f2_(p->get<RealType>("Fiber 2 maximum damage", 1.0)),
00030 saturation_f2_(p->get<RealType>("Fiber 2 damage saturation", 0.0)),
00031 volume_fraction_m_(p->get<RealType>("Matrix volume fraction", 1.0)),
00032 max_damage_m_(p->get<RealType>("Matrix maximum damage", 1.0)),
00033 saturation_m_(p->get<RealType>("Matrix damage saturation", 0.0)),
00034 direction_f1_(p->get<Teuchos::Array<RealType> >("Fiber 1 Orientation Vector")
00035 .toVector()),
00036 direction_f2_(p->get<Teuchos::Array<RealType> >("Fiber 2 Orientation Vector")
00037 .toVector())
00038 {
00039 std::string F_string = (*field_name_map_)["F"];
00040 std::string J_string = (*field_name_map_)["J"];
00041
00042
00043 this->dep_field_map_.insert(std::make_pair(F_string, dl->qp_tensor));
00044 this->dep_field_map_.insert(std::make_pair(J_string, dl->qp_scalar));
00045 this->dep_field_map_.insert(std::make_pair("Poissons Ratio", dl->qp_scalar));
00046 this->dep_field_map_.insert(std::make_pair("Elastic Modulus", dl->qp_scalar));
00047
00048
00049 std::string cauchy_string = (*field_name_map_)["Cauchy_Stress"];
00050 std::string matrix_energy_string = (*field_name_map_)["Matrix_Energy"];
00051 std::string f1_energy_string = (*field_name_map_)["F1_Energy"];
00052 std::string f2_energy_string = (*field_name_map_)["F2_Energy"];
00053 std::string matrix_damage_string = (*field_name_map_)["Matrix_Damage"];
00054 std::string f1_damage_string = (*field_name_map_)["F1_Damage"];
00055 std::string f2_damage_string = (*field_name_map_)["F2_Damage"];
00056
00057
00058 this->eval_field_map_.insert(std::make_pair(cauchy_string, dl->qp_tensor));
00059 this->eval_field_map_.insert(
00060 std::make_pair(matrix_energy_string, dl->qp_scalar));
00061 this->eval_field_map_.insert(std::make_pair(f1_energy_string, dl->qp_scalar));
00062 this->eval_field_map_.insert(std::make_pair(f2_energy_string, dl->qp_scalar));
00063 this->eval_field_map_.insert(
00064 std::make_pair(matrix_damage_string, dl->qp_scalar));
00065 this->eval_field_map_.insert(std::make_pair(f1_damage_string, dl->qp_scalar));
00066 this->eval_field_map_.insert(std::make_pair(f2_damage_string, dl->qp_scalar));
00067
00068
00069
00070
00071 this->num_state_variables_++;
00072 this->state_var_names_.push_back(cauchy_string);
00073 this->state_var_layouts_.push_back(dl->qp_tensor);
00074 this->state_var_init_types_.push_back("scalar");
00075 this->state_var_init_values_.push_back(0.0);
00076 this->state_var_old_state_flags_.push_back(false);
00077 this->state_var_output_flags_.push_back(p->get<bool>("Output Cauchy Stress", false));
00078
00079
00080 this->num_state_variables_++;
00081 this->state_var_names_.push_back(matrix_energy_string);
00082 this->state_var_layouts_.push_back(dl->qp_scalar);
00083 this->state_var_init_types_.push_back("scalar");
00084 this->state_var_init_values_.push_back(0.0);
00085 this->state_var_old_state_flags_.push_back(true);
00086 this->state_var_output_flags_.push_back(p->get<bool>("Output Matrix Energy", false));
00087
00088
00089 this->num_state_variables_++;
00090 this->state_var_names_.push_back(f1_energy_string);
00091 this->state_var_layouts_.push_back(dl->qp_scalar);
00092 this->state_var_init_types_.push_back("scalar");
00093 this->state_var_init_values_.push_back(0.0);
00094 this->state_var_old_state_flags_.push_back(true);
00095 this->state_var_output_flags_.push_back(p->get<bool>("Output Fiber 1 Energy", false));
00096
00097
00098 this->num_state_variables_++;
00099 this->state_var_names_.push_back(f2_energy_string);
00100 this->state_var_layouts_.push_back(dl->qp_scalar);
00101 this->state_var_init_types_.push_back("scalar");
00102 this->state_var_init_values_.push_back(0.0);
00103 this->state_var_old_state_flags_.push_back(true);
00104 this->state_var_output_flags_.push_back(p->get<bool>("Output Fiber 2 Energy", false));
00105
00106
00107 this->num_state_variables_++;
00108 this->state_var_names_.push_back(matrix_damage_string);
00109 this->state_var_layouts_.push_back(dl->qp_scalar);
00110 this->state_var_init_types_.push_back("scalar");
00111 this->state_var_init_values_.push_back(0.0);
00112 this->state_var_old_state_flags_.push_back(true);
00113 this->state_var_output_flags_.push_back(p->get<bool>("Output Matrix Damage", false));
00114
00115
00116 this->num_state_variables_++;
00117 this->state_var_names_.push_back(f1_damage_string);
00118 this->state_var_layouts_.push_back(dl->qp_scalar);
00119 this->state_var_init_types_.push_back("scalar");
00120 this->state_var_init_values_.push_back(0.0);
00121 this->state_var_old_state_flags_.push_back(true);
00122 this->state_var_output_flags_.push_back(p->get<bool>("Output Fiber 1 Damage", false));
00123
00124
00125 this->num_state_variables_++;
00126 this->state_var_names_.push_back(f2_damage_string);
00127 this->state_var_layouts_.push_back(dl->qp_scalar);
00128 this->state_var_init_types_.push_back("scalar");
00129 this->state_var_init_values_.push_back(0.0);
00130 this->state_var_old_state_flags_.push_back(true);
00131 this->state_var_output_flags_.push_back(p->get<bool>("Output Fiber 2 Damage", false));
00132
00133 }
00134
00135 template<typename EvalT, typename Traits>
00136 void AnisotropicHyperelasticDamageModel<EvalT, Traits>::
00137 computeState(typename Traits::EvalData workset,
00138 std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > dep_fields,
00139 std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > eval_fields)
00140 {
00141 bool print = false;
00142
00143
00144
00145
00146 std::string F_string = (*field_name_map_)["F"];
00147 std::string J_string = (*field_name_map_)["J"];
00148 std::string cauchy_string = (*field_name_map_)["Cauchy_Stress"];
00149 std::string matrix_energy_string = (*field_name_map_)["Matrix_Energy"];
00150 std::string f1_energy_string = (*field_name_map_)["F1_Energy"];
00151 std::string f2_energy_string = (*field_name_map_)["F2_Energy"];
00152 std::string matrix_damage_string = (*field_name_map_)["Matrix_Damage"];
00153 std::string f1_damage_string = (*field_name_map_)["F1_Damage"];
00154 std::string f2_damage_string = (*field_name_map_)["F2_Damage"];
00155
00156
00157 PHX::MDField<ScalarT> def_grad = *dep_fields[F_string];
00158 PHX::MDField<ScalarT> J = *dep_fields[J_string];
00159 PHX::MDField<ScalarT> poissons_ratio = *dep_fields["Poissons Ratio"];
00160 PHX::MDField<ScalarT> elastic_modulus = *dep_fields["Elastic Modulus"];
00161
00162
00163 PHX::MDField<ScalarT> stress = *eval_fields[cauchy_string];
00164 PHX::MDField<ScalarT> energy_m = *eval_fields[matrix_energy_string];
00165 PHX::MDField<ScalarT> energy_f1 = *eval_fields[f1_energy_string];
00166 PHX::MDField<ScalarT> energy_f2 = *eval_fields[f2_energy_string];
00167 PHX::MDField<ScalarT> damage_m = *eval_fields[matrix_damage_string];
00168 PHX::MDField<ScalarT> damage_f1 = *eval_fields[f1_damage_string];
00169 PHX::MDField<ScalarT> damage_f2 = *eval_fields[f2_damage_string];
00170
00171
00172 Albany::MDArray energy_m_old =
00173 (*workset.stateArrayPtr)[matrix_energy_string + "_old"];
00174 Albany::MDArray energy_f1_old =
00175 (*workset.stateArrayPtr)[f1_energy_string + "_old"];
00176 Albany::MDArray energy_f2_old =
00177 (*workset.stateArrayPtr)[f2_energy_string + "_old"];
00178
00179 ScalarT kappa, mu, Jm53, Jm23, p, I4_f1, I4_f2;
00180 ScalarT alpha_f1, alpha_f2, alpha_m;
00181
00182
00183 Intrepid::Tensor<ScalarT> I(Intrepid::eye<ScalarT>(num_dims_));
00184 Intrepid::Tensor<ScalarT> F(num_dims_), s(num_dims_), b(num_dims_), C(
00185 num_dims_);
00186 Intrepid::Tensor<ScalarT> sigma_m(num_dims_), sigma_f1(num_dims_), sigma_f2(
00187 num_dims_);
00188 Intrepid::Tensor<ScalarT> M1dyadM1(num_dims_), M2dyadM2(num_dims_);
00189 Intrepid::Tensor<ScalarT> S0_f1(num_dims_), S0_f2(num_dims_);
00190
00191 Intrepid::Vector<ScalarT> M1(num_dims_), M2(num_dims_);
00192
00193 for (std::size_t cell = 0; cell < workset.numCells; ++cell) {
00194 for (std::size_t pt = 0; pt < num_pts_; ++pt) {
00195
00196 kappa = elastic_modulus(cell, pt)
00197 / (3. * (1. - 2. * poissons_ratio(cell, pt)));
00198 mu = elastic_modulus(cell, pt) / (2. * (1. + poissons_ratio(cell, pt)));
00199 Jm53 = std::pow(J(cell, pt), -5. / 3.);
00200 Jm23 = std::pow(J(cell, pt), -2. / 3.);
00201 F.fill(&def_grad(cell, pt, 0, 0));
00202
00203
00204 b = F * Intrepid::transpose(F);
00205 s = mu * Jm53 * Intrepid::dev(b);
00206
00207 p = 0.5 * kappa * (J(cell, pt) - 1. / (J(cell, pt)));
00208
00209 sigma_m = s + p * I;
00210
00211
00212 energy_m(cell, pt) = 0.5 * kappa
00213 * (0.5 * (J(cell, pt) * J(cell, pt) - 1.0) - std::log(J(cell, pt)))
00214 + 0.5 * mu * (Jm23 * Intrepid::trace(b) - 3.0);
00215
00216
00217 alpha_m = energy_m_old(cell, pt);
00218 if (energy_m(cell, pt) > alpha_m) alpha_m = energy_m(cell, pt);
00219
00220 damage_m(cell, pt) = max_damage_m_
00221 * (1 - std::exp(-alpha_m / saturation_m_));
00222
00223
00224
00225
00226 C = Intrepid::transpose(F) * F;
00227
00228
00229
00230
00231 for (std::size_t i = 0; i < num_dims_; ++i) {
00232 M1(i) = direction_f1_[i];
00233 }
00234 M1 = M1 / norm(M1);
00235
00236
00237 for (std::size_t i = 0; i < num_dims_; ++i) {
00238 M2(i) = direction_f2_[i];
00239 }
00240 M2 = M2 / norm(M2);
00241
00242
00243 I4_f1 = Intrepid::dot(M1, Intrepid::dot(C, M1));
00244 I4_f2 = Intrepid::dot(M2, Intrepid::dot(C, M2));
00245 M1dyadM1 = Intrepid::dyad(M1, M1);
00246 M2dyadM2 = Intrepid::dyad(M2, M2);
00247
00248
00249 S0_f1 = (4.0 * k_f1_ * (I4_f1 - 1.0)
00250 * std::exp(q_f1_ * (I4_f1 - 1) * (I4_f1 - 1))) * M1dyadM1;
00251 S0_f2 = (4.0 * k_f2_ * (I4_f2 - 1.0)
00252 * std::exp(q_f2_ * (I4_f2 - 1) * (I4_f2 - 1))) * M2dyadM2;
00253
00254
00255 energy_f1(cell, pt) = k_f1_
00256 * (std::exp(q_f1_ * (I4_f1 - 1) * (I4_f1 - 1)) - 1) / q_f1_;
00257 energy_f2(cell, pt) = k_f2_
00258 * (std::exp(q_f2_ * (I4_f2 - 1) * (I4_f2 - 1)) - 1) / q_f2_;
00259
00260
00261 sigma_f1 = (1.0 / J(cell, pt))
00262 * Intrepid::dot(F, Intrepid::dot(S0_f1, Intrepid::transpose(F)));
00263 sigma_f2 = (1.0 / J(cell, pt))
00264 * Intrepid::dot(F, Intrepid::dot(S0_f2, Intrepid::transpose(F)));
00265
00266
00267 alpha_f1 = energy_f1_old(cell, pt);
00268 alpha_f2 = energy_f2_old(cell, pt);
00269
00270 if (energy_f1(cell, pt) > alpha_f1) alpha_f1 = energy_f1(cell, pt);
00271
00272 if (energy_f2(cell, pt) > alpha_f2) alpha_f2 = energy_f2(cell, pt);
00273
00274
00275 damage_f1(cell, pt) = max_damage_f1_
00276 * (1 - std::exp(-alpha_f1 / saturation_f1_));
00277 damage_f2(cell, pt) = max_damage_f2_
00278 * (1 - std::exp(-alpha_f2 / saturation_f2_));
00279
00280
00281 for (std::size_t i(0); i < num_dims_; ++i) {
00282 for (std::size_t j(0); j < num_dims_; ++j) {
00283 stress(cell, pt, i, j) =
00284 volume_fraction_m_ * (1 - damage_m(cell, pt)) * sigma_m(i, j)
00285 + volume_fraction_f1_ * (1 - damage_f1(cell, pt))
00286 * sigma_f1(i, j)
00287 + volume_fraction_f2_ * (1 - damage_f2(cell, pt))
00288 * sigma_f2(i, j);
00289 }
00290 }
00291
00292 if (print) {
00293 std::cout << " matrix damage: " << damage_m(cell, pt) << std::endl;
00294 std::cout << " matrix energy: " << energy_m(cell, pt) << std::endl;
00295 std::cout << " fiber1 damage: " << damage_f1(cell, pt) << std::endl;
00296 std::cout << " fiber1 energy: " << energy_f1(cell, pt) << std::endl;
00297 std::cout << " fiber2 damage: " << damage_f2(cell, pt) << std::endl;
00298 std::cout << " fiber2 energy: " << energy_f2(cell, pt) << std::endl;
00299 }
00300 }
00301 }
00302 }
00303
00304 }
00305