Point contact
Click here to return back to Solid Mechanics
Contents
[hide]Point contact on a 2D half-plane
A starting point to solve problems in contact mechanics is to understand the effect of a point-load applied to a homogeneous, linear elastic, isotropic half-plane. This problem may be defined either as plane stress or plain strain (for solution with FreeFem++ we have choosen the latter). The traction boundary conditions for this problem are: \begin{equation}\label{eq:bc} \sigma_{xy}(x,0) = 0, \quad \sigma_{yy}(x,y) = -P\delta(x,y) \end{equation}
The analytical relations for the stresses can be found from the Flamant solution (stress distributions in a linear elastic wedge loaded by point forces a the tip. When the "wedge" is flat we get a half-plane. The derivation uses polar coordinates.) and are given as: \begin{equation} \sigma_{xx} = -\frac{2P}{\pi} \frac{x^2y}{\left(x^2+y^2\right)^2}, \end{equation}
Numerical solution with FreeFem++
Due to the known analytical solution the point-contact problem can be used for benchmarking numerical PDE solvers in terms of accuracy (as well as computational efficiency). The purpose of this section is to compare the numerical solution obtained by FreeFem++ with the analytical solution, as well as provide a reference numerical solution for the C++ library developed in our laboratory.
For purposes of simplicity we limit ourselves to the domain (x,y) \in \Omega = [-1,1] \times[-1,-0.1] and prescribe Dirichlet displacement on the boundaries \Gamma_D from the known analytical solution (\ref{eq:dispx}, \ref{eq:dispy}). This way we avoid having to deal with the Dirac delta traction boundary condition (\ref{eq:bc}). The problem can be described as find \boldsymbol{u(\boldsymbol{x})} that satisfies \begin{equation} \boldsymbol{\nabla}\cdot\boldsymbol{\sigma}= 0 \qquad \text{on }\Omega \end{equation}
To solve the point-contact problem in FreeFem++ we must first provide the weak form of the balance equation: \begin{equation*} \boldsymbol{\nabla}\cdot\boldsymbol{\sigma} + \boldsymbol{b} = 0. \end{equation*}
Equation (\ref{eq:weak}) can be handed to FreeFem++ with the help of Voigt or Mandel notation, that reduces the symmetric tensors \boldsymbol{\sigma} and \boldsymbol{\varepsilon} to vectors. The benefit of Mandel notation is that it allows the tensor scalar product to be performed as a scalar product of two vectors. For this reason we create the following macros:
macro u [ux,uy] // displacements macro v [vx,vy] // test function macro b [bx,by] // body forces macro e(u) [dx(u[0]),dy(u[1]),(dx(u[1])+dy(u[0]))/2] // strain (for post-processing) macro em(u) [dx(u[0]),dy(u[1]),sqrt(2)*(dx(u[1])+dy(u[0]))/2] // strain in Mandel notation macro A [[2*mu+lambda,mu,0],[mu,2*mu+lambda,0],[0,0,2*mu]] // stress-strain matrix
The weak form (\ref{eq:weak}) can then be expressed naturally in FreeFem++ syntax as
int2d(Th)((A*em(u))'*em(v)) - int2d(Th)(b'*v)
Numerical solution with Meshless Local Strong Form Method (MLSM)
1 RectangleDomain<vec_t> domain(O.domain_lo, O.domain_hi);
2 domain.fillUniformInteriorWithStep(O.d_space);
3 domain.fillUniformBoundaryWithStep(O.d_space);
4 domain.findSupport(O.n, internal);
5 domain.findSupport(O.n, boundary, internal, true); //search only among internal nodes + itself
6
7 EngineMLS<vec_t, Gaussians, Gaussians> mls(
8 {pow(domain.characteristicDistance()*O.sigmaB,2), O.m}, //basis functionsa
9 domain.positions[domain.support[0]],
10 pow(domain.characteristicDistance()*O.sigmaW,2)); //weight function */
11
12 auto mlsm = make_mlsm(domain, mls, internal);
13 for (auto& i : boundary) u_3[i] = u_anal(i);
14 /// [MAIN TEMPORAL LOOP]
15 for (size_t step = 0; step * O.dt < O.time; step++) {
16 int i;
17 ///[NAVIER-CAUCHY EQUATION :: explicit Plane Stress]
18 #pragma omp parallel for private(i) schedule(static)
19 for (i=0;i<internal.size();++i){
20 u_3[i] = O.dt * O.dt / O.rho * (
21 O.mu * mlsm.lap(u_2,i) + O.E/(2-2*O.v ) * mlsm.graddiv(u_2,i) +
22 force[i] - O.dampCoef * (u_2[i] - u_1[i])/O.dt
23 ) /// navier part
24 + 2 * u_2[i] - u_1[i];
25 }
26 ///[STEP FORWARD]
27 u_1 = u_2;
28 u_2 = u_3;
29 }
Results
Purely Dirichlet case
For starters we check the solution of the problem with Dirichlet BC. The conditions are obtained from closed form solution
Comparison of MLS based MLSM, FEM and analytical solution of displacements at different cross-sections.
The L^2 error norm is used to measure the "difference" between solutions. Since the displacements are the variable that obtain from FreeFem++, we use the displacement magnitude |\boldsymbol{u}| = \sqrt{u_x^2+u_y^2} to define our L^2-error norm. The exact equation we have used is \begin{equation} L^2\text{-norm} = \sqrt{\frac{\int_\Omega (|\boldsymbol{u_{\text{numerical}}}|-|\boldsymbol{u_{\text{analytical}}}|)^2d\Omega}{\int_\Omega|\boldsymbol{u_{\text{analytical}}}|^2d\Omega}}. \end{equation}