00001
00002
00003
00004
00005
00006
00007 #include <Intrepid_MiniTensor.h>
00008 #include "Teuchos_TestForException.hpp"
00009 #include "Phalanx_DataLayout.hpp"
00010
00011 #include "LocalNonlinearSolver.hpp"
00012
00013 namespace LCM
00014 {
00015
00016
00017 template<typename EvalT, typename Traits>
00018 GursonModel<EvalT, Traits>::
00019 GursonModel(Teuchos::ParameterList* p,
00020 const Teuchos::RCP<Albany::Layouts>& dl) :
00021 LCM::ConstitutiveModel<EvalT, Traits>(p, dl),
00022 sat_mod_(p->get<RealType>("Saturation Modulus", 0.0)),
00023 sat_exp_(p->get<RealType>("Saturation Exponent", 0.0)),
00024 f0_(p->get<RealType>("Initial Void Volume", 0.0)),
00025 kw_(p->get<RealType>("Shear Damage Parameter", 0.0)),
00026 eN_(p->get<RealType>("Void Nucleation Parameter eN", 0.0)),
00027 sN_(p->get<RealType>("Void Nucleation Parameter sN", 0.1)),
00028 fN_(p->get<RealType>("Void Nucleation Parameter fN", 0.0)),
00029 fc_(p->get<RealType>("Critical Void Volume", 1.0)),
00030 ff_(p->get<RealType>("Failure Void Volume", 1.0)),
00031 q1_(p->get<RealType>("Yield Parameter q1", 1.0)),
00032 q2_(p->get<RealType>("Yield Parameter q2", 1.0)),
00033 q3_(p->get<RealType>("Yield Parameter q3", 1.0))
00034 {
00035
00036 this->dep_field_map_.insert(std::make_pair("F", dl->qp_tensor));
00037 this->dep_field_map_.insert(std::make_pair("J", dl->qp_scalar));
00038 this->dep_field_map_.insert(std::make_pair("Poissons Ratio", dl->qp_scalar));
00039 this->dep_field_map_.insert(std::make_pair("Elastic Modulus", dl->qp_scalar));
00040 this->dep_field_map_.insert(std::make_pair("Yield Strength", dl->qp_scalar));
00041 this->dep_field_map_.insert(
00042 std::make_pair("Hardening Modulus", dl->qp_scalar));
00043
00044
00045 std::string cauchy_string = (*field_name_map_)["Cauchy_Stress"];
00046 std::string Fp_string = (*field_name_map_)["Fp"];
00047 std::string eqps_string = (*field_name_map_)["eqps"];
00048 std::string void_string = (*field_name_map_)["Void_Volume"];
00049
00050
00051 this->eval_field_map_.insert(std::make_pair(cauchy_string, dl->qp_tensor));
00052 this->eval_field_map_.insert(std::make_pair(Fp_string, dl->qp_tensor));
00053 this->eval_field_map_.insert(std::make_pair(eqps_string, dl->qp_scalar));
00054 this->eval_field_map_.insert(std::make_pair(void_string, dl->qp_scalar));
00055
00056
00057
00058
00059 this->num_state_variables_++;
00060 this->state_var_names_.push_back(cauchy_string);
00061 this->state_var_layouts_.push_back(dl->qp_tensor);
00062 this->state_var_init_types_.push_back("scalar");
00063 this->state_var_init_values_.push_back(0.0);
00064 this->state_var_old_state_flags_.push_back(false);
00065 this->state_var_output_flags_.push_back(true);
00066
00067
00068 this->num_state_variables_++;
00069 this->state_var_names_.push_back(Fp_string);
00070 this->state_var_layouts_.push_back(dl->qp_tensor);
00071 this->state_var_init_types_.push_back("identity");
00072 this->state_var_init_values_.push_back(1.0);
00073 this->state_var_old_state_flags_.push_back(true);
00074 this->state_var_output_flags_.push_back(false);
00075
00076
00077 this->num_state_variables_++;
00078 this->state_var_names_.push_back(eqps_string);
00079 this->state_var_layouts_.push_back(dl->qp_scalar);
00080 this->state_var_init_types_.push_back("scalar");
00081 this->state_var_init_values_.push_back(0.0);
00082 this->state_var_old_state_flags_.push_back(true);
00083 this->state_var_output_flags_.push_back(true);
00084
00085
00086 this->num_state_variables_++;
00087 this->state_var_names_.push_back(void_string);
00088 this->state_var_layouts_.push_back(dl->qp_scalar);
00089 this->state_var_init_types_.push_back("scalar");
00090 this->state_var_init_values_.push_back(f0_);
00091 this->state_var_old_state_flags_.push_back(true);
00092 this->state_var_output_flags_.push_back(true);
00093 }
00094
00095 template<typename EvalT, typename Traits>
00096 void GursonModel<EvalT, Traits>::
00097 computeState(typename Traits::EvalData workset,
00098 std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > dep_fields,
00099 std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > eval_fields)
00100 {
00101
00102 PHX::MDField<ScalarT> def_grad = *dep_fields["F"];
00103 PHX::MDField<ScalarT> J = *dep_fields["J"];
00104 PHX::MDField<ScalarT> poissons_ratio = *dep_fields["Poissons Ratio"];
00105 PHX::MDField<ScalarT> elastic_modulus = *dep_fields["Elastic Modulus"];
00106 PHX::MDField<ScalarT> yield_strength = *dep_fields["Yield Strength"];
00107 PHX::MDField<ScalarT> hardening_modulus = *dep_fields["Hardening Modulus"];
00108
00109
00110 std::string cauchy_string = (*field_name_map_)["Cauchy_Stress"];
00111 std::string Fp_string = (*field_name_map_)["Fp"];
00112 std::string eqps_string = (*field_name_map_)["eqps"];
00113 std::string void_string = (*field_name_map_)["Void_Volume"];
00114
00115
00116 PHX::MDField<ScalarT> stress = *eval_fields[cauchy_string];
00117 PHX::MDField<ScalarT> Fp = *eval_fields[Fp_string];
00118 PHX::MDField<ScalarT> eqps = *eval_fields[eqps_string];
00119 PHX::MDField<ScalarT> void_volume = *eval_fields[void_string];
00120
00121
00122 Albany::MDArray Fp_old =
00123 (*workset.stateArrayPtr)[Fp_string + "_old"];
00124 Albany::MDArray eqps_old =
00125 (*workset.stateArrayPtr)[eqps_string + "_old"];
00126 Albany::MDArray void_volume_old =
00127 (*workset.stateArrayPtr)[void_string + "_old"];
00128
00129 Intrepid::Tensor<ScalarT> F(num_dims_), be(num_dims_), logbe(num_dims_);
00130 Intrepid::Tensor<ScalarT> s(num_dims_), sigma(num_dims_), N(num_dims_);
00131 Intrepid::Tensor<ScalarT> A(num_dims_), expA(num_dims_), Fpnew(num_dims_);
00132 Intrepid::Tensor<ScalarT> Fpn(num_dims_), Fpinv(num_dims_), Cpinv(num_dims_);
00133 Intrepid::Tensor<ScalarT> dPhi(num_dims_);
00134 Intrepid::Tensor<ScalarT> I(Intrepid::eye<ScalarT>(num_dims_));
00135
00136 ScalarT kappa, mu, K, Y;
00137 ScalarT p, trlogbeby3, detbe;
00138 ScalarT fvoid, fvoid_star, eq, Phi, dgam, Ybar;
00139
00140
00141 std::vector<ScalarT> X(4);
00142 std::vector<ScalarT> R(4);
00143 std::vector<ScalarT> dRdX(16);
00144
00145 for (std::size_t cell(0); cell < workset.numCells; ++cell) {
00146 for (std::size_t pt(0); pt < num_pts_; ++pt) {
00147 kappa = elastic_modulus(cell, pt)
00148 / (3.0 * (1.0 - 2.0 * poissons_ratio(cell, pt)));
00149 mu = elastic_modulus(cell, pt)
00150 / (2.0 * (1.0 + poissons_ratio(cell, pt)));
00151 K = hardening_modulus(cell, pt);
00152 Y = yield_strength(cell, pt);
00153
00154
00155 F.fill(&def_grad(cell, pt, 0, 0));
00156
00157 for (std::size_t i(0); i < num_dims_; ++i) {
00158 for (std::size_t j(0); j < num_dims_; ++j) {
00159 Fpn(i, j) = static_cast<ScalarT>(Fp_old(cell, pt, i, j));
00160 }
00161 }
00162
00163
00164 Fpinv = Intrepid::inverse(Fpn);
00165 Cpinv = Fpinv * Intrepid::transpose(Fpinv);
00166 be = F * Cpinv * Intrepid::transpose(F);
00167 logbe = Intrepid::log_sym<ScalarT>(be);
00168 trlogbeby3 = Intrepid::trace(logbe) / 3.0;
00169 detbe = Intrepid::det<ScalarT>(be);
00170 s = mu * (logbe - trlogbeby3 * I);
00171 p = 0.5 * kappa * std::log(detbe);
00172 fvoid = void_volume_old(cell, pt);
00173 eq = eqps_old(cell, pt);
00174
00175
00176 Phi = YieldFunction(s, p, fvoid, eq, K, Y, J(cell, pt),
00177 elastic_modulus(cell, pt));
00178
00179 dgam = 0.0;
00180 if (Phi > 0.0) {
00181
00182
00183 X[0] = dgam;
00184 X[1] = p;
00185 X[2] = fvoid;
00186 X[3] = eq;
00187
00188 LocalNonlinearSolver<EvalT, Traits> solver;
00189
00190 int iter = 0;
00191 ScalarT norm_residual0(0.0), norm_residual(0.0), relative_residual(0.0);
00192
00193
00194 while (true) {
00195
00196 ResidualJacobian(X, R, dRdX, p, fvoid, eq, s, mu, kappa, K, Y,
00197 J(cell, pt));
00198
00199 norm_residual = 0.0;
00200 for (int i = 0; i < 4; i++)
00201 norm_residual += R[i] * R[i];
00202
00203 norm_residual = std::sqrt(norm_residual);
00204
00205 if (iter == 0)
00206 norm_residual0 = norm_residual;
00207
00208 if (norm_residual0 != 0)
00209 relative_residual = norm_residual / norm_residual0;
00210 else
00211 relative_residual = norm_residual0;
00212
00213
00214
00215
00216 if (relative_residual < 1.0e-11 || norm_residual < 1.0e-11)
00217 break;
00218
00219 if (iter > 20)
00220 break;
00221
00222
00223 solver.solve(dRdX, X, R);
00224
00225 iter++;
00226 }
00227
00228
00229
00230 solver.computeFadInfo(dRdX, X, R);
00231
00232
00233 dgam = X[0];
00234 p = X[1];
00235 fvoid = X[2];
00236 eq = X[3];
00237
00238
00239 fvoid_star = fvoid;
00240 if ((fvoid > fc_) && (fvoid < ff_)) {
00241 if ((ff_ - fc_) != 0.0) {
00242 fvoid_star = fc_ + (fvoid - fc_) * (1.0 / q1_ - fc_) / (ff_ - fc_);
00243 }
00244 }
00245 else if (fvoid >= ff_) {
00246 fvoid_star = 1.0 / q1_;
00247 if (fvoid_star > 1.0)
00248 fvoid_star = 1.0;
00249 }
00250
00251
00252 s = (1.0 / (1.0 + 2.0 * mu * dgam)) * s;
00253
00254
00255 Ybar = Y + sat_mod_ * (1.0 - std::exp(-sat_exp_ * eq)) + K * eq;
00256
00257
00258 Ybar = Ybar * J(cell, pt);
00259
00260
00261 ScalarT tmp = 1.5 * q2_ * p / Ybar;
00262 dPhi =
00263 s + 1.0 / 3.0 * q1_ * q2_ * Ybar * fvoid_star * std::sinh(tmp) * I;
00264
00265 expA = Intrepid::exp(dgam * dPhi);
00266
00267 for (std::size_t i(0); i < num_dims_; ++i) {
00268 for (std::size_t j(0); j < num_dims_; ++j) {
00269 Fp(cell, pt, i, j) = 0.0;
00270 for (std::size_t k(0); k < num_dims_; ++k) {
00271 Fp(cell, pt, i, j) += expA(i, k) * Fpn(k, j);
00272 }
00273 }
00274 }
00275
00276 eqps(cell, pt) = eq;
00277 void_volume(cell, pt) = fvoid;
00278
00279 }
00280 else {
00281
00282 eqps(cell, pt) = eqps_old(cell, pt);
00283 void_volume(cell, pt) = void_volume_old(cell, pt);
00284
00285 for (std::size_t i(0); i < num_dims_; ++i) {
00286 for (std::size_t j(0); j < num_dims_; ++j) {
00287 Fp(cell, pt, i, j) = Fp_old(cell, pt, i, j);
00288 }
00289 }
00290
00291 }
00292
00293
00294
00295
00296 for (std::size_t i(0); i < num_dims_; ++i) {
00297 for (std::size_t j(0); j < num_dims_; ++j) {
00298 stress(cell, pt, i, j) = s(i, j) / J(cell, pt);
00299 }
00300 stress(cell, pt, i, i) += p / J(cell, pt);
00301 }
00302
00303 }
00304 }
00305
00306 }
00307
00308
00309
00310 template<typename EvalT, typename Traits>
00311 typename EvalT::ScalarT
00312 GursonModel<EvalT, Traits>::YieldFunction(Intrepid::Tensor<ScalarT> const & s,
00313 ScalarT const & p, ScalarT const & fvoid, ScalarT const & eq,
00314 ScalarT const & K, ScalarT const & Y, ScalarT const & jacobian,
00315 ScalarT const & E)
00316 {
00317
00318 ScalarT Ybar = Y + sat_mod_ * (1.0 - std::exp(-sat_exp_ * eq)) + K * eq;
00319
00320
00321 Ybar = Ybar * jacobian;
00322
00323 ScalarT tmp = 1.5 * q2_ * p / Ybar;
00324
00325
00326 ScalarT fvoid_star = fvoid;
00327 if ((fvoid > fc_) && (fvoid < ff_)) {
00328 if ((ff_ - fc_) != 0.0) {
00329 fvoid_star = fc_ + (fvoid - fc_) * (1. / q1_ - fc_) / (ff_ - fc_);
00330 }
00331 }
00332 else if (fvoid >= ff_) {
00333 fvoid_star = 1.0 / q1_;
00334 if (fvoid_star > 1.0)
00335 fvoid_star = 1.0;
00336 }
00337
00338 ScalarT psi = 1.0 + q3_ * fvoid_star * fvoid_star
00339 - 2.0 * q1_ * fvoid_star * std::cosh(tmp);
00340
00341
00342 ScalarT Phi = 0.5 * Intrepid::dotdot(s, s) - psi * Ybar * Ybar / 3.0;
00343
00344
00345
00346
00347
00348
00349
00350 return Phi;
00351 }
00352
00353 template<typename EvalT, typename Traits>
00354 void
00355 GursonModel<EvalT, Traits>::ResidualJacobian(std::vector<ScalarT> & X,
00356 std::vector<ScalarT> & R, std::vector<ScalarT> & dRdX, const ScalarT & p,
00357 const ScalarT & fvoid, const ScalarT & eq, Intrepid::Tensor<ScalarT> & s,
00358 const ScalarT & mu, const ScalarT & kappa, const ScalarT & K,
00359 const ScalarT & Y, const ScalarT & jacobian)
00360 {
00361 ScalarT sq32 = std::sqrt(3.0 / 2.0);
00362 ScalarT sq23 = std::sqrt(2.0 / 3.0);
00363 std::vector<DFadType> Rfad(4);
00364 std::vector<DFadType> Xfad(4);
00365
00366
00367
00368
00369
00370 std::vector<ScalarT> Xval(4);
00371 for (std::size_t i = 0; i < 4; ++i) {
00372 Xval[i] = Sacado::ScalarValue<ScalarT>::eval(X[i]);
00373 Xfad[i] = DFadType(4, i, Xval[i]);
00374 }
00375
00376 DFadType dgam = Xfad[0];
00377 DFadType pFad = Xfad[1];
00378 DFadType fvoidFad = Xfad[2];
00379 DFadType eqFad = Xfad[3];
00380
00381
00382 DFadType fvoidFad_star = fvoidFad;
00383
00384 if ((fvoidFad > fc_) && (fvoidFad < ff_)) {
00385 if ((ff_ - fc_) != 0.0) {
00386 fvoidFad_star = fc_ + (fvoidFad - fc_) * (1. / q1_ - fc_) / (ff_ - fc_);
00387 }
00388 }
00389 else if (fvoidFad >= ff_) {
00390 fvoidFad_star = 1.0 / q1_;
00391 if (fvoidFad_star > 1.0)
00392 fvoidFad_star = 1.0;
00393 }
00394
00395
00396 DFadType Ybar =
00397 Y + sat_mod_ * (1.0 - std::exp(-sat_exp_ * eqFad)) + K * eqFad;
00398
00399
00400 Ybar = Ybar * jacobian;
00401
00402 DFadType tmp = 1.5 * q2_ * pFad / Ybar;
00403
00404 DFadType psi =
00405 1.0 + q3_ * fvoidFad_star * fvoidFad_star
00406 - 2.0 * q1_ * fvoidFad_star * std::cosh(tmp);
00407
00408 DFadType factor = 1.0 / (1.0 + (2.0 * (mu * dgam)));
00409
00410
00411 Intrepid::Tensor<DFadType> sfad(num_dims_);
00412 for (std::size_t i = 0; i < num_dims_; ++i) {
00413 for (std::size_t j = 0; j < num_dims_; ++j) {
00414 sfad(i, j) = factor * s(i, j);
00415 }
00416 }
00417
00418
00419
00420
00421
00422 DFadType omega(0.0), J3(0.0), taue(0.0), smag2, smag;
00423 J3 = Intrepid::det(sfad);
00424 smag2 = Intrepid::dotdot(sfad, sfad);
00425 if (smag2 > 0.0) {
00426 smag = std::sqrt(smag2);
00427 taue = sq32 * smag;
00428 }
00429
00430 if (taue > 0.0)
00431 omega = 1.0
00432 - (27.0 * J3 / 2.0 / taue / taue / taue)
00433 * (27.0 * J3 / 2.0 / taue / taue / taue);
00434
00435 DFadType deq(0.0);
00436 if (smag != 0.0) {
00437 deq = dgam
00438 * (smag2 + q1_ * q2_ * pFad * Ybar * fvoidFad_star * std::sinh(tmp))
00439 / (1.0 - fvoidFad) / Ybar;
00440 }
00441 else {
00442 deq = dgam * (q1_ * q2_ * pFad * Ybar * fvoidFad_star * std::sinh(tmp))
00443 / (1.0 - fvoidFad) / Ybar;
00444 }
00445
00446
00447 DFadType dfn(0.0);
00448 DFadType An(0.0), eratio(0.0);
00449 eratio = -0.5 * (eqFad - eN_) * (eqFad - eN_) / sN_ / sN_;
00450
00451 const double pi = acos(-1.0);
00452 if (pFad >= 0.0) {
00453 An = fN_ / sN_ / (std::sqrt(2.0 * pi)) * std::exp(eratio);
00454 }
00455
00456 dfn = An * deq;
00457
00458
00459
00460 DFadType dfg(0.0);
00461 if (taue > 0.0) {
00462 dfg = dgam * q1_ * q2_ * (1.0 - fvoidFad) * fvoidFad_star * Ybar
00463 * std::sinh(tmp) + sq23 * dgam * kw_ * fvoidFad * omega * smag;
00464 }
00465 else {
00466 dfg = dgam * q1_ * q2_ * (1.0 - fvoidFad)
00467 * fvoidFad_star * Ybar * std::sinh(tmp);
00468 }
00469
00470 DFadType Phi;
00471 Phi = 0.5 * smag2 - psi * Ybar * Ybar / 3.0;
00472
00473
00474 Rfad[0] = Phi;
00475 Rfad[1] = pFad - p
00476 + dgam * q1_ * q2_ * kappa * Ybar * fvoidFad_star * std::sinh(tmp);
00477 Rfad[2] = fvoidFad - fvoid - dfg - dfn;
00478 Rfad[3] = eqFad - eq - deq;
00479
00480
00481 for (int i = 0; i < 4; i++)
00482 R[i] = Rfad[i].val();
00483
00484
00485 for (int i = 0; i < 4; i++)
00486 for (int j = 0; j < 4; j++)
00487 dRdX[i + 4 * j] = Rfad[i].dx(j);
00488
00489 }
00490
00491 }
00492