13 #if !defined(KRATOS_HEXAHEDRA_INTERFACE_3D_8_H_INCLUDED )
14 #define KRATOS_HEXAHEDRA_INTERFACE_3D_8_H_INCLUDED
287 typename PointType::Pointer pPoint2,
288 typename PointType::Pointer pPoint3,
289 typename PointType::Pointer pPoint4,
290 typename PointType::Pointer pPoint5,
291 typename PointType::Pointer pPoint6,
292 typename PointType::Pointer pPoint7,
293 typename PointType::Pointer pPoint8 )
298 noalias(vx) += - *pPoint1 - *pPoint5;
299 noalias(vx) += - *pPoint4 - *pPoint8;
300 noalias(vx) += *pPoint2 + *pPoint6;
301 noalias(vx) += *pPoint3 + *pPoint7;
306 noalias(vy) += *pPoint4 + *pPoint3;
307 noalias(vy) += *pPoint7 + *pPoint8;
308 noalias(vy) += - *pPoint1 - *pPoint2;
309 noalias(vy) += - *pPoint6 - *pPoint5;
314 noalias(vz) += *pPoint5 + *pPoint6;
315 noalias(vz) += *pPoint7 + *pPoint8;
316 noalias(vz) += - *pPoint1 - *pPoint2;
317 noalias(vz) += - *pPoint3 - *pPoint4;
382 :
BaseType( ThisPoints, &msGeometryData )
464 template<
class TOtherPo
intType>
503 return std::sqrt(
Area());
550 double pos = 0.5 + 0.5 / std::sqrt(3.0);
553 double C1 = pos*(p1z - p2z + p3z - p4z);
554 double C2 = pos*(p1y - p2y + p3y - p4y);
555 double C3 = pos*(p1x - p2x + p3x - p4x);
556 double C4 = C1 - p1z + p2z;
557 double C5 = C1 + p1z - p2z;
558 double C6 = C2 + p1y - p2y;
559 double C7 = C2 - p1y + p2y;
560 double C8 = C3 - p1x + p2x;
561 double C9 = C3 + p1x - p2x;
562 double C10 = C1 - p1z + p4z;
563 double C11 = C2 - p1y + p4y;
564 double C12 = C3 - p1x + p4x;
565 double C13 = C1 + p1z - p4z;
566 double C14 = C2 + p1y - p4y;
567 double C15 = C3 + p1x - p4x;
570 std::sqrt( std::pow(C4*C11 - C7*C10, 2) + std::pow(C4*C12 - C8*C10, 2) + std::pow(C7*C12 - C8*C11, 2)) +
571 std::sqrt( std::pow(C5*C11 - C6*C10, 2) + std::pow(C5*C12 - C9*C10, 2) + std::pow(C6*C12 - C9*C11, 2)) +
572 std::sqrt( std::pow(C4*C14 - C7*C13, 2) + std::pow(C4*C15 - C8*C13, 2) + std::pow(C7*C15 - C8*C14, 2)) +
573 std::sqrt( std::pow(C5*C14 - C6*C13, 2) + std::pow(C5*C15 - C9*C13, 2) + std::pow(C6*C15 - C9*C14, 2))
608 if ( rResult.size1() != 8 || rResult.size2() != 3 )
609 rResult.
resize( 8, 3 ,
false);
611 rResult( 0, 0 ) = -1.0;
612 rResult( 0, 1 ) = -1.0;
613 rResult( 0, 2 ) = -1.0;
615 rResult( 1, 0 ) = 1.0;
616 rResult( 1, 1 ) = -1.0;
617 rResult( 1, 2 ) = -1.0;
619 rResult( 2, 0 ) = 1.0;
620 rResult( 2, 1 ) = 1.0;
621 rResult( 2, 2 ) = -1.0;
623 rResult( 3, 0 ) = -1.0;
624 rResult( 3, 1 ) = 1.0;
625 rResult( 3, 2 ) = -1.0;
627 rResult( 4, 0 ) = -1.0;
628 rResult( 4, 1 ) = -1.0;
629 rResult( 4, 2 ) = 1.0;
631 rResult( 5, 0 ) = 1.0;
632 rResult( 5, 1 ) = -1.0;
633 rResult( 5, 2 ) = 1.0;
635 rResult( 6, 0 ) = 1.0;
636 rResult( 6, 1 ) = 1.0;
637 rResult( 6, 2 ) = 1.0;
639 rResult( 7, 0 ) = -1.0;
640 rResult( 7, 1 ) = 1.0;
641 rResult( 7, 2 ) = 1.0;
657 const double Tolerance = std::numeric_limits<double>::epsilon()
662 if ( std::abs( rResult[0] ) <= (1.0 + Tolerance) )
664 if ( std::abs( rResult[1] ) <= (1.0 + Tolerance) )
666 if ( std::abs( rResult[2] ) <= (1.0 + Tolerance) )
689 for(
unsigned int i=0;
i<4;
i++)
709 for (
int k = 0;
k < maxiter;
k++ )
714 noalias( CurrentGlobalCoords ) = rPoint - CurrentGlobalCoords;
718 Matrix shape_functions_gradients;
726 double det_j =
J( 0, 0 ) *
J( 1, 1 ) -
J( 0, 1 ) *
J( 1, 0 );
729 invJ( 0, 0 ) = (
J( 1, 1 ) ) / ( det_j );
730 invJ( 1, 0 ) = -(
J( 1, 0 ) ) / ( det_j );
731 invJ( 0, 1 ) = -(
J( 0, 1 ) ) / ( det_j );
732 invJ( 1, 1 ) = (
J( 0, 0 ) ) / ( det_j );
735 DeltaXi( 0 ) = invJ( 0, 0 ) *
res[0] + invJ( 0, 1 ) *
res[1];
736 DeltaXi( 1 ) = invJ( 1, 0 ) *
res[0] + invJ( 1, 1 ) *
res[1];
738 rResult[0] += DeltaXi[0];
739 rResult[1] += DeltaXi[1];
742 if (
norm_2( DeltaXi ) > 300 )
746 std::cout <<
"detJ =" << det_j <<
"DeltaX = " << DeltaXi <<
" stopping calculation and assigning the baricenter" << std::endl;
789 rResult.
resize( 3, 2,
false );
793 CalculateShapeFunctionsIntegrationPointsLocalGradients( ThisMethod );
794 Matrix& ShapeFunctionsGradientInIntegrationPoint = shape_functions_gradients( IntegrationPointIndex );
802 ( this->
GetPoint(
i ).X() ) * ( ShapeFunctionsGradientInIntegrationPoint(
i, 0 ) );
804 ( this->
GetPoint(
i ).X() ) * ( ShapeFunctionsGradientInIntegrationPoint(
i, 1 ) );
806 ( this->
GetPoint(
i ).Y() ) * ( ShapeFunctionsGradientInIntegrationPoint(
i, 0 ) );
808 ( this->
GetPoint(
i ).Y() ) * ( ShapeFunctionsGradientInIntegrationPoint(
i, 1 ) );
810 ( this->
GetPoint(
i ).Z() ) * ( ShapeFunctionsGradientInIntegrationPoint(
i, 0 ) );
812 ( this->
GetPoint(
i ).Z() ) * ( ShapeFunctionsGradientInIntegrationPoint(
i, 1 ) );
841 const Matrix& DeltaPosition )
const override
844 rResult.
resize( 3, 2 ,
false);
848 CalculateShapeFunctionsIntegrationPointsLocalGradients( ThisMethod );
849 Matrix& ShapeFunctionsGradientInIntegrationPoint = shape_functions_gradients( IntegrationPointIndex );
857 ( this->
GetPoint(
i ).X() - DeltaPosition(
i,0) ) * ( ShapeFunctionsGradientInIntegrationPoint(
i, 0 ) );
859 ( this->
GetPoint(
i ).X() - DeltaPosition(
i,0) ) * ( ShapeFunctionsGradientInIntegrationPoint(
i, 1 ) );
861 ( this->
GetPoint(
i ).Y() - DeltaPosition(
i,1) ) * ( ShapeFunctionsGradientInIntegrationPoint(
i, 0 ) );
863 ( this->
GetPoint(
i ).Y() - DeltaPosition(
i,1) ) * ( ShapeFunctionsGradientInIntegrationPoint(
i, 1 ) );
865 ( this->
GetPoint(
i ).Z() - DeltaPosition(
i,2) ) * ( ShapeFunctionsGradientInIntegrationPoint(
i, 0 ) );
867 ( this->
GetPoint(
i ).Z() - DeltaPosition(
i,2) ) * ( ShapeFunctionsGradientInIntegrationPoint(
i, 1 ) );
891 rResult.
resize( 3, 2 ,
false );
894 Matrix shape_functions_gradients;
896 shape_functions_gradients, rPoint );
902 rResult( 0, 0 ) += ( this->
GetPoint(
i ).X() ) * ( shape_functions_gradients(
i, 0 ) );
903 rResult( 0, 1 ) += ( this->
GetPoint(
i ).X() ) * ( shape_functions_gradients(
i, 1 ) );
904 rResult( 1, 0 ) += ( this->
GetPoint(
i ).Y() ) * ( shape_functions_gradients(
i, 0 ) );
905 rResult( 1, 1 ) += ( this->
GetPoint(
i ).Y() ) * ( shape_functions_gradients(
i, 1 ) );
906 rResult( 2, 0 ) += ( this->
GetPoint(
i ).Z() ) * ( shape_functions_gradients(
i, 0 ) );
907 rResult( 2, 1 ) += ( this->
GetPoint(
i ).Z() ) * ( shape_functions_gradients(
i, 1 ) );
931 if( rResult.size() != this->IntegrationPointsNumber( ThisMethod ) )
941 this->
Jacobian( J, pnt, ThisMethod);
943 Tangent0[0] =
J(0,0);
944 Tangent0[1] =
J(1,0);
945 Tangent0[2] =
J(2,0);
947 Tangent1[0] =
J(0,1);
948 Tangent1[1] =
J(1,1);
949 Tangent1[2] =
J(2,1);
984 this->
Jacobian( J, IntegrationPointIndex, ThisMethod);
990 Tangent0[0] =
J(0,0);
991 Tangent0[1] =
J(1,0);
992 Tangent0[2] =
J(2,0);
994 Tangent1[0] =
J(0,1);
995 Tangent1[1] =
J(1,1);
996 Tangent1[2] =
J(2,1);
1232 switch ( ShapeFunctionIndex )
1235 return( 0.125*( 1.0 - rPoint[0] )*( 1.0 - rPoint[1] )*( 1.0 - rPoint[2] ) );
1237 return( 0.125*( 1.0 + rPoint[0] )*( 1.0 - rPoint[1] )*( 1.0 - rPoint[2] ) );
1239 return( 0.125*( 1.0 + rPoint[0] )*( 1.0 + rPoint[1] )*( 1.0 - rPoint[2] ) );
1241 return( 0.125*( 1.0 - rPoint[0] )*( 1.0 + rPoint[1] )*( 1.0 - rPoint[2] ) );
1243 return( 0.125*( 1.0 - rPoint[0] )*( 1.0 - rPoint[1] )*( 1.0 + rPoint[2] ) );
1245 return( 0.125*( 1.0 + rPoint[0] )*( 1.0 - rPoint[1] )*( 1.0 + rPoint[2] ) );
1247 return( 0.125*( 1.0 + rPoint[0] )*( 1.0 + rPoint[1] )*( 1.0 + rPoint[2] ) );
1249 return( 0.125*( 1.0 - rPoint[0] )*( 1.0 + rPoint[1] )*( 1.0 + rPoint[2] ) );
1251 KRATOS_ERROR <<
"Wrong index of shape function!" << *
this << std::endl;
1270 if(rResult.size() != 8) rResult.
resize(8,
false);
1271 rResult[0] = 0.125*( 1.0 - rCoordinates[0] )*( 1.0 - rCoordinates[1] )*( 1.0 - rCoordinates[2] ) ;
1272 rResult[1] = 0.125*( 1.0 + rCoordinates[0] )*( 1.0 - rCoordinates[1] )*( 1.0 - rCoordinates[2] ) ;
1273 rResult[2] = 0.125*( 1.0 + rCoordinates[0] )*( 1.0 + rCoordinates[1] )*( 1.0 - rCoordinates[2] ) ;
1274 rResult[3] = 0.125*( 1.0 - rCoordinates[0] )*( 1.0 + rCoordinates[1] )*( 1.0 - rCoordinates[2] ) ;
1275 rResult[4] = 0.125*( 1.0 - rCoordinates[0] )*( 1.0 - rCoordinates[1] )*( 1.0 + rCoordinates[2] ) ;
1276 rResult[5] = 0.125*( 1.0 + rCoordinates[0] )*( 1.0 - rCoordinates[1] )*( 1.0 + rCoordinates[2] ) ;
1277 rResult[6] = 0.125*( 1.0 + rCoordinates[0] )*( 1.0 + rCoordinates[1] )*( 1.0 + rCoordinates[2] ) ;
1278 rResult[7] = 0.125*( 1.0 - rCoordinates[0] )*( 1.0 + rCoordinates[1] )*( 1.0 + rCoordinates[2] ) ;
1292 if ( rResult.size1() != 8 || rResult.size2() != 3 )
1293 rResult.
resize( 8, 3 ,
false);
1297 rResult( 0, 0 ) = -0.125 * ( 1.0 - rPoint[1] ) * ( 1.0 - rPoint[2] );
1298 rResult( 0, 1 ) = -0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 - rPoint[2] );
1299 rResult( 0, 2 ) = -0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 - rPoint[1] );
1300 rResult( 1, 0 ) = 0.125 * ( 1.0 - rPoint[1] ) * ( 1.0 - rPoint[2] );
1301 rResult( 1, 1 ) = -0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 - rPoint[2] );
1302 rResult( 1, 2 ) = -0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 - rPoint[1] );
1303 rResult( 2, 0 ) = 0.125 * ( 1.0 + rPoint[1] ) * ( 1.0 - rPoint[2] );
1304 rResult( 2, 1 ) = 0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 - rPoint[2] );
1305 rResult( 2, 2 ) = -0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 + rPoint[1] );
1306 rResult( 3, 0 ) = -0.125 * ( 1.0 + rPoint[1] ) * ( 1.0 - rPoint[2] );
1307 rResult( 3, 1 ) = 0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 - rPoint[2] );
1308 rResult( 3, 2 ) = -0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 + rPoint[1] );
1309 rResult( 4, 0 ) = -0.125 * ( 1.0 - rPoint[1] ) * ( 1.0 + rPoint[2] );
1310 rResult( 4, 1 ) = -0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 + rPoint[2] );
1311 rResult( 4, 2 ) = 0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 - rPoint[1] );
1312 rResult( 5, 0 ) = 0.125 * ( 1.0 - rPoint[1] ) * ( 1.0 + rPoint[2] );
1313 rResult( 5, 1 ) = -0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 + rPoint[2] );
1314 rResult( 5, 2 ) = 0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 - rPoint[1] );
1315 rResult( 6, 0 ) = 0.125 * ( 1.0 + rPoint[1] ) * ( 1.0 + rPoint[2] );
1316 rResult( 6, 1 ) = 0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 + rPoint[2] );
1317 rResult( 6, 2 ) = 0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 + rPoint[1] );
1318 rResult( 7, 0 ) = -0.125 * ( 1.0 + rPoint[1] ) * ( 1.0 + rPoint[2] );
1319 rResult( 7, 1 ) = 0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 + rPoint[2] );
1320 rResult( 7, 2 ) = 0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 + rPoint[1] );
1344 const unsigned int integration_points_number =
1347 if ( integration_points_number == 0 )
1348 KRATOS_ERROR <<
"This integration method is not supported" << *
this << std::endl;
1351 if ( rResult.size() != integration_points_number )
1361 CalculateShapeFunctionsIntegrationPointsLocalGradients( ThisMethod );
1369 for (
unsigned int pnt = 0; pnt < integration_points_number; pnt++ )
1371 rResult[pnt].
resize( 8, 3,
false );
1373 for (
int i = 0;
i < 8;
i++ )
1375 for (
int j = 0;
j < 3;
j++ )
1377 rResult[pnt](
i,
j ) =
1378 ( locG[pnt](
i, 0 ) * invJ[pnt](
j, 0 ) )
1379 + ( locG[pnt](
i, 1 ) * invJ[pnt](
j, 1 ) )
1380 + ( locG[pnt](
i, 2 ) * invJ[pnt](
j, 2 ) );
1389 Vector& determinants_of_jacobian,
1395 if ( integration_points_number == 0 )
1396 KRATOS_ERROR <<
"This integration method is not supported" << *
this << std::endl;
1399 if ( rResult.size() != integration_points_number )
1407 if ( determinants_of_jacobian.size() != integration_points_number)
1408 determinants_of_jacobian.
resize(integration_points_number,
false);
1412 CalculateShapeFunctionsIntegrationPointsLocalGradients( ThisMethod );
1419 for (
unsigned int pnt = 0; pnt < integration_points_number; pnt++ )
1431 determinants_of_jacobian[pnt] = DetJ;
1433 rResult[pnt].
resize( 4, 3 ,
false);
1435 for (
int i = 0;
i < 4;
i++ )
1437 for (
int j = 0;
j < 3;
j++ )
1439 rResult[pnt](
i,
j ) =
1443 ( locG[pnt](
i, 0 ) * invJ( 0,
j ) )
1444 + ( locG[pnt](
i, 1 ) * invJ( 1,
j ) )
1445 + ( locG[pnt](
i, 2 ) * invJ( 2,
j ) );
1465 return "3 dimensional hexahedra with eight nodes in 3D space";
1477 rOStream <<
"3 dimensional hexahedra with eight nodes in 3D space";
1493 std::cout << std::endl;
1499 rOStream <<
" Jacobian\t : " << jacobian;
1517 void save(
Serializer& rSerializer )
const override
1547 result( 0, 0 ) = -0.125 * ( 1.0 - rPoint[1] ) * ( 1.0 - rPoint[2] );
1548 result( 0, 1 ) = -0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 - rPoint[2] );
1549 result( 0, 2 ) = -0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 - rPoint[1] );
1550 result( 1, 0 ) = 0.125 * ( 1.0 - rPoint[1] ) * ( 1.0 - rPoint[2] );
1551 result( 1, 1 ) = -0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 - rPoint[2] );
1552 result( 1, 2 ) = -0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 - rPoint[1] );
1553 result( 2, 0 ) = 0.125 * ( 1.0 + rPoint[1] ) * ( 1.0 - rPoint[2] );
1554 result( 2, 1 ) = 0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 - rPoint[2] );
1555 result( 2, 2 ) = -0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 + rPoint[1] );
1556 result( 3, 0 ) = -0.125 * ( 1.0 + rPoint[1] ) * ( 1.0 - rPoint[2] );
1557 result( 3, 1 ) = 0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 - rPoint[2] );
1558 result( 3, 2 ) = -0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 + rPoint[1] );
1559 result( 4, 0 ) = -0.125 * ( 1.0 - rPoint[1] ) * ( 1.0 + rPoint[2] );
1560 result( 4, 1 ) = -0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 + rPoint[2] );
1561 result( 4, 2 ) = 0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 - rPoint[1] );
1562 result( 5, 0 ) = 0.125 * ( 1.0 - rPoint[1] ) * ( 1.0 + rPoint[2] );
1563 result( 5, 1 ) = -0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 + rPoint[2] );
1564 result( 5, 2 ) = 0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 - rPoint[1] );
1565 result( 6, 0 ) = 0.125 * ( 1.0 + rPoint[1] ) * ( 1.0 + rPoint[2] );
1566 result( 6, 1 ) = 0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 + rPoint[2] );
1567 result( 6, 2 ) = 0.125 * ( 1.0 + rPoint[0] ) * ( 1.0 + rPoint[1] );
1568 result( 7, 0 ) = -0.125 * ( 1.0 + rPoint[1] ) * ( 1.0 + rPoint[2] );
1569 result( 7, 1 ) = 0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 + rPoint[2] );
1570 result( 7, 2 ) = 0.125 * ( 1.0 - rPoint[0] ) * ( 1.0 + rPoint[1] );
1586 static Matrix CalculateShapeFunctionsIntegrationPointsValues(
1593 const int integration_points_number = integration_points.size();
1595 const int points_number = 8;
1597 Matrix shape_function_values( integration_points_number, points_number );
1600 for (
int pnt = 0; pnt < integration_points_number; pnt++ )
1602 shape_function_values( pnt, 0 ) =
1603 0.125 * ( 1.0 - integration_points[pnt].X() )
1604 * ( 1.0 - integration_points[pnt].
Y() )
1605 * ( 1.0 - integration_points[pnt].
Z() );
1606 shape_function_values( pnt, 1 ) =
1607 0.125 * ( 1.0 + integration_points[pnt].X() )
1608 * ( 1.0 - integration_points[pnt].
Y() )
1609 * ( 1.0 - integration_points[pnt].
Z() );
1610 shape_function_values( pnt, 2 ) =
1611 0.125 * ( 1.0 + integration_points[pnt].X() )
1612 * ( 1.0 + integration_points[pnt].
Y() )
1613 * ( 1.0 - integration_points[pnt].
Z() );
1614 shape_function_values( pnt, 3 ) =
1615 0.125 * ( 1.0 - integration_points[pnt].X() )
1616 * ( 1.0 + integration_points[pnt].
Y() )
1617 * ( 1.0 - integration_points[pnt].
Z() );
1618 shape_function_values( pnt, 4 ) =
1619 0.125 * ( 1.0 - integration_points[pnt].X() )
1620 * ( 1.0 - integration_points[pnt].
Y() )
1621 * ( 1.0 + integration_points[pnt].
Z() );
1622 shape_function_values( pnt, 5 ) =
1623 0.125 * ( 1.0 + integration_points[pnt].X() )
1624 * ( 1.0 - integration_points[pnt].
Y() )
1625 * ( 1.0 + integration_points[pnt].
Z() );
1626 shape_function_values( pnt, 6 ) =
1627 0.125 * ( 1.0 + integration_points[pnt].X() )
1628 * ( 1.0 + integration_points[pnt].
Y() )
1629 * ( 1.0 + integration_points[pnt].
Z() );
1630 shape_function_values( pnt, 7 ) =
1631 0.125 * ( 1.0 - integration_points[pnt].X() )
1632 * ( 1.0 + integration_points[pnt].
Y() )
1633 * ( 1.0 + integration_points[pnt].
Z() );
1636 return shape_function_values;
1655 AllIntegrationPoints();
1658 const int integration_points_number = integration_points.size();
1663 for (
int pnt = 0; pnt < integration_points_number; pnt++ )
1665 Matrix& result = d_shape_f_values[pnt];
1668 -0.125 * ( 1.0 - integration_points[pnt].Y() )
1669 * ( 1.0 - integration_points[pnt].
Z() );
1671 -0.125 * ( 1.0 - integration_points[pnt].X() )
1672 * ( 1.0 - integration_points[pnt].
Z() );
1674 -0.125 * ( 1.0 - integration_points[pnt].X() )
1675 * ( 1.0 - integration_points[pnt].
Y() );
1677 0.125 * ( 1.0 - integration_points[pnt].Y() )
1678 * ( 1.0 - integration_points[pnt].
Z() );
1680 -0.125 * ( 1.0 + integration_points[pnt].X() )
1681 * ( 1.0 - integration_points[pnt].
Z() );
1683 -0.125 * ( 1.0 + integration_points[pnt].X() )
1684 * ( 1.0 - integration_points[pnt].
Y() );
1686 0.125 * ( 1.0 + integration_points[pnt].Y() )
1687 * ( 1.0 - integration_points[pnt].
Z() );
1689 0.125 * ( 1.0 + integration_points[pnt].X() )
1690 * ( 1.0 - integration_points[pnt].
Z() );
1692 -0.125 * ( 1.0 + integration_points[pnt].X() )
1693 * ( 1.0 + integration_points[pnt].
Y() );
1695 -0.125 * ( 1.0 + integration_points[pnt].Y() )
1696 * ( 1.0 - integration_points[pnt].
Z() );
1698 0.125 * ( 1.0 - integration_points[pnt].X() )
1699 * ( 1.0 - integration_points[pnt].
Z() );
1701 -0.125 * ( 1.0 - integration_points[pnt].X() )
1702 * ( 1.0 + integration_points[pnt].
Y() );
1704 -0.125 * ( 1.0 - integration_points[pnt].Y() )
1705 * ( 1.0 + integration_points[pnt].
Z() );
1707 -0.125 * ( 1.0 - integration_points[pnt].X() )
1708 * ( 1.0 + integration_points[pnt].
Z() );
1710 0.125 * ( 1.0 - integration_points[pnt].X() )
1711 * ( 1.0 - integration_points[pnt].
Y() );
1713 0.125 * ( 1.0 - integration_points[pnt].Y() )
1714 * ( 1.0 + integration_points[pnt].
Z() );
1716 -0.125 * ( 1.0 + integration_points[pnt].X() )
1717 * ( 1.0 + integration_points[pnt].
Z() );
1719 0.125 * ( 1.0 + integration_points[pnt].X() )
1720 * ( 1.0 - integration_points[pnt].
Y() );
1722 0.125 * ( 1.0 + integration_points[pnt].Y() )
1723 * ( 1.0 + integration_points[pnt].
Z() );
1725 0.125 * ( 1.0 + integration_points[pnt].X() )
1726 * ( 1.0 + integration_points[pnt].
Z() );
1728 0.125 * ( 1.0 + integration_points[pnt].X() )
1729 * ( 1.0 + integration_points[pnt].
Y() );
1731 -0.125 * ( 1.0 + integration_points[pnt].Y() )
1732 * ( 1.0 + integration_points[pnt].
Z() );
1734 0.125 * ( 1.0 - integration_points[pnt].X() )
1735 * ( 1.0 + integration_points[pnt].
Z() );
1737 0.125 * ( 1.0 - integration_points[pnt].X() )
1738 * ( 1.0 + integration_points[pnt].
Y() );
1741 return d_shape_f_values;
1749 Quadrature < HexahedronGaussLobattoIntegrationPoints1,
1750 3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
1751 Quadrature < HexahedronGaussLobattoIntegrationPoints2,
1752 3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
1757 return integration_points;
1765 HexahedraInterface3D8<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1767 HexahedraInterface3D8<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1773 return shape_functions_values;
1780 AllShapeFunctionsLocalGradients()
1785 HexahedraInterface3D8<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1787 HexahedraInterface3D8<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1793 return shape_functions_local_gradients;
1828 rOStream << std::endl;
1835 template<
class TPo
intType>
const
1836 GeometryData HexahedraInterface3D8<TPointType>::msGeometryData(
1839 HexahedraInterface3D8<TPointType>::AllIntegrationPoints(),
1840 HexahedraInterface3D8<TPointType>::AllShapeFunctionsValues(),
1841 AllShapeFunctionsLocalGradients()
1844 template<
class TPo
intType>
Definition: geometry_data.h:60
KratosGeometryType
Definition: geometry_data.h:110
SizeType IntegrationPointsNumber() const
Definition: geometry_data.h:430
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
virtual CoordinatesArrayType & GlobalCoordinates(CoordinatesArrayType &rResult, CoordinatesArrayType const &LocalCoordinates) const
Definition: geometry.h:2377
std::vector< IntegrationPointType > IntegrationPointsArrayType
Definition: geometry.h:161
virtual array_1d< double, 3 > Normal(const CoordinatesArrayType &rPointLocalCoordinates) const
It returns a vector that is normal to its corresponding geometry in the given local point.
Definition: geometry.h:1543
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
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
const ShapeFunctionsGradientsType & ShapeFunctionsLocalGradients() const
Definition: geometry.h:3539
JacobiansType & Jacobian(JacobiansType &rResult) const
Definition: geometry.h:2901
TPointType const & GetPoint(const int Index) const
Definition: geometry.h:1816
GeometryData::ShapeFunctionsLocalGradientsContainerType ShapeFunctionsLocalGradientsContainerType
Definition: geometry.h:177
SizeType IntegrationPointsNumber() const
Definition: geometry.h:2257
GeometryData::IntegrationMethod IntegrationMethod
Definition: geometry.h:122
Definition: hexahedra_interface_3d_8.h:36
BaseType::IndexType IndexType
Definition: hexahedra_interface_3d_8.h:79
Matrix & PointsLocalCoordinates(Matrix &rResult) const override
Definition: hexahedra_interface_3d_8.h:606
HexahedraInterface3D8(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)
Definition: hexahedra_interface_3d_8.h:286
double Volume() const override
This method calculate and return volume of this geometry.
Definition: hexahedra_interface_3d_8.h:577
SizeType EdgesNumber() const override
This method gives you number of all edges of this geometry.
Definition: hexahedra_interface_3d_8.h:1098
double ShapeFunctionValue(IndexType ShapeFunctionIndex, const CoordinatesArrayType &rPoint) const override
Definition: hexahedra_interface_3d_8.h:1229
double DeterminantOfJacobian(IndexType IntegrationPointIndex, IntegrationMethod ThisMethod) const override
Definition: hexahedra_interface_3d_8.h:980
BaseType::CoordinatesArrayType CoordinatesArrayType
Definition: hexahedra_interface_3d_8.h:153
BaseType::ShapeFunctionsLocalGradientsContainerType ShapeFunctionsLocalGradientsContainerType
Definition: hexahedra_interface_3d_8.h:129
double DomainSize() const override
Definition: hexahedra_interface_3d_8.h:596
BaseType::IntegrationPointsContainerType IntegrationPointsContainerType
Definition: hexahedra_interface_3d_8.h:115
BaseType::Pointer Create(PointsArrayType const &ThisPoints) const override
Definition: hexahedra_interface_3d_8.h:477
BaseType::SizeType SizeType
Definition: hexahedra_interface_3d_8.h:87
GeometryData::KratosGeometryFamily GetGeometryFamily() const override
Definition: hexahedra_interface_3d_8.h:422
BaseType::GeometriesArrayType GeometriesArrayType
Definition: hexahedra_interface_3d_8.h:67
void PrintInfo(std::ostream &rOStream) const override
Definition: hexahedra_interface_3d_8.h:1475
double Area() const override
Definition: hexahedra_interface_3d_8.h:520
SizeType FacesNumber() const override
Returns the number of faces of the current geometry.
Definition: hexahedra_interface_3d_8.h:1165
Matrix & InverseOfJacobian(Matrix &rResult, const CoordinatesArrayType &rPoint) const override
Definition: hexahedra_interface_3d_8.h:1077
BaseType::IntegrationPointsArrayType IntegrationPointsArrayType
Definition: hexahedra_interface_3d_8.h:108
HexahedraInterface3D8 & operator=(HexahedraInterface3D8< TOtherPointType > const &rOther)
Definition: hexahedra_interface_3d_8.h:465
Matrix & ShapeFunctionsLocalGradients(Matrix &rResult, const CoordinatesArrayType &rPoint) const override
Definition: hexahedra_interface_3d_8.h:1290
HexahedraInterface3D8(const PointsArrayType &ThisPoints)
Definition: hexahedra_interface_3d_8.h:381
HexahedraInterface3D8(HexahedraInterface3D8 const &rOther)
Definition: hexahedra_interface_3d_8.h:397
BaseType::PointsArrayType PointsArrayType
Definition: hexahedra_interface_3d_8.h:93
virtual ~HexahedraInterface3D8()
Destructor. Does nothing!!!
Definition: hexahedra_interface_3d_8.h:420
CoordinatesArrayType & PointLocalCoordinates(CoordinatesArrayType &rResult, const CoordinatesArrayType &rPoint) const override
Returns the local coordinates of a given arbitrary point.
Definition: hexahedra_interface_3d_8.h:682
GeometriesArrayType GenerateFaces() const override
Returns all faces of the current geometry.
Definition: hexahedra_interface_3d_8.h:1178
Vector & DeterminantOfJacobian(Vector &rResult, IntegrationMethod ThisMethod) const override
Definition: hexahedra_interface_3d_8.h:928
TPointType PointType
Definition: hexahedra_interface_3d_8.h:72
void ShapeFunctionsIntegrationPointsGradients(ShapeFunctionsGradientsType &rResult, IntegrationMethod ThisMethod) const override
Definition: hexahedra_interface_3d_8.h:1340
Geometry< TPointType > BaseType
Definition: hexahedra_interface_3d_8.h:45
Matrix & Jacobian(Matrix &rResult, const CoordinatesArrayType &rPoint) const override
Definition: hexahedra_interface_3d_8.h:888
BaseType::IntegrationPointType IntegrationPointType
Definition: hexahedra_interface_3d_8.h:100
GeometriesArrayType GenerateEdges() const override
This method gives you all edges of this geometry.
Definition: hexahedra_interface_3d_8.h:1111
GeometryData::IntegrationMethod IntegrationMethod
Definition: hexahedra_interface_3d_8.h:61
HexahedraInterface3D8 & operator=(const HexahedraInterface3D8 &rOther)
Definition: hexahedra_interface_3d_8.h:447
Line3D2< TPointType > EdgeType
Definition: hexahedra_interface_3d_8.h:50
Vector & ShapeFunctionsValues(Vector &rResult, const CoordinatesArrayType &rCoordinates) const override
Definition: hexahedra_interface_3d_8.h:1268
std::string Info() const override
Definition: hexahedra_interface_3d_8.h:1463
BaseType::ShapeFunctionsValuesContainerType ShapeFunctionsValuesContainerType
Definition: hexahedra_interface_3d_8.h:122
BaseType::JacobiansType JacobiansType
Definition: hexahedra_interface_3d_8.h:136
bool IsInside(const CoordinatesArrayType &rPoint, CoordinatesArrayType &rResult, const double Tolerance=std::numeric_limits< double >::epsilon()) const override
Definition: hexahedra_interface_3d_8.h:654
KRATOS_CLASS_POINTER_DEFINITION(HexahedraInterface3D8)
double Length() const override
Definition: hexahedra_interface_3d_8.h:501
Matrix & Jacobian(Matrix &rResult, IndexType IntegrationPointIndex, IntegrationMethod ThisMethod, const Matrix &DeltaPosition) const override
Definition: hexahedra_interface_3d_8.h:838
Matrix MatrixType
Definition: hexahedra_interface_3d_8.h:158
JacobiansType & InverseOfJacobian(JacobiansType &rResult, IntegrationMethod ThisMethod) const override
Definition: hexahedra_interface_3d_8.h:1022
Matrix & InverseOfJacobian(Matrix &rResult, IndexType IntegrationPointIndex, IntegrationMethod ThisMethod) const override
Definition: hexahedra_interface_3d_8.h:1053
Quadrilateral3D4< TPointType > FaceType
Definition: hexahedra_interface_3d_8.h:51
BaseType::NormalType NormalType
Definition: hexahedra_interface_3d_8.h:148
BaseType::ShapeFunctionsGradientsType ShapeFunctionsGradientsType
Definition: hexahedra_interface_3d_8.h:143
void ShapeFunctionsIntegrationPointsGradients(ShapeFunctionsGradientsType &rResult, Vector &determinants_of_jacobian, IntegrationMethod ThisMethod) const override
Definition: hexahedra_interface_3d_8.h:1387
GeometryData::KratosGeometryType GetGeometryType() const override
Definition: hexahedra_interface_3d_8.h:427
friend class HexahedraInterface3D8
Definition: hexahedra_interface_3d_8.h:1801
void PrintData(std::ostream &rOStream) const override
Definition: hexahedra_interface_3d_8.h:1489
HexahedraInterface3D8(HexahedraInterface3D8< TOtherPointType > const &rOther)
Definition: hexahedra_interface_3d_8.h:414
Matrix & Jacobian(Matrix &rResult, IndexType IntegrationPointIndex, IntegrationMethod ThisMethod) const override
Definition: hexahedra_interface_3d_8.h:784
Short class definition.
Definition: integration_point.h:52
Definition: amatrix_interface.h:41
void swap(Matrix &Other)
Definition: amatrix_interface.h:289
void resize(std::size_t NewSize1, std::size_t NewSize2, bool preserve=0)
Definition: amatrix_interface.h:224
An two node 3D line geometry with linear shape functions.
Definition: line_3d_2.h:64
static void InvertMatrix3(const TMatrix1 &rInputMatrix, TMatrix2 &rInvertedMatrix, double &rInputMatrixDet)
It inverts matrices of order 3.
Definition: math_utils.h:449
static double Norm3(const TVectorType &a)
Calculates the norm of vector "a" which is assumed to be of size 3.
Definition: math_utils.h:691
static T CrossProduct(const T &a, const T &b)
Performs the vector product of the two input vectors a,b.
Definition: math_utils.h:762
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
void push_back(const TPointerType &x)
Definition: pointer_vector.h:270
A four node 3D quadrilateral geometry with bi-linear shape functions.
Definition: quadrilateral_3d_4.h:76
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
Short class definition.
Definition: array_1d.h:61
BOOST_UBLAS_INLINE void clear()
Definition: array_1d.h:325
#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
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
TExpressionType::data_type norm_2(AMatrix::MatrixExpression< TExpressionType, TCategory > const &TheExpression)
Definition: amatrix_interface.h:625
KratosZeroMatrix< double > ZeroMatrix
Definition: amatrix_interface.h:559
const GeometryData HexahedraInterface3D8< TPointType >::msGeometryData & msGeometryDimension(), HexahedraInterface3D8< TPointType >::AllShapeFunctionsValues(), AllShapeFunctionsLocalGradients()
Definition: brep_curve.h:483
Internals::Matrix< double, AMatrix::dynamic, AMatrix::dynamic > Matrix
Definition: amatrix_interface.h:470
AMatrix::MatrixProductExpression< TExpression1Type, TExpression2Type > prod(AMatrix::MatrixExpression< TExpression1Type, TCategory1 > const &First, AMatrix::MatrixExpression< TExpression2Type, TCategory2 > const &Second)
Definition: amatrix_interface.h:568
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
const GeometryData BrepCurve< TContainerPointType, TContainerPointEmbeddedType >::msGeometryData & msGeometryDimension
Definition: brep_curve.h:483
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
AMatrix::TransposeMatrix< const T > trans(const T &TheMatrix)
Definition: amatrix_interface.h:486
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
res
Definition: generate_convection_diffusion_explicit_element.py:211
w
Definition: generate_convection_diffusion_explicit_element.py:108
DN
Definition: generate_convection_diffusion_explicit_element.py:98
int tol
Definition: hinsberg_optimization.py:138
def load(f)
Definition: ode_solve.py:307
int k
Definition: quadrature.py:595
int j
Definition: quadrature.py:648
float temp
Definition: rotating_cone.py:85
J
Definition: sensitivityMatrix.py:58
integer i
Definition: TensorModule.f:17