void MultiPlasticityLinearSystem::nrStep(const RankTwoTensor & stress, const std::vector<Real> & intnl_old, const std::vector<Real> & intnl, const std::vector<Real> & pm, const RankFourTensor & E_inv, const RankTwoTensor & delta_dp, RankTwoTensor & dstress, std::vector<Real> & dpm, std::vector<Real> & dintnl, const std::vector<bool> & active, std::vector<bool> & deactivated_due_to_ld) { // Calculate RHS and Jacobian std::vector<Real> rhs; calculateRHS(stress, intnl_old, intnl, pm, delta_dp, rhs, active, true, deactivated_due_to_ld); std::vector<std::vector<Real> > jac; calculateJacobian(stress, intnl, pm, E_inv, active, deactivated_due_to_ld, jac); // prepare for LAPACKgesv_ routine provided by PETSc int system_size = rhs.size(); std::vector<double> a(system_size*system_size); // Fill in the a "matrix" by going down columns unsigned ind = 0; for (int col = 0 ; col < system_size ; ++col) for (int row = 0 ; row < system_size ; ++row) a[ind++] = jac[row][col]; int nrhs = 1; std::vector<int> ipiv(system_size); int info; LAPACKgesv_(&system_size, &nrhs, &a[0], &system_size, &ipiv[0], &rhs[0], &system_size, &info); if (info != 0) mooseError("In solving the linear system in a Newton-Raphson process, the PETSC LAPACK gsev routine returned with error code " << info); // Extract the results back to dstress, dpm and dintnl std::vector<bool> active_not_deact(_num_surfaces); for (unsigned surface = 0 ; surface < _num_surfaces ; ++surface) active_not_deact[surface] = (active[surface] && !deactivated_due_to_ld[surface]); unsigned int dim = 3; ind = 0; for (unsigned i = 0 ; i < dim ; ++i) for (unsigned j = 0 ; j <= i ; ++j) dstress(i, j) = dstress(j, i) = rhs[ind++]; dpm.assign(_num_surfaces, 0); for (unsigned surface = 0 ; surface < _num_surfaces ; ++surface) if (active_not_deact[surface]) dpm[surface] = rhs[ind++]; dintnl.assign(_num_models, 0); for (unsigned model = 0 ; model < _num_models ; ++model) if (anyActiveSurfaces(model, active_not_deact)) dintnl[model] = rhs[ind++]; mooseAssert(static_cast<int>(ind) == system_size, "Incorrect extracting of changes from NR solution in nrStep"); }
void MultiPlasticityLinearSystem::calculateRHS(const RankTwoTensor & stress, const std::vector<Real> & intnl_old, const std::vector<Real> & intnl, const std::vector<Real> & pm, const RankTwoTensor & delta_dp, std::vector<Real> & rhs, const std::vector<bool> & active, bool eliminate_ld, std::vector<bool> & deactivated_due_to_ld) { // see comments at the start of .h file mooseAssert(intnl_old.size() == _num_models, "Size of intnl_old is " << intnl_old.size() << " which is incorrect in calculateRHS"); mooseAssert(intnl.size() == _num_models, "Size of intnl is " << intnl.size() << " which is incorrect in calculateRHS"); mooseAssert(pm.size() == _num_surfaces, "Size of pm is " << pm.size() << " which is incorrect in calculateRHS"); mooseAssert(active.size() == _num_surfaces, "Size of active is " << active.size() << " which is incorrect in calculateRHS"); std::vector<Real> f; // the yield functions RankTwoTensor epp; // the plastic-strain constraint ("direction constraint") std::vector<Real> ic; // the "internal constraints" std::vector<RankTwoTensor> r; calculateConstraints(stress, intnl_old, intnl, pm, delta_dp, f, r, epp, ic, active); if (eliminate_ld) eliminateLinearDependence(stress, intnl, f, r, active, deactivated_due_to_ld); else deactivated_due_to_ld.assign(_num_surfaces, false); std::vector<bool> active_not_deact(_num_surfaces); for (unsigned surface = 0 ; surface < _num_surfaces ; ++surface) active_not_deact[surface] = (active[surface] && !deactivated_due_to_ld[surface]); unsigned num_active_f = 0; for (unsigned surface = 0 ; surface < _num_surfaces ; ++surface) if (active_not_deact[surface]) num_active_f++; unsigned num_active_ic = 0; for (unsigned model = 0 ; model < _num_models ; ++model) if (anyActiveSurfaces(model, active_not_deact)) num_active_ic++; unsigned int dim = 3; unsigned int system_size = 6 + num_active_f + num_active_ic; // "6" comes from symmeterizing epp, num_active_f comes from "f", num_active_f comes from "ic" rhs.resize(system_size); unsigned ind = 0; for (unsigned i = 0 ; i < dim ; ++i) for (unsigned j = 0 ; j <= i ; ++j) rhs[ind++] = -epp(i, j); unsigned active_surface = 0; for (unsigned surface = 0 ; surface < _num_surfaces ; ++surface) if (active[surface]) { if (!deactivated_due_to_ld[surface]) rhs[ind++] = -f[active_surface]; active_surface++; } unsigned active_model = 0; for (unsigned model = 0 ; model < _num_models ; ++model) if (anyActiveSurfaces(model, active)) { if (anyActiveSurfaces(model, active_not_deact)) rhs[ind++] = -ic[active_model]; active_model++; } mooseAssert(ind == system_size, "Incorrect filling of the rhs in calculateRHS"); }
void MultiPlasticityDebugger::checkSolution(const RankFourTensor & E_inv) { Moose::err << "\n\n+++++++++++++++++++++\nChecking the Solution\n"; Moose::err << "(Ie, checking Ax = b)\n+++++++++++++++++++++\n"; outputAndCheckDebugParameters(); std::vector<bool> act; act.assign(_num_surfaces, true); std::vector<bool> deactivated_due_to_ld; deactivated_due_to_ld.assign(_num_surfaces, false); RankTwoTensor delta_dp = -E_inv * _fspb_debug_stress; std::vector<Real> orig_rhs; calculateRHS(_fspb_debug_stress, _fspb_debug_intnl, _fspb_debug_intnl, _fspb_debug_pm, delta_dp, orig_rhs, act, true, deactivated_due_to_ld); Moose::err << "\nb = "; for (unsigned i = 0; i < orig_rhs.size(); ++i) Moose::err << orig_rhs[i] << " "; Moose::err << "\n\n"; std::vector<std::vector<Real>> jac_coded; calculateJacobian(_fspb_debug_stress, _fspb_debug_intnl, _fspb_debug_pm, E_inv, act, deactivated_due_to_ld, jac_coded); Moose::err << "Before checking Ax=b is correct, check that the Jacobians given below are equal.\n"; Moose::err << "The hand-coded Jacobian is used in calculating the solution 'x', given 'b' above.\n"; Moose::err << "Note that this only includes degrees of freedom that aren't deactivated due to " "linear dependence.\n"; Moose::err << "Hand-coded Jacobian:\n"; for (unsigned row = 0; row < jac_coded.size(); ++row) { for (unsigned col = 0; col < jac_coded.size(); ++col) Moose::err << jac_coded[row][col] << " "; Moose::err << "\n"; } deactivated_due_to_ld.assign(_num_surfaces, false); // this potentially gets changed by nrStep, below RankTwoTensor dstress; std::vector<Real> dpm; std::vector<Real> dintnl; nrStep(_fspb_debug_stress, _fspb_debug_intnl, _fspb_debug_intnl, _fspb_debug_pm, E_inv, delta_dp, dstress, dpm, dintnl, act, deactivated_due_to_ld); std::vector<bool> active_not_deact(_num_surfaces); for (unsigned surface = 0; surface < _num_surfaces; ++surface) active_not_deact[surface] = !deactivated_due_to_ld[surface]; std::vector<Real> x; x.assign(orig_rhs.size(), 0); unsigned ind = 0; for (unsigned i = 0; i < 3; ++i) for (unsigned j = 0; j <= i; ++j) x[ind++] = dstress(i, j); for (unsigned surface = 0; surface < _num_surfaces; ++surface) if (active_not_deact[surface]) x[ind++] = dpm[surface]; for (unsigned model = 0; model < _num_models; ++model) if (anyActiveSurfaces(model, active_not_deact)) x[ind++] = dintnl[model]; mooseAssert(ind == orig_rhs.size(), "Incorrect extracting of changes from NR solution in the " "finite-difference checking of nrStep"); Moose::err << "\nThis yields x ="; for (unsigned i = 0; i < orig_rhs.size(); ++i) Moose::err << x[i] << " "; Moose::err << "\n"; std::vector<std::vector<Real>> jac_fd; fdJacobian(_fspb_debug_stress, _fspb_debug_intnl, _fspb_debug_intnl, _fspb_debug_pm, delta_dp, E_inv, true, jac_fd); Moose::err << "\nThe finite-difference Jacobian is used to multiply by this 'x',\n"; Moose::err << "in order to check that the solution is correct\n"; Moose::err << "Finite-difference Jacobian:\n"; for (unsigned row = 0; row < jac_fd.size(); ++row) { for (unsigned col = 0; col < jac_fd.size(); ++col) Moose::err << jac_fd[row][col] << " "; Moose::err << "\n"; } Real L2_numer = 0; Real L2_denom = 0; for (unsigned row = 0; row < jac_coded.size(); ++row) for (unsigned col = 0; col < jac_coded.size(); ++col) { L2_numer += Utility::pow<2>(jac_coded[row][col] - jac_fd[row][col]); L2_denom += Utility::pow<2>(jac_coded[row][col] + jac_fd[row][col]); } Moose::err << "Relative L2norm of the hand-coded and finite-difference Jacobian is " << std::sqrt(L2_numer / L2_denom) / 0.5 << "\n"; std::vector<Real> fd_times_x; fd_times_x.assign(orig_rhs.size(), 0); for (unsigned row = 0; row < orig_rhs.size(); ++row) for (unsigned col = 0; col < orig_rhs.size(); ++col) fd_times_x[row] += jac_fd[row][col] * x[col]; Moose::err << "\n(Finite-difference Jacobian)*x =\n"; for (unsigned i = 0; i < orig_rhs.size(); ++i) Moose::err << fd_times_x[i] << " "; Moose::err << "\n"; Moose::err << "Recall that b = \n"; for (unsigned i = 0; i < orig_rhs.size(); ++i) Moose::err << orig_rhs[i] << " "; Moose::err << "\n"; L2_numer = 0; L2_denom = 0; for (unsigned i = 0; i < orig_rhs.size(); ++i) { L2_numer += Utility::pow<2>(orig_rhs[i] - fd_times_x[i]); L2_denom += Utility::pow<2>(orig_rhs[i] + fd_times_x[i]); } Moose::err << "\nRelative L2norm of these is " << std::sqrt(L2_numer / L2_denom) / 0.5 << "\n"; }