Real PenetrationAux::computeValue() { const Node * current_node = NULL; if (_nodal) current_node = _current_node; else current_node = _mesh.getQuadratureNode(_current_elem, _current_side, _qp); PenetrationInfo * pinfo = _penetration_locator._penetration_info[current_node->id()]; Real retVal(NotPenetrated); if (pinfo) switch (_quantity) { case PA_DISTANCE: retVal = pinfo->_distance; break; case PA_TANG_DISTANCE: retVal = pinfo->_tangential_distance; break; case PA_NORMAL_X: retVal = pinfo->_normal(0); break; case PA_NORMAL_Y: retVal = pinfo->_normal(1); break; case PA_NORMAL_Z: retVal = pinfo->_normal(2); break; case PA_CLOSEST_POINT_X: retVal = pinfo->_closest_point(0); break; case PA_CLOSEST_POINT_Y: retVal = pinfo->_closest_point(1); break; case PA_CLOSEST_POINT_Z: retVal = pinfo->_closest_point(2); break; case PA_ELEM_ID: retVal = static_cast<Real>(pinfo->_elem->id() + 1); break; case PA_SIDE: retVal = static_cast<Real>(pinfo->_side_num); break; case PA_INCREMENTAL_SLIP_MAG: retVal = pinfo->isCaptured() ? pinfo->_incremental_slip.norm() : 0; break; case PA_INCREMENTAL_SLIP_X: retVal = pinfo->isCaptured() ? pinfo->_incremental_slip(0) : 0; break; case PA_INCREMENTAL_SLIP_Y: retVal = pinfo->isCaptured() ? pinfo->_incremental_slip(1) : 0; break; case PA_INCREMENTAL_SLIP_Z: retVal = pinfo->isCaptured() ? pinfo->_incremental_slip(2) : 0; break; case PA_ACCUMULATED_SLIP: retVal = pinfo->_accumulated_slip; break; case PA_FORCE_X: retVal = pinfo->_contact_force(0); break; case PA_FORCE_Y: retVal = pinfo->_contact_force(1); break; case PA_FORCE_Z: retVal = pinfo->_contact_force(2); break; case PA_NORMAL_FORCE_MAG: retVal = -pinfo->_contact_force * pinfo->_normal; break; case PA_NORMAL_FORCE_X: retVal = (pinfo->_contact_force * pinfo->_normal) * pinfo->_normal(0); break; case PA_NORMAL_FORCE_Y: retVal = (pinfo->_contact_force * pinfo->_normal) * pinfo->_normal(1); break; case PA_NORMAL_FORCE_Z: retVal = (pinfo->_contact_force * pinfo->_normal) * pinfo->_normal(2); break; case PA_TANGENTIAL_FORCE_MAG: { RealVectorValue contact_force_normal((pinfo->_contact_force * pinfo->_normal) * pinfo->_normal); RealVectorValue contact_force_tangential(pinfo->_contact_force - contact_force_normal); retVal = contact_force_tangential.norm(); break; } case PA_TANGENTIAL_FORCE_X: retVal = pinfo->_contact_force(0) - (pinfo->_contact_force * pinfo->_normal) * pinfo->_normal(0); break; case PA_TANGENTIAL_FORCE_Y: retVal = pinfo->_contact_force(1) - (pinfo->_contact_force * pinfo->_normal) * pinfo->_normal(1); break; case PA_TANGENTIAL_FORCE_Z: retVal = pinfo->_contact_force(2) - (pinfo->_contact_force * pinfo->_normal) * pinfo->_normal(2); break; case PA_FRICTIONAL_ENERGY: retVal = pinfo->_frictional_energy; break; case PA_LAGRANGE_MULTIPLIER: retVal = pinfo->_lagrange_multiplier; break; case PA_MECH_STATUS: retVal = pinfo->_mech_status; break; default: mooseError("Unknown PA_ENUM"); } // switch return retVal; }
void ContactMaster::computeContactForce(PenetrationInfo * pinfo) { std::map<unsigned int, Real> & lagrange_multiplier = _penetration_locator._lagrange_multiplier; const Node * node = pinfo->_node; RealVectorValue res_vec; // Build up residual vector for (unsigned int i=0; i<_mesh_dimension; ++i) { long int dof_number = node->dof_number(0, _vars(i), 0); res_vec(i) = _residual_copy(dof_number); } const Real area = nodalArea(*pinfo); RealVectorValue distance_vec(_mesh.node(node->id()) - pinfo->_closest_point); RealVectorValue pen_force(_penalty * distance_vec); if (_normalize_penalty) pen_force *= area; RealVectorValue tan_residual(0,0,0); if (_model == CM_FRICTIONLESS) { switch (_formulation) { case CF_DEFAULT: pinfo->_contact_force = -pinfo->_normal * (pinfo->_normal * res_vec); break; case CF_PENALTY: pinfo->_contact_force = pinfo->_normal * (pinfo->_normal * pen_force); break; case CF_AUGMENTED_LAGRANGE: pinfo->_contact_force = (pinfo->_normal * (pinfo->_normal * //( pen_force + (lagrange_multiplier[node->id()]/distance_vec.size())*distance_vec))); ( pen_force + lagrange_multiplier[node->id()] * pinfo->_normal))); break; default: mooseError("Invalid contact formulation"); break; } pinfo->_mech_status=PenetrationInfo::MS_SLIPPING; } else if (_model == CM_COULOMB && _formulation == CF_PENALTY) { distance_vec = pinfo->_incremental_slip + (pinfo->_normal * (_mesh.node(node->id()) - pinfo->_closest_point)) * pinfo->_normal; pen_force = _penalty * distance_vec; if (_normalize_penalty) pen_force *= area; // Frictional capacity // const Real capacity( _friction_coefficient * (pen_force * pinfo->_normal < 0 ? -pen_force * pinfo->_normal : 0) ); const Real capacity( _friction_coefficient * (res_vec * pinfo->_normal > 0 ? res_vec * pinfo->_normal : 0) ); // Elastic predictor pinfo->_contact_force = pen_force + (pinfo->_contact_force_old - pinfo->_normal*(pinfo->_normal*pinfo->_contact_force_old)); RealVectorValue contact_force_normal( (pinfo->_contact_force*pinfo->_normal) * pinfo->_normal ); RealVectorValue contact_force_tangential( pinfo->_contact_force - contact_force_normal ); // Tangential magnitude of elastic predictor const Real tan_mag( contact_force_tangential.size() ); if ( tan_mag > capacity ) { pinfo->_contact_force = contact_force_normal + capacity * contact_force_tangential / tan_mag; pinfo->_mech_status=PenetrationInfo::MS_SLIPPING; } else { pinfo->_mech_status=PenetrationInfo::MS_STICKING; } } else if (_model == CM_GLUED || (_model == CM_COULOMB && _formulation == CF_DEFAULT)) { switch (_formulation) { case CF_DEFAULT: pinfo->_contact_force = -res_vec; break; case CF_PENALTY: pinfo->_contact_force = pen_force; break; case CF_AUGMENTED_LAGRANGE: pinfo->_contact_force = pen_force + lagrange_multiplier[node->id()]*distance_vec/distance_vec.size(); break; default: mooseError("Invalid contact formulation"); break; } pinfo->_mech_status=PenetrationInfo::MS_STICKING; } else { mooseError("Invalid or unavailable contact model"); } }
void MechanicalContactConstraint::computeContactForce(PenetrationInfo * pinfo) { const Node * node = pinfo->_node; RealVectorValue res_vec; // Build up residual vector for (unsigned int i=0; i<_mesh_dimension; ++i) { dof_id_type dof_number = node->dof_number(0, _vars(i), 0); res_vec(i) = _residual_copy(dof_number); } RealVectorValue distance_vec(_mesh.node(node->id()) - pinfo->_closest_point); const Real penalty = getPenalty(*pinfo); RealVectorValue pen_force(penalty * distance_vec); switch (_model) { case CM_FRICTIONLESS: switch (_formulation) { case CF_KINEMATIC: pinfo->_contact_force = -pinfo->_normal * (pinfo->_normal * res_vec); break; case CF_PENALTY: pinfo->_contact_force = pinfo->_normal * (pinfo->_normal * pen_force); break; case CF_AUGMENTED_LAGRANGE: pinfo->_contact_force = (pinfo->_normal * (pinfo->_normal * ( pen_force + pinfo->_lagrange_multiplier * pinfo->_normal))); break; default: mooseError("Invalid contact formulation"); break; } pinfo->_mech_status=PenetrationInfo::MS_SLIPPING; break; case CM_COULOMB: switch (_formulation) { case CF_KINEMATIC: pinfo->_contact_force = -res_vec; break; case CF_PENALTY: { distance_vec = pinfo->_incremental_slip + (pinfo->_normal * (_mesh.node(node->id()) - pinfo->_closest_point)) * pinfo->_normal; pen_force = penalty * distance_vec; // Frictional capacity // const Real capacity( _friction_coefficient * (pen_force * pinfo->_normal < 0 ? -pen_force * pinfo->_normal : 0) ); const Real capacity( _friction_coefficient * (res_vec * pinfo->_normal > 0 ? res_vec * pinfo->_normal : 0) ); // Elastic predictor pinfo->_contact_force = pen_force + (pinfo->_contact_force_old - pinfo->_normal*(pinfo->_normal*pinfo->_contact_force_old)); RealVectorValue contact_force_normal( (pinfo->_contact_force*pinfo->_normal) * pinfo->_normal ); RealVectorValue contact_force_tangential( pinfo->_contact_force - contact_force_normal ); // Tangential magnitude of elastic predictor const Real tan_mag( contact_force_tangential.size() ); if ( tan_mag > capacity ) { pinfo->_contact_force = contact_force_normal + capacity * contact_force_tangential / tan_mag; pinfo->_mech_status=PenetrationInfo::MS_SLIPPING; } else pinfo->_mech_status=PenetrationInfo::MS_STICKING; break; } case CF_AUGMENTED_LAGRANGE: pinfo->_contact_force = pen_force + pinfo->_lagrange_multiplier*distance_vec/distance_vec.size(); break; default: mooseError("Invalid contact formulation"); break; } break; case CM_GLUED: switch (_formulation) { case CF_KINEMATIC: pinfo->_contact_force = -res_vec; break; case CF_PENALTY: pinfo->_contact_force = pen_force; break; case CF_AUGMENTED_LAGRANGE: pinfo->_contact_force = pen_force + pinfo->_lagrange_multiplier*distance_vec/distance_vec.size(); break; default: mooseError("Invalid contact formulation"); break; } pinfo->_mech_status=PenetrationInfo::MS_STICKING; break; default: mooseError("Invalid or unavailable contact model"); break; } }