232 typename PointType::Pointer pPoint2,
233 typename PointType::Pointer pPoint3,
234 typename PointType::Pointer pPoint4,
235 typename PointType::Pointer pPoint5,
236 typename PointType::Pointer pPoint6,
237 typename PointType::Pointer pPoint7,
238 typename PointType::Pointer pPoint8,
239 typename PointType::Pointer pPoint9,
240 typename PointType::Pointer pPoint10,
241 typename PointType::Pointer pPoint11,
242 typename PointType::Pointer pPoint12,
243 typename PointType::Pointer pPoint13,
244 typename PointType::Pointer pPoint14,
245 typename PointType::Pointer pPoint15,
246 typename PointType::Pointer pPoint16,
247 typename PointType::Pointer pPoint17,
248 typename PointType::Pointer pPoint18,
249 typename PointType::Pointer pPoint19,
250 typename PointType::Pointer pPoint20,
251 typename PointType::Pointer pPoint21,
252 typename PointType::Pointer pPoint22,
253 typename PointType::Pointer pPoint23,
254 typename PointType::Pointer pPoint24,
255 typename PointType::Pointer pPoint25,
256 typename PointType::Pointer pPoint26,
257 typename PointType::Pointer pPoint27 )
290 :
BaseType( ThisPoints, &msGeometryData )
300 ) :
BaseType( GeometryId, rThisPoints, &msGeometryData)
307 const std::string& rGeometryName,
309 ) :
BaseType( rGeometryName, rThisPoints, &msGeometryData)
390 template<
class TOtherPo
intType>
414 return typename BaseType::Pointer(
new Hexahedra3D27( NewId, rThisPoints ) );
428 auto p_geometry =
typename BaseType::Pointer(
new Hexahedra3D27( NewGeometryId, rGeometry.
Points() ) );
429 p_geometry->SetData(rGeometry.
GetData());
519 const CoordinatesArrayType& rPoint,
521 const double Tolerance = std::numeric_limits<double>::epsilon()
526 if ( std::abs( rResult[0] ) <= (1.0 + Tolerance) )
528 if ( std::abs( rResult[1] ) <= (1.0 + Tolerance) )
530 if ( std::abs( rResult[2] ) <= (1.0 + Tolerance) )
635 constexpr
double w = 1.0/12.0;
637 return std::accumulate(
641 [](
double sum,
const auto& rEdge) {return sum + rEdge.Length();}
755 double fx1 = 0.5 * ( rPoint[0] - 1.0 ) * ( rPoint[0] );
756 double fx2 = 0.5 * ( rPoint[0] + 1.0 ) * ( rPoint[0] );
757 double fx3 = 1.0 - ( rPoint[0] * rPoint[0] );
758 double fy1 = 0.5 * ( rPoint[1] - 1.0 ) * ( rPoint[1] );
759 double fy2 = 0.5 * ( rPoint[1] + 1.0 ) * ( rPoint[1] );
760 double fy3 = 1.0 - ( rPoint[1] * rPoint[1] );
761 double fz1 = 0.5 * ( rPoint[2] - 1.0 ) * ( rPoint[2] );
762 double fz2 = 0.5 * ( rPoint[2] + 1.0 ) * ( rPoint[2] );
763 double fz3 = 1.0 - ( rPoint[2] * rPoint[2] );
765 switch ( ShapeFunctionIndex )
768 return( fx1*fy1*fz1 );
770 return( fx2*fy1*fz1 );
772 return( fx2*fy2*fz1 );
774 return( fx1*fy2*fz1 );
776 return( fx1*fy1*fz2 );
778 return( fx2*fy1*fz2 );
780 return( fx2*fy2*fz2 );
782 return( fx1*fy2*fz2 );
784 return( fx3*fy1*fz1 );
786 return( fx2*fy3*fz1 );
788 return( fx3*fy2*fz1 );
790 return( fx1*fy3*fz1 );
792 return( fx1*fy1*fz3 );
794 return( fx2*fy1*fz3 );
796 return( fx2*fy2*fz3 );
798 return( fx1*fy2*fz3 );
800 return( fx3*fy1*fz2 );
802 return( fx2*fy3*fz2 );
804 return( fx3*fy2*fz2 );
806 return( fx1*fy3*fz2 );
808 return( fx3*fy3*fz1 );
810 return( fx3*fy1*fz3 );
812 return( fx2*fy3*fz3 );
814 return( fx3*fy2*fz3 );
816 return( fx1*fy3*fz3 );
818 return( fx3*fy3*fz2 );
820 return( fx3*fy3*fz3 );
823 KRATOS_ERROR <<
"Wrong index of shape function!" << *
this << std::endl;
842 if(rResult.size() != 27) rResult.
resize(27,
false);
844 double fx1 = 0.5 * ( rCoordinates[0] - 1.0 ) * ( rCoordinates[0] );
845 double fx2 = 0.5 * ( rCoordinates[0] + 1.0 ) * ( rCoordinates[0] );
846 double fx3 = 1.0 - ( rCoordinates[0] * rCoordinates[0] );
847 double fy1 = 0.5 * ( rCoordinates[1] - 1.0 ) * ( rCoordinates[1] );
848 double fy2 = 0.5 * ( rCoordinates[1] + 1.0 ) * ( rCoordinates[1] );
849 double fy3 = 1.0 - ( rCoordinates[1] * rCoordinates[1] );
850 double fz1 = 0.5 * ( rCoordinates[2] - 1.0 ) * ( rCoordinates[2] );
851 double fz2 = 0.5 * ( rCoordinates[2] + 1.0 ) * ( rCoordinates[2] );
852 double fz3 = 1.0 - ( rCoordinates[2] * rCoordinates[2] );
854 rResult[0] = ( fx1*fy1*fz1 );
855 rResult[1] = ( fx2*fy1*fz1 );
856 rResult[2] = ( fx2*fy2*fz1 );
857 rResult[3] = ( fx1*fy2*fz1 );
858 rResult[4] = ( fx1*fy1*fz2 );
859 rResult[5] = ( fx2*fy1*fz2 );
860 rResult[6] = ( fx2*fy2*fz2 );
861 rResult[7] = ( fx1*fy2*fz2 );
862 rResult[8] = ( fx3*fy1*fz1 );
863 rResult[9] = ( fx2*fy3*fz1 );
864 rResult[10] = ( fx3*fy2*fz1 );
865 rResult[11] = ( fx1*fy3*fz1 );
866 rResult[12] = ( fx1*fy1*fz3 );
867 rResult[13] = ( fx2*fy1*fz3 );
868 rResult[14] = ( fx2*fy2*fz3 );
869 rResult[15] = ( fx1*fy2*fz3 );
870 rResult[16] = ( fx3*fy1*fz2 );
871 rResult[17] = ( fx2*fy3*fz2 );
872 rResult[18] = ( fx3*fy2*fz2 );
873 rResult[19] = ( fx1*fy3*fz2 );
874 rResult[20] = ( fx3*fy3*fz1 );
875 rResult[21] = ( fx3*fy1*fz3 );
876 rResult[22] = ( fx2*fy3*fz3 );
877 rResult[23] = ( fx3*fy2*fz3 );
878 rResult[24] = ( fx1*fy3*fz3 );
879 rResult[25] = ( fx3*fy3*fz2 );
880 rResult[26] = ( fx3*fy3*fz3 );
887 constexpr std::array<std::array<std::size_t, 3>, 48> triangle_faces{{
888 { 0, 11, 8}, {11, 20, 8}, { 8, 20, 1}, {20, 9, 1}, {11, 3, 20}, { 3, 10, 20}, {20, 10, 9}, {10, 2, 9},
889 { 0, 12, 24}, {24, 11, 0}, {11, 24, 15}, {15, 3, 11}, {12, 4, 19}, {19, 24, 12}, {24, 19, 7}, { 7, 15, 24},
890 { 3, 15, 23}, {23, 10, 3}, {10, 23, 14}, {14, 2, 10}, {15, 7, 18}, {18, 23, 15}, {23, 18, 6}, { 6, 14, 23},
891 { 2, 14, 22}, {22, 9, 2}, { 9, 22, 13}, {13, 1, 9}, {14, 6, 17}, {17, 22, 14}, {22, 17, 5}, { 5, 13, 22},
892 {12, 0, 8}, {21, 12, 8}, {21, 8, 1}, {13, 21, 1}, { 4, 12, 21}, {21, 16, 4}, {16, 21, 13}, {13, 5, 16},
893 { 4, 16, 19}, {16, 25, 19}, {16, 5, 25}, { 5, 17, 25}, {19, 25, 7}, {25, 18, 7}, {25, 17, 18}, {17, 6, 18}
896 for (
const auto& r_nodes: triangle_faces) {
899 if (face.HasIntersection(rLowPoint, rHighPoint))
return true;
904 return IsInside(rLowPoint,local_coordinates);
918 std::string
Info()
const override
920 return "3 dimensional hexahedra with 27 nodes and quadratic shape functions in 3D space";
932 rOStream <<
"3 dimensional hexahedra with 27 nodes and quadratic shape functions in 3D space";
948 std::cout << std::endl;
954 rOStream <<
" Jacobian in the origin\t : " << jacobian;
972 double fx1 = 0.5 * ( rPoint[0] - 1.0 ) * ( rPoint[0] );
973 double fx2 = 0.5 * ( rPoint[0] + 1.0 ) * ( rPoint[0] );
974 double fx3 = 1.0 - ( rPoint[0] * rPoint[0] );
975 double fy1 = 0.5 * ( rPoint[1] - 1.0 ) * ( rPoint[1] );
976 double fy2 = 0.5 * ( rPoint[1] + 1.0 ) * ( rPoint[1] );
977 double fy3 = 1.0 - ( rPoint[1] * rPoint[1] );
978 double fz1 = 0.5 * ( rPoint[2] - 1.0 ) * ( rPoint[2] );
979 double fz2 = 0.5 * ( rPoint[2] + 1.0 ) * ( rPoint[2] );
980 double fz3 = 1.0 - ( rPoint[2] * rPoint[2] );
982 double gx1 = 0.5 * ( 2.0 * rPoint[0] - 1.0 );
983 double gx2 = 0.5 * ( 2.0 * rPoint[0] + 1.0 );
984 double gx3 = -2 * rPoint[0];
985 double gy1 = 0.5 * ( 2.0 * rPoint[1] - 1.0 );
986 double gy2 = 0.5 * ( 2.0 * rPoint[1] + 1.0 );
987 double gy3 = -2 * rPoint[1];
988 double gz1 = 0.5 * ( 2.0 * rPoint[2] - 1.0 );
989 double gz2 = 0.5 * ( 2.0 * rPoint[2] + 1.0 );
990 double gz3 = -2 * rPoint[2];
993 if ( result.size1() != 27 || result.size2() != 3 )
994 result.
resize( 27, 3,
false );
996 result( 0, 0 ) = gx1 * fy1 * fz1;
998 result( 0, 1 ) = fx1 * gy1 * fz1;
1000 result( 0, 2 ) = fx1 * fy1 * gz1;
1002 result( 1, 0 ) = gx2 * fy1 * fz1;
1004 result( 1, 1 ) = fx2 * gy1 * fz1;
1006 result( 1, 2 ) = fx2 * fy1 * gz1;
1008 result( 2, 0 ) = gx2 * fy2 * fz1;
1010 result( 2, 1 ) = fx2 * gy2 * fz1;
1012 result( 2, 2 ) = fx2 * fy2 * gz1;
1014 result( 3, 0 ) = gx1 * fy2 * fz1;
1016 result( 3, 1 ) = fx1 * gy2 * fz1;
1018 result( 3, 2 ) = fx1 * fy2 * gz1;
1020 result( 4, 0 ) = gx1 * fy1 * fz2;
1022 result( 4, 1 ) = fx1 * gy1 * fz2;
1024 result( 4, 2 ) = fx1 * fy1 * gz2;
1026 result( 5, 0 ) = gx2 * fy1 * fz2;
1028 result( 5, 1 ) = fx2 * gy1 * fz2;
1030 result( 5, 2 ) = fx2 * fy1 * gz2;
1032 result( 6, 0 ) = gx2 * fy2 * fz2;
1034 result( 6, 1 ) = fx2 * gy2 * fz2;
1036 result( 6, 2 ) = fx2 * fy2 * gz2;
1038 result( 7, 0 ) = gx1 * fy2 * fz2;
1040 result( 7, 1 ) = fx1 * gy2 * fz2;
1042 result( 7, 2 ) = fx1 * fy2 * gz2;
1044 result( 8, 0 ) = gx3 * fy1 * fz1;
1046 result( 8, 1 ) = fx3 * gy1 * fz1;
1048 result( 8, 2 ) = fx3 * fy1 * gz1;
1050 result( 9, 0 ) = gx2 * fy3 * fz1;
1052 result( 9, 1 ) = fx2 * gy3 * fz1;
1054 result( 9, 2 ) = fx2 * fy3 * gz1;
1056 result( 10, 0 ) = gx3 * fy2 * fz1;
1058 result( 10, 1 ) = fx3 * gy2 * fz1;
1060 result( 10, 2 ) = fx3 * fy2 * gz1;
1062 result( 11, 0 ) = gx1 * fy3 * fz1;
1064 result( 11, 1 ) = fx1 * gy3 * fz1;
1066 result( 11, 2 ) = fx1 * fy3 * gz1;
1068 result( 12, 0 ) = gx1 * fy1 * fz3;
1070 result( 12, 1 ) = fx1 * gy1 * fz3;
1072 result( 12, 2 ) = fx1 * fy1 * gz3;
1074 result( 13, 0 ) = gx2 * fy1 * fz3;
1076 result( 13, 1 ) = fx2 * gy1 * fz3;
1078 result( 13, 2 ) = fx2 * fy1 * gz3;
1080 result( 14, 0 ) = gx2 * fy2 * fz3;
1082 result( 14, 1 ) = fx2 * gy2 * fz3;
1084 result( 14, 2 ) = fx2 * fy2 * gz3;
1086 result( 15, 0 ) = gx1 * fy2 * fz3;
1088 result( 15, 1 ) = fx1 * gy2 * fz3;
1090 result( 15, 2 ) = fx1 * fy2 * gz3;
1092 result( 16, 0 ) = gx3 * fy1 * fz2;
1094 result( 16, 1 ) = fx3 * gy1 * fz2;
1096 result( 16, 2 ) = fx3 * fy1 * gz2;
1098 result( 17, 0 ) = gx2 * fy3 * fz2;
1100 result( 17, 1 ) = fx2 * gy3 * fz2;
1102 result( 17, 2 ) = fx2 * fy3 * gz2;
1104 result( 18, 0 ) = gx3 * fy2 * fz2;
1106 result( 18, 1 ) = fx3 * gy2 * fz2;
1108 result( 18, 2 ) = fx3 * fy2 * gz2;
1110 result( 19, 0 ) = gx1 * fy3 * fz2;
1112 result( 19, 1 ) = fx1 * gy3 * fz2;
1114 result( 19, 2 ) = fx1 * fy3 * gz2;
1116 result( 20, 0 ) = gx3 * fy3 * fz1;
1118 result( 20, 1 ) = fx3 * gy3 * fz1;
1120 result( 20, 2 ) = fx3 * fy3 * gz1;
1122 result( 21, 0 ) = gx3 * fy1 * fz3;
1124 result( 21, 1 ) = fx3 * gy1 * fz3;
1126 result( 21, 2 ) = fx3 * fy1 * gz3;
1128 result( 22, 0 ) = gx2 * fy3 * fz3;
1130 result( 22, 1 ) = fx2 * gy3 * fz3;
1132 result( 22, 2 ) = fx2 * fy3 * gz3;
1134 result( 23, 0 ) = gx3 * fy2 * fz3;
1136 result( 23, 1 ) = fx3 * gy2 * fz3;
1138 result( 23, 2 ) = fx3 * fy2 * gz3;
1140 result( 24, 0 ) = gx1 * fy3 * fz3;
1142 result( 24, 1 ) = fx1 * gy3 * fz3;
1144 result( 24, 2 ) = fx1 * fy3 * gz3;
1146 result( 25, 0 ) = gx3 * fy3 * fz2;
1148 result( 25, 1 ) = fx3 * gy3 * fz2;
1150 result( 25, 2 ) = fx3 * fy3 * gz2;
1152 result( 26, 0 ) = gx3 * fy3 * fz3;
1154 result( 26, 1 ) = fx3 * gy3 * fz3;
1156 result( 26, 2 ) = fx3 * fy3 * gz3;
1175 void save(
Serializer& rSerializer )
const override
1205 static Matrix CalculateShapeFunctionsIntegrationPointsValues(
1209 AllIntegrationPoints();
1212 const int integration_points_number = integration_points.size();
1214 const int points_number = 27;
1216 Matrix shape_function_values( integration_points_number, points_number );
1219 for (
int pnt = 0; pnt < integration_points_number; pnt++ )
1221 double fx1 = 0.5 * ( integration_points[pnt].X() - 1.0 ) * ( integration_points[pnt].
X() );
1222 double fx2 = 0.5 * ( integration_points[pnt].X() + 1.0 ) * ( integration_points[pnt].
X() );
1223 double fx3 = 1.0 - ( integration_points[pnt].X() * integration_points[pnt].X() );
1224 double fy1 = 0.5 * ( integration_points[pnt].Y() - 1.0 ) * ( integration_points[pnt].
Y() );
1225 double fy2 = 0.5 * ( integration_points[pnt].Y() + 1.0 ) * ( integration_points[pnt].
Y() );
1226 double fy3 = 1.0 - ( integration_points[pnt].Y() * integration_points[pnt].Y() );
1227 double fz1 = 0.5 * ( integration_points[pnt].Z() - 1.0 ) * ( integration_points[pnt].
Z() );
1228 double fz2 = 0.5 * ( integration_points[pnt].Z() + 1.0 ) * ( integration_points[pnt].
Z() );
1229 double fz3 = 1.0 - ( integration_points[pnt].Z() * integration_points[pnt].Z() );
1231 shape_function_values( pnt, 0 ) = ( fx1 * fy1 * fz1 );
1232 shape_function_values( pnt, 1 ) = ( fx2 * fy1 * fz1 );
1233 shape_function_values( pnt, 2 ) = ( fx2 * fy2 * fz1 );
1234 shape_function_values( pnt, 3 ) = ( fx1 * fy2 * fz1 );
1235 shape_function_values( pnt, 4 ) = ( fx1 * fy1 * fz2 );
1236 shape_function_values( pnt, 5 ) = ( fx2 * fy1 * fz2 );
1237 shape_function_values( pnt, 6 ) = ( fx2 * fy2 * fz2 );
1238 shape_function_values( pnt, 7 ) = ( fx1 * fy2 * fz2 );
1239 shape_function_values( pnt, 8 ) = ( fx3 * fy1 * fz1 );
1240 shape_function_values( pnt, 9 ) = ( fx2 * fy3 * fz1 );
1241 shape_function_values( pnt, 10 ) = ( fx3 * fy2 * fz1 );
1242 shape_function_values( pnt, 11 ) = ( fx1 * fy3 * fz1 );
1243 shape_function_values( pnt, 12 ) = ( fx1 * fy1 * fz3 );
1244 shape_function_values( pnt, 13 ) = ( fx2 * fy1 * fz3 );
1245 shape_function_values( pnt, 14 ) = ( fx2 * fy2 * fz3 );
1246 shape_function_values( pnt, 15 ) = ( fx1 * fy2 * fz3 );
1247 shape_function_values( pnt, 16 ) = ( fx3 * fy1 * fz2 );
1248 shape_function_values( pnt, 17 ) = ( fx2 * fy3 * fz2 );
1249 shape_function_values( pnt, 18 ) = ( fx3 * fy2 * fz2 );
1250 shape_function_values( pnt, 19 ) = ( fx1 * fy3 * fz2 );
1251 shape_function_values( pnt, 20 ) = ( fx3 * fy3 * fz1 );
1252 shape_function_values( pnt, 21 ) = ( fx3 * fy1 * fz3 );
1253 shape_function_values( pnt, 22 ) = ( fx2 * fy3 * fz3 );
1254 shape_function_values( pnt, 23 ) = ( fx3 * fy2 * fz3 );
1255 shape_function_values( pnt, 24 ) = ( fx1 * fy3 * fz3 );
1256 shape_function_values( pnt, 25 ) = ( fx3 * fy3 * fz2 );
1257 shape_function_values( pnt, 26 ) = ( fx3 * fy3 * fz3 );
1260 return shape_function_values;
1275 CalculateShapeFunctionsIntegrationPointsLocalGradients(
1279 AllIntegrationPoints();
1282 const int integration_points_number = integration_points.size();
1286 for (
int pnt = 0; pnt < integration_points_number; pnt++ )
1288 double fx1 = 0.5 * ( integration_points[pnt].X() - 1.0 ) * ( integration_points[pnt].
X() );
1289 double fx2 = 0.5 * ( integration_points[pnt].X() + 1.0 ) * ( integration_points[pnt].
X() );
1290 double fx3 = 1.0 - ( integration_points[pnt].X() * integration_points[pnt].X() );
1291 double fy1 = 0.5 * ( integration_points[pnt].Y() - 1.0 ) * ( integration_points[pnt].
Y() );
1292 double fy2 = 0.5 * ( integration_points[pnt].Y() + 1.0 ) * ( integration_points[pnt].
Y() );
1293 double fy3 = 1.0 - ( integration_points[pnt].Y() * integration_points[pnt].Y() );
1294 double fz1 = 0.5 * ( integration_points[pnt].Z() - 1.0 ) * ( integration_points[pnt].
Z() );
1295 double fz2 = 0.5 * ( integration_points[pnt].Z() + 1.0 ) * ( integration_points[pnt].
Z() );
1296 double fz3 = 1.0 - ( integration_points[pnt].Z() * integration_points[pnt].Z() );
1298 double gx1 = 0.5 * ( 2.0 * integration_points[pnt].X() - 1.0 );
1299 double gx2 = 0.5 * ( 2.0 * integration_points[pnt].X() + 1.0 );
1300 double gx3 = -2 * integration_points[pnt].X();
1301 double gy1 = 0.5 * ( 2.0 * integration_points[pnt].Y() - 1.0 );
1302 double gy2 = 0.5 * ( 2.0 * integration_points[pnt].Y() + 1.0 );
1303 double gy3 = -2 * integration_points[pnt].Y();
1304 double gz1 = 0.5 * ( 2.0 * integration_points[pnt].Z() - 1.0 );
1305 double gz2 = 0.5 * ( 2.0 * integration_points[pnt].Z() + 1.0 );
1306 double gz3 = -2 * integration_points[pnt].Z();
1310 result( 0, 0 ) = gx1 * fy1 * fz1;
1311 result( 0, 1 ) = fx1 * gy1 * fz1;
1312 result( 0, 2 ) = fx1 * fy1 * gz1;
1314 result( 1, 0 ) = gx2 * fy1 * fz1;
1315 result( 1, 1 ) = fx2 * gy1 * fz1;
1316 result( 1, 2 ) = fx2 * fy1 * gz1;
1318 result( 2, 0 ) = gx2 * fy2 * fz1;
1319 result( 2, 1 ) = fx2 * gy2 * fz1;
1320 result( 2, 2 ) = fx2 * fy2 * gz1;
1322 result( 3, 0 ) = gx1 * fy2 * fz1;
1323 result( 3, 1 ) = fx1 * gy2 * fz1;
1324 result( 3, 2 ) = fx1 * fy2 * gz1;
1326 result( 4, 0 ) = gx1 * fy1 * fz2;
1327 result( 4, 1 ) = fx1 * gy1 * fz2;
1328 result( 4, 2 ) = fx1 * fy1 * gz2;
1330 result( 5, 0 ) = gx2 * fy1 * fz2;
1331 result( 5, 1 ) = fx2 * gy1 * fz2;
1332 result( 5, 2 ) = fx2 * fy1 * gz2;
1334 result( 6, 0 ) = gx2 * fy2 * fz2;
1335 result( 6, 1 ) = fx2 * gy2 * fz2;
1336 result( 6, 2 ) = fx2 * fy2 * gz2;
1338 result( 7, 0 ) = gx1 * fy2 * fz2;
1339 result( 7, 1 ) = fx1 * gy2 * fz2;
1340 result( 7, 2 ) = fx1 * fy2 * gz2;
1342 result( 8, 0 ) = gx3 * fy1 * fz1;
1343 result( 8, 1 ) = fx3 * gy1 * fz1;
1344 result( 8, 2 ) = fx3 * fy1 * gz1;
1346 result( 9, 0 ) = gx2 * fy3 * fz1;
1347 result( 9, 1 ) = fx2 * gy3 * fz1;
1348 result( 9, 2 ) = fx2 * fy3 * gz1;
1350 result( 10, 0 ) = gx3 * fy2 * fz1;
1351 result( 10, 1 ) = fx3 * gy2 * fz1;
1352 result( 10, 2 ) = fx3 * fy2 * gz1;
1354 result( 11, 0 ) = gx1 * fy3 * fz1;
1355 result( 11, 1 ) = fx1 * gy3 * fz1;
1356 result( 11, 2 ) = fx1 * fy3 * gz1;
1358 result( 12, 0 ) = gx1 * fy1 * fz3;
1359 result( 12, 1 ) = fx1 * gy1 * fz3;
1360 result( 12, 2 ) = fx1 * fy1 * gz3;
1362 result( 13, 0 ) = gx2 * fy1 * fz3;
1363 result( 13, 1 ) = fx2 * gy1 * fz3;
1364 result( 13, 2 ) = fx2 * fy1 * gz3;
1366 result( 14, 0 ) = gx2 * fy2 * fz3;
1367 result( 14, 1 ) = fx2 * gy2 * fz3;
1368 result( 14, 2 ) = fx2 * fy2 * gz3;
1370 result( 15, 0 ) = gx1 * fy2 * fz3;
1371 result( 15, 1 ) = fx1 * gy2 * fz3;
1372 result( 15, 2 ) = fx1 * fy2 * gz3;
1374 result( 16, 0 ) = gx3 * fy1 * fz2;
1375 result( 16, 1 ) = fx3 * gy1 * fz2;
1376 result( 16, 2 ) = fx3 * fy1 * gz2;
1378 result( 17, 0 ) = gx2 * fy3 * fz2;
1379 result( 17, 1 ) = fx2 * gy3 * fz2;
1380 result( 17, 2 ) = fx2 * fy3 * gz2;
1382 result( 18, 0 ) = gx3 * fy2 * fz2;
1383 result( 18, 1 ) = fx3 * gy2 * fz2;
1384 result( 18, 2 ) = fx3 * fy2 * gz2;
1386 result( 19, 0 ) = gx1 * fy3 * fz2;
1387 result( 19, 1 ) = fx1 * gy3 * fz2;
1388 result( 19, 2 ) = fx1 * fy3 * gz2;
1390 result( 20, 0 ) = gx3 * fy3 * fz1;
1391 result( 20, 1 ) = fx3 * gy3 * fz1;
1392 result( 20, 2 ) = fx3 * fy3 * gz1;
1394 result( 21, 0 ) = gx3 * fy1 * fz3;
1395 result( 21, 1 ) = fx3 * gy1 * fz3;
1396 result( 21, 2 ) = fx3 * fy1 * gz3;
1398 result( 22, 0 ) = gx2 * fy3 * fz3;
1399 result( 22, 1 ) = fx2 * gy3 * fz3;
1400 result( 22, 2 ) = fx2 * fy3 * gz3;
1402 result( 23, 0 ) = gx3 * fy2 * fz3;
1403 result( 23, 1 ) = fx3 * gy2 * fz3;
1404 result( 23, 2 ) = fx3 * fy2 * gz3;
1406 result( 24, 0 ) = gx1 * fy3 * fz3;
1407 result( 24, 1 ) = fx1 * gy3 * fz3;
1408 result( 24, 2 ) = fx1 * fy3 * gz3;
1410 result( 25, 0 ) = gx3 * fy3 * fz2;
1411 result( 25, 1 ) = fx3 * gy3 * fz2;
1412 result( 25, 2 ) = fx3 * fy3 * gz2;
1414 result( 26, 0 ) = gx3 * fy3 * fz3;
1415 result( 26, 1 ) = fx3 * gy3 * fz3;
1416 result( 26, 2 ) = fx3 * fy3 * gz3;
1418 d_shape_f_values[pnt] = result;
1421 return d_shape_f_values;
1432 if ( rResult.size() != this->PointsNumber() ) {
1438 rResult[
i].resize( 3, 3,
false );
1441 const double fx1 = 0.5 * ( rPoint[0] - 1 ) * rPoint[0];
1442 const double fx2 = 0.5 * ( rPoint[0] + 1 ) * rPoint[0];
1443 const double fx3 = 1 - rPoint[0] * rPoint[0];
1444 const double fy1 = 0.5 * ( rPoint[1] - 1 ) * rPoint[1];
1445 const double fy2 = 0.5 * ( rPoint[1] + 1 ) * rPoint[1];
1446 const double fy3 = 1 - rPoint[1] * rPoint[1];
1447 const double fz1 = 0.5 * ( rPoint[2] - 1 ) * rPoint[2];
1448 const double fz2 = 0.5 * ( rPoint[2] + 1 ) * rPoint[2];
1449 const double fz3 = 1 - rPoint[2] * rPoint[2];
1451 const double gx1 = 0.5 * ( 2 * rPoint[0] - 1 );
1452 const double gx2 = 0.5 * ( 2 * rPoint[0] + 1 );
1453 const double gx3 = -2.0 * rPoint[0];
1454 const double gy1 = 0.5 * ( 2 * rPoint[1] - 1 );
1455 const double gy2 = 0.5 * ( 2 * rPoint[1] + 1 );
1456 const double gy3 = -2.0 * rPoint[1];
1457 const double gz1 = 0.5 * ( 2 * rPoint[2] - 1 );
1458 const double gz2 = 0.5 * ( 2 * rPoint[2] + 1 );
1459 const double gz3 = -2.0 * rPoint[2];
1461 const double hx1 = 1.0;
1462 const double hx2 = 1.0;
1463 const double hx3 = -2.0;
1464 const double hy1 = 1.0;
1465 const double hy2 = 1.0;
1466 const double hy3 = -2.0;
1467 const double hz1 = 1.0;
1468 const double hz2 = 1.0;
1469 const double hz3 = -2.0;
1471 rResult[0]( 0, 0 ) = hx1 * fy1 * fz1;
1472 rResult[0]( 0, 1 ) = gx1 * gy1 * fz1;
1473 rResult[0]( 0, 2 ) = gx1 * fy1 * gz1;
1474 rResult[0]( 1, 0 ) = gx1 * gy1 * fz1;
1475 rResult[0]( 1, 1 ) = fx1 * hy1 * fz1;
1476 rResult[0]( 1, 2 ) = fx1 * gy1 * gz1;
1477 rResult[0]( 2, 0 ) = gx1 * fy1 * gz1;
1478 rResult[0]( 2, 1 ) = fx1 * gy1 * gz1;
1479 rResult[0]( 2, 2 ) = fx1 * fy1 * hz1;
1481 rResult[1]( 0, 0 ) = hx2 * fy1 * fz1;
1482 rResult[1]( 0, 1 ) = gx2 * gy1 * fz1;
1483 rResult[1]( 0, 2 ) = gx2 * fy1 * gz1;
1484 rResult[1]( 1, 0 ) = gx2 * gy1 * fz1;
1485 rResult[1]( 1, 1 ) = fx2 * hy1 * fz1;
1486 rResult[1]( 1, 2 ) = fx2 * gy1 * gz1;
1487 rResult[1]( 2, 0 ) = gx2 * fy1 * gz1;
1488 rResult[1]( 2, 1 ) = fx2 * gy1 * gz1;
1489 rResult[1]( 2, 2 ) = fx2 * fy1 * hz1;
1491 rResult[2]( 0, 0 ) = hx2 * fy2 * fz1;
1492 rResult[2]( 0, 1 ) = gx2 * gy2 * fz1;
1493 rResult[2]( 0, 2 ) = gx2 * fy2 * gz1;
1494 rResult[2]( 1, 0 ) = gx2 * gy2 * fz1;
1495 rResult[2]( 1, 1 ) = fx2 * hy2 * fz1;
1496 rResult[2]( 1, 2 ) = fx2 * gy2 * gz1;
1497 rResult[2]( 2, 0 ) = gx2 * fy2 * gz1;
1498 rResult[2]( 2, 1 ) = fx2 * gy2 * gz1;
1499 rResult[2]( 2, 2 ) = fx2 * fy2 * hz1;
1501 rResult[3]( 0, 0 ) = hx1 * fy2 * fz1;
1502 rResult[3]( 0, 1 ) = gx1 * gy2 * fz1;
1503 rResult[3]( 0, 2 ) = gx1 * fy2 * gz1;
1504 rResult[3]( 1, 0 ) = gx1 * gy2 * fz1;
1505 rResult[3]( 1, 1 ) = fx1 * hy2 * fz1;
1506 rResult[3]( 1, 2 ) = fx1 * gy2 * gz1;
1507 rResult[3]( 2, 0 ) = gx1 * fy2 * gz1;
1508 rResult[3]( 2, 1 ) = fx1 * gy2 * gz1;
1509 rResult[3]( 2, 2 ) = fx1 * fy2 * hz1;
1511 rResult[4]( 0, 0 ) = hx1 * fy1 * fz2;
1512 rResult[4]( 0, 1 ) = gx1 * gy1 * fz2;
1513 rResult[4]( 0, 2 ) = gx1 * fy1 * gz2;
1514 rResult[4]( 1, 0 ) = gx1 * gy1 * fz2;
1515 rResult[4]( 1, 1 ) = fx1 * hy1 * fz2;
1516 rResult[4]( 1, 2 ) = fx1 * gy1 * gz2;
1517 rResult[4]( 2, 0 ) = gx1 * fy1 * gz2;
1518 rResult[4]( 2, 1 ) = fx1 * gy1 * gz2;
1519 rResult[4]( 2, 2 ) = fx1 * fy1 * hz2;
1521 rResult[5]( 0, 0 ) = hx2 * fy1 * fz2;
1522 rResult[5]( 0, 1 ) = gx2 * gy1 * fz2;
1523 rResult[5]( 0, 2 ) = gx2 * fy1 * gz2;
1524 rResult[5]( 1, 0 ) = gx2 * gy1 * fz2;
1525 rResult[5]( 1, 1 ) = fx2 * hy1 * fz2;
1526 rResult[5]( 1, 2 ) = fx2 * gy1 * gz2;
1527 rResult[5]( 2, 0 ) = gx2 * fy1 * gz2;
1528 rResult[5]( 2, 1 ) = fx2 * gy1 * gz2;
1529 rResult[5]( 2, 2 ) = fx2 * fy1 * hz2;
1531 rResult[6]( 0, 0 ) = hx2 * fy2 * fz2;
1532 rResult[6]( 0, 1 ) = gx2 * gy2 * fz2;
1533 rResult[6]( 0, 2 ) = gx2 * fy2 * gz2;
1534 rResult[6]( 1, 0 ) = gx2 * gy2 * fz2;
1535 rResult[6]( 1, 1 ) = fx2 * hy2 * fz2;
1536 rResult[6]( 1, 2 ) = fx2 * gy2 * gz2;
1537 rResult[6]( 2, 0 ) = gx2 * fy2 * gz2;
1538 rResult[6]( 2, 1 ) = fx2 * gy2 * gz2;
1539 rResult[6]( 2, 2 ) = fx2 * fy2 * hz2;
1541 rResult[7]( 0, 0 ) = hx1 * fy2 * fz2;
1542 rResult[7]( 0, 1 ) = gx1 * gy2 * fz2;
1543 rResult[7]( 0, 2 ) = gx1 * fy2 * gz2;
1544 rResult[7]( 1, 0 ) = gx1 * gy2 * fz2;
1545 rResult[7]( 1, 1 ) = fx1 * hy2 * fz2;
1546 rResult[7]( 1, 2 ) = fx1 * gy2 * gz2;
1547 rResult[7]( 2, 0 ) = gx1 * fy2 * gz2;
1548 rResult[7]( 2, 1 ) = fx1 * gy2 * gz2;
1549 rResult[7]( 2, 2 ) = fx1 * fy2 * hz2;
1551 rResult[8]( 0, 0 ) = hx3 * fy1 * fz1;
1552 rResult[8]( 0, 1 ) = gx3 * gy1 * fz1;
1553 rResult[8]( 0, 2 ) = gx3 * fy1 * gz1;
1554 rResult[8]( 1, 0 ) = gx3 * gy1 * fz1;
1555 rResult[8]( 1, 1 ) = fx3 * hy1 * fz1;
1556 rResult[8]( 1, 2 ) = fx3 * gy1 * gz1;
1557 rResult[8]( 2, 0 ) = gx3 * fy1 * gz1;
1558 rResult[8]( 2, 1 ) = fx3 * gy1 * gz1;
1559 rResult[8]( 2, 2 ) = fx3 * fy1 * hz1;
1561 rResult[9]( 0, 0 ) = hx2 * fy3 * fz1;
1562 rResult[9]( 0, 1 ) = gx2 * gy3 * fz1;
1563 rResult[9]( 0, 2 ) = gx2 * fy3 * gz1;
1564 rResult[9]( 1, 0 ) = gx2 * gy3 * fz1;
1565 rResult[9]( 1, 1 ) = fx2 * hy3 * fz1;
1566 rResult[9]( 1, 2 ) = fx2 * gy3 * gz1;
1567 rResult[9]( 2, 0 ) = gx2 * fy3 * gz1;
1568 rResult[9]( 2, 1 ) = fx2 * gy3 * gz1;
1569 rResult[9]( 2, 2 ) = fx2 * fy3 * hz1;
1571 rResult[10]( 0, 0 ) = hx3 * fy2 * fz1;
1572 rResult[10]( 0, 1 ) = gx3 * gy2 * fz1;
1573 rResult[10]( 0, 2 ) = gx3 * fy2 * gz1;
1574 rResult[10]( 1, 0 ) = gx3 * gy2 * fz1;
1575 rResult[10]( 1, 1 ) = fx3 * hy2 * fz1;
1576 rResult[10]( 1, 2 ) = fx3 * gy2 * gz1;
1577 rResult[10]( 2, 0 ) = gx3 * fy2 * gz1;
1578 rResult[10]( 2, 1 ) = fx3 * gy2 * gz1;
1579 rResult[10]( 2, 2 ) = fx3 * fy2 * hz1;
1581 rResult[11]( 0, 0 ) = hx1 * fy3 * fz1;
1582 rResult[11]( 0, 1 ) = gx1 * gy3 * fz1;
1583 rResult[11]( 0, 2 ) = gx1 * fy3 * gz1;
1584 rResult[11]( 1, 0 ) = gx1 * gy3 * fz1;
1585 rResult[11]( 1, 1 ) = fx1 * hy3 * fz1;
1586 rResult[11]( 1, 2 ) = fx1 * gy3 * gz1;
1587 rResult[11]( 2, 0 ) = gx1 * fy3 * gz1;
1588 rResult[11]( 2, 1 ) = fx1 * gy3 * gz1;
1589 rResult[11]( 2, 2 ) = fx1 * fy3 * hz1;
1591 rResult[12]( 0, 0 ) = hx1 * fy1 * fz3;
1592 rResult[12]( 0, 1 ) = gx1 * gy1 * fz3;
1593 rResult[12]( 0, 2 ) = gx1 * fy1 * gz3;
1594 rResult[12]( 1, 0 ) = gx1 * gy1 * fz3;
1595 rResult[12]( 1, 1 ) = fx1 * hy1 * fz3;
1596 rResult[12]( 1, 2 ) = fx1 * gy1 * gz3;
1597 rResult[12]( 2, 0 ) = gx1 * fy1 * gz3;
1598 rResult[12]( 2, 1 ) = fx1 * gy1 * gz3;
1599 rResult[12]( 2, 2 ) = fx1 * fy1 * hz3;
1601 rResult[13]( 0, 0 ) = hx2 * fy1 * fz3;
1602 rResult[13]( 0, 1 ) = gx2 * gy1 * fz3;
1603 rResult[13]( 0, 2 ) = gx2 * fy1 * gz3;
1604 rResult[13]( 1, 0 ) = gx2 * gy1 * fz3;
1605 rResult[13]( 1, 1 ) = fx2 * hy1 * fz3;
1606 rResult[13]( 1, 2 ) = fx2 * gy1 * gz3;
1607 rResult[13]( 2, 0 ) = gx2 * fy1 * gz3;
1608 rResult[13]( 2, 1 ) = fx2 * gy1 * gz3;
1609 rResult[13]( 2, 2 ) = fx2 * fy1 * hz3;
1611 rResult[14]( 0, 0 ) = hx2 * fy2 * fz3;
1612 rResult[14]( 0, 1 ) = gx2 * gy2 * fz3;
1613 rResult[14]( 0, 2 ) = gx2 * fy2 * gz3;
1614 rResult[14]( 1, 0 ) = gx2 * gy2 * fz3;
1615 rResult[14]( 1, 1 ) = fx2 * hy2 * fz3;
1616 rResult[14]( 1, 2 ) = fx2 * gy2 * gz3;
1617 rResult[14]( 2, 0 ) = gx2 * fy2 * gz3;
1618 rResult[14]( 2, 1 ) = fx2 * gy2 * gz3;
1619 rResult[14]( 2, 2 ) = fx2 * fy2 * hz3;
1621 rResult[15]( 0, 0 ) = hx1 * fy2 * fz3;
1622 rResult[15]( 0, 1 ) = gx1 * gy2 * fz3;
1623 rResult[15]( 0, 2 ) = gx1 * fy2 * gz3;
1624 rResult[15]( 1, 0 ) = gx1 * gy2 * fz3;
1625 rResult[15]( 1, 1 ) = fx1 * hy2 * fz3;
1626 rResult[15]( 1, 2 ) = fx1 * gy2 * gz3;
1627 rResult[15]( 2, 0 ) = gx1 * fy2 * gz3;
1628 rResult[15]( 2, 1 ) = fx1 * gy2 * gz3;
1629 rResult[15]( 2, 2 ) = fx1 * fy2 * hz3;
1631 rResult[16]( 0, 0 ) = hx3 * fy1 * fz2;
1632 rResult[16]( 0, 1 ) = gx3 * gy1 * fz2;
1633 rResult[16]( 0, 2 ) = gx3 * fy1 * gz2;
1634 rResult[16]( 1, 0 ) = gx3 * gy1 * fz2;
1635 rResult[16]( 1, 1 ) = fx3 * hy1 * fz2;
1636 rResult[16]( 1, 2 ) = fx3 * gy1 * gz2;
1637 rResult[16]( 2, 0 ) = gx3 * fy1 * gz2;
1638 rResult[16]( 2, 1 ) = fx3 * gy1 * gz2;
1639 rResult[16]( 2, 2 ) = fx3 * fy1 * hz2;
1641 rResult[17]( 0, 0 ) = hx2 * fy3 * fz2;
1642 rResult[17]( 0, 1 ) = gx2 * gy3 * fz2;
1643 rResult[17]( 0, 2 ) = gx2 * fy3 * gz2;
1644 rResult[17]( 1, 0 ) = gx2 * gy3 * fz2;
1645 rResult[17]( 1, 1 ) = fx2 * hy3 * fz2;
1646 rResult[17]( 1, 2 ) = fx2 * gy3 * gz2;
1647 rResult[17]( 2, 0 ) = gx2 * fy3 * gz2;
1648 rResult[17]( 2, 1 ) = fx2 * gy3 * gz2;
1649 rResult[17]( 2, 2 ) = fx2 * fy3 * hz2;
1651 rResult[18]( 0, 0 ) = hx3 * fy2 * fz2;
1652 rResult[18]( 0, 1 ) = gx3 * gy2 * fz2;
1653 rResult[18]( 0, 2 ) = gx3 * fy2 * gz2;
1654 rResult[18]( 1, 0 ) = gx3 * gy2 * fz2;
1655 rResult[18]( 1, 1 ) = fx3 * hy2 * fz2;
1656 rResult[18]( 1, 2 ) = fx3 * gy2 * gz2;
1657 rResult[18]( 2, 0 ) = gx3 * fy2 * gz2;
1658 rResult[18]( 2, 1 ) = fx3 * gy2 * gz2;
1659 rResult[18]( 2, 2 ) = fx3 * fy2 * hz2;
1661 rResult[19]( 0, 0 ) = hx1 * fy3 * fz2;
1662 rResult[19]( 0, 1 ) = gx1 * gy3 * fz2;
1663 rResult[19]( 0, 2 ) = gx1 * fy3 * gz2;
1664 rResult[19]( 1, 0 ) = gx1 * gy3 * fz2;
1665 rResult[19]( 1, 1 ) = fx1 * hy3 * fz2;
1666 rResult[19]( 1, 2 ) = fx1 * gy3 * gz2;
1667 rResult[19]( 2, 0 ) = gx1 * fy3 * gz2;
1668 rResult[19]( 2, 1 ) = fx1 * gy3 * gz2;
1669 rResult[19]( 2, 2 ) = fx1 * fy3 * hz2;
1671 rResult[20]( 0, 0 ) = hx3 * fy3 * fz1;
1672 rResult[20]( 0, 1 ) = gx3 * gy3 * fz1;
1673 rResult[20]( 0, 2 ) = gx3 * fy3 * gz1;
1674 rResult[20]( 1, 0 ) = gx3 * gy3 * fz1;
1675 rResult[20]( 1, 1 ) = fx3 * hy3 * fz1;
1676 rResult[20]( 1, 2 ) = fx3 * gy3 * gz1;
1677 rResult[20]( 2, 0 ) = gx3 * fy3 * gz1;
1678 rResult[20]( 2, 1 ) = fx3 * gy3 * gz1;
1679 rResult[20]( 2, 2 ) = fx3 * fy3 * hz1;
1681 rResult[21]( 0, 0 ) = hx3 * fy1 * fz3;
1682 rResult[21]( 0, 1 ) = gx3 * gy1 * fz3;
1683 rResult[21]( 0, 2 ) = gx3 * fy1 * gz3;
1684 rResult[21]( 1, 0 ) = gx3 * gy1 * fz3;
1685 rResult[21]( 1, 1 ) = fx3 * hy1 * fz3;
1686 rResult[21]( 1, 2 ) = fx3 * gy1 * gz3;
1687 rResult[21]( 2, 0 ) = gx3 * fy1 * gz3;
1688 rResult[21]( 2, 1 ) = fx3 * gy1 * gz3;
1689 rResult[21]( 2, 2 ) = fx3 * fy1 * hz3;
1691 rResult[22]( 0, 0 ) = hx2 * fy3 * fz3;
1692 rResult[22]( 0, 1 ) = gx2 * gy3 * fz3;
1693 rResult[22]( 0, 2 ) = gx2 * fy3 * gz3;
1694 rResult[22]( 1, 0 ) = gx2 * gy3 * fz3;
1695 rResult[22]( 1, 1 ) = fx2 * hy3 * fz3;
1696 rResult[22]( 1, 2 ) = fx2 * gy3 * gz3;
1697 rResult[22]( 2, 0 ) = gx2 * fy3 * gz3;
1698 rResult[22]( 2, 1 ) = fx2 * gy3 * gz3;
1699 rResult[22]( 2, 2 ) = fx2 * fy3 * hz3;
1701 rResult[23]( 0, 0 ) = hx3 * fy2 * fz3;
1702 rResult[23]( 0, 1 ) = gx3 * gy2 * fz3;
1703 rResult[23]( 0, 2 ) = gx3 * fy2 * gz3;
1704 rResult[23]( 1, 0 ) = gx3 * gy2 * fz3;
1705 rResult[23]( 1, 1 ) = fx3 * hy2 * fz3;
1706 rResult[23]( 1, 2 ) = fx3 * gy2 * gz3;
1707 rResult[23]( 2, 0 ) = gx3 * fy2 * gz3;
1708 rResult[23]( 2, 1 ) = fx3 * gy2 * gz3;
1709 rResult[23]( 2, 2 ) = fx3 * fy2 * hz3;
1711 rResult[24]( 0, 0 ) = hx1 * fy3 * fz3;
1712 rResult[24]( 0, 1 ) = gx1 * gy3 * fz3;
1713 rResult[24]( 0, 2 ) = gx1 * fy3 * gz3;
1714 rResult[24]( 1, 0 ) = gx1 * gy3 * fz3;
1715 rResult[24]( 1, 1 ) = fx1 * hy3 * fz3;
1716 rResult[24]( 1, 2 ) = fx1 * gy3 * gz3;
1717 rResult[24]( 2, 0 ) = gx1 * fy3 * gz3;
1718 rResult[24]( 2, 1 ) = fx1 * gy3 * gz3;
1719 rResult[24]( 2, 2 ) = fx1 * fy3 * hz3;
1721 rResult[25]( 0, 0 ) = hx3 * fy3 * fz2;
1722 rResult[25]( 0, 1 ) = gx3 * gy3 * fz2;
1723 rResult[25]( 0, 2 ) = gx3 * fy3 * gz2;
1724 rResult[25]( 1, 0 ) = gx3 * gy3 * fz2;
1725 rResult[25]( 1, 1 ) = fx3 * hy3 * fz2;
1726 rResult[25]( 1, 2 ) = fx3 * gy3 * gz2;
1727 rResult[25]( 2, 0 ) = gx3 * fy3 * gz2;
1728 rResult[25]( 2, 1 ) = fx3 * gy3 * gz2;
1729 rResult[25]( 2, 2 ) = fx3 * fy3 * hz2;
1731 rResult[26]( 0, 0 ) = hx3 * fy3 * fz3;
1732 rResult[26]( 0, 1 ) = gx3 * gy3 * fz3;
1733 rResult[26]( 0, 2 ) = gx3 * fy3 * gz3;
1734 rResult[26]( 1, 0 ) = gx3 * gy3 * fz3;
1735 rResult[26]( 1, 1 ) = fx3 * hy3 * fz3;
1736 rResult[26]( 1, 2 ) = fx3 * gy3 * gz3;
1737 rResult[26]( 2, 0 ) = gx3 * fy3 * gz3;
1738 rResult[26]( 2, 1 ) = fx3 * gy3 * gz3;
1739 rResult[26]( 2, 2 ) = fx3 * fy3 * hz3;
1753 Quadrature < HexahedronGaussLegendreIntegrationPoints1,
1754 3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
1755 Quadrature < HexahedronGaussLegendreIntegrationPoints2,
1756 3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
1757 Quadrature < HexahedronGaussLegendreIntegrationPoints3,
1758 3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
1759 Quadrature < HexahedronGaussLegendreIntegrationPoints4,
1760 3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
1761 Quadrature < HexahedronGaussLegendreIntegrationPoints5,
1762 3, IntegrationPoint<3> >::GenerateIntegrationPoints()
1765 return integration_points;
1773 Hexahedra3D27<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1775 Hexahedra3D27<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1777 Hexahedra3D27<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1779 Hexahedra3D27<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1781 Hexahedra3D27<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1785 return shape_functions_values;
1789 AllShapeFunctionsLocalGradients()
1794 Hexahedra3D27<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1796 Hexahedra3D27<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1798 Hexahedra3D27<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1800 Hexahedra3D27<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1802 Hexahedra3D27<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1806 return shape_functions_local_gradients;
1841 rOStream << std::endl;
1847 template<
class TPo
intType>
const
1848 GeometryData Hexahedra3D27<TPointType>::msGeometryData(
1851 Hexahedra3D27<TPointType>::AllIntegrationPoints(),
1852 Hexahedra3D27<TPointType>::AllShapeFunctionsValues(),
1853 AllShapeFunctionsLocalGradients()
1856 template<
class TPo
intType>
Definition: geometry_data.h:60
KratosGeometryType
Definition: geometry_data.h:110
IntegrationMethod
Definition: geometry_data.h:76
KratosGeometryFamily
Definition: geometry_data.h:91
Definition: geometry_dimension.h:42
Geometry base class.
Definition: geometry.h:71
SizeType PointsNumber() const
Definition: geometry.h:528
virtual void PrintData(std::ostream &rOStream) const
Print object's data.
Definition: geometry.h:3834
Geometry & operator=(const Geometry &rOther)
Definition: geometry.h:400
std::vector< IntegrationPointType > IntegrationPointsArrayType
Definition: geometry.h:161
DataValueContainer & GetData()
Definition: geometry.h:591
std::size_t SizeType
Definition: geometry.h:144
const TPointType::Pointer pGetPoint(const int Index) const
Definition: geometry.h:1790
std::size_t IndexType
Definition: geometry.h:137
std::array< Matrix, static_cast< int >GeometryData::IntegrationMethod::NumberOfIntegrationMethods)> ShapeFunctionsValuesContainerType
Definition: geometry.h:172
Vector & DeterminantOfJacobian(Vector &rResult) const
Definition: geometry.h:3154
virtual CoordinatesArrayType & PointLocalCoordinates(CoordinatesArrayType &rResult, const CoordinatesArrayType &rPoint) const
Returns the local coordinates of a given arbitrary point.
Definition: geometry.h:1854
const PointsArrayType & Points() const
Definition: geometry.h:1768
bool AllPointsAreValid() const
Checks if the geometry points are valid Checks if the geometry points are valid from the pointer valu...
Definition: geometry.h:4093
std::array< IntegrationPointsArrayType, static_cast< int >GeometryData::IntegrationMethod::NumberOfIntegrationMethods)> IntegrationPointsContainerType
Definition: geometry.h:167
JacobiansType & Jacobian(JacobiansType &rResult) const
Definition: geometry.h:2901
GeometryData::ShapeFunctionsLocalGradientsContainerType ShapeFunctionsLocalGradientsContainerType
Definition: geometry.h:177
GeometryData::IntegrationMethod IntegrationMethod
Definition: geometry.h:122
A twenty-seven node hexahedra geometry with second order shape functions.
Definition: hexahedra_3d_27.h:53
double Volume() const override
This method calculate and return volume of this geometry.
Definition: hexahedra_3d_27.h:486
BaseType::IndexType IndexType
Definition: hexahedra_3d_27.h:96
TPointType PointType
Definition: hexahedra_3d_27.h:89
GeometryData::KratosGeometryFamily GetGeometryFamily() const override
Definition: hexahedra_3d_27.h:348
Hexahedra3D27(const std::string &rGeometryName, const PointsArrayType &rThisPoints)
Constructor with Geometry Name.
Definition: hexahedra_3d_27.h:306
BaseType::ShapeFunctionsValuesContainerType ShapeFunctionsValuesContainerType
Definition: hexahedra_3d_27.h:139
SizeType EdgesNumber() const override
This method gives you number of all edges of this geometry.
Definition: hexahedra_3d_27.h:555
BaseType::Pointer Create(IndexType NewId, PointsArrayType const &rThisPoints) const override
Creates a new geometry pointer.
Definition: hexahedra_3d_27.h:409
BaseType::Pointer Create(const IndexType NewGeometryId, const BaseType &rGeometry) const override
Creates a new geometry pointer.
Definition: hexahedra_3d_27.h:423
Hexahedra3D27(Hexahedra3D27< TOtherPointType > const &rOther)
Definition: hexahedra_3d_27.h:340
GeometryData::IntegrationMethod IntegrationMethod
Definition: hexahedra_3d_27.h:78
void PrintData(std::ostream &rOStream) const override
Definition: hexahedra_3d_27.h:944
Hexahedra3D27(typename PointType::Pointer pPoint1, typename PointType::Pointer pPoint2, typename PointType::Pointer pPoint3, typename PointType::Pointer pPoint4, typename PointType::Pointer pPoint5, typename PointType::Pointer pPoint6, typename PointType::Pointer pPoint7, typename PointType::Pointer pPoint8, typename PointType::Pointer pPoint9, typename PointType::Pointer pPoint10, typename PointType::Pointer pPoint11, typename PointType::Pointer pPoint12, typename PointType::Pointer pPoint13, typename PointType::Pointer pPoint14, typename PointType::Pointer pPoint15, typename PointType::Pointer pPoint16, typename PointType::Pointer pPoint17, typename PointType::Pointer pPoint18, typename PointType::Pointer pPoint19, typename PointType::Pointer pPoint20, typename PointType::Pointer pPoint21, typename PointType::Pointer pPoint22, typename PointType::Pointer pPoint23, typename PointType::Pointer pPoint24, typename PointType::Pointer pPoint25, typename PointType::Pointer pPoint26, typename PointType::Pointer pPoint27)
Definition: hexahedra_3d_27.h:231
Matrix & ShapeFunctionsLocalGradients(Matrix &result, const CoordinatesArrayType &rPoint) const override
Definition: hexahedra_3d_27.h:969
Hexahedra3D27(const PointType &Point1, const PointType &Point2, const PointType &Point3, const PointType &Point4, const PointType &Point5, const PointType &Point6, const PointType &Point7, const PointType &Point8, const PointType &Point9, const PointType &Point10, const PointType &Point11, const PointType &Point12, const PointType &Point13, const PointType &Point14, const PointType &Point15, const PointType &Point16, const PointType &Point17, const PointType &Point18, const PointType &Point19, const PointType &Point20, const PointType &Point21, const PointType &Point22, const PointType &Point23, const PointType &Point24, const PointType &Point25, const PointType &Point26, const PointType &Point27)
Definition: hexahedra_3d_27.h:190
BaseType::CoordinatesArrayType CoordinatesArrayType
Definition: hexahedra_3d_27.h:178
double DomainSize() const override
Definition: hexahedra_3d_27.h:505
GeometriesArrayType GenerateEdges() const override
This method gives you all edges of this geometry.
Definition: hexahedra_3d_27.h:568
double ShapeFunctionValue(IndexType ShapeFunctionIndex, const CoordinatesArrayType &rPoint) const override
Definition: hexahedra_3d_27.h:752
Quadrilateral3D9< TPointType > FaceType
Definition: hexahedra_3d_27.h:68
KRATOS_CLASS_POINTER_DEFINITION(Hexahedra3D27)
void PrintInfo(std::ostream &rOStream) const override
Definition: hexahedra_3d_27.h:930
BaseType::PointsArrayType PointsArrayType
Definition: hexahedra_3d_27.h:110
BaseType::JacobiansType JacobiansType
Definition: hexahedra_3d_27.h:153
GeometriesArrayType GenerateFaces() const override
Returns all faces of the current geometry.
Definition: hexahedra_3d_27.h:669
~Hexahedra3D27() override
Destructor. Does nothing!!!
Definition: hexahedra_3d_27.h:346
BaseType::IntegrationPointsContainerType IntegrationPointsContainerType
Definition: hexahedra_3d_27.h:132
BaseType::ShapeFunctionsGradientsType ShapeFunctionsGradientsType
Definition: hexahedra_3d_27.h:160
Line3D3< TPointType > EdgeType
Definition: hexahedra_3d_27.h:67
bool IsInside(const CoordinatesArrayType &rPoint, CoordinatesArrayType &rResult, const double Tolerance=std::numeric_limits< double >::epsilon()) const override
Returns whether given arbitrary point is inside the Geometry and the respective local point for the g...
Definition: hexahedra_3d_27.h:518
BaseType::NormalType NormalType
Definition: hexahedra_3d_27.h:173
BaseType::IntegrationPointType IntegrationPointType
Definition: hexahedra_3d_27.h:117
Geometry< TPointType > BaseType
Definition: hexahedra_3d_27.h:62
Hexahedra3D27(Hexahedra3D27 const &rOther)
Definition: hexahedra_3d_27.h:323
SizeType FacesNumber() const override
Returns the number of faces of the current geometry.
Definition: hexahedra_3d_27.h:656
Hexahedra3D27(const IndexType GeometryId, const PointsArrayType &rThisPoints)
Constructor with Geometry Id.
Definition: hexahedra_3d_27.h:297
Hexahedra3D27 & operator=(const Hexahedra3D27 &rOther)
Definition: hexahedra_3d_27.h:373
double Area() const override
Definition: hexahedra_3d_27.h:473
double AverageEdgeLength() const override
Definition: hexahedra_3d_27.h:633
GeometryData::KratosGeometryType GetGeometryType() const override
Definition: hexahedra_3d_27.h:353
Hexahedra3D27(const PointsArrayType &ThisPoints)
Definition: hexahedra_3d_27.h:289
Hexahedra3D27 & operator=(Hexahedra3D27< TOtherPointType > const &rOther)
Definition: hexahedra_3d_27.h:391
Matrix MatrixType
Definition: hexahedra_3d_27.h:183
BaseType::IntegrationPointsArrayType IntegrationPointsArrayType
Definition: hexahedra_3d_27.h:125
BaseType::ShapeFunctionsSecondDerivativesType ShapeFunctionsSecondDerivativesType
Definition: hexahedra_3d_27.h:168
std::string Info() const override
Definition: hexahedra_3d_27.h:918
friend class Hexahedra3D27
Definition: hexahedra_3d_27.h:1814
BaseType::GeometriesArrayType GeometriesArrayType
Definition: hexahedra_3d_27.h:84
bool HasIntersection(const Point &rLowPoint, const Point &rHighPoint) const override
Definition: hexahedra_3d_27.h:885
double Length() const override
Definition: hexahedra_3d_27.h:453
BaseType::ShapeFunctionsLocalGradientsContainerType ShapeFunctionsLocalGradientsContainerType
Definition: hexahedra_3d_27.h:146
BaseType::SizeType SizeType
Definition: hexahedra_3d_27.h:104
Vector & ShapeFunctionsValues(Vector &rResult, const CoordinatesArrayType &rCoordinates) const override
Definition: hexahedra_3d_27.h:840
Short class definition.
Definition: integration_point.h:52
static double ComputeVolume3DGeometry(const Geometry< TPointType > &rGeometry)
This method calculates and returns the volume of the geometry from a 3D geometry.
Definition: integration_utilities.h:168
Definition: amatrix_interface.h:41
void resize(std::size_t NewSize1, std::size_t NewSize2, bool preserve=0)
Definition: amatrix_interface.h:224
An three node 3D line geometry with quadratic shape functions.
Definition: line_3d_3.h:66
This class defines the node.
Definition: node.h:65
Point class.
Definition: point.h:59
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
iterator end()
Definition: pointer_vector.h:177
iterator begin()
Definition: pointer_vector.h:169
void push_back(const TPointerType &x)
Definition: pointer_vector.h:270
A nine node 3D quadrilateral geometry with quadratic shape functions.
Definition: quadrilateral_3d_9.h:48
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
A three node 3D triangle geometry with linear shape functions.
Definition: triangle_3d_3.h:77
Short class definition.
Definition: array_1d.h:61
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
#define KRATOS_ERROR
Definition: exception.h:161
#define KRATOS_ERROR_IF(conditional)
Definition: exception.h:162
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroMatrix< double > ZeroMatrix
Definition: amatrix_interface.h:559
Internals::Matrix< double, AMatrix::dynamic, AMatrix::dynamic > Matrix
Definition: amatrix_interface.h:470
TExpressionType::data_type sum(AMatrix::MatrixExpression< TExpressionType, TCategory > const &TheExpression)
Definition: amatrix_interface.h:675
const GeometryData Hexahedra3D27< TPointType >::msGeometryData & msGeometryDimension(), Hexahedra3D27< TPointType >::AllShapeFunctionsValues(), AllShapeFunctionsLocalGradients()
Definition: brep_curve.h:483
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
const GeometryData BrepCurve< TContainerPointType, TContainerPointEmbeddedType >::msGeometryData & msGeometryDimension
Definition: brep_curve.h:483
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
w
Definition: generate_convection_diffusion_explicit_element.py:108
def load(f)
Definition: ode_solve.py:307
integer i
Definition: TensorModule.f:17