29 #ifndef OOMPH_SOLID_TRACTION_ELEMENTS_HEADER
30 #define OOMPH_SOLID_TRACTION_ELEMENTS_HEADER
34 #include <oomph-lib-config.h>
39 #include "../generic/Qelements.h"
40 #include "../generic/hermite_elements.h"
48 namespace SolidTractionElementHelper
58 unsigned n_dim = load.size();
59 for (
unsigned i = 0;
i < n_dim;
i++)
75 template<
class ELEMENT>
120 const bool& called_from_refineable_constructor =
false)
133 if (!called_from_refineable_constructor)
135 if (element_pt->
dim() == 3)
145 "This face element will not work correctly if nodes are "
146 "hanging.\nUse the refineable version instead. ",
147 OOMPH_CURRENT_FUNCTION,
148 OOMPH_EXCEPTION_LOCATION);
199 void output(std::ostream& outfile,
const unsigned& n_plot)
212 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
233 for (
unsigned i = 0;
i < n_dim;
i++)
235 outfile << x[
i] <<
" ";
239 for (
unsigned i = 0;
i < n_dim;
i++)
245 for (
unsigned i = 0;
i < n_dim;
i++)
247 outfile << unit_normal[
i] <<
" ";
250 outfile << std::endl;
264 void output(FILE* file_pt,
const unsigned& n_plot)
285 template<
class ELEMENT>
289 unsigned n_dim = this->nodal_dimension();
293 interpolated_x(
s, x);
298 this->interpolated_xi(
s, xi);
302 outer_unit_normal(
s, unit_normal);
308 get_traction(ipt, xi, x, unit_normal, traction);
315 template<
class ELEMENT>
320 unsigned n_node = nnode();
323 unsigned n_position_type = this->nnodal_position_type();
326 unsigned n_dim = this->nodal_dimension();
334 Shape psi(n_node, n_position_type);
335 DShape dpsids(n_node, n_position_type, n_dim - 1);
338 unsigned n_intpt = integral_pt()->nweight();
341 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
344 double w = integral_pt()->weight(ipt);
347 dshape_local_at_knot(ipt, psi, dpsids);
357 for (
unsigned l = 0; l < n_node; l++)
360 for (
unsigned k = 0; k < n_position_type; k++)
363 for (
unsigned i = 0;
i < n_dim;
i++)
367 nodal_position_gen(l, bulk_position_type(k),
i) * psi(l, k);
369 interpolated_xi[
i] +=
370 this->lagrangian_position_gen(l, bulk_position_type(k),
i) *
375 for (
unsigned j = 0; j < n_dim - 1; j++)
377 interpolated_A(j,
i) +=
378 nodal_position_gen(l, bulk_position_type(k),
i) *
387 for (
unsigned i = 0;
i < n_dim - 1;
i++)
389 for (
unsigned j = 0; j < n_dim - 1; j++)
394 for (
unsigned k = 0; k < n_dim; k++)
396 A(
i, j) += interpolated_A(
i, k) * interpolated_A(j, k);
403 outer_unit_normal(ipt, interpolated_normal);
413 Adet = A(0, 0) * A(1, 1) - A(0, 1) * A(1, 0);
417 "Wrong dimension in SolidTractionElement",
418 "SolidTractionElement::fill_in_contribution_to_residuals()",
419 OOMPH_EXCEPTION_LOCATION);
424 double W = w * sqrt(Adet);
429 ipt, interpolated_xi, interpolated_x, interpolated_normal, traction);
434 for (
unsigned l = 0; l < n_node; l++)
437 for (
unsigned k = 0; k < n_position_type; k++)
440 for (
unsigned i = 0;
i < n_dim;
i++)
442 local_eqn = this->position_local_eqn(l, bulk_position_type(k),
i);
447 residuals[local_eqn] -= traction[
i] * psi(l, k) *
W;
470 template<
class ELEMENT>
537 template<
class ELEMENT>
543 unsigned n_node = nnode();
546 unsigned n_position_type = this->nnodal_position_type();
550 if (n_position_type != 1)
553 "RefineableSolidTractionElement only works for n_position_type=1",
554 OOMPH_CURRENT_FUNCTION,
555 OOMPH_EXCEPTION_LOCATION);
561 unsigned n_dim = this->nodal_dimension();
569 Shape psi(n_node, n_position_type);
570 DShape dpsids(n_node, n_position_type, n_dim - 1);
573 unsigned n_intpt = integral_pt()->nweight();
576 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
579 double w = integral_pt()->weight(ipt);
582 dshape_local_at_knot(ipt, psi, dpsids);
592 for (
unsigned l = 0; l < n_node; l++)
595 for (
unsigned k = 0; k < n_position_type; k++)
598 for (
unsigned i = 0;
i < n_dim;
i++)
602 nodal_position_gen(l, this->bulk_position_type(k),
i) * psi(l, k);
604 interpolated_xi[
i] +=
605 this->lagrangian_position_gen(l, this->bulk_position_type(k),
i) *
610 for (
unsigned j = 0; j < n_dim - 1; j++)
612 interpolated_A(j,
i) +=
613 nodal_position_gen(l, this->bulk_position_type(k),
i) *
622 for (
unsigned i = 0;
i < n_dim - 1;
i++)
624 for (
unsigned j = 0; j < n_dim - 1; j++)
629 for (
unsigned k = 0; k < n_dim; k++)
631 A(
i, j) += interpolated_A(
i, k) * interpolated_A(j, k);
638 this->outer_unit_normal(ipt, interpolated_normal);
648 Adet = A(0, 0) * A(1, 1) - A(0, 1) * A(1, 0);
652 "Wrong dimension in RefineableSolidTractionElement",
653 OOMPH_CURRENT_FUNCTION,
654 OOMPH_EXCEPTION_LOCATION);
659 double W = w * sqrt(Adet);
664 ipt, interpolated_xi, interpolated_x, interpolated_normal, traction);
667 unsigned n_master = 1;
668 double hang_weight = 1.0;
673 for (
unsigned l = 0; l < n_node; l++)
676 Node* local_node_pt = node_pt(l);
679 bool is_hanging = local_node_pt->
is_hanging();
698 for (
unsigned m = 0; m < n_master; m++)
703 position_local_eqn_at_node = local_position_hang_eqn(
712 for (
unsigned k = 0; k < n_position_type; k++)
715 for (
unsigned i = 0;
i < n_dim;
i++)
717 position_local_eqn_at_node(k,
i) =
718 this->position_local_eqn(l, k,
i);
727 for (
unsigned k = 0; k < n_position_type; k++)
730 for (
unsigned i = 0;
i < n_dim;
i++)
732 local_eqn = position_local_eqn_at_node(k,
i);
754 residuals[local_eqn] -=
755 traction[
i] * psi(l, k) *
W * hang_weight;
777 template<
class ELEMENT,
unsigned DIM>
798 const bool& called_from_refineable_constructor =
false)
800 element_pt,
face_index, called_from_refineable_constructor),
803 unsigned n_lagr = DIM - 1;
804 unsigned n_dim = DIM;
832 OOMPH_CURRENT_FUNCTION,
833 OOMPH_EXCEPTION_LOCATION);
867 for (
unsigned i = 0;
i < DIM;
i++)
878 void output(std::ostream& outfile,
const unsigned& n_plot)
881 outfile <<
"ZONE" << std::endl;
888 unsigned el_dim = this->
dim();
895 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
899 for (
unsigned i = 0;
i < el_dim;
i++)
925 for (
unsigned i = 0;
i < DIM;
i++)
927 outfile << r[
i] <<
" ";
929 for (
unsigned i = 0;
i < DIM;
i++)
933 for (
unsigned i = 0;
i < DIM;
i++)
935 outfile << unit_normal[
i] <<
" ";
937 outfile << std::endl;
949 throw OomphLibError(
"It doesn't make sense to specify an external "
950 "traction in an FSI context",
951 OOMPH_CURRENT_FUNCTION,
952 OOMPH_EXCEPTION_LOCATION);
972 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const;
984 template<
class ELEMENT,
unsigned DIM>
986 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const
989 std::pair<unsigned, unsigned> dof_lookup;
992 const unsigned n_node = this->nnode();
995 const unsigned n_position_type = nnodal_position_type();
996 const unsigned nodal_dim = nodal_dimension();
999 int local_unknown = 0;
1002 for (
unsigned n = 0; n < n_node; n++)
1005 for (
unsigned k = 0; k < n_position_type; k++)
1008 for (
unsigned i = 0;
i < nodal_dim;
i++)
1011 local_unknown = position_local_eqn(n, k,
i);
1014 if (local_unknown >= 0)
1018 dof_lookup.first = this->eqn_number(local_unknown);
1019 dof_lookup.second = 0;
1022 dof_lookup_list.push_front(dof_lookup);
1041 template<
class ELEMENT,
unsigned DIM>
1076 residuals, jacobian);
1103 template<
class ELEMENT>
1117 const unsigned&
id = 0,
1118 const bool& called_from_refineable_constructor =
false)
1140 if (!called_from_refineable_constructor)
1142 if (element_pt->
dim() == 3)
1152 "This face element will not work correctly if nodes are "
1153 "hanging\nUse the refineable version instead. ",
1154 OOMPH_CURRENT_FUNCTION,
1155 OOMPH_EXCEPTION_LOCATION);
1168 "ImposeDisplacementByLagrangeMultiplierElement cannot (currently) "
1169 "be used with elements that have generalised positional dofs",
1170 OOMPH_CURRENT_FUNCTION,
1171 OOMPH_EXCEPTION_LOCATION);
1177 unsigned dim = element_pt->
dim();
1236 std::ostringstream error_message;
1237 error_message <<
"About to wipe external data for "
1238 <<
"ImposeDisplacementByLagrangeMultiplierElement.\n"
1239 <<
"I noted that N_assigned_geom_data = "
1243 <<
"so we're going to wipe some external data that\n"
1244 <<
"is not geometric Data of the GeomObject that\n"
1245 <<
"specifies the desired boundary shape.\n"
1248 OOMPH_CURRENT_FUNCTION,
1249 OOMPH_EXCEPTION_LOCATION);
1253 for (
unsigned i = 0;
i < n_geom_data;
i++)
1268 unsigned n_node =
nnode();
1274 unsigned dim_el =
dim();
1284 std::ostringstream error_message;
1285 error_message <<
"About to wipe external data for "
1286 <<
"ImposeDisplacementByLagrangeMultiplierElement.\n"
1287 <<
"I noted that N_assigned_geom_data = "
1291 <<
"so we're going to wipe some external data that\n"
1292 <<
"is not geometric Data of the GeomObject that\n"
1293 <<
"specifies the desired boundary shape.\n"
1296 OOMPH_CURRENT_FUNCTION,
1297 OOMPH_EXCEPTION_LOCATION);
1314 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
1324 for (
unsigned j = 0; j < n_node; j++)
1326 for (
unsigned k = 0; k < n_position_type; k++)
1329 for (
unsigned i = 0;
i < dim_el;
i++)
1343 for (
unsigned i = 0;
i < n_geom_data;
i++)
1369 residuals, jacobian, 1);
1391 const unsigned n_dof = this->
ndof();
1392 for (
unsigned i = 0;
i < n_dof;
i++)
1394 residuals[
i] *= -1.0;
1395 for (
unsigned j = 0; j < n_dof; j++)
1397 jacobian(
i, j) *= -1.0;
1404 void output(std::ostream& outfile,
const unsigned& n_plot)
1407 unsigned dim_el =
dim();
1413 if (n_position_type != 1)
1416 "ImposeDisplacementByLagrangeMultiplierElement cannot (currently) be "
1417 "used with elements that have generalised positional dofs",
1418 OOMPH_CURRENT_FUNCTION,
1419 OOMPH_EXCEPTION_LOCATION);
1428 unsigned n_node =
nnode();
1429 Shape psi(n_node, n_position_type);
1436 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
1448 for (
unsigned j = 0; j < n_node; j++)
1459 unsigned first_index =
1463 for (
unsigned i = 0;
i < dim_el + 1;
i++)
1468 lambda[
i] += nod_pt->
value(first_index +
i) * psi(j, 0);
1471 for (
unsigned i = 0;
i < dim_el;
i++)
1474 for (
unsigned k = 0; k < n_position_type; k++)
1486 for (
unsigned i = 0;
i < dim_el + 1;
i++)
1488 outfile << x[
i] <<
" ";
1490 for (
unsigned i = 0;
i < dim_el + 1;
i++)
1492 outfile << -lambda[
i] <<
" ";
1494 for (
unsigned i = 0;
i < dim_el + 1;
i++)
1496 outfile << r_prescribed[
i] <<
" ";
1502 outfile << std::endl;
1510 unsigned n_plot = 5;
1523 if (n_position_type != 1)
1526 "ImposeDisplacementByLagrangeMultiplierElement cannot (currently) be "
1527 "used with elements that have generalised positional dofs",
1528 OOMPH_CURRENT_FUNCTION,
1529 OOMPH_EXCEPTION_LOCATION);
1534 unsigned n_node =
nnode();
1537 unsigned dim_el =
dim();
1541 DShape dpsids(n_node, dim_el);
1548 double squared_error = 0.0;
1551 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
1566 for (
unsigned j = 0; j < n_node; j++)
1576 unsigned first_index =
1580 for (
unsigned i = 0;
i < dim_el + 1;
i++)
1583 lambda[
i] += nod_pt->
value(first_index +
i) * psi(j);
1584 for (
unsigned ii = 0; ii < dim_el; ii++)
1586 interpolated_a(ii,
i) +=
1592 for (
unsigned k = 0; k < n_position_type; k++)
1595 for (
unsigned i = 0;
i < dim_el;
i++)
1608 for (
unsigned i = 0;
i < dim_el;
i++)
1610 for (
unsigned j = 0; j < dim_el; j++)
1615 for (
unsigned k = 0; k < dim_el + 1; k++)
1617 a(
i, j) += interpolated_a(
i, k) * interpolated_a(j, k);
1632 adet = a(0, 0) * a(1, 1) - a(0, 1) * a(1, 0);
1638 "fill_in_generic_contribution_to_residuals_displ_lagr_multiplier",
1639 "ImposeDisplacementByLagrangeMultiplierElement::fill_in_generic_"
1640 "contribution_to_residuals_displ_lagr_multiplier()",
1641 OOMPH_EXCEPTION_LOCATION);
1657 double W = w * sqrt(adet);
1662 for (
unsigned i = 0;
i < dim_el + 1;
i++)
1665 (x[
i] - r_prescribed[
i]) * (x[
i] - r_prescribed[
i]) *
W;
1669 return squared_error;
1679 const unsigned& flag)
1685 if (n_position_type != 1)
1688 "ImposeDisplacementByLagrangeMultiplierElement cannot (currently) be "
1689 "used with elements that have generalised positional dofs",
1690 OOMPH_CURRENT_FUNCTION,
1691 OOMPH_EXCEPTION_LOCATION);
1696 unsigned n_node =
nnode();
1699 unsigned dim_el =
dim();
1703 DShape dpsids(n_node, dim_el);
1709 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
1724 for (
unsigned j = 0; j < n_node; j++)
1734 unsigned first_index =
1738 for (
unsigned i = 0;
i < dim_el + 1;
i++)
1741 lambda[
i] += nod_pt->
value(first_index +
i) * psi(j);
1742 for (
unsigned ii = 0; ii < dim_el; ii++)
1744 interpolated_a(ii,
i) +=
1750 for (
unsigned k = 0; k < n_position_type; k++)
1753 for (
unsigned i = 0;
i < dim_el;
i++)
1766 for (
unsigned i = 0;
i < dim_el;
i++)
1768 for (
unsigned j = 0; j < dim_el; j++)
1773 for (
unsigned k = 0; k < dim_el + 1; k++)
1775 a(
i, j) += interpolated_a(
i, k) * interpolated_a(j, k);
1790 adet = a(0, 0) * a(1, 1) - a(0, 1) * a(1, 0);
1796 "fill_in_generic_contribution_to_residuals_displ_lagr_multiplier",
1797 "ImposeDisplacementByLagrangeMultiplierElement::fill_in_generic_"
1798 "contribution_to_residuals_displ_lagr_multiplier()",
1799 OOMPH_EXCEPTION_LOCATION);
1815 double W = w * sqrt(adet);
1820 for (
unsigned i = 0;
i < dim_el + 1;
i++)
1823 for (
unsigned j = 0; j < n_node; j++)
1839 residuals[local_eqn] += (x[
i] - r_prescribed[
i]) * psi(j) *
W;
1846 for (
unsigned jj = 0; jj < n_node; jj++)
1849 if (local_unknown >= 0)
1851 jacobian(local_eqn, local_unknown) += psi(jj) * psi(j) *
W;
1865 residuals[local_eqn] += lambda[
i] * psi(j) *
W;
1872 for (
unsigned jj = 0; jj < n_node; jj++)
1883 if (local_unknown >= 0)
1885 jacobian(local_eqn, local_unknown) += psi(jj) * psi(j) *
W;
1903 return this->
dim() + 1;
1914 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const
1917 std::pair<unsigned, unsigned> dof_lookup;
1920 const unsigned n_node = this->
nnode();
1923 unsigned dim_el = this->
dim();
1924 for (
unsigned i = 0;
i < dim_el + 1;
i++)
1927 for (
unsigned j = 0; j < n_node; j++)
1940 dof_lookup.first = this->
eqn_number(local_eqn);
1941 dof_lookup.second =
i;
1944 dof_lookup_list.push_front(dof_lookup);
2008 template<
class ELEMENT>
2023 const unsigned&
id = 0)
2052 residuals, jacobian, 1);
2066 const unsigned& flag)
2072 if (n_position_type != 1)
2075 "RefineableImposeDisplacementByLagrangeMultiplierElement \n cannot "
2076 "(currently) be used with elements that have generalised\n "
2077 "positional dofs. Upgrade should be straightforward though the code "
2078 "is in a \n bit of state, with generalised degrees of freedom "
2079 "sometimes half taken into \n account, sometimes completely "
2081 OOMPH_CURRENT_FUNCTION,
2082 OOMPH_EXCEPTION_LOCATION);
2087 unsigned n_node = this->
nnode();
2090 unsigned dim_el = this->
dim();
2094 DShape dpsids(n_node, dim_el);
2102 int local_unknown = 0;
2105 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
2122 for (
unsigned j = 0; j < n_node; j++)
2131 unsigned first_index =
2135 for (
unsigned i = 0;
i < dim_el + 1;
i++)
2138 lambda[
i] += nod_pt->
value(first_index +
i) * psi(j);
2139 for (
unsigned ii = 0; ii < dim_el; ii++)
2141 interpolated_a(ii,
i) +=
2147 for (
unsigned k = 0; k < n_position_type; k++)
2150 for (
unsigned i = 0;
i < dim_el;
i++)
2163 for (
unsigned i = 0;
i < dim_el;
i++)
2165 for (
unsigned j = 0; j < dim_el; j++)
2170 for (
unsigned k = 0; k < dim_el + 1; k++)
2172 a(
i, j) += interpolated_a(
i, k) * interpolated_a(j, k);
2187 adet = a(0, 0) * a(1, 1) - a(0, 1) * a(1, 0);
2193 "refineable_fill_in_generic_contribution_to_residuals_displ_lagr_"
2195 "RefineableImposeDisplacementByLagrangeMultiplierElement::"
2196 "refineablefill_in_generic_contribution_to_residuals_displ_lagr_"
2198 OOMPH_EXCEPTION_LOCATION);
2214 double W = w * sqrt(adet);
2221 unsigned n_master = 1;
2222 unsigned n_master2 = 1;
2223 double hang_weight = 1.0;
2224 double hang_weight2 = 1.0;
2232 for (
unsigned j = 0; j < n_node; j++)
2238 bool is_node_hanging = local_node_pt->
is_hanging();
2241 if (is_node_hanging)
2246 n_master = hang_info_pt->
nmaster();
2255 for (
unsigned m = 0; m < n_master; m++)
2258 for (
unsigned i = 0;
i < dim_el + 1;
i++)
2263 if (is_node_hanging)
2300 residuals[local_eqn] +=
2301 (x[
i] - r_prescribed[
i]) * psi(j) *
W * hang_weight;
2308 for (
unsigned jj = 0; jj < n_node; jj++)
2314 bool is_node_hanging2 = local_node2_pt->
is_hanging();
2317 if (is_node_hanging2)
2319 hang_info2_pt = local_node2_pt->
hanging_pt();
2322 n_master2 = hang_info2_pt->
nmaster();
2331 for (
unsigned m2 = 0; m2 < n_master2; m2++)
2336 n_position_type, dim_el + 1);
2338 if (is_node_hanging2)
2362 position_local_eqn_at_node2(k2, i2) =
2371 local_unknown = position_local_eqn_at_node2(k2,
i);
2372 if (local_unknown >= 0)
2374 jacobian(local_eqn, local_unknown) +=
2375 psi(jj) * psi(j) *
W * hang_weight * hang_weight2;
2390 if (is_node_hanging)
2409 position_local_eqn_at_node(k2, i2) =
2414 local_eqn = position_local_eqn_at_node(k,
i);
2420 residuals[local_eqn] += lambda[
i] * psi(j) *
W * hang_weight;
2427 for (
unsigned jj = 0; jj < n_node; jj++)
2433 bool is_node_hanging2 = local_node2_pt->
is_hanging();
2436 if (is_node_hanging2)
2438 hang_info2_pt = local_node2_pt->
hanging_pt();
2441 n_master2 = hang_info2_pt->
nmaster();
2450 for (
unsigned m2 = 0; m2 < n_master2; m2++)
2455 if (is_node_hanging2)
2484 ->index_of_first_value_assigned_by_face_element(
2492 if (local_unknown >= 0)
2494 jacobian(local_eqn, local_unknown) +=
2495 psi(jj) * psi(j) *
W * hang_weight * hang_weight2;
2533 template<
class ELEMENT>
2563 const unsigned&
id = 0,
2564 const bool& called_from_refineable_constructor =
false)
2580 if (!called_from_refineable_constructor)
2582 if (element_pt->
dim() == 3)
2592 "This face element will not work correctly if nodes are "
2593 "hanging\nUse the refineable version instead. ",
2594 OOMPH_CURRENT_FUNCTION,
2595 OOMPH_EXCEPTION_LOCATION);
2607 throw OomphLibError(
"FSIImposeDisplacementByLagrangeMultiplierElement"
2608 " cannot (currently) be used with elements that "
2609 "have generalised positional dofs",
2610 OOMPH_CURRENT_FUNCTION,
2611 OOMPH_EXCEPTION_LOCATION);
2617 unsigned dim = element_pt->
dim();
2645 residuals, jacobian, 1);
2665 const unsigned n_dof = this->
ndof();
2666 for (
unsigned i = 0;
i < n_dof;
i++)
2668 residuals[
i] *= -1.0;
2669 for (
unsigned j = 0; j < n_dof; j++)
2671 jacobian(
i, j) *= -1.0;
2678 void output(std::ostream& outfile,
const unsigned& n_plot)
2681 unsigned dim_el =
dim();
2687 if (n_position_type != 1)
2690 "FSIImposeDisplacementByLagrangeMultiplierElement cannot (currently) "
2691 "be used with elements that have generalised positional dofs",
2692 OOMPH_CURRENT_FUNCTION,
2693 OOMPH_EXCEPTION_LOCATION);
2702 unsigned n_node =
nnode();
2703 Shape psi(n_node, n_position_type);
2710 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
2722 for (
unsigned j = 0; j < n_node; j++)
2733 unsigned first_index =
2737 for (
unsigned i = 0;
i < dim_el + 1;
i++)
2742 lambda[
i] += nod_pt->
value(first_index +
i) * psi(j, 0);
2745 for (
unsigned i = 0;
i < dim_el;
i++)
2748 for (
unsigned k = 0; k < n_position_type; k++)
2756 for (
unsigned i = 0;
i < dim_el + 1;
i++)
2758 outfile << x[
i] <<
" ";
2760 for (
unsigned i = 0;
i < dim_el + 1;
i++)
2762 outfile << -lambda[
i] <<
" ";
2764 outfile << std::endl;
2772 unsigned n_plot = 5;
2783 const unsigned& flag)
2789 if (n_position_type != 1)
2792 "FSIImposeDisplacementByLagrangeMultiplierElement cannot (currently) "
2793 "be used with elements that have generalised positional dofs",
2794 OOMPH_CURRENT_FUNCTION,
2795 OOMPH_EXCEPTION_LOCATION);
2800 unsigned n_node =
nnode();
2803 unsigned dim_el =
dim();
2807 DShape dpsids(n_node, dim_el);
2813 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
2828 for (
unsigned j = 0; j < n_node; j++)
2838 unsigned first_index =
2842 for (
unsigned i = 0;
i < dim_el + 1;
i++)
2845 lambda[
i] += nod_pt->
value(first_index +
i) * psi(j);
2846 for (
unsigned ii = 0; ii < dim_el; ii++)
2848 interpolated_a(ii,
i) +=
2856 for (
unsigned i = 0;
i < dim_el;
i++)
2858 for (
unsigned j = 0; j < dim_el; j++)
2863 for (
unsigned k = 0; k < dim_el + 1; k++)
2865 a(
i, j) += interpolated_a(
i, k) * interpolated_a(j, k);
2880 adet = a(0, 0) * a(1, 1) - a(0, 1) * a(1, 0);
2886 "fill_in_generic_contribution_to_residuals_displ_lagr_multiplier",
2887 "FSIImposeDisplacementByLagrangeMultiplierElement::fill_in_"
2888 "generic_contribution_to_residuals_fsi_displ_lagr_multiplier()",
2889 OOMPH_EXCEPTION_LOCATION);
2905 double W = w * sqrt(adet);
2910 for (
unsigned i = 0;
i < dim_el + 1;
i++)
2913 for (
unsigned j = 0; j < n_node; j++)
2929 residuals[local_eqn] += (x[
i] - r_prescribed[
i]) * psi(j) *
W;
2936 for (
unsigned jj = 0; jj < n_node; jj++)
2939 if (local_unknown >= 0)
2941 jacobian(local_eqn, local_unknown) += psi(jj) * psi(j) *
W;
2955 residuals[local_eqn] += lambda[
i] * psi(j) *
W;
2962 for (
unsigned jj = 0; jj < n_node; jj++)
2973 if (local_unknown >= 0)
2975 jacobian(local_eqn, local_unknown) += psi(jj) * psi(j) *
W;
2992 return this->
dim() + 1;
3003 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const
3006 std::pair<unsigned, unsigned> dof_lookup;
3009 const unsigned n_node = this->
nnode();
3012 unsigned dim_el = this->
dim();
3013 for (
unsigned i = 0;
i < dim_el + 1;
i++)
3016 for (
unsigned j = 0; j < n_node; j++)
3029 dof_lookup.first = this->
eqn_number(local_eqn);
3030 dof_lookup.second =
i;
3033 dof_lookup_list.push_front(dof_lookup);
3068 template<
class ELEMENT>
3094 const unsigned&
id = 0)
3123 residuals, jacobian, 1);
3136 const unsigned& flag)
3142 if (n_position_type != 1)
3145 "RefineableImposeDisplacementByLagrangeMultiplierElement \n cannot "
3146 "(currently) be used with elements that have generalised\n "
3147 "positional dofs. Upgrade should be straightforward though the code "
3148 "is in a \n bit of state, with generalised degrees of freedom "
3149 "sometimes half taken into \n account, sometimes completely "
3151 OOMPH_CURRENT_FUNCTION,
3152 OOMPH_EXCEPTION_LOCATION);
3157 unsigned n_node = this->
nnode();
3160 unsigned dim_el = this->
dim();
3164 DShape dpsids(n_node, dim_el);
3172 int local_unknown = 0;
3175 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
3192 for (
unsigned j = 0; j < n_node; j++)
3201 unsigned first_index =
3205 for (
unsigned i = 0;
i < dim_el + 1;
i++)
3208 lambda[
i] += nod_pt->
value(first_index +
i) * psi(j);
3209 for (
unsigned ii = 0; ii < dim_el; ii++)
3211 interpolated_a(ii,
i) +=
3220 for (
unsigned i = 0;
i < dim_el;
i++)
3222 for (
unsigned j = 0; j < dim_el; j++)
3227 for (
unsigned k = 0; k < dim_el + 1; k++)
3229 a(
i, j) += interpolated_a(
i, k) * interpolated_a(j, k);
3244 adet = a(0, 0) * a(1, 1) - a(0, 1) * a(1, 0);
3250 "refineable_fill_in_generic_contribution_to_residuals_displ_lagr_"
3252 "RefineableImposeDisplacementByLagrangeMultiplierElement::"
3253 "refineablefill_in_generic_contribution_to_residuals_displ_lagr_"
3255 OOMPH_EXCEPTION_LOCATION);
3272 double W = w * sqrt(adet);
3279 unsigned n_master = 1;
3280 unsigned n_master2 = 1;
3281 double hang_weight = 1.0;
3282 double hang_weight2 = 1.0;
3290 for (
unsigned j = 0; j < n_node; j++)
3296 bool is_node_hanging = local_node_pt->
is_hanging();
3299 if (is_node_hanging)
3304 n_master = hang_info_pt->
nmaster();
3313 for (
unsigned m = 0; m < n_master; m++)
3316 for (
unsigned i = 0;
i < dim_el + 1;
i++)
3321 if (is_node_hanging)
3358 residuals[local_eqn] +=
3359 (x[
i] - r_prescribed[
i]) * psi(j) *
W * hang_weight;
3366 for (
unsigned jj = 0; jj < n_node; jj++)
3372 bool is_node_hanging2 = local_node2_pt->
is_hanging();
3375 if (is_node_hanging2)
3377 hang_info2_pt = local_node2_pt->
hanging_pt();
3380 n_master2 = hang_info2_pt->
nmaster();
3389 for (
unsigned m2 = 0; m2 < n_master2; m2++)
3394 n_position_type, dim_el + 1);
3396 if (is_node_hanging2)
3420 position_local_eqn_at_node2(k2, i2) =
3429 local_unknown = position_local_eqn_at_node2(k2,
i);
3430 if (local_unknown >= 0)
3432 jacobian(local_eqn, local_unknown) +=
3433 psi(jj) * psi(j) *
W * hang_weight * hang_weight2;
3448 if (is_node_hanging)
3467 position_local_eqn_at_node(k2, i2) =
3472 local_eqn = position_local_eqn_at_node(k,
i);
3478 residuals[local_eqn] += lambda[
i] * psi(j) *
W * hang_weight;
3485 for (
unsigned jj = 0; jj < n_node; jj++)
3491 bool is_node_hanging2 = local_node2_pt->
is_hanging();
3494 if (is_node_hanging2)
3496 hang_info2_pt = local_node2_pt->
hanging_pt();
3499 n_master2 = hang_info2_pt->
nmaster();
3508 for (
unsigned m2 = 0; m2 < n_master2; m2++)
3513 if (is_node_hanging2)
3542 ->index_of_first_value_assigned_by_face_element(
3550 if (local_unknown >= 0)
3552 jacobian(local_eqn, local_unknown) +=
3553 psi(jj) * psi(j) *
W * hang_weight * hang_weight2;
A class that contains the information required by Nodes that are located on Mesh boundaries....
unsigned index_of_first_value_assigned_by_face_element(const unsigned &face_id=0) const
Return the index of the first value associated with the i-th face element value. If no argument is sp...
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
This is a base class for all elements that require external sources (e.g. FSI, multi-domain problems ...
Vector< double > & external_element_local_coord(const unsigned &interaction_index, const unsigned &ipt)
Access function to get source element's local coords for specified interaction index at specified int...
void set_ninteraction(const unsigned &n_interaction)
Set the number of interactions in the element This function is usually called in the specific element...
void describe_local_dofs(std::ostream &out, const std::string &curr_string) const
Function to describe the local dofs of the element. The ostream specifies the output stream to which ...
FiniteElement *& external_element_pt(const unsigned &interaction_index, const unsigned &ipt)
Access function to source element for specified interaction index at specified integration point.
void fill_in_jacobian_from_external_interaction_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the contributions to the jacobian from all external interaction degrees of freedom (geometr...
///////////////////////////////////////////////////////////////////////// ///////////////////////////...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
FSIImposeDisplacementByLagrangeMultiplierElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0, const bool &called_from_refineable_constructor=false)
Constructor takes a "bulk" element and the index that identifies which face the FaceElement is suppos...
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Fill in contribution to Mass matrix and Jacobian. There is no contributiont to mass matrix so simply ...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals.
void fill_in_generic_contribution_to_residuals_fsi_displ_lagr_multiplier(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Helper function to compute the residuals and, if flag==1, the Jacobian.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: Just the soli...
void describe_local_dofs(std::ostream &out, const std::string ¤t_string) const
Function to describe the local dofs of the elements. The ostream specifies the output stream to which...
void output(std::ostream &outfile)
Output function.
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
////////////////////////////////////////////////////////////////////// //////////////////////////////...
void set_normal_pointing_into_fluid()
Set the normal computed by underlying FaceElement to point into the fluid.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: Just the soli...
void dposition_dlagrangian_at_local_coordinate(const Vector< double > &s, DenseMatrix< double > &drdxi) const
Derivative of position vector w.r.t. the SolidFiniteElement's Lagrangian coordinates; evaluated at cu...
FSISolidTractionElement(FiniteElement *const &element_pt, const int &face_index, const bool &called_from_refineable_constructor=false)
Constructor: Create element as FSIWallElement (and thus, by inheritance, a GeomObject) with DIM-1 Lag...
~FSISolidTractionElement()
Destructor: empty.
virtual void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Final overload... Forwards to the version in the FSIWallElement.
virtual void get_traction(const unsigned &intpt, const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Get the load vector: Pass number of the integration point, Lagr. coordinate, Eulerian coordinate (nei...
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: Note we can only output the traction at Gauss points so n_plot is actually ignored.
virtual void(*&)(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction) traction_fct_pt()
Broken overloaded reference to the traction function pointer. It doesn't make sense to specify an ext...
void set_normal_pointing_out_of_fluid()
Set the normal computed by underlying FaceElement to point out of the fluid.
bool Normal_points_into_fluid
Boolean flag to indicate whether the normal is directed into the fluid.
////////////////////////////////////////////////////////////////////////// //////////////////////////...
void setup_fsi_wall_element(const unsigned &nlagr_solid, const unsigned &ndim_fluid)
Setup: Assign storage – pass the Eulerian dimension of the "adjacent" fluid elements and the number o...
void fluid_load_vector(const unsigned &intpt, const Vector< double > &N, Vector< double > &load)
Get FE Jacobian by systematic finite differencing w.r.t. nodal positition Data, internal and external...
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
bool Boundary_number_in_bulk_mesh_has_been_set
Has the Boundary_number_in_bulk_mesh been set? Only included if compiled with PARANOID switched on.
unsigned Boundary_number_in_bulk_mesh
The boundary number in the bulk mesh to which this element is attached.
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
void add_additional_values(const Vector< unsigned > &nadditional_values, const unsigned &id)
Helper function adding additional values for the unknowns associated with the FaceElement....
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
const unsigned & boundary_number_in_bulk_mesh() const
Broken assignment operator.
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s. Overloaded to get information from bulk...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
A general Finite Element class.
virtual void dshape_local_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsids) const
Return the geometric shape function and its derivative w.r.t. the local coordinates at the ipt-th int...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
unsigned nnodal_position_type() const
Return the number of coordinate types that the element requires to interpolate the geometry between t...
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
unsigned nnode() const
Return the number of nodes.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
double nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. If the node is hanging, the appropriate interpolation is ...
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Function for building a lower dimensional FaceElement on the specified face of the FiniteElement....
virtual void describe_nodal_local_dofs(std::ostream &out, const std::string ¤t_string) const
Function to describe the local dofs of the element[s]. The ostream specifies the output stream to whi...
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Add tecplot zone "footer" to output stream (when plotting nplot points in each "coordinate direction"...
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Return the geometric shape function at the ipt-th integration point.
bool has_hanging_nodes() const
Return boolean to indicate if any of the element's nodes are geometrically hanging.
void fill_in_jacobian_from_external_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
Calculate the contributions to the jacobian from the external degrees of freedom using finite differe...
unsigned nexternal_data() const
Return the number of external data objects.
bool is_halo() const
Is this element a halo?
unsigned ndof() const
Return the number of equations/dofs in the element.
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
void flush_external_data()
Flush all external data.
unsigned add_external_data(Data *const &data_pt, const bool &fd=true)
Add a (pointer to an) external data object to the element and return its index (i....
/////////////////////////////////////////////////////////////////////
virtual void position(const Vector< double > &zeta, Vector< double > &r) const =0
Parametrised position on object at current time: r(zeta).
virtual unsigned ngeom_data() const
How many items of Data does the shape of the object depend on? This is implemented as a broken virtua...
virtual void locate_zeta(const Vector< double > &zeta, GeomObject *&sub_geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
A geometric object may be composed of may sub-objects (e.g. a finite-element representation of a boun...
virtual Data * geom_data_pt(const unsigned &j)
Return pointer to the j-th Data item that the object's shape depends on. This is implemented as a bro...
Class that contains data for hanging nodes.
double const & master_weight(const unsigned &i) const
Return weight for dofs on i-th master node.
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
unsigned nmaster() const
Return the number of master nodes.
///////////////////////////////////////////////////////////////////////// ///////////////////////////...
Vector< Vector< double > > Zeta_sub_geom_object
Storage for local coordinates in sub-GeomObjects at integration points.
unsigned N_assigned_geom_data
Bool to record the number of geom Data that has been assigned to external data (we're keeping a recor...
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals.
Vector< GeomObject * > Sub_geom_object_pt
Storage for sub-GeomObject at the integration points.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: We only label...
void output(std::ostream &outfile)
Output function.
void set_boundary_shape_geom_object_pt(GeomObject *boundary_shape_geom_object_pt, const unsigned &boundary_number_in_bulk_mesh)
Set GeomObject that specifies the prescribed boundary displacement; GeomObject is assumed to be param...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
double square_of_l2_norm_of_error()
Compute square of L2 norm of error between prescribed and actual displacement.
GeomObject * Boundary_shape_geom_object_pt
GeomObject that specifies the prescribed boundary displacement; GeomObject is assumed to be parametri...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
GeomObject * boundary_shape_geom_object_pt() const
Access to GeomObject that specifies the prescribed boundary displacement; GeomObject is assumed to be...
bool Sparsify
Boolean flag to indicate if we're using geometric Data of sub-GeomObjects that make up the (possibly ...
ImposeDisplacementByLagrangeMultiplierElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0, const bool &called_from_refineable_constructor=false)
Constructor takes a "bulk" element and the index that identifies which face the FaceElement is suppos...
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Fill in contribution to Mass matrix and Jacobian. There is no contributiont to mass matrix so simply ...
void fill_in_generic_contribution_to_residuals_displ_lagr_multiplier(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Helper function to compute the residuals and, if flag==1, the Jacobian.
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
HangInfo *const & hanging_pt() const
Return pointer to hanging node data (this refers to the geometric hanging node status) (const version...
bool is_hanging() const
Test whether the node is geometrically hanging.
double value(const unsigned &i) const
Return i-th value (dofs or pinned) at this node either directly or via hanging node representation....
/////////////////////////////////////////////////////////////////// /////////////////////////////////...
An OomphLibError object which should be thrown when an run-time error is encountered....
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
int local_hang_eqn(Node *const &node_pt, const unsigned &i)
Access function that returns the local equation number for the hanging node variables (values stored ...
///////////////////////////////////////////////////////////////////////// ///////////////////////////...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
RefineableFSIImposeDisplacementByLagrangeMultiplierElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Constructor takes a "bulk" element and the index that identifies which face the FaceElement is suppos...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals.
void refineable_fill_in_generic_contribution_to_residuals_fsi_displ_lagr_multiplier(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Helper function to compute the residuals and, if flag==1, the Jacobian.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: Same for all nodes since it's a refineable element.
////////////////////////////////////////////////////////////////////// //////////////////////////////...
virtual void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Final overload. Get contributions from refineable solid traction element and derivatives from externa...
~RefineableFSISolidTractionElement()
Destructor: empty.
RefineableFSISolidTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor: Create element as FSIWallElement (and thus, by inheritance, a GeomObject) with DIM-1 Lag...
///////////////////////////////////////////////////////////////////////// ///////////////////////////...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals.
void refineable_fill_in_generic_contribution_to_residuals_displ_lagr_multiplier(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Helper function to compute the residuals and, if flag==1, the Jacobian.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: Same for all nodes since it's a refineable element.
RefineableImposeDisplacementByLagrangeMultiplierElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Constructor takes a "bulk" element and the index that identifies which face the FaceElement is suppos...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
DenseMatrix< int > & local_position_hang_eqn(Node *const &node_pt)
Access the local equation number of of hanging node variables associated with nodal positions....
////////////////////////////////////////////////////////////////////// //////////////////////////////...
~RefineableSolidTractionElement()
Destructor should not delete anything.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
RefineableSolidTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the face index.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values are the same as those in the bulk element.
void refineable_fill_in_contribution_to_residuals_solid_traction(Vector< double > &residuals)
This function returns the residuals for the traction function.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the Jacobian.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
SolidFaceElements combine FaceElements and SolidFiniteElements and overload various functions so they...
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
The "global" intrinsic coordinate of the element when viewed as part of a geometric object should be ...
double interpolated_xi(const Vector< double > &s, const unsigned &i) const
Return i-th FE-interpolated Lagrangian coordinate xi[i] at local coordinate s. Overloaded from SolidF...
double lagrangian_position(const unsigned &n, const unsigned &i) const
Return i-th Lagrangian coordinate at local node n.
Data * geom_data_pt(const unsigned &j)
Return pointer to the j-th Data item that the object's shape depends on. (Redirects to the nodes' pos...
int position_local_eqn(const unsigned &n, const unsigned &k, const unsigned &j) const
Access function that returns the local equation number that corresponds to the j-th coordinate of the...
A class for elements that allow the imposition of an applied traction in the principle of virtual dis...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
void fill_in_contribution_to_residuals_solid_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
SolidTractionElement(FiniteElement *const &element_pt, const int &face_index, const bool &called_from_refineable_constructor=false)
Constructor, which takes a "bulk" element and the value of the index and its limit.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
void(*&)(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction) traction_fct_pt()
Reference to the traction function pointer.
void output(std::ostream &outfile)
Output function.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
void traction(const Vector< double > &s, Vector< double > &traction)
Compute traction vector at specified local coordinate Should only be used for post-processing; ignore...
virtual void get_traction(const unsigned &intpt, const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Get the traction vector: Pass number of integration point (dummy), Lagr. coordinate and normal vector...
void output(FILE *file_pt)
C_style output function.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
void(* Traction_fct_pt)(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Pointer to an imposed traction function. Arguments: Lagrangian coordinate; Eulerian coordinate; outer...
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
void Zero_traction_fct(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Default load function (zero traction)
//////////////////////////////////////////////////////////////////// ////////////////////////////////...