25 #include "getfem/getfem_contact_and_friction_large_sliding.h" 
   37 #define FRICTION_LAW 1 
   42   template <
typename VEC, 
typename VEC2, 
typename VECR>
 
   43   void aug_friction(
const VEC &lambda, scalar_type g, 
const VEC &Vs,
 
   44                     const VEC &n, scalar_type r, 
const VEC2 &f, VECR &F) {
 
   47     scalar_type lambdan_aug = gmm::neg(lambdan + r * g);
 
   49     scalar_type tau = ((i >= 3) ? f[2] : scalar_type(0)) + f[0]*lambdan_aug;
 
   50     if (i >= 2) tau = std::min(tau, f[1]);
 
   52     if (tau > scalar_type(0)) {
 
   53       gmm::add(lambda, gmm::scaled(Vs, -r), F);
 
   57       if (norm > tau) gmm::scale(F, tau / norm);
 
   60     gmm::add(gmm::scaled(n, -lambdan_aug/nn), F);
 
   63   template <
typename VEC, 
typename VEC2, 
typename VECR, 
typename MAT>
 
   64   void aug_friction_grad(
const VEC &lambda, scalar_type g, 
const VEC &Vs,
 
   65                          const VEC &n, scalar_type r, 
const VEC2 &f, VECR &F,
 
   66                          MAT &dlambda, VECR &dg, MAT &dn, MAT &dVs) {
 
   70     scalar_type lambdan_aug = gmm::neg(lambdan + r * g);
 
   72     scalar_type tau = ((i >= 3) ? f[2] : scalar_type(0)) + f[0]*lambdan_aug;
 
   73     if (i >= 2) tau = std::min(tau, f[1]);
 
   76     if (tau > scalar_type(0)) {
 
   77       gmm::add(lambda, gmm::scaled(Vs, -r), F);
 
   82       gmm::scale(dn, -mu/nn);
 
   83       gmm::rank_one_update(dn, gmm::scaled(n, mu/(nn*nn*nn)), n);
 
   84       gmm::rank_one_update(dn, gmm::scaled(n, scalar_type(-1)/(nn*nn)), F);
 
   86       gmm::rank_one_update(dVs, n, gmm::scaled(n, scalar_type(-1)/(nn*nn)));
 
   89         gmm::rank_one_update(dVs, F,
 
   90                              gmm::scaled(F, scalar_type(-1)/(norm*norm)));
 
   91         gmm::scale(dVs, tau / norm);
 
   92         gmm::copy(gmm::scaled(F, scalar_type(1)/norm), dg);
 
   93         gmm::rank_one_update(dn, gmm::scaled(F, mu/(norm*norm*nn)), F);
 
   94         gmm::scale(dn, tau / norm);
 
   95         gmm::scale(F, tau / norm);
 
  103     if (norm > tau && ((i <= 1) || tau < f[1]) && ((i <= 2) || tau > f[2])) {
 
  104       gmm::rank_one_update(dn, dg, gmm::scaled(lambda, -f[0]/nn));
 
  105       gmm::rank_one_update(dn, dg, gmm::scaled(n, f[0]*lambdan/(nn*nn)));
 
  106       gmm::rank_one_update(dlambda, dg, gmm::scaled(n, -f[0]/nn));
 
  107       gmm::scale(dg, -f[0]*r);
 
  109     if (lambdan_aug > scalar_type(0)) {
 
  111       gmm::rank_one_update(dlambda, n, gmm::scaled(n, scalar_type(1)/(nn*nn)));
 
  112       gmm::rank_one_update(dn, gmm::scaled(n, scalar_type(1)/(nn*nn)), lambda);
 
  113       gmm::rank_one_update(dn,
 
  114                            gmm::scaled(n,(lambdan_aug-lambdan)/(nn*nn*nn)), n);
 
  115       for (
size_type j = 0; j < N; ++j) dn(j,j) -= lambdan_aug/nn;
 
  117     gmm::add(gmm::scaled(n, -lambdan_aug/nn), F);
 
  122 #elif FRICTION_LAW == 2  
  124   template <
typename VEC, 
typename VEC2, 
typename VECR>
 
  125   void aug_friction(
const VEC &lambda, scalar_type g, 
const VEC &,
 
  126                     const VEC &n, scalar_type r, 
const VEC2 &, VECR &F) {
 
  129     scalar_type lambdan_aug = gmm::neg(lambdan + r * g);
 
  130     gmm::copy(gmm::scaled(n, -lambdan_aug/nn), F);
 
  133   template <
typename VEC, 
typename VEC2, 
typename VECR, 
typename MAT>
 
  134   void aug_friction_grad(
const VEC &lambda, scalar_type g, 
const VEC &,
 
  135                          const VEC &n, scalar_type r, 
const VEC2 &, VECR &F,
 
  136                          MAT &dlambda, VECR &dg, MAT &dn, MAT &dVs) {
 
  140     scalar_type lambdan_aug = gmm::neg(lambdan + r * g);
 
  147     if (lambdan_aug > scalar_type(0)) {
 
  149       gmm::rank_one_update(dlambda, n, gmm::scaled(n, scalar_type(1)/(nn*nn)));
 
  150       gmm::rank_one_update(dn, gmm::scaled(n, scalar_type(1)/(nn*nn)), lambda);
 
  151       gmm::rank_one_update(dn,
 
  152                            gmm::scaled(n,(lambdan_aug-lambdan)/(nn*nn*nn)), n);
 
  153       for (
size_type j = 0; j < N; ++j) dn(j,j) -= lambdan_aug/nn;
 
  155     gmm::add(gmm::scaled(n, -lambdan_aug/nn), F);
 
  162 #elif FRICTION_LAW == 3  
  164   template <
typename VEC, 
typename VEC2, 
typename VECR>
 
  165   void aug_friction(
const VEC &lambda, scalar_type g, 
const VEC &Vs,
 
  166                     const VEC &n, scalar_type r, 
const VEC2 &f, VECR &F) {
 
  167     gmm::copy(gmm::scaled(lambda, g*r*f[0]), F); 
 
  173   template <
typename VEC, 
typename VEC2, 
typename VECR, 
typename MAT>
 
  174   void aug_friction_grad(
const VEC &lambda, scalar_type g, 
const VEC &Vs,
 
  175                          const VEC &n, scalar_type r, 
const VEC2 &f, VECR &F,
 
  176                          MAT &dlambda, VECR &dg, MAT &dn, MAT &dVs) {
 
  177     gmm::copy(gmm::scaled(lambda, g*r*f[0]), F); 
 
  187 #elif FRICTION_LAW == 4  
  189   template <
typename VEC, 
typename VEC2, 
typename VECR>
 
  190   void aug_friction(
const VEC &lambda, scalar_type g, 
const VEC &Vs,
 
  191                     const VEC &n, scalar_type r, 
const VEC2 &f, VECR &F) {
 
  192     gmm::copy(gmm::scaled(lambda, g*r*f[0]*n[0]*Vs[0]), F); 
 
  196   template <
typename VEC, 
typename VEC2, 
typename VECR, 
typename MAT>
 
  197   void aug_friction_grad(
const VEC &lambda, scalar_type g, 
const VEC &Vs,
 
  198                          const VEC &n, scalar_type r, 
const VEC2 &f, VECR &F,
 
  199                          MAT &dlambda, VECR &dg, MAT &dn, MAT &dVs) {
 
  200     gmm::copy(gmm::scaled(lambda, g*r*f[0]*n[0]*Vs[0]), F); 
 
  205     gmm::copy(gmm::identity_matrix(), dlambda);
 
  208 #elif FRICTION_LAW == 5  
  210   template <
typename VEC, 
typename VEC2, 
typename VECR>
 
  211   void aug_friction(
const VEC &lambda, scalar_type g, 
const VEC &Vs,
 
  212                     const VEC &n, scalar_type r, 
const VEC2 &f, VECR &F) {
 
  213     gmm::copy(gmm::scaled(lambda, g*r*f[0]*n[0]*Vs[0]), F); 
 
  217   template <
typename VEC, 
typename VEC2, 
typename VECR, 
typename MAT>
 
  218   void aug_friction_grad(
const VEC &lambda, scalar_type g, 
const VEC &Vs,
 
  219                          const VEC &n, scalar_type r, 
const VEC2 &f, VECR &F,
 
  220                          MAT &dlambda, VECR &dg, MAT &dn, MAT &dVs) {
 
  221     gmm::copy(gmm::scaled(lambda, g*r*f[0]*n[0]*Vs[0]), F); 
 
  243   struct integral_large_sliding_contact_brick : 
public virtual_brick {
 
  245     multi_contact_frame &mcf;
 
  249     virtual void asm_real_tangent_terms(
const model &md, 
size_type ,
 
  250                                         const model::varnamelist &vl,
 
  251                                         const model::varnamelist &dl,
 
  252                                         const model::mimlist &mims,
 
  253                                         model::real_matlist &matl,
 
  254                                         model::real_veclist &vecl,
 
  255                                         model::real_veclist &,
 
  257                                         build_version version) 
const;
 
  259     integral_large_sliding_contact_brick(multi_contact_frame &mcff,
 
  261       : mcf(mcff), with_friction(with_fric) {
 
  262       set_flags(
"Integral large sliding contact brick",
 
  271   struct gauss_point_precomp {
 
  273     fem_precomp_pool fppool;
 
  274     const multi_contact_frame &mcf;
 
  276     const multi_contact_frame::contact_pair *cp;
 
  278     const base_node &x(
void)
 const { 
return cp->slave_point; }
 
  279     const base_node &nx(
void)
 const { 
return cp->slave_n; }
 
  280     const base_node &y(
void)
 const { 
return cp->master_point; }
 
  281     const base_node &y_ref(
void)
 const { 
return cp->master_point_ref; }
 
  282     const base_node &ny(
void)
 const { 
return cp->master_n; }
 
  283     scalar_type g(
void)
 const { 
return cp->signed_dist; }
 
  286     bool I_nxnx_computed;
 
  287     const base_matrix &I_nxnx(
void) {
 
  288       if (!I_nxnx_computed) {
 
  289         gmm::copy(gmm::identity_matrix(), I_nxnx_);
 
  290         gmm::rank_one_update(I_nxnx_, nx(), gmm::scaled(nx(),scalar_type(-1)));
 
  291         I_nxnx_computed = 
true;
 
  297     bool I_nyny_computed;
 
  298     const base_matrix &I_nyny(
void) {
 
  299       if (!I_nyny_computed) {
 
  300         gmm::copy(gmm::identity_matrix(), I_nyny_);
 
  301         gmm::rank_one_update(I_nyny_, ny(), gmm::scaled(ny(),scalar_type(-1)));
 
  302         I_nyny_computed = 
true;
 
  308     bool I_nxny_computed;
 
  309     const base_matrix &I_nxny(
void) {
 
  310       if (!I_nxny_computed) {
 
  311         gmm::copy(gmm::identity_matrix(), I_nxny_);
 
  312         gmm::rank_one_update(I_nxny_, nx(),
 
  313                              gmm::scaled(ny(),scalar_type(-1)/nxny));
 
  314         I_nxny_computed = 
true;
 
  320     scalar_type nxdotny(
void)
 const { 
return nxny; }
 
  323     bool isrigid(
void) { 
return isrigid_; }
 
  325     base_tensor base_ux, base_uy, base_lx, base_ly;
 
  326     base_matrix vbase_ux_, vbase_uy_, vbase_lx_, vbase_ly_;
 
  327     bool vbase_ux_init, vbase_uy_init, vbase_lx_init, vbase_ly_init;
 
  328     base_tensor grad_base_ux, grad_base_uy, vgrad_base_ux_, vgrad_base_uy_;
 
  329     bool vgrad_base_ux_init, vgrad_base_uy_init;
 
  330     bool have_lx, have_ly;
 
  332     fem_interpolation_context ctx_ux_, ctx_uy_, ctx_lx_, ctx_ly_;
 
  333     bool ctx_ux_init, ctx_uy_init, ctx_lx_init, ctx_ly_init;
 
  335     const mesh_fem *mf_ux_, *mf_uy_, *mf_lx_, *mf_ly_;
 
  336     gmm::sub_interval I_ux_, I_uy_, I_lx_, I_ly_;
 
  337     pfem pf_ux, pf_uy, pf_lx, pf_ly;
 
  338     size_type ndof_ux_, qdim_ux, ndof_uy_, qdim_uy, ndof_lx_, qdim_lx;
 
  339     size_type ndof_ly_, qdim_ly, cvx_, cvy_, ibx, iby;
 
  343     pintegration_method pim;
 
  346     scalar_type weight(
void) { 
return weight_; }
 
  348     const mesh &meshx(
void)
 const { 
return mf_ux_->linked_mesh(); }
 
  349     const mesh &meshy(
void)
 const { 
return mf_uy_->linked_mesh(); }
 
  350     const mesh_fem *mf_ux(
void)
 const { 
return mf_ux_; }
 
  351     const mesh_fem *mf_uy(
void)
 const { 
return mf_uy_; }
 
  352     const mesh_fem *mf_lx(
void)
 const { 
return mf_lx_; }
 
  353     const mesh_fem *mf_ly(
void)
 const { 
return mf_ly_; }
 
  354     size_type ndof_ux(
void)
 const { 
return ndof_ux_; }
 
  355     size_type ndof_uy(
void)
 const { 
return ndof_uy_; }
 
  356     size_type ndof_lx(
void)
 const { 
return ndof_lx_; }
 
  357     size_type ndof_ly(
void)
 const { 
return ndof_ly_; }
 
  358     size_type cvx(
void)
 const { 
return cvx_; }
 
  359     size_type cvy(
void)
 const { 
return cvy_; }
 
  360     const gmm::sub_interval I_ux(
void)
 const { 
return I_ux_; }
 
  361     const gmm::sub_interval I_uy(
void)
 const { 
return I_uy_; }
 
  362     const gmm::sub_interval I_lx(
void)
 const { 
return I_lx_; }
 
  363     const gmm::sub_interval I_ly(
void)
 const { 
return I_ly_; }
 
  366     fem_interpolation_context &ctx_ux(
void) {
 
  368         bgeot::vectors_to_base_matrix(Gx, meshx().points_of_convex(cvx_));
 
  370           = fppool(pf_ux, pim->approx_method()->pintegration_points());
 
  371         ctx_ux_ = fem_interpolation_context(pgtx, pfp_ux, cp->slave_ind_pt,
 
  378     fem_interpolation_context &ctx_lx(
void) {
 
  379       GMM_ASSERT1(have_lx, 
"No multiplier defined on the slave surface");
 
  382           = fppool(pf_lx, pim->approx_method()->pintegration_points());
 
  383         ctx_lx_ = fem_interpolation_context(pgtx, pfp_lx, cp->slave_ind_pt,
 
  384                                             ctx_ux().G(), cvx_, fx);
 
  390     fem_interpolation_context &ctx_uy(
void) {
 
  391       GMM_ASSERT1(!isrigid(), 
"Rigid obstacle master node: no fem defined");
 
  393         bgeot::vectors_to_base_matrix(Gy, meshy().points_of_convex(cvy_));
 
  394         ctx_uy_ = fem_interpolation_context(pgty, pf_uy, y_ref(), Gy, cvy_, fy);
 
  400     fem_interpolation_context &ctx_ly(
void) {
 
  401       GMM_ASSERT1(have_ly, 
"No multiplier defined on the master surface");
 
  403         ctx_ly_ = fem_interpolation_context(pgty, pf_ly, y_ref(),
 
  404                                             ctx_uy().G(), cvy_, fy);
 
  410     const base_matrix &vbase_ux(
void) {
 
  411       if (!vbase_ux_init) {
 
  412         ctx_ux().base_value(base_ux);
 
  413         vectorize_base_tensor(base_ux, vbase_ux_, ndof_ux_, qdim_ux, N);
 
  414         vbase_ux_init = 
true;
 
  419     const base_matrix &vbase_uy(
void) {
 
  420       if (!vbase_uy_init) {
 
  421         ctx_uy().base_value(base_uy);
 
  422         vectorize_base_tensor(base_uy, vbase_uy_, ndof_uy_, qdim_uy, N);
 
  423         vbase_uy_init = 
true;
 
  428     const base_matrix &vbase_lx(
void) {
 
  429       if (!vbase_lx_init) {
 
  430         ctx_lx().base_value(base_lx);
 
  431         vectorize_base_tensor(base_lx, vbase_lx_, ndof_lx_, qdim_lx, N);
 
  432         vbase_lx_init = 
true;
 
  437     const base_matrix &vbase_ly(
void) {
 
  438       if (!vbase_ly_init) {
 
  439         ctx_ly().base_value(base_ly);
 
  440         vectorize_base_tensor(base_ly, vbase_ly_, ndof_ly_, qdim_ly, N);
 
  441         vbase_ly_init = 
true;
 
  446     const base_tensor &vgrad_base_ux(
void) {
 
  447       if (!vgrad_base_ux_init) {
 
  448         ctx_ux().grad_base_value(grad_base_ux);
 
  449         vectorize_grad_base_tensor(grad_base_ux, vgrad_base_ux_, ndof_ux_,
 
  451         vgrad_base_ux_init = 
true;
 
  453       return vgrad_base_ux_;
 
  456     const base_tensor &vgrad_base_uy(
void) {
 
  457       if (!vgrad_base_uy_init) {
 
  458         ctx_uy().grad_base_value(grad_base_uy);
 
  459         vectorize_grad_base_tensor(grad_base_uy, vgrad_base_uy_, ndof_uy_,
 
  461         vgrad_base_uy_init = 
true;
 
  463       return vgrad_base_uy_;
 
  466     base_small_vector lambda_x_, lambda_y_;
 
  467     bool lambda_x_init, lambda_y_init;
 
  470     const base_small_vector &lx(
void) {
 
  471       if (!lambda_x_init) {
 
  472         pfem pf = ctx_lx().pf();
 
  475         pf->interpolation(ctx_lx(), coeff, lambda_x_, dim_type(N));
 
  476         lambda_x_init = 
true;
 
  481     const base_small_vector &ly(
void) {
 
  482       if (!lambda_y_init) {
 
  483         pfem pf = ctx_ly().pf();
 
  486         pf->interpolation(ctx_ly(), coeff, lambda_y_, dim_type(N));
 
  487         lambda_y_init = 
true;
 
  492     base_matrix grad_phix_, grad_phix_inv_, grad_phiy_, grad_phiy_inv_;
 
  493     bool grad_phix_init, grad_phix_inv_init;
 
  494     bool grad_phiy_init, grad_phiy_inv_init;
 
  496     const base_matrix &grad_phix(
void) {
 
  497       if (!grad_phix_init) {
 
  498         pfem pf = ctx_ux().pf();
 
  501         pf->interpolation_grad(ctx_ux(), coeff, grad_phix_, dim_type(N));
 
  502         gmm::add(gmm::identity_matrix(), grad_phix_);
 
  503         grad_phix_init = 
true;
 
  508     const base_matrix &grad_phix_inv(
void) {
 
  509       if (!grad_phix_inv_init) {
 
  513         grad_phix_inv_init = 
true;
 
  515       return grad_phix_inv_;
 
  518     const base_matrix &grad_phiy(
void) {
 
  519       if (!grad_phiy_init) {
 
  520         pfem pf = ctx_uy().pf();
 
  523         pf->interpolation_grad(ctx_uy(), coeff, grad_phiy_, dim_type(N));
 
  524         gmm::add(gmm::identity_matrix(), grad_phiy_);
 
  525         grad_phiy_init = 
true;
 
  530     const base_matrix &grad_phiy_inv(
void) {
 
  531       if (!grad_phiy_inv_init) {
 
  535         grad_phiy_inv_init = 
true;
 
  537       return grad_phiy_inv_;
 
  541     base_small_vector x0_, y0_, nx0_, Vs_;
 
  542     bool x0_init, y0_init, nx0_init, Vs_init;
 
  543     base_matrix grad_phiy0_;
 
  544     bool grad_phiy0_init;
 
  546     const base_small_vector &x0(
void) {
 
  548         const model_real_plain_vector &all_x0 = mcf.disp0_of_boundary(ibx);
 
  550           pfem pf = ctx_ux().pf();
 
  552           pf->interpolation(ctx_ux(), coeff, x0_, dim_type(N));
 
  560     const base_small_vector &y0(
void) {
 
  563           const model_real_plain_vector &all_y0 = mcf.disp0_of_boundary(iby);
 
  565             pfem pf = ctx_uy().pf();
 
  567             pf->interpolation(ctx_uy(), coeff, y0_, dim_type(N));
 
  576     const base_small_vector &nx0(
void) {
 
  578         const model_real_plain_vector &all_x0 = mcf.disp0_of_boundary(ibx);
 
  580           pfem pf = ctx_ux().pf();
 
  582           base_small_vector nx00_(N);
 
  583           base_matrix grad_phix0_(N,N);
 
  584           compute_normal(ctx_ux(), fx, 
false, coeff, nx00_, nx0_, grad_phix0_);
 
  592     const base_small_vector &Vs(
void) { 
 
  594         if (alpha != scalar_type(0)) {
 
  595 #ifdef CONSIDER_FRAME_INDIFFERENCE 
  597             gmm::add(y0(), gmm::scaled(x0(), scalar_type(-1)), Vs_);
 
  598             gmm::add(gmm::scaled(nx0(), -g()), Vs_);
 
  602             gmm::add(x(), gmm::scaled(y(), scalar_type(-1)), Vs_);
 
  603             gmm::add(gmm::scaled(x0(), scalar_type(-1)), Vs_);
 
  606           gmm::scale(Vs_, alpha);
 
  613     const base_matrix &grad_phiy0(
void) { 
 
  615       if (!grad_phiy0_init) {
 
  616         const model_real_plain_vector &all_y0 = mcf.disp0_of_boundary(iby);
 
  617         if (!isrigid() && all_y0.size()) {
 
  618           pfem pf = ctx_uy().pf();
 
  620           pf->interpolation_grad(ctx_uy(), coeff, grad_phiy0_, dim_type(N));
 
  621           gmm::add(gmm::identity_matrix(), grad_phiy0_);
 
  622         } 
else gmm::copy(gmm::identity_matrix(), grad_phiy0_);
 
  623         grad_phiy0_init = 
true;
 
  628     base_small_vector un;
 
  630     void set_pair(
const multi_contact_frame::contact_pair &cp_) {
 
  632       I_nxnx_computed = I_nyny_computed = I_nxny_computed = 
false;
 
  633       ctx_ux_init = ctx_uy_init = ctx_lx_init = ctx_ly_init = 
false;
 
  634       vbase_ux_init = vbase_uy_init = vbase_lx_init = vbase_ly_init = 
false;
 
  635       vgrad_base_ux_init = vgrad_base_uy_init = 
false;
 
  636       lambda_x_init = lambda_y_init = 
false;
 
  637       have_lx = have_ly = 
false;
 
  638       grad_phix_init = grad_phiy_init = 
false;
 
  639       grad_phix_inv_init = grad_phiy_inv_init = 
false;
 
  640       x0_init = y0_init = Vs_init = grad_phiy0_init = 
false;
 
  642       isrigid_ = (cp->irigid_obstacle != 
size_type(-1));
 
  644       cvx_ = cp->slave_ind_element;
 
  645       ibx = cp->slave_ind_boundary;
 
  646       mf_ux_ = &(mcf.mfdisp_of_boundary(ibx));
 
  647       pf_ux = mf_ux_->fem_of_element(cvx_);
 
  648       qdim_ux = pf_ux->target_dim();
 
  649       ndof_ux_ = pf_ux->nb_dof(cvx_) * N / qdim_ux;
 
  650       fx = cp->slave_ind_face;
 
  651       pgtx = meshx().trans_of_convex(cvx_);
 
  652       mim = &(mcf.mim_of_boundary(ibx));
 
  653       pim = mim->int_method_of_element(cvx_);
 
  654       weight_ = pim->approx_method()->coeff(cp->slave_ind_pt) * ctx_ux().J();
 
  655       gmm::mult(ctx_ux().B(), pgtx->normals()[fx], un);
 
  657       const std::string &name_ux = mcf.varname_of_boundary(ibx);
 
  658       I_ux_ = md.interval_of_variable(name_ux);
 
  660       const std::string &name_lx = mcf.multname_of_boundary(ibx);
 
  661       have_lx = (name_lx.size() > 0);
 
  663         mf_lx_ = &(mcf.mfmult_of_boundary(ibx));
 
  664         I_lx_ = md.interval_of_variable(name_lx);
 
  665         pf_lx = mf_lx_->fem_of_element(cvx_);
 
  666         qdim_lx = pf_lx->target_dim();
 
  667         ndof_lx_ = pf_lx->nb_dof(cvx_) * N / qdim_lx;
 
  671         cvy_ = cp->master_ind_element;
 
  672         iby = cp->master_ind_boundary;
 
  673         fy = cp->master_ind_face;
 
  674         mf_uy_ = &(mcf.mfdisp_of_boundary(iby));
 
  675         pf_uy = mf_uy_->fem_of_element(cvy_);
 
  676         qdim_uy = pf_uy->target_dim();
 
  677         ndof_uy_ = pf_uy->nb_dof(cvy_) * N / qdim_uy;
 
  678         pgty = meshy().trans_of_convex(cvy_);
 
  680         const std::string &name_uy = mcf.varname_of_boundary(iby);
 
  681         I_uy_ = md.interval_of_variable(name_uy);
 
  682         const std::string &name_ly = mcf.multname_of_boundary(iby);
 
  683         have_ly = (name_ly.size() > 0);
 
  685           mf_ly_ = &(mcf.mfmult_of_boundary(iby));
 
  686           I_ly_ = md.interval_of_variable(name_ly);
 
  687           pf_ly = mf_ly_->fem_of_element(cvy_);
 
  688           qdim_ly = pf_ly->target_dim();
 
  689           ndof_ly_ = pf_ly->nb_dof(cvy_) * N / qdim_ly;
 
  694     gauss_point_precomp(
size_type N_, 
const model &md_,
 
  695                         const multi_contact_frame &mcf_, scalar_type alpha_) :
 
  696       N(N_), mcf(mcf_), md(md_),
 
  697       I_nxnx_(N,N), I_nyny_(N,N), I_nxny_(N,N),
 
  698       lambda_x_(N), lambda_y_(N),
 
  699       grad_phix_(N, N), grad_phix_inv_(N, N),
 
  700       grad_phiy_(N, N), grad_phiy_inv_(N, N), 
alpha(alpha_),
 
  701       x0_(N), y0_(N), nx0_(N), Vs_(N), grad_phiy0_(N, N), un(N) {}
 
  779   void integral_large_sliding_contact_brick::asm_real_tangent_terms
 
  780   (
const model &md, 
size_type , 
const model::varnamelist &vl,
 
  781    const model::varnamelist &dl, 
const model::mimlist &,
 
  782    model::real_matlist &matl, model::real_veclist &vecl,
 
  784    build_version version)
 const {
 
  787     GMM_ASSERT1((with_friction && dl.size() >= 2 && dl.size() <= 3)
 
  788                 || (!with_friction && dl.size() >= 1 && dl.size() <= 2),
 
  789             "Wrong number of data for integral large sliding contact brick");
 
  791     GMM_ASSERT1(vl.size() == mcf.nb_variables() + mcf.nb_multipliers(),
 
  792                 "For the moment, it is not allowed to add boundaries to " 
  793                 "the multi contact frame object after a model brick has " 
  796     const model_real_plain_vector &vr = md.real_variable(dl[0]);
 
  797     GMM_ASSERT1(gmm::vect_size(vr) == 1, 
"Large sliding contact " 
  798                     "brick: parameter r should be a scalar");
 
  799     scalar_type r = vr[0];
 
  801     model_real_plain_vector f_coeff;
 
  803       f_coeff = md.real_variable(dl[1]);
 
  804       GMM_ASSERT1(gmm::vect_size(f_coeff) <= 3,
 
  805                   "Large sliding contact " 
  806                   "brick: the friction law has less than 3 parameters");
 
  808     if (gmm::vect_size(f_coeff) == 0) 
 
  809       { f_coeff.resize(1); f_coeff[0] = scalar_type(0); }
 
  811     scalar_type 
alpha(0);
 
  813     if (dl.size() >= ind+1) {
 
  814       GMM_ASSERT1(md.real_variable(dl[ind]).size() == 1,
 
  815                   "Large sliding contact " 
  816                   "brick: parameter alpha should be a scalar");
 
  817       alpha = md.real_variable(dl[ind])[0];
 
  820     GMM_ASSERT1(matl.size() == 1,
 
  821                 "Large sliding contact brick should have only one term");
 
  822     model_real_sparse_matrix &M = matl[0]; 
gmm::clear(M);
 
  823     model_real_plain_vector &V = vecl[0]; 
gmm::clear(V);
 
  825     mcf.set_raytrace(
true);
 
  826     mcf.set_nodes_mode(0);
 
  827     mcf.compute_contact_pairs();
 
  831     base_matrix dlambdaF(N, N), dnF(N, N), dVsF(N, N);
 
  832     base_small_vector F(N), dgF(N);
 
  834     scalar_type FMULT = 1.;
 
  837     for (
size_type i = 0; i < mcf.nb_boundaries(); ++i)
 
  838       if (mcf.is_self_contact() || mcf.is_slave_boundary(i)) {
 
  839         size_type region = mcf.region_of_boundary(i);
 
  840         const std::string &name_lx = mcf.multname_of_boundary(i);
 
  841         GMM_ASSERT1(name_lx.size() > 0, 
"This brick need " 
  842                     "multipliers defined on the multi_contact_frame object");
 
  843         const mesh_fem &mflambda = mcf.mfmult_of_boundary(i);
 
  844         const mesh_im &mim = mcf.mim_of_boundary(i);
 
  845         const gmm::sub_interval &I = md.interval_of_variable(name_lx);
 
  847         if (version & model::BUILD_MATRIX) { 
 
  848           model_real_sparse_matrix M1(mflambda.nb_dof(), mflambda.nb_dof());
 
  850           gmm::add(gmm::scaled(M1, FMULT/r), gmm::sub_matrix(M, I, I));
 
  853         if (version & model::BUILD_RHS) { 
 
  854           model_real_plain_vector V1(mflambda.nb_dof());
 
  856             (V1, mim, mflambda, mflambda,
 
  857              md.real_variable(mcf.multname_of_boundary(i)), region);
 
  858           gmm::add(gmm::scaled(V1, -FMULT/r), gmm::sub_vector(V, I));
 
  862     gauss_point_precomp gpp(N, md, mcf, alpha);
 
  866     base_matrix auxNN1(N, N), auxNN2(N, N);
 
  867     base_small_vector auxN1(N), auxN2(N);
 
  870     for (
size_type icp = 0; icp < mcf.nb_contact_pairs(); ++icp) {
 
  871       const multi_contact_frame::contact_pair &cp = mcf.get_contact_pair(icp);
 
  873       const base_small_vector &nx = gpp.nx(), &ny = gpp.ny();
 
  874       const mesh_fem *mf_ux = gpp.mf_ux(), *mf_lx = gpp.mf_lx(), *mf_uy(0);
 
  875       size_type ndof_ux = gpp.ndof_ux(), ndof_uy(0), ndof_lx = gpp.ndof_lx();
 
  877       const gmm::sub_interval &I_ux = gpp.I_ux(), &I_lx = gpp.I_lx();
 
  878       gmm::sub_interval I_uy;
 
  879       bool isrigid = gpp.isrigid();
 
  881         ndof_uy = gpp.ndof_uy(); I_uy = gpp.I_uy();
 
  882         mf_uy = gpp.mf_uy(); cvy =  gpp.cvy();
 
  884       scalar_type weight = gpp.weight(), g = gpp.g();
 
  885       const base_small_vector &lambda = gpp.lx();
 
  887       base_vector auxUX(ndof_ux), auxUY(ndof_uy);
 
  889       if (version & model::BUILD_MATRIX) {
 
  891         base_matrix auxUYN(ndof_uy, N);
 
  892         base_matrix auxLXN1(ndof_lx, N), auxLXN2(ndof_lx, N);
 
  894         aug_friction_grad(lambda, g, gpp.Vs(), nx, r, f_coeff, F, dlambdaF,
 
  898         const base_tensor &vgrad_base_ux = gpp.vgrad_base_ux();
 
  899         base_matrix graddeltaunx(ndof_ux, N);
 
  903               graddeltaunx(i, j) += nx[k] * vgrad_base_ux(i, k, j);
 
  905 #define CONSIDER_TERM1 
  906 #define CONSIDER_TERM2 
  907 #define CONSIDER_TERM3 
  910 #ifdef CONSIDER_TERM1 
  913         gmm::mult(gpp.vbase_ux(), gmm::transposed(gpp.vbase_lx()), Melem);
 
  914         gmm::scale(Melem, -weight);
 
  915         mat_elem_assembly(M, I_ux, I_lx, Melem, *mf_ux, cvx, *mf_lx, cvx);
 
  918 #ifdef CONSIDER_TERM2 
  923           gmm::mult(gpp.vbase_uy(), gmm::transposed(gpp.vbase_lx()), Melem);
 
  924           gmm::scale(Melem, weight);
 
  925           mat_elem_assembly(M, I_uy, I_lx, Melem, *mf_uy, cvy, *mf_lx, cvx);
 
  929           const base_tensor &vgrad_base_uy = gpp.vgrad_base_uy();
 
  933                 auxUYN(i, j) += lambda[k] * vgrad_base_uy(i, k, j);
 
  934           base_matrix lgraddeltavgradphiyinv(ndof_uy, N);
 
  935           gmm::mult(auxUYN, gpp.grad_phiy_inv(), lgraddeltavgradphiyinv);
 
  939           gmm::mult(lgraddeltavgradphiyinv, gpp.I_nxny(), auxUYN);
 
  940           gmm::mult(auxUYN, gmm::transposed(gpp.vbase_uy()), Melem);
 
  942           gmm::scale(Melem, -weight);
 
  943           mat_elem_assembly(M, I_uy, I_uy, Melem, *mf_uy, cvy, *mf_uy, cvy);
 
  949           gmm::mult(auxUYN, gmm::transposed(gpp.vbase_ux()), Melem);
 
  952           base_matrix auxUYUX(ndof_uy, ndof_ux);
 
  953           gmm::mult(gpp.I_nxny(), gmm::transposed(gpp.grad_phix_inv()), auxNN1);
 
  954           gmm::mult(lgraddeltavgradphiyinv, auxNN1, auxUYN);
 
  955           gmm::mult(auxUYN, gmm::transposed(graddeltaunx), auxUYUX);
 
  956           gmm::scale(auxUYUX, -g);
 
  958           gmm::scale(Melem, weight);
 
  959           mat_elem_assembly(M, I_uy, I_ux, Melem, *mf_uy, cvy, *mf_ux, cvx);
 
  965 #ifdef CONSIDER_TERM3 
  971         gmm::copy(gmm::scaled(dlambdaF, scalar_type(-1)/r), auxNN1);
 
  972         gmm::mult(gpp.vbase_lx(), auxNN1, auxLXN1);
 
  973         gmm::mult(auxLXN1, gmm::transposed(gpp.vbase_lx()), Melem);
 
  974         gmm::scale(Melem, weight*FMULT);
 
  975         mat_elem_assembly(M, I_lx, I_lx, Melem, *mf_lx, cvx, *mf_lx, cvx);
 
  981         gmm::mult(auxLXN1, gpp.I_nxnx(), auxLXN2);
 
  982         gmm::mult(auxLXN2,  gmm::transposed(gpp.grad_phix_inv()), auxLXN1);
 
  983         gmm::mult(auxLXN1, gmm::transposed(graddeltaunx), Melem);
 
  984         gmm::scale(Melem, scalar_type(1)/r);
 
  988         base_vector deltamudgF(ndof_lx);
 
  990                   gmm::scaled(dgF, scalar_type(1)/(r*gpp.nxdotny())),
 
  997         gmm::mult(gpp.I_nxnx(), gmm::scaled(ny, -g), auxN1);
 
  998         gmm::mult(gpp.grad_phix_inv(), auxN1, auxN2);
 
 1000         gmm::rank_one_update(Melem, deltamudgF,  auxUX);
 
 1001         gmm::scale(Melem, weight*FMULT);
 
 1002         mat_elem_assembly(M, I_lx, I_ux, Melem, *mf_lx, cvx, *mf_ux, cvx);
 
 1009           gmm::rank_one_update(Melem, deltamudgF, auxUY);
 
 1010           gmm::scale(Melem, -weight*FMULT);
 
 1011           mat_elem_assembly(M, I_lx, I_uy, Melem, *mf_lx, cvx, *mf_uy, cvy);
 
 1014         if (alpha != scalar_type(0)) {
 
 1018 #ifdef CONSIDER_FRAME_INDIFFERENCE 
 1019             base_matrix gphiy0gphiyinv(N, N);
 
 1020             gmm::mult(gpp.grad_phiy0(), gpp.grad_phiy_inv(), gphiy0gphiyinv);
 
 1021             gmm::mult(gphiy0gphiyinv, gpp.I_nxny(), auxNN1);
 
 1022             gmm::rank_one_update(auxNN1, gpp.nx0(),
 
 1023                                  gmm::scaled(gpp.ny(),scalar_type(1)/gpp.nxdotny()));
 
 1030             gmm::mult(gpp.vbase_lx(), gmm::transposed(auxNN2), auxLXN1);
 
 1032             gmm::mult(auxLXN1, gmm::transposed(gpp.vbase_ux()), Melem);
 
 1035             base_matrix auxLXUX(ndof_lx, ndof_ux);
 
 1036             gmm::mult(auxNN2, gpp.I_nxnx(), auxNN1);
 
 1037             gmm::mult(auxNN1, gmm::transposed(gpp.grad_phix_inv()), auxNN2);
 
 1038             gmm::mult(gpp.vbase_lx(), gmm::transposed(auxNN2), auxLXN2);
 
 1039             gmm::mult(auxLXN2, gmm::transposed(graddeltaunx), auxLXUX);
 
 1040             gmm::scale(auxLXUX, -g);
 
 1042             gmm::scale(Melem, -weight*alpha*FMULT/r);
 
 1043             mat_elem_assembly(M, I_lx, I_ux, Melem, *mf_lx, cvx, *mf_ux, cvx);
 
 1049             gmm::mult(auxLXN1, gmm::transposed(gpp.vbase_uy()), Melem);
 
 1050             gmm::scale(Melem, weight*alpha*FMULT/r);
 
 1051             mat_elem_assembly(M, I_lx, I_uy, Melem, *mf_lx, cvx, *mf_uy, cvy);
 
 1053             base_matrix I_gphiy0gphiyinv(N, N);
 
 1054             gmm::mult(gmm::scaled(gpp.grad_phiy0(), scalar_type(-1)),
 
 1055                       gpp.grad_phiy_inv(), I_gphiy0gphiyinv);
 
 1056             gmm::add(gmm::identity_matrix(), I_gphiy0gphiyinv);
 
 1061             gmm::mult(I_gphiy0gphiyinv, gpp.I_nxny(), auxNN1);
 
 1062             for (
size_type j = 0; j < N; ++j) auxNN1(j,j) -= scalar_type(1);
 
 1064             gmm::mult(gpp.vbase_lx(), gmm::transposed(auxNN2), auxLXN2);
 
 1066             gmm::mult(auxLXN2, gmm::transposed(gpp.vbase_ux()), Melem);
 
 1069             base_matrix auxLXUX(ndof_lx, ndof_ux);
 
 1070             gmm::mult(dVsF, I_gphiy0gphiyinv, auxNN1);
 
 1071             gmm::mult(auxNN1, gpp.I_nxny(), auxNN2);
 
 1072             gmm::mult(auxNN2, gmm::transposed(gpp.grad_phix_inv()), auxNN1);
 
 1073             gmm::mult(gpp.vbase_lx(), gmm::transposed(auxNN1), auxLXN1);
 
 1074             gmm::mult(auxLXN1, gmm::transposed(graddeltaunx), auxLXUX);
 
 1075             gmm::scale(auxLXUX, -g);
 
 1077             gmm::scale(Melem, weight*alpha*FMULT/r);
 
 1078             mat_elem_assembly(M, I_lx, I_ux, Melem, *mf_lx, cvx, *mf_ux, cvx);
 
 1084             gmm::mult(auxLXN2, gmm::transposed(gpp.vbase_uy()), Melem);
 
 1085             gmm::scale(Melem, -weight*alpha*FMULT/r);
 
 1086             mat_elem_assembly(M, I_lx, I_uy, Melem, *mf_lx, cvx, *mf_uy, cvy);
 
 1090             gmm::mult(gpp.vbase_lx(), gmm::transposed(dVsF), auxLXN1);
 
 1091             gmm::mult(auxLXN1, gmm::transposed(gpp.vbase_ux()), Melem);
 
 1092             gmm::scale(Melem, -weight*alpha*FMULT/r);
 
 1093             mat_elem_assembly(M, I_lx, I_ux, Melem, *mf_lx, cvx, *mf_ux, cvx);
 
 1099       if (version & model::BUILD_RHS) {
 
 1101         if (!(version & model::BUILD_MATRIX))
 
 1102           aug_friction(lambda, g, gpp.Vs(), nx, r, f_coeff, F);
 
 1104 #ifdef CONSIDER_TERM1 
 1107         gmm::mult(gpp.vbase_ux(), lambda, auxUX);
 
 1108         gmm::scale(auxUX, weight);
 
 1109         vec_elem_assembly(V, I_ux, auxUX, *mf_ux, cvx);
 
 1112 #ifdef CONSIDER_TERM2 
 1116           gmm::mult(gpp.vbase_uy(), lambda, auxUY);
 
 1117           gmm::scale(auxUY, -weight);
 
 1118           vec_elem_assembly(V, I_uy, auxUY, *mf_uy, cvy);
 
 1122 #ifdef CONSIDER_TERM3 
 1126         base_vector auxLX(ndof_lx);
 
 1127         gmm::mult(gpp.vbase_lx(), gmm::scaled(F, weight*FMULT/r), auxLX);
 
 1128         vec_elem_assembly(V, I_lx, auxLX, *mf_lx, cvx);
 
 1137   (
model &md, multi_contact_frame &mcf,
 
 1138    const std::string &dataname_r, 
const std::string &dataname_friction_coeff,
 
 1139    const std::string &dataname_alpha) {
 
 1141     bool with_friction = (dataname_friction_coeff.size() > 0);
 
 1143       = std::make_shared<integral_large_sliding_contact_brick>(mcf,
 
 1147     tl.push_back(model::term_description(
true, 
false));
 
 1149     model::varnamelist dl(1, dataname_r);
 
 1150     if (with_friction) dl.push_back(dataname_friction_coeff);
 
 1151     if (dataname_alpha.size()) dl.push_back(dataname_alpha);
 
 1153     model::varnamelist vl;
 
 1155     bool selfcontact = mcf.is_self_contact();
 
 1157     dal::bit_vector uvar, mvar;
 
 1158     for (
size_type i = 0; i < mcf.nb_boundaries(); ++i) {
 
 1159       size_type ind_u = mcf.ind_varname_of_boundary(i);
 
 1160       if (!(uvar.is_in(ind_u))) {
 
 1161         vl.push_back(mcf.varname(ind_u));
 
 1164       size_type ind_lambda = mcf.ind_multname_of_boundary(i);
 
 1166       if (selfcontact || mcf.is_slave_boundary(i))
 
 1167         GMM_ASSERT1(ind_lambda != 
size_type(-1), 
"Large sliding contact " 
 1168                     "brick: a multiplier should be associated to each slave " 
 1169                     "boundary in the multi_contact_frame object.");
 
 1170       if (ind_lambda != 
size_type(-1) && !(mvar.is_in(ind_lambda))) {
 
 1171         vl.push_back(mcf.multname(ind_lambda));
 
 1190 #if GETFEM_HAVE_MUPARSER_MUPARSER_H 
 1191 #include <muParser/muParser.h> 
 1192 #elif GETFEM_HAVE_MUPARSER_H 
 1193 #include <muParser.h> 
 1201   struct contact_frame {
 
 1204     scalar_type friction_coef;
 
 1205     std::vector<const model_real_plain_vector *> Us;
 
 1206     std::vector<model_real_plain_vector> ext_Us;
 
 1207     std::vector<const model_real_plain_vector *> lambdas;
 
 1208     std::vector<model_real_plain_vector> ext_lambdas;
 
 1209     struct contact_boundary {
 
 1216     std::vector<contact_boundary> contact_boundaries;
 
 1218     gmm::dense_matrix< model_real_sparse_matrix * > UU;
 
 1219     gmm::dense_matrix< model_real_sparse_matrix * > UL;
 
 1220     gmm::dense_matrix< model_real_sparse_matrix * > LU;
 
 1221     gmm::dense_matrix< model_real_sparse_matrix * > LL;
 
 1223     std::vector< model_real_plain_vector *> Urhs;
 
 1224     std::vector< model_real_plain_vector *> Lrhs;
 
 1228     std::vector<std::string> coordinates;
 
 1230 #if GETFEM_HAVE_MUPARSER_MUPARSER_H || GETFEM_HAVE_MUPARSER_H 
 1231     std::vector<mu::Parser> obstacles_parsers;
 
 1233     std::vector<std::string> obstacles;
 
 1234     std::vector<std::string> obstacles_velocities;
 
 1237                     const model_real_plain_vector &U) {
 
 1239       for (; i < Us.size(); ++i) 
if (Us[i] == &U) 
return i;
 
 1242       mfu.extend_vector(U, ext_U);
 
 1243       ext_Us.push_back(ext_U);
 
 1248                          const model_real_plain_vector &l) {
 
 1250       for (; i < lambdas.size(); ++i) 
if (lambdas[i] == &l) 
return i;
 
 1251       lambdas.push_back(&l);
 
 1253       mfl.extend_vector(l, ext_l);
 
 1254       ext_lambdas.push_back(ext_l);
 
 1258     void extend_vectors(
void) {
 
 1259       for (
size_type i = 0; i < contact_boundaries.size(); ++i) {
 
 1260         size_type ind_U = contact_boundaries[i].ind_U;
 
 1261         contact_boundaries[i].mfu->extend_vector(*(Us[ind_U]), ext_Us[ind_U]);
 
 1262         size_type ind_lambda = contact_boundaries[i].ind_lambda;
 
 1263         contact_boundaries[i].mflambda->extend_vector(*(lambdas[ind_lambda]),
 
 1264                                                       ext_lambdas[ind_lambda]);
 
 1270     { 
return *(contact_boundaries[n].mfu); }
 
 1272     { 
return *(contact_boundaries[n].mflambda); }
 
 1273     const model_real_plain_vector &disp_of_boundary(
size_type n)
 const 
 1274     { 
return ext_Us[contact_boundaries[n].ind_U]; }
 
 1275     const model_real_plain_vector &lambda_of_boundary(
size_type n)
 const 
 1276     { 
return ext_lambdas[contact_boundaries[n].ind_lambda]; }
 
 1278     { 
return contact_boundaries[n].region; }
 
 1280     { 
return *(UU(contact_boundaries[n].ind_U, contact_boundaries[m].ind_U)); }
 
 1282       return *(LU(contact_boundaries[n].ind_lambda,
 
 1283                   contact_boundaries[m].ind_U));
 
 1286       return *(UL(contact_boundaries[n].ind_U,
 
 1287                   contact_boundaries[m].ind_lambda));
 
 1290       return *(LL(contact_boundaries[n].ind_lambda,
 
 1291                   contact_boundaries[m].ind_lambda));
 
 1293     model_real_plain_vector &U_vector(
size_type n)
 const 
 1294     { 
return *(Urhs[contact_boundaries[n].ind_U]); }
 
 1295     model_real_plain_vector &L_vector(
size_type n)
 const 
 1296     { 
return *(Lrhs[contact_boundaries[n].ind_lambda]); }
 
 1298     contact_frame(
size_type NN) : N(NN), coordinates(N), pt_eval(N) {
 
 1299       if (N > 0) coordinates[0] = 
"x";
 
 1300       if (N > 1) coordinates[1] = 
"y";
 
 1301       if (N > 2) coordinates[2] = 
"z";
 
 1302       if (N > 3) coordinates[3] = 
"w";
 
 1303       GMM_ASSERT1(N <= 4, 
"Complete the definition for contact in " 
 1304                   "dimension greater than 4");
 
 1307     size_type add_obstacle(
const std::string &obs) {
 
 1309       obstacles.push_back(obs);
 
 1310       obstacles_velocities.push_back(
"");
 
 1311 #if GETFEM_HAVE_MUPARSER_MUPARSER_H || GETFEM_HAVE_MUPARSER_H 
 1314       obstacles_parsers.push_back(mu);
 
 1315       obstacles_parsers[ind].SetExpr(obstacles[ind]);
 
 1317         obstacles_parsers[ind].DefineVar(coordinates[k], &pt_eval[k]);
 
 1319       GMM_ASSERT1(
false, 
"You have to link muparser with getfem to deal " 
 1320                   "with rigid body obstacles");
 
 1326                            const model_real_plain_vector &U,
 
 1328                            const model_real_plain_vector &l,
 
 1330       contact_boundary cb;
 
 1334       cb.ind_U = add_U(mfu, U);
 
 1335       cb.ind_lambda = add_lambda(mfl, l);
 
 1336       size_type ind = contact_boundaries.size();
 
 1337       contact_boundaries.push_back(cb);
 
 1354   struct contact_elements {
 
 1361     std::vector<size_type> boundary_of_elements;
 
 1362     std::vector<size_type> ind_of_elements;
 
 1363     std::vector<size_type> face_of_elements;
 
 1364     std::vector<base_node> unit_normal_of_elements;
 
 1366     contact_elements(contact_frame &ccf) : cf(ccf) {}
 
 1368     bool add_point_contribution(
size_type boundary_num,
 
 1371                                 scalar_type weight, scalar_type f_coeff,
 
 1372                                 scalar_type r, model::build_version version);
 
 1376   void contact_elements::init(
void) {
 
 1377     fem_precomp_pool fppool;
 
 1380     element_boxes.clear();
 
 1381     unit_normal_of_elements.resize(0);
 
 1382     boundary_of_elements.resize(0);
 
 1383     ind_of_elements.resize(0);
 
 1384     face_of_elements.resize(0);
 
 1388     model_real_plain_vector coeff;
 
 1389     cf.extend_vectors();
 
 1390     for (
size_type i = 0; i < cf.contact_boundaries.size(); ++i) {
 
 1391       size_type bnum = cf.region_of_boundary(i);
 
 1392       const mesh_fem &mfu = cf.mfu_of_boundary(i);
 
 1393       const model_real_plain_vector &U = cf.disp_of_boundary(i);
 
 1394       const mesh &m = mfu.linked_mesh();
 
 1395       if (i == 0) N = m.dim();
 
 1396       GMM_ASSERT1(m.dim() == N,
 
 1397                   "Meshes are of mixed dimensions, cannot deal with that");
 
 1398       base_node val(N), bmin(N), bmax(N), n0(N), n(N), n_mean(N);
 
 1399       base_matrix grad(N,N);
 
 1400       mesh_region region = m.region(bnum);
 
 1401       GMM_ASSERT1(mfu.get_qdim() == N,
 
 1402                   "Wrong mesh_fem qdim to compute contact pairs");
 
 1404       dal::bit_vector points_already_interpolated;
 
 1405       std::vector<base_node> transformed_points(m.nb_max_points());
 
 1409         pfem pf_s = mfu.fem_of_element(cv);
 
 1412         bgeot::vectors_to_base_matrix
 
 1413           (G, mfu.linked_mesh().points_of_convex(cv));
 
 1415         pfem_precomp pfp = fppool(pf_s, &(pgt->geometric_nodes()));
 
 1416         fem_interpolation_context ctx(pgt, pfp, 
size_type(-1), G, cv);
 
 1421           size_type ind = m.ind_points_of_convex(cv)[ip];
 
 1424           if (!(points_already_interpolated.is_in(ind))) {
 
 1426             pf_s->interpolation(ctx, coeff, val, dim_type(N));
 
 1428             transformed_points[ind] = val;
 
 1429             points_already_interpolated.add(ind);
 
 1431             val = transformed_points[ind];
 
 1434           bool is_on_face = 
false;
 
 1436           for (
size_type k = 0; k < cvs->nb_points_of_face(v.f()); ++k)
 
 1437             if (cvs->ind_points_of_face(v.f())[k] == ip) is_on_face = 
true;
 
 1441             pf_s->interpolation_grad(ctx, coeff, grad, dim_type(N));
 
 1442             gmm::add(gmm::identity_matrix(), grad);
 
 1443             scalar_type J = bgeot::lu_inverse(&(*(grad.begin())), N);
 
 1444             if (J <= scalar_type(0)) GMM_WARNING1(
"Inverted element ! " << J);
 
 1445             gmm::mult(gmm::transposed(grad), n0, n);
 
 1455               bmin[k] = std::min(bmin[k], val[k]);
 
 1456               bmax[k] = std::max(bmax[k], val[k]);
 
 1461         GMM_ASSERT1(nb_pt_on_face,
 
 1462                     "This element has not vertex on considered face !");
 
 1466         scalar_type h = bmax[0] - bmin[0];
 
 1468           h = std::max(h, bmax[k] - bmin[k]);
 
 1470           { bmin[k] -= h; bmax[k] += h; }
 
 1473         element_boxes.add_box(bmin, bmax, unit_normal_of_elements.size());
 
 1475         unit_normal_of_elements.push_back(n_mean);
 
 1476         boundary_of_elements.push_back(i);
 
 1477         ind_of_elements.push_back(cv);
 
 1478         face_of_elements.push_back(v.f());
 
 1481     element_boxes.build_tree();
 
 1486   bool contact_elements::add_point_contribution
 
 1489    scalar_type , scalar_type r, model::build_version version) {
 
 1490     const mesh_fem &mfu = cf.mfu_of_boundary(boundary_num);
 
 1491     const mesh_fem &mfl = cf.mflambda_of_boundary(boundary_num);
 
 1492     const model_real_plain_vector &U = cf.disp_of_boundary(boundary_num);
 
 1493     const model_real_plain_vector &L = cf.lambda_of_boundary(boundary_num);
 
 1495     base_node x0 = ctxu.
xreal();
 
 1506     base_small_vector n(N), val(N), h(N);
 
 1507     base_matrix gradinv(N,N), grad(N,N), gradtot(N,N), G;
 
 1508     size_type cvnbdofu = mfu.nb_basic_dof_of_element(cv);
 
 1509     size_type cvnbdofl = mfl.nb_basic_dof_of_element(cv);
 
 1510     base_vector coeff(cvnbdofu);
 
 1512     ctxu.
pf()->interpolation(ctxu, coeff, val, dim_type(N));
 
 1513     base_node x = x0 + val;
 
 1515     ctxu.
pf()->interpolation_grad(ctxu, coeff, gradinv, dim_type(N));
 
 1516     gmm::add(gmm::identity_matrix(), gradinv);
 
 1517     scalar_type J = bgeot::lu_inverse(&(*(grad_inv.begin())), N); 
 
 1518     if (J <= scalar_type(0)) {
 
 1519       GMM_WARNING1(
"Inverted element !");
 
 1521       GMM_ASSERT1(!(version & model::BUILD_MATRIX), 
"Impossible to build " 
 1522                   "tangent matrix for large sliding contact");
 
 1523       if (version & model::BUILD_RHS) {
 
 1524         base_vector Velem(cvnbdofl);
 
 1525         for (
size_type i = 0; i < cvnbdofl; ++i) Velem[i] = 1E200;
 
 1526         vec_elem_assembly(cf.L_vector(boundary_num), Velem, mfl, cv);
 
 1531     gmm::mult(gmm::transposed(gradinv), n0, n);
 
 1538     bgeot::rtree::pbox_set bset;
 
 1539     element_boxes.find_boxes_at_point(x, bset);
 
 1541     if (noisy) cout << 
"Number of boxes found : " << bset.size() << endl;
 
 1548     bgeot::rtree::pbox_set::iterator it = bset.begin(), itnext;
 
 1549     for (; it != bset.end(); it = itnext) {
 
 1550       itnext = it; ++itnext;
 
 1551       if (gmm::vect_sp(unit_normal_of_elements[(*it)->id], n)
 
 1552           >= -scalar_type(1)/scalar_type(20)) bset.erase(it);
 
 1556       cout << 
"Number of boxes satisfying the unit normal criterion : " 
 1557            << bset.size() << endl;
 
 1567     std::vector<base_node> y0s;
 
 1568     std::vector<base_small_vector> n0_y0s;
 
 1569     std::vector<scalar_type> d0s;
 
 1570     std::vector<scalar_type> d1s;
 
 1571     std::vector<size_type> elt_nums;
 
 1572     std::vector<fem_interpolation_context> ctx_y0s;
 
 1573     for (; it != bset.end(); ++it) {
 
 1574       size_type boundary_num_y0 = boundary_of_elements[(*it)->id];
 
 1575       size_type cv_y0 = ind_of_elements[(*it)->id];
 
 1577       const mesh_fem &mfu_y0 = cf.mfu_of_boundary(boundary_num_y0);
 
 1578       pfem pf_s_y0 = mfu_y0.fem_of_element(cv_y0);
 
 1579       const model_real_plain_vector &U_y0
 
 1580         = cf.disp_of_boundary(boundary_num_y0);
 
 1581       const mesh &m_y0 = mfu_y0.linked_mesh();
 
 1588       for (; ind_dep_point < cvs_y0->nb_points(); ++ind_dep_point) {
 
 1589         bool is_on_face = 
false;
 
 1591              k < cvs_y0->nb_points_of_face(face_y0); ++k)
 
 1592           if (cvs_y0->ind_points_of_face(face_y0)[k]
 
 1593               == ind_dep_point) is_on_face = 
true;
 
 1594         if (!is_on_face) 
break;
 
 1596       GMM_ASSERT1(ind_dep_point < cvs_y0->nb_points(),
 
 1597                   "No interior point found !");
 
 1599       base_node y0_ref = pgt_y0->convex_ref()->points()[ind_dep_point];
 
 1603       bgeot::vectors_to_base_matrix(G, m_y0.points_of_convex(cv_y0));
 
 1605       fem_interpolation_context ctx_y0(pgt_y0, pf_s_y0, y0_ref, G, cv_y0);
 
 1610         pf_s_y0->interpolation(ctx_y0, coeff, val, dim_type(N));
 
 1611         val += ctx_y0.xreal() - x;
 
 1614         if (init_res < 1E-12) 
break;
 
 1615         if (newton_iter > 100) {
 
 1616           GMM_WARNING1(
"Newton has failed to invert transformation"); 
 
 1617           GMM_ASSERT1(!(version & model::BUILD_MATRIX), 
"Impossible to build " 
 1618                       "tangent matrix for large sliding contact");
 
 1619           if (version & model::BUILD_RHS) {
 
 1620             base_vector Velem(cvnbdofl);
 
 1621             for (
size_type i = 0; i < cvnbdofl; ++i) Velem[i] = 1E200;
 
 1622             vec_elem_assembly(cf.L_vector(boundary_num), Velem, mfl, cv);
 
 1627         pf_s_y0->interpolation_grad(ctx_y0, coeff, grad, dim_type(N));
 
 1628         gmm::add(gmm::identity_matrix(), grad);
 
 1631         std::vector<long> ipvt(N);
 
 1633         GMM_ASSERT1(!info, 
"Singular system, pivot = " << info); 
 
 1639         for (alpha = 1; 
alpha >= 1E-5; 
alpha/=scalar_type(2)) {
 
 1641           ctx_y0.set_xref(y0_ref - alpha*h);
 
 1642           pf_s_y0->interpolation(ctx_y0, coeff, val, dim_type(N));
 
 1643           val += ctx_y0.xreal() - x;
 
 1645           if (gmm::vect_norm2(val) < init_res) { ok = 
true; 
break; }
 
 1648           GMM_WARNING1(
"Line search has failed to invert transformation");
 
 1650         ctx_y0.set_xref(y0_ref);
 
 1654       base_node y0 = ctx_y0.xreal();
 
 1656       scalar_type d0_ref = pgt_y0->convex_ref()->is_in_face(face_y0, y0_ref);
 
 1660       scalar_type d1 = d0_ref; 
 
 1663       for (
short_type k = 0; k <  pgt_y0->structure()->nb_faces(); ++k) {
 
 1664         scalar_type dd = pgt_y0->convex_ref()->is_in_face(k, y0_ref);
 
 1665         if (dd > scalar_type(0) && dd > gmm::abs(d1)) { d1 = dd; ifd = k; }
 
 1670         if (gmm::abs(d1) < gmm::abs(d0)) d1 = d0;
 
 1677       if (noisy) cout << 
"gmm::vect_norm2(n0_y0) = " << 
gmm::vect_norm2(n0_y0) << endl;
 
 1679       if (noisy) cout << 
"autocontact status : x0 = " << x0 << 
" y0 = " << y0 << 
"  " <<  
gmm::vect_dist2(y0, x0) << 
" : " << d0*0.75 << 
" : " << d1*0.75 << endl;
 
 1680       if (noisy) cout << 
"n = " << n << 
" unit_normal_of_elements[(*it)->id] = " << unit_normal_of_elements[(*it)->id] << endl;
 
 1682       if (d0 < scalar_type(0)
 
 1684                && (gmm::vect_dist2(y0, x0) < gmm::abs(d1)*scalar_type(3)/scalar_type(4)))
 
 1685               || gmm::abs(d1) > 0.05)) {
 
 1686         if (noisy)  cout << 
"Eliminated x0 = " << x0 << 
" y0 = " << y0
 
 1687                          << 
" d0 = " << d0 << endl;
 
 1699       y0s.push_back(ctx_y0.xreal()); 
 
 1700       elt_nums.push_back((*it)->id);
 
 1703       ctx_y0s.push_back(ctx_y0);
 
 1705       n0_y0s.push_back(n0_y0);
 
 1707       if (noisy) cout << 
"dist0 = " << d0 << 
" dist0 * area = " 
 1708                       << pgt_y0->convex_ref()->is_in(y0_ref) << endl;
 
 1717     scalar_type d0 = 1E100, d1 = 1E100;
 
 1718     base_small_vector grad_obs(N);
 
 1721     for (
size_type k = 0; k < y0s.size(); ++k)
 
 1722       if (d1s[k] < d1) { d0 = d0s[k]; d1 = d1s[k]; ibound = k; state = 1; }
 
 1726 #if GETFEM_HAVE_MUPARSER_MUPARSER_H || GETFEM_HAVE_MUPARSER_H 
 1728     for (
size_type i = 0; i < cf.obstacles.size(); ++i) {
 
 1729       scalar_type d0_o = scalar_type(cf.obstacles_parsers[i].Eval());
 
 1730       if (d0_o < d0) { d0 = d0_o; irigid_obstacle = i; state = 2; }
 
 1733       scalar_type EPS = face_factor * 1E-9;
 
 1735         cf.pt_eval[k] += EPS;
 
 1737           (scalar_type(cf.obstacles_parsers[irigid_obstacle].Eval())-d0)/EPS;
 
 1738         cf.pt_eval[k] -= EPS;
 
 1743     if (cf.obstacles.size() > 0)
 
 1744       GMM_WARNING1(
"Rigid obstacles are ignored. Recompile with " 
 1745                    "muParser to account for rigid obstacles");
 
 1754     if (noisy && state == 1) {
 
 1755       cout  << 
"Point : " << x0 << 
" of boundary " << boundary_num
 
 1756             << 
" and element " << cv << 
" state = " << int(state);
 
 1757       if (version & model::BUILD_RHS) cout << 
" RHS";
 
 1758       if (version & model::BUILD_MATRIX) cout << 
" MATRIX";
 
 1761       size_type boundary_num_y0 = boundary_of_elements[elt_nums[ibound]];
 
 1762       const mesh_fem &mfu_y0 = cf.mfu_of_boundary(boundary_num_y0);
 
 1763       const mesh &m_y0 = mfu_y0.linked_mesh();
 
 1764       size_type cv_y0 = ind_of_elements[elt_nums[ibound]];
 
 1766       if (noisy) cout << 
" y0 = " << y0s[ibound] << 
" of element " 
 1767                             << cv_y0  << 
" of boundary " << boundary_num_y0 << endl;
 
 1768       for (
size_type k = 0; k < m_y0.nb_points_of_convex(cv_y0); ++k)
 
 1769         if (noisy) cout << 
"point " << k << 
" : " 
 1770                         << m_y0.points()[m_y0.ind_points_of_convex(cv_y0)[k]] << endl;
 
 1771       if (boundary_num_y0 == 0 && boundary_num == 0 && d0 < 0.0 && (version & model::BUILD_MATRIX)) GMM_ASSERT1(
false, 
"oups");
 
 1773     if (noisy) cout << 
" d0 = " << d0 << endl;
 
 1779     GMM_ASSERT1(ctxu.
pf()->target_dim() == 1 && ctxl.
pf()->target_dim() == 1,
 
 1780                 "Large sliding contact assembly procedure has to be adapted " 
 1781                 "to intrinsic vectorial elements. To be done.");
 
 1793     base_small_vector lambda(N);
 
 1795     ctxl.
pf()->interpolation(ctxl, coeff, lambda, dim_type(N));
 
 1796     GMM_ASSERT1(!(isnan(lambda[0])), 
"internal error");
 
 1801     scalar_type aux1, aux2;
 
 1806       base_small_vector zeta(N);
 
 1807       gmm::add(lambda, gmm::scaled(n, r*d0), zeta);
 
 1814       size_type boundary_num_y0 = 0, cv_y0 = 0, cvnbdofu_y0 = 0;
 
 1816         ctx_y0s[ibound].base_value(tu_y0);
 
 1817         boundary_num_y0 = boundary_of_elements[elt_nums[ibound]];
 
 1818         cv_y0 = ind_of_elements[elt_nums[ibound]];
 
 1819         cvnbdofu_y0 = cf.mfu_of_boundary(boundary_num_y0).nb_basic_dof_of_element(cv_y0);
 
 1821       const mesh_fem &mfu_y0 = (state == 1) ?
 
 1822                                cf.mfu_of_boundary(boundary_num_y0) : mfu;
 
 1824       if (version & model::BUILD_RHS) {
 
 1829         base_small_vector vecaux(N);
 
 1831         De_Saxce_projection(vecaux, n, scalar_type(0));
 
 1832         gmm::scale(vecaux, -scalar_type(1));
 
 1834         for (
size_type i = 0; i < cvnbdofl; ++i)
 
 1835           Velem[i] = tl[i/N] * vecaux[i%N] * weight/r;
 
 1836         vec_elem_assembly(cf.L_vector(boundary_num), Velem, mfl, cv);
 
 1840         for (
size_type i = 0; i < cvnbdofu; ++i)
 
 1841           Velem[i] = tu[i/N] * lambda[i%N] * weight;
 
 1842         vec_elem_assembly(cf.U_vector(boundary_num), Velem, mfu, cv);
 
 1846           for (
size_type i = 0; i < cvnbdofu_y0; ++i)
 
 1847             Velem[i] = -tu_y0[i/N] * lambda[i%N] * weight;
 
 1848           vec_elem_assembly(cf.U_vector(boundary_num_y0), Velem, mfu_y0, cv_y0);
 
 1852       if (version & model::BUILD_MATRIX) {
 
 1854         base_small_vector gradinv_n(N);
 
 1858         base_matrix pgrad(N,N), pgradn(N,N);
 
 1859         De_Saxce_projection_grad(zeta, n, scalar_type(0), pgrad);
 
 1860         De_Saxce_projection_gradn(zeta, n, scalar_type(0), pgradn);
 
 1862         base_small_vector pgrad_n(N), pgradn_n(N);
 
 1865         base_matrix gradinv_pgrad(N,N), gradinv_pgradn(N,N);
 
 1866         gmm::mult(gradinv, gmm::transposed(pgrad), gradinv_pgrad);
 
 1867         gmm::mult(gradinv, gmm::transposed(pgradn), gradinv_pgradn);
 
 1872         for (
size_type i = 0; i < cvnbdofl; i += N) {
 
 1873           aux1 = -tl[i/N] * weight/r;
 
 1874           for (
size_type j = 0; j < cvnbdofl; j += N) {
 
 1875             aux2 = aux1 * tl[j/N];
 
 1876             for (
size_type k = 0; k < N; k++) Melem(i+k,j+k) = aux2;
 
 1880         for (
size_type i = 0, ii = 0; i < cvnbdofl; ++i, ii = i%N)
 
 1881           for (
size_type j = 0, jj = 0; j < cvnbdofl; ++j, jj = j%N)
 
 1882             Melem(i,j) += tl[i/N] * tl[j/N] * pgrad(ii,jj) * weight/r;
 
 1883         mat_elem_assembly(cf.LL_matrix(boundary_num, boundary_num),
 
 1884                           Melem, mfl, cv, mfl, cv);
 
 1889         for (
size_type i = 0; i < cvnbdofu; i += N) {
 
 1890           aux1 = -tu[i/N] * weight;
 
 1891           for (
size_type j = 0; j < cvnbdofl; j += N) {
 
 1892             aux2 = aux1 * tl[j/N];
 
 1893             for (
size_type k = 0; k < N; k++) Melem(i+k,j+k) = aux2;
 
 1896         mat_elem_assembly(cf.UL_matrix(boundary_num, boundary_num),
 
 1897                           Melem, mfu, cv, mfl, cv);
 
 1903         for (
size_type i = 0, ii = 0; i < cvnbdofl; ++i, ii = i%N)
 
 1904           for (
size_type j = 0, jj = 0; j < cvnbdofu; ++j, jj = j%N) {
 
 1905             aux1 = aux2 = scalar_type(0);
 
 1907               aux1 += tgradu[j/N+N*k] * gradinv_n[k];
 
 1908               aux2 += tgradu[j/N+N*k] * gradinv_pgrad(k,ii);
 
 1910             Melem(i,j) = d0 * tl[i/N] * (pgrad_n[ii] * aux1 - aux2) * n[jj] * weight;
 
 1916         for (
size_type i = 0, ii = 0; i < cvnbdofl; ++i, ii = i%N)
 
 1917           for (
size_type j = 0, jj = 0; j < cvnbdofu; ++j, jj = j%N) {
 
 1918             aux1 = aux2 = scalar_type(0);
 
 1920               aux1 += tgradu[j/N+N*k] * gradinv_n[k];
 
 1921               aux2 += tgradu[j/N+N*k] * gradinv_pgradn(k,ii);
 
 1923             Melem(i,j) += tl[i/N] * (pgradn_n[ii] * aux1 - aux2) * n[jj] * weight / r;
 
 1925         mat_elem_assembly(cf.LU_matrix(boundary_num, boundary_num),
 
 1926                           Melem, mfl, cv, mfu, cv);
 
 1931           base_tensor tgradu_y0;
 
 1932           ctx_y0s[ibound].grad_base_value(tgradu_y0);
 
 1934           base_matrix gradinv_y0(N,N);
 
 1935           base_small_vector ntilde_y0(N);
 
 1937             base_matrix grad_y0(N,N);
 
 1938             base_vector coeff_y0(cvnbdofu_y0);
 
 1939             const model_real_plain_vector &U_y0
 
 1940               = cf.disp_of_boundary(boundary_num_y0);
 
 1942             ctx_y0s[ibound].pf()->interpolation_grad(ctx_y0s[ibound], coeff_y0,
 
 1943                                                    grad_y0, dim_type(N));
 
 1944             gmm::add(gmm::identity_matrix(), grad_y0);
 
 1948             gmm::mult(gmm::transposed(gradinv_y0), n0_y0s[ibound], ntilde_y0); 
 
 1953           for (
size_type i = 0; i < cvnbdofu_y0; i += N) {
 
 1954             aux1 = tu_y0[i/N] * weight;
 
 1955             for (
size_type j = 0; j < cvnbdofl; j += N) {
 
 1956               aux2 = aux1 * tl[j/N];
 
 1957               for (
size_type k = 0; k < N; k++) Melem(i+k,j+k) = aux2;
 
 1960           mat_elem_assembly(cf.UL_matrix(boundary_num_y0, boundary_num),
 
 1961                             Melem, mfu_y0, cv_y0, mfl, cv);
 
 1969           for (
size_type i = 0, ii = 0; i < cvnbdofu_y0; ++i, ii = i%N)
 
 1970             for (
size_type j = 0, jj = 0; j < cvnbdofu; ++j, jj = j%N) {
 
 1971               aux1 = scalar_type(0);
 
 1973                 aux1 += tgradu_y0[i/N+N*k]* gradinv_y0(k,jj);
 
 1974               Melem(i,j) = lambda[ii] * aux1 * tu[j/N] * weight;
 
 1976           mat_elem_assembly(cf.UU_matrix(boundary_num_y0, boundary_num),
 
 1977                             Melem, mfu_y0, cv_y0, mfu, cv);
 
 1982           for (
size_type i = 0, ii = 0; i < cvnbdofu_y0; ++i, ii = i%N)
 
 1983             for (
size_type j = 0, jj = 0; j < cvnbdofu_y0; ++j, jj = j%N) {
 
 1984               aux1 = scalar_type(0);
 
 1986                 aux1 += tgradu_y0[i/N+N*k] * gradinv_y0(k,jj);
 
 1987               Melem(i,j) = - lambda[ii] * aux1 * tu_y0[j/N] * weight;
 
 1989           mat_elem_assembly(cf.UU_matrix(boundary_num_y0, boundary_num_y0),
 
 1990                             Melem, mfu_y0, cv_y0, mfu_y0, cv_y0);
 
 1995           for (
size_type i = 0; i < cvnbdofl; ++i) {
 
 1996             aux1 = tl[i/N] * pgrad_n[i%N] * weight;
 
 1997             for (
size_type j = 0; j < cvnbdofu_y0; ++j)
 
 1998               Melem(i,j) = - aux1 * tu_y0[j/N] * ntilde_y0[j%N];
 
 2000           mat_elem_assembly(cf.LU_matrix(boundary_num, boundary_num_y0),
 
 2001                             Melem, mfl, cv, mfu_y0, cv_y0);
 
 2006           for (
size_type i = 0; i < cvnbdofl; ++i) {
 
 2007             aux1 = tl[i/N] * pgrad_n[i%N] * weight;
 
 2008             for (
size_type j = 0; j < cvnbdofu; ++j)
 
 2009               Melem(i,j) = aux1 * tu[j/N] * ntilde_y0[j%N];
 
 2016           for (
size_type i = 0; i < cvnbdofl; ++i) {
 
 2017             aux1 = tl[i/N] * pgrad_n[i%N] * weight;
 
 2018             for (
size_type j = 0; j < cvnbdofu; ++j)
 
 2019               Melem(i,j) = aux1 * tu[j/N] * grad_obs[j%N];
 
 2022         mat_elem_assembly(cf.LU_matrix(boundary_num, boundary_num),
 
 2023                           Melem, mfl, cv, mfu, cv);
 
 2030       if (version & model::BUILD_RHS) {
 
 2032         for (
size_type i = 0; i < cvnbdofl; ++i)
 
 2033           Velem[i] = tl[i/N] * lambda[i%N] * weight/r;
 
 2034         vec_elem_assembly(cf.L_vector(boundary_num), Velem, mfl, cv);
 
 2038       if (version & model::BUILD_MATRIX) {
 
 2040         for (
size_type i = 0; i < cvnbdofl; i += N) {
 
 2041           aux1 = -tl[i/N] * weight/r;
 
 2042           for (
size_type j = 0; j < cvnbdofl; j += N) {
 
 2043             aux2 = aux1 * tl[j/N];
 
 2044             for (
size_type k = 0; k < N; k++) Melem(i+k,j+k) = aux2;
 
 2047         mat_elem_assembly(cf.LL_matrix(boundary_num, boundary_num),
 
 2048                           Melem, mfl, cv, mfl, cv);
 
 2059   struct integral_large_sliding_contact_brick_field_extension : 
public virtual_brick {
 
 2062     struct contact_boundary {
 
 2064       std::string varname;
 
 2065       std::string multname;
 
 2069     std::vector<contact_boundary> boundaries;
 
 2070     std::vector<std::string> obstacles;
 
 2072     void add_boundary(
const std::string &varn, 
const std::string &multn,
 
 2074       contact_boundary cb;
 
 2075       cb.region = region; cb.varname = varn; cb.multname = multn; cb.mim=&mim;
 
 2076       boundaries.push_back(cb);
 
 2079     void add_obstacle(
const std::string &obs)
 
 2080     { obstacles.push_back(obs); }
 
 2082     void build_contact_frame(
const model &md, contact_frame &cf)
 const {
 
 2083       for (
size_type i = 0; i < boundaries.size(); ++i) {
 
 2084         const contact_boundary &cb = boundaries[i];
 
 2085         cf.add_boundary(md.mesh_fem_of_variable(cb.varname),
 
 2086                         md.real_variable(cb.varname),
 
 2087                         md.mesh_fem_of_variable(cb.multname),
 
 2088                         md.real_variable(cb.multname), cb.region);
 
 2090       for (
size_type i = 0; i < obstacles.size(); ++i)
 
 2091         cf.add_obstacle(obstacles[i]);
 
 2095     virtual void asm_real_tangent_terms(
const model &md, 
size_type ,
 
 2096                                         const model::varnamelist &vl,
 
 2097                                         const model::varnamelist &dl,
 
 2098                                         const model::mimlist &mims,
 
 2099                                         model::real_matlist &matl,
 
 2100                                         model::real_veclist &vecl,
 
 2101                                         model::real_veclist &,
 
 2103                                         build_version version) 
const;
 
 2105     integral_large_sliding_contact_brick_field_extension() {
 
 2106       set_flags(
"Integral large sliding contact brick",
 
 2117   void integral_large_sliding_contact_brick_field_extension::asm_real_tangent_terms
 
 2118   (
const model &md, 
size_type , 
const model::varnamelist &vl,
 
 2119    const model::varnamelist &dl, 
const model::mimlist &,
 
 2120    model::real_matlist &matl, model::real_veclist &vecl,
 
 2122    build_version version)
 const {
 
 2124     fem_precomp_pool fppool;
 
 2126     size_type N = md.mesh_fem_of_variable(vl[0]).linked_mesh().dim();
 
 2127     contact_frame cf(N);
 
 2128     build_contact_frame(md, cf);
 
 2130     size_type Nvar = vl.size(), Nu = cf.Urhs.size(), Nl = cf.Lrhs.size();
 
 2131     GMM_ASSERT1(Nvar == Nu+Nl, 
"Wrong size of variable list for integral " 
 2132                 "large sliding contact brick");
 
 2133     GMM_ASSERT1(matl.size() == Nvar*Nvar, 
"Wrong size of terms for " 
 2134                 "integral large sliding contact brick");
 
 2136     if (version & model::BUILD_MATRIX) {
 
 2140           if (i <  Nu && j <  Nu) cf.UU(i,j)       = &(matl[i*Nvar+j]);
 
 2141           if (i >= Nu && j <  Nu) cf.LU(i-Nu,j)    = &(matl[i*Nvar+j]);
 
 2142           if (i <  Nu && j >= Nu) cf.UL(i,j-Nu)    = &(matl[i*Nvar+j]);
 
 2143           if (i >= Nu && j >= Nu) cf.LL(i-Nu,j-Nu) = &(matl[i*Nvar+j]);
 
 2146     if (version & model::BUILD_RHS) {
 
 2147       for (
size_type i = 0; i < vl.size(); ++i) {
 
 2148         if (i < Nu) cf.Urhs[i] = &(vecl[i*Nvar]);
 
 2149         else cf.Lrhs[i-Nu] = &(vecl[i*Nvar]);
 
 2154     GMM_ASSERT1(dl.size() == 2, 
"Wrong number of data for integral large " 
 2155                 "sliding contact brick");
 
 2157     const model_real_plain_vector &vr = md.real_variable(dl[0]);
 
 2158     GMM_ASSERT1(gmm::vect_size(vr) == 1, 
"Parameter r should be a scalar");
 
 2160     const model_real_plain_vector &f_coeff = md.real_variable(dl[1]);
 
 2161     GMM_ASSERT1(gmm::vect_size(f_coeff) == 1,
 
 2162                 "Friction coefficient should be a scalar");
 
 2164     contact_elements ce(cf);
 
 2167     for (
size_type bnum = 0; bnum < boundaries.size(); ++bnum) {
 
 2168       mesh_region rg(boundaries[bnum].region);
 
 2169       const mesh_fem &mfu=md.mesh_fem_of_variable(boundaries[bnum].varname);
 
 2170       const mesh_fem &mfl=md.mesh_fem_of_variable(boundaries[bnum].multname);
 
 2171       const mesh_im &mim = *(boundaries[bnum].mim);
 
 2172       const mesh &m = mfu.linked_mesh();
 
 2173       mfu.linked_mesh().intersect_with_mpi_region(rg);
 
 2179         pfem pf_s = mfu.fem_of_element(cv);
 
 2180         pfem pf_sl = mfl.fem_of_element(cv);
 
 2181         pintegration_method pim = mim.int_method_of_element(cv);
 
 2182         bgeot::vectors_to_base_matrix(G, m.points_of_convex(cv));
 
 2185           = fppool(pf_s, pim->approx_method()->pintegration_points());
 
 2187           = fppool(pf_sl, pim->approx_method()->pintegration_points());
 
 2188         fem_interpolation_context ctxu(pgt, pfpu, 
size_type(-1), G, cv, v.f());
 
 2189         fem_interpolation_context ctxl(pgt, pfpl, 
size_type(-1), G, cv, v.f());
 
 2192              k < pim->approx_method()->nb_points_on_face(v.f()); ++k) {
 
 2194             = pim->approx_method()->ind_first_point_on_face(v.f()) + k;
 
 2197           if (!(ce.add_point_contribution
 
 2198                (bnum, ctxu, ctxl,pim->approx_method()->coeff(ind),
 
 2199                 f_coeff[0], vr[0], version))) 
return;
 
 2209   size_type add_integral_large_sliding_contact_brick_field_extension
 
 2210   (model &md, 
const mesh_im &mim, 
const std::string &varname_u,
 
 2211    const std::string &multname, 
const std::string &dataname_r,
 
 2212    const std::string &dataname_friction_coeff, 
size_type region) {
 
 2215       =std::make_shared<integral_large_sliding_contact_brick_field_extension>();
 
 2217     pbr->add_boundary(varname_u, multname, mim, region);
 
 2220     tl.push_back(model::term_description(varname_u, varname_u, 
false));
 
 2221     tl.push_back(model::term_description(varname_u, multname,  
false));
 
 2222     tl.push_back(model::term_description(multname,  varname_u, 
false));
 
 2223     tl.push_back(model::term_description(multname,  multname,  
false));
 
 2225     model::varnamelist dl(1, dataname_r);
 
 2226     dl.push_back(dataname_friction_coeff);
 
 2228     model::varnamelist vl(1, varname_u);
 
 2229     vl.push_back(multname);
 
 2231     return md.add_brick(pbr, vl, dl, tl, model::mimlist(1, &mim), region);
 
 2235   void add_boundary_to_large_sliding_contact_brick
 
 2236   (model &md, 
size_type indbrick, 
const mesh_im &mim,
 
 2237    const std::string &varname_u, 
const std::string &multname,
 
 2239     dim_type N = md.mesh_fem_of_variable(varname_u).linked_mesh().dim();
 
 2240     pbrick pbr = md.brick_pointer(indbrick);
 
 2241     md.touch_brick(indbrick);
 
 2242     integral_large_sliding_contact_brick_field_extension *p
 
 2243       = 
dynamic_cast<integral_large_sliding_contact_brick_field_extension *
> 
 2244       (
const_cast<virtual_brick *
>(pbr.get()));
 
 2245     GMM_ASSERT1(p, 
"Wrong type of brick");
 
 2246     p->add_boundary(varname_u, multname, mim, region);
 
 2247     md.add_mim_to_brick(indbrick, mim);
 
 2249     contact_frame cf(N);
 
 2250     p->build_contact_frame(md, cf);
 
 2252     model::varnamelist vl;
 
 2254     for (
size_type i = 0; i < cf.contact_boundaries.size(); ++i)
 
 2255       if (cf.contact_boundaries[i].ind_U >= nvaru)
 
 2256         { vl.push_back(p->boundaries[i].varname); ++nvaru; }
 
 2259     for (
size_type i = 0; i < cf.contact_boundaries.size(); ++i)
 
 2260       if (cf.contact_boundaries[i].ind_lambda >= nvarl)
 
 2261         { vl.push_back(p->boundaries[i].multname); ++nvarl; }
 
 2262     md.change_variables_of_brick(indbrick, vl);
 
 2265     for (
size_type i = 0; i < vl.size(); ++i)
 
 2266       for (
size_type j = 0; j < vl.size(); ++j)
 
 2267         tl.push_back(model::term_description(vl[i], vl[j], 
false));
 
 2269     md.change_terms_of_brick(indbrick, tl);
 
 2273   (model &md, 
size_type indbrick, 
const std::string &obs) { 
 
 2274     pbrick pbr = md.brick_pointer(indbrick);
 
 2275     md.touch_brick(indbrick);
 
 2276     integral_large_sliding_contact_brick_field_extension *p
 
 2277       = 
dynamic_cast<integral_large_sliding_contact_brick_field_extension *
> 
 2278       (
const_cast<virtual_brick *
>(pbr.get()));
 
 2279     GMM_ASSERT1(p, 
"Wrong type of brick");
 
 2280     p->add_obstacle(obs);
 
 2292   struct intergral_large_sliding_contact_brick_raytracing
 
 2293     : 
public virtual_brick {
 
 2295     struct contact_boundary {
 
 2297       std::string varname_u;
 
 2298       std::string varname_lambda;
 
 2299       std::string varname_w;
 
 2307     std::vector<contact_boundary> boundaries;
 
 2308     std::string transformation_name;
 
 2309     std::string u_group;
 
 2310     std::string w_group;
 
 2311     std::string friction_coeff;
 
 2313     std::string augmentation_param;
 
 2314     model::varnamelist vl, dl;
 
 2317     bool sym_version, frame_indifferent;
 
 2319     void add_contact_boundary(model &md, 
const mesh_im &mim, 
size_type region,
 
 2320                               bool is_master, 
bool is_slave,
 
 2321                               const std::string &u,
 
 2322                               const std::string &lambda,
 
 2323                               const std::string &w = 
"") {
 
 2324       std::string test_u = 
"Test_" + sup_previous_and_dot_to_varname(u);
 
 2325       std::string test_u_group = 
"Test_" + sup_previous_and_dot_to_varname(u_group);
 
 2326       std::string test_lambda = 
"Test_" + sup_previous_and_dot_to_varname(lambda);
 
 2327       GMM_ASSERT1(is_slave || is_master, 
"The contact boundary should be " 
 2328                   "either master, slave or both");
 
 2329       const mesh_fem *mf = md.pmesh_fem_of_variable(u);
 
 2330       GMM_ASSERT1(mf, 
"The displacement variable should be a f.e.m. one");
 
 2331       GMM_ASSERT1(&(mf->linked_mesh()) == &(mim.linked_mesh()),
 
 2332                   "The displacement variable and the integration method " 
 2333                   "should share the same mesh");
 
 2335         const mesh_fem *mf_l = md.pmesh_fem_of_variable(lambda);
 
 2336         const im_data *mimd_l = md.pim_data_of_variable(lambda);
 
 2337         GMM_ASSERT1(mf_l || mimd_l,
 
 2338                     "The multiplier variable should be defined on a " 
 2339                     "mesh_fem or an im_data");
 
 2340         GMM_ASSERT1(&(mf_l ? mf_l->linked_mesh() : mimd_l->linked_mesh())
 
 2341                     == &(mim.linked_mesh()),
 
 2342                     "The displacement and the multiplier fields " 
 2343                     "should be defined on the same mesh");
 
 2347         const mesh_fem *mf2 =  md.pmesh_fem_of_variable(w);
 
 2348         GMM_ASSERT1(!mf2 || &(mf2->linked_mesh()) == &(mf->linked_mesh()),
 
 2349                     "The data for the sliding velocity should be defined on " 
 2350                     " the same mesh as the displacement variable");
 
 2353       for (
size_type i = 0; i < boundaries.size(); ++i) {
 
 2354         const contact_boundary &cb = boundaries[i];
 
 2355         if (&(md.mesh_fem_of_variable(cb.varname_u).linked_mesh())
 
 2356             == &(mf->linked_mesh()) && cb.region == region)
 
 2357           GMM_ASSERT1(
false, 
"This contact boundary has already been added");
 
 2361           (md, transformation_name, mf->linked_mesh(), u_group, region);
 
 2364           (md, transformation_name, mf->linked_mesh(), u_group, region);
 
 2366       boundaries.push_back(contact_boundary());
 
 2367       contact_boundary &cb = boundaries.back();
 
 2370       if (is_slave) cb.varname_lambda = lambda;
 
 2372       cb.is_master = is_master;
 
 2373       cb.is_slave = is_slave;
 
 2376         std::string n, n0, Vs, g, Y;
 
 2377         n = 
"Transformed_unit_vector(Grad_"+u+
", Normal)";
 
 2378         n0 = 
"Transformed_unit_vector(Grad_"+w+
", Normal)";
 
 2385         Y = 
"Interpolate(X,"+transformation_name+
")";
 
 2386         g = 
"("+Y+
"+Interpolate("+u_group+
","+transformation_name+
")-X-"+u+
")."+n;
 
 2387         Vs = 
"("+u+
"-Interpolate("+u_group+
","+transformation_name+
")";
 
 2389           Vs += 
"-"+w+
"+Interpolate("+w_group+
","+transformation_name+
")";
 
 2390           if (frame_indifferent)
 
 2391             Vs += 
"+("+g+
")*("+n+
"-"+n0+
")";
 
 2395         std::string coupled_projection_def =
 
 2396           "Coulomb_friction_coupled_projection(" 
 2397           + lambda+
","+n+
","+Vs+
","+g+
","+friction_coeff+
"," 
 2398           + augmentation_param+
")";
 
 2404         g = 
"(Interpolate(X,"+transformation_name+
")-X-"+u+
")."+n;
 
 2405         if (frame_indifferent && w.size())
 
 2406           Vs = 
"("+u+
"-"+w+
"+"+g+
"*("+n+
"-"+n0+
"))*"+
alpha;
 
 2408           Vs = 
"("+u+(w.size() ? (
"-"+w):
"")+
")*"+
alpha;
 
 2410         std::string coupled_projection_rig =
 
 2411           "Coulomb_friction_coupled_projection(" 
 2412           + lambda+
","+n+
","+Vs+
","+g+
","+ friction_coeff+
"," 
 2413           + augmentation_param+
")";
 
 2417           (sym_version ? 
"" : (
"-"+lambda+
"."+test_u))
 
 2420           + (sym_version ? (
"+ Interpolate_filter("+transformation_name+
",-" 
 2421                             +coupled_projection_def+
"."+test_u+
",1)") : 
"")
 
 2422           + (sym_version ? (
"+ Interpolate_filter("+transformation_name+
",-" 
 2423                             +coupled_projection_rig+
"."+test_u+
",2)") : 
"")
 
 2429           + 
"+ Interpolate_filter("+transformation_name+
"," 
 2430           + (sym_version ? coupled_projection_def : lambda)
 
 2431           + 
".Interpolate("+test_u_group+
"," + transformation_name+
"), 1)" 
 2433           + 
"-(1/"+augmentation_param+
")*"+lambda+
"."+test_lambda
 
 2436           + 
"+ Interpolate_filter("+transformation_name+
"," 
 2437           + 
"(1/"+augmentation_param+
")*"+ coupled_projection_rig
 
 2438           + 
"."+test_lambda+
", 2)" 
 2441           + 
"+ Interpolate_filter("+transformation_name+
"," 
 2442           + 
"(1/"+augmentation_param+
")*" + coupled_projection_def + 
"." 
 2443           + test_lambda+
", 1)";
 
 2447     virtual void asm_real_tangent_terms(
const model &md, 
size_type ,
 
 2448                                         const model::varnamelist &,
 
 2449                                         const model::varnamelist &,
 
 2450                                         const model::mimlist &,
 
 2451                                         model::real_matlist &,
 
 2452                                         model::real_veclist &,
 
 2453                                         model::real_veclist &,
 
 2455                                         build_version)
 const {
 
 2460       for (
const contact_boundary &cb : boundaries) {
 
 2462           md.add_generic_expression(cb.expr, *(cb.mim), cb.region);
 
 2467     intergral_large_sliding_contact_brick_raytracing
 
 2468     (
const std::string &r,
 
 2469      const std::string &f_coeff,
const std::string &ug,
 
 2470      const std::string &wg, 
const std::string &tr,
 
 2471      const std::string &alpha_ = 
"1", 
bool sym_v = 
false,
 
 2472      bool frame_indiff = 
false) {
 
 2473       transformation_name = tr;
 
 2474       u_group = ug; w_group = wg;
 
 2475       friction_coeff = f_coeff;
 
 2477       augmentation_param = r;
 
 2478       sym_version = sym_v;
 
 2479       frame_indifferent = frame_indiff;
 
 2481       set_flags(
"Integral large sliding contact bick raytracing",
 
 2492     pbrick pbr = md.brick_pointer(indbrick);
 
 2493     intergral_large_sliding_contact_brick_raytracing *p
 
 2494       = 
dynamic_cast<intergral_large_sliding_contact_brick_raytracing *
> 
 2496     GMM_ASSERT1(p, 
"Wrong type of brick");
 
 2497     return p->transformation_name;
 
 2502     pbrick pbr = md.brick_pointer(indbrick);
 
 2503     intergral_large_sliding_contact_brick_raytracing *p
 
 2504       = 
dynamic_cast<intergral_large_sliding_contact_brick_raytracing *
> 
 2506     GMM_ASSERT1(p, 
"Wrong type of brick");
 
 2512     pbrick pbr = md.brick_pointer(indbrick);
 
 2513     intergral_large_sliding_contact_brick_raytracing *p
 
 2514       = 
dynamic_cast<intergral_large_sliding_contact_brick_raytracing *
> 
 2516     GMM_ASSERT1(p, 
"Wrong type of brick");
 
 2522     pbrick pbr = md.brick_pointer(indbrick);
 
 2523     intergral_large_sliding_contact_brick_raytracing *p
 
 2524       = 
dynamic_cast<intergral_large_sliding_contact_brick_raytracing *
> 
 2526     GMM_ASSERT1(p, 
"Wrong type of brick");
 
 2528       (md, p->transformation_name, expr, N);
 
 2533    bool is_master, 
bool is_slave, 
const std::string &u,
 
 2534    const std::string &lambda, 
const std::string &w) {
 
 2536     pbrick pbr = md.brick_pointer(indbrick);
 
 2537     intergral_large_sliding_contact_brick_raytracing *p
 
 2538       = 
dynamic_cast<intergral_large_sliding_contact_brick_raytracing *
> 
 2540     GMM_ASSERT1(p, 
"Wrong type of brick");
 
 2542     bool found_u = 
false, found_lambda = 
false;
 
 2543     for (
const auto & v : p->vl) {
 
 2544       if (v.compare(u) == 0) found_u = 
true;
 
 2545       if (v.compare(lambda) == 0) found_lambda = 
true;
 
 2547     if (!found_u) p->vl.push_back(u);
 
 2548     GMM_ASSERT1(!is_slave || lambda.size(),
 
 2549                 "You should define a multiplier on each slave boundary");
 
 2550     if (is_slave && !found_lambda) p->vl.push_back(lambda);
 
 2551     if (!found_u || (is_slave && !found_lambda))
 
 2554     std::vector<std::string> ug = md.variable_group(p->u_group);
 
 2556     for (
const auto &uu : ug)
 
 2557       if (uu.compare(u) == 0) { found_u = 
true; 
break; }
 
 2560       md.define_variable_group(p->u_group, ug);
 
 2564       bool found_w = 
false;
 
 2565       for (
const auto &ww : p->dl)
 
 2566         if (ww.compare(w) == 0) { found_w = 
true; 
break; }
 
 2571       std::vector<std::string> wg = md.variable_group(p->w_group);
 
 2573       for (
const auto &ww : wg)
 
 2574         if (ww.compare(w) == 0) { found_w = 
true; 
break; }
 
 2577         md.define_variable_group(p->w_group, wg);
 
 2581     bool found_mim = 
false;
 
 2582     for (
const auto &pmim : p->ml)
 
 2583       if (pmim == &mim) { found_mim = 
true; 
break; }
 
 2585       p->ml.push_back(&mim);
 
 2589     p->add_contact_boundary(md, mim, region, is_master, is_slave,
 
 2594   (
model &md, 
const std::string &augm_param,
 
 2595    scalar_type release_distance, 
const std::string &f_coeff, 
 
 2596    const std::string &alpha, 
bool sym_v, 
bool frame_indifferent) {
 
 2598     char ugroupname[50], wgroupname[50], transname[50];
 
 2599     for (
int i = 0; i < 10000; ++i) {
 
 2600       snprintf(ugroupname, 49, 
"disp__group_raytracing_%d", i);
 
 2601       if (!(md.variable_group_exists(ugroupname))
 
 2605     md.define_variable_group(ugroupname, std::vector<std::string>());
 
 2607     for (
int i = 0; i < 10000; ++i) {
 
 2608       snprintf(wgroupname, 49, 
"w__group_raytracing_%d", i);
 
 2609       if (!(md.variable_group_exists(wgroupname))
 
 2613     md.define_variable_group(wgroupname, std::vector<std::string>());
 
 2615     for (
int i = 0; i < 10000; ++i) {
 
 2616       snprintf(transname, 49, 
"trans__raytracing_%d", i);
 
 2622     model::varnamelist vl, dl;
 
 2627     auto p = std::make_shared<intergral_large_sliding_contact_brick_raytracing>
 
 2628       (augm_param, f_coeff, ugroupname, wgroupname, transname, alpha,
 
 2629        sym_v, frame_indifferent);
 
 2633     return md.
add_brick(pbr, p->vl, p->dl, model::termlist(), model::mimlist(),
 
 2639  struct Nitsche_large_sliding_contact_brick_raytracing
 
 2640     : 
public virtual_brick {
 
 2642     struct contact_boundary {
 
 2644       std::string varname_u;
 
 2645       std::string sigma_u;
 
 2646       std::string varname_w;
 
 2655     std::vector<contact_boundary> boundaries;
 
 2656     std::string transformation_name;
 
 2657     std::string u_group;
 
 2658     std::string w_group;
 
 2659     std::string friction_coeff;
 
 2661     std::string Nitsche_param;
 
 2662     model::varnamelist vl, dl;
 
 2665     bool  sym_version, frame_indifferent, unbiased;
 
 2667     void add_contact_boundary(model &md, 
const mesh_im &mim, 
size_type region,
 
 2668                               bool is_master, 
bool is_slave, 
bool is_unbiased,
 
 2669                               const std::string &u,
 
 2670                               const std::string &sigma_u,
 
 2671                               const std::string &w = 
"") {
 
 2672       std::string test_u = 
"Test_" + sup_previous_and_dot_to_varname(u);
 
 2673       std::string test_u_group = 
"Test_" + sup_previous_and_dot_to_varname(u_group);
 
 2674       GMM_ASSERT1(is_slave || is_master, 
"The contact boundary should be " 
 2675                   "either master, slave or both");
 
 2677                         GMM_ASSERT1((is_slave && is_master), 
"The contact boundary should be " 
 2678                           "both master and slave for the unbiased version");
 
 2679                         is_slave=
true; is_master=
true;
 
 2681       const mesh_fem *mf = md.pmesh_fem_of_variable(u);
 
 2682       GMM_ASSERT1(mf, 
"The displacement variable should be a f.e.m. one");
 
 2683       GMM_ASSERT1(&(mf->linked_mesh()) == &(mim.linked_mesh()),
 
 2684                   "The displacement variable and the integration method " 
 2685                   "should share the same mesh");
 
 2688         const mesh_fem *mf2 =  md.pmesh_fem_of_variable(w);
 
 2689         GMM_ASSERT1(!mf2 || &(mf2->linked_mesh()) == &(mf->linked_mesh()),
 
 2690                     "The data for the sliding velocity should be defined on " 
 2691                     " the same mesh as the displacement variable");
 
 2694       for (
size_type i = 0; i < boundaries.size(); ++i) {
 
 2695         const contact_boundary &cb = boundaries[i];
 
 2696         if (&(md.mesh_fem_of_variable(cb.varname_u).linked_mesh())
 
 2697             == &(mf->linked_mesh()) && cb.region == region)
 
 2698           GMM_ASSERT1(
false, 
"This contact boundary has already been added");
 
 2702           (md, transformation_name, mf->linked_mesh(), u_group, region);
 
 2705           (md, transformation_name, mf->linked_mesh(), u_group, region);
 
 2707       boundaries.push_back(contact_boundary());
 
 2708       contact_boundary &cb = boundaries.back();
 
 2711       if (is_slave) cb.sigma_u = sigma_u;
 
 2713       cb.is_master = is_master;
 
 2714       cb.is_slave = is_slave;
 
 2715       cb.is_unbiased= is_unbiased;
 
 2718         std::string n, n0, Vs, g, Y, gamma;
 
 2720         gamma =
"("+Nitsche_param+
"/element_size)";
 
 2721         n = 
"Transformed_unit_vector(Grad_"+u+
", Normal)";
 
 2722         n0 = 
"Transformed_unit_vector(Grad_"+w+
", Normal)";
 
 2730         Y = 
"Interpolate(X,"+transformation_name+
")";
 
 2731         g = 
"("+Y+
"+Interpolate("+u_group+
","+transformation_name+
")-X-"+u+
")."+n;
 
 2732         Vs = 
"("+u+
"-Interpolate("+u_group+
","+transformation_name+
")";
 
 2734           Vs += 
"-"+w+
"+Interpolate("+w_group+
","+transformation_name+
")";
 
 2735           if (frame_indifferent)
 
 2736             Vs += 
"+("+g+
")*("+n+
"-"+n0+
")";
 
 2740         std::string coupled_projection_def =
 
 2741           "Coulomb_friction_coupled_projection(" 
 2742           + sigma_u +
","+n+
","+Vs+
","+g+
","+friction_coeff+
"," 
 2750         g = 
"(Interpolate(X,"+transformation_name+
")-X-"+u+
")."+n;
 
 2751         if (frame_indifferent && w.size())
 
 2752           Vs = 
"("+u+
"-"+w+
"+"+g+
"*("+n+
"-"+n0+
"))*"+
alpha;
 
 2754           Vs = 
"("+u+(w.size() ? (
"-"+w):
"")+
")*"+
alpha;
 
 2756         std::string coupled_projection_rig =
 
 2757           "Coulomb_friction_coupled_projection(" 
 2758           + sigma_u +
","+n+
","+Vs+
","+g+
","+ friction_coeff+
"," 
 2763           (is_unbiased ? 
"-0.5*" : 
"-")
 
 2765           + (
"Interpolate_filter("+transformation_name+
"," 
 2766                             +coupled_projection_def+
"."+test_u+
",1) ") 
 
 2767           +(is_unbiased ? 
"":
"-Interpolate_filter("+transformation_name+
"," 
 2768                             +coupled_projection_rig+
"."+test_u+
",2) ")
 
 2774           + (is_unbiased ? 
"+ 0.5*" : 
"+ ")
 
 2775           +
"Interpolate_filter("+transformation_name+
"," 
 2776                             +coupled_projection_def +
".Interpolate("+test_u_group+
"," + transformation_name+
"), 1)" 
 2777           +(is_unbiased ? 
"":
"+ Interpolate_filter("+transformation_name+
"," 
 2778                             +coupled_projection_rig +
".Interpolate("+test_u_group+
"," + transformation_name+
"), 2)");                  
 
 2782     virtual void asm_real_tangent_terms(
const model &md, 
size_type ,
 
 2783                                         const model::varnamelist &,
 
 2784                                         const model::varnamelist &,
 
 2785                                         const model::mimlist &,
 
 2786                                         model::real_matlist &,
 
 2787                                         model::real_veclist &,
 
 2788                                         model::real_veclist &,
 
 2790                                         build_version)
 const {
 
 2795       for (
const contact_boundary &cb : boundaries) {
 
 2797           md.add_generic_expression(cb.expr, *(cb.mim), cb.region);
 
 2802     Nitsche_large_sliding_contact_brick_raytracing
 
 2804       const std::string &Nitsche_parameter,
 
 2805      const std::string &f_coeff,
const std::string &ug,
 
 2806      const std::string &wg, 
const std::string &tr,
 
 2807      const std::string &alpha_ = 
"1", 
bool sym_v = 
false,
 
 2808      bool frame_indiff = 
false) {
 
 2809       transformation_name = tr;
 
 2810       u_group = ug; w_group = wg;
 
 2811       friction_coeff = f_coeff;
 
 2813       Nitsche_param = Nitsche_parameter;
 
 2814       sym_version = sym_v;
 
 2815       frame_indifferent = frame_indiff;
 
 2818       set_flags(
"Integral large sliding contact bick raytracing",
 
 2829     pbrick pbr = md.brick_pointer(indbrick);
 
 2830     Nitsche_large_sliding_contact_brick_raytracing *p
 
 2831       = 
dynamic_cast<Nitsche_large_sliding_contact_brick_raytracing *
> 
 2833     GMM_ASSERT1(p, 
"Wrong type of brick");
 
 2834     return p->transformation_name;
 
 2839     pbrick pbr = md.brick_pointer(indbrick);
 
 2840     Nitsche_large_sliding_contact_brick_raytracing *p
 
 2841       = 
dynamic_cast<Nitsche_large_sliding_contact_brick_raytracing *
> 
 2843     GMM_ASSERT1(p, 
"Wrong type of brick");
 
 2849     pbrick pbr = md.brick_pointer(indbrick);
 
 2850     Nitsche_large_sliding_contact_brick_raytracing *p
 
 2851       = 
dynamic_cast<Nitsche_large_sliding_contact_brick_raytracing *
> 
 2853     GMM_ASSERT1(p, 
"Wrong type of brick");
 
 2859     pbrick pbr = md.brick_pointer(indbrick);
 
 2860     Nitsche_large_sliding_contact_brick_raytracing *p
 
 2861       = 
dynamic_cast<Nitsche_large_sliding_contact_brick_raytracing *
> 
 2863     GMM_ASSERT1(p, 
"Wrong type of brick");
 
 2865       (md, p->transformation_name, expr, N);
 
 2870    bool is_master, 
bool is_slave, 
bool is_unbiased, 
const std::string &u,
 
 2871    const std::string &sigma_u, 
const std::string &w) {
 
 2873     pbrick pbr = md.brick_pointer(indbrick);
 
 2874     Nitsche_large_sliding_contact_brick_raytracing *p
 
 2875       = 
dynamic_cast<Nitsche_large_sliding_contact_brick_raytracing *
> 
 2877     GMM_ASSERT1(p, 
"Wrong type of brick");
 
 2879     bool found_u = 
false, found_sigma = 
false;
 
 2880     for (
const auto & v : p->vl) {
 
 2881       if (v.compare(u) == 0) found_u = 
true;
 
 2882       if (v.compare(sigma_u) == 0) found_sigma = 
true;
 
 2884     if (!found_u) p->vl.push_back(u);
 
 2885     GMM_ASSERT1(!is_slave || sigma_u.size(),
 
 2886                 "You should define an expression of sigma(u) on each slave boundary");
 
 2891     std::vector<std::string> ug = md.variable_group(p->u_group);
 
 2893     for (
const auto &uu : ug)
 
 2894       if (uu.compare(u) == 0) { found_u = 
true; 
break; }
 
 2897       md.define_variable_group(p->u_group, ug);
 
 2901       bool found_w = 
false;
 
 2902       for (
const auto &ww : p->dl)
 
 2903         if (ww.compare(w) == 0) { found_w = 
true; 
break; }
 
 2908       std::vector<std::string> wg = md.variable_group(p->w_group);
 
 2910       for (
const auto &ww : wg)
 
 2911         if (ww.compare(w) == 0) { found_w = 
true; 
break; }
 
 2914         md.define_variable_group(p->w_group, wg);
 
 2918     bool found_mim = 
false;
 
 2919     for (
const auto &pmim : p->ml)
 
 2920       if (pmim == &mim) { found_mim = 
true; 
break; }
 
 2922       p->ml.push_back(&mim);
 
 2926     p->add_contact_boundary(md, mim, region, is_master, is_slave,is_unbiased,
 
 2931   (
model &md, 
bool is_unbiased, 
const std::string &Nitsche_param,
 
 2932    scalar_type release_distance, 
const std::string &f_coeff, 
 
 2933    const std::string &alpha, 
bool sym_v, 
bool frame_indifferent) {
 
 2935     char ugroupname[50], wgroupname[50], transname[50];
 
 2936     for (
int i = 0; i < 10000; ++i) {
 
 2937       snprintf(ugroupname, 49, 
"disp__group_raytracing_%d", i);
 
 2938       if (!(md.variable_group_exists(ugroupname))
 
 2942     md.define_variable_group(ugroupname, std::vector<std::string>());
 
 2944     for (
int i = 0; i < 10000; ++i) {
 
 2945       snprintf(wgroupname, 49, 
"w__group_raytracing_%d", i);
 
 2946       if (!(md.variable_group_exists(wgroupname))
 
 2950     md.define_variable_group(wgroupname, std::vector<std::string>());
 
 2952     for (
int i = 0; i < 10000; ++i) {
 
 2953       snprintf(transname, 49, 
"trans__raytracing_%d", i);
 
 2959     model::varnamelist vl, dl;
 
 2964     auto p = std::make_shared<Nitsche_large_sliding_contact_brick_raytracing>
 
 2965       (is_unbiased, Nitsche_param , f_coeff, ugroupname, wgroupname, transname, alpha,
 
 2966        sym_v, frame_indifferent);
 
 2970     return md.
add_brick(pbr, p->vl, p->dl, model::termlist(), model::mimlist(),
 
region-tree for window/point search on a set of rectangles.
 
const base_node & xreal() const
coordinates of the current point, in the real convex.
 
void set_ii(size_type ii__)
change the current point (assuming a geotrans_precomp_ is used)
 
Balanced tree of n-dimensional rectangles.
 
structure passed as the argument of fem interpolation functions.
 
Describe a finite element method linked to a mesh.
 
virtual size_type nb_basic_dof() const
Return the total number of basic degrees of freedom (before the optional reduction).
 
Describe an integration method linked to a mesh.
 
"iterator" class for regions.
 
`‘Model’' variables store the variables, the data and the description of a model.
 
bool interpolate_transformation_exists(const std::string &name) const
Tests if name corresponds to an interpolate transformation.
 
void change_mims_of_brick(size_type ib, const mimlist &ml)
Change the mim list of a brick.
 
size_type add_brick(pbrick pbr, const varnamelist &varnames, const varnamelist &datanames, const termlist &terms, const mimlist &mims, size_type region)
Add a brick to the model.
 
void change_data_of_brick(size_type ib, const varnamelist &vl)
Change the data list of a brick.
 
void change_variables_of_brick(size_type ib, const varnamelist &vl)
Change the variable list of a brick.
 
bool variable_exists(const std::string &name) const
States if a name corresponds to a declared variable.
 
The virtual brick has to be derived to describe real model bricks.
 
Miscelleanous assembly routines for common terms. Use the low-level generic assembly....
 
void mult_add(const L1 &l1, const L2 &l2, L3 &l3)
*/
 
void copy(const L1 &l1, L2 &l2)
*/
 
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2(const V &v)
Euclidean norm of a vector.
 
void clear(L &l)
clear (fill with zeros) a vector or matrix.
 
number_traits< typename linalg_traits< V1 >::value_type >::magnitude_type vect_dist2(const V1 &v1, const V2 &v2)
Euclidean distance between two vectors.
 
void resize(V &v, size_type n)
*/
 
void mult(const L1 &l1, const L2 &l2, L3 &l3)
*/
 
strongest_value_type< V1, V2 >::value_type vect_sp(const V1 &v1, const V2 &v2)
*/
 
void add(const L1 &l1, L2 &l2)
*/
 
computation of the condition number of dense matrices.
 
void lu_inverse(const DenseMatrixLU &LU, const Pvector &pvector, const DenseMatrix &AInv_)
Given an LU factored matrix, build the inverse of the matrix.
 
size_type lu_factor(DenseMatrix &A, Pvector &ipvt)
LU Factorization of a general (dense) matrix (real or complex).
 
void lu_solve(const DenseMatrix &LU, const Pvector &pvector, VectorX &x, const VectorB &b)
LU Solve : Solve equation Ax=b, given an LU factored matrix.
 
void asm_mass_matrix(const MAT &M, const mesh_im &mim, const mesh_fem &mf1, const mesh_region &rg=mesh_region::all_convexes())
generic mass matrix assembly (on the whole mesh or on the specified convex set or boundary)
 
void asm_source_term(const VECT1 &B, const mesh_im &mim, const mesh_fem &mf, const mesh_fem &mf_data, const VECT2 &F, const mesh_region &rg=mesh_region::all_convexes())
source term (for both volumic sources and boundary (Neumann) sources).
 
size_type convex_num() const
get the current convex number
 
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
 
void grad_base_value(base_tensor &t, bool withM=true) const
fill the tensor with the gradient of the base functions (taken at point this->xref())
 
void base_value(base_tensor &t, bool withM=true) const
fill the tensor with the values of the base functions (taken at point this->xref())
 
short_type face_num() const
get the current face number
 
const pfem pf() const
get the current FEM descriptor
 
gmm::uint16_type short_type
used as the common short type integer in the library
 
std::shared_ptr< const convex_structure > pconvex_structure
Pointer on a convex structure description.
 
base_small_vector compute_normal(const geotrans_interpolation_context &c, size_type face)
norm of returned vector is the ratio between the face surface on the real element and the face surfac...
 
size_t size_type
used as the common size type in the library
 
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
 
size_type alpha(short_type n, short_type d)
Return the value of  which is the number of monomials of a polynomial of  variables and degree .
 
GEneric Tool for Finite Element Methods.
 
void add_rigid_obstacle_to_raytracing_transformation(model &md, const std::string &transname, const std::string &expr, size_type N)
Add a rigid obstacle whose geometry corresponds to the zero level-set of the high-level generic assem...
 
const std::string & transformation_name_of_large_sliding_contact_brick(model &md, size_type indbrick)
Gives the name of the raytracing interpolate transformation for an existing large sliding contact bri...
 
const std::string & sliding_data_group_name_of_Nitsche_large_sliding_contact_brick(model &md, size_type indbrick)
Gives the name of the group of variables corresponding to the sliding data for an existing large slid...
 
size_type add_integral_large_sliding_contact_brick_raytracing(model &md, const std::string &augm_param, scalar_type release_distance, const std::string &f_coeff="0", const std::string &alpha="1", bool sym_v=false, bool frame_indifferent=false)
Adds a large sliding contact with friction brick to the model.
 
size_type add_integral_large_sliding_contact_brick_raytrace(model &md, multi_contact_frame &mcf, const std::string &dataname_r, const std::string &dataname_friction_coeff=std::string(), const std::string &dataname_alpha=std::string())
Adds a large sliding contact with friction brick to the model.
 
void add_raytracing_transformation(model &md, const std::string &transname, scalar_type release_distance)
Add a raytracing interpolate transformation called 'transname' to a model to be used by the generic a...
 
size_type add_Nitsche_large_sliding_contact_brick_raytracing(model &md, bool is_unbiased, const std::string &Nitsche_param, scalar_type release_distance, const std::string &f_coeff="0", const std::string &alpha="1", bool sym_v=false, bool frame_indifferent=false)
Adds a large sliding contact with friction brick to the model based on Nitsche's method.
 
void add_rigid_obstacle_to_Nitsche_large_sliding_contact_brick(model &md, size_type indbrick, const std::string &expr, size_type N)
Adds a rigid obstacle to an existing large sliding contact with friction brick.
 
const std::string & sliding_data_group_name_of_large_sliding_contact_brick(model &md, size_type indbrick)
Gives the name of the group of variables corresponding to the sliding data for an existing large slid...
 
void add_contact_boundary_to_Nitsche_large_sliding_contact_brick(model &md, size_type indbrick, const mesh_im &mim, size_type region, bool is_master, bool is_slave, bool is_unbiased, const std::string &u, const std::string &lambda="", const std::string &w="")
Adds a contact boundary to an existing large sliding contact with friction brick.
 
void add_rigid_obstacle_to_large_sliding_contact_brick(model &md, size_type indbrick, const std::string &expr, size_type N)
Adds a rigid obstacle to an existing large sliding contact with friction brick.
 
void add_master_contact_boundary_to_raytracing_transformation(model &md, const std::string &transname, const mesh &m, const std::string &dispname, size_type region)
Add a master boundary with corresponding displacement variable 'dispname' on a specific boundary 'reg...
 
void slice_vector_on_basic_dof_of_element(const mesh_fem &mf, const VEC1 &vec, size_type cv, VEC2 &coeff, size_type qmult1=size_type(-1), size_type qmult2=size_type(-1))
Given a mesh_fem.
 
std::shared_ptr< const virtual_brick > pbrick
type of pointer on a brick
 
const std::string & transformation_name_of_Nitsche_large_sliding_contact_brick(model &md, size_type indbrick)
Gives the name of the raytracing interpolate transformation for an existing large sliding contact bri...
 
const std::string & displacement_group_name_of_large_sliding_contact_brick(model &md, size_type indbrick)
Gives the name of the group of variables corresponding to the displacement for an existing large slid...
 
const std::string & displacement_group_name_of_Nitsche_large_sliding_contact_brick(model &md, size_type indbrick)
Gives the name of the group of variables corresponding to the displacement for an existing large slid...
 
void add_slave_contact_boundary_to_raytracing_transformation(model &md, const std::string &transname, const mesh &m, const std::string &dispname, size_type region)
Add a slave boundary with corresponding displacement variable 'dispname' on a specific boundary 'regi...
 
void add_contact_boundary_to_large_sliding_contact_brick(model &md, size_type indbrick, const mesh_im &mim, size_type region, bool is_master, bool is_slave, const std::string &u, const std::string &lambda="", const std::string &w="")
Adds a contact boundary to an existing large sliding contact with friction brick.