27 #ifndef OOMPH_TELEMENT_HEADER
28 #define OOMPH_TELEMENT_HEADER
32 #include <oomph-lib-config.h>
65 std::ostringstream error_stream;
66 error_stream <<
"TFace cannot have two identical vertex nodes\n";
68 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
73 std::set<Node*> nodes;
77 std::set<Node*>::iterator it = nodes.begin();
180 std::set<unsigned> set1;
181 std::set<unsigned>* set1_pt = &set1;
183 std::set<unsigned> set2;
184 std::set<unsigned>* set2_pt = &set2;
186 std::set<unsigned> set3;
187 std::set<unsigned>* set3_pt = &set3;
189 std::set<unsigned> aux_set;
190 set_intersection((*set1_pt).begin(),
194 inserter(aux_set, aux_set.begin()));
195 set_intersection(aux_set.begin(),
199 inserter((*boundaries_pt), (*boundaries_pt).begin()));
227 template<
unsigned DIM,
unsigned NNODE_1D>
256 std::ostringstream error_message;
258 <<
"Element only has two nodes; called with node number " << j
262 OOMPH_CURRENT_FUNCTION,
263 OOMPH_EXCEPTION_LOCATION);
300 this->dshape_local(
s, psi, dpsids);
333 std::ostringstream error_message;
335 <<
"Element only has three nodes; called with node number " << j
339 OOMPH_CURRENT_FUNCTION,
340 OOMPH_EXCEPTION_LOCATION);
350 psi[0] = 2.0 * (
s[0] - 1.0) * (
s[0] - 0.5);
351 psi[1] = 4.0 * (1.0 -
s[0]) *
s[0];
352 psi[2] = 2.0 * (
s[0] - 0.5) *
s[0];
364 dpsids(0, 0) = 4.0 *
s[0] - 3.0;
365 dpsids(1, 0) = 4.0 - 8.0 *
s[0];
366 dpsids(2, 0) = 4.0 *
s[0] - 1.0;
380 this->dshape_local(
s, psi, dpsids);
383 d2psids(1, 0) = -8.0;
417 std::ostringstream error_message;
419 <<
"Element only has four nodes; called with node number " << j
423 OOMPH_CURRENT_FUNCTION,
424 OOMPH_EXCEPTION_LOCATION);
434 psi[0] = 0.5 * (1.0 -
s[0]) * (3.0 *
s[0] - 2.0) * (3.0 *
s[0] - 1.0);
435 psi[1] = -4.5 *
s[0] * (1.0 -
s[0]) * (3.0 *
s[0] - 2.0);
436 psi[2] = 4.5 *
s[0] * (1.0 -
s[0]) * (3.0 *
s[0] - 1.0);
437 psi[3] = 0.5 *
s[0] * (3.0 *
s[0] - 2.0) * (3.0 *
s[0] - 1.0);
448 dpsids(0, 0) = -13.5 *
s[0] *
s[0] + 18.0 *
s[0] - 5.5;
449 dpsids(1, 0) = 40.5 *
s[0] *
s[0] - 45.0 *
s[0] + 9.0;
450 dpsids(2, 0) = -40.5 *
s[0] *
s[0] + 36.0 *
s[0] - 4.5;
451 dpsids(3, 0) = 13.5 *
s[0] *
s[0] - 9.0 *
s[0] + 1.0;
464 "Not checked yet", OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
467 this->dshape_local(
s, psi, dpsids);
506 std::ostringstream error_message;
508 <<
"Element only has three nodes; called with node number " << j
512 OOMPH_CURRENT_FUNCTION,
513 OOMPH_EXCEPTION_LOCATION);
525 psi[2] = 1.0 -
s[0] -
s[1];
557 this->dshape_local(
s, psi, dpsids);
559 for (
unsigned i = 0;
i < 3;
i++)
612 std::ostringstream error_message;
614 <<
"Element only has six nodes; called with node number " << j
618 OOMPH_CURRENT_FUNCTION,
619 OOMPH_EXCEPTION_LOCATION);
630 double s_2 = 1.0 -
s[0] -
s[1];
635 psi[0] = 2.0 *
s[0] * (
s[0] - 0.5);
636 psi[1] = 2.0 *
s[1] * (
s[1] - 0.5);
637 psi[2] = 2.0 * s_2 * (s_2 - 0.5);
638 psi[3] = 4.0 *
s[0] *
s[1];
639 psi[4] = 4.0 *
s[1] * s_2;
640 psi[5] = 4.0 * s_2 *
s[0];
652 dpsids(0, 0) = 4.0 *
s[0] - 1.0;
655 dpsids(1, 1) = 4.0 *
s[1] - 1.0;
656 dpsids(2, 0) = 2.0 * (2.0 *
s[0] - 1.5 + 2.0 *
s[1]);
657 dpsids(2, 1) = 2.0 * (2.0 *
s[0] - 1.5 + 2.0 *
s[1]);
658 dpsids(3, 0) = 4.0 *
s[1];
659 dpsids(3, 1) = 4.0 *
s[0];
660 dpsids(4, 0) = -4.0 *
s[1];
661 dpsids(4, 1) = 4.0 * (1.0 -
s[0] - 2.0 *
s[1]);
662 dpsids(5, 0) = 4.0 * (1.0 - 2.0 *
s[0] -
s[1]);
663 dpsids(5, 1) = -4.0 *
s[0];
679 this->dshape_local(
s, psi, dpsids);
698 d2psids(4, 1) = -8.0;
699 d2psids(4, 2) = -4.0;
701 d2psids(5, 0) = -8.0;
703 d2psids(5, 2) = -4.0;
771 std::ostringstream error_message;
773 <<
"Element only has ten nodes; called with node number " << j
777 OOMPH_CURRENT_FUNCTION,
778 OOMPH_EXCEPTION_LOCATION);
788 psi[0] = 0.5 *
s[0] * (3.0 *
s[0] - 2.0) * (3.0 *
s[0] - 1.0);
789 psi[1] = 0.5 *
s[1] * (3.0 *
s[1] - 2.0) * (3.0 *
s[1] - 1.0);
790 psi[2] = 0.5 * (1.0 -
s[0] -
s[1]) * (1.0 - 3.0 *
s[0] - 3.0 *
s[1]) *
791 (2.0 - 3.0 *
s[0] - 3.0 *
s[1]);
792 psi[3] = 4.5 *
s[0] *
s[1] * (3.0 *
s[0] - 1.0);
793 psi[4] = 4.5 *
s[0] *
s[1] * (3.0 *
s[1] - 1.0);
794 psi[5] = 4.5 *
s[1] * (1.0 -
s[0] -
s[1]) * (3.0 *
s[1] - 1.0);
796 4.5 *
s[1] * (1.0 -
s[0] -
s[1]) * (3.0 * (1.0 -
s[0] -
s[1]) - 1.0);
797 psi[7] = 4.5 *
s[0] * (1.0 -
s[0] -
s[1]) * (2.0 - 3 *
s[0] - 3 *
s[1]);
798 psi[8] = 4.5 *
s[0] * (1.0 -
s[0] -
s[1]) * (3.0 *
s[0] - 1.0);
799 psi[9] = 27.0 *
s[0] *
s[1] * (1.0 -
s[0] -
s[1]);
810 dpsids(0, 0) = 13.5 *
s[0] *
s[0] - 9.0 *
s[0] + 1.0;
813 dpsids(1, 1) = 13.5 *
s[1] *
s[1] - 9.0 *
s[1] + 1.0;
814 dpsids(2, 0) = 0.5 * (36.0 *
s[0] + 36.0 *
s[1] - 27.0 *
s[0] *
s[0] -
815 27.0 *
s[1] *
s[1] - 54.0 *
s[0] *
s[1] - 11.0);
816 dpsids(2, 1) = 0.5 * (36.0 *
s[0] + 36.0 *
s[1] - 27.0 *
s[0] *
s[0] -
817 27.0 *
s[1] *
s[1] - 54.0 *
s[0] *
s[1] - 11.0);
818 dpsids(3, 0) = 27.0 *
s[0] *
s[1] - 4.5 *
s[1];
819 dpsids(3, 1) = 4.5 *
s[0] * (3.0 *
s[0] - 1.0);
820 dpsids(4, 0) = 4.5 *
s[1] * (3.0 *
s[1] - 1.0);
821 dpsids(4, 1) = 27.0 *
s[0] *
s[1] - 4.5 *
s[0];
822 dpsids(5, 0) = 4.5 * (
s[1] - 3.0 *
s[1] *
s[1]);
824 4.5 * (
s[0] - 6.0 *
s[0] *
s[1] - 9.0 *
s[1] *
s[1] + 8 *
s[1] - 1.0);
825 dpsids(6, 0) = 4.5 * (6.0 *
s[0] *
s[1] - 5.0 *
s[1] + 6.0 *
s[1] *
s[1]);
827 4.5 * (2.0 - 5.0 *
s[0] + 3.0 *
s[0] *
s[0] + 12.0 *
s[0] *
s[1] -
828 10.0 *
s[1] + 9.0 *
s[1] *
s[1]);
830 4.5 * (2.0 - 10.0 *
s[0] + 9.0 *
s[0] *
s[0] + 12.0 *
s[0] *
s[1] -
831 5.0 *
s[1] + 3.0 *
s[1] *
s[1]);
832 dpsids(7, 1) = 4.5 * (6.0 *
s[0] *
s[0] - 5.0 *
s[0] + 6.0 *
s[0] *
s[1]);
834 4.5 * (
s[1] - 6.0 *
s[0] *
s[1] - 9.0 *
s[0] *
s[0] + 8 *
s[0] - 1.0);
835 dpsids(8, 1) = 4.5 * (
s[0] - 3.0 *
s[0] *
s[0]);
836 dpsids(9, 0) = 27.0 *
s[1] - 54.0 *
s[0] *
s[1] - 27.0 *
s[1] *
s[1];
837 dpsids(9, 1) = 27.0 *
s[0] - 54.0 *
s[0] *
s[1] - 27.0 *
s[0] *
s[0];
852 "Not checked yet", OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
855 this->dshape_local(
s, psi, dpsids);
857 d2psids(0, 0) = 9.0 * (3.0 *
s[0] - 1.0);
862 d2psids(1, 2) = 9.0 * (3.0 *
s[1] - 1.0);
863 d2psids(2, 0) = 9.0 * (2.0 - 3.0 *
s[0] - 3.0 *
s[1]);
864 d2psids(2, 1) = 9.0 * (2.0 - 3.0 *
s[0] - 3.0 *
s[1]);
865 d2psids(2, 2) = 9.0 * (2.0 - 3.0 *
s[0] - 3.0 *
s[1]);
866 d2psids(3, 0) = 27.0 *
s[1];
868 d2psids(3, 2) = 27.0 *
s[0] - 4.5;
870 d2psids(4, 1) = 27.0 *
s[0];
871 d2psids(4, 2) = 27.0 *
s[1] - 4.5;
873 d2psids(5, 1) = 9.0 * (4.0 - 3.0 *
s[0] - 9.0 *
s[1]);
874 d2psids(5, 2) = 4.5 * (1.0 - 6.0 *
s[1]);
875 d2psids(6, 0) = 27.0 *
s[1];
876 d2psids(6, 1) = 9.0 * (6.0 *
s[0] + 9.0 *
s[1] - 5.0);
877 d2psids(6, 2) = 4.5 * (6.0 *
s[0] + 12.0 *
s[1] - 5.0);
878 d2psids(8, 0) = 9.0 * (4.0 - 9.0 *
s[0] - 3.0 *
s[1]);
880 d2psids(8, 2) = 4.5 * (1.0 - 6.0 *
s[0]);
881 d2psids(7, 0) = 9.0 * (9.0 *
s[0] + 6.0 *
s[1] - 5.0);
882 d2psids(7, 1) = 27.0 *
s[0];
883 d2psids(7, 2) = 4.5 * (12.0 *
s[0] + 6.0 *
s[1] - 5.0);
884 d2psids(9, 0) = -54.0 *
s[1];
885 d2psids(9, 1) = -54.0 *
s[0];
886 d2psids(9, 2) = 27.0 - 54.0 *
s[0] - 54.0 *
s[1];
899 template<
unsigned DIM,
unsigned NNODE_1D>
972 std::ostringstream error_message;
974 <<
"Element only has seven nodes; called with node number " << j
978 OOMPH_CURRENT_FUNCTION,
979 OOMPH_EXCEPTION_LOCATION);
990 const double s_2 = 1.0 -
s[0] -
s[1];
993 const double cubic_bubble =
s[0] *
s[1] * s_2;
1001 psi[0] = 2.0 *
s[0] * (
s[0] - 0.5) + 3.0 * cubic_bubble;
1002 psi[1] = 2.0 *
s[1] * (
s[1] - 0.5) + 3.0 * cubic_bubble;
1003 psi[2] = 2.0 * s_2 * (s_2 - 0.5) + 3.0 * cubic_bubble;
1004 psi[3] = 4.0 *
s[0] *
s[1] - 12.0 * cubic_bubble;
1005 psi[4] = 4.0 *
s[1] * s_2 - 12.0 * cubic_bubble;
1006 psi[5] = 4.0 * s_2 *
s[0] - 12.0 * cubic_bubble;
1008 psi[6] = 27.0 * cubic_bubble;
1018 this->
shape(s, psi);
1021 const double d_bubble_ds0 =
s[1] * (1.0 -
s[1] - 2.0 *
s[0]);
1022 const double d_bubble_ds1 =
s[0] * (1.0 -
s[0] - 2.0 *
s[1]);
1025 dpsids(0, 0) = 4.0 *
s[0] - 1.0 + 3.0 * d_bubble_ds0;
1026 dpsids(0, 1) = 0.0 + 3.0 * d_bubble_ds1;
1027 dpsids(1, 0) = 0.0 + 3.0 * d_bubble_ds0;
1028 dpsids(1, 1) = 4.0 *
s[1] - 1.0 + 3.0 * d_bubble_ds1;
1029 dpsids(2, 0) = 2.0 * (2.0 *
s[0] - 1.5 + 2.0 *
s[1]) + 3.0 * d_bubble_ds0;
1030 dpsids(2, 1) = 2.0 * (2.0 *
s[0] - 1.5 + 2.0 *
s[1]) + 3.0 * d_bubble_ds1;
1031 dpsids(3, 0) = 4.0 *
s[1] - 12.0 * d_bubble_ds0;
1032 dpsids(3, 1) = 4.0 *
s[0] - 12.0 * d_bubble_ds1;
1033 dpsids(4, 0) = -4.0 *
s[1] - 12.0 * d_bubble_ds0;
1034 dpsids(4, 1) = 4.0 * (1.0 -
s[0] - 2.0 *
s[1]) - 12.0 * d_bubble_ds1;
1035 dpsids(5, 0) = 4.0 * (1.0 - 2.0 *
s[0] -
s[1]) - 12.0 * d_bubble_ds0;
1036 dpsids(5, 1) = -4.0 *
s[0] - 12.0 * d_bubble_ds1;
1037 dpsids(6, 0) = 27.0 * d_bubble_ds0;
1038 dpsids(6, 1) = 27.0 * d_bubble_ds1;
1055 this->dshape_local(
s, psi, dpsids);
1058 const double d2_bubble_ds0 = -2.0 *
s[1];
1059 const double d2_bubble_ds1 = -2.0 *
s[0];
1060 const double d2_bubble_ds2 = 1.0 - 2.0 *
s[0] - 2.0 *
s[1];
1062 d2psids(0, 0) = 4.0 + 3.0 * d2_bubble_ds0;
1063 d2psids(0, 1) = 0.0 + 3.0 * d2_bubble_ds1;
1064 d2psids(0, 2) = 0.0 + 3.0 * d2_bubble_ds2;
1066 d2psids(1, 0) = 0.0 + 3.0 * d2_bubble_ds0;
1067 d2psids(1, 1) = 4.0 + 3.0 * d2_bubble_ds1;
1068 d2psids(1, 2) = 0.0 + 3.0 * d2_bubble_ds2;
1070 d2psids(2, 0) = 4.0 + 3.0 * d2_bubble_ds0;
1071 d2psids(2, 1) = 4.0 + 3.0 * d2_bubble_ds1;
1072 d2psids(2, 2) = 4.0 + 3.0 * d2_bubble_ds2;
1074 d2psids(3, 0) = 0.0 - 12.0 * d2_bubble_ds0;
1075 d2psids(3, 1) = 0.0 - 12.0 * d2_bubble_ds1;
1076 d2psids(3, 2) = 4.0 - 12.0 * d2_bubble_ds2;
1078 d2psids(4, 0) = 0.0 - 12.0 * d2_bubble_ds0;
1079 d2psids(4, 1) = -8.0 - 12.0 * d2_bubble_ds1;
1080 d2psids(4, 2) = -4.0 - 12.0 * d2_bubble_ds2;
1082 d2psids(5, 0) = -8.0 - 12.0 * d2_bubble_ds0;
1083 d2psids(5, 1) = 0.0 - 12.0 * d2_bubble_ds1;
1084 d2psids(5, 2) = -4.0 - 12.0 * d2_bubble_ds2;
1086 d2psids(6, 0) = 27.0 * d2_bubble_ds0;
1087 d2psids(6, 1) = 27.0 * d2_bubble_ds1;
1088 d2psids(6, 2) = 27.0 * d2_bubble_ds2;
1151 unsigned ncoord =
dim();
1153 for (
unsigned i = 0;
i < ncoord;
i++)
1178 unsigned ncoord =
dim();
1180 for (
unsigned i = 0;
i < ncoord;
i++)
1183 if (
s[
i] < 0.0)
s[
i] = 0.0;
1188 double excess = sum - 1.0;
1192 double sub = excess / double(ncoord);
1193 for (
unsigned i = 0;
i < ncoord;
i++)
1206 template<
unsigned DIM,
unsigned NNODE_1D>
1218 template<
unsigned NNODE_1D>
1244 "One-dimensional TElements are currently only implemented for\n";
1245 error_message +=
"three and six nodes, i.e. NNODE_1D=2 , 3 , 4\n";
1248 error_message, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1252 unsigned n_node = NNODE_1D;
1253 this->set_n_node(n_node);
1259 set_integration_scheme(&Default_integration_scheme);
1299 return node_pt(NNODE_1D - 1);
1304 std::ostringstream error_message;
1306 <<
"Element only has two vertex nodes; called with node number "
1309 OOMPH_CURRENT_FUNCTION,
1310 OOMPH_EXCEPTION_LOCATION);
1348 return FiniteElement::invert_jacobian<1>(jacobian, inverse_jacobian);
1388 const unsigned& nplot,
1389 unsigned& counter)
const
1392 unsigned plot = nsub_elements_paraview(nplot);
1395 for (
unsigned i = 0;
i < plot;
i++)
1397 file_out <<
i + counter <<
" " <<
i + 1 + counter << std::endl;
1399 counter += nplot_points_paraview(nplot);
1406 const unsigned& nplot)
const
1408 unsigned local_loop = nsub_elements_paraview(nplot);
1409 for (
unsigned i = 0;
i < local_loop;
i++)
1411 file_out <<
"3" << std::endl;
1419 const unsigned& nplot,
1420 unsigned& offset_sum)
const
1424 unsigned local_loop = nsub_elements_paraview(nplot);
1425 for (
unsigned i = 0;
i < local_loop;
i++)
1428 file_out << offset_sum << std::endl;
1436 void output(std::ostream& outfile,
const unsigned& nplot);
1439 void output(FILE* file_pt);
1442 void output(FILE* file_pt,
const unsigned& n_plot);
1448 const unsigned& nplot,
1450 const bool& use_equally_spaced_interior_sample_points =
false)
const
1454 s[0] = double(
i) / double(nplot - 1);
1456 if (use_equally_spaced_interior_sample_points)
1459 double dx_new = range / double(nplot);
1460 double range_new = double(nplot - 1) * dx_new;
1461 s[0] = 0.5 * dx_new + range_new *
s[0] / range;
1474 std::ostringstream header;
1475 header <<
"ZONE I=" << nplot <<
"\n";
1476 return header.str();
1491 void build_face_element(
const int& face_index,
1502 template<
unsigned NNODE_1D>
1532 "Triangles are currently only implemented for\n";
1533 error_message +=
"three and six nodes, i.e. NNODE_1D=2 , 3 , 4\n";
1536 error_message, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1540 unsigned n_node = (NNODE_1D * (NNODE_1D + 1)) / 2;
1541 this->set_n_node(n_node);
1547 set_integration_scheme(&Default_integration_scheme);
1554 if (!allow_high_order)
1562 unsigned n_node = (NNODE_1D * (NNODE_1D + 1)) / 2;
1563 this->set_n_node(n_node);
1569 set_integration_scheme(&Default_integration_scheme);
1599 const unsigned&
i)
const
1611 std::ostringstream error_message;
1613 <<
"Element only has three vertex nodes; called with node number "
1616 OOMPH_CURRENT_FUNCTION,
1617 OOMPH_EXCEPTION_LOCATION);
1657 return FiniteElement::invert_jacobian<2>(jacobian, inverse_jacobian);
1683 unsigned node_sum = 0;
1684 for (
unsigned i = 1;
i <= nplot;
i++)
1695 unsigned local_sum = 0;
1696 for (
unsigned i = 1;
i < nplot;
i++)
1698 local_sum += 2 * (nplot -
i - 1) + 1;
1707 const unsigned& nplot,
1708 unsigned& counter)
const
1713 unsigned node_count = 0;
1714 for (
unsigned i = 0;
i < nplot - 1;
i++)
1716 for (
unsigned j = 0; j < nplot -
i - 1; j++)
1718 file_out << j + node_count + counter <<
" "
1719 << j + node_count + 1 + counter <<
" "
1720 << j + nplot + node_count -
i + counter << std::endl;
1722 if (j < nplot -
i - 2)
1724 file_out << j + node_count + 1 + counter <<
" "
1725 << j + nplot + node_count -
i + 1 + counter <<
" "
1726 << j + nplot + node_count -
i + counter << std::endl;
1729 node_count += (nplot -
i);
1731 counter += nplot_points_paraview(nplot);
1738 const unsigned& nplot)
const
1740 unsigned local_loop = nsub_elements_paraview(nplot);
1743 for (
unsigned i = 0;
i < local_loop;
i++)
1745 file_out <<
"5" << std::endl;
1753 const unsigned& nplot,
1754 unsigned& offset_sum)
const
1756 unsigned local_loop = nsub_elements_paraview(nplot);
1760 for (
unsigned i = 0;
i < local_loop;
i++)
1763 file_out << offset_sum << std::endl;
1771 void output(std::ostream& outfile,
const unsigned& nplot);
1774 void output(FILE* file_pt);
1777 void output(FILE* file_pt,
const unsigned& n_plot);
1782 const unsigned& iplot,
1783 const unsigned& nplot,
1785 const bool& use_equally_spaced_interior_sample_points =
false)
const
1789 unsigned np = 0,
i, j;
1790 for (
i = 0;
i < nplot; ++
i)
1792 for (j = 0; j < nplot -
i; ++j)
1796 s[0] = double(j) / double(nplot - 1);
1797 s[1] = double(
i) / double(nplot - 1);
1798 if (use_equally_spaced_interior_sample_points)
1801 double dx_new = range / (double(nplot) + 0.5);
1802 double range_new = double(nplot - 1) * dx_new;
1803 s[0] = 0.5 * dx_new + range_new *
s[0] / range;
1804 s[1] = 0.5 * dx_new + range_new *
s[1] / range;
1823 std::ostringstream header;
1825 for (
unsigned i = 1;
i < nplot;
i++)
1829 header <<
"ZONE N=" << nplot_points(nplot) <<
", E=" << nel
1830 <<
", F=FEPOINT, ET=TRIANGLE\n";
1831 return header.str();
1839 const unsigned& nplot)
const
1843 unsigned nod_count = 1;
1844 for (
unsigned i = 0;
i < nplot;
i++)
1846 for (
unsigned j = 0; j < nplot -
i; j++)
1848 if (j < nplot -
i - 1)
1850 outfile << nod_count <<
" " << nod_count + 1 <<
" "
1851 << nod_count + nplot -
i << std::endl;
1852 if (j < nplot -
i - 2)
1854 outfile << nod_count + 1 <<
" " << nod_count + nplot -
i + 1
1855 <<
" " << nod_count + nplot -
i << std::endl;
1871 unsigned nod_count = 1;
1872 for (
unsigned i = 0;
i < nplot;
i++)
1874 for (
unsigned j = 0; j < nplot -
i; j++)
1876 if (j < nplot -
i - 1)
1882 nod_count + nplot -
i);
1883 if (j < nplot -
i - 2)
1888 nod_count + nplot -
i + 1,
1889 nod_count + nplot -
i);
1902 for (
unsigned i = 1;
i <= nplot;
i++)
1915 void build_face_element(
const int& face_index,
1963 std::ostringstream error_message;
1965 <<
"Element only has four nodes; called with node number " << j
1969 OOMPH_CURRENT_FUNCTION,
1970 OOMPH_EXCEPTION_LOCATION);
1983 psi[3] = 1.0 -
s[0] -
s[1] -
s[2];
1993 this->
shape(s, psi);
2008 dpsids(3, 0) = -1.0;
2009 dpsids(3, 1) = -1.0;
2010 dpsids(3, 2) = -1.0;
2029 this->dshape_local(
s, psi, dpsids);
2031 for (
unsigned i = 0;
i < 4;
i++)
2033 d2psids(
i, 0) = 0.0;
2034 d2psids(
i, 1) = 0.0;
2035 d2psids(
i, 2) = 0.0;
2036 d2psids(
i, 3) = 0.0;
2037 d2psids(
i, 4) = 0.0;
2038 d2psids(
i, 5) = 0.0;
2118 std::ostringstream error_message;
2120 <<
"Element only has ten nodes; called with node number " << j
2124 OOMPH_CURRENT_FUNCTION,
2125 OOMPH_EXCEPTION_LOCATION);
2135 double s3 = 1.0 -
s[0] -
s[1] -
s[2];
2136 psi[0] = (2.0 *
s[0] - 1.0) *
s[0];
2137 psi[1] = (2.0 *
s[1] - 1.0) *
s[1];
2138 psi[2] = (2.0 *
s[2] - 1.0) *
s[2];
2139 psi[3] = (2.0 * s3 - 1.0) * s3;
2140 psi[4] = 4.0 *
s[0] *
s[1];
2141 psi[5] = 4.0 *
s[0] *
s[2];
2142 psi[6] = 4.0 *
s[0] * s3;
2143 psi[7] = 4.0 *
s[1] *
s[2];
2144 psi[8] = 4.0 *
s[2] * s3;
2145 psi[9] = 4.0 *
s[1] * s3;
2155 this->
shape(s, psi);
2158 double s3 = 1.0 -
s[0] -
s[1] -
s[2];
2160 dpsids(0, 0) = 4.0 *
s[0] - 1.0;
2165 dpsids(1, 1) = 4.0 *
s[1] - 1.0;
2170 dpsids(2, 2) = 4.0 *
s[2] - 1.0;
2172 dpsids(3, 0) = -4.0 * s3 + 1.0;
2173 dpsids(3, 1) = -4.0 * s3 + 1.0;
2174 dpsids(3, 2) = -4.0 * s3 + 1.0;
2176 dpsids(4, 0) = 4.0 *
s[1];
2177 dpsids(4, 1) = 4.0 *
s[0];
2180 dpsids(5, 0) = 4.0 *
s[2];
2182 dpsids(5, 2) = 4.0 *
s[0];
2184 dpsids(6, 0) = 4.0 * (s3 -
s[0]);
2185 dpsids(6, 1) = -4.0 *
s[0];
2186 dpsids(6, 2) = -4.0 *
s[0];
2189 dpsids(7, 1) = 4.0 *
s[2];
2190 dpsids(7, 2) = 4.0 *
s[1];
2192 dpsids(8, 0) = -4.0 *
s[2];
2193 dpsids(8, 1) = -4.0 *
s[2];
2194 dpsids(8, 2) = 4.0 * (s3 -
s[2]);
2196 dpsids(9, 0) = -4.0 *
s[1];
2197 dpsids(9, 1) = 4.0 * (s3 -
s[1]);
2198 dpsids(9, 2) = -4.0 *
s[1];
2217 this->dshape_local(
s, psi, dpsids);
2223 d2psids(0, 0) = 4.0;
2224 d2psids(0, 1) = 0.0;
2225 d2psids(0, 2) = 0.0;
2226 d2psids(0, 3) = 0.0;
2227 d2psids(0, 4) = 0.0;
2228 d2psids(0, 5) = 0.0;
2231 d2psids(1, 0) = 0.0;
2232 d2psids(1, 1) = 4.0;
2233 d2psids(1, 2) = 0.0;
2234 d2psids(1, 3) = 0.0;
2235 d2psids(1, 4) = 0.0;
2236 d2psids(1, 5) = 0.0;
2238 d2psids(2, 0) = 0.0;
2239 d2psids(2, 1) = 0.0;
2240 d2psids(2, 2) = 4.0;
2241 d2psids(2, 3) = 0.0;
2242 d2psids(2, 4) = 0.0;
2243 d2psids(2, 5) = 0.0;
2245 d2psids(3, 0) = 4.0;
2246 d2psids(3, 1) = 4.0;
2247 d2psids(3, 2) = 4.0;
2248 d2psids(3, 3) = 4.0;
2249 d2psids(3, 4) = 4.0;
2250 d2psids(3, 5) = 4.0;
2252 d2psids(4, 0) = 0.0;
2253 d2psids(4, 1) = 0.0;
2254 d2psids(4, 2) = 0.0;
2255 d2psids(4, 3) = 4.0;
2256 d2psids(4, 4) = 0.0;
2257 d2psids(4, 5) = 0.0;
2259 d2psids(5, 0) = 0.0;
2260 d2psids(5, 1) = 0.0;
2261 d2psids(5, 2) = 0.0;
2262 d2psids(5, 3) = 0.0;
2263 d2psids(5, 4) = 4.0;
2264 d2psids(5, 5) = 0.0;
2266 d2psids(6, 0) = -8.0;
2267 d2psids(6, 1) = 0.0;
2268 d2psids(6, 2) = 0.0;
2269 d2psids(6, 3) = -4.0;
2270 d2psids(6, 4) = -4.0;
2271 d2psids(6, 5) = 0.0;
2273 d2psids(7, 0) = 0.0;
2274 d2psids(7, 1) = 0.0;
2275 d2psids(7, 2) = 0.0;
2276 d2psids(7, 3) = 0.0;
2277 d2psids(7, 4) = 0.0;
2278 d2psids(7, 5) = 4.0;
2280 d2psids(8, 0) = 0.0;
2281 d2psids(8, 1) = 0.0;
2282 d2psids(8, 2) = -8.0;
2283 d2psids(8, 3) = 0.0;
2284 d2psids(8, 4) = -4.0;
2285 d2psids(8, 5) = -4.0;
2287 d2psids(9, 0) = 0.0;
2288 d2psids(9, 1) = -8.0;
2289 d2psids(9, 2) = 0.0;
2290 d2psids(9, 3) = -4.0;
2291 d2psids(9, 4) = 0.0;
2292 d2psids(9, 5) = -4.0;
2423 std::ostringstream error_message;
2425 <<
"Element only has fifteen nodes; called with node number " << j
2429 OOMPH_CURRENT_FUNCTION,
2430 OOMPH_EXCEPTION_LOCATION);
2441 const double s3 = 1.0 -
s[0] -
s[1] -
s[2];
2443 const double quartic_bubble =
s[0] *
s[1] *
s[2] * s3;
2444 const double cubic_bubble012 =
s[0] *
s[1] *
s[2];
2445 const double cubic_bubble013 =
s[0] *
s[1] * s3;
2446 const double cubic_bubble023 =
s[0] *
s[2] * s3;
2447 const double cubic_bubble123 =
s[1] *
s[2] * s3;
2454 psi[0] = (2.0 *
s[0] - 1.0) *
s[0] +
2455 3.0 * (cubic_bubble012 + cubic_bubble013 + cubic_bubble023) -
2456 4.0 * quartic_bubble;
2457 psi[1] = (2.0 *
s[1] - 1.0) *
s[1] +
2458 3.0 * (cubic_bubble012 + cubic_bubble013 + cubic_bubble123) -
2459 4.0 * quartic_bubble;
2460 psi[2] = (2.0 *
s[2] - 1.0) *
s[2] +
2461 3.0 * (cubic_bubble012 + cubic_bubble023 + cubic_bubble123) -
2462 4.0 * quartic_bubble;
2463 psi[3] = (2.0 * s3 - 1.0) * s3 +
2464 3.0 * (cubic_bubble013 + cubic_bubble023 + cubic_bubble123) -
2465 4.0 * quartic_bubble;
2466 psi[4] = 4.0 *
s[0] *
s[1] - 12.0 * (cubic_bubble012 + cubic_bubble013) +
2467 32.0 * quartic_bubble;
2468 psi[5] = 4.0 *
s[0] *
s[2] - 12.0 * (cubic_bubble012 + cubic_bubble023) +
2469 32.0 * quartic_bubble;
2470 psi[6] = 4.0 *
s[0] * s3 - 12.0 * (cubic_bubble013 + cubic_bubble023) +
2471 32.0 * quartic_bubble;
2472 psi[7] = 4.0 *
s[1] *
s[2] - 12.0 * (cubic_bubble012 + cubic_bubble123) +
2473 32.0 * quartic_bubble;
2474 psi[8] = 4.0 *
s[2] * s3 - 12.0 * (cubic_bubble023 + cubic_bubble123) +
2475 32.0 * quartic_bubble;
2476 psi[9] = 4.0 *
s[1] * s3 - 12.0 * (cubic_bubble013 + cubic_bubble123) +
2477 32.0 * quartic_bubble;
2479 psi[10] = 27.0 * cubic_bubble013 - 108.0 * quartic_bubble;
2481 psi[11] = 27.0 * cubic_bubble012 - 108.0 * quartic_bubble;
2483 psi[12] = 27.0 * cubic_bubble023 - 108.0 * quartic_bubble;
2485 psi[13] = 27.0 * cubic_bubble123 - 108.0 * quartic_bubble;
2487 psi[14] = 256.0 * quartic_bubble;
2497 this->
shape(s, psi);
2500 const double s3 = 1.0 -
s[0] -
s[1] -
s[2];
2503 const double d_quartic_bubble_ds0 =
2504 s[1] *
s[2] * (1.0 -
s[1] -
s[2] - 2.0 *
s[0]);
2505 const double d_quartic_bubble_ds1 =
2506 s[0] *
s[2] * (1.0 -
s[0] -
s[2] - 2.0 *
s[1]);
2507 const double d_quartic_bubble_ds2 =
2508 s[0] *
s[1] * (1.0 -
s[0] -
s[1] - 2.0 *
s[2]);
2510 const double d_cubic_bubble012_ds0 =
s[1] *
s[2];
2511 const double d_cubic_bubble012_ds1 =
s[0] *
s[2];
2512 const double d_cubic_bubble012_ds2 =
s[0] *
s[1];
2514 const double d_cubic_bubble013_ds0 =
s[1] * (s3 -
s[0]);
2515 const double d_cubic_bubble013_ds1 =
s[0] * (s3 -
s[1]);
2516 const double d_cubic_bubble013_ds2 = -
s[0] *
s[1];
2518 const double d_cubic_bubble023_ds0 =
s[2] * (s3 -
s[0]);
2519 const double d_cubic_bubble023_ds1 = -
s[0] *
s[2];
2520 const double d_cubic_bubble023_ds2 =
s[0] * (s3 -
s[2]);
2522 const double d_cubic_bubble123_ds0 = -
s[1] *
s[2];
2523 const double d_cubic_bubble123_ds1 =
s[2] * (s3 -
s[1]);
2524 const double d_cubic_bubble123_ds2 =
s[1] * (s3 -
s[2]);
2529 dpsids(0, 0) = 4.0 *
s[0] - 1.0 +
2530 3.0 * (d_cubic_bubble012_ds0 + d_cubic_bubble013_ds0 +
2531 d_cubic_bubble023_ds0) -
2532 4.0 * d_quartic_bubble_ds0;
2533 dpsids(0, 1) = 0.0 +
2534 3.0 * (d_cubic_bubble012_ds1 + d_cubic_bubble013_ds1 +
2535 d_cubic_bubble023_ds1) -
2536 4.0 * d_quartic_bubble_ds1;
2537 dpsids(0, 2) = 0.0 +
2538 3.0 * (d_cubic_bubble012_ds2 + d_cubic_bubble013_ds2 +
2539 d_cubic_bubble023_ds2) -
2540 4.0 * d_quartic_bubble_ds2;
2542 dpsids(1, 0) = 0.0 +
2543 3.0 * (d_cubic_bubble012_ds0 + d_cubic_bubble013_ds0 +
2544 d_cubic_bubble123_ds0) -
2545 4.0 * d_quartic_bubble_ds0;
2546 dpsids(1, 1) = 4.0 *
s[1] - 1.0 +
2547 3.0 * (d_cubic_bubble012_ds1 + d_cubic_bubble013_ds1 +
2548 d_cubic_bubble123_ds1) -
2549 4.0 * d_quartic_bubble_ds1;
2550 dpsids(1, 2) = 0.0 +
2551 3.0 * (d_cubic_bubble012_ds2 + d_cubic_bubble013_ds2 +
2552 d_cubic_bubble123_ds2) -
2553 4.0 * d_quartic_bubble_ds2;
2555 dpsids(2, 0) = 0.0 +
2556 3.0 * (d_cubic_bubble012_ds0 + d_cubic_bubble023_ds0 +
2557 d_cubic_bubble123_ds0) -
2558 4.0 * d_quartic_bubble_ds0;
2559 dpsids(2, 1) = 0.0 +
2560 3.0 * (d_cubic_bubble012_ds1 + d_cubic_bubble023_ds1 +
2561 d_cubic_bubble123_ds1) -
2562 4.0 * d_quartic_bubble_ds1;
2563 dpsids(2, 2) = 4.0 *
s[2] - 1.0 +
2564 3.0 * (d_cubic_bubble012_ds2 + d_cubic_bubble023_ds2 +
2565 d_cubic_bubble123_ds2) -
2566 4.0 * d_quartic_bubble_ds2;
2568 dpsids(3, 0) = -4.0 * s3 + 1.0 +
2569 3.0 * (d_cubic_bubble013_ds0 + d_cubic_bubble023_ds0 +
2570 d_cubic_bubble123_ds0) -
2571 4.0 * d_quartic_bubble_ds0;
2572 dpsids(3, 1) = -4.0 * s3 + 1.0 +
2573 3.0 * (d_cubic_bubble013_ds1 + d_cubic_bubble023_ds1 +
2574 d_cubic_bubble123_ds1) -
2575 4.0 * d_quartic_bubble_ds1;
2576 dpsids(3, 2) = -4.0 * s3 + 1.0 +
2577 3.0 * (d_cubic_bubble013_ds2 + d_cubic_bubble023_ds2 +
2578 d_cubic_bubble123_ds2) -
2579 4.0 * d_quartic_bubble_ds2;
2581 dpsids(4, 0) = 4.0 *
s[1] -
2582 12.0 * (d_cubic_bubble012_ds0 + d_cubic_bubble013_ds0) +
2583 32.0 * d_quartic_bubble_ds0;
2584 dpsids(4, 1) = 4.0 *
s[0] -
2585 12.0 * (d_cubic_bubble012_ds1 + d_cubic_bubble013_ds1) +
2586 32.0 * d_quartic_bubble_ds1;
2587 dpsids(4, 2) = 0.0 -
2588 12.0 * (d_cubic_bubble012_ds2 + d_cubic_bubble013_ds2) +
2589 32.0 * d_quartic_bubble_ds2;
2591 dpsids(5, 0) = 4.0 *
s[2] -
2592 12.0 * (d_cubic_bubble012_ds0 + d_cubic_bubble023_ds0) +
2593 32.0 * d_quartic_bubble_ds0;
2594 dpsids(5, 1) = 0.0 -
2595 12.0 * (d_cubic_bubble012_ds1 + d_cubic_bubble023_ds1) +
2596 32.0 * d_quartic_bubble_ds1;
2597 dpsids(5, 2) = 4.0 *
s[0] -
2598 12.0 * (d_cubic_bubble012_ds2 + d_cubic_bubble023_ds2) +
2599 32.0 * d_quartic_bubble_ds2;
2601 dpsids(6, 0) = 4.0 * (s3 -
s[0]) -
2602 12.0 * (d_cubic_bubble013_ds0 + d_cubic_bubble023_ds0) +
2603 32.0 * d_quartic_bubble_ds0;
2604 dpsids(6, 1) = -4.0 *
s[0] -
2605 12.0 * (d_cubic_bubble013_ds1 + d_cubic_bubble023_ds1) +
2606 32.0 * d_quartic_bubble_ds1;
2607 dpsids(6, 2) = -4.0 *
s[0] -
2608 12.0 * (d_cubic_bubble013_ds2 + d_cubic_bubble023_ds2) +
2609 32.0 * d_quartic_bubble_ds2;
2611 dpsids(7, 0) = 0.0 -
2612 12.0 * (d_cubic_bubble012_ds0 + d_cubic_bubble123_ds0) +
2613 32.0 * d_quartic_bubble_ds0;
2614 dpsids(7, 1) = 4.0 *
s[2] -
2615 12.0 * (d_cubic_bubble012_ds1 + d_cubic_bubble123_ds1) +
2616 32.0 * d_quartic_bubble_ds1;
2617 dpsids(7, 2) = 4.0 *
s[1] -
2618 12.0 * (d_cubic_bubble012_ds2 + d_cubic_bubble123_ds2) +
2619 32.0 * d_quartic_bubble_ds2;
2621 dpsids(8, 0) = -4.0 *
s[2] -
2622 12.0 * (d_cubic_bubble023_ds0 + d_cubic_bubble123_ds0) +
2623 32.0 * d_quartic_bubble_ds0;
2624 dpsids(8, 1) = -4.0 *
s[2] -
2625 12.0 * (d_cubic_bubble023_ds1 + d_cubic_bubble123_ds1) +
2626 32.0 * d_quartic_bubble_ds1;
2627 dpsids(8, 2) = 4.0 * (s3 -
s[2]) -
2628 12.0 * (d_cubic_bubble023_ds2 + d_cubic_bubble123_ds2) +
2629 32.0 * d_quartic_bubble_ds2;
2631 dpsids(9, 0) = -4.0 *
s[1] -
2632 12.0 * (d_cubic_bubble013_ds0 + d_cubic_bubble123_ds0) +
2633 32.0 * d_quartic_bubble_ds0;
2634 dpsids(9, 1) = 4.0 * (s3 -
s[1]) -
2635 12.0 * (d_cubic_bubble013_ds1 + d_cubic_bubble123_ds1) +
2636 32.0 * d_quartic_bubble_ds1;
2637 dpsids(9, 2) = -4.0 *
s[1] -
2638 12.0 * (d_cubic_bubble013_ds2 + d_cubic_bubble123_ds2) +
2639 32.0 * d_quartic_bubble_ds2;
2643 27.0 * d_cubic_bubble013_ds0 - 108.0 * d_quartic_bubble_ds0;
2645 27.0 * d_cubic_bubble013_ds1 - 108.0 * d_quartic_bubble_ds1;
2647 27.0 * d_cubic_bubble013_ds2 - 108.0 * d_quartic_bubble_ds2;
2651 27.0 * d_cubic_bubble012_ds0 - 108.0 * d_quartic_bubble_ds0;
2653 27.0 * d_cubic_bubble012_ds1 - 108.0 * d_quartic_bubble_ds1;
2655 27.0 * d_cubic_bubble012_ds2 - 108.0 * d_quartic_bubble_ds2;
2659 27.0 * d_cubic_bubble023_ds0 - 108.0 * d_quartic_bubble_ds0;
2661 27.0 * d_cubic_bubble023_ds1 - 108.0 * d_quartic_bubble_ds1;
2663 27.0 * d_cubic_bubble023_ds2 - 108.0 * d_quartic_bubble_ds2;
2667 27.0 * d_cubic_bubble123_ds0 - 108.0 * d_quartic_bubble_ds0;
2669 27.0 * d_cubic_bubble123_ds1 - 108.0 * d_quartic_bubble_ds1;
2671 27.0 * d_cubic_bubble123_ds2 - 108.0 * d_quartic_bubble_ds2;
2674 dpsids(14, 0) = 256.0 * d_quartic_bubble_ds0;
2675 dpsids(14, 1) = 256.0 * d_quartic_bubble_ds1;
2676 dpsids(14, 2) = 256.0 * d_quartic_bubble_ds2;
2695 this->dshape_local(
s, psi, dpsids);
2698 const double s3 = 1.0 -
s[0] -
s[1] -
s[2];
2706 const double d2_quartic_bubble_ds0 = -2.0 *
s[1] *
s[2];
2707 const double d2_quartic_bubble_ds1 = -2.0 *
s[0] *
s[2];
2708 const double d2_quartic_bubble_ds2 = -2.0 *
s[0] *
s[1];
2709 const double d2_quartic_bubble_ds3 =
2710 s[2] * (1.0 - 2.0 *
s[0] - 2.0 *
s[1] -
s[2]);
2711 const double d2_quartic_bubble_ds4 =
2712 s[1] * (1.0 - 2.0 *
s[0] - 2.0 *
s[2] -
s[1]);
2713 const double d2_quartic_bubble_ds5 =
2714 s[0] * (1.0 - 2.0 *
s[1] - 2.0 *
s[2] -
s[0]);
2716 const double d2_cubic_bubble012_ds0 = 0.0;
2717 const double d2_cubic_bubble012_ds1 = 0.0;
2718 const double d2_cubic_bubble012_ds2 = 0.0;
2719 const double d2_cubic_bubble012_ds3 =
s[2];
2720 const double d2_cubic_bubble012_ds4 =
s[1];
2721 const double d2_cubic_bubble012_ds5 =
s[0];
2723 const double d2_cubic_bubble013_ds0 = -2.0 *
s[1];
2724 const double d2_cubic_bubble013_ds1 = -2.0 *
s[0];
2725 const double d2_cubic_bubble013_ds2 = 0.0;
2726 const double d2_cubic_bubble013_ds3 = s3 -
s[0] -
s[1];
2727 const double d2_cubic_bubble013_ds4 = -
s[1];
2728 const double d2_cubic_bubble013_ds5 = -
s[0];
2730 const double d2_cubic_bubble023_ds0 = -2.0 *
s[2];
2731 const double d2_cubic_bubble023_ds1 = 0.0;
2732 const double d2_cubic_bubble023_ds2 = -2.0 *
s[0];
2733 const double d2_cubic_bubble023_ds3 = -
s[2];
2734 const double d2_cubic_bubble023_ds4 = s3 -
s[0] -
s[2];
2735 const double d2_cubic_bubble023_ds5 = -
s[0];
2737 const double d2_cubic_bubble123_ds0 = 0.0;
2738 const double d2_cubic_bubble123_ds1 = -2.0 *
s[2];
2739 const double d2_cubic_bubble123_ds2 = -2.0 *
s[1];
2740 const double d2_cubic_bubble123_ds3 = -
s[2];
2741 const double d2_cubic_bubble123_ds4 = -
s[1];
2742 const double d2_cubic_bubble123_ds5 = s3 -
s[1] -
s[2];
2745 d2psids(0, 0) = 4.0 +
2746 3.0 * (d2_cubic_bubble012_ds0 + d2_cubic_bubble013_ds0 +
2747 d2_cubic_bubble023_ds0) -
2748 4.0 * d2_quartic_bubble_ds0;
2749 d2psids(0, 1) = 0.0 +
2750 3.0 * (d2_cubic_bubble012_ds1 + d2_cubic_bubble013_ds1 +
2751 d2_cubic_bubble023_ds1) -
2752 4.0 * d2_quartic_bubble_ds1;
2753 d2psids(0, 2) = 0.0 +
2754 3.0 * (d2_cubic_bubble012_ds2 + d2_cubic_bubble013_ds2 +
2755 d2_cubic_bubble023_ds2) -
2756 4.0 * d2_quartic_bubble_ds2;
2757 d2psids(0, 3) = 0.0 +
2758 3.0 * (d2_cubic_bubble012_ds3 + d2_cubic_bubble013_ds3 +
2759 d2_cubic_bubble023_ds3) -
2760 4.0 * d2_quartic_bubble_ds3;
2761 d2psids(0, 4) = 0.0 +
2762 3.0 * (d2_cubic_bubble012_ds4 + d2_cubic_bubble013_ds4 +
2763 d2_cubic_bubble023_ds4) -
2764 4.0 * d2_quartic_bubble_ds4;
2765 d2psids(0, 5) = 0.0 +
2766 3.0 * (d2_cubic_bubble012_ds5 + d2_cubic_bubble013_ds5 +
2767 d2_cubic_bubble023_ds5) -
2768 4.0 * d2_quartic_bubble_ds5;
2771 d2psids(1, 0) = 0.0 +
2772 3.0 * (d2_cubic_bubble012_ds0 + d2_cubic_bubble013_ds0 +
2773 d2_cubic_bubble123_ds0) -
2774 4.0 * d2_quartic_bubble_ds0;
2775 d2psids(1, 1) = 4.0 +
2776 3.0 * (d2_cubic_bubble012_ds1 + d2_cubic_bubble013_ds1 +
2777 d2_cubic_bubble123_ds1) -
2778 4.0 * d2_quartic_bubble_ds1;
2779 d2psids(1, 2) = 0.0 +
2780 3.0 * (d2_cubic_bubble012_ds2 + d2_cubic_bubble013_ds2 +
2781 d2_cubic_bubble123_ds2) -
2782 4.0 * d2_quartic_bubble_ds2;
2783 d2psids(1, 3) = 0.0 +
2784 3.0 * (d2_cubic_bubble012_ds3 + d2_cubic_bubble013_ds3 +
2785 d2_cubic_bubble123_ds3) -
2786 4.0 * d2_quartic_bubble_ds3;
2787 d2psids(1, 4) = 0.0 +
2788 3.0 * (d2_cubic_bubble012_ds4 + d2_cubic_bubble013_ds4 +
2789 d2_cubic_bubble123_ds4) -
2790 4.0 * d2_quartic_bubble_ds4;
2791 d2psids(1, 5) = 0.0 +
2792 3.0 * (d2_cubic_bubble012_ds5 + d2_cubic_bubble013_ds5 +
2793 d2_cubic_bubble123_ds5) -
2794 4.0 * d2_quartic_bubble_ds5;
2797 d2psids(2, 0) = 0.0 +
2798 3.0 * (d2_cubic_bubble012_ds0 + d2_cubic_bubble023_ds0 +
2799 d2_cubic_bubble123_ds0) -
2800 4.0 * d2_quartic_bubble_ds0;
2801 d2psids(2, 1) = 0.0 +
2802 3.0 * (d2_cubic_bubble012_ds1 + d2_cubic_bubble023_ds1 +
2803 d2_cubic_bubble123_ds1) -
2804 4.0 * d2_quartic_bubble_ds1;
2805 d2psids(2, 2) = 4.0 +
2806 3.0 * (d2_cubic_bubble012_ds2 + d2_cubic_bubble023_ds2 +
2807 d2_cubic_bubble123_ds2) -
2808 4.0 * d2_quartic_bubble_ds2;
2809 d2psids(2, 3) = 0.0 +
2810 3.0 * (d2_cubic_bubble012_ds3 + d2_cubic_bubble023_ds3 +
2811 d2_cubic_bubble123_ds3) -
2812 4.0 * d2_quartic_bubble_ds3;
2813 d2psids(2, 4) = 0.0 +
2814 3.0 * (d2_cubic_bubble012_ds4 + d2_cubic_bubble023_ds4 +
2815 d2_cubic_bubble123_ds4) -
2816 4.0 * d2_quartic_bubble_ds4;
2817 d2psids(2, 5) = 0.0 +
2818 3.0 * (d2_cubic_bubble012_ds5 + d2_cubic_bubble023_ds5 +
2819 d2_cubic_bubble123_ds5) -
2820 4.0 * d2_quartic_bubble_ds5;
2823 d2psids(3, 0) = 4.0 +
2824 3.0 * (d2_cubic_bubble013_ds0 + d2_cubic_bubble023_ds0 +
2825 d2_cubic_bubble123_ds0) -
2826 4.0 * d2_quartic_bubble_ds0;
2827 d2psids(3, 1) = 4.0 +
2828 3.0 * (d2_cubic_bubble013_ds1 + d2_cubic_bubble023_ds1 +
2829 d2_cubic_bubble123_ds1) -
2830 4.0 * d2_quartic_bubble_ds1;
2831 d2psids(3, 2) = 4.0 +
2832 3.0 * (d2_cubic_bubble013_ds2 + d2_cubic_bubble023_ds2 +
2833 d2_cubic_bubble123_ds2) -
2834 4.0 * d2_quartic_bubble_ds2;
2835 d2psids(3, 3) = 4.0 +
2836 3.0 * (d2_cubic_bubble013_ds3 + d2_cubic_bubble023_ds3 +
2837 d2_cubic_bubble123_ds3) -
2838 4.0 * d2_quartic_bubble_ds3;
2839 d2psids(3, 4) = 4.0 +
2840 3.0 * (d2_cubic_bubble013_ds4 + d2_cubic_bubble023_ds4 +
2841 d2_cubic_bubble123_ds4) -
2842 4.0 * d2_quartic_bubble_ds4;
2843 d2psids(3, 5) = 4.0 +
2844 3.0 * (d2_cubic_bubble013_ds5 + d2_cubic_bubble023_ds5 +
2845 d2_cubic_bubble123_ds5) -
2846 4.0 * d2_quartic_bubble_ds5;
2849 d2psids(4, 0) = 0.0 -
2850 12.0 * (d2_cubic_bubble012_ds0 + d2_cubic_bubble013_ds0) +
2851 32.0 * d2_quartic_bubble_ds0;
2852 d2psids(4, 1) = 0.0 -
2853 12.0 * (d2_cubic_bubble012_ds1 + d2_cubic_bubble013_ds1) +
2854 32.0 * d2_quartic_bubble_ds1;
2855 d2psids(4, 2) = 0.0 -
2856 12.0 * (d2_cubic_bubble012_ds2 + d2_cubic_bubble013_ds2) +
2857 32.0 * d2_quartic_bubble_ds2;
2858 d2psids(4, 3) = 4.0 -
2859 12.0 * (d2_cubic_bubble012_ds3 + d2_cubic_bubble013_ds3) +
2860 32.0 * d2_quartic_bubble_ds3;
2861 d2psids(4, 4) = 0.0 -
2862 12.0 * (d2_cubic_bubble012_ds4 + d2_cubic_bubble013_ds4) +
2863 32.0 * d2_quartic_bubble_ds4;
2864 d2psids(4, 5) = 0.0 -
2865 12.0 * (d2_cubic_bubble012_ds5 + d2_cubic_bubble013_ds5) +
2866 32.0 * d2_quartic_bubble_ds5;
2869 d2psids(5, 0) = 0.0 -
2870 12.0 * (d2_cubic_bubble012_ds0 + d2_cubic_bubble023_ds0) +
2871 32.0 * d2_quartic_bubble_ds0;
2872 d2psids(5, 1) = 0.0 -
2873 12.0 * (d2_cubic_bubble012_ds1 + d2_cubic_bubble023_ds1) +
2874 32.0 * d2_quartic_bubble_ds1;
2875 d2psids(5, 2) = 0.0 -
2876 12.0 * (d2_cubic_bubble012_ds2 + d2_cubic_bubble023_ds2) +
2877 32.0 * d2_quartic_bubble_ds2;
2878 d2psids(5, 3) = 0.0 -
2879 12.0 * (d2_cubic_bubble012_ds3 + d2_cubic_bubble023_ds3) +
2880 32.0 * d2_quartic_bubble_ds3;
2881 d2psids(5, 4) = 4.0 -
2882 12.0 * (d2_cubic_bubble012_ds4 + d2_cubic_bubble023_ds4) +
2883 32.0 * d2_quartic_bubble_ds4;
2884 d2psids(5, 5) = 0.0 -
2885 12.0 * (d2_cubic_bubble012_ds5 + d2_cubic_bubble023_ds5) +
2886 32.0 * d2_quartic_bubble_ds5;
2889 d2psids(6, 0) = -8.0 -
2890 12.0 * (d2_cubic_bubble013_ds0 + d2_cubic_bubble023_ds0) +
2891 32.0 * d2_quartic_bubble_ds0;
2892 d2psids(6, 1) = 0.0 -
2893 12.0 * (d2_cubic_bubble013_ds1 + d2_cubic_bubble023_ds1) +
2894 32.0 * d2_quartic_bubble_ds1;
2895 d2psids(6, 2) = 0.0 -
2896 12.0 * (d2_cubic_bubble013_ds2 + d2_cubic_bubble023_ds2) +
2897 32.0 * d2_quartic_bubble_ds2;
2898 d2psids(6, 3) = -4.0 -
2899 12.0 * (d2_cubic_bubble013_ds3 + d2_cubic_bubble023_ds3) +
2900 32.0 * d2_quartic_bubble_ds3;
2901 d2psids(6, 4) = -4.0 -
2902 12.0 * (d2_cubic_bubble013_ds4 + d2_cubic_bubble023_ds4) +
2903 32.0 * d2_quartic_bubble_ds4;
2904 d2psids(6, 5) = 0.0 -
2905 12.0 * (d2_cubic_bubble013_ds5 + d2_cubic_bubble023_ds5) +
2906 32.0 * d2_quartic_bubble_ds5;
2908 d2psids(7, 0) = 0.0 -
2909 12.0 * (d2_cubic_bubble012_ds0 + d2_cubic_bubble123_ds0) +
2910 32.0 * d2_quartic_bubble_ds0;
2911 d2psids(7, 1) = 0.0 -
2912 12.0 * (d2_cubic_bubble012_ds1 + d2_cubic_bubble123_ds1) +
2913 32.0 * d2_quartic_bubble_ds1;
2914 d2psids(7, 2) = 0.0 -
2915 12.0 * (d2_cubic_bubble012_ds2 + d2_cubic_bubble123_ds2) +
2916 32.0 * d2_quartic_bubble_ds2;
2917 d2psids(7, 3) = 0.0 -
2918 12.0 * (d2_cubic_bubble012_ds3 + d2_cubic_bubble123_ds3) +
2919 32.0 * d2_quartic_bubble_ds3;
2920 d2psids(7, 4) = 0.0 -
2921 12.0 * (d2_cubic_bubble012_ds4 + d2_cubic_bubble123_ds4) +
2922 32.0 * d2_quartic_bubble_ds4;
2923 d2psids(7, 5) = 4.0 -
2924 12.0 * (d2_cubic_bubble012_ds5 + d2_cubic_bubble123_ds5) +
2925 32.0 * d2_quartic_bubble_ds5;
2927 d2psids(8, 0) = 0.0 -
2928 12.0 * (d2_cubic_bubble023_ds0 + d2_cubic_bubble123_ds0) +
2929 32.0 * d2_quartic_bubble_ds0;
2930 d2psids(8, 1) = 0.0 -
2931 12.0 * (d2_cubic_bubble023_ds1 + d2_cubic_bubble123_ds1) +
2932 32.0 * d2_quartic_bubble_ds1;
2933 d2psids(8, 2) = -8.0 -
2934 12.0 * (d2_cubic_bubble023_ds2 + d2_cubic_bubble123_ds2) +
2935 32.0 * d2_quartic_bubble_ds2;
2936 d2psids(8, 3) = 0.0 -
2937 12.0 * (d2_cubic_bubble023_ds3 + d2_cubic_bubble123_ds3) +
2938 32.0 * d2_quartic_bubble_ds3;
2939 d2psids(8, 4) = -4.0 -
2940 12.0 * (d2_cubic_bubble023_ds4 + d2_cubic_bubble123_ds4) +
2941 32.0 * d2_quartic_bubble_ds4;
2942 d2psids(8, 5) = -4.0 -
2943 12.0 * (d2_cubic_bubble023_ds5 + d2_cubic_bubble123_ds5) +
2944 32.0 * d2_quartic_bubble_ds5;
2946 d2psids(9, 0) = 0.0 -
2947 12.0 * (d2_cubic_bubble013_ds0 + d2_cubic_bubble123_ds0) +
2948 32.0 * d2_quartic_bubble_ds0;
2949 d2psids(9, 1) = -8.0 -
2950 12.0 * (d2_cubic_bubble013_ds1 + d2_cubic_bubble123_ds1) +
2951 32.0 * d2_quartic_bubble_ds1;
2952 d2psids(9, 2) = 0.0 -
2953 12.0 * (d2_cubic_bubble013_ds2 + d2_cubic_bubble123_ds2) +
2954 32.0 * d2_quartic_bubble_ds3;
2955 d2psids(9, 3) = -4.0 -
2956 12.0 * (d2_cubic_bubble013_ds3 + d2_cubic_bubble123_ds3) +
2957 32.0 * d2_quartic_bubble_ds3;
2958 d2psids(9, 4) = 0.0 -
2959 12.0 * (d2_cubic_bubble013_ds4 + d2_cubic_bubble123_ds4) +
2960 32.0 * d2_quartic_bubble_ds4;
2961 d2psids(9, 5) = -4.0 -
2962 12.0 * (d2_cubic_bubble013_ds5 + d2_cubic_bubble123_ds5) +
2963 32.0 * d2_quartic_bubble_ds5;
2967 27.0 * d2_cubic_bubble013_ds0 - 108.0 * d2_quartic_bubble_ds0;
2969 27.0 * d2_cubic_bubble013_ds1 - 108.0 * d2_quartic_bubble_ds1;
2971 27.0 * d2_cubic_bubble013_ds2 - 108.0 * d2_quartic_bubble_ds2;
2973 27.0 * d2_cubic_bubble013_ds3 - 108.0 * d2_quartic_bubble_ds3;
2975 27.0 * d2_cubic_bubble013_ds4 - 108.0 * d2_quartic_bubble_ds4;
2977 27.0 * d2_cubic_bubble013_ds5 - 108.0 * d2_quartic_bubble_ds5;
2981 27.0 * d2_cubic_bubble012_ds0 - 108.0 * d2_quartic_bubble_ds0;
2983 27.0 * d2_cubic_bubble012_ds1 - 108.0 * d2_quartic_bubble_ds1;
2985 27.0 * d2_cubic_bubble012_ds2 - 108.0 * d2_quartic_bubble_ds2;
2987 27.0 * d2_cubic_bubble012_ds3 - 108.0 * d2_quartic_bubble_ds3;
2989 27.0 * d2_cubic_bubble012_ds4 - 108.0 * d2_quartic_bubble_ds4;
2991 27.0 * d2_cubic_bubble012_ds5 - 108.0 * d2_quartic_bubble_ds5;
2995 27.0 * d2_cubic_bubble023_ds0 - 108.0 * d2_quartic_bubble_ds0;
2997 27.0 * d2_cubic_bubble023_ds1 - 108.0 * d2_quartic_bubble_ds1;
2999 27.0 * d2_cubic_bubble023_ds2 - 108.0 * d2_quartic_bubble_ds2;
3001 27.0 * d2_cubic_bubble023_ds3 - 108.0 * d2_quartic_bubble_ds3;
3003 27.0 * d2_cubic_bubble023_ds4 - 108.0 * d2_quartic_bubble_ds4;
3005 27.0 * d2_cubic_bubble023_ds5 - 108.0 * d2_quartic_bubble_ds5;
3009 27.0 * d2_cubic_bubble123_ds0 - 108.0 * d2_quartic_bubble_ds0;
3011 27.0 * d2_cubic_bubble123_ds1 - 108.0 * d2_quartic_bubble_ds1;
3013 27.0 * d2_cubic_bubble123_ds2 - 108.0 * d2_quartic_bubble_ds2;
3015 27.0 * d2_cubic_bubble123_ds3 - 108.0 * d2_quartic_bubble_ds3;
3017 27.0 * d2_cubic_bubble123_ds4 - 108.0 * d2_quartic_bubble_ds4;
3019 27.0 * d2_cubic_bubble123_ds5 - 108.0 * d2_quartic_bubble_ds5;
3022 d2psids(14, 0) = 256.0 * d2_quartic_bubble_ds0;
3023 d2psids(14, 1) = 256.0 * d2_quartic_bubble_ds1;
3024 d2psids(14, 2) = 256.0 * d2_quartic_bubble_ds2;
3025 d2psids(14, 3) = 256.0 * d2_quartic_bubble_ds3;
3026 d2psids(14, 4) = 256.0 * d2_quartic_bubble_ds4;
3027 d2psids(14, 5) = 256.0 * d2_quartic_bubble_ds5;
3039 template<
unsigned NNODE_1D>
3070 "Tets are currently only implemented for\n";
3071 error_message +=
"four and ten nodes, i.e. NNODE_1D=2 , 3 \n";
3074 error_message, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
3080 (NNODE_1D * (NNODE_1D + 1)) / 2 + 1 + 3 * (NNODE_1D - 2);
3081 this->set_n_node(n_node);
3087 set_integration_scheme(&Default_integration_scheme);
3117 const unsigned&
i)
const
3129 std::ostringstream error_message;
3131 <<
"Element only has four vertex nodes; called with node number " << j
3134 OOMPH_CURRENT_FUNCTION,
3135 OOMPH_EXCEPTION_LOCATION);
3178 return FiniteElement::invert_jacobian<3>(jacobian, inverse_jacobian);
3204 unsigned node_sum = 0;
3205 for (
unsigned j = 1; j <= nplot; j++)
3207 for (
unsigned i = 1;
i <= j;
i++)
3219 return (nplot - 1) * (nplot - 1) * (nplot - 1);
3226 const unsigned& nplot,
3227 unsigned& counter)
const
3232 unsigned paraview_fix = counter - 1;
3233 unsigned nod_count = 1;
3235 for (
unsigned i = 0;
i < nplot;
i++)
3237 for (
unsigned j = 0; j < nplot -
i; j++)
3239 for (
unsigned k = 0; k < nplot -
i - j; k++)
3241 if (k < nplot -
i - j - 1)
3243 file_out << nod_count + paraview_fix <<
" "
3244 << nod_count + 1 + paraview_fix <<
" "
3245 << nod_count + nplot -
i - j + paraview_fix <<
" "
3246 << nod_count + nplot -
i - j +
3247 ((nplot - 1 -
i) * (nplot -
i) / 2) + paraview_fix
3249 if (k < nplot -
i - j - 2)
3251 file_out << nod_count + 1 + paraview_fix <<
" "
3252 << nod_count + nplot -
i - j + paraview_fix <<
" "
3253 << nod_count + nplot -
i - j +
3254 ((nplot - 1 -
i) * (nplot -
i) / 2) + paraview_fix
3256 << nod_count + 2 * (nplot -
i - j) - 1 +
3257 ((nplot - 1 -
i) * (nplot -
i) / 2) + paraview_fix
3259 file_out << nod_count + 1 + paraview_fix <<
" "
3260 << nod_count + nplot -
i - j + paraview_fix <<
" "
3261 << nod_count + nplot -
i - j + 1 + paraview_fix <<
" "
3262 << nod_count + 2 * (nplot -
i - j) - 1 +
3263 ((nplot - 1 -
i) * (nplot -
i) / 2) + paraview_fix
3265 file_out << nod_count + 1 + paraview_fix <<
" "
3266 << nod_count + nplot -
i - j +
3267 ((nplot - 1 -
i) * (nplot -
i) / 2) + paraview_fix
3269 << nod_count + nplot -
i - j +
3270 ((nplot - 1 -
i) * (nplot -
i) / 2) + 1 +
3273 << nod_count + 2 * (nplot -
i - j) - 1 +
3274 ((nplot - 1 -
i) * (nplot -
i) / 2) + paraview_fix
3276 file_out << nod_count + 1 + paraview_fix <<
" "
3277 << nod_count + nplot -
i - j + 1 + paraview_fix <<
" "
3278 << nod_count + nplot -
i - j +
3279 ((nplot - 1 -
i) * (nplot -
i) / 2) + 1 +
3282 << nod_count + 2 * (nplot -
i - j) - 1 +
3283 ((nplot - 1 -
i) * (nplot -
i) / 2) + paraview_fix
3288 file_out << nod_count + nplot -
i - j - 1 + paraview_fix <<
" "
3289 << nod_count + nplot -
i - j +
3290 ((nplot - 1 -
i) * (nplot -
i) / 2) - 1 +
3293 << nod_count + 2 * (nplot -
i - j - 1) +
3294 ((nplot - 1 -
i) * (nplot -
i) / 2) - 1 +
3297 << nod_count + 2 * (nplot -
i - j - 1) +
3298 ((nplot - 1 -
i) * (nplot -
i) / 2) + paraview_fix
3308 counter += nplot_points_paraview(nplot);
3316 const unsigned& nplot)
const
3318 unsigned local_loop = nsub_elements_paraview(nplot);
3319 for (
unsigned i = 0;
i < local_loop;
i++)
3321 file_out <<
"10" << std::endl;
3329 const unsigned& nplot,
3330 unsigned& offset_sum)
const
3332 unsigned local_loop = nsub_elements_paraview(nplot);
3333 for (
unsigned i = 0;
i < local_loop;
i++)
3336 file_out << offset_sum << std::endl;
3344 void output(std::ostream& outfile,
const unsigned& nplot);
3347 void output(FILE* file_pt);
3350 void output(FILE* file_pt,
const unsigned& n_plot);
3355 const unsigned& iplot,
3356 const unsigned& nplot,
3358 const bool& use_equally_spaced_interior_sample_points =
false)
const
3363 for (
unsigned i = 0;
i < nplot; ++
i)
3365 for (
unsigned j = 0; j < nplot -
i; ++j)
3367 for (
unsigned k = 0; k < nplot -
i - j; ++k)
3372 s[0] = double(j) / double(nplot - 1);
3373 s[1] = double(
i) / double(nplot - 1);
3374 s[2] = double(k) / double(nplot - 1);
3375 if (use_equally_spaced_interior_sample_points)
3378 double dx_new = range / double(nplot + 1);
3379 double range_new = double(nplot - 1) * dx_new;
3380 s[0] = 0.5 * dx_new + range_new *
s[0] / range;
3381 s[1] = 0.5 * dx_new + range_new *
s[1] / range;
3382 s[2] = 0.5 * dx_new + range_new *
s[2] / range;
3404 std::ostringstream header;
3406 nel = (nplot - 1) * (nplot - 1) * (nplot - 1);
3407 header <<
"ZONE N=" << nplot_points(nplot) <<
", E=" << nel
3408 <<
", F=FEPOINT, ET=TETRAHEDRON\n";
3409 return header.str();
3417 const unsigned& nplot)
const
3421 unsigned nod_count = 1;
3422 for (
unsigned i = 0;
i < nplot;
i++)
3424 for (
unsigned j = 0; j < nplot -
i; j++)
3426 for (
unsigned k = 0; k < nplot -
i - j; k++)
3428 if (k < nplot -
i - j - 1)
3430 outfile << nod_count <<
" " << nod_count + 1 <<
" "
3431 << nod_count + nplot -
i - j <<
" "
3432 << nod_count + nplot -
i - j +
3433 ((nplot - 1 -
i) * (nplot -
i) / 2)
3435 if (k < nplot -
i - j - 2)
3437 outfile << nod_count + 1 <<
" " << nod_count + nplot -
i - j
3439 << nod_count + nplot -
i - j +
3440 ((nplot - 1 -
i) * (nplot -
i) / 2)
3442 << nod_count + 2 * (nplot -
i - j) - 1 +
3443 ((nplot - 1 -
i) * (nplot -
i) / 2)
3445 outfile << nod_count + 1 <<
" " << nod_count + nplot -
i - j
3446 <<
" " << nod_count + nplot -
i - j + 1 <<
" "
3447 << nod_count + 2 * (nplot -
i - j) - 1 +
3448 ((nplot - 1 -
i) * (nplot -
i) / 2)
3450 outfile << nod_count + 1 <<
" "
3451 << nod_count + nplot -
i - j +
3452 ((nplot - 1 -
i) * (nplot -
i) / 2)
3454 << nod_count + nplot -
i - j +
3455 ((nplot - 1 -
i) * (nplot -
i) / 2) + 1
3457 << nod_count + 2 * (nplot -
i - j) - 1 +
3458 ((nplot - 1 -
i) * (nplot -
i) / 2)
3460 outfile << nod_count + 1 <<
" " << nod_count + nplot -
i - j + 1
3462 << nod_count + nplot -
i - j +
3463 ((nplot - 1 -
i) * (nplot -
i) / 2) + 1
3465 << nod_count + 2 * (nplot -
i - j) - 1 +
3466 ((nplot - 1 -
i) * (nplot -
i) / 2)
3471 outfile << nod_count + nplot -
i - j - 1 <<
" "
3472 << nod_count + nplot -
i - j +
3473 ((nplot - 1 -
i) * (nplot -
i) / 2) - 1
3475 << nod_count + 2 * (nplot -
i - j - 1) +
3476 ((nplot - 1 -
i) * (nplot -
i) / 2) - 1
3478 << nod_count + 2 * (nplot -
i - j - 1) +
3479 ((nplot - 1 -
i) * (nplot -
i) / 2)
3498 unsigned nod_count = 1;
3499 for (
unsigned i = 0;
i < nplot;
i++)
3501 for (
unsigned j = 0; j < nplot -
i; j++)
3503 for (
unsigned k = 0; k < nplot -
i - j; k++)
3505 if (j < nplot -
i - 1)
3511 nod_count + nplot -
i);
3512 if (j < nplot -
i - 2)
3517 nod_count + nplot -
i + 1,
3518 nod_count + nplot -
i);
3535 for (
unsigned i = 3;
i <= nplot;
i++)
3537 res += (
i * (
i + 1) / 2);
3552 void build_face_element(
const int& face_index,
3568 template<
unsigned DIM,
unsigned NNODE_1D>
3579 template<
unsigned DIM,
unsigned NPTS_1D>
3620 template<
unsigned DIM>
3638 unsigned n_node = this->nnode();
3639 this->set_n_node(n_node +
3642 this->set_integration_scheme(&Default_enriched_integration_scheme);
3681 s, psi, dpsids, d2psids);
3726 template<
unsigned DIM,
unsigned NNODE_1D>
3735 template<
unsigned NNODE_1D>
3744 set_lagrangian_dimension(1);
3758 inline void build_face_element(
const int& face_index,
3779 template<
unsigned NNODE_1D>
3781 const int& face_index,
FaceElement* face_element_pt)
3788 ->set_lagrangian_dimension(
3796 template<
unsigned NNODE_1D>
3809 set_lagrangian_dimension(2);
3823 inline void build_face_element(
const int& face_index,
3839 template<
unsigned NNODE_1D>
3841 const int& face_index,
FaceElement* face_element_pt)
3848 ->set_lagrangian_dimension(
3856 template<
unsigned NNODE_1D>
3869 set_lagrangian_dimension(3);
3886 inline void build_face_element(
const int& face_index,
3900 template<
unsigned NNODE_1D>
3902 const int& face_index,
FaceElement* face_element_pt)
3909 ->set_lagrangian_dimension(
3924 template<
unsigned DIM,
unsigned NNODE_1D>
3933 template<
unsigned DIM>
3967 template<
unsigned DIM,
unsigned NNODE_1D>
3969 :
public virtual TElement<DIM - 1, NNODE_1D>
3981 template<
unsigned NNODE_1D>
4003 template<
unsigned NNODE_1D>
4005 :
public virtual TElement<1, NNODE_1D>
4021 template<
unsigned NNODE_1D>
4043 template<
unsigned DIM,
unsigned NNODE_1D>
4057 template<
unsigned NNODE_1D>
4080 template<
unsigned NNODE_1D>
4098 template<
unsigned NNODE_1D>
A class that contains the information required by Nodes that are located on Mesh boundaries....
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....
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
A general Finite Element class.
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
virtual bool is_on_boundary() const
Test whether the Node lies on a boundary. The "bulk" Node cannot lie on a boundary,...
virtual void get_boundaries_pt(std::set< unsigned > *&boundaries_pt)
Return a pointer to set of mesh boundaries that this node occupies; this will be overloaded by Bounda...
An OomphLibError object which should be thrown when an run-time error is encountered....
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
////////////////////////////////////////////////////////////////////// //////////////////////////////...
A Class for nodes that deform elastically (i.e. position is an unknown in the problem)....
unsigned nlagrangian() const
Return number of lagrangian coordinates.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
SolidTBubbleEnrichedElement(const SolidTBubbleEnrichedElement &)=delete
Broken copy constructor.
SolidTBubbleEnrichedElement()
Constructor.
void build_face_element(const int &face_index, FaceElement *face_element_pt)
Build the lower-dimensional FaceElement Need to put in a final override here.
~SolidTBubbleEnrichedElement()
Broken assignment operator.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
SolidTElement(const SolidTElement &)=delete
Broken copy constructor.
SolidTElement()
Constructor.
SolidTElement()
Constructor.
SolidTElement(const SolidTElement &)=delete
Broken copy constructor.
SolidTElement()
Constructor.
SolidTElement(const SolidTElement &)=delete
Broken copy constructor.
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
unsigned n_enriched_nodes()
Return the number of nodes required for enrichement.
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Second derivatives of shape functions for specific TBubbleElement<2,3>: d2psids(i,...
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Derivatives of shape functions for specific TBubbleElement<2,3>
void shape(const Vector< double > &s, Shape &psi) const
Shape function for specific TBubbleEnrichedElement<2,3>
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Return local coordinates of node j.
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Second derivatives of shape functions for specific TElement<3,3>: d2psids(i,0) = d2psids(i,...
unsigned n_enriched_nodes()
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Derivatives of shape functions for specific TElement<3,3>
void shape(const Vector< double > &s, Shape &psi) const
Shape function for specific TBubbleEnrichedElement<3,3>
A class for those member functions that must be fully specialised for Telements that are enriched by ...
TBubbleEnrichedElement(const TBubbleEnrichedElement &)=delete
Broken copy constructor.
void build_face_element(const int &face_index, FaceElement *face_element_pt)
Build the lower-dimensional FaceElement.
void shape(const Vector< double > &s, Shape &psi) const
Calculate the geometric shape functions at local coordinate s.
~TBubbleEnrichedElement()
Broken assignment operator.
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Compute the geometric shape functions and derivatives w.r.t. local coordinates at local coordinate s.
static TBubbleEnrichedGauss< DIM, 3 > Default_enriched_integration_scheme
TBubbleEnrichedElement()
Constructor.
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Return local coordinates of node j.
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Compute the geometric shape functions, derivatives and second derivatives w.r.t local coordinates at ...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
const unsigned Central_node_on_face[3]
Define integration schemes that are required to exactly integrate the mass matrices of the bubble-enr...
/////////////////////////////////////////////////////////////////// /////////////////////////////////...
void move_local_coord_back_into_element(Vector< double > &s) const
Adjust local coordinates so that they're located inside the element.
TElementBase()
Empty default constructor.
ElementGeometry::ElementGeometry element_geometry() const
Broken assignment operator.
TElementBase(const TElementBase &)=delete
Broken copy constructor.
bool local_coord_is_valid(const Vector< double > &s)
Check whether the local coordinates are valid or not.
/////////////////////////////////////////////////////////////////// /////////////////////////////////...
TElementGeometricBase(const TElementGeometricBase &)=delete
Broken copy constructor.
TElementGeometricBase()
Empty default constructor.
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Return local coordinates of node j.
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Derivatives of shape functions for specific TElement<2,2>
void shape(const Vector< double > &s, Shape &psi) const
Shape function for specific TElement<1,2>
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Second derivatives of shape functions for specific TElement<1,2>: d2psids(i,0) = .
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Return local coordinates of node j.
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Derivatives of shape functions for specific TElement<1,3>
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Second derivatives of shape functions for specific TElement<1,3>: d2psids(i,0) = .
void shape(const Vector< double > &s, Shape &psi) const
Shape function for specific TElement<1,3>
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Return local coordinates of node j.
void shape(const Vector< double > &s, Shape &psi) const
Shape function for specific TElement<1,4>
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Second derivatives of shape functions for specific TElement<2,4>: d2psids(i,0) = .
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Derivatives of shape functions for specific TElement<1,4>
void shape(const Vector< double > &s, Shape &psi) const
Shape function for specific TElement<2,2>
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Return local coordinates of node j.
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Second derivatives of shape functions for specific TElement<2,2>: d2psids(i,0) = d2psids(i,...
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Derivatives of shape functions for specific TElement<2,2>
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Return local coordinates of node j.
void shape(const Vector< double > &s, Shape &psi) const
Shape function for specific TElement<2,3>
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Derivatives of shape functions for specific TElement<2,3>
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Second derivatives of shape functions for specific TElement<2,3>: d2psids(i,0) = d2psids(i,...
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Derivatives of shape functions for specific TElement<2,4>
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Return local coordinates of node j.
void shape(const Vector< double > &s, Shape &psi) const
Shape function for specific TElement<2,4>
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Second derivatives of shape functions for specific TElement<2,4>: d2psids(i,0) = d2psids(i,...
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Derivatives of shape functions for specific TElement<3,2>
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Second derivatives of shape functions for specific TElement<3,2>: d2psids(i,0) = d2psids(i,...
void shape(const Vector< double > &s, Shape &psi) const
Shape function for specific TElement<3,2>
void shape(const Vector< double > &s, Shape &psi) const
Shape function for specific TElement<3,3>
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Second derivatives of shape functions for specific TElement<3,3>: d2psids(i,0) = d2psids(i,...
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Derivatives of shape functions for specific TElement<3,3>
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Compute the geometric shape functions and derivatives w.r.t. local coordinates at local coordinate s.
std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction)
unsigned nvertex_node() const
Number of vertex nodes in the element: One more than spatial dimension.
unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction)
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Compute the geometric shape functions, derivatives and second derivatives w.r.t local coordinates at ...
void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &use_equally_spaced_interior_sample_points=false) const
Get vector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
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...
~TElement()
Broken assignment operator.
double s_max() const
Max. value of local coordinate.
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...
unsigned nsub_elements_paraview(const unsigned &nplot) const
Return the number of local sub-elements for paraview plot with parameter nplot.
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....
TElement(const TElement &)=delete
Broken copy constructor.
unsigned nnode_1d() const
Number of nodes along each element edge.
void shape(const Vector< double > &s, Shape &psi) const
Calculate the geometric shape functions at local coordinate s.
static TGauss< 1, NNODE_1D > Default_integration_scheme
Default integration rule: Gaussian integration of same 'order' as the element.
double s_min() const
Min. value of local coordinate.
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Return local coordinates of node j.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
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 nplot_points_paraview(const unsigned &nplot) const
Return the number of actual plot points for paraview plot with parameter nplot.
General TElement class specialised to two spatial dimensions Ordering of nodes as in Zienkiwizc sketc...
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...
unsigned nvertex_node() const
Number of vertex nodes in the element: One more than spatial dimension.
double s_min() const
Min. value of local coordinate.
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)...
void shape(const Vector< double > &s, Shape &psi) const
Calculate the geometric shape functions at local coordinate s.
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...
std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction)
unsigned nplot_points_paraview(const unsigned &nplot) const
Return the number of actual plot points for paraview plot with parameter nplot.
double s_max() const
Max. value of local coordinate.
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...
~TElement()
Broken assignment operator.
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 get_s_plot(const unsigned &iplot, const unsigned &nplot, Vector< double > &s, const bool &use_equally_spaced_interior_sample_points=false) const
Get vector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Compute the geometric shape functions and derivatives w.r.t. local coordinates at local coordinate s.
static TGauss< 2, NNODE_1D > Default_integration_scheme
Default integration rule: Gaussian integration of same 'order' as the element.
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Compute the geometric shape functions, derivatives and second derivatives w.r.t local coordinates at ...
unsigned nsub_elements_paraview(const unsigned &nplot) const
Return the number of local sub-elements for paraview plot with parameter nplot.
unsigned nnode_1d() const
Number of nodes along each element edge.
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Return local coordinates of node j.
unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction)
void write_tecplot_zone_footer(FILE *file_pt, const unsigned &nplot) const
Add tecplot zone "footer" to C-style output. (when plotting nplot points in each "coordinate directio...
unsigned get_bulk_node_number(const int &face_index, const unsigned &i) const
Public access function for Node_on_face.
TElement(const TElement &)=delete
Broken copy constructor.
TElement(const bool &allow_high_order)
Alternative constructor.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
unsigned get_bulk_node_number(const int &face_index, const unsigned &i) const
Public access function for Node_on_face.
TElement(const TElement &)=delete
Broken copy constructor.
unsigned nnode_1d() const
Number of nodes along each element edge.
double s_max() const
Max. value of local coordinate.
static TGauss< 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 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)...
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Compute the geometric shape functions, derivatives and second derivatives w.r.t local coordinates at ...
~TElement()
Broken assignment operator.
unsigned nplot_points_paraview(const unsigned &nplot) const
Return the number of actual plot points for paraview plot with parameter nplot.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
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_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...
void write_tecplot_zone_footer(FILE *file_pt, const unsigned &nplot) const
Add tecplot zone "footer" to C-style output. (when plotting nplot points in each "coordinate directio...
void shape(const Vector< double > &s, Shape &psi) const
Calculate the geometric shape functions at local coordinate s.
unsigned nvertex_node() const
Number of vertex nodes in the element: One more than spatial dimension.
unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction)
void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Return local coordinates of node j.
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 s_min() const
Min. 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...
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Compute the geometric shape functions and derivatives w.r.t. local coordinates at local coordinate s.
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 get_s_plot(const unsigned &iplot, const unsigned &nplot, Vector< double > &s, const bool &use_equally_spaced_interior_sample_points=false) const
Get vector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
const unsigned Node_on_face[3][2]
Assign the nodal translation schemes.
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
bool is_on_boundary() const
Test whether the face lies on a boundary. Relatively simple test, based on all vertices lying on (som...
Node * Node2_pt
Second vertex node.
Node * node2_pt() const
Access to the second vertex node.
Node * Node3_pt
Third vertex node.
void get_boundaries_pt(std::set< unsigned > *&boundaries_pt)
Access to pointer to set of mesh boundaries that this face occupies; NULL if the node is not on any b...
TFace(Node *node1_pt, Node *node2_pt, Node *node3_pt)
Constructor: Pass in the three vertex nodes.
Node * node1_pt() const
Access to the first vertex node.
Node * Node1_pt
First vertex node.
bool operator<(const TFace &other) const
Less-than operator.
bool is_boundary_face() const
Test whether the face is a boundary face, i.e. does it connnect three boundary nodes?
Node * node3_pt() const
Access to the third vertex node.
bool operator==(const TFace &other) const
Comparison operator.
Base class for Solid Telements.
TSolidElementBase()
Constructor: Empty.
TSolidElementBase(const TSolidElementBase &)=delete
Broken copy constructor.
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...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...