27 #ifndef GETFEM_HAVE_QHULL_QHULL_H
38 static void orthonormal_basis_to_unit_vec(
size_type d,
const base_node &un,
41 for (
size_type k = 0; k <= d && n < d; ++k) {
44 ut[n][k] = scalar_type(1);
48 ut[n] -= gmm::vect_sp(ut[nn], ut[n]) * ut[nn];
50 if (gmm::vect_norm2(ut[n]) < 1e-3)
continue;
54 GMM_ASSERT1(n == d,
"Gram-Schmidt algorithm to find an "
55 "orthonormal basis for the tangential displacement failed");
65 std::vector<size_type> cvs;
66 std::vector<short_type> fcs;
68 contact_node() : mf(0), cvs(0), fcs(0) {}
69 contact_node(
const mesh_fem &mf_) {mf = &mf_;}
73 struct contact_node_pair {
74 contact_node cn_s, cn_m;
77 contact_node_pair(scalar_type threshold=10.) : cn_s(), cn_m()
78 {dist2 = threshold * threshold; is_active =
false;}
82 class contact_node_pair_list :
public std::vector<contact_node_pair> {
84 void contact_node_list_from_region
85 (
const mesh_fem &mf,
size_type contact_region,
86 std::vector<contact_node> &cnl) {
89 const mesh &m = mf.linked_mesh();
91 std::map<size_type, size_type> dof_to_cnid;
93 dal::bit_vector dofs = mf.basic_dof_on_region(contact_region);
94 for (dal::bv_visitor dof(dofs); !dof.finished(); ++dof)
95 if ( dof % qdim == 0) {
96 dof_to_cnid[dof] = cnid++;
97 contact_node new_cn(mf);
99 cnl.push_back(new_cn);
101 for (mr_visitor face(m.region(contact_region));
102 !face.finished(); ++face) {
103 assert(face.is_face());
104 mesh_fem::ind_dof_face_ct
105 face_dofs = mf.ind_basic_dof_of_face_of_element(face.cv(),face.f());
106 for (
size_type it=0; it < face_dofs.size(); it += qdim ) {
108 cnid = dof_to_cnid[dof];
109 cnl[cnid].cvs.push_back(face.cv());
110 cnl[cnid].fcs.push_back(face.f());
116 contact_node_pair_list() : std::vector<contact_node_pair>() {}
118 void append_min_dist_cn_pairs(
const mesh_fem &mf1,
const mesh_fem &mf2,
120 bool slave1=
true,
bool slave2=
false) {
122 std::vector<contact_node> cnl1(0), cnl2(0);
123 contact_node_list_from_region(mf1, rg1, cnl1);
124 contact_node_list_from_region(mf2, rg2, cnl2);
128 size_type size1 = slave1 ? cnl1.size() : 0;
129 size_type size2 = slave2 ? cnl2.size() : 0;
130 this->
resize( size0 + size1 + size2 );
131 # ifndef GETFEM_HAVE_QHULL_QHULL_H
133 for (
size_type i1 = 0; i1 < cnl1.size(); ++i1) {
134 contact_node *cn1 = &cnl1[i1];
137 for (
size_type i2 = 0; i2 < cnl2.size(); ++i2) {
138 contact_node *cn2 = &cnl2[i2];
143 for (
size_type i1 = 0; i1 < cnl1.size(); ++i1, ++ii1) {
144 contact_node *cn1 = &cnl1[i1];
145 base_node node1 = cn1->mf->point_of_basic_dof(cn1->dof);
147 scalar_type dist2 = tree2.nearest_neighbor(ipt, node1);
148 if (ipt.i !=
size_type(-1) && dist2 < (*
this)[ii1].dist2) {
149 (*this)[ii1].cn_s = *cn1;
150 (*this)[ii1].cn_m = cnl2[ipt.i];
151 (*this)[ii1].dist2 = dist2;
152 (*this)[ii1].is_active =
true;
158 for (
size_type i2 = 0; i2 < cnl2.size(); ++i2, ++ii2) {
159 contact_node *cn2 = &cnl2[i2];
160 base_node node2 = cn2->mf->point_of_basic_dof(cn2->dof);
162 scalar_type dist2 = tree1.nearest_neighbor(ipt, node2);
163 if (ipt.i !=
size_type(-1) && dist2 < (*
this)[ii2].dist2) {
164 (*this)[ii2].cn_s = *cn2;
165 (*this)[ii2].cn_m = cnl1[ipt.i];
166 (*this)[ii2].dist2 = dist2;
167 (*this)[ii2].is_active =
true;
172 std::vector<base_node> pts;
173 for (
size_type i1 = 0; i1 < cnl1.size(); ++i1) {
174 contact_node *cn1 = &cnl1[i1];
175 pts.push_back(cn1->mf->point_of_basic_dof(cn1->dof));
177 for (
size_type i2 = 0; i2 < cnl2.size(); ++i2) {
178 contact_node *cn2 = &cnl2[i2];
179 pts.push_back(cn2->mf->point_of_basic_dof(cn2->dof));
181 gmm::dense_matrix<size_type> simplexes;
183 getfem::delaunay(pts, simplexes);
185 size_type nb_vertices = gmm::mat_nrows(simplexes);
186 std::vector<size_type> facet_vertices(nb_vertices);
187 std::vector< std::vector<size_type> > pt1_neighbors(size1);
188 for (
size_type i = 0; i < gmm::mat_ncols(simplexes); ++i) {
189 gmm::copy(gmm::mat_col(simplexes, i), facet_vertices);
190 for (
size_type iv1 = 0; iv1 < nb_vertices-1; ++iv1) {
192 bool v1_on_surface1 = (v1 < size1);
193 for (
size_type iv2 = iv1 + 1; iv2 < nb_vertices; ++iv2) {
195 bool v2_on_surface1 = (v2 < size1);
196 if (v1_on_surface1 ^ v2_on_surface1) {
197 bool already_in =
false;
198 size_type vv1 = (v1_on_surface1 ? v1 : v2);
199 size_type vv2 = (v2_on_surface1 ? v1 : v2);
200 for (
size_type j = 0; j < pt1_neighbors[vv1].size(); ++j)
201 if (pt1_neighbors[vv1][j] == vv2) {
205 if (!already_in) pt1_neighbors[vv1].push_back(vv2);
212 for (
size_type j = 0; j < pt1_neighbors[i1].size(); ++j) {
213 size_type i2 = pt1_neighbors[i1][j] - size1;
216 contact_node *cn1 = &cnl1[i1];
217 base_node node1 = cn1->mf->point_of_basic_dof(cn1->dof);
218 contact_node *cn2 = &cnl2[i2];
219 base_node node2 = cn2->mf->point_of_basic_dof(cn2->dof);
221 if (slave1 && dist2 < (*
this)[ii1].dist2) {
222 (*this)[ii1].cn_s = *cn1;
223 (*this)[ii1].cn_m = *cn2;
224 (*this)[ii1].dist2 = dist2;
225 (*this)[ii1].is_active =
true;
227 if (slave2 && dist2 < (*
this)[ii2].dist2) {
228 (*this)[ii2].cn_s = *cn2;
229 (*this)[ii2].cn_m = *cn1;
230 (*this)[ii2].dist2 = dist2;
231 (*this)[ii2].is_active =
true;
237 void append_min_dist_cn_pairs(
const mesh_fem &mf,
239 bool slave1=
true,
bool slave2=
false) {
240 append_min_dist_cn_pairs(mf, mf, rg1, rg2, slave1, slave2);
244 scalar_type projection_on_convex_face
246 const base_node &master_node,
const base_node &slave_node,
247 base_node &un, base_node &proj_node, base_node &proj_node_ref) {
251 if (pgt->is_linear()) {
253 un = m.normal_of_face_of_convex(cv,fc);
256 gmm::add(master_node, gmm::scaled(slave_node, -1.), proj_node);
257 gmm::copy(gmm::scaled(un, gmm::vect_sp(proj_node, un)), proj_node);
261 gic.init(m.points_of_convex(cv), pgt);
262 gic.
invert(proj_node, proj_node_ref);
263 return pgt->convex_ref()->is_in(proj_node_ref);
271 size_type nb_pts_fc = pgt->structure()->nb_points_of_face(fc);
273 bgeot::convex_ind_ct ind_pts_fc = pgt->structure()->ind_points_of_face(fc);
274 ref_mesh_face_pt_ct pts_fc = m.points_of_face_of_convex(cv, fc);
277 base_matrix base_ref_fc(P-1,N);
280 GMM_ASSERT1( dref_pts_fc.size() == P,
"Dimensions mismatch");
281 base_node vec(dref_pts_fc[0].size());
283 vec = dref_pts_fc[i+1] - dref_pts_fc[0];
284 gmm::copy(vec,gmm::mat_row(base_ref_fc,i));
288 GMM_ASSERT1( slave_node.size() == N,
"Dimensions mismatch");
289 const base_node &xx = slave_node;
290 base_node &xxp = proj_node; xxp.resize(N);
291 base_node &xp = proj_node_ref; xp.resize(P);
296 xp = gmm::mean_value(pgt->convex_ref()->points_of_face(fc));
299 base_vector val(nb_pts_fc);
300 pgt->poly_vector_val(xp, ind_pts_fc, val);
301 for (
size_type l = 0; l < nb_pts_fc; ++l)
302 gmm::add(gmm::scaled(pts_fc[l], val[l] ), xxp);
304 base_matrix G(N, nb_pts_fc);
305 vectors_to_base_matrix(G, pts_fc);
307 base_matrix K(N,P-1);
309 base_matrix grad_fc(nb_pts_fc, P);
310 base_matrix grad_fc1(nb_pts_fc, P-1);
311 base_matrix B(N,P-1), BB(N,P), CS(P-1,P-1);
313 scalar_type EPS = 10E-12;
315 while (res > EPS && --cnt) {
317 pgt->poly_vector_grad(xp, ind_pts_fc, grad_fc);
318 gmm::mult(grad_fc, gmm::transposed(base_ref_fc), grad_fc1);
321 bgeot::lu_inverse(&(*(CS.begin())), P-1);
328 pgt->poly_vector_val(xp, ind_pts_fc, val);
329 for (
size_type l = 0; l < nb_pts_fc; ++l)
330 gmm::add(gmm::scaled(pts_fc[l], val[l]), xxp);
332 gmm::mult(gmm::transposed(BB), xx - xxp, vres);
335 GMM_ASSERT1( res <= EPS,
336 "Iterative pojection on convex face did not converge");
338 pgt->poly_vector_grad(xp, ind_pts_fc, grad_fc);
339 gmm::mult(grad_fc, gmm::transposed(base_ref_fc), grad_fc1);
351 base_matrix grad_cv(nb_pts_cv, P);
352 pgt->poly_vector_grad(xp, grad_cv);
354 base_matrix GG(N, nb_pts_cv);
355 m.points_of_convex(cv, GG);
360 base_matrix bases_product(P-1, P);
361 gmm::mult(gmm::transposed(K), KK, bases_product);
364 std::vector<size_type> ind(0);
366 if (j != i ) ind.push_back(j);
367 gmm::sub_index SUBI(ind);
370 gmm::sub_interval(0, P-1), SUBI));
371 gmm::add(gmm::scaled(gmm::mat_col(KK, i), (i % 2) ? -det : +det ),
376 gmm::scale(un, 1/gmm::vect_norm2(un));
378 if (gmm::vect_sp(un, gmm::mean_value(pts_fc) -
379 gmm::mean_value(m.points_of_convex(cv))) < 0)
380 gmm::scale(un,scalar_type(-1));
382 return pgt->convex_ref()->is_in(proj_node_ref);
386 void compute_contact_matrices
387 (
const mesh_fem &mf_disp1,
const mesh_fem &mf_disp2,
388 contact_node_pair_list &cnpl, model_real_plain_vector &gap,
389 CONTACT_B_MATRIX *BN1, CONTACT_B_MATRIX *BN2 = 0,
390 CONTACT_B_MATRIX *BT1 = 0, CONTACT_B_MATRIX *BT2 = 0) {
392 GMM_ASSERT1(gmm::vect_size(gap) == cnpl.size(),
393 "Wrong number of contact node pairs or wrong size of gap");
395 GMM_ASSERT1( gmm::mat_nrows(*BN1) == cnpl.size(),
"Wrong size of BN1");
398 GMM_ASSERT1( gmm::mat_nrows(*BN2) == cnpl.size(),
"Wrong size of BN2");
400 dim_type qdim = mf_disp1.get_qdim();
404 GMM_ASSERT1( gmm::mat_nrows(*BT1) == cnpl.size() * d,
"Wrong size of BT1");
408 GMM_ASSERT1( gmm::mat_nrows(*BT2) == cnpl.size() * d,
"Wrong size of BT2");
411 for (
size_type row = 0; row < cnpl.size(); ++row) {
412 contact_node_pair *cnp = &cnpl[row];
413 if (cnp->is_active) {
414 contact_node *cn_s = &cnp->cn_s;
415 contact_node *cn_m = &cnp->cn_m;
416 const mesh &mesh_m = cn_m->mf->linked_mesh();
417 base_node slave_node = cn_s->mf->point_of_basic_dof(cn_s->dof);
418 base_node master_node = cn_m->mf->point_of_basic_dof(cn_m->dof);
419 GMM_ASSERT1(slave_node.size() == qdim && master_node.size() == qdim,
421 base_node un_sel(qdim), proj_node_sel(qdim), proj_node_ref_sel(qdim);
422 scalar_type is_in_min = 1e5;
425 std::vector<size_type>::iterator cv;
426 std::vector<short_type>::iterator fc;
427 for (cv = cn_m->cvs.begin(), fc = cn_m->fcs.begin();
428 cv != cn_m->cvs.end() && fc != cn_m->fcs.end(); cv++, fc++) {
429 base_node un(qdim), proj_node(qdim), proj_node_ref(qdim);
430 scalar_type is_in = projection_on_convex_face
431 (mesh_m, *cv, *fc, master_node, slave_node, un, proj_node, proj_node_ref);
432 if (is_in < is_in_min) {
437 proj_node_sel = proj_node;
438 proj_node_ref_sel = proj_node_ref;
441 if (is_in_min < 0.05) {
442 gap[row] =
gmm::vect_sp(slave_node-proj_node_sel, un_sel);
444 std::vector<base_node> ut(d);
445 if (BT1) orthonormal_basis_to_unit_vec(d, un_sel, &(ut[0]));
447 CONTACT_B_MATRIX *BN = 0;
448 CONTACT_B_MATRIX *BT = 0;
449 if (cn_s->mf == &mf_disp1) {
452 }
else if (cn_s->mf == &mf_disp2) {
458 (*BN)(row, cn_s->dof + k) -= un_sel[k];
462 (*BT)(row * d + n, cn_s->dof + k) -= ut[n][k];
465 const mesh_fem *mf_disp = 0;
466 if (cn_m->mf == &mf_disp1) {
470 }
else if (cn_m->mf == &mf_disp2) {
477 base_matrix M(qdim, mf_disp->nb_basic_dof_of_element(cv_sel));
478 mesh_m.points_of_convex(cv_sel, G);
479 pfem pf = mf_disp->fem_of_element(cv_sel);
481 fem_interpolation_context
482 ctx(pgt, pf, proj_node_ref_sel, G, cv_sel, fc_sel);
483 pf->interpolation (ctx, M,
int(qdim));
486 master_dofs = mf_disp->ind_basic_dof_of_element(cv_sel);
488 model_real_plain_vector MT_u(mf_disp->nb_basic_dof_of_element(cv_sel));
489 gmm::mult(gmm::transposed(M), un_sel, MT_u);
490 for (
size_type j = 0; j < master_dofs.size(); ++j)
491 (*BN)(row, master_dofs[j]) += MT_u[j];
495 gmm::mult(gmm::transposed(M), ut[n], MT_u);
496 for (
size_type j = 0; j < master_dofs.size(); ++j)
497 (*BT)(row * d + n, master_dofs[j]) += MT_u[j];
516 struct Coulomb_friction_brick :
public virtual_brick {
518 mutable CONTACT_B_MATRIX BN1, BT1, BN2, BT2;
519 mutable CONTACT_B_MATRIX DN, DDN, DT, DDT;
520 mutable CONTACT_B_MATRIX BBN1, BBT1, BBN2, BBT2;
521 mutable model_real_plain_vector gap, threshold, friction_coeff,
alpha;
522 mutable model_real_plain_vector RLN, RLT;
523 mutable scalar_type r, gamma;
524 mutable bool is_init;
525 bool Tresca_version, contact_only;
526 bool really_stationary, friction_dynamic_term;
527 bool two_variables, Hughes_stabilized;
528 int augmentation_version;
533 void init_BBN_BBT(
void)
const {
534 gmm::resize(BBN1, gmm::mat_nrows(BN1), gmm::mat_ncols(BN1));
536 if (Hughes_stabilized) {
537 gmm::resize(DDN, gmm::mat_nrows(DN), gmm::mat_ncols(DN));
541 gmm::resize(BBN2, gmm::mat_nrows(BN2), gmm::mat_ncols(BN2));
545 if (Hughes_stabilized) {
546 gmm::resize(DDT, gmm::mat_nrows(DT), gmm::mat_ncols(DT));
549 gmm::resize(BBT1, gmm::mat_nrows(BT1), gmm::mat_ncols(BT1));
552 gmm::resize(BBT2, gmm::mat_nrows(BT2), gmm::mat_ncols(BT2));
559 gmm::scale(gmm::mat_row(BBN1, i), alpha[i]);
560 if (Hughes_stabilized) gmm::scale(gmm::mat_row(DDN, i), alpha[i]);
562 gmm::scale(gmm::mat_row(BBN2, i), alpha[i]);
565 if (Hughes_stabilized)
566 gmm::scale(gmm::mat_row(DDT, d*i+k), alpha[i]);
567 gmm::scale(gmm::mat_row(BBT1, d*i+k), alpha[i]);
569 gmm::scale(gmm::mat_row(BBT2, d*i+k), alpha[i]);
575 void set_stationary(
void) { really_stationary =
true; }
577 void precomp(
const model_real_plain_vector &u1,
578 const model_real_plain_vector &u2,
579 const model_real_plain_vector &lambda_n,
580 const model_real_plain_vector &lambda_t,
581 const model_real_plain_vector &wt1,
582 const model_real_plain_vector &wt2)
const {
584 for (
size_type i = 0; i < gmm::mat_nrows(BN1); ++i) RLN[i] *= alpha[i];
587 if (Hughes_stabilized)
589 if (two_variables)
gmm::mult_add(BBN2, gmm::scaled(u2, -r), RLN);
592 if (friction_dynamic_term) {
597 if (!really_stationary) {
599 if (two_variables)
gmm::mult_add(BBT2, gmm::scaled(u2, -r), RLT);
601 if (Hughes_stabilized)
607 void basic_asm_real_tangent_terms(
const model_real_plain_vector &u1,
608 const model_real_plain_vector &u2,
609 const model_real_plain_vector &lambda_n,
610 const model_real_plain_vector &lambda_t,
611 const model_real_plain_vector &wt1,
612 const model_real_plain_vector &wt2,
613 model::real_matlist &matl,
614 model::real_veclist &vecl,
615 build_version version)
const {
616 size_type nbt = 4 + (contact_only ? 0 : 4) + (two_variables ? 3 : 0)
617 + (two_variables && !contact_only ? 2 : 0);
618 GMM_ASSERT1(matl.size() == nbt,
619 "Wrong number of terms for the contact brick");
621 const scalar_type vt1 = scalar_type(1), vt0 = scalar_type(0);
627 model_real_sparse_matrix &T_u1_u1 = matl[nt++], &T_u2_u2 = matl[nt++];
628 if (!two_variables) nt--;
629 model_real_sparse_matrix &T_u1_n = matl[nt++], &T_n_u1 = matl[nt++];
630 if (!two_variables) nt -= 2;
631 model_real_sparse_matrix &T_u2_n = matl[nt++], &T_n_u2 = matl[nt++];
633 model_real_sparse_matrix &T_n_n = matl[nt++];
634 if (contact_only) nt -= 2;
635 model_real_sparse_matrix &T_u1_t = matl[nt++], &T_t_u1 = matl[nt++];
636 if (contact_only || !two_variables) nt -= 2;
637 model_real_sparse_matrix &T_u2_t = matl[nt++], &T_t_u2 = matl[nt++];
638 if (contact_only) nt -= 2;
640 model_real_sparse_matrix &T_t_t = matl[nt++], &T_t_n = matl[nt++];
643 model_real_plain_vector &ru1 = vecl[0];
644 model_real_plain_vector &ru2 = vecl[1];
645 model_real_plain_vector &rlambda_n = vecl[nvec_lambda_n];
646 model_real_plain_vector &rlambda_t = vecl[nvec_lambda_t];
649 if (!is_init) init_BBN_BBT();
652 if (augmentation_version <= 2)
653 precomp(u1, u2, lambda_n, lambda_t, wt1, wt2);
655 if (version & model::BUILD_MATRIX) {
656 base_matrix pg(d, d);
668 switch (augmentation_version) {
670 gmm::copy(gmm::scaled(gmm::transposed(BN1), -vt1), T_u1_n);
672 gmm::copy(gmm::scaled(gmm::transposed(BN2), -vt1), T_u2_n);
676 if (two_variables)
gmm::clear(gmm::mat_col(T_u2_n, i));
677 T_n_n(i, i) = -vt1/(r*
alpha[i]);
679 if (Hughes_stabilized && RLN[i] <= vt0)
680 gmm::copy(gmm::scaled(gmm::mat_row(DN, i), -vt1),
681 gmm::mat_col(T_n_n, i));
683 if (Hughes_stabilized) {
684 model_real_sparse_matrix aux(nbc, nbc);
688 gmm::copy(gmm::transposed(T_u1_n), T_n_u1);
689 if (two_variables)
gmm::copy(gmm::transposed(T_u2_n), T_n_u2);
692 gmm::sub_interval SUBI(i*d, d);
693 scalar_type th = Tresca_version ? threshold[i]
694 : - (std::min(vt0, RLN[i])) * friction_coeff[i];
695 ball_projection_grad(gmm::sub_vector(RLT, SUBI), th, pg);
696 if (!really_stationary)
699 gmm::add(gmm::scaled(gmm::mat_row(BT1,i*d+k1),-pg(k2,k1)),
700 gmm::mat_col(T_u1_t, i*d+k2));
702 gmm::add(gmm::scaled(gmm::mat_row(BT2,i*d+k1),
704 gmm::mat_col(T_u2_t, i*d+k2));
707 if (!Tresca_version) {
709 ball_projection_grad_r(gmm::sub_vector(RLT, SUBI), th, vg);
711 T_t_n(i*d+k, i) = - friction_coeff[i] * vg[k]/(r*
alpha[i]);
714 for (
size_type k = 0; k < d; ++k) pg(k,k) -= vt1;
715 gmm::copy(gmm::scaled(pg, vt1/(r*alpha[i])),
716 gmm::sub_matrix(T_t_t,SUBI));
717 if (Hughes_stabilized) {
720 gmm::add(gmm::scaled(gmm::mat_row(DT, d*i+l), -pg(k,l)),
721 gmm::mat_col(T_t_t, d*i+k));
726 if (Hughes_stabilized) {
727 model_real_sparse_matrix aux(gmm::mat_nrows(T_t_t),
728 gmm::mat_nrows(T_t_t));
732 gmm::copy(gmm::transposed(T_u1_t), T_t_u1);
733 if (two_variables)
gmm::copy(gmm::transposed(T_u2_t), T_t_u2);
736 if (augmentation_version == 1) {
737 gmm::copy(gmm::scaled(gmm::transposed(BN1), -vt1), T_u1_n);
739 gmm::copy(gmm::scaled(gmm::transposed(BN2), -vt1), T_u2_n);
741 gmm::copy(gmm::scaled(gmm::transposed(BT1), -vt1), T_u1_t);
743 gmm::copy(gmm::scaled(gmm::transposed(BT2), -vt1), T_u2_t);
746 model_real_sparse_matrix tmp1(gmm::mat_ncols(BN1),
747 gmm::mat_ncols(BN1));
748 gmm::mult(gmm::transposed(gmm::scaled(BBN1,-r)), T_n_u1, tmp1);
751 gmm::mult(gmm::transposed(gmm::scaled(BBN2,-r)), T_n_u2, tmp1);
756 gmm::mult(gmm::transposed(gmm::scaled(BBT1,-r)), T_t_u1, tmp1);
759 gmm::mult(gmm::transposed(gmm::scaled(BBT2,-r)), T_t_u2, tmp1);
765 if (!contact_only && !Tresca_version) {
767 model_real_sparse_matrix tmp5(gmm::mat_ncols(BT1),
768 gmm::mat_ncols(BT1));
769 model_real_sparse_matrix tmp6(gmm::mat_ncols(BT1),
770 gmm::mat_ncols(BT1));
771 model_real_sparse_matrix tmp7(gmm::mat_ncols(BT2),
772 gmm::mat_ncols(BT2));
773 model_real_sparse_matrix tmp8(gmm::mat_ncols(BT2),
774 gmm::mat_ncols(BT2));
775 model_real_sparse_matrix tmp3(gmm::mat_ncols(T_t_u1),
776 gmm::mat_nrows(T_t_u1));
777 model_real_sparse_matrix tmp4(gmm::mat_ncols(T_t_u2),
778 gmm::mat_nrows(T_t_u2));
781 gmm::sub_interval SUBI(i*d, d);
782 scalar_type th = - (std::min(vt0, RLN[i])) * friction_coeff[i];
784 ball_projection_grad_r(gmm::sub_vector(RLT, SUBI), th, vg);
786 gmm::add(gmm::scaled(gmm::mat_row(BN1,i),
787 vg[k]*friction_coeff[i]), gmm::mat_col(tmp3, i*d+k));
789 gmm::add(gmm::scaled(gmm::mat_row(BN2,i),
790 vg[k]*friction_coeff[i]), gmm::mat_col(tmp4, i*d+k));
792 if (augmentation_version == 2) {
794 gmm::add(gmm::scaled(gmm::mat_row(BT1,i*d+k),
795 vg[k]*friction_coeff[i]), gmm::mat_col(T_u1_n, i));
797 gmm::add(gmm::scaled(gmm::mat_row(BT2,i*d+k),
798 vg[k]*friction_coeff[i]), gmm::mat_col(T_u2_n, i));
800 gmm::copy(gmm::scaled(gmm::mat_row(BBT1,i*d+k),
801 -r*friction_coeff[i]*vg[k]),
802 gmm::mat_col(tmp5, i*d+k));
804 gmm::mat_col(tmp6, i*d+k));
806 gmm::copy(gmm::scaled(gmm::mat_row(BBT2,i*d+k),
807 -r*friction_coeff[i]*vg[k]),
808 gmm::mat_col(tmp7, i*d+k));
810 gmm::mat_col(tmp8, i*d+k));
817 gmm::add(gmm::transposed(tmp3), T_t_u1);
819 gmm::add(gmm::transposed(tmp4), T_t_u2);
821 if (augmentation_version == 2) {
822 model_real_sparse_matrix tmp1(gmm::mat_ncols(BN1),
823 gmm::mat_ncols(BN1));
824 gmm::mult(tmp5, gmm::transposed(tmp6), gmm::transposed(tmp1));
825 gmm::add(gmm::transposed(tmp1), T_u1_u1);
827 gmm::mult(tmp7, gmm::transposed(tmp8),gmm::transposed(tmp1));
828 gmm::add(gmm::transposed(tmp1), T_u2_u2);
835 gmm::copy(gmm::scaled(gmm::transposed(BN1), -vt1), T_u1_n);
837 gmm::copy(gmm::scaled(gmm::transposed(BN2), -vt1), T_u2_n);
839 if (lambda_n[i] > vt0) {
841 if (two_variables)
gmm::clear(gmm::mat_col(T_u2_n, i));
842 T_n_n(i, i) = -vt1/r;
845 gmm::copy(gmm::scaled(BN1, -vt1), T_n_u1);
846 if (two_variables)
gmm::copy(gmm::scaled(BN2, -r), T_n_u2);
849 gmm::sub_interval SUBI(i*d, d);
850 scalar_type th = Tresca_version ? threshold[i]
851 : gmm::neg(lambda_n[i]) * friction_coeff[i];
852 ball_projection_grad(gmm::sub_vector(lambda_t, SUBI), th, pg);
853 if (!really_stationary)
856 gmm::add(gmm::scaled(gmm::mat_row(BT1,i*d+k1),-pg(k2,k1)),
857 gmm::mat_col(T_u1_t, i*d+k2));
859 gmm::add(gmm::scaled(gmm::mat_row(BT2,i*d+k1),
861 gmm::mat_col(T_u2_t, i*d+k2));
863 if (!Tresca_version) {
864 ball_projection_grad_r(gmm::sub_vector(lambda_t, SUBI),th,vg);
866 gmm::add(gmm::scaled(gmm::mat_row(BT1,i*d+k1),
867 friction_coeff[i]*vg[k1]),
868 gmm::mat_col(T_u1_n, i));
870 gmm::add(gmm::scaled(gmm::mat_row(BT2,i*d+k1),
871 friction_coeff[i]*vg[k1]),
872 gmm::mat_col(T_u2_n, i));
873 T_t_n(i*d+k1, i) = friction_coeff[i] * vg[k1] / (r*
alpha[i]);
876 for (
size_type k = 0; k < d; ++k) pg(k,k) -= vt1;
878 gmm::copy(gmm::scaled(pg, vt1/(r*alpha[i])),
879 gmm::sub_matrix(T_t_t, SUBI));
882 gmm::copy(gmm::scaled(BT1, -vt1), T_t_u1);
883 if (two_variables)
gmm::copy(gmm::scaled(BT2, -r), T_t_u2);
888 base_small_vector x(d+1), n(d+1), u(d); n[0] = vt1;
889 base_matrix g(d+1, d+1);
890 model_real_sparse_matrix T_n_u1_transp(gmm::mat_ncols(T_n_u1), nbc);
891 model_real_sparse_matrix T_n_u2_transp(gmm::mat_ncols(T_n_u2), nbc);
898 for (
size_type j=0; j < d; ++j) x[1+j] = lambda_t[i*d+j];
899 De_Saxce_projection_grad(x, n, friction_coeff[i], g);
901 gmm::add(gmm::scaled(gmm::mat_row(BN1, i), -g(0,0)),
902 gmm::mat_col(T_u1_n, i));
904 gmm::add(gmm::scaled(gmm::mat_row(BN2, i), -g(0,0)),
905 gmm::mat_col(T_u2_n, i));
906 T_n_n(i, i) = (g(0,0) - vt1)/(r*alpha[i]);
908 gmm::copy(gmm::sub_vector(RLT, gmm::sub_interval(i*d,d)), u);
912 gmm::add(gmm::scaled(gmm::mat_row(BT1, i*d+j),
913 friction_coeff[i] * u[j] / nu),
914 gmm::mat_col(T_n_u1_transp, i));
916 gmm::add(gmm::scaled(gmm::mat_row(BT2, i*d+j),
917 friction_coeff[i] * u[j] / nu),
918 gmm::mat_col(T_n_u2_transp, i));
922 gmm::add(gmm::scaled(gmm::mat_row(BT1, i*d+j), -g(0,j+1)),
923 gmm::mat_col(T_u1_n, i));
925 gmm::add(gmm::scaled(gmm::mat_row(BT2, i*d+j), -g(0,j+1)),
926 gmm::mat_col(T_u2_n, i));
928 gmm::add(gmm::scaled(gmm::mat_row(BN1, i), -g(j+1,0)),
929 gmm::mat_col(T_u1_t, i*d+j));
931 gmm::add(gmm::scaled(gmm::mat_row(BN2, i), -g(j+1,0)),
932 gmm::mat_col(T_u2_t, i*d+j));
935 gmm::add(gmm::scaled(gmm::mat_row(BT1, i*d+k), -g(1+j,1+k)),
936 gmm::mat_col(T_u1_t, i*d+j));
938 gmm::add(gmm::scaled(gmm::mat_row(BT2, i*d+k), -g(1+j,1+k)),
939 gmm::mat_col(T_u2_t, i*d+j));
940 T_t_t(i*d+j, i*d+k) = g(1+j, 1+k)/r;
942 T_t_t(i*d+j, i*d+j) -= vt1/(r*
alpha[i]);
943 T_t_n(i*d+j, i) = g(1+j,0)/(r*
alpha[i]);
947 gmm::copy(gmm::scaled(BN1, -vt1), T_n_u1);
948 if (two_variables)
gmm::copy(gmm::scaled(BN2, -vt1), T_n_u2);
949 gmm::add(gmm::transposed(T_n_u1_transp), T_n_u1);
950 if (two_variables)
gmm::add(gmm::transposed(T_n_u2_transp), T_n_u2);
951 gmm::copy(gmm::scaled(BT1, -vt1), T_t_u1);
952 if (two_variables)
gmm::copy(gmm::scaled(BT2, -vt1), T_t_u2);
957 if (version & model::BUILD_RHS) {
959 switch (augmentation_version) {
962 RLN[i] = std::min(scalar_type(0), RLN[i]);
964 scalar_type radius = Tresca_version ? threshold[i]
965 : -friction_coeff[i]*RLN[i];
967 (gmm::sub_vector(RLT, gmm::sub_interval(i*d,d)), radius);
979 rlambda_n[i] = (lambda_n[i] - RLN[i]) / (r * alpha[i]);
983 = (lambda_t[i*d+k] - RLT[i*d+k]) / (r *
alpha[i]);
988 RLN[i] = std::min(vt0, RLN[i]);
990 scalar_type radius = Tresca_version ? threshold[i]
991 : -friction_coeff[i]*RLN[i];
993 (gmm::sub_vector(RLT, gmm::sub_interval(i*d,d)), radius);
997 if (two_variables)
gmm::mult_add(gmm::transposed(BN2), RLN, ru2);
1000 if (two_variables)
gmm::mult_add(gmm::transposed(BT2), RLT, ru2);
1003 rlambda_n[i] = (lambda_n[i] - RLN[i]) / (r * alpha[i]);
1007 = (lambda_t[i*d+k] - RLT[i*d+k]) / (r *
alpha[i]);
1011 if (!contact_only)
gmm::copy(lambda_t, RLT);
1013 RLN[i] = -gmm::neg(lambda_n[i]);
1014 rlambda_n[i] = gmm::pos(lambda_n[i])/r -
alpha[i]*gap[i];
1016 if (!contact_only) {
1017 scalar_type radius = Tresca_version ? threshold[i]
1018 : friction_coeff[i]*gmm::neg(lambda_n[i]);
1020 (gmm::sub_vector(RLT, gmm::sub_interval(i*d,d)), radius);
1023 gmm::mult(gmm::transposed(BN1), RLN, ru1);
1024 if (two_variables)
gmm::mult(gmm::transposed(BN2), RLN, ru2);
1027 if (!contact_only) {
1029 if (two_variables)
gmm::mult_add(gmm::transposed(BT2), RLT, ru2);
1030 gmm::add(gmm::scaled(lambda_t, vt1/r), gmm::scaled(RLT,-vt1/r),
1036 rlambda_n[i] /=
alpha[i];
1038 for (
size_type k = 0; k < d; ++k) rlambda_t[i*d+k] /= alpha[i];
1042 base_small_vector x(d+1), n(d+1);
1044 GMM_ASSERT1(!Tresca_version,
1045 "Augmentation version incompatible with Tresca friction law");
1051 gmm::copy(gmm::sub_vector(lambda_t, gmm::sub_interval(i*d,d)),
1052 gmm::sub_vector(x, gmm::sub_interval(1, d)));
1053 De_Saxce_projection(x, n, friction_coeff[i]);
1055 gmm::copy(gmm::sub_vector(x, gmm::sub_interval(1, d)),
1056 gmm::sub_vector(RLT, gmm::sub_interval(i*d,d)));
1057 rlambda_n[i] = lambda_n[i]/r - x[0]/r -
alpha[i]*gap[i]
1059 gmm::sub_interval(i*d,d)));
1062 if (two_variables)
gmm::mult_add(gmm::transposed(BT2), RLT, ru2);
1064 if (two_variables)
gmm::mult_add(gmm::transposed(BN2), RLN, ru2);
1065 gmm::add(gmm::scaled(lambda_t, vt1/r), rlambda_t);
1066 gmm::add(gmm::scaled(RLT, -vt1/r), rlambda_t);
1070 rlambda_n[i] /=
alpha[i];
1072 for (
size_type k = 0; k < d; ++k) rlambda_t[i*d+k] /= alpha[i];
1080 virtual void asm_real_tangent_terms(
const model &md,
size_type ib,
1081 const model::varnamelist &vl,
1082 const model::varnamelist &dl,
1083 const model::mimlist &mims,
1084 model::real_matlist &matl,
1085 model::real_veclist &vecl,
1086 model::real_veclist &,
1088 build_version version)
const {
1089 if (MPI_IS_MASTER()) {
1090 GMM_ASSERT1(mims.size() == 0,
"Contact brick need no mesh_im");
1091 size_type nbvar = 2 + (contact_only ? 0 : 1) + (two_variables ? 1 : 0);
1092 GMM_ASSERT1(vl.size() == nbvar,
1093 "Wrong number of variables for contact brick");
1094 size_type nbdl = 3 + (contact_only ? 0 : 1) + (Tresca_version ? 1 : 0)
1095 + (friction_dynamic_term ? 2 : 0);
1096 GMM_ASSERT1(dl.size() == nbdl,
"Wrong number of data for contact brick, "
1097 << dl.size() <<
" should be " << nbdl);
1107 const model_real_plain_vector &u1 = md.real_variable(vl[nv++]);
1108 const model_real_plain_vector &u2 = md.real_variable(vl[nv++]);
1109 if (!two_variables) nv--;
1110 const model_real_plain_vector &lambda_n = md.real_variable(vl[nv++]);
1111 if (contact_only) nv--;
1112 const model_real_plain_vector &lambda_t = md.real_variable(vl[nv]);
1116 size_type np = 0, np_wt1 = 0, np_wt2 = 0, np_alpha = 0;
1117 const model_real_plain_vector &vr = md.real_variable(dl[np++]);
1118 GMM_ASSERT1(gmm::vect_size(vr) == 1,
"Parameter r should be a scalar");
1120 const model_real_plain_vector &vgap = md.real_variable(dl[np++]);
1121 GMM_ASSERT1(gmm::vect_size(vgap) == 1 || gmm::vect_size(vgap) == nbc,
1122 "Parameter gap has a wrong size");
1124 if (gmm::vect_size(vgap) == 1)
1129 const model_real_plain_vector &valpha = md.real_variable(dl[np_alpha]);
1130 GMM_ASSERT1(gmm::vect_size(valpha)== 1 || gmm::vect_size(valpha) == nbc,
1131 "Parameter alpha has a wrong size");
1133 if (gmm::vect_size(valpha) == 1)
1137 if (!contact_only) {
1138 const model_real_plain_vector &vfr = md.real_variable(dl[np++]);
1139 GMM_ASSERT1(gmm::vect_size(vfr)==1 || gmm::vect_size(vfr) == nbc,
1140 "Parameter friction_coeff has a wrong size");
1142 if (gmm::vect_size(vfr) == 1)
1146 if (friction_dynamic_term) {
1147 const model_real_plain_vector &vg = md.real_variable(dl[np++]);
1148 GMM_ASSERT1(gmm::vect_size(vg) == 1,
1149 "Parameter gamma should be a scalar");
1152 if (two_variables) np_wt2 = np++;
1154 if (Tresca_version) {
1155 const model_real_plain_vector &vth = md.real_variable(dl[np++]);
1156 GMM_ASSERT1(gmm::vect_size(vth) == 1 || gmm::vect_size(vth) == nbc,
1157 "Parameter threshold has a wrong size");
1159 if (gmm::vect_size(vth) == 1)
1166 if (md.is_var_newer_than_brick(dl[np_alpha], ib)) is_init =
false;
1168 basic_asm_real_tangent_terms
1169 (u1, u2, lambda_n, lambda_t, md.real_variable(dl[np_wt1]),
1170 md.real_variable(dl[np_wt2]), matl, vecl, version);
1175 Coulomb_friction_brick(
int aug_version,
bool contact_only_,
1176 bool two_variables_=
false,
1177 bool Tresca_version_=
false,
1178 bool Hughes_stabilized_=
false,
1179 bool friction_dynamic_term_=
false) {
1181 if (aug_version == 4 && contact_only_) aug_version = 3;
1182 augmentation_version = aug_version;
1183 GMM_ASSERT1(aug_version >= 1 && aug_version <= 4,
1184 "Wrong augmentation version");
1185 GMM_ASSERT1(!Hughes_stabilized_ || aug_version <= 2,
1186 "The Hughes stabilized version is only for Alart-Curnier "
1188 contact_only = contact_only_;
1190 Tresca_version = Tresca_version_;
1191 really_stationary =
false;
1192 friction_dynamic_term = friction_dynamic_term_;
1193 two_variables = two_variables_;
1194 Hughes_stabilized = Hughes_stabilized_;
1195 set_flags(
"Coulomb friction brick",
false ,
1197 (augmentation_version == 2) && (contact_only||Tresca_version),
1202 void set_BN1(CONTACT_B_MATRIX &BN1_) {
1203 gmm::resize(BN1, gmm::mat_nrows(BN1_), gmm::mat_ncols(BN1_));
1208 void set_BN2(CONTACT_B_MATRIX &BN2_) {
1209 gmm::resize(BN2, gmm::mat_nrows(BN2_), gmm::mat_ncols(BN2_));
1214 void set_DN(CONTACT_B_MATRIX &DN_) {
1215 gmm::resize(DN, gmm::mat_nrows(DN_), gmm::mat_ncols(DN_));
1220 void set_DT(CONTACT_B_MATRIX &DT_) {
1221 gmm::resize(DT, gmm::mat_nrows(DT_), gmm::mat_ncols(DT_));
1226 void set_BT1(CONTACT_B_MATRIX &BT1_) {
1227 gmm::resize(BT1, gmm::mat_nrows(BT1_), gmm::mat_ncols(BT1_));
1232 CONTACT_B_MATRIX &get_BN1(
void) {
return BN1; }
1233 CONTACT_B_MATRIX &get_DN(
void) {
return DN; }
1234 CONTACT_B_MATRIX &get_DT(
void) {
return DT; }
1235 CONTACT_B_MATRIX &get_BT1(
void) {
return BT1; }
1236 const CONTACT_B_MATRIX &get_BN1(
void)
const {
return BN1; }
1237 const CONTACT_B_MATRIX &get_DN(
void)
const {
return DN; }
1238 const CONTACT_B_MATRIX &get_BT1(
void)
const {
return BT1; }
1244 pbrick pbr = md.brick_pointer(indbrick);
1246 Coulomb_friction_brick *p =
dynamic_cast<Coulomb_friction_brick *
>
1248 GMM_ASSERT1(p,
"Wrong type of brick");
1249 p->set_stationary();
1254 pbrick pbr = md.brick_pointer(indbrick);
1256 Coulomb_friction_brick *p =
dynamic_cast<Coulomb_friction_brick *
>
1258 GMM_ASSERT1(p,
"Wrong type of brick");
1259 return p->get_BN1();
1265 pbrick pbr = md.brick_pointer(indbrick);
1267 Coulomb_friction_brick *p =
dynamic_cast<Coulomb_friction_brick *
>
1269 GMM_ASSERT1(p,
"Wrong type of brick");
1275 pbrick pbr = md.brick_pointer(indbrick);
1277 Coulomb_friction_brick *p =
dynamic_cast<Coulomb_friction_brick *
>
1279 GMM_ASSERT1(p,
"Wrong type of brick");
1286 pbrick pbr = md.brick_pointer(indbrick);
1288 Coulomb_friction_brick *p =
dynamic_cast<Coulomb_friction_brick *
>
1290 GMM_ASSERT1(p,
"Wrong type of brick");
1291 return p->get_BT1();
1299 (
model &md,
const std::string &varname_u,
const std::string &multname_n,
1300 const std::string &dataname_r, CONTACT_B_MATRIX &BN,
1301 std::string dataname_gap, std::string dataname_alpha,
1302 int aug_version,
bool Hughes_stabilized) {
1304 auto pbr_ = std::make_shared<Coulomb_friction_brick>
1305 (aug_version,
true,
false,
false, Hughes_stabilized);
1310 tl.push_back(model::term_description(varname_u, varname_u,
false));
1311 tl.push_back(model::term_description(varname_u, multname_n,
false));
1312 tl.push_back(model::term_description(multname_n, varname_u,
false));
1313 tl.push_back(model::term_description(multname_n, multname_n,
false));
1314 model::varnamelist dl(1, dataname_r);
1316 if (dataname_gap.size() == 0) {
1317 dataname_gap = md.
new_name(
"contact_gap_on_" + varname_u);
1319 (dataname_gap, model_real_plain_vector(1, scalar_type(0)));
1321 dl.push_back(dataname_gap);
1323 if (dataname_alpha.size() == 0) {
1324 dataname_alpha = md.
new_name(
"contact_parameter_alpha_on_"+ multname_n);
1326 (dataname_alpha, model_real_plain_vector(1, scalar_type(1)));
1328 dl.push_back(dataname_alpha);
1330 model::varnamelist vl(1, varname_u);
1331 vl.push_back(multname_n);
1342 (
model &md,
const std::string &varname_u1,
const std::string &varname_u2,
1343 const std::string &multname_n,
1344 const std::string &dataname_r, CONTACT_B_MATRIX &BN1, CONTACT_B_MATRIX &BN2,
1345 std::string dataname_gap, std::string dataname_alpha,
1346 int aug_version,
bool Hughes_stabilized) {
1348 auto pbr_ = std::make_shared<Coulomb_friction_brick>
1349 (aug_version,
true,
true,
false, Hughes_stabilized);
1355 tl.push_back(model::term_description(varname_u1, varname_u1,
false));
1356 tl.push_back(model::term_description(varname_u2, varname_u2,
false));
1357 tl.push_back(model::term_description(varname_u1, multname_n,
false));
1358 tl.push_back(model::term_description(multname_n, varname_u1,
false));
1359 tl.push_back(model::term_description(varname_u2, multname_n,
false));
1360 tl.push_back(model::term_description(multname_n, varname_u2,
false));
1361 tl.push_back(model::term_description(multname_n, multname_n,
false));
1362 model::varnamelist dl(1, dataname_r);
1364 if (dataname_gap.size() == 0) {
1365 dataname_gap = md.
new_name(
"contact_gap_on_" + varname_u1);
1367 (dataname_gap, model_real_plain_vector(1, scalar_type(0)));
1369 dl.push_back(dataname_gap);
1371 if (dataname_alpha.size() == 0) {
1372 dataname_alpha = md.
new_name(
"contact_parameter_alpha_on_"+ multname_n);
1374 (dataname_alpha, model_real_plain_vector(1, scalar_type(1)));
1376 dl.push_back(dataname_alpha);
1378 model::varnamelist vl(1, varname_u1);
1379 vl.push_back(varname_u2);
1380 vl.push_back(multname_n);
1391 (
model &md,
const std::string &varname_u,
const std::string &multname_n,
1392 const std::string &multname_t,
const std::string &dataname_r,
1393 CONTACT_B_MATRIX &BN, CONTACT_B_MATRIX &BT,
1394 std::string dataname_friction_coeff,
1395 std::string dataname_gap, std::string dataname_alpha,
1396 int aug_version,
bool Tresca_version,
const std::string dataname_threshold,
1397 std::string dataname_gamma, std::string dataname_wt,
bool Hughes_stabilized) {
1399 bool dynamic_terms = (dataname_gamma.size() > 0);
1401 auto pbr_ = std::make_shared<Coulomb_friction_brick>
1402 (aug_version,
false,
false, Tresca_version, Hughes_stabilized,
1409 tl.push_back(model::term_description(varname_u, varname_u,
false));
1410 tl.push_back(model::term_description(varname_u, multname_n,
false));
1411 tl.push_back(model::term_description(multname_n, varname_u,
false));
1412 tl.push_back(model::term_description(multname_n, multname_n,
false));
1413 tl.push_back(model::term_description(varname_u, multname_t,
false));
1414 tl.push_back(model::term_description(multname_t, varname_u,
false));
1415 tl.push_back(model::term_description(multname_t, multname_t,
false));
1416 tl.push_back(model::term_description(multname_t, multname_n,
1417 (aug_version == 4)));
1418 model::varnamelist dl(1, dataname_r);
1419 if (dataname_gap.size() == 0) {
1420 dataname_gap = md.
new_name(
"contact_gap_on_" + varname_u);
1422 (dataname_gap, model_real_plain_vector(1, scalar_type(0)));
1424 dl.push_back(dataname_gap);
1426 if (dataname_alpha.size() == 0) {
1427 dataname_alpha = md.
new_name(
"contact_parameter_alpha_on_"+ multname_n);
1429 (dataname_alpha, model_real_plain_vector(1, scalar_type(1)));
1431 dl.push_back(dataname_alpha);
1432 dl.push_back(dataname_friction_coeff);
1433 if (dataname_gamma.size()) {
1434 dl.push_back(dataname_gamma);
1435 dl.push_back(dataname_wt);
1438 dl.push_back(dataname_threshold);
1440 model::varnamelist vl(1, varname_u);
1441 vl.push_back(multname_n);
1442 vl.push_back(multname_t);
1455 struct Coulomb_friction_brick_rigid_obstacle
1456 :
public Coulomb_friction_brick {
1458 std::string obstacle;
1462 virtual void asm_real_tangent_terms(
const model &md,
size_type ib,
1463 const model::varnamelist &vl,
1464 const model::varnamelist &dl,
1465 const model::mimlist &mims,
1466 model::real_matlist &matl,
1467 model::real_veclist &vecl,
1468 model::real_veclist &,
1470 build_version version)
const {
1471 if (MPI_IS_MASTER()) {
1472 GMM_ASSERT1(mims.size() == 1,
"This contact brick needs one mesh_im");
1473 size_type nbvar = 2 + (contact_only ? 0 : 1);
1474 GMM_ASSERT1(vl.size() == nbvar,
1475 "Wrong number of variables for contact brick: "
1476 << vl.size() <<
" should be " << nbvar);
1477 size_type nbdl = 1 + (contact_only ? 0 : 1) + (Tresca_version ? 1 : 0)
1478 + (friction_dynamic_term ? 1 : 0);
1479 GMM_ASSERT1(dl.size() == nbdl,
1480 "Wrong number of data for contact brick: "
1481 << dl.size() <<
" should be " << nbdl);
1482 GMM_ASSERT1(!two_variables,
"internal error");
1483 const mesh_im &mim = *mims[0];
1489 const model_real_plain_vector &u1 = md.real_variable(vl[nv++]);
1490 const mesh_fem &mf_u1 = md.mesh_fem_of_variable(vl[0]);
1491 const model_real_plain_vector &lambda_n = md.real_variable(vl[nv++]);
1492 if (contact_only) nv--;
1493 const model_real_plain_vector &lambda_t = md.real_variable(vl[nv]);
1498 const model_real_plain_vector &vr = md.real_variable(dl[np++]);
1499 GMM_ASSERT1(gmm::vect_size(vr) == 1,
"Parameter r should be a scalar");
1503 if (md.is_var_mf_newer_than_brick(vl[0], ib)) {
1506 GMM_ASSERT1(!(mf_u1.is_reduced()),
1507 "This contact brick works only for pure Lagrange fems");
1508 dal::bit_vector dofs = mf_u1.basic_dof_on_region(region);
1509 for (dal::bv_visitor
id(dofs); !
id.finished(); ++id) {
1510 size_type cv = mf_u1.first_convex_of_basic_dof(
id);
1511 GMM_ASSERT1(mf_u1.fem_of_element(cv)->is_lagrange(),
1512 "This contact brick works only for pure Lagrange fems");
1514 size_type d = mf_u1.get_qdim() - 1, i = 0, j = 0;
1515 nbc = dofs.card() / (d+1);
1518 base_node Pmin, Pmax;
1519 mf_u1.linked_mesh().bounding_box(Pmin, Pmax);
1520 scalar_type l = scalar_type(0);
1521 for (i = 0; i < Pmin.size(); ++i)
1522 l = std::max(l, gmm::abs(Pmax[i] - Pmin[i]));
1524 CONTACT_B_MATRIX MM(mf_u1.nb_dof(), mf_u1.nb_dof());
1528 for (dal::bv_visitor
id(dofs); !
id.finished(); ++id, ++i)
1529 if ((i % (d+1)) == 0)
alpha[j++] = MM(
id,
id) / l;
1531 getfem::ga_workspace gw;
1533 getfem::model_real_plain_vector pt(N);
1534 gw.add_fixed_size_constant(
"X", pt);
1535 if (N >= 1) gw.add_macro(
"x",
"X(1)");
1536 if (N >= 2) gw.add_macro(
"y",
"X(2)");
1537 if (N >= 3) gw.add_macro(
"z",
"X(3)");
1538 if (N >= 4) gw.add_macro(
"w",
"X(4)");
1540 getfem::ga_function f(gw, obstacle);
1546 if (!contact_only) {
1550 base_node grad(d+1), ut[3];
1553 for (dal::bv_visitor
id(dofs); !
id.finished(); ++id, ++i) {
1554 if ((i % (d+1)) == 0) {
1555 gmm::copy(mf_u1.point_of_basic_dof(
id), pt);
1558 gap[j] = (f.eval())[0];
1561 size_type cv = mf_u1.first_convex_of_basic_dof(
id);
1563 = mf_u1.linked_mesh().convex_radius_estimate(cv) * 1E-3;
1566 grad[k] = ((f.eval())[0] - gap[j]) / eps;
1573 BN1(j,
id + k) = un[k];
1576 if (!contact_only) {
1578 orthonormal_basis_to_unit_vec(d, un, ut);
1582 BT1(j*d+nn,
id + k) = ut[nn][k];
1589 GMM_ASSERT1(gmm::vect_size(md.real_variable(vl[1])) == nbc,
1590 "Wrong size of multiplier for the contact condition");
1593 GMM_ASSERT1(gmm::vect_size(md.real_variable(vl[2])) == nbc*d,
1594 "Wrong size of multiplier for the friction condition");
1599 nbc = gmm::mat_nrows(BN1);
1601 if (!contact_only) {
1602 const model_real_plain_vector &vfr = md.real_variable(dl[np++]);
1603 GMM_ASSERT1(gmm::vect_size(vfr)==1 || gmm::vect_size(vfr) == nbc,
1604 "Parameter friction_coeff has a wrong size");
1606 if (gmm::vect_size(vfr) == 1)
1610 if (friction_dynamic_term) {
1611 const model_real_plain_vector &vg = md.real_variable(dl[np++]);
1612 GMM_ASSERT1(gmm::vect_size(vg) == 1,
1613 "Parameter gamma should be a scalar");
1617 if (Tresca_version) {
1618 const model_real_plain_vector &vth = md.real_variable(dl[np++]);
1619 GMM_ASSERT1(gmm::vect_size(vth) == 1 || gmm::vect_size(vth) == nbc,
1620 "Parameter threshold has a wrong size");
1622 if (gmm::vect_size(vth) == 1)
1629 basic_asm_real_tangent_terms
1630 (u1, u1, lambda_n, lambda_t, md.real_variable(dl[np_wt1]),
1631 md.real_variable(dl[np_wt1]), matl, vecl, version);
1636 Coulomb_friction_brick_rigid_obstacle
1637 (
int aug_version,
bool contact_only_,
const std::string &obs)
1638 : Coulomb_friction_brick(aug_version, contact_only_), obstacle(obs) {}
1649 (
model &md,
const mesh_im &mim,
const std::string &varname_u,
1650 const std::string &multname_n,
const std::string &dataname_r,
1651 size_type region,
const std::string &obstacle,
int aug_version) {
1652 pbrick pbr = std::make_shared<Coulomb_friction_brick_rigid_obstacle>
1653 (aug_version,
true, obstacle);
1656 tl.push_back(model::term_description(varname_u, varname_u,
false));
1657 tl.push_back(model::term_description(varname_u, multname_n,
false));
1658 tl.push_back(model::term_description(multname_n, varname_u,
false));
1659 tl.push_back(model::term_description(multname_n, multname_n,
false));
1660 model::varnamelist dl(1, dataname_r);
1662 model::varnamelist vl(1, varname_u);
1663 vl.push_back(multname_n);
1665 return md.
add_brick(pbr, vl, dl, tl, model::mimlist(1, &mim), region);
1675 (
model &md,
const mesh_im &mim,
const std::string &varname_u,
1676 const std::string &multname_n,
const std::string &multname_t,
1677 const std::string &dataname_r,
const std::string &dataname_friction_coeff,
1678 size_type region,
const std::string &obstacle,
int aug_version) {
1679 pbrick pbr = std::make_shared<Coulomb_friction_brick_rigid_obstacle>
1680 (aug_version,
false, obstacle);
1683 tl.push_back(model::term_description(varname_u, varname_u,
false));
1684 tl.push_back(model::term_description(varname_u, multname_n,
false));
1685 tl.push_back(model::term_description(multname_n, varname_u,
false));
1686 tl.push_back(model::term_description(multname_n, multname_n,
false));
1687 tl.push_back(model::term_description(varname_u, multname_t,
false));
1688 tl.push_back(model::term_description(multname_t, varname_u,
false));
1689 tl.push_back(model::term_description(multname_t, multname_t,
false));
1690 tl.push_back(model::term_description(multname_t, multname_n,
1691 (aug_version == 4)));
1692 model::varnamelist dl(1, dataname_r);
1693 dl.push_back(dataname_friction_coeff);
1695 model::varnamelist vl(1, varname_u);
1696 vl.push_back(multname_n);
1697 vl.push_back(multname_t);
1699 return md.
add_brick(pbr, vl, dl, tl, model::mimlist(1, &mim), region);
1712 struct Coulomb_friction_brick_nonmatching_meshes
1713 :
public Coulomb_friction_brick {
1715 std::vector<size_type> rg1, rg2;
1719 bool slave1, slave2;
1723 virtual void asm_real_tangent_terms(
const model &md,
size_type ib,
1724 const model::varnamelist &vl,
1725 const model::varnamelist &dl,
1726 const model::mimlist &mims,
1727 model::real_matlist &matl,
1728 model::real_veclist &vecl,
1729 model::real_veclist &,
1731 build_version version)
const {
1732 if (MPI_IS_MASTER()) {
1733 GMM_ASSERT1(mims.size() == 2,
"This contact brick needs two mesh_im");
1734 const mesh_im &mim1 = *mims[0];
1735 const mesh_im &mim2 = *mims[1];
1743 std::string varname_u1 = vl[nv];
1744 const model_real_plain_vector &u1 = md.real_variable(varname_u1);
1745 const mesh_fem &mf_u1 = md.mesh_fem_of_variable(varname_u1);
1746 if (two_variables) nv++;
1747 std::string varname_u2 = vl[nv++];
1748 const model_real_plain_vector &u2 = md.real_variable(varname_u2);
1749 const mesh_fem &mf_u2 = md.mesh_fem_of_variable(varname_u2);
1750 const model_real_plain_vector &lambda_n = md.real_variable(vl[nv]);
1751 if (!contact_only) nv++;
1752 const model_real_plain_vector &lambda_t = md.real_variable(vl[nv]);
1757 const model_real_plain_vector &vr = md.real_variable(dl[0]);
1758 GMM_ASSERT1(gmm::vect_size(vr) == 1,
"Parameter r should be a scalar");
1760 if (!contact_only) {
1761 const model_real_plain_vector &vfr = md.real_variable(dl[1]);
1762 GMM_ASSERT1(gmm::vect_size(vfr)==1 || gmm::vect_size(vfr) == nbc,
1763 "Parameter friction_coeff has a wrong size");
1764 gmm::resize(friction_coeff, nbc);
1765 if (gmm::vect_size(vfr) == 1)
1766 gmm::fill(friction_coeff, vfr[0]);
1768 gmm::copy(vfr, friction_coeff);
1772 if ( md.is_var_mf_newer_than_brick(varname_u1, ib)
1773 || md.is_var_mf_newer_than_brick(varname_u2, ib)) {
1776 const mesh_fem &mf_u = i ? mf_u2 : mf_u1;
1778 GMM_ASSERT1(!(mf_u.is_reduced()),
1779 "This contact brick works only for pure Lagrange fems");
1780 dal::bit_vector dofs = mf_u.basic_dof_on_region(region);
1781 for (dal::bv_visitor
id(dofs); !
id.finished(); ++id) {
1782 size_type cv = mf_u.first_convex_of_basic_dof(
id);
1783 GMM_ASSERT1(mf_u.fem_of_element(cv)->is_lagrange(),
1784 "This contact brick works only for pure Lagrange fems");
1788 contact_node_pair_list cnpl;
1789 for (
size_type it = 0; it < rg1.size() && it < rg2.size(); ++it)
1790 cnpl.append_min_dist_cn_pairs
1791 (mf_u1, mf_u2, rg1[it], rg2[it], slave1, slave2);
1797 if (!two_variables) {
1798 compute_contact_matrices(mf_u1, mf_u2, cnpl, gap, &BN1);
1801 compute_contact_matrices(mf_u1, mf_u2, cnpl, gap, &BN1, &BN2);
1806 if (!two_variables) {
1807 compute_contact_matrices(mf_u1, mf_u2, cnpl, gap, &BN1, 0, &BT1);
1812 compute_contact_matrices(mf_u1, mf_u2, cnpl, gap, &BN1, &BN2, &BT1, &BT2);
1817 scalar_type l = scalar_type(0);
1819 const mesh_fem &mf_u = i ? mf_u2 : mf_u1;
1820 base_node Pmin, Pmax;
1821 mf_u.linked_mesh().bounding_box(Pmin, Pmax);
1822 for (
size_type j = 0; j < Pmin.size(); ++j)
1823 l = std::max(l, gmm::abs(Pmax[j] - Pmin[j]));
1827 for (
size_type it = 0; it < rg1.size() && it < rg2.size(); ++it) {
1828 for (
size_type swap = 0; swap <= 1; ++swap) {
1829 if (swap ? slave2 : slave1) {
1830 size_type rg = swap ? rg2[it] : rg1[it];
1831 const mesh_fem &mf_u = swap ? mf_u2 : mf_u1;
1832 const mesh_im &mim = swap ? mim2 : mim1;
1833 CONTACT_B_MATRIX MM(mf_u.nb_dof(), mf_u.nb_dof());
1836 dal::bit_vector rg_dofs = mf_u.basic_dof_on_region(rg);
1837 for (dal::bv_visitor
id(rg_dofs); !
id.finished(); ++id)
1838 if (
id % qdim == 0)
alpha[mult_id++] = MM(
id,
id) / l;
1844 const model_real_plain_vector dummy_wt;
1845 basic_asm_real_tangent_terms
1846 (u1, u2, lambda_n, lambda_t, dummy_wt, dummy_wt, matl, vecl, version);
1850 Coulomb_friction_brick_nonmatching_meshes
1851 (
int aug_version,
bool contact_only_,
bool two_variables_,
1852 const std::vector<size_type> &rg1_,
const std::vector<size_type> &rg2_,
1853 bool slave1_=
true,
bool slave2_=
false)
1854 : Coulomb_friction_brick(aug_version, contact_only_, two_variables_),
1855 rg1(rg1_), rg2(rg2_), slave1(slave1_), slave2(slave2_) {}
1867 const std::string &varname_u1,
const std::string &varname_u2,
1868 std::string &multname_n,
const std::string &dataname_r,
1869 const std::vector<size_type> &rg1,
const std::vector<size_type> &rg2,
1870 bool slave1,
bool slave2,
int aug_version) {
1872 bool two_variables = (varname_u1.compare(varname_u2) != 0);
1874 pbrick pbr = std::make_shared<Coulomb_friction_brick_nonmatching_meshes>
1875 (aug_version,
true, two_variables, rg1, rg2, slave1, slave2);
1881 for (
size_type it = 0; it < rg1.size() && it < rg2.size(); ++it) {
1882 for (
size_type swap = 0; swap <= 1; ++swap) {
1883 if (swap ? slave2 : slave1) {
1884 const mesh_fem &mf = swap ? mf_u2 : mf_u1;
1885 size_type rg = swap ? rg2[it] : rg1[it];
1887 nbc += rg_dofs.card() / mf.
get_qdim();
1892 if (multname_n.size() == 0)
1893 multname_n = md.
new_name(
"contact_multiplier");
1895 GMM_ASSERT1(multname_n.compare(md.
new_name(multname_n)) == 0,
1896 "The given name for the multiplier is already reserved in the model");
1900 tl.push_back(model::term_description(varname_u1, varname_u1,
false));
1901 if (two_variables) {
1902 tl.push_back(model::term_description(varname_u2, varname_u2,
false));
1904 tl.push_back(model::term_description(varname_u1, multname_n,
false));
1905 tl.push_back(model::term_description(multname_n, varname_u1,
false));
1906 if (two_variables) {
1907 tl.push_back(model::term_description(varname_u2, multname_n,
false));
1908 tl.push_back(model::term_description(multname_n, varname_u2,
false));
1910 tl.push_back(model::term_description(multname_n, multname_n,
false));
1913 model::varnamelist vl;
1914 vl.push_back(varname_u1);
1915 if (two_variables) vl.push_back(varname_u2);
1916 vl.push_back(multname_n);
1919 model::varnamelist dl;
1920 dl.push_back(dataname_r);
1923 ml.push_back(&mim1);
1924 ml.push_back(&mim2);
1937 const std::string &varname_u1,
const std::string &varname_u2,
1938 std::string &multname_n, std::string &multname_t,
1939 const std::string &dataname_r,
const std::string &dataname_friction_coeff,
1940 const std::vector<size_type> &rg1,
const std::vector<size_type> &rg2,
1941 bool slave1,
bool slave2,
int aug_version) {
1943 bool two_variables = (varname_u1.compare(varname_u2) != 0);
1945 pbrick pbr = std::make_shared<Coulomb_friction_brick_nonmatching_meshes>
1946 (aug_version,
false, two_variables, rg1, rg2, slave1, slave2);
1952 for (
size_type it = 0; it < rg1.size() && it < rg2.size(); ++it) {
1953 for (
size_type swap = 0; swap <= 1; ++swap) {
1954 if (swap ? slave2 : slave1) {
1955 const mesh_fem &mf = swap ? mf_u2 : mf_u1;
1956 size_type rg = swap ? rg2[it] : rg1[it];
1958 nbc += rg_dofs.card() / mf.
get_qdim();
1963 if (multname_n.size() == 0)
1964 multname_n = md.
new_name(
"contact_normal_multiplier");
1966 GMM_ASSERT1(multname_n.compare(md.
new_name(multname_n)) == 0,
1967 "The given name for the multiplier is already reserved in the model");
1969 if (multname_t.size() == 0)
1970 multname_t = md.
new_name(
"contact_tangent_multiplier");
1972 GMM_ASSERT1(multname_t.compare(md.
new_name(multname_t)) == 0,
1973 "The given name for the multiplier is already reserved in the model");
1977 tl.push_back(model::term_description(varname_u1, varname_u1,
false));
1978 if (two_variables) {
1979 tl.push_back(model::term_description(varname_u2, varname_u2,
false));
1982 tl.push_back(model::term_description(varname_u1, multname_n,
false));
1983 tl.push_back(model::term_description(multname_n, varname_u1,
false));
1984 if (two_variables) {
1985 tl.push_back(model::term_description(varname_u2, multname_n,
false));
1986 tl.push_back(model::term_description(multname_n, varname_u2,
false));
1988 tl.push_back(model::term_description(multname_n, multname_n,
false));
1990 tl.push_back(model::term_description(varname_u1, multname_t,
false));
1991 tl.push_back(model::term_description(multname_t, varname_u1,
false));
1992 if (two_variables) {
1993 tl.push_back(model::term_description(varname_u2, multname_t,
false));
1994 tl.push_back(model::term_description(multname_t, varname_u2,
false));
1996 tl.push_back(model::term_description(multname_t, multname_t,
false));
1997 tl.push_back(model::term_description(multname_t, multname_n,
1998 (aug_version == 4)));
2001 model::varnamelist vl;
2002 vl.push_back(varname_u1);
2003 if (two_variables) vl.push_back(varname_u2);
2004 vl.push_back(multname_n);
2005 vl.push_back(multname_t);
2008 model::varnamelist dl;
2009 dl.push_back(dataname_r);
2010 dl.push_back(dataname_friction_coeff);
2013 ml.push_back(&mim1);
2014 ml.push_back(&mim2);
Simple implementation of a KD-tree.
dref_convex_pt_ct dir_points_of_face(short_type i) const
Direct points for a given face.
does the inversion of the geometric transformation for a given convex
bool invert(const base_node &n, base_node &n_ref, scalar_type IN_EPS=1e-12, bool project_into_element=false)
given the node on the real element, returns the node on the reference element (even if it is outside ...
Balanced tree over a set of points.
void add_point_with_id(const base_node &n, size_type i)
insert a new point, with an associated number.
Describe a finite element method linked to a mesh.
virtual dim_type get_qdim() const
Return the Q dimension.
virtual dal::bit_vector basic_dof_on_region(const mesh_region &b) const
Get a list of basic dof lying on a given mesh_region.
Describe an integration method linked to a mesh.
`‘Model’' variables store the variables, the data and the description of a model.
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.
const mesh_fem & mesh_fem_of_variable(const std::string &name) const
Gives the access to the mesh_fem of a variable if any.
void add_initialized_fixed_size_data(const std::string &name, const VECT &v)
Add a fixed size data (assumed to be a vector) to the model and initialized with v.
void add_fixed_size_variable(const std::string &name, size_type size, size_type niter=1)
Add a fixed size variable to the model assumed to be a vector.
std::string new_name(const std::string &name)
Gives a non already existing variable name begining by name.
void touch_brick(size_type ib)
Force the re-computation of a brick for the next assembly.
The virtual brick has to be derived to describe real model bricks.
indexed array reference (given a container X, and a set of indexes I, this class provides a pseudo-co...
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 fill(L &l, typename gmm::linalg_traits< L >::value_type x)
*/
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2_sqr(const V &v)
squared Euclidean norm of a vector.
void clear(L &l)
clear (fill with zeros) a vector or matrix.
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)
*/
linalg_traits< DenseMatrixLU >::value_type lu_det(const DenseMatrixLU &LU, const Pvector &pvector)
Compute the matrix determinant (via a LU factorization)
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)
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
gmm::uint16_type short_type
used as the common short type integer in the library
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.
size_type add_basic_contact_brick(model &md, const std::string &varname_u, const std::string &multname_n, const std::string &dataname_r, CONTACT_B_MATRIX &BN, std::string dataname_gap="", std::string dataname_alpha="", int aug_version=1, bool Hughes_stabilized=false)
Add a frictionless contact condition to the model.
size_type add_basic_contact_brick_two_deformable_bodies(model &md, const std::string &varname_u1, const std::string &varname_u2, const std::string &multname_n, const std::string &dataname_r, CONTACT_B_MATRIX &BN1, CONTACT_B_MATRIX &BN2, std::string dataname_gap="", std::string dataname_alpha="", int aug_version=1, bool Hughes_stabilized=false)
Add a frictionless contact condition to the model between two deformable bodies.
CONTACT_B_MATRIX & contact_brick_set_BN(model &md, size_type indbrick)
Can be used to change the matrix BN of a basic contact/friction brick.
void contact_brick_set_stationary(model &md, size_type indbrick)
Can be used to set the stationary option.
CONTACT_B_MATRIX & contact_brick_set_DT(model &md, size_type indbrick)
Can be used to change the matrix DT of a basic contact/friction brick.
CONTACT_B_MATRIX & contact_brick_set_DN(model &md, size_type indbrick)
Can be used to change the matrix DN of a basic contact/friction brick.
std::shared_ptr< const virtual_brick > pbrick
type of pointer on a brick
size_type add_nodal_contact_between_nonmatching_meshes_brick(model &md, const mesh_im &mim1, const mesh_im &mim2, const std::string &varname_u1, const std::string &varname_u2, std::string &multname_n, const std::string &dataname_r, const std::vector< size_type > &rg1, const std::vector< size_type > &rg2, bool slave1=true, bool slave2=false, int aug_version=1)
Add a frictionless contact condition between two faces of one or two elastic bodies.
size_type add_nodal_contact_with_rigid_obstacle_brick(model &md, const mesh_im &mim, const std::string &varname_u, const std::string &multname_n, const std::string &dataname_r, size_type region, const std::string &obstacle, int aug_version=1)
Add a frictionless contact condition with a rigid obstacle to the model.
CONTACT_B_MATRIX & contact_brick_set_BT(model &md, size_type indbrick)
Can be used to change the matrix BT of a basic contact/friction brick.
store a point and the associated index for the kdtree.