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. ",
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()",
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",
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",
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;
867 for (
unsigned i = 0;
i <
DIM;
i++)
878 void output(std::ostream& outfile,
const unsigned& n_plot)
881 outfile <<
"ZONE" << std::endl;
895 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
918 if (!this->is_halo())
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",
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)
1137 N_assigned_geom_data = 0;
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. ",
1168 "ImposeDisplacementByLagrangeMultiplierElement cannot (currently) "
1169 "be used with elements that have generalised positional dofs",
1177 unsigned dim = element_pt->
dim();
1215 Boundary_number_in_bulk_mesh_has_been_set =
true;
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 = "
1240 << 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"
1253 for (
unsigned i = 0;
i < n_geom_data;
i++)
1258 N_assigned_geom_data = n_geom_data;
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 = "
1288 << 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"
1305 N_assigned_geom_data = 0;
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++)
1348 N_assigned_geom_data += n_geom_data;
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",
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++)
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++)
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",
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++)
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()",
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",
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++)
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()",
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);
1961 unsigned N_assigned_geom_data;
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 "
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++)
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_"
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. ",
2607 throw OomphLibError(
"FSIImposeDisplacementByLagrangeMultiplierElement"
2608 " cannot (currently) be used with elements that "
2609 "have generalised positional dofs",
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",
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++)
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++)
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",
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++)
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()",
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 "
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++)
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_"
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;
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
cout<< "The eigenvalues of A are:"<< endl<< ces.eigenvalues()<< endl;cout<< "The matrix of eigenvectors, V, is:"<< endl<< ces.eigenvectors()<< endl<< endl;complex< float > lambda
Definition: ComplexEigenSolver_compute.cpp:9
RowVector3d w
Definition: Matrix_resize_int.cpp:3
void load(Archive &ar, ParticleHandler &handl)
Definition: Particles.h:21
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
unsigned index_of_first_value_assigned_by_face_element(const unsigned &face_id=0) const
Definition: nodes.h:2061
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
Definition: nodes.h:483
Definition: element_with_external_element.h:56
Vector< double > & external_element_local_coord(const unsigned &interaction_index, const unsigned &ipt)
Definition: element_with_external_element.h:136
void set_ninteraction(const unsigned &n_interaction)
Definition: element_with_external_element.h:178
void describe_local_dofs(std::ostream &out, const std::string &curr_string) const
Definition: element_with_external_element.cc:205
FiniteElement *& external_element_pt(const unsigned &interaction_index, const unsigned &ipt)
Definition: element_with_external_element.h:107
void fill_in_jacobian_from_external_interaction_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: element_with_external_element.h:345
Definition: solid_traction_elements.h:2538
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
Definition: solid_traction_elements.h:2678
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
Definition: solid_traction_elements.h:2640
FSIImposeDisplacementByLagrangeMultiplierElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0, const bool &called_from_refineable_constructor=false)
Definition: solid_traction_elements.h:2560
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Definition: solid_traction_elements.h:2656
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals.
Definition: solid_traction_elements.h:2631
void fill_in_generic_contribution_to_residuals_fsi_displ_lagr_multiplier(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: solid_traction_elements.h:2780
unsigned ndof_types() const
Definition: solid_traction_elements.h:2990
void describe_local_dofs(std::ostream &out, const std::string ¤t_string) const
Definition: solid_traction_elements.h:2548
unsigned Id
Lagrange Id.
Definition: solid_traction_elements.h:3040
void output(std::ostream &outfile)
Output function.
Definition: solid_traction_elements.h:2770
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Definition: solid_traction_elements.h:3002
Definition: solid_traction_elements.h:780
void set_normal_pointing_into_fluid()
Definition: solid_traction_elements.h:814
unsigned ndof_types() const
Definition: solid_traction_elements.h:960
void dposition_dlagrangian_at_local_coordinate(const Vector< double > &s, DenseMatrix< double > &drdxi) const
Definition: solid_traction_elements.h:828
FSISolidTractionElement(FiniteElement *const &element_pt, const int &face_index, const bool &called_from_refineable_constructor=false)
Definition: solid_traction_elements.h:795
~FSISolidTractionElement()
Destructor: empty.
Definition: solid_traction_elements.h:809
virtual void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Final overload... Forwards to the version in the FSIWallElement.
Definition: solid_traction_elements.h:838
virtual void get_traction(const unsigned &intpt, const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Definition: solid_traction_elements.h:854
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Definition: solid_traction_elements.h:985
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: solid_traction_elements.h:878
virtual void(*&)(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction) traction_fct_pt()
Definition: solid_traction_elements.h:944
void set_normal_pointing_out_of_fluid()
Definition: solid_traction_elements.h:821
bool Normal_points_into_fluid
Boolean flag to indicate whether the normal is directed into the fluid.
Definition: solid_traction_elements.h:783
void setup_fsi_wall_element(const unsigned &nlagr_solid, const unsigned &ndim_fluid)
Definition: fsi.cc:121
void fluid_load_vector(const unsigned &intpt, const Vector< double > &N, Vector< double > &load)
Definition: fsi.cc:162
Definition: elements.h:4338
unsigned Boundary_number_in_bulk_mesh
The boundary number in the bulk mesh to which this element is attached.
Definition: elements.h:4388
int & face_index()
Definition: elements.h:4626
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
Definition: elements.cc:6006
void add_additional_values(const Vector< unsigned > &nadditional_values, const unsigned &id)
Definition: elements.h:4428
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
Definition: elements.h:4735
const unsigned & boundary_number_in_bulk_mesh() const
Broken assignment operator.
Definition: elements.h:4475
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Definition: elements.h:4528
Definition: elements.h:4998
Definition: elements.h:1313
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
virtual void output(std::ostream &outfile)
Definition: elements.h:3050
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
unsigned nnodal_position_type() const
Definition: elements.h:2463
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Definition: elements.cc:5132
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
virtual void dshape_local_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsids) const
Definition: elements.cc:3239
virtual void shape(const Vector< double > &s, Shape &psi) const =0
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned dim() const
Definition: elements.h:2611
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186
double nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.h:2317
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
Definition: elements.h:2484
virtual void describe_nodal_local_dofs(std::ostream &out, const std::string ¤t_string) const
Definition: elements.cc:1727
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Definition: elements.h:3174
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Definition: elements.cc:3220
bool has_hanging_nodes() const
Definition: elements.h:2470
unsigned nexternal_data() const
Return the number of external data objects.
Definition: elements.h:829
void flush_external_data()
Flush all external data.
Definition: elements.cc:387
unsigned ndof() const
Return the number of equations/dofs in the element.
Definition: elements.h:835
unsigned long eqn_number(const unsigned &ieqn_local) const
Definition: elements.h:704
unsigned add_external_data(Data *const &data_pt, const bool &fd=true)
Definition: elements.cc:307
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
void fill_in_jacobian_from_external_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
Definition: elements.cc:1199
Definition: geom_objects.h:101
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
Definition: geom_objects.h:209
virtual void locate_zeta(const Vector< double > &zeta, GeomObject *&sub_geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
Definition: geom_objects.h:381
virtual Data * geom_data_pt(const unsigned &j)
Definition: geom_objects.h:233
double const & master_weight(const unsigned &i) const
Return weight for dofs on i-th master node.
Definition: nodes.h:808
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
Definition: nodes.h:791
unsigned nmaster() const
Return the number of master nodes.
Definition: nodes.h:785
Definition: solid_traction_elements.h:1107
Vector< Vector< double > > Zeta_sub_geom_object
Definition: solid_traction_elements.h:1977
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Definition: solid_traction_elements.h:1913
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals.
Definition: solid_traction_elements.h:1355
Vector< GeomObject * > Sub_geom_object_pt
Storage for sub-GeomObject at the integration points.
Definition: solid_traction_elements.h:1973
unsigned ndof_types() const
Definition: solid_traction_elements.h:1901
unsigned Id
Lagrange Id.
Definition: solid_traction_elements.h:1952
void output(std::ostream &outfile)
Output function.
Definition: solid_traction_elements.h:1508
void set_boundary_shape_geom_object_pt(GeomObject *boundary_shape_geom_object_pt, const unsigned &boundary_number_in_bulk_mesh)
Definition: solid_traction_elements.h:1209
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
Definition: solid_traction_elements.h:1364
double square_of_l2_norm_of_error()
Definition: solid_traction_elements.h:1517
GeomObject * Boundary_shape_geom_object_pt
Definition: solid_traction_elements.h:1970
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
Definition: solid_traction_elements.h:1404
GeomObject * boundary_shape_geom_object_pt() const
Definition: solid_traction_elements.h:1195
bool Sparsify
Definition: solid_traction_elements.h:1982
ImposeDisplacementByLagrangeMultiplierElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0, const bool &called_from_refineable_constructor=false)
Definition: solid_traction_elements.h:1114
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Definition: solid_traction_elements.h:1382
void fill_in_generic_contribution_to_residuals_displ_lagr_multiplier(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: solid_traction_elements.h:1676
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.
HangInfo *const & hanging_pt() const
Definition: nodes.h:1228
bool is_hanging() const
Test whether the node is geometrically hanging.
Definition: nodes.h:1285
double value(const unsigned &i) const
Definition: nodes.cc:2408
Definition: refineable_elements.h:1042
Definition: oomph_definitions.h:222
Definition: refineable_elements.h:97
int local_hang_eqn(Node *const &node_pt, const unsigned &i)
Definition: refineable_elements.h:278
Definition: solid_traction_elements.h:3073
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
Definition: solid_traction_elements.h:3118
RefineableFSIImposeDisplacementByLagrangeMultiplierElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Definition: solid_traction_elements.h:3091
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals.
Definition: solid_traction_elements.h:3109
void refineable_fill_in_generic_contribution_to_residuals_fsi_displ_lagr_multiplier(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: solid_traction_elements.h:3133
unsigned ncont_interpolated_values() const
Definition: solid_traction_elements.h:3103
Definition: solid_traction_elements.h:1046
virtual void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: solid_traction_elements.h:1071
~RefineableFSISolidTractionElement()
Destructor: empty.
Definition: solid_traction_elements.h:1066
RefineableFSISolidTractionElement(FiniteElement *const &element_pt, const int &face_index)
Definition: solid_traction_elements.h:1057
Definition: solid_traction_elements.h:2013
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals.
Definition: solid_traction_elements.h:2038
void refineable_fill_in_generic_contribution_to_residuals_displ_lagr_multiplier(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: solid_traction_elements.h:2063
unsigned ncont_interpolated_values() const
Definition: solid_traction_elements.h:2032
RefineableImposeDisplacementByLagrangeMultiplierElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Definition: solid_traction_elements.h:2020
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
Definition: solid_traction_elements.h:2047
DenseMatrix< int > & local_position_hang_eqn(Node *const &node_pt)
Definition: refineable_elements.h:1005
Definition: solid_traction_elements.h:474
~RefineableSolidTractionElement()
Destructor should not delete anything.
Definition: solid_traction_elements.h:485
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
Definition: solid_traction_elements.h:496
RefineableSolidTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the face index.
Definition: solid_traction_elements.h:477
unsigned ncont_interpolated_values() const
Definition: solid_traction_elements.h:489
void refineable_fill_in_contribution_to_residuals_solid_traction(Vector< double > &residuals)
Definition: solid_traction_elements.h:539
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the Jacobian.
Definition: solid_traction_elements.h:504
Definition: elements.h:4914
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: elements.h:4921
double interpolated_xi(const Vector< double > &s, const unsigned &i) const
Definition: elements.h:4937
double lagrangian_position(const unsigned &n, const unsigned &i) const
Return i-th Lagrangian coordinate at local node n.
Definition: elements.h:3905
Data * geom_data_pt(const unsigned &j)
Definition: elements.h:3615
int position_local_eqn(const unsigned &n, const unsigned &k, const unsigned &j) const
Definition: elements.h:4137
Definition: solid_traction_elements.h:78
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
Definition: solid_traction_elements.h:199
void fill_in_contribution_to_residuals_solid_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
Definition: solid_traction_elements.h:317
SolidTractionElement(FiniteElement *const &element_pt, const int &face_index, const bool &called_from_refineable_constructor=false)
Definition: solid_traction_elements.h:118
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
Definition: solid_traction_elements.h:176
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.
Definition: solid_traction_elements.h:159
void output(std::ostream &outfile)
Output function.
Definition: solid_traction_elements.h:191
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
Definition: solid_traction_elements.h:169
void traction(const Vector< double > &s, Vector< double > &traction)
Definition: solid_traction_elements.h:286
virtual void get_traction(const unsigned &intpt, const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Definition: solid_traction_elements.h:97
void output(FILE *file_pt)
C_style output function.
Definition: solid_traction_elements.h:258
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
Definition: solid_traction_elements.h:264
void(* Traction_fct_pt)(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Definition: solid_traction_elements.h:85
@ N
Definition: constructor.cpp:22
RealScalar s
Definition: level1_cplx_impl.h:130
const Scalar * a
Definition: level2_cplx_impl.h:32
int * m
Definition: level2_cplx_impl.h:294
char char char int int * k
Definition: level2_impl.h:374
#define DIM
Definition: linearised_navier_stokes_elements.h:44
EIGEN_STRONG_INLINE const Eigen::CwiseBinaryOp< Eigen::internal::scalar_zeta_op< typename DerivedX::Scalar >, const DerivedX, const DerivedQ > zeta(const Eigen::ArrayBase< DerivedX > &x, const Eigen::ArrayBase< DerivedQ > &q)
Definition: SpecialFunctionsArrayAPI.h:152
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
@ W
Definition: quadtree.h:63
void Zero_traction_fct(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Default load function (zero traction)
Definition: solid_traction_elements.h:53
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
list x
Definition: plotDoE.py:28
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
unsigned el_dim
dimension
Definition: overloaded_cartesian_element_body.h:30
std::ofstream out("Result.txt")
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2