//======================================================================== //======================================================================== // // NAME: double bessel_2(int n, double arg) // // DESC: Calculates the Bessel function of the second kind (Yn). // // INPUT: // int n:: Order of Bessel function // double arg: Bessel function argument // // OUTPUT: // Yn:: Bessel function of the second kind of order n // // NOTES: 1) The Bessel function of the second kind is defined as: // Yn(x) = (2.0*Jn(x)/pi)*(ln(x/2) + gamma_e) // + sum{m=0, m=inf}[(((-1)^(m-1))*(hs(m) + // hs(m+n)) + x^(2*m))/(((2.0^(2*m+n))*m!*(m+n)!) // - sum{m=0, m=(n-1)}[((n-m-1)!*x^(2*m))/((2.0^(2*m-n))*m!) // // 2) Y-n = ((-1)^n)*Yn // //======================================================================== //======================================================================== double bessel_2(int n, double arg) { // get a tolerance for the "infinite" sum double tol = 100.0*depsilon(); int m, mmax; double Jn, Yn, sum1, sum1_prev, sum2; sum1 = sum1_prev = sum2 = 0.0; // Make sure that we calculate the positive order Bessel function int k = abs(n); // GET THE BESSEL FUNCTION OF THE FIRST KIND Jn = bessel_1(k, arg); // !!! my concern with this do loop is that if a certain term contributes 0 then the // loop may inappropriately exit mmax = static_cast<int>(2.0*arg) + 1; m = 0; // GET TERM SUM 1 do { sum1_prev = sum1; sum1 += pow((0.0-1.0), (m-1))*(h_s(m) + h_s(m+k))*pow(arg, (2*m))/(pow(2.0, (2*m+k))*dfactorial(m)*dfactorial(m+k)); m++; } while((fabs(sum1 - sum1_prev) > tol) || (m < mmax)); sum1 = pow(arg, k)*sum1/(1.0*PI); // GET TERM SUM 2 for(m = 0; m <= (k-1); ++m) { sum2 = sum2 + dfactorial(k-m-1)*pow(arg, (2*m))/(pow(2.0, (2*m-k))*dfactorial(m)); } sum2 = pow(arg, (0-k))*sum2/(1.0*PI); // NOW GET Yn Yn = (2.0*Jn/(1.0*PI))*(log(arg/2.0) + EULERC) + sum1 - sum2; // NOW WE WILL MAKE USE OF THE BESSEL FUNCTION RELATION FOR NEGATIVE n: // Y(n-1) = (-1)^n*Yn if(n < 0) { Yn = pow((0.0-1.0), k)*Yn; } return Yn; }
//======================================================================== //======================================================================== // // NAME: complex<double> bessel_1_complex(int n, complex<double> arg) // // DESC: Calculates the Bessel function of the first kind (Jn) for // complex arguments. // // INPUT: // int n:: Order of Bessel function // complex<double> arg: Bessel function argument // // OUTPUT: // Jn:: Bessel function of the first kind of order n // // NOTES: 1) The Bessel function of the first kind is defined as: // Jn(x) = ((x/2.0)^n)*sum{m=0, m=inf}[(((-x^2)/4.0)^m)/(m!*(n+m)!) // // 2) J-n = ((-1)^n)*Jn // // 3) The infinite series representing the Bessel function // converges for all arguments given that enough terms are taken: // k >> |arg| // //======================================================================== //======================================================================== complex<double> bessel_1_complex(int n, complex<double> arg) { complex<double> Jn(0.0, 0.0); complex<double> Jn_series(0.0, 0.0); complex<double> Jn_seriesm1(0.0, 0.0); double tolerance = 100.0*depsilon(); // Make sure that we calculate the positive order Bessel function int m = abs(n); // !!! I don't know if this will be large enough int k_max_r, k_max_im; k_max_r = static_cast<int>(real(arg)); k_max_im = static_cast<int>(imag(arg)); /* int k_max = 10*imax(k_max_r, k_max_im); for(int k = 0; k <= k_max; k++) { Jn_series = Jn_series + pow((0.0 - arg*arg/4.0), k)/(dfactorial(k)*dfactorial(k+m)); } */ // !!! my concern with this do loop is that if a certain term contributes 0 then the // loop may inappropriately exit int k_max2 = static_cast<int>(2.0*static_cast<double>(imax(k_max_r, k_max_im))) + 1; int k = 0; do { Jn_seriesm1 = Jn_series; Jn_series = Jn_series + pow((0.0 - arg*arg/4.0), k)/(dfactorial(k)*dfactorial(k+m)); k = k + 1; } while ((sqrt(real(Jn_series*conj(Jn_series) - Jn_seriesm1*conj(Jn_seriesm1))) > tolerance) || (k < k_max2)); Jn = pow((arg/2.0), m)*Jn_series; // NOW WE WILL MAKE USE OF THE BESSEL FUNCTION RELATION FOR NEGATIVE n: // J(n-1) = (-1)^n*Jn if(n < 0) { Jn = pow((0.0 - 1.0), m)*Jn; } return Jn; }
//======================================================================== //======================================================================== // // NAME: double bessel_1(int n, double arg) // // DESC: Calculates the Bessel function of the first kind (Jn) of order n. // // INPUT: // int n:: Order of Bessel function // double arg: Bessel function argument // // OUTPUT: // Jn:: Bessel function of the first kind of order n // // NOTES: 1) The Bessel function of the first kind is defined as: // Jn(x) = ((x/2.0)^n)*sum{m=0, m=inf}[(((-x^2)/4.0)^m)/(m!*(n+m)!) // // 2) J-n = ((-1)^n)*Jn // // 3) The infinite series representing the Bessel function // converges for all arguments given that enough terms are taken: // k >> |arg| // //======================================================================== //======================================================================== double bessel_1(int n, double arg) { double tolerance = 100.0*depsilon(); double Jn, Jn1; Jn = Jn1 = 0.0; // make sure we calculate the positive order Bessel function int m = abs(n); // !!! my concern with this do loop is that if a certain term contributes 0 then the // loop may inappropriately exit, and I don't know if this will be large enough int k_max2 = int(2.0*arg + 1.0); int k = 0; do { Jn1 = Jn; Jn += pow(-arg*arg/4.0, k)/(dfactorial(k)*dfactorial(k+m)); k++; } while ((fabs(Jn - Jn1) > tolerance) || (k < k_max2)); Jn *= pow((arg/2.0), m); // NOW WE WILL MAKE USE OF THE BESSEL FUNCTION RELATION FOR NEGATIVE n // J(n-1) = (-1)^n*Jn if(n < 0) { Jn *= pow(-1.0, m); } return Jn; }
//======================================================================== //======================================================================== // // NAME: complex<double> bessel_2_complex(int n, complex<double> arg) // // DESC: Calculates the Bessel function of the second kind (Yn) for // complex arguments. // // INPUT: // int n:: Order of Bessel function // complex<double> arg: Bessel function argument // // OUTPUT: // complex<double> Yn:: Bessel function of the second kind of // order n // // NOTES: 1) The Bessel function of the second kind is defined as: // Yn(x) = (2.0*Jn(x)/pi)*(ln(x/2) + gamma_e) // + sum{m=0, m=inf}[(((-1)^(m-1))*(hs(m) + // hs(m+n)) + x^(2*m))/(((2.0^(2*m+n))*m!*(m+n)!) // - sum{m=0, m=(n-1)}[((n-m-1)!*x^(2*m))/((2.0^(2*m-n))*m!) // // 2) Y-n = ((-1)^n)*Yn // //======================================================================== //======================================================================== complex<double> bessel_2_complex(int n, complex<double> arg) { complex<double> Yn(0.0, 0.0); // get a tolerance for the "infinite" sum double tolerance = 100.0*depsilon(); // make sure we calculate the positive order Bessel function int k = abs(n); // GET BESSEL FUNCTIONS OF THE FIRST KIND complex<double> Jn = bessel_1_complex(k, arg); // GET TERM SUM 1 complex<double> sum1(0.0,0.0), sum11(0.0,0.0); /* for(int i = 0; i <= 20; i++) { sum1 = sum1 + pow((0.0-1.0), (i-1))*(h_s(i) + h_s(i+k))*pow(arg, (2*i))/(pow(2.0, (2*i+k))*dfactorial(i)*dfactorial(i+k)); } */ // !!! my concern with this do loop is that if a certain term contributes 0 then the // loop may inappropriately exit int mm_max = static_cast<int>(2.0*real(arg)) + 1; int mm = 0; do { sum11 = sum1; sum1 = sum1 + pow((0.0-1.0), (mm-1))*(h_s(mm) + h_s(mm+k))*pow(arg, (2*mm))/(pow(2.0, (2*mm+k))*dfactorial(mm)*dfactorial(mm+k)); mm = mm + 1; } while((fabs(real(sqrt(sum1*conj(sum1) - sum11*conj(sum11)))) > tolerance) || (mm < mm_max)); sum1 = pow(arg, k)*sum1/(1.0*PI); // GET TERM SUM 2 complex<double> sum2(0.0,0.0); for(int m = 0; m <= (k-1); m++) { sum2 = sum2 + dfactorial((k-m-1))*pow(arg, (2*m))/(pow(2.0, (2*m-k))*dfactorial(m)); } sum2 = pow(arg, (0-k))*sum2/(1.0*PI); // NOW GET Yn Yn = (2.0*Jn/(1.0*PI))*(log(arg/2.0) + EULERC) + sum1 - sum2; // NOW WE WILL MAKE USE OF THE BESSEL FUNCTION RELATION FOR NEGATIVE n: // Y(n-1) = (-1)^n*Yn if(n < 0) { Yn = pow((0.0-1.0), k)*Yn; } return Yn; }
void DruckerPragerModel<EvalT, Traits>:: computeState(typename Traits::EvalData workset, std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > dep_fields, std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > eval_fields) { // extract dependent MDFields PHX::MDField<ScalarT> strain = *dep_fields["Strain"]; PHX::MDField<ScalarT> poissons_ratio = *dep_fields["Poissons Ratio"]; PHX::MDField<ScalarT> elastic_modulus = *dep_fields["Elastic Modulus"]; // retrieve appropriate field name strings std::string cauchy_string = (*field_name_map_)["Cauchy_Stress"]; std::string strain_string = (*field_name_map_)["Strain"]; std::string eqps_string = (*field_name_map_)["eqps"]; std::string friction_string = (*field_name_map_)["Friction_Parameter"]; // extract evaluated MDFields PHX::MDField<ScalarT> stress = *eval_fields[cauchy_string]; PHX::MDField<ScalarT> eqps = *eval_fields[eqps_string]; PHX::MDField<ScalarT> friction = *eval_fields[friction_string]; PHX::MDField<ScalarT> tangent = *eval_fields["Material Tangent"]; // get State Variables Albany::MDArray strainold = (*workset.stateArrayPtr)[strain_string + "_old"]; Albany::MDArray stressold = (*workset.stateArrayPtr)[cauchy_string + "_old"]; Albany::MDArray eqpsold = (*workset.stateArrayPtr)[eqps_string + "_old"]; Albany::MDArray frictionold = (*workset.stateArrayPtr)[friction_string + "_old"]; Intrepid::Tensor<ScalarT> id(Intrepid::eye<ScalarT>(num_dims_)); Intrepid::Tensor4<ScalarT> id1(Intrepid::identity_1<ScalarT>(num_dims_)); Intrepid::Tensor4<ScalarT> id2(Intrepid::identity_2<ScalarT>(num_dims_)); Intrepid::Tensor4<ScalarT> id3(Intrepid::identity_3<ScalarT>(num_dims_)); Intrepid::Tensor4<ScalarT> Celastic (num_dims_); Intrepid::Tensor<ScalarT> sigma(num_dims_), sigmaN(num_dims_), s(num_dims_); Intrepid::Tensor<ScalarT> epsilon(num_dims_), epsilonN(num_dims_); Intrepid::Tensor<ScalarT> depsilon(num_dims_); Intrepid::Tensor<ScalarT> nhat(num_dims_); ScalarT lambda, mu, kappa; ScalarT alpha, alphaN; ScalarT p, q, ptr, qtr; ScalarT eq, eqN, deq; ScalarT snorm; ScalarT Phi; //local unknowns and residual vectors std::vector<ScalarT> X(4); std::vector<ScalarT> R(4); std::vector<ScalarT> dRdX(16); for (std::size_t cell(0); cell < workset.numCells; ++cell) { for (std::size_t pt(0); pt < num_pts_; ++pt) { lambda = ( elastic_modulus(cell,pt) * poissons_ratio(cell,pt) ) / ( ( 1 + poissons_ratio(cell,pt) ) * ( 1 - 2 * poissons_ratio(cell,pt) ) ); mu = elastic_modulus(cell,pt) / ( 2 * ( 1 + poissons_ratio(cell,pt) ) ); kappa = lambda + 2.0 * mu / 3.0; // 4-th order elasticity tensor Celastic = lambda * id3 + mu * (id1 + id2); // previous state (the fill doesn't work for state virable) //sigmaN.fill( &stressold(cell,pt,0,0) ); //epsilonN.fill( &strainold(cell,pt,0,0) ); for (std::size_t i(0); i < num_dims_; ++i) { for (std::size_t j(0); j < num_dims_; ++j) { sigmaN(i, j) = stressold(cell, pt, i, j); epsilonN(i, j) = strainold(cell, pt, i, j); //epsilon(i,j) = strain(cell,pt,i,j); } } epsilon.fill( &strain(cell,pt,0,0) ); depsilon = epsilon - epsilonN; alphaN = frictionold(cell,pt); eqN = eqpsold(cell,pt); // trial state sigma = sigmaN + Intrepid::dotdot(Celastic,depsilon); ptr = Intrepid::trace(sigma) / 3.0; s = sigma - ptr * id; snorm = Intrepid::dotdot(s,s); if(snorm > 0) snorm = std::sqrt(snorm); qtr = sqrt(3.0/2.0) * snorm; // unit deviatoric tensor if (snorm > 0){ nhat = s / snorm; } else{ nhat = id; } // check yielding Phi = qtr + alphaN * ptr - Cf_; alpha = alphaN; p = ptr; q = qtr; deq = 0.0; if(Phi > 1.0e-12) {// plastic yielding // initialize local unknown vector X[0] = ptr; X[1] = qtr; X[2] = alpha; X[3] = deq; LocalNonlinearSolver<EvalT, Traits> solver; int iter = 0; ScalarT norm_residual0(0.0), norm_residual(0.0), relative_residual(0.0); // local N-R loop while (true) { ResidualJacobian(X, R, dRdX, ptr, qtr, eqN, mu, kappa); norm_residual = 0.0; for (int i = 0; i < 4; i++) norm_residual += R[i] * R[i]; norm_residual = std::sqrt(norm_residual); if (iter == 0) norm_residual0 = norm_residual; if (norm_residual0 != 0) relative_residual = norm_residual / norm_residual0; else relative_residual = norm_residual0; //std::cout << iter << " " //<< Sacado::ScalarValue<ScalarT>::eval(norm_residual) //<< " " << Sacado::ScalarValue<ScalarT>::eval(relative_residual) //<< std::endl; if (relative_residual < 1.0e-11 || norm_residual < 1.0e-11) break; if (iter > 20) break; // call local nonlinear solver solver.solve(dRdX, X, R); iter++; } // end of local N-R loop // compute sensitivity information w.r.t. system parameters // and pack the sensitivity back to X solver.computeFadInfo(dRdX, X, R); // update p = X[0]; q = X[1]; alpha = X[2]; deq = X[3]; }//end plastic yielding eq = eqN + deq; s = sqrt(2.0/3.0) * q * nhat; sigma = s + p * id; eqps(cell, pt) = eq; friction(cell,pt) = alpha; for (std::size_t i(0); i < num_dims_; ++i) { for (std::size_t j(0); j < num_dims_; ++j) { stress(cell,pt,i,j) = sigma(i, j); } } }// end loop over pt } // end loop over cell }