28 #ifndef OOMPH_QELEMENT_HEADER
29 #define OOMPH_QELEMENT_HEADER
33 #include <oomph-lib-config.h>
115 unsigned ncoord =
dim();
116 for (
unsigned i = 0;
i < ncoord;
i++)
131 unsigned ncoord =
dim();
132 for (
unsigned i = 0;
i < ncoord;
i++)
146 unsigned n_dim =
dim();
173 for (
unsigned i = 0;
i < n_dim;
i++)
192 OOMPH_CURRENT_FUNCTION,
193 OOMPH_EXCEPTION_LOCATION);
208 OOMPH_CURRENT_FUNCTION,
209 OOMPH_EXCEPTION_LOCATION);
224 OOMPH_CURRENT_FUNCTION,
225 OOMPH_EXCEPTION_LOCATION);
240 OOMPH_CURRENT_FUNCTION,
241 OOMPH_EXCEPTION_LOCATION);
255 throw OomphLibError(
"Macro Element pointer not set in this element\n",
256 OOMPH_CURRENT_FUNCTION,
257 OOMPH_EXCEPTION_LOCATION);
261 unsigned el_dim =
dim();
263 for (
unsigned i = 0;
i < el_dim;
i++)
281 throw OomphLibError(
"Macro Element pointer not set in this element\n",
282 OOMPH_CURRENT_FUNCTION,
283 OOMPH_EXCEPTION_LOCATION);
287 unsigned el_dim =
dim();
289 for (
unsigned i = 0;
i < el_dim;
i++)
302 return static_cast<unsigned>(std::pow(
static_cast<double>(
nnode_1d()),
303 static_cast<double>(
dim() - 1)));
387 unsigned n_xi = xi_fe.size();
388 for (
unsigned i = 0;
i < n_xi;
i++)
396 unsigned n_xi = xi.size();
397 for (
unsigned i = 0;
i < n_xi;
i++)
405 unsigned el_dim =
dim();
407 for (
unsigned i = 0;
i < el_dim;
i++)
417 unsigned n_x = x_fe.size();
418 for (
unsigned i = 0;
i < n_x;
i++)
426 for (
unsigned i = 0;
i < n_x;
i++)
434 unsigned el_dim =
dim();
436 for (
unsigned i = 0;
i < el_dim;
i++)
457 template<
unsigned DIM,
unsigned NNODE_1D>
481 template<
unsigned NNODE_1D>
497 this->set_n_node(NNODE_1D);
499 this->set_dimension(1);
501 this->set_integration_scheme(&Default_integration_scheme);
533 return FiniteElement::invert_jacobian<1>(jacobian, inverse_jacobian);
557 unsigned n_node_1d = nnode_1d();
565 nod_pt = node_pt(n_node_1d - 1);
568 std::ostringstream error_message;
569 error_message <<
"Vertex node number is " << j
570 <<
" but must be from 0 to 1\n";
573 OOMPH_CURRENT_FUNCTION,
574 OOMPH_EXCEPTION_LOCATION);
583 s[0] = this->s_min() +
584 double(j) / double(NNODE_1D - 1) * (this->s_max() - this->s_min());
591 s_fraction.resize(1);
592 s_fraction[0] = double(j) / double(NNODE_1D - 1);
602 return double(n1d) / double(NNODE_1D - 1);
632 const unsigned& nplot,
633 unsigned& counter)
const
636 unsigned plot = nsub_elements_paraview(nplot);
639 for (
unsigned i = 0;
i < plot;
i++)
641 file_out <<
i + counter <<
" " <<
i + 1 + counter << std::endl;
643 counter += nplot_points_paraview(nplot);
651 const unsigned& nplot)
const
653 unsigned local_loop = nsub_elements_paraview(nplot);
654 for (
unsigned i = 0;
i < local_loop;
i++)
656 file_out <<
"3" << std::endl;
664 const unsigned& nplot,
665 unsigned& offset_sum)
const
669 unsigned local_loop = nsub_elements_paraview(nplot);
670 for (
unsigned i = 0;
i < local_loop;
i++)
673 file_out << offset_sum << std::endl;
678 void output(std::ostream& outfile);
681 void output(std::ostream& outfile,
const unsigned& n_plot);
684 void output(FILE* file_pt);
687 void output(FILE* file_pt,
const unsigned& n_plot);
694 const unsigned& nplot,
696 const bool& use_equally_spaced_interior_sample_points =
false)
const
700 s[0] = -1.0 + 2.0 * double(
i) / double(nplot - 1);
701 if (use_equally_spaced_interior_sample_points)
704 double dx_new = range / double(nplot);
705 double range_new = double(nplot - 1) * dx_new;
706 s[0] = -1.0 + 0.5 * dx_new + range_new * (1.0 +
s[0]) / range;
719 std::ostringstream header;
720 header <<
"ZONE I=" << nplot <<
"\n";
730 for (
unsigned i = 0;
i < DIM;
i++)
740 const unsigned&
i)
const
742 face_node_number_error_check(
i);
744 if (face_index == -1)
748 else if (face_index == +1)
750 return nnode_1d() - 1;
754 std::string err =
"Face index should be in {-1, +1}.";
756 err, OOMPH_EXCEPTION_LOCATION, OOMPH_CURRENT_FUNCTION);
764 if (std::abs(face_index) != 1)
766 std::string err =
"Face index should be in {-1, +1}.";
768 err, OOMPH_EXCEPTION_LOCATION, OOMPH_CURRENT_FUNCTION);
777 const int& face_index)
const
783 else if (face_index == -1)
789 std::string err =
"Face index should be in {-1, +1}.";
791 err, OOMPH_EXCEPTION_LOCATION, OOMPH_CURRENT_FUNCTION);
798 const int& face_index)
const
804 else if (face_index == -1)
810 std::string err =
"Face index should be in {-1, +1}.";
812 err, OOMPH_EXCEPTION_LOCATION, OOMPH_CURRENT_FUNCTION);
836 template<
unsigned NNODE_1D>
852 this->set_n_node(NNODE_1D * NNODE_1D);
856 set_integration_scheme(&Default_integration_scheme);
890 return FiniteElement::invert_jacobian<2>(jacobian, inverse_jacobian);
915 unsigned n_node_1d = nnode_1d();
923 nod_pt = node_pt(n_node_1d - 1);
926 nod_pt = node_pt(n_node_1d * (n_node_1d - 1));
929 nod_pt = node_pt(n_node_1d * n_node_1d - 1);
932 std::ostringstream error_message;
933 error_message <<
"Vertex node number is " << j
934 <<
" but must be from 0 to 3\n";
937 OOMPH_CURRENT_FUNCTION,
938 OOMPH_EXCEPTION_LOCATION);
948 unsigned j0 = j % NNODE_1D;
949 unsigned j1 = j / NNODE_1D;
950 const double S_min = this->s_min();
951 const double S_range = this->s_max() - S_min;
952 s[0] = S_min + double(j0) / double(NNODE_1D - 1) * S_range;
953 s[1] = S_min + double(j1) / double(NNODE_1D - 1) * S_range;
959 s_fraction.resize(2);
960 unsigned j0 = j % NNODE_1D;
961 unsigned j1 = j / NNODE_1D;
962 s_fraction[0] = double(j0) / double(NNODE_1D - 1);
963 s_fraction[1] = double(j1) / double(NNODE_1D - 1);
972 return double(n1d) / double(NNODE_1D - 1);
988 return nplot * nplot;
995 return (nplot - 1) * (nplot - 1);
1002 const unsigned& nplot,
1003 unsigned& counter)
const
1006 unsigned plot = nsub_elements_paraview(nplot);
1009 for (
unsigned i = 0;
i < plot;
i++)
1011 unsigned d = (
i - (
i % (nplot - 1))) / (nplot - 1);
1013 file_out <<
i % (nplot - 1) + d * nplot + counter <<
" "
1014 <<
i % (nplot - 1) + 1 + d * nplot + counter <<
" "
1015 <<
i % (nplot - 1) + 1 + (d + 1) * nplot + counter <<
" "
1016 <<
i % (nplot - 1) + (d + 1) * nplot + counter << std::endl;
1018 counter += nplot_points_paraview(nplot);
1026 const unsigned& nplot)
const
1028 unsigned local_loop = nsub_elements_paraview(nplot);
1029 for (
unsigned i = 0;
i < local_loop;
i++)
1031 file_out <<
"9" << std::endl;
1039 const unsigned& nplot,
1040 unsigned& offset_sum)
const
1044 unsigned local_loop = nsub_elements_paraview(nplot);
1045 for (
unsigned i = 0;
i < local_loop;
i++)
1048 file_out << offset_sum << std::endl;
1053 void output(std::ostream& outfile);
1056 void output(std::ostream& outfile,
const unsigned& n_plot);
1059 void output(FILE* file_pt);
1062 void output(FILE* file_pt,
const unsigned& n_plot);
1069 const unsigned& nplot,
1071 const bool& use_equally_spaced_interior_sample_points =
false)
const
1075 unsigned i0 =
i % nplot;
1076 unsigned i1 = (
i - i0) / nplot;
1078 s[0] = -1.0 + 2.0 * double(i0) / double(nplot - 1);
1079 s[1] = -1.0 + 2.0 * double(i1) / double(nplot - 1);
1080 if (use_equally_spaced_interior_sample_points)
1083 double dx_new = range / double(nplot);
1084 double range_new = double(nplot - 1) * dx_new;
1085 s[0] = -1.0 + 0.5 * dx_new + range_new * (1.0 +
s[0]) / range;
1086 s[1] = -1.0 + 0.5 * dx_new + range_new * (1.0 +
s[1]) / range;
1100 std::ostringstream header;
1101 header <<
"ZONE I=" << nplot <<
", J=" << nplot <<
"\n";
1102 return header.str();
1111 for (
unsigned i = 0;
i < DIM;
i++)
1121 const unsigned&
i)
const
1123 face_node_number_error_check(
i);
1125 const unsigned nn1d = nnode_1d();
1127 if (face_index == -1)
1131 else if (face_index == +1)
1133 return nn1d *
i + nn1d - 1;
1135 else if (face_index == -2)
1139 else if (face_index == +2)
1141 return nn1d * (nn1d - 1) +
i;
1145 std::string err =
"Face index should be in {-1, +1, -2, +2}.";
1147 err, OOMPH_EXCEPTION_LOCATION, OOMPH_CURRENT_FUNCTION);
1158 else if (face_index > 0)
1164 std::string err =
"Face index should be one of {-1, +1, -2, +2}.";
1166 err, OOMPH_EXCEPTION_LOCATION, OOMPH_CURRENT_FUNCTION);
1173 const int& face_index)
const
1175 if (face_index == 1)
1179 else if (face_index == -1)
1183 else if (face_index == -2)
1187 else if (face_index == 2)
1193 std::string err =
"Face index should be in {-1, +1}.";
1195 err, OOMPH_EXCEPTION_LOCATION, OOMPH_CURRENT_FUNCTION);
1202 const int& face_index)
const
1204 if (face_index == 1)
1208 else if (face_index == -1)
1212 else if (face_index == -2)
1216 else if (face_index == 2)
1222 std::string err =
"Face index should be in {-1, +1}.";
1224 err, OOMPH_EXCEPTION_LOCATION, OOMPH_CURRENT_FUNCTION);
1248 template<
unsigned NNODE_1D>
1264 this->set_n_node(NNODE_1D * NNODE_1D * NNODE_1D);
1268 set_integration_scheme(&Default_integration_scheme);
1307 return FiniteElement::invert_jacobian<3>(jacobian, inverse_jacobian);
1331 unsigned N = nnode_1d();
1336 nod_pt = node_pt(0);
1339 nod_pt = node_pt(
N - 1);
1342 nod_pt = node_pt(
N * (
N - 1));
1345 nod_pt = node_pt(
N *
N - 1);
1348 nod_pt = node_pt(
N *
N * (
N - 1));
1351 nod_pt = node_pt(
N *
N * (
N - 1) + (
N - 1));
1354 nod_pt = node_pt(
N *
N *
N -
N);
1357 nod_pt = node_pt(
N *
N *
N - 1);
1360 std::ostringstream error_message;
1361 error_message <<
"Vertex node number is " << j
1362 <<
" but must be from 0 to 7\n";
1365 OOMPH_CURRENT_FUNCTION,
1366 OOMPH_EXCEPTION_LOCATION);
1375 unsigned j0 = j % NNODE_1D;
1376 unsigned j1 = (j / NNODE_1D) % NNODE_1D;
1377 unsigned j2 = j / (NNODE_1D * NNODE_1D);
1378 const double S_min = this->s_min();
1379 const double S_range = this->s_max() - S_min;
1381 s[0] = S_min + double(j0) / double(NNODE_1D - 1) * S_range;
1382 s[1] = S_min + double(j1) / double(NNODE_1D - 1) * S_range;
1383 s[2] = S_min + double(j2) / double(NNODE_1D - 1) * S_range;
1389 s_fraction.resize(3);
1390 unsigned j0 = j % NNODE_1D;
1391 unsigned j1 = (j / NNODE_1D) % NNODE_1D;
1392 unsigned j2 = j / (NNODE_1D * NNODE_1D);
1393 s_fraction[0] = double(j0) / double(NNODE_1D - 1);
1394 s_fraction[1] = double(j1) / double(NNODE_1D - 1);
1395 s_fraction[2] = double(j2) / double(NNODE_1D - 1);
1404 return double(n1d) / double(NNODE_1D - 1);
1420 return nplot * nplot * nplot;
1427 return (nplot - 1) * (nplot - 1) * (nplot - 1);
1434 const unsigned& nplot,
1435 unsigned& counter)
const
1438 unsigned plot = nsub_elements_paraview(nplot);
1440 for (
unsigned j = 0; j < plot; j += (nplot - 1) * (nplot - 1) + 1)
1443 unsigned r = ((j - (j % ((nplot - 1) * (nplot - 1)))) /
1444 ((nplot - 1) * (nplot - 1)));
1447 unsigned sub_plot = (nplot - 1) * (nplot - 1);
1450 for (
unsigned i = 0;
i < sub_plot;
i++)
1452 unsigned d = ((
i - (
i % (nplot - 1))) / (nplot - 1));
1457 <<
i % (nplot - 1) + d * nplot + r * nplot * nplot + counter <<
" "
1458 <<
i % (nplot - 1) + 1 + d * nplot + r * nplot * nplot + counter
1460 <<
i % (nplot - 1) + 1 + (d + 1) * nplot + r * nplot * nplot +
1463 <<
i % (nplot - 1) + (d + 1) * nplot + r * nplot * nplot + counter
1467 <<
i % (nplot - 1) + d * nplot + (r + 1) * nplot * nplot + counter
1469 <<
i % (nplot - 1) + 1 + d * nplot + (r + 1) * nplot * nplot +
1472 <<
i % (nplot - 1) + 1 + (d + 1) * nplot + (r + 1) * nplot * nplot +
1475 <<
i % (nplot - 1) + (d + 1) * nplot + (r + 1) * nplot * nplot +
1480 counter += nplot_points_paraview(nplot);
1488 const unsigned& nplot)
const
1490 unsigned local_loop = nsub_elements_paraview(nplot);
1491 for (
unsigned i = 0;
i < local_loop;
i++)
1493 file_out <<
"12" << std::endl;
1501 const unsigned& nplot,
1502 unsigned& offset_sum)
const
1504 unsigned local_loop = nsub_elements_paraview(nplot);
1505 for (
unsigned i = 0;
i < local_loop;
i++)
1508 file_out << offset_sum << std::endl;
1513 void output(std::ostream& outfile);
1516 void output(std::ostream& outfile,
const unsigned& n_plot);
1519 void output(FILE* file_pt);
1522 void output(FILE* file_pt,
const unsigned& n_plot);
1528 const unsigned& nplot,
1530 const bool& use_equally_spaced_interior_sample_points =
false)
const
1534 unsigned i01 =
i % (nplot * nplot);
1535 unsigned i0 = i01 % nplot;
1536 unsigned i1 = (i01 - i0) / nplot;
1537 unsigned i2 = (
i - i01) / (nplot * nplot);
1539 s[0] = -1.0 + 2.0 * double(i0) / double(nplot - 1);
1540 s[1] = -1.0 + 2.0 * double(i1) / double(nplot - 1);
1541 s[2] = -1.0 + 2.0 * double(i2) / double(nplot - 1);
1542 if (use_equally_spaced_interior_sample_points)
1545 double dx_new = range / double(nplot);
1546 double range_new = double(nplot - 1) * dx_new;
1547 s[0] = -1.0 + 0.5 * dx_new + range_new * (1.0 +
s[0]) / range;
1548 s[1] = -1.0 + 0.5 * dx_new + range_new * (1.0 +
s[1]) / range;
1549 s[2] = -1.0 + 0.5 * dx_new + range_new * (1.0 +
s[2]) / range;
1564 std::ostringstream header;
1565 header <<
"ZONE I=" << nplot <<
", J=" << nplot <<
", K=" << nplot
1567 return header.str();
1576 for (
unsigned i = 0;
i < DIM;
i++)
1586 const unsigned&
i)
const
1588 face_node_number_error_check(
i);
1590 const unsigned nn1d = nnode_1d();
1592 if (face_index == -1)
1596 else if (face_index == +1)
1598 return i * nn1d + (nn1d - 1);
1600 else if (face_index == -2)
1602 return i % nn1d + (
i / nn1d) * nn1d * nn1d;
1604 else if (face_index == +2)
1606 return i % nn1d + (
i / nn1d) * nn1d * nn1d + (nn1d * (nn1d - 1));
1608 else if (face_index == -3)
1612 else if (face_index == +3)
1614 return i + (nn1d * nn1d) * (nn1d - 1);
1618 std::string err =
"Face index should be in {-1, +1, -2, +2, -3, +3}.";
1620 err, OOMPH_EXCEPTION_LOCATION, OOMPH_CURRENT_FUNCTION);
1627 if (face_index == -3)
1631 else if (face_index == +3)
1635 else if (face_index == -2)
1639 else if (face_index == 2)
1643 else if (face_index == -1)
1647 else if (face_index == 1)
1653 std::string err =
"Face index should be in {-1, +1, -2, +2, -3, +3}.";
1655 err, OOMPH_EXCEPTION_LOCATION, OOMPH_CURRENT_FUNCTION);
1662 const int& face_index)
const
1664 if (face_index == 1)
1668 else if (face_index == -1)
1672 else if (face_index == -2)
1676 else if (face_index == 2)
1680 else if (face_index == -3)
1684 else if (face_index == 3)
1690 std::string err =
"Face index should be in {-1, +1}.";
1692 err, OOMPH_EXCEPTION_LOCATION, OOMPH_CURRENT_FUNCTION);
1699 const int& face_index)
const
1701 if (face_index == 1)
1705 else if (face_index == -1)
1709 else if (face_index == -2)
1713 else if (face_index == 2)
1717 else if (face_index == -3)
1721 else if (face_index == 3)
1727 std::string err =
"Face index should be in {-1, +1}.";
1729 err, OOMPH_EXCEPTION_LOCATION, OOMPH_CURRENT_FUNCTION);
1740 template<
unsigned DIM,
unsigned NNODE_1D>
1749 template<
unsigned NNODE_1D>
1758 set_lagrangian_dimension(1);
1774 inline void output(std::ostream& outfile,
const unsigned& n_plot);
1783 inline void output(FILE* file_pt,
const unsigned& n_plot);
1791 inline void build_face_element(
const int& face_index,
1812 template<
unsigned NNODE_1D>
1814 const unsigned& n_plot)
1820 outfile <<
"ZONE I=" << n_plot << std::endl;
1823 unsigned n_dim = this->nodal_dimension();
1826 unsigned n_lagr =
static_cast<SolidNode*
>(node_pt(0))->nlagrangian();
1829 for (
unsigned l = 0; l < n_plot; l++)
1831 s[0] = -1.0 + l * 2.0 / (n_plot - 1);
1834 for (
unsigned i = 0;
i < n_dim;
i++)
1836 outfile << QElement<1, NNODE_1D>::interpolated_x(
s,
i) <<
" ";
1839 for (
unsigned i = 0;
i < n_lagr;
i++)
1841 outfile << SolidQElement<1, NNODE_1D>::interpolated_xi(
s,
i) <<
" ";
1843 outfile << std::endl;
1845 outfile << std::endl;
1852 template<
unsigned NNODE_1D>
1859 fprintf(file_pt,
"ZONE I=%i\n", n_plot);
1862 unsigned n_dim = this->nodal_dimension();
1865 unsigned n_lagr =
static_cast<SolidNode*
>(node_pt(0))->nlagrangian();
1868 for (
unsigned l = 0; l < n_plot; l++)
1870 s[0] = -1.0 + l * 2.0 / (n_plot - 1);
1873 for (
unsigned i = 0;
i < n_dim;
i++)
1878 for (
unsigned i = 0;
i < n_lagr;
i++)
1883 fprintf(file_pt,
"\n");
1885 fprintf(file_pt,
"\n");
1893 template<
unsigned NNODE_1D>
1895 const int& face_index,
FaceElement* face_element_pt)
1902 ->set_lagrangian_dimension(
1910 template<
unsigned NNODE_1D>
1923 set_lagrangian_dimension(2);
1939 inline void output(std::ostream& outfile,
const unsigned& n_plot);
1948 inline void output(FILE* file_pt,
const unsigned& n_plot);
1958 inline void build_face_element(
const int& face_index,
1970 template<
unsigned NNODE_1D>
1972 const unsigned& n_plot)
1978 outfile <<
"ZONE I=" << n_plot <<
", J=" << n_plot << std::endl;
1981 unsigned n_dim = this->nodal_dimension();
1984 unsigned n_lagr =
static_cast<SolidNode*
>(node_pt(0))->nlagrangian();
1987 for (
unsigned l2 = 0; l2 < n_plot; l2++)
1989 s[1] = -1.0 + l2 * 2.0 / (n_plot - 1);
1990 for (
unsigned l1 = 0; l1 < n_plot; l1++)
1992 s[0] = -1.0 + l1 * 2.0 / (n_plot - 1);
1995 for (
unsigned i = 0;
i < n_dim;
i++)
1997 outfile << QElement<2, NNODE_1D>::interpolated_x(
s,
i) <<
" ";
2000 for (
unsigned i = 0;
i < n_lagr;
i++)
2002 outfile << SolidQElement<2, NNODE_1D>::interpolated_xi(
s,
i) <<
" ";
2004 outfile << std::endl;
2007 outfile << std::endl;
2014 template<
unsigned NNODE_1D>
2021 fprintf(file_pt,
"ZONE I=%i, J=%i\n", n_plot, n_plot);
2024 unsigned n_dim = this->nodal_dimension();
2027 unsigned n_lagr =
static_cast<SolidNode*
>(node_pt(0))->nlagrangian();
2030 for (
unsigned l2 = 0; l2 < n_plot; l2++)
2032 s[1] = -1.0 + l2 * 2.0 / (n_plot - 1);
2033 for (
unsigned l1 = 0; l1 < n_plot; l1++)
2035 s[0] = -1.0 + l1 * 2.0 / (n_plot - 1);
2038 for (
unsigned i = 0;
i < n_dim;
i++)
2043 for (
unsigned i = 0;
i < n_lagr;
i++)
2048 fprintf(file_pt,
"\n");
2051 fprintf(file_pt,
"\n");
2059 template<
unsigned NNODE_1D>
2061 const int& face_index,
FaceElement* face_element_pt)
2068 ->set_lagrangian_dimension(
2076 template<
unsigned NNODE_1D>
2089 set_lagrangian_dimension(3);
2105 inline void output(std::ostream& outfile,
const unsigned& n_plot);
2114 inline void output(FILE* file_pt,
const unsigned& n_plot);
2127 inline void build_face_element(
const int& face_index,
2139 template<
unsigned NNODE_1D>
2141 const unsigned& n_plot)
2147 outfile <<
"ZONE I=" << n_plot <<
", J=" << n_plot <<
", K=" << n_plot
2151 unsigned n_dim = this->nodal_dimension();
2154 unsigned n_lagr =
static_cast<SolidNode*
>(node_pt(0))->nlagrangian();
2157 for (
unsigned l3 = 0; l3 < n_plot; l3++)
2159 s[2] = -1.0 + l3 * 2.0 / (n_plot - 1);
2160 for (
unsigned l2 = 0; l2 < n_plot; l2++)
2162 s[1] = -1.0 + l2 * 2.0 / (n_plot - 1);
2163 for (
unsigned l1 = 0; l1 < n_plot; l1++)
2165 s[0] = -1.0 + l1 * 2.0 / (n_plot - 1);
2168 for (
unsigned i = 0;
i < n_dim;
i++)
2170 outfile << QElement<3, NNODE_1D>::interpolated_x(
s,
i) <<
" ";
2173 for (
unsigned i = 0;
i < n_lagr;
i++)
2175 outfile << SolidQElement<3, NNODE_1D>::interpolated_xi(
s,
i) <<
" ";
2177 outfile << std::endl;
2181 outfile << std::endl;
2188 template<
unsigned NNODE_1D>
2195 fprintf(file_pt,
"ZONE I=%i, J=%i, K=%i\n", n_plot, n_plot, n_plot);
2198 unsigned n_dim = this->nodal_dimension();
2201 unsigned n_lagr =
static_cast<SolidNode*
>(node_pt(0))->nlagrangian();
2204 for (
unsigned l3 = 0; l3 < n_plot; l3++)
2206 s[2] = -1.0 + l3 * 2.0 / (n_plot - 1);
2207 for (
unsigned l2 = 0; l2 < n_plot; l2++)
2209 s[1] = -1.0 + l2 * 2.0 / (n_plot - 1);
2210 for (
unsigned l1 = 0; l1 < n_plot; l1++)
2212 s[0] = -1.0 + l1 * 2.0 / (n_plot - 1);
2215 for (
unsigned i = 0;
i < n_dim;
i++)
2221 for (
unsigned i = 0;
i < n_lagr;
i++)
2227 fprintf(file_pt,
"\n");
2231 fprintf(file_pt,
"\n");
2239 template<
unsigned NNODE_1D>
2241 const int& face_index,
FaceElement* face_element_pt)
2248 ->set_lagrangian_dimension(
2257 template<
unsigned DIM>
2272 template<
unsigned DIM,
unsigned INITIAL_NNODE_1D = 2>
2284 template<
unsigned DIM>
Base class for all brick elements.
virtual Node * vertex_node_pt(const unsigned &j) const =0
Pointer to the j-th vertex node in the element.
virtual unsigned nvertex_node() const =0
Number of vertex nodes in the element.
BrickElementBase()
Constructor. Empty.
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
A general Finite Element class.
virtual void set_macro_elem_pt(MacroElement *macro_elem_pt)
Set pointer to macro element – can be overloaded in derived elements to perform additional tasks.
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
virtual double s_min() const
Min value of local coordinate.
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
MacroElement * Macro_elem_pt
Pointer to the element's macro element (NULL by default)
virtual double s_max() const
Max. value of local coordinate.
MacroElement * macro_elem_pt()
Access function to pointer to macro element.
virtual unsigned nnode_1d() const
Return the number of nodes along one edge of the element Default is to return zero — must be overload...
Base class for all line elements.
LineElementBase()
Constructor. Empty.
virtual Node * vertex_node_pt(const unsigned &j) const =0
Pointer to the j-th vertex node in the element.
virtual unsigned nvertex_node() const =0
Number of vertex nodes in the element.
Base class for MacroElement s that are used during mesh refinement in domains with curvlinear and/or ...
void macro_map(const Vector< double > &s, Vector< double > &r)
The mapping from local to global coordinates at the current time : r(s)
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
An OomphLibError object which should be thrown when an run-time error is encountered....
A class that is used to template the p-refineable Q elements by dimension. It's really nothing more t...
PRefineableQElement()
Empty constuctor.
/////////////////////////////////////////////////////////////////// /////////////////////////////////...
virtual void set_macro_elem_pt(MacroElement *macro_elem_pt)
Set pointer to macro element also sets up storage for the reference coordinates and initialises them.
QElementBase(const QElementBase &)=delete
Broken copy constructor.
Vector< double > * S_macro_ur_pt
Pointer to vector of upper right vertex coords. in macro element.
QElementBase()
Constructor: Initialise pointers to macro element reference coords.
bool local_coord_is_valid(const Vector< double > &s)
Check whether the local coordinate are valid or not.
ElementGeometry::ElementGeometry element_geometry() const
It's a Q element!
void get_x_from_macro_element(const unsigned &t, const Vector< double > &s, Vector< double > &x)
Global coordinates as function of local coordinates at previous time "level" t (t=0: present; t>0: pr...
double & s_macro_ll(const unsigned &i)
Access fct to the i-th coordinate of the element's "lower left" vertex in the associated MacroElement...
void get_x_from_macro_element(const Vector< double > &s, Vector< double > &x) const
Global coordinates as function of local coordinates. using the macro element representation.
void move_local_coord_back_into_element(Vector< double > &s) const
Adjust local coordinates so that they're located inside the element.
double s_macro_ll(const unsigned &i) const
Access fct to the i-th coordinate of the element's "lower left" vertex in the associated MacroElement...
unsigned nnode_on_face() const
Return number of nodes on one face of the element. Always nnode_1d^(el_dim - 1).
Vector< double > * S_macro_ll_pt
Pointer to vector of lower left vertex coords. in macro element.
double s_macro_ur(const unsigned &i) const
Access fct to the i-th coordinate of the element's "upper right" vertex in the associated MacroElemen...
double & s_macro_ur(const unsigned &i)
Access fct to the i-th coordinate of the element's "upper right" vertex in the associated MacroElemen...
virtual ~QElementBase()
Broken assignment operator.
/////////////////////////////////////////////////////////////////// /////////////////////////////////...
QElementGeometricBase()
Empty default constructor.
QElementGeometricBase(const QElementGeometricBase &)=delete
Broken copy constructor.
std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction)
BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt(const int &face_index) const
Get a pointer to the derivative of the mapping from face to bulk coordinates.
void write_paraview_output_offset_information(std::ofstream &file_out, const unsigned &nplot, unsigned &counter) const
Fill in the offset information for paraview plot. Needs to be implemented for each new geometric elem...
QElement(const QElement &)=delete
Broken copy constructor.
void local_fraction_of_node(const unsigned &j, Vector< double > &s_fraction)
Get the local fraction of node j in the element.
unsigned get_bulk_node_number(const int &face_index, const unsigned &i) const
Get the number of the ith node on face face_index in the bulk node vector.
unsigned nvertex_node() const
Number of vertex nodes in the element.
unsigned nplot_points_paraview(const unsigned &nplot) const
Return the number of actual plot points for paraview plot with parameter nplot.
static Gauss< 1, NNODE_1D > Default_integration_scheme
Default integration rule: Gaussian integration of same 'order' as the element.
double invert_jacobian_mapping(const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
Overload the template-free interface for the calculation of inverse jacobian matrix....
CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt(const int &face_index) const
Get a pointer to the function mapping face coordinates to bulk coordinates.
double s_min() const
Min. value of local coordinate.
unsigned nnode_1d() const
Number of nodes along each element edge.
double local_one_d_fraction_of_node(const unsigned &n1d, const unsigned &i)
This function returns the local fraction of all nodes at the n-th position in a one dimensional expan...
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
unsigned nsub_elements_paraview(const unsigned &nplot) const
Return the number of local sub-elements for paraview plot with parameter nplot.
void write_paraview_type(std::ofstream &file_out, const unsigned &nplot) const
Return the paraview element type. Needs to be implemented for each new geometric element type; see ht...
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Get local coordinates of node j in the element; vector sets its own size.
unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction)
void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &use_equally_spaced_interior_sample_points=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
int face_outer_unit_normal_sign(const int &face_index) const
Get the sign of the outer unit normal on the face given by face_index.
double s_max() const
Max. value of local coordinate.
void write_paraview_offsets(std::ofstream &file_out, const unsigned &nplot, unsigned &offset_sum) const
Return the offsets for the paraview sub-elements. Needs to be implemented for each new geometric elem...
CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt(const int &face_index) const
Get a pointer to the function mapping face coordinates to bulk coordinates.
double s_min() const
Min. value of local coordinate.
unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction)
void write_paraview_type(std::ofstream &file_out, const unsigned &nplot) const
Return the paraview element type. Needs to be implemented for each new geometric element type; see ht...
double invert_jacobian_mapping(const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
Overload the template-free interface for the calculation of inverse jacobian matrix....
void write_paraview_offsets(std::ofstream &file_out, const unsigned &nplot, unsigned &offset_sum) const
Return the offsets for the paraview sub-elements. Needs to be implemented for each new geometric elem...
unsigned nvertex_node() const
Number of vertex nodes in the element.
BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt(const int &face_index) const
Get a pointer to the derivative of the mapping from face to bulk coordinates.
void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &use_equally_spaced_interior_sample_points=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
void local_fraction_of_node(const unsigned &j, Vector< double > &s_fraction)
Get the local fraction of node j in the element.
unsigned nnode_1d() const
Number of nodes along each element edge.
unsigned nplot_points_paraview(const unsigned &nplot) const
Return the number of actual plot points for paraview plot with parameter nplot.
void write_paraview_output_offset_information(std::ofstream &file_out, const unsigned &nplot, unsigned &counter) const
Fill in the offset information for paraview plot. Needs to be implemented for each new geometric elem...
static Gauss< 2, NNODE_1D > Default_integration_scheme
Default integration rule: Gaussian integration of same 'order' as the element.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
QElement(const QElement &)=delete
Broken copy constructor.
std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction)
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Get local coordinates of node j in the element; vector sets its own size.
unsigned nsub_elements_paraview(const unsigned &nplot) const
Return the number of local sub-elements for paraview plot with parameter nplot.
double local_one_d_fraction_of_node(const unsigned &n1d, const unsigned &i)
This function returns the local fraction of ant nodes in the n-th positoin in a one dimensional expan...
double s_max() const
Max. value of local coordinate.
int face_outer_unit_normal_sign(const int &face_index) const
Get the sign of the outer unit normal on the face given by face_index.
unsigned get_bulk_node_number(const int &face_index, const unsigned &i) const
Get the number of the ith node on face face_index (in the bulk node vector).
unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction)
double local_one_d_fraction_of_node(const unsigned &n1d, const unsigned &i)
This function returns the local fraction of any nodes in the n-th positoin in a one dimensional expan...
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
int face_outer_unit_normal_sign(const int &face_index) const
Get the sign of the outer unit normal on the face given by face_index.
BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt(const int &face_index) const
Get a pointer to the derivative of the mapping from face to bulk coordinates.
double s_max() const
Max. value of local coordinate.
unsigned nnode_1d() const
Number of nodes along each element edge.
void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &use_equally_spaced_interior_sample_points=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
double s_min() const
Min. value of local coordinate.
unsigned nplot_points_paraview(const unsigned &nplot) const
Return the number of actual plot points for paraview plot with parameter nplot.
CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt(const int &face_index) const
Get a pointer to the function mapping face coordinates to bulk coordinates.
void local_fraction_of_node(const unsigned &j, Vector< double > &s_fraction)
Get the local fraction of node j in the element.
unsigned nvertex_node() const
Number of vertex nodes in the element.
unsigned get_bulk_node_number(const int &face_index, const unsigned &i) const
Get the number of the ith node on face face_index in the bulk node vector.
std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction)
void write_paraview_offsets(std::ofstream &file_out, const unsigned &nplot, unsigned &offset_sum) const
Return the offsets for the paraview sub-elements. Needs to be implemented for each new geometric elem...
double invert_jacobian_mapping(const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
Overload the template-free interface for the calculation of the inverse jacobian mapping....
QElement(const QElement &)=delete
Broken copy constructor.
void write_paraview_type(std::ofstream &file_out, const unsigned &nplot) const
Return the paraview element type. Needs to be implemented for each new geometric element type; see ht...
static Gauss< 3, NNODE_1D > Default_integration_scheme
Default integration rule: Gaussian integration of same 'order' as the element.
unsigned nsub_elements_paraview(const unsigned &nplot) const
Return the number of local sub-elements for paraview plot with parameter nplot.
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Get local coordinates of node j in the element; vector sets its own size.
void write_paraview_output_offset_information(std::ofstream &file_out, const unsigned &nplot, unsigned &counter) const
Fill in the offset information for paraview plot. Needs to be implemented for each new geometric elem...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
virtual void set_macro_elem_pt(MacroElement *macro_elem_pt)
Broken assignment operator.
void get_x_and_xi(const Vector< double > &s, Vector< double > &x_fe, Vector< double > &x, Vector< double > &xi_fe, Vector< double > &xi) const
Eulerian and Lagrangian coordinates as function of the local coordinates: The Eulerian position is re...
QSolidElementBase()
Constructor: Empty.
virtual void set_macro_elem_pt(MacroElement *macro_elem_pt, MacroElement *undeformed_macro_elem_pt)
Set pointers to "current" and "undeformed" MacroElements.
QSolidElementBase(const QSolidElementBase &)=delete
Broken copy constructor.
Base class for all quad elements.
QuadElementBase()
Constructor. Empty.
virtual Node * vertex_node_pt(const unsigned &j) const =0
Pointer to the j-th vertex node in the element.
virtual unsigned nvertex_node() const =0
Number of vertex nodes in the element.
A class that is used to template the refineable Q elements by dimension. It's really nothing more tha...
RefineableQElement()
Empty constuctor.
A class that is used to template the solid refineable Q elements by dimension. It's really nothing mo...
RefineableSolidQElement()
Empty constuctor.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
////////////////////////////////////////////////////////////////////// //////////////////////////////...
void set_undeformed_macro_elem_pt(MacroElement *undeformed_macro_elem_pt)
Set pointer to "undeformed" macro element. Can be overloaded in derived classes to perform additional...
MacroElement * undeformed_macro_elem_pt()
Access function to pointer to "undeformed" macro element.
MacroElement * Undeformed_macro_elem_pt
Pointer to the element's "undeformed" macro element (NULL by default)
virtual double interpolated_xi(const Vector< double > &s, const unsigned &i) const
Return i-th FE-interpolated Lagrangian coordinate xi[i] at local coordinate s.
A Class for nodes that deform elastically (i.e. position is an unknown in the problem)....
unsigned nlagrangian() const
Return number of lagrangian coordinates.
void output(FILE *file_pt)
C-style output.
SolidQElement()
Constructor.
void output(std::ostream &outfile)
Broken assignment operator.
SolidQElement(const SolidQElement &)=delete
Broken copy constructor.
void output(FILE *file_pt)
C-style output.
SolidQElement(const SolidQElement &)=delete
Broken copy constructor.
void output(std::ostream &outfile)
Broken assignment operator.
SolidQElement()
Constructor.
void output(std::ostream &outfile)
Broken assignment operator.
SolidQElement()
Constructor.
void output(FILE *file_pt)
C-style output.
SolidQElement(const SolidQElement &)=delete
Broken copy constructor.
SolidQElement elements are quadrilateral elements whose derivatives also include those based upon the...
void output()
Doc the command line arguments.
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
void shape(const double &s, double *Psi)
Definition for 1D Lagrange shape functions. The value of all the shape functions at the local coordin...
void faces0(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for both faces – the bulk coordinate is fixed on both.
void face1(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the face s0 = 1.0.
void face0(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the face s0 = -1.0.
void faces1(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for the north and south faces, along which s1 is fixed.
void faces0(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for the east and west faces, along which s0 is fixed.
void face2(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the east face (s0 = 1.0)
void face1(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the south face (s1 = -1.0)
void face0(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the west face (s0 = -1.0)
void face3(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the north face (s1 = 1.0)
void faces2(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for the left and right faces, along which s2 is fixed.
void faces0(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for the back and front faces, along which s0 is fixed.
void faces1(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for the up and down faces, along which s1 is fixed.
void face4(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the up face (s1 = 1.0)
void face1(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the down face (s1 = -1.0)
void face2(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the back face (s2 = -1.0)
void face0(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the left face (s0 = -1.0)
void face5(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the front face (s2 = 1.0)
void face3(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the right face (s0 = 1.0)
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void(* BulkCoordinateDerivativesFctPt)(const Vector< double > &s, DenseMatrix< double > &ds_bulk_dsface, unsigned &interior_direction)
Typedef for the function that returns the partial derivative of the local coordinates in the bulk ele...
void(* CoordinateMappingFctPt)(const Vector< double > &s, Vector< double > &s_bulk)
Typedef for the function that translates the face coordinate to the coordinate in the bulk element.