Ersatz Method

class ersatz_method.ErsatzMethod(level_set, V_ls, space_displacement, ds, bc, bc_velocity, parameters, shift)[source]

This is the Ersatz class.

Some details about the initialization of linear elasticity problem with Ersatz method.

Definition of the linear elasticity PDE

Linear elasticity problem is given by: Find \(u:\Omega \rightarrow \mathbb{R}^{d}\)

\[\begin{split}\begin{align} \begin{cases} -\text{div} \sigma(u) & \!\!\!\!=0 \text{ in }\Omega\\ u& \!\!\!\!=0\text{ on }\Gamma_{D}\\ \sigma(u)\cdot n & \!\!\!\!=g\text{ on }\Gamma_{N} \end{cases} \end{align}\end{split}\]

Where \(d\) the dimension of the problem. We assume small deformations and zero volumetric forces.

The fictitious material method uses a soft material in the domain \(D \ \overline{\Omega}\). Thus, we denote the adapted Lamé coefficients as follows: \(\widetilde{\lambda}\) and \(\widetilde{\mu}\). For more details on this method, we refer you to the explanation page CutFEM for Immersed geometry discretization in section Ersatz material method.

The following weak formulation: Find \(u \in V\), such that for all \(v \in V\) we have

\[a\left(u,v\right)=l\left(v\right)\]

with :

\[\begin{split}\begin{align} a\left(u,v\right)&=2\widetilde{\mu}\left(\varepsilon(u),\varepsilon\left(v\right)\right)_{L^{2}(\Omega)}+\widetilde{\lambda}\left(\nabla\cdot u,\nabla\cdot v\right)_{L^{2}(\Omega)}\\ l\left(v\right)&=\left(g,v\right)_{L^{2}\left(\Gamma_{N}\right)}, \end{align}\end{split}\]

Bilinear form:

import ufl

self.u = ufl.TrialFunction(self.space_displacement)
self.v = ufl.TestFunction(self.space_displacement)

self.uh = fem.Function(self.space_displacement)
self.ph = fem.Function(self.space_displacement)


self.lame_mu_fic = self.compute_lame_mu(self.xsi)
self.lame_lambda_fic = self.compute_lame_lambda(self.xsi)

self.a = 2.0*self.lame_mu_fic*ufl.inner(mechanics_tool.strain(self.u), mechanics_tool.strain(self.v)) * ufl.dx +\
    self.lame_lambda_fic * ufl.inner(ufl.nabla_div(self.u), ufl.nabla_div(self.v)) * ufl.dx

Linear form:

self.L = ufl.dot(self.shift,self.v) * self.ds(2)
advection(level_set, dt, velocity_field)[source]

Resolution of advection equation without any stabilization.

Parameters:
  • level_set (fem.Function) – The level set field which defined implicitly the domain \(\Omega\).

  • dt (float) – The :math: dt time parameters.

  • velocity_field (ufl.Expression) – The value of advection velocity_field, in normal direction of the interface \(\partial\Omega\).

Returns:

The values of the advected level set.

Return type:

fem.Function

compute_lame_lambda(xsi)[source]

Returns the Lamé coefficient \(\lambda\) according to the smoothing xsi.

Parameters:

xsi (fem.Expression) – The heaviside function \(\chi\).

Returns:

The fictitious \(\lambda\) Lamé coefficient update with the smoothing function.

Return type:

fem.Expression

compute_lame_mu(xsi)[source]

Returns the Lamé coefficient \(\mu\) according to the smoothing xsi.

Parameters:

xsi (fem.Expression) – The heaviside function \(\chi\).

Returns:

The fictitious \(\mu\) Lamé coefficient update with the smoothing function.

Return type:

fem.Expression

descent_direction(level_set, parameters, rest_constraint, constraint_integrand, cost_integrand, xsi_temp)[source]

Determine the descent direction by solving the following equation:

Find \(v'_{\text{reg}}\in H_{\Gamma_{D}}^{1}=\left\{ v\in H^{1}\left(D\right)\text{ such that }v=0\text{ on }\Gamma_{D}\right\}\) such that \(\forall w\in H_{\Gamma_{D}}^{1}\)

\[\alpha\left(\nabla v'_{\text{reg}},\nabla w\right)_{L^{2}\left(D\right)}+\left(v'_{\text{reg}},w\right)_{L^{2}\left(D\right)}=-J'(\Omega)\left(w\right)\]

with \(J\) the cost function and \(\alpha>0\) is a smoothing parameter instantiated in the Parameter class.

Parameters:
  • level_set (fem.Function) – The level set field which defined implicitly the domain \(\Omega\).

  • parameters (Parameters) – The object parameters.

  • rest_constraint (float) – The value of the constraint function \(C(\Omega)\).

  • constraint_integrand (fem.Expression) – The integrand of the constraint function.

  • cost_integrand (fem.Expression) – The integrand of the cost function.

Returns:

The velocity field, defined in D.

Return type:

fem.Function

ersatz_solver(level_set, parameters, adjoint=0)[source]

Resolution of the primal and dual problem with Ersatz method.

Parameters:
  • level_set (fem.Function) – The level set field which defined implicitly the domain \(\Omega\).

  • parameters (Parameters) – The object parameters.

  • adjoint (ufl.Expression) – The adjoint operator if needed.

Returns:

The values of the primal and dual solution.

Return type:

fem.Function, fem.Function

euclidean_norm_grad(func)[source]

Calculation of the integrand of the L2-norm of the gradient of the function provided, given by the following equality:

\[\left|\nabla\phi\right| = \sqrt{\nabla\phi\cdot\nabla\phi}\]
Parameters:

func (fem.Function) – Function field.

Returns:

The values of integrand of the L2-norm of the gradient.

Return type:

fem.Expression

heaviside(level_set)[source]

Returns the smoothing function \(\chi\).

Parameters:

level_set (fem.Expression) – The level_set function \(\phi\).

Returns:

The fictitious \(\mu\) Lamé coefficient update with the smoothing function.

Return type:

fem.Expression

primal_problem(level_set, parameters)[source]

Resolution of the primal problem with the Ersatz method.

Parameters:
  • level_set (fem.Function) – The level set field which defined implicitly the domain \(\Omega\).

  • parameters (Parameters) – The object parameters.

Returns:

The primal solution.

Return type:

fem.Function

set_bilin_form()[source]

Returns the update of the primal bilinear form with the smoothed value of the Lamé coefficients.

Returns:

The bilinear form of the primal and dual problems.

Return type:

fem.form