KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
prism_3d_15.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Riccardo Rossi
11 // Janosch Stascheit
12 // Felix Nagel
13 // contributors: Hoang Giang Bui
14 // Josep Maria Carbonell
15 // Vicente Mataix Ferrandiz
16 //
17 
18 #pragma once
19 
20 // System includes
21 
22 // External includes
23 
24 // Project includes
29 
30 namespace Kratos
31 {
56 template<class TPointType>
57 class Prism3D15
58  : public Geometry<TPointType>
59 {
60 public:
69 
76 
81 
86 
92 
96  typedef TPointType PointType;
97 
103  typedef typename BaseType::IndexType IndexType;
104 
105 
111  typedef typename BaseType::SizeType SizeType;
112 
118 
125 
133 
140 
147 
154 
161 
168 
173 
178 
183 
184 
188  Prism3D15( typename PointType::Pointer pPoint1,
189  typename PointType::Pointer pPoint2,
190  typename PointType::Pointer pPoint3,
191  typename PointType::Pointer pPoint4,
192  typename PointType::Pointer pPoint5,
193  typename PointType::Pointer pPoint6,
194  typename PointType::Pointer pPoint7,
195  typename PointType::Pointer pPoint8,
196  typename PointType::Pointer pPoint9,
197  typename PointType::Pointer pPoint10,
198  typename PointType::Pointer pPoint11,
199  typename PointType::Pointer pPoint12,
200  typename PointType::Pointer pPoint13,
201  typename PointType::Pointer pPoint14,
202  typename PointType::Pointer pPoint15 )
203  : BaseType( PointsArrayType(), &msGeometryData )
204  {
205  this->Points().reserve( 15 );
206  this->Points().push_back( pPoint1 );
207  this->Points().push_back( pPoint2 );
208  this->Points().push_back( pPoint3 );
209  this->Points().push_back( pPoint4 );
210  this->Points().push_back( pPoint5 );
211  this->Points().push_back( pPoint6 );
212  this->Points().push_back( pPoint7 );
213  this->Points().push_back( pPoint8 );
214  this->Points().push_back( pPoint9 );
215  this->Points().push_back( pPoint10 );
216  this->Points().push_back( pPoint11 );
217  this->Points().push_back( pPoint12 );
218  this->Points().push_back( pPoint13 );
219  this->Points().push_back( pPoint14 );
220  this->Points().push_back( pPoint15 );
221  }
222 
223  Prism3D15( const PointsArrayType& rThisPoints )
224  : BaseType( rThisPoints, &msGeometryData )
225  {
226  KRATOS_ERROR_IF( this->PointsNumber() != 15 ) << "Invalid points number. Expected 15, given " << this->PointsNumber() << std::endl;
227  }
228 
230  explicit Prism3D15(
231  const IndexType GeometryId,
232  const PointsArrayType& rThisPoints
233  ) : BaseType(GeometryId, rThisPoints, &msGeometryData)
234  {
235  KRATOS_ERROR_IF( this->PointsNumber() != 15 ) << "Invalid points number. Expected 15, given " << this->PointsNumber() << std::endl;
236  }
237 
239  explicit Prism3D15(
240  const std::string& rGeometryName,
241  const PointsArrayType& rThisPoints
242  ) : BaseType(rGeometryName, rThisPoints, &msGeometryData)
243  {
244  KRATOS_ERROR_IF(this->PointsNumber() != 15) << "Invalid points number. Expected 15, given " << this->PointsNumber() << std::endl;
245  }
246 
256  Prism3D15( Prism3D15 const& rOther )
257  : BaseType( rOther )
258  {
259  }
260 
273  template<class TOtherPointType> Prism3D15( Prism3D15<TOtherPointType> const& rOther )
274  : BaseType( rOther )
275  {
276  }
277 
279  ~Prism3D15() override {}
280 
282  {
284  }
285 
287  {
289  }
290 
306  Prism3D15& operator=( const Prism3D15& rOther )
307  {
308  BaseType::operator=( rOther );
309  return *this;
310  }
311 
323  template<class TOtherPointType>
325  {
326  BaseType::operator=( rOther );
327 
328  return *this;
329  }
330 
331 
342  typename BaseType::Pointer Create(
343  const IndexType NewGeometryId,
344  PointsArrayType const& rThisPoints
345  ) const override
346  {
347  return typename BaseType::Pointer( new Prism3D15( NewGeometryId, rThisPoints ) );
348  }
349 
356  typename BaseType::Pointer Create(
357  const IndexType NewGeometryId,
358  const BaseType& rGeometry
359  ) const override
360  {
361  auto p_geometry = typename BaseType::Pointer( new Prism3D15( NewGeometryId, rGeometry.Points() ) );
362  p_geometry->SetData(rGeometry.GetData());
363  return p_geometry;
364  }
365 
378  double DomainSize() const override
379  {
380  return this->Volume();
381  }
382 
391  double Volume() const override
392  {
394  }
395 
396  double Area() const override
397  {
398  return std::abs( this->DeterminantOfJacobian( PointType() ) ) * 0.5;
399  }
400 
401  double Length() const override
402  {
403  const double volume = Volume();
404  return std::pow(volume, 1.0/3.0)/3.0;
405  }
406 
407 
413  Matrix& PointsLocalCoordinates( Matrix& rResult ) const override
414  {
415  if ( rResult.size1() != 15 || rResult.size2() != 3 )
416  rResult.resize( 15, 3, false );
417 
418  rResult( 0, 0 ) = 0.0; //corners of bottom base
419  rResult( 0, 1 ) = 0.0;
420  rResult( 0, 2 ) = 0.0;
421 
422  rResult( 1, 0 ) = 1.0;
423  rResult( 1, 1 ) = 0.0;
424  rResult( 1, 2 ) = 0.0;
425 
426  rResult( 2, 0 ) = 0.0;
427  rResult( 2, 1 ) = 1.0;
428  rResult( 2, 2 ) = 0.0;
429 
430  rResult( 3, 0 ) = 0.0; //corners of top base
431  rResult( 3, 1 ) = 0.0;
432  rResult( 3, 2 ) = 1.0;
433 
434  rResult( 4, 0 ) = 1.0;
435  rResult( 4, 1 ) = 0.0;
436  rResult( 4, 2 ) = 1.0;
437 
438  rResult( 5, 0 ) = 0.0;
439  rResult( 5, 1 ) = 1.0;
440  rResult( 5, 2 ) = 1.0;
441 
442  rResult( 6, 0 ) = 0.5; //mid point of bottom base
443  rResult( 6, 1 ) = 0.0;
444  rResult( 6, 2 ) = -1.0;
445 
446  rResult( 7, 0 ) = 0.5;
447  rResult( 7, 1 ) = 0.5;
448  rResult( 7, 2 ) = -1.0;
449 
450  rResult( 8, 0 ) = 0.0;
451  rResult( 8, 1 ) = 0.5;
452  rResult( 8, 2 ) = -1.0;
453 
454  rResult( 9, 0 ) = 0.0; //veritices of mid height
455  rResult( 9, 1 ) = 0.0;
456  rResult( 9, 2 ) = 0.5;
457 
458  rResult( 10, 0 ) = 1.0;
459  rResult( 10, 1 ) = 0.0;
460  rResult( 10, 2 ) = 0.5;
461 
462  rResult( 11, 0 ) = 0.0;
463  rResult( 11, 1 ) = 1.0;
464  rResult( 11, 2 ) = 0.5;
465 
466  rResult( 12, 0 ) = 0.5; //mid point of top base
467  rResult( 12, 1 ) = 0.0;
468  rResult( 12, 2 ) = 1.0;
469 
470  rResult( 13, 0 ) = 0.5;
471  rResult( 13, 1 ) = 0.5;
472  rResult( 13, 2 ) = 1.0;
473 
474  rResult( 14, 0 ) = 0.0;
475  rResult( 14, 1 ) = 0.5;
476  rResult( 14, 2 ) = 1.0;
477 
478  return rResult;
479  }
480 
489  bool IsInside(
490  const CoordinatesArrayType& rPoint,
491  CoordinatesArrayType& rResult,
492  const double Tolerance = std::numeric_limits<double>::epsilon()
493  ) const override
494  {
495  this->PointLocalCoordinates( rResult, rPoint );
496 
497  if ( (rResult[0] >= (0.0 - Tolerance)) && (rResult[0] <= (1.0 + Tolerance)) )
498  if ( (rResult[1] >= (0.0 - Tolerance)) && (rResult[1] <= (1.0 + Tolerance)) )
499  if ( (rResult[2] >= (-1.0 - Tolerance)) && (rResult[2] <= (1.0 + Tolerance)) )
500  if ((( 1.0 - ( rResult[0] + rResult[1] ) ) >= (0.0 - Tolerance) ) && (( 1.0 - ( rResult[0] + rResult[1] ) ) <= (1.0 + Tolerance) ) )
501  return true;
502 
503  return false;
504  }
505 
509 
521  SizeType EdgesNumber() const override
522  {
523  return 9;
524  }
525 
535  {
537  typedef typename Geometry<TPointType>::Pointer EdgePointerType;
538  edges.push_back( EdgePointerType( new EdgeType(
539  this->pGetPoint( 0 ),
540  this->pGetPoint( 1 ),
541  this->pGetPoint( 6 ) ) ) );
542  edges.push_back( EdgePointerType( new EdgeType(
543  this->pGetPoint( 1 ),
544  this->pGetPoint( 2 ),
545  this->pGetPoint( 7 ) ) ) );
546  edges.push_back( EdgePointerType( new EdgeType(
547  this->pGetPoint( 2 ),
548  this->pGetPoint( 0 ),
549  this->pGetPoint( 8 ) ) ) );
550  edges.push_back( EdgePointerType( new EdgeType(
551  this->pGetPoint( 3 ),
552  this->pGetPoint( 4 ),
553  this->pGetPoint( 12 ) ) ) );
554  edges.push_back( EdgePointerType( new EdgeType(
555  this->pGetPoint( 4 ),
556  this->pGetPoint( 5 ),
557  this->pGetPoint( 13 ) ) ) );
558  edges.push_back( EdgePointerType( new EdgeType(
559  this->pGetPoint( 5 ),
560  this->pGetPoint( 3 ),
561  this->pGetPoint( 14 ) ) ) );
562  edges.push_back( EdgePointerType( new EdgeType(
563  this->pGetPoint( 0 ),
564  this->pGetPoint( 3 ),
565  this->pGetPoint( 9 ) ) ) );
566  edges.push_back( EdgePointerType( new EdgeType(
567  this->pGetPoint( 1 ),
568  this->pGetPoint( 4 ),
569  this->pGetPoint( 10 ) ) ) );
570  edges.push_back( EdgePointerType( new EdgeType(
571  this->pGetPoint( 2 ),
572  this->pGetPoint( 5 ),
573  this->pGetPoint( 11 ) ) ) );
574  return edges;
575  }
576 
580 
588  SizeType FacesNumber() const override
589  {
590  return 5;
591  }
592 
602  {
604  typedef typename Geometry<TPointType>::Pointer FacePointerType;
605  faces.push_back( FacePointerType( new FaceType1(
606  this->pGetPoint( 0 ),
607  this->pGetPoint( 2 ),
608  this->pGetPoint( 1 ),
609  this->pGetPoint( 8 ),
610  this->pGetPoint( 7 ),
611  this->pGetPoint( 6 ) ) ) );
612  faces.push_back( FacePointerType( new FaceType1(
613  this->pGetPoint( 3 ),
614  this->pGetPoint( 4 ),
615  this->pGetPoint( 5 ),
616  this->pGetPoint( 12 ),
617  this->pGetPoint( 13 ),
618  this->pGetPoint( 14 ) ) ) );
619  faces.push_back( FacePointerType( new FaceType2(
620  this->pGetPoint( 0 ),
621  this->pGetPoint( 1 ),
622  this->pGetPoint( 4 ),
623  this->pGetPoint( 3 ),
624  this->pGetPoint( 6 ),
625  this->pGetPoint( 10 ),
626  this->pGetPoint( 12 ),
627  this->pGetPoint( 9 ) ) ) );
628  faces.push_back( FacePointerType( new FaceType2(
629  this->pGetPoint( 2 ),
630  this->pGetPoint( 0 ),
631  this->pGetPoint( 3 ),
632  this->pGetPoint( 5 ),
633  this->pGetPoint( 8 ),
634  this->pGetPoint( 9 ),
635  this->pGetPoint( 14 ),
636  this->pGetPoint( 11 ) ) ) );
637  faces.push_back( FacePointerType( new FaceType2(
638  this->pGetPoint( 1 ),
639  this->pGetPoint( 2 ),
640  this->pGetPoint( 5 ),
641  this->pGetPoint( 4 ),
642  this->pGetPoint( 7 ),
643  this->pGetPoint( 11 ),
644  this->pGetPoint( 13 ),
645  this->pGetPoint( 10 ) ) ) );
646  return faces;
647  }
648 
663  IndexType ShapeFunctionIndex,
664  const CoordinatesArrayType& rPoint
665  ) const override
666  {
667  return CalculateShapeFunctionValue(ShapeFunctionIndex, rPoint);
668  }
669 
671 
682  Matrix& ShapeFunctionsLocalGradients( Matrix& rResult, const CoordinatesArrayType& rPoint ) const override
683  {
684  return CalculateShapeFunctionsLocalGradients( rResult, rPoint );
685  }
686 
698  std::string Info() const override
699  {
700  return "3 dimensional prism with fiftheen nodes in 3D space";
701  }
702 
710  void PrintInfo( std::ostream& rOStream ) const override
711  {
712  rOStream << "3 dimensional prism with fifthteen nodes in 3D space";
713  }
714 
724  void PrintData( std::ostream& rOStream ) const override
725  {
726  // Base Geometry class PrintData call
727  BaseType::PrintData( rOStream );
728  std::cout << std::endl;
729 
730  // If the geometry has valid points, calculate and output its data
731  if (this->AllPointsAreValid()) {
732  Matrix jacobian;
733  this->Jacobian( jacobian, PointType() );
734  rOStream << " Jacobian in the origin\t : " << jacobian;
735  }
736  }
737 
738 private:
741 
742  static const GeometryData msGeometryData;
743 
744  static const GeometryDimension msGeometryDimension;
745 
749 
750  friend class Serializer;
751 
752  void save( Serializer& rSerializer ) const override
753  {
755  }
756 
757  void load( Serializer& rSerializer ) override
758  {
760  }
761 
762  Prism3D15(): BaseType( PointsArrayType(), &msGeometryData ) {}
763 
764 
778  static double CalculateShapeFunctionValue(
779  const IndexType ShapeFunctionIndex,
780  const CoordinatesArrayType& rPoint
781  )
782  {
783  const double x = rPoint[0];
784  const double y = rPoint[1];
785  const double z = rPoint[2];
786  switch ( ShapeFunctionIndex )
787  {
788  case 0 : return (1.0/2.0)*(2*z - 2.0)*(2*z - 1)*(-2.0*x - 2.0*y + 1.0)*(-x - y + 1.0) ;
789  case 1 : return (1.0/2.0)*x*(2.0*x - 1.0)*(2*z - 2.0)*(2*z - 1) ;
790  case 2 : return (1.0/2.0)*y*(2.0*y - 1.0)*(2*z - 2.0)*(2*z - 1) ;
791  case 3 : return z*(2*z - 1)*(-2.0*x - 2.0*y + 1.0)*(-x - y + 1.0) ;
792  case 4 : return x*z*(2.0*x - 1.0)*(2*z - 1) ;
793  case 5 : return y*z*(2.0*y - 1.0)*(2*z - 1) ;
794  case 6 : return (1.0/2.0)*x*(2*z - 2.0)*(2*z - 1)*(-4.0*x - 4.0*y + 4.0) ;
795  case 7 : return 2.0*x*y*(2*z - 2.0)*(2*z - 1) ;
796  case 8 : return 2.0*y*(2*z - 2.0)*(2*z - 1)*(-x - y + 1.0) ;
797  case 9 : return (1.0 - std::pow(2*z - 1, 2))*(-x - y + 1.0) ;
798  case 10 : return x*(1.0 - std::pow(2*z - 1, 2)) ;
799  case 11 : return y*(1.0 - std::pow(2*z - 1, 2)) ;
800  case 12 : return x*z*(2*z - 1)*(-4.0*x - 4.0*y + 4.0) ;
801  case 13 : return 4.0*x*y*z*(2*z - 1) ;
802  case 14 : return 4.0*y*z*(2*z - 1)*(-x - y + 1.0) ;
803  default:
804  KRATOS_ERROR << "Wrong index of shape function!" << ShapeFunctionIndex << std::endl;
805  }
806  }
807 
817  Vector& rResult,
818  const CoordinatesArrayType& rCoordinates
819  ) const override
820  {
821  const double x = rCoordinates[0];
822  const double y = rCoordinates[1];
823  const double z = rCoordinates[2];
824 
825  rResult( 0 ) = (1.0/2.0)*(2*z - 2.0)*(2*z - 1)*(-2.0*x - 2.0*y + 1.0)*(-x - y + 1.0) ;
826  rResult( 1 ) = (1.0/2.0)*x*(2.0*x - 1.0)*(2*z - 2.0)*(2*z - 1) ;
827  rResult( 2 ) = (1.0/2.0)*y*(2.0*y - 1.0)*(2*z - 2.0)*(2*z - 1) ;
828  rResult( 3 ) = z*(2*z - 1)*(-2.0*x - 2.0*y + 1.0)*(-x - y + 1.0) ;
829  rResult( 4 ) = x*z*(2.0*x - 1.0)*(2*z - 1) ;
830  rResult( 5 ) = y*z*(2.0*y - 1.0)*(2*z - 1) ;
831  rResult( 6 ) = (1.0/2.0)*x*(2*z - 2.0)*(2*z - 1)*(-4.0*x - 4.0*y + 4.0) ;
832  rResult( 7 ) = 2.0*x*y*(2*z - 2.0)*(2*z - 1) ;
833  rResult( 8 ) = 2.0*y*(2*z - 2.0)*(2*z - 1)*(-x - y + 1.0) ;
834  rResult( 9 ) = (1.0 - std::pow(2*z - 1, 2))*(-x - y + 1.0) ;
835  rResult( 10 ) = x*(1.0 - std::pow(2*z - 1, 2)) ;
836  rResult( 11 ) = y*(1.0 - std::pow(2*z - 1, 2)) ;
837  rResult( 12 ) = x*z*(2*z - 1)*(-4.0*x - 4.0*y + 4.0) ;
838  rResult( 13 ) = 4.0*x*y*z*(2*z - 1) ;
839  rResult( 14 ) = 4.0*y*z*(2*z - 1)*(-x - y + 1.0) ;
840 
841  return rResult;
842  }
843 
852  static Matrix& CalculateShapeFunctionsLocalGradients( Matrix& DN, const CoordinatesArrayType& rPoint )
853  {
854  const double x = rPoint[0];
855  const double y = rPoint[1];
856  const double z = rPoint[2];
857  DN.resize(15,3,false);
858 
859  DN( 0 , 0 )= (1.0/2.0)*(2*z - 2.0)*(2*z - 1)*(4.0*x + 4.0*y - 3.0) ;
860  DN( 0 , 1 )= (1.0/2.0)*(2*z - 2.0)*(2*z - 1)*(4.0*x + 4.0*y - 3.0) ;
861  DN( 0 , 2 )= (4*z - 3.0)*(x + y - 1.0)*(2.0*x + 2.0*y - 1.0) ;
862  DN( 1 , 0 )= (1.0/2.0)*(4.0*x - 1.0)*(2*z - 2.0)*(2*z - 1) ;
863  DN( 1 , 1 )= 0 ;
864  DN( 1 , 2 )= x*(2.0*x - 1.0)*(4*z - 3.0) ;
865  DN( 2 , 0 )= 0 ;
866  DN( 2 , 1 )= (1.0/2.0)*(4.0*y - 1.0)*(2*z - 2.0)*(2*z - 1) ;
867  DN( 2 , 2 )= y*(2.0*y - 1.0)*(4*z - 3.0) ;
868  DN( 3 , 0 )= z*(2*z - 1)*(4.0*x + 4.0*y - 3.0) ;
869  DN( 3 , 1 )= z*(2*z - 1)*(4.0*x + 4.0*y - 3.0) ;
870  DN( 3 , 2 )= (4*z - 1)*(x + y - 1.0)*(2.0*x + 2.0*y - 1.0) ;
871  DN( 4 , 0 )= z*(4.0*x - 1.0)*(2*z - 1) ;
872  DN( 4 , 1 )= 0 ;
873  DN( 4 , 2 )= x*(2.0*x - 1.0)*(4*z - 1) ;
874  DN( 5 , 0 )= 0 ;
875  DN( 5 , 1 )= z*(4.0*y - 1.0)*(2*z - 1) ;
876  DN( 5 , 2 )= y*(2.0*y - 1.0)*(4*z - 1) ;
877  DN( 6 , 0 )= 2.0*(2*z - 2.0)*(2*z - 1)*(-2*x - y + 1) ;
878  DN( 6 , 1 )= x*(-8.0*std::pow(z, 2) + 12.0*z - 4.0) ;
879  DN( 6 , 2 )= 4.0*x*(3.0 - 4*z)*(x + y - 1) ;
880  DN( 7 , 0 )= y*(8.0*std::pow(z, 2) - 12.0*z + 4.0) ;
881  DN( 7 , 1 )= x*(8.0*std::pow(z, 2) - 12.0*z + 4.0) ;
882  DN( 7 , 2 )= x*y*(16.0*z - 12.0) ;
883  DN( 8 , 0 )= y*(-8.0*std::pow(z, 2) + 12.0*z - 4.0) ;
884  DN( 8 , 1 )= -(2*z - 2.0)*(2.0*y*(2*z - 1) + (4.0*z - 2.0)*(x + y - 1.0)) ;
885  DN( 8 , 2 )= 4.0*y*(3.0 - 4*z)*(x + y - 1.0) ;
886  DN( 9 , 0 )= 4*z*(z - 1) ;
887  DN( 9 , 1 )= 4*z*(z - 1) ;
888  DN( 9 , 2 )= 4*(2*z - 1)*(x + y - 1.0) ;
889  DN( 10 , 0 )= 4*z*(1 - z) ;
890  DN( 10 , 1 )= 0 ;
891  DN( 10 , 2 )= 4*x*(1 - 2*z) ;
892  DN( 11 , 0 )= 0 ;
893  DN( 11 , 1 )= 4*z*(1 - z) ;
894  DN( 11 , 2 )= 4*y*(1 - 2*z) ;
895  DN( 12 , 0 )= 4.0*z*(2*z - 1)*(-2*x - y + 1) ;
896  DN( 12 , 1 )= x*z*(4.0 - 8.0*z) ;
897  DN( 12 , 2 )= x*(4.0 - 16.0*z)*(x + y - 1) ;
898  DN( 13 , 0 )= y*z*(8.0*z - 4.0) ;
899  DN( 13 , 1 )= x*z*(8.0*z - 4.0) ;
900  DN( 13 , 2 )= x*y*(16.0*z - 4.0) ;
901  DN( 14 , 0 )= y*z*(4.0 - 8.0*z) ;
902  DN( 14 , 1 )= 4.0*z*(2*z - 1)*(-x - 2*y + 1.0) ;
903  DN( 14 , 2 )= y*(4.0 - 16.0*z)*(x + y - 1.0) ;
904 
905  return DN;
906  }
907 
916  static Matrix CalculateShapeFunctionsIntegrationPointsValues(typename BaseType::IntegrationMethod ThisMethod )
917  {
918  IntegrationPointsContainerType all_integration_points = AllIntegrationPoints();
919  IntegrationPointsArrayType integration_points = all_integration_points[static_cast<int>(ThisMethod)];
920 
921  // Number of integration points
922  const std::size_t integration_points_number = integration_points.size();
923 
924  //Setting up return matrix
925  Matrix shape_function_values( integration_points_number, 15 );
926 
927  // Loop over all integration points
928  double x, y, z;
929  for ( std::size_t pnt = 0; pnt < integration_points_number; pnt++ ) {
930  const auto& r_point = integration_points[pnt];
931  x = r_point[0];
932  y = r_point[1];
933  z = r_point[2];
934 
935  shape_function_values( pnt, 0 ) = (1.0/2.0)*(2*z - 2.0)*(2*z - 1)*(-2.0*x - 2.0*y + 1.0)*(-x - y + 1.0) ;
936  shape_function_values( pnt, 1 ) = (1.0/2.0)*x*(2.0*x - 1.0)*(2*z - 2.0)*(2*z - 1) ;
937  shape_function_values( pnt, 2 ) = (1.0/2.0)*y*(2.0*y - 1.0)*(2*z - 2.0)*(2*z - 1) ;
938  shape_function_values( pnt, 3 ) = z*(2*z - 1)*(-2.0*x - 2.0*y + 1.0)*(-x - y + 1.0) ;
939  shape_function_values( pnt, 4 ) = x*z*(2.0*x - 1.0)*(2*z - 1) ;
940  shape_function_values( pnt, 5 ) = y*z*(2.0*y - 1.0)*(2*z - 1) ;
941  shape_function_values( pnt, 6 ) = (1.0/2.0)*x*(2*z - 2.0)*(2*z - 1)*(-4.0*x - 4.0*y + 4.0) ;
942  shape_function_values( pnt, 7 ) = 2.0*x*y*(2*z - 2.0)*(2*z - 1) ;
943  shape_function_values( pnt, 8 ) = 2.0*y*(2*z - 2.0)*(2*z - 1)*(-x - y + 1.0) ;
944  shape_function_values( pnt, 9 ) = (1.0 - std::pow(2*z - 1, 2))*(-x - y + 1.0) ;
945  shape_function_values( pnt, 10 ) = x*(1.0 - std::pow(2*z - 1, 2)) ;
946  shape_function_values( pnt, 11 ) = y*(1.0 - std::pow(2*z - 1, 2)) ;
947  shape_function_values( pnt, 12 ) = x*z*(2*z - 1)*(-4.0*x - 4.0*y + 4.0) ;
948  shape_function_values( pnt, 13 ) = 4.0*x*y*z*(2*z - 1) ;
949  shape_function_values( pnt, 14 ) = 4.0*y*z*(2*z - 1)*(-x - y + 1.0) ;
950  }
951 
952  return shape_function_values;
953  }
954 
964  static ShapeFunctionsGradientsType CalculateShapeFunctionsIntegrationPointsLocalGradients(typename BaseType::IntegrationMethod ThisMethod )
965  {
966  IntegrationPointsContainerType all_integration_points = AllIntegrationPoints();
967  IntegrationPointsArrayType integration_points = all_integration_points[static_cast<int>(ThisMethod)];
968 
969  // Number of integration points
970  const std::size_t integration_points_number = integration_points.size();
971  ShapeFunctionsGradientsType d_shape_f_values( integration_points_number );
972 
973  // Initialising container
974  Matrix result = ZeroMatrix( 15, 3 );
975 
976  // Loop over all integration points
977  for ( std::size_t pnt = 0; pnt < integration_points_number; ++pnt ) {
978  CalculateShapeFunctionsLocalGradients(result, integration_points[pnt]);
979  d_shape_f_values[pnt] = result;
980  }
981 
982  return d_shape_f_values;
983  }
984 
985  static const IntegrationPointsContainerType AllIntegrationPoints()
986  {
987  IntegrationPointsContainerType integration_points =
988  {
989  {
990  Quadrature < PrismGaussLegendreIntegrationPoints1,
991  3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
992  Quadrature < PrismGaussLegendreIntegrationPoints2,
993  3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
994  Quadrature < PrismGaussLegendreIntegrationPoints3,
995  3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
996  Quadrature < PrismGaussLegendreIntegrationPoints4,
997  3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
998  Quadrature < PrismGaussLegendreIntegrationPoints5,
999  3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
1000  Quadrature < PrismGaussLegendreIntegrationPointsExt1,
1001  3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
1002  Quadrature < PrismGaussLegendreIntegrationPointsExt2,
1003  3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
1004  Quadrature < PrismGaussLegendreIntegrationPointsExt3,
1005  3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
1006  Quadrature < PrismGaussLegendreIntegrationPointsExt4,
1007  3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
1008  Quadrature < PrismGaussLegendreIntegrationPointsExt5,
1009  3, IntegrationPoint<3> >::GenerateIntegrationPoints(),
1010  }
1011  };
1012  return integration_points;
1013  }
1014 
1015  static const ShapeFunctionsValuesContainerType AllShapeFunctionsValues()
1016  {
1017  ShapeFunctionsValuesContainerType shape_functions_values =
1018  {
1019  {
1020  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1022  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1024  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1026  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1028  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1030  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1032  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1034  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1036  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1038  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(
1040  }
1041  };
1042  return shape_functions_values;
1043  }
1044 
1045  static const ShapeFunctionsLocalGradientsContainerType AllShapeFunctionsLocalGradients()
1046  {
1047  ShapeFunctionsLocalGradientsContainerType shape_functions_local_gradients =
1048  {
1049  {
1050  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1052  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1054  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1056  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1058  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1060  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1062  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1064  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1066  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1068  Prism3D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(
1070  }
1071  };
1072  return shape_functions_local_gradients;
1073  }
1074 
1079  template<class TOtherPointType> friend class Prism3D15;
1080 
1081 
1086 };// Class Prism3D15
1087 
1088 
1094 template<class TPointType> inline std::istream& operator >> (
1095  std::istream& rIStream, Prism3D15<TPointType>& rThis );
1096 
1098 template<class TPointType> inline std::ostream& operator << (
1099  std::ostream& rOStream, const Prism3D15<TPointType>& rThis )
1100 {
1101  rThis.PrintInfo( rOStream );
1102  rOStream << std::endl;
1103  rThis.PrintData( rOStream );
1104 
1105  return rOStream;
1106 }
1107 
1108 
1109 template<class TPointType> const
1110 GeometryData Prism3D15<TPointType>::msGeometryData(
1113  Prism3D15<TPointType>::AllIntegrationPoints(),
1114  Prism3D15<TPointType>::AllShapeFunctionsValues(),
1115  AllShapeFunctionsLocalGradients()
1116 );
1117 
1118 template<class TPointType> const
1119 GeometryDimension Prism3D15<TPointType>::msGeometryDimension(3, 3);
1120 
1121 }// namespace Kratos.
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
const Matrix & ShapeFunctionsValues() const
Definition: geometry.h:3393
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
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
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
void reserve(size_type dim)
Definition: pointer_vector.h:319
void push_back(const TPointerType &x)
Definition: pointer_vector.h:270
A fifteen node prism geometry with quadratic shape functions.
Definition: prism_3d_15.h:59
bool IsInside(const CoordinatesArrayType &rPoint, CoordinatesArrayType &rResult, const double Tolerance=std::numeric_limits< double >::epsilon()) const override
Definition: prism_3d_15.h:489
GeometryData::KratosGeometryFamily GetGeometryFamily() const override
Definition: prism_3d_15.h:281
SizeType FacesNumber() const override
Returns the number of faces of the current geometry.
Definition: prism_3d_15.h:588
SizeType EdgesNumber() const override
This method gives you number of all edges of this geometry.
Definition: prism_3d_15.h:521
Matrix & PointsLocalCoordinates(Matrix &rResult) const override
Definition: prism_3d_15.h:413
TPointType PointType
Definition: prism_3d_15.h:96
GeometryData::KratosGeometryType GetGeometryType() const override
Definition: prism_3d_15.h:286
GeometryData::IntegrationMethod IntegrationMethod
Definition: prism_3d_15.h:85
GeometriesArrayType GenerateFaces() const override
Returns all faces of the current geometry.
Definition: prism_3d_15.h:601
Matrix MatrixType
Definition: prism_3d_15.h:182
double Area() const override
This method calculate and return area or surface area of this geometry depending to it's dimension.
Definition: prism_3d_15.h:396
BaseType::ShapeFunctionsValuesContainerType ShapeFunctionsValuesContainerType
Definition: prism_3d_15.h:146
KRATOS_CLASS_POINTER_DEFINITION(Prism3D15)
Quadrilateral3D8< TPointType > FaceType2
Definition: prism_3d_15.h:75
BaseType::Pointer Create(const IndexType NewGeometryId, const BaseType &rGeometry) const override
Creates a new geometry pointer.
Definition: prism_3d_15.h:356
Prism3D15(const PointsArrayType &rThisPoints)
Definition: prism_3d_15.h:223
BaseType::IntegrationPointsArrayType IntegrationPointsArrayType
Definition: prism_3d_15.h:132
Matrix & ShapeFunctionsLocalGradients(Matrix &rResult, const CoordinatesArrayType &rPoint) const override
Definition: prism_3d_15.h:682
~Prism3D15() override
Destructor. Does nothing!!!
Definition: prism_3d_15.h:279
Prism3D15(Prism3D15 const &rOther)
Definition: prism_3d_15.h:256
BaseType::NormalType NormalType
Definition: prism_3d_15.h:172
BaseType::PointsArrayType PointsArrayType
Definition: prism_3d_15.h:117
GeometriesArrayType GenerateEdges() const override
This method gives you all edges of this geometry.
Definition: prism_3d_15.h:534
BaseType::JacobiansType JacobiansType
Definition: prism_3d_15.h:160
BaseType::ShapeFunctionsLocalGradientsContainerType ShapeFunctionsLocalGradientsContainerType
Definition: prism_3d_15.h:153
double ShapeFunctionValue(IndexType ShapeFunctionIndex, const CoordinatesArrayType &rPoint) const override
Definition: prism_3d_15.h:662
BaseType::IntegrationPointsContainerType IntegrationPointsContainerType
Definition: prism_3d_15.h:139
BaseType::Pointer Create(const IndexType NewGeometryId, PointsArrayType const &rThisPoints) const override
Creates a new geometry pointer.
Definition: prism_3d_15.h:342
BaseType::IndexType IndexType
Definition: prism_3d_15.h:103
Prism3D15(const IndexType GeometryId, const PointsArrayType &rThisPoints)
Constructor with Geometry Id.
Definition: prism_3d_15.h:230
BaseType::IntegrationPointType IntegrationPointType
Definition: prism_3d_15.h:124
Prism3D15(const std::string &rGeometryName, const PointsArrayType &rThisPoints)
Constructor with Geometry Name.
Definition: prism_3d_15.h:239
BaseType::ShapeFunctionsGradientsType ShapeFunctionsGradientsType
Definition: prism_3d_15.h:167
Prism3D15 & operator=(const Prism3D15 &rOther)
Definition: prism_3d_15.h:306
Prism3D15(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)
Definition: prism_3d_15.h:188
void PrintData(std::ostream &rOStream) const override
Definition: prism_3d_15.h:724
double Volume() const override
This method calculate and return volume of this geometry.
Definition: prism_3d_15.h:391
BaseType::GeometriesArrayType GeometriesArrayType
Definition: prism_3d_15.h:91
Geometry< TPointType > BaseType
Definition: prism_3d_15.h:68
Prism3D15 & operator=(Prism3D15< TOtherPointType > const &rOther)
Definition: prism_3d_15.h:324
Line3D3< TPointType > EdgeType
Definition: prism_3d_15.h:73
BaseType::SizeType SizeType
Definition: prism_3d_15.h:111
BaseType::CoordinatesArrayType CoordinatesArrayType
Definition: prism_3d_15.h:177
Prism3D15(Prism3D15< TOtherPointType > const &rOther)
Definition: prism_3d_15.h:273
Triangle3D6< TPointType > FaceType1
Definition: prism_3d_15.h:74
friend class Prism3D15
Definition: prism_3d_15.h:1079
double Length() const override
Definition: prism_3d_15.h:401
double DomainSize() const override
Definition: prism_3d_15.h:378
void PrintInfo(std::ostream &rOStream) const override
Definition: prism_3d_15.h:710
std::string Info() const override
Definition: prism_3d_15.h:698
A eight node 3D quadrilateral geometry with quadratic shape functions.
Definition: quadrilateral_3d_8.h:49
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
A six node 3D triangular geometry with quadratic shape functions.
Definition: triangle_3d_6.h:68
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
z
Definition: GenerateWind.py:163
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
Internals::Matrix< double, AMatrix::dynamic, 1 > Vector
Definition: amatrix_interface.h:472
KratosZeroMatrix< double > ZeroMatrix
Definition: amatrix_interface.h:559
Internals::Matrix< double, AMatrix::dynamic, AMatrix::dynamic > Matrix
Definition: amatrix_interface.h:470
const GeometryData Prism3D15< TPointType >::msGeometryData & msGeometryDimension(), Prism3D15< 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
y
Other simbols definition.
Definition: generate_axisymmetric_navier_stokes_element.py:54
DN
Definition: generate_convection_diffusion_explicit_element.py:98
def load(f)
Definition: ode_solve.py:307
x
Definition: sensitivityMatrix.py:49
volume
Definition: sp_statistics.py:15