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.
triangle_2d_15.h
Go to the documentation of this file.
1 #pragma once
2 // | / |
3 // ' / __| _` | __| _ \ __|
4 // . \ | ( | | ( |\__ `
5 // _|\_\_| \__,_|\__|\___/ ____/
6 // Multi-Physics
7 //
8 // License: BSD License
9 // Kratos default license: kratos/license.txt
10 //
11 // Main authors: Mohamed Nabi
12 //
13 //
14 // Contributors:
15 //
16 //
17 
18 #pragma once
19 
20 // System includes
21 
22 // External includes
23 
24 // Project includes
25 #include "geometries/line_2d_5.h"
27 
28 namespace Kratos
29 {
32 
36 
40 
44 
48 
68  template<class TPointType>
69  class Triangle2D15 : public Geometry<TPointType>
70  {
71  public:
75 
80 
85 
90 
95 
101 
105  typedef TPointType PointType;
106 
113  typedef typename BaseType::IndexType IndexType;
114 
120  typedef typename BaseType::SizeType SizeType;
121 
127 
132 
138 
147 
154 
160 
166 
173 
180 
188 
196 
201 
205 
206  Triangle2D15(const PointType& Point01, const PointType& Point02, const PointType& Point03,
207  const PointType& Point04, const PointType& Point05, const PointType& Point06,
208  const PointType& Point07, const PointType& Point08, const PointType& Point09,
209  const PointType& Point10, const PointType& Point11, const PointType& Point12,
210  const PointType& Point13, const PointType& Point14, const PointType& Point15)
211  : BaseType(PointsArrayType(), &msGeometryData)
212  {
213  auto& r_points = this->Points();
214  r_points.push_back(typename PointType::Pointer(new PointType(Point01)));
215  r_points.push_back(typename PointType::Pointer(new PointType(Point02)));
216  r_points.push_back(typename PointType::Pointer(new PointType(Point03)));
217  r_points.push_back(typename PointType::Pointer(new PointType(Point04)));
218  r_points.push_back(typename PointType::Pointer(new PointType(Point05)));
219  r_points.push_back(typename PointType::Pointer(new PointType(Point06)));
220  r_points.push_back(typename PointType::Pointer(new PointType(Point07)));
221  r_points.push_back(typename PointType::Pointer(new PointType(Point08)));
222  r_points.push_back(typename PointType::Pointer(new PointType(Point09)));
223  r_points.push_back(typename PointType::Pointer(new PointType(Point10)));
224  r_points.push_back(typename PointType::Pointer(new PointType(Point11)));
225  r_points.push_back(typename PointType::Pointer(new PointType(Point12)));
226  r_points.push_back(typename PointType::Pointer(new PointType(Point13)));
227  r_points.push_back(typename PointType::Pointer(new PointType(Point14)));
228  r_points.push_back(typename PointType::Pointer(new PointType(Point15)));
229  }
230 
231  Triangle2D15(typename PointType::Pointer pPoint01, typename PointType::Pointer pPoint02,
232  typename PointType::Pointer pPoint03, typename PointType::Pointer pPoint04,
233  typename PointType::Pointer pPoint05, typename PointType::Pointer pPoint06,
234  typename PointType::Pointer pPoint07, typename PointType::Pointer pPoint08,
235  typename PointType::Pointer pPoint09, typename PointType::Pointer pPoint10,
236  typename PointType::Pointer pPoint11, typename PointType::Pointer pPoint12,
237  typename PointType::Pointer pPoint13, typename PointType::Pointer pPoint14,
238  typename PointType::Pointer pPoint15) : BaseType(PointsArrayType(), &msGeometryData)
239  {
240  auto& r_points = this->Points();
241  r_points.push_back(pPoint01);
242  r_points.push_back(pPoint02);
243  r_points.push_back(pPoint03);
244  r_points.push_back(pPoint04);
245  r_points.push_back(pPoint05);
246  r_points.push_back(pPoint06);
247  r_points.push_back(pPoint07);
248  r_points.push_back(pPoint08);
249  r_points.push_back(pPoint09);
250  r_points.push_back(pPoint10);
251  r_points.push_back(pPoint11);
252  r_points.push_back(pPoint12);
253  r_points.push_back(pPoint13);
254  r_points.push_back(pPoint14);
255  r_points.push_back(pPoint15);
256  }
257 
258  Triangle2D15(const PointsArrayType& ThisPoints) : BaseType(ThisPoints, &msGeometryData)
259  {
260  if (this->PointsNumber() != 15) {
261  KRATOS_ERROR << "Invalid points number. Expected 15, given " << this->PointsNumber()
262  << std::endl;
263  }
264  }
265 
267  explicit Triangle2D15(const IndexType GeometryId, const PointsArrayType& rThisPoints)
268  : BaseType(GeometryId, rThisPoints, &msGeometryData)
269  {
270  KRATOS_ERROR_IF(this->PointsNumber() != 15) << "Invalid points number. Expected 15, given "
271  << this->PointsNumber() << std::endl;
272  }
273 
275  explicit Triangle2D15(const std::string& rGeometryName, const PointsArrayType& rThisPoints)
276  : BaseType(rGeometryName, rThisPoints, &msGeometryData)
277  {
278  KRATOS_ERROR_IF(this->PointsNumber() != 15) << "Invalid points number. Expected 15, given "
279  << this->PointsNumber() << std::endl;
280  }
281 
291  Triangle2D15(Triangle2D15 const& rOther) : BaseType(rOther)
292  {
293  }
294 
307  template<class TOtherPointType> Triangle2D15(Triangle2D15<TOtherPointType> const& rOther)
308  : BaseType(rOther)
309  {
310  }
311 
315  ~Triangle2D15() override {}
316 
318  {
320  }
321 
323  {
325  }
326 
330 
343  {
344  BaseType::operator=(rOther);
345  return *this;
346  }
347 
359  template<class TOtherPointType>
361  {
362  BaseType::operator=(rOther);
363  return *this;
364  }
365 
369 
375  typename BaseType::Pointer Create(PointsArrayType const& rThisPoints) const override
376  {
377  return typename BaseType::Pointer(new Triangle2D15(rThisPoints));
378  }
379 
386  typename BaseType::Pointer Create(const IndexType NewGeometryId,
387  PointsArrayType const& rThisPoints) const override
388  {
389  return typename BaseType::Pointer(new Triangle2D15(NewGeometryId, rThisPoints));
390  }
391 
397  typename BaseType::Pointer Create(const BaseType& rGeometry) const override
398  {
399  auto p_geometry = typename BaseType::Pointer(new Triangle2D15(rGeometry.Points()));
400  p_geometry->SetData(rGeometry.GetData());
401  return p_geometry;
402  }
403 
410  typename BaseType::Pointer Create(const IndexType NewGeometryId, const BaseType& rGeometry) const override
411  {
412  auto p_geometry = typename BaseType::Pointer(new Triangle2D15(NewGeometryId, rGeometry.Points()));
413  p_geometry->SetData(rGeometry.GetData());
414  return p_geometry;
415  }
416 
422  Matrix& PointsLocalCoordinates(Matrix& rResult) const override
423  {
424  rResult.resize(15, 2, false);
425  //
426  rResult(0, 0) = 0.00;
427  rResult(0, 1) = 0.00;
428  //
429  rResult(1, 0) = 1.00;
430  rResult(1, 1) = 0.00;
431  //
432  rResult(2, 0) = 0.00;
433  rResult(2, 1) = 1.00;
434  //
435  rResult(3, 0) = 0.25;
436  rResult(3, 1) = 0.00;
437  //
438  rResult(4, 0) = 0.50;
439  rResult(4, 1) = 0.00;
440  //
441  rResult(5, 0) = 0.75;
442  rResult(5, 1) = 0.00;
443  //
444  rResult(6, 0) = 0.75;
445  rResult(6, 1) = 0.25;
446  //
447  rResult(7, 0) = 0.50;
448  rResult(7, 1) = 0.50;
449  //
450  rResult(8, 0) = 0.25;
451  rResult(8, 1) = 0.75;
452  //
453  rResult(9, 0) = 0.00;
454  rResult(9, 1) = 0.75;
455  //
456  rResult(10, 0) = 0.00;
457  rResult(10, 1) = 0.50;
458  //
459  rResult(11, 0) = 0.00;
460  rResult(11, 1) = 0.25;
461  //
462  rResult(12, 0) = 0.25;
463  rResult(12, 1) = 0.25;
464  //
465  rResult(13, 0) = 0.50;
466  rResult(13, 1) = 0.25;
467  //
468  rResult(14, 0) = 0.25;
469  rResult(14, 1) = 0.50;
470  //
471  return rResult;
472  }
473 
477 
493  double Length() const override
494  {
495  return std::sqrt(std::abs(Area()));
496  }
497 
509  double Area() const override
510  {
511  Vector temp;
512  this->DeterminantOfJacobian(temp, msGeometryData.DefaultIntegrationMethod());
513  const IntegrationPointsArrayType& integration_points = this->IntegrationPoints(msGeometryData.DefaultIntegrationMethod());
514  double area = 0.00;
515  for (unsigned int i = 0; i < integration_points.size(); ++i)
516  {
517  area += temp[i] * integration_points[i].Weight();
518  }
519  return area;
520  }
521 
532  double DomainSize() const override
533  {
534  return Area();
535  }
536 
545  bool IsInside(const CoordinatesArrayType& rPoint, CoordinatesArrayType& rResult,
546  const double Tolerance = std::numeric_limits<double>::epsilon()) const override
547  {
548  this->PointLocalCoordinates(rResult, rPoint);
549  if ((rResult[0] >= (0.0 - Tolerance)) && (rResult[0] <= (1.0 + Tolerance))) {
550  if ((rResult[1] >= (0.0 - Tolerance)) && (rResult[1] <= (1.0 + Tolerance))) {
551  if ((rResult[0] + rResult[1]) <= (1.0 + Tolerance)) {
552  return true;
553  }
554  }
555  }
556  return false;
557  }
558 
562 
569  Vector& ShapeFunctionsValues(Vector& rResult, const CoordinatesArrayType& rCoordinates) const override
570  {
571  if (rResult.size() != 15) rResult.resize(15, false);
572  //
573  const double cof1 = 128.0 / 3.0;
574  const double cof2 = 32.0 / 3.0;
575  const double xi = rCoordinates[0];
576  const double et = rCoordinates[1];
577  const double zt = 1.0 - xi - et;
578  //
579  rResult[0] = zt * (zt - 0.25) * (zt - 0.5) * (zt - 0.75) * cof2;
580  rResult[1] = xi * (xi - 0.25) * (xi - 0.5) * (xi - 0.75) * cof2;
581  rResult[2] = et * (et - 0.25) * (et - 0.5) * (et - 0.75) * cof2;
582  rResult[3] = xi * zt * (zt - 0.25) * (zt - 0.50) * cof1;
583  rResult[4] = xi * zt * (zt - 0.25) * (xi - 0.25) * 64.0;
584  rResult[5] = xi * zt * (xi - 0.25) * (xi - 0.50) * cof1;
585  rResult[6] = xi * et * (xi - 0.25) * (xi - 0.50) * cof1;
586  rResult[7] = xi * et * (xi - 0.25) * (et - 0.25) * 64.0;
587  rResult[8] = xi * et * (et - 0.25) * (et - 0.50) * cof1;
588  rResult[9] = et * zt * (et - 0.25) * (et - 0.50) * cof1;
589  rResult[10] = et * zt * (et - 0.25) * (zt - 0.25) * 64.0;
590  rResult[11] = et * zt * (zt - 0.25) * (zt - 0.50) * cof1;
591  rResult[12] = xi * et * zt * (zt - 0.25) * 128.0;
592  rResult[13] = xi * et * zt * (xi - 0.25) * 128.0;
593  rResult[14] = xi * et * zt * (et - 0.25) * 128.0;
594  //
595  return rResult;
596  }
597 
607  double ShapeFunctionValue(IndexType ShapeFunctionIndex, const CoordinatesArrayType& rPoint) const override
608  {
609  const double cof1 = 128.0 / 3.0;
610  const double cof2 = 32.0 / 3.0;
611  const double xi = rPoint[0];
612  const double et = rPoint[1];
613  const double zt = 1.0 - xi - et;
614  double shape = 0.0;
615  //
616  switch (ShapeFunctionIndex) {
617  case 0:
618  shape = zt * (zt - 0.25) * (zt - 0.5) * (zt - 0.75) * cof2;
619  break;
620  case 1:
621  shape = xi * (xi - 0.25) * (xi - 0.5) * (xi - 0.75) * cof2;
622  break;
623  case 2:
624  shape = et * (et - 0.25) * (et - 0.5) * (et - 0.75) * cof2;
625  break;
626  case 3:
627  shape = xi * zt * (zt - 0.25) * (zt - 0.5) * cof1;
628  break;
629  case 4:
630  shape = xi * zt * (zt - 0.25) * (xi - 0.25) * 64.0;
631  break;
632  case 5:
633  shape = xi * zt * (xi - 0.25) * (xi - 0.5) * cof1;
634  break;
635  case 6:
636  shape = xi * et * (xi - 0.25) * (xi - 0.5) * cof1;
637  break;
638  case 7:
639  shape = xi * et * (xi - 0.25) * (et - 0.25) * 64.0;
640  break;
641  case 8:
642  shape = xi * et * (et - 0.25) * (et - 0.5) * cof1;
643  break;
644  case 9:
645  shape = et * zt * (et - 0.25) * (et - 0.5) * cof1;
646  break;
647  case 10:
648  shape = et * zt * (et - 0.25) * (zt - 0.25) * 64.0;
649  break;
650  case 11:
651  shape = et * zt * (zt - 0.25) * (zt - 0.5) * cof1;
652  break;
653  case 12:
654  shape = xi * et * zt * (zt - 0.25) * 128.0;
655  break;
656  case 13:
657  shape = xi * et * zt * (xi - 0.25) * 128.0;
658  break;
659  case 14:
660  shape = xi * et * zt * (et - 0.25) * 128.0;
661  break;
662  default:
663  KRATOS_ERROR << "Wrong index of shape function!" << *this << std::endl;
664  break;
665  }
666  //
667  return shape;
668  }
669 
673 
681  std::string Info() const override
682  {
683  return "2 dimensional triangle with fifteen nodes in 2D space";
684  }
685 
692  void PrintInfo(std::ostream& rOStream) const override
693  {
694  rOStream << "2 dimensional triangle with fifteen nodes in 2D space";
695  }
696 
707  void PrintData( std::ostream& rOStream ) const override
708  {
709  // Base Geometry class PrintData call
710  BaseType::PrintData( rOStream );
711  std::cout << std::endl;
712 
713  // If the geometry has valid points, calculate and output its data
714  if (this->AllPointsAreValid()) {
715  Matrix jacobian;
716  this->Jacobian( jacobian, PointType() );
717  rOStream << " Jacobian in the origin\t : " << jacobian;
718  }
719  }
720 
724 
736  SizeType EdgesNumber() const override
737  {
738  return 3;
739  }
740 
750  {
752  edges.push_back(Kratos::make_shared<EdgeType>(this->pGetPoint(0), this->pGetPoint(1), this->pGetPoint(3), this->pGetPoint(4), this->pGetPoint(5)));
753  edges.push_back(Kratos::make_shared<EdgeType>(this->pGetPoint(1), this->pGetPoint(2), this->pGetPoint(6), this->pGetPoint(7), this->pGetPoint(8)));
754  edges.push_back(Kratos::make_shared<EdgeType>(this->pGetPoint(2), this->pGetPoint(0), this->pGetPoint(9), this->pGetPoint(10), this->pGetPoint(11)));
755  return edges;
756  }
757 
758  SizeType FacesNumber() const override
759  {
760  return 3;
761  }
762 
763  //Connectivities of faces required
765  {
766  if (NumberNodesInFaces.size() != 3) NumberNodesInFaces.resize(3, false);
767  NumberNodesInFaces[0] = 5;
768  NumberNodesInFaces[1] = 5;
769  NumberNodesInFaces[2] = 5;
770  }
771 
773  {
774  // faces in columns
775  if (NodesInFaces.size1() != 6 || NodesInFaces.size2() != 3)
776  NodesInFaces.resize(6, 3, false);
777 
778  //compatible to Edges() function. Triangle library considers a different order
779  NodesInFaces(0, 0) = 0;//face or master node
780  NodesInFaces(1, 0) = 1;
781  NodesInFaces(2, 0) = 6;
782  NodesInFaces(3, 0) = 7;
783  NodesInFaces(4, 0) = 8;
784  NodesInFaces(5, 0) = 2;
785  //
786  NodesInFaces(0, 1) = 1;//face or master node
787  NodesInFaces(1, 1) = 2;
788  NodesInFaces(2, 1) = 9;
789  NodesInFaces(3, 1) = 10;
790  NodesInFaces(4, 1) = 11;
791  NodesInFaces(5, 1) = 0;
792  //
793  NodesInFaces(0, 2) = 2;//face or master node
794  NodesInFaces(1, 2) = 0;
795  NodesInFaces(2, 2) = 3;
796  NodesInFaces(3, 2) = 4;
797  NodesInFaces(4, 2) = 5;
798  NodesInFaces(5, 2) = 1;
799  }
800 
806  {
807  ShapeFunctionsGradientsType localGradients
808  = CalculateShapeFunctionsIntegrationPointsLocalGradients(ThisMethod);
809  const int integration_points_number = msGeometryData.IntegrationPointsNumber(ThisMethod);
810  ShapeFunctionsGradientsType Result(integration_points_number);
811  //
812  for (int pnt = 0; pnt < integration_points_number; ++pnt)
813  {
814  Result[pnt] = localGradients[pnt];
815  }
816  return Result;
817  }
818 
824  {
825  IntegrationMethod ThisMethod = msGeometryData.DefaultIntegrationMethod();
826  ShapeFunctionsGradientsType localGradients
827  = CalculateShapeFunctionsIntegrationPointsLocalGradients(ThisMethod);
828  const int integration_points_number = msGeometryData.IntegrationPointsNumber(ThisMethod);
829  ShapeFunctionsGradientsType Result(integration_points_number);
830  //
831  for (int pnt = 0; pnt < integration_points_number; ++pnt)
832  {
833  Result[pnt] = localGradients[pnt];
834  }
835  return Result;
836  }
837 
847  Matrix& ShapeFunctionsLocalGradients(Matrix& rResult, const CoordinatesArrayType& rPoint) const override
848  {
849  rResult.resize(15, 2, false);
850  const double xi1 = rPoint[0];
851  const double xi2 = xi1 * xi1;
852  const double xi3 = xi2 * xi1;
853  const double et1 = rPoint[1];
854  const double et2 = et1 * et1;
855  const double et3 = et2 * et1;
856  const double zt1 = 1.0 - xi1 - et1;
857  const double zt2 = zt1 * zt1;
858  const double zt3 = zt2 * zt1;
859  //
860  noalias(rResult) = ZeroMatrix(15, 2);
861  //
862  rResult(0, 0) = -(128.0 * zt3 - 144.0 * zt2 + 44.0 * zt1 - 3.0) / 3.0;
863  rResult(0, 1) = -(128.0 * zt3 - 144.0 * zt2 + 44.0 * zt1 - 3.0) / 3.0;
864  //
865  rResult(1, 0) = (128.0 * xi3 - 144.0 * xi2 + 44.0 * xi1 - 3.0) / 3.0;
866  rResult(1, 1) = 0.0;
867  //
868  rResult(2, 0) = 0.0;
869  rResult(2, 1) = (128.0 * et3 - 144.0 * et2 + 44.0 * et1 - 3.0) / 3.0;
870  //
871  rResult(3, 0) = -128.0 * (zt2 - 0.5 * zt1 + 1.0 / 24.0) * xi1 + (128.0 * zt3 - 96.0 * zt2 + 16.0 * zt1) / 3.0;
872  rResult(3, 1) = -16.0 * xi1 * (24.0 * zt2 - 12.0 * zt1 + 1.0) / 3.0;
873  //
874  rResult(4, 0) = -128.0 * (xi1 - 0.25) * (zt1 - 0.125) * xi1 + 128.0 * (xi1 - 0.125) * (zt1 - 0.25) * zt1;
875  rResult(4, 1) = -4.0 * xi1 * (4.0 * xi1 - 1.0) * (8.0 * zt1 - 1.0);
876  //
877  rResult(5, 0) = -(128.0 * xi3 - 96.0 * xi2 + 16.0 * xi1) / 3.0 + 128.0 * (xi2 - 0.5 * xi1 + 1.0 / 24.0) * zt1;
878  rResult(5, 1) = -16.0 * xi1 * (8.0 * xi2 - 6.0 * xi1 + 1.0) / 3.0;
879  //
880  rResult(6, 0) = 16.0 * et1 * (24.0 * xi2 - 12.0 * xi1 + 1.0) / 3.0;
881  rResult(6, 1) = (128.0 * xi3 - 96.0 * xi2 + 16.0 * xi1) / 3.0;
882  //
883  rResult(7, 0) = 4.0 * (8.0 * xi1 - 1.0) * (4.0 * et1 - 1.0) * et1;
884  rResult(7, 1) = 4.0 * (4.0 * xi1 - 1.0) * (8.0 * et1 - 1.0) * xi1;
885  //
886  rResult(8, 0) = (128.0 * et3 - 96.0 * et2 + 16.0 * et1) / 3.0;
887  rResult(8, 1) = 16.0 * xi1 * (24.0 * et2 - 12.0 * et1 + 1.0) / 3.0;
888  //
889  rResult(9, 0) = -16.0 * et1 * (8.0 * et2 - 6.0 * et1 + 1.0) / 3.0;
890  rResult(9, 1) = -(128.0 * et3 - 96.0 * et2 + 16.0 * et1) / 3.0 + 128.0 * (et2 - 0.5 * et1 + 1.0 / 24.0) * zt1;
891  //
892  rResult(10, 0) = -4.0 * et1 * (4.0 * et1 - 1.0) * (8.0 * zt1 - 1.0);
893  rResult(10, 1) = -128.0 * (et1 - 0.25) * (zt1 - 0.125) * et1 + 128.0 * (et1 - 0.125) * zt1 * (zt1 - 0.25);
894  //
895  rResult(11, 0) = -16.0 * et1 * (24.0 * zt2 - 12.0 * zt1 + 1.0) / 3.0;
896  rResult(11, 1) = -128.0 * (zt2 - 0.5 * zt1 + 1.0 / 24.0) * et1 + (128.0 * zt3 - 96.0 * zt2 + 16.0 * zt1) / 3.0;
897  //
898  rResult(12, 0) = 256.0 * et1 * (-xi1 * (zt1 - 0.125) + 0.5 * zt2 - 0.125 * zt1);
899  rResult(12, 1) = 256.0 * xi1 * (-et1 * (zt1 - 0.125) + 0.5 * zt2 - 0.125 * zt1);
900  //
901  rResult(13, 0) = -32.0 * et1 * (4.0 * xi2 - xi1) + 256.0 * (xi1 - 0.125) * et1 * zt1;
902  rResult(13, 1) = 128.0 * (xi1 - 0.25) * (-et1 + zt1) * xi1;
903  //
904  rResult(14, 0) = 128.0 * (et1 - 0.25) * et1 * (-xi1 + zt1);
905  rResult(14, 1) = -32.0 * xi1 * (4.0 * et2 - et1) + 256.0 * (et1 - 0.125) * zt1 * xi1;
906  //
907  return rResult;
908  }
909 
919  virtual Matrix& ShapeFunctionsGradients(Matrix& rResult, const CoordinatesArrayType& rPoint)
920  {
921  rResult.resize(15, 2, false);
922  const double xi1 = rPoint[0];
923  const double xi2 = xi1 * xi1;
924  const double xi3 = xi2 * xi1;
925  const double et1 = rPoint[1];
926  const double et2 = et1 * et1;
927  const double et3 = et2 * et1;
928  const double zt1 = 1.0 - xi1 - et1;
929  const double zt2 = zt1 * zt1;
930  const double zt3 = zt2 * zt1;
931  //
932  noalias(rResult) = ZeroMatrix(15, 2);
933  //
934  rResult(0, 0) = -(128.0 * zt3 - 144.0 * zt2 + 44.0 * zt1 - 3.0) / 3.0;
935  rResult(0, 1) = -(128.0 * zt3 - 144.0 * zt2 + 44.0 * zt1 - 3.0) / 3.0;
936  //
937  rResult(1, 0) = (128.0 * xi3 - 144.0 * xi2 + 44.0 * xi1 - 3.0) / 3.0;
938  rResult(1, 1) = 0.0;
939  //
940  rResult(2, 0) = 0.0;
941  rResult(2, 1) = (128.0 * et3 - 144.0 * et2 + 44.0 * et1 - 3.0) / 3.0;
942  //
943  rResult(3, 0) = -128.0 * (zt2 - 0.5 * zt1 + 1.0 / 24.0) * xi1 + (128.0 * zt3 - 96.0 * zt2 + 16.0 * zt1) / 3.0;
944  rResult(3, 1) = -16.0 * xi1 * (24.0 * zt2 - 12.0 * zt1 + 1.0) / 3.0;
945  //
946  rResult(4, 0) = -128.0 * (xi1 - 0.25) * (zt1 - 0.125) * xi1 + 128.0 * (xi1 - 0.125) * (zt1 - 0.25) * zt1;
947  rResult(4, 1) = -4.0 * xi1 * (4.0 * xi1 - 1.0) * (8.0 * zt1 - 1.0);
948  //
949  rResult(5, 0) = -(128.0 * xi3 - 96.0 * xi2 + 16.0 * xi1) / 3.0 + 128.0 * (xi2 - 0.5 * xi1 + 1.0 / 24.0) * zt1;
950  rResult(5, 1) = -16.0 * xi1 * (8.0 * xi2 - 6.0 * xi1 + 1.0) / 3.0;
951  //
952  rResult(6, 0) = 16.0 * et1 * (24.0 * xi2 - 12.0 * xi1 + 1.0) / 3.0;
953  rResult(6, 1) = (128.0 * xi3 - 96.0 * xi2 + 16.0 * xi1) / 3.0;
954  //
955  rResult(7, 0) = 4.0 * (8.0 * xi1 - 1.0) * (4.0 * et1 - 1.0) * et1;
956  rResult(7, 1) = 4.0 * (4.0 * xi1 - 1.0) * (8.0 * et1 - 1.0) * xi1;
957  //
958  rResult(8, 0) = (128.0 * et3 - 96.0 * et2 + 16.0 * et1) / 3.0;
959  rResult(8, 1) = 16.0 * xi1 * (24.0 * et2 - 12.0 * et1 + 1.0) / 3.0;
960  //
961  rResult(9, 0) = -16.0 * et1 * (8.0 * et2 - 6.0 * et1 + 1.0) / 3.0;
962  rResult(9, 1) = -(128.0 * et3 - 96.0 * et2 + 16.0 * et1) / 3.0 + 128.0 * (et2 - 0.5 * et1 + 1.0 / 24.0) * zt1;
963  //
964  rResult(10, 0) = -4.0 * et1 * (4.0 * et1 - 1.0) * (8.0 * zt1 - 1.0);
965  rResult(10, 1) = -128.0 * (et1 - 0.25) * (zt1 - 0.125) * et1 + 128.0 * (et1 - 0.125) * zt1 * (zt1 - 0.25);
966  //
967  rResult(11, 0) = -16.0 * et1 * (24.0 * zt2 - 12.0 * zt1 + 1.0) / 3.0;
968  rResult(11, 1) = -128.0 * (zt2 - 0.5 * zt1 + 1.0 / 24.0) * et1 + (128.0 * zt3 - 96.0 * zt2 + 16.0 * zt1) / 3.0;
969  //
970  rResult(12, 0) = 256.0 * et1 * (-xi1 * (zt1 - 0.125) + 0.5 * zt2 - 0.125 * zt1);
971  rResult(12, 1) = 256.0 * xi1 * (-et1 * (zt1 - 0.125) + 0.5 * zt2 - 0.125 * zt1);
972  //
973  rResult(13, 0) = -32.0 * et1 * (4.0 * xi2 - xi1) + 256.0 * (xi1 - 0.125) * et1 * zt1;
974  rResult(13, 1) = 128.0 * (xi1 - 0.25) * (-et1 + zt1) * xi1;
975  //
976  rResult(14, 0) = 128.0 * (et1 - 0.25) * et1 * (-xi1 + zt1);
977  rResult(14, 1) = -32.0 * xi1 * (4.0 * et2 - et1) + 256.0 * (et1 - 0.125) * zt1 * xi1;
978  //
979  return rResult;
980  }
981 
989  ShapeFunctionsSecondDerivativesType& rResult, const CoordinatesArrayType& rPoint) const override
990  {
991  if (rResult.size() != this->PointsNumber())
992  {
994  rResult.swap(temp);
995  }
996  const double xi1 = rPoint[0];
997  const double xi2 = xi1 * xi1;
998  const double et1 = rPoint[1];
999  const double et2 = et1 * et1;
1000  //
1001  rResult[0].resize(2, 2, false);
1002  rResult[1].resize(2, 2, false);
1003  rResult[2].resize(2, 2, false);
1004  rResult[3].resize(2, 2, false);
1005  rResult[4].resize(2, 2, false);
1006  rResult[5].resize(2, 2, false);
1007  rResult[6].resize(2, 2, false);
1008  rResult[7].resize(2, 2, false);
1009  rResult[8].resize(2, 2, false);
1010  rResult[9].resize(2, 2, false);
1011  rResult[10].resize(2, 2, false);
1012  rResult[11].resize(2, 2, false);
1013  rResult[12].resize(2, 2, false);
1014  rResult[13].resize(2, 2, false);
1015  rResult[14].resize(2, 2, false);
1016  //
1017  rResult[0](0, 0) = 128.0 * (xi2 + et2) - 160.0 * (xi1 + et1) + 256.0 * xi1 * et1 + 140.0 / 3.0;
1018  rResult[0](0, 1) = 128.0 * (xi2 + et2) - 160.0 * (xi1 + et1) + 256.0 * xi1 * et1 + 140.0 / 3.0;
1019  rResult[0](1, 0) = 128.0 * (xi2 + et2) - 160.0 * (xi1 + et1) + 256.0 * xi1 * et1 + 140.0 / 3.0;
1020  rResult[0](1, 1) = 128.0 * (xi2 + et2) - 160.0 * (xi1 + et1) + 256.0 * xi1 * et1 + 140.0 / 3.0;
1021  //
1022  rResult[1](0, 0) = 128.0 * xi2 - 96.0 * xi1 + 44.0 / 3.0;
1023  rResult[1](0, 1) = 0.0;
1024  rResult[1](1, 0) = 0.0;
1025  rResult[1](1, 1) = 0.0;
1026  //
1027  rResult[2](0, 0) = 0.0;
1028  rResult[2](0, 1) = 0.0;
1029  rResult[2](1, 0) = 0.0;
1030  rResult[2](1, 1) = 128.0 * et2 - 96.0 * et1 + 44.0 / 3.0;
1031  //
1032  rResult[3](0, 0) = -512.0 * xi2 - 768.0 * xi1 * et1 - 256.0 * et2 + 576.0 * xi1 + 384.0 * et1 - 416.0 / 3.0;
1033  rResult[3](0, 1) = -384.0 * xi2 - 512.0 * xi1 * et1 - 128.0 * et2 + 384.0 * xi1 + 192.0 * et1 - 208.0 / 3.0;
1034  rResult[3](1, 0) = -384.0 * xi2 - 512.0 * xi1 * et1 - 128.0 * et2 + 384.0 * xi1 + 192.0 * et1 - 208.0 / 3.0;
1035  rResult[3](1, 1) = 64.0 * xi1 * (3.0 - 4.0 * xi1 - 4.0 * et1);
1036  //
1037  rResult[4](0, 0) = 768.0 * xi2 + 768.0 * (et1 - 1.0) * xi1 + 128.0 * et2 - 288.0 * et1 + 152.0;
1038  rResult[4](0, 1) = 384.0 * xi2 + (256.0 * et1 - 288.0) * xi1 - 32.0 * et1 + 28.0;
1039  rResult[4](1, 0) = 384.0 * xi2 + (256.0 * et1 - 288.0) * xi1 - 32.0 * et1 + 28.0;
1040  rResult[4](1, 1) = 32.0 * xi1 * (4.0 * xi1 - 1.0);
1041  //
1042  rResult[5](0, 0) = -512.0 * xi2 + 448.0 * xi1 - 256.0 * xi1 * et1 + 64.0 * et1 - 224.0 / 3.0;
1043  rResult[5](0, 1) = -128.0 * xi2 + 64.0 * xi1 - 16.0 / 3.0;
1044  rResult[5](1, 0) = -128.0 * xi2 + 64.0 * xi1 - 16.0 / 3.0;
1045  rResult[5](1, 1) = 0.0;
1046  //
1047  rResult[6](0, 0) = 64.0 * et1 * (4.0 * xi1 - 1.0);
1048  rResult[6](0, 1) = 128.0 * xi2 - 64.0 * xi1 + 16.0 / 3.0;
1049  rResult[6](1, 0) = 128.0 * xi2 - 64.0 * xi1 + 16.0 / 3.0;
1050  rResult[6](1, 1) = 0.0;
1051  //
1052  rResult[7](0, 0) = 32.0 * et1 * (4.0 * et1 - 1.0);
1053  rResult[7](0, 1) = 32.0 * xi1 * (8.0 * et1 - 1.0) - 32.0 * et1 + 4.0;
1054  rResult[7](1, 0) = 32.0 * xi1 * (8.0 * et1 - 1.0) - 32.0 * et1 + 4.0;
1055  rResult[7](1, 1) = 32.0 * xi1 * (4.0 * xi1 - 1.0);
1056  //
1057  rResult[8](0, 0) = 0.0;
1058  rResult[8](0, 1) = 128.0 * et2 - 64.0 * et1 + 16.0 / 3.0;
1059  rResult[8](1, 0) = 128.0 * et2 - 64.0 * et1 + 16.0 / 3.0;
1060  rResult[8](1, 1) = (256.0 * et1 - 64.0) * xi1;
1061  //
1062  rResult[9](0, 0) = 0.0;
1063  rResult[9](0, 1) = -128.0 * et2 + 64.0 * et1 - 16.0 / 3.0;
1064  rResult[9](1, 0) = -128.0 * et2 + 64.0 * et1 - 16.0 / 3.0;;
1065  rResult[9](1, 1) = 448.0 * et1 - 224.0 / 3.0 - 256.0 * xi1 * et1 + 64.0 * xi1 - 512.0 * et2;
1066  //
1067  rResult[10](0, 0) = 128.0 * et2 - 32.0 * et1;
1068  rResult[10](0, 1) = 384.0 * et2 + (256.0 * xi1 - 288.0) * et1 - 32.0 * xi1 + 28.0;
1069  rResult[10](1, 0) = 384.0 * et2 + (256.0 * xi1 - 288.0) * et1 - 32.0 * xi1 + 28.0;
1070  rResult[10](1, 1) = 128.0 * xi2 + (768.0 * et1 - 288.0) * xi1 + 768.0 * et2 - 768.0 * et1 + 152.0;
1071  //
1072  rResult[11](0, 0) = -64.0 * et1 * (-3.0 + 4.0 * xi1 + 4.0 * et1);
1073  rResult[11](0, 1) = -384.0 * et2 + (-512.0 * xi1 + 384.0) * et1 - 128.0 * xi2 + 192.0 * xi1 - 208.0 / 3.0;
1074  rResult[11](1, 0) = -384.0 * et2 + (-512.0 * xi1 + 384.0) * et1 - 128.0 * xi2 + 192.0 * xi1 - 208.0 / 3.0;
1075  rResult[11](1, 1) = -512.0 * et2 + (-768.0 * xi1 + 576.0) * et1 - 256.0 * xi2 + 384.0 * xi1 - 416.0 / 3.0;
1076  //
1077  rResult[12](0, 0) = 64.0 * et1 * (12.0 * xi1 + 8.0 * et1 - 7.0);
1078  rResult[12](0, 1) = 384.0 * xi2 + (1024.0 * et1 - 448.0) * xi1 + 384.0 * et2 - 448.0 * et1 + 96.0;
1079  rResult[12](1, 0) = 384.0 * xi2 + (1024.0 * et1 - 448.0) * xi1 + 384.0 * et2 - 448.0 * et1 + 96.0;
1080  rResult[12](1, 1) = 64.0 * xi1 * (8.0 * xi1 + 12.0 * et1 - 7.0);
1081  //
1082  rResult[13](0, 0) = -64.0 * et1 * (12.0 * xi1 + 4.0 * et1 - 5.0);
1083  rResult[13](0, 1) = -384.0 * xi2 + (-512.0 * et1 + 320.0) * xi1 + 64.0 * et1 - 32.0;
1084  rResult[13](1, 0) = -384.0 * xi2 + (-512.0 * et1 + 320.0) * xi1 + 64.0 * et1 - 32.0;
1085  rResult[13](1, 1) = -256.0 * xi2 + 64.0 * xi1;
1086  //
1087  rResult[14](0, 0) = -256.0 * et2 + 64.0 * et1;
1088  rResult[14](0, 1) = -384.0 * et2 + (-512.0 * xi1 + 320.0) * et1 + 64.0 * xi1 - 32.0;
1089  rResult[14](1, 0) = -384.0 * et2 + (-512.0 * xi1 + 320.0) * et1 + 64.0 * xi1 - 32.0;
1090  rResult[14](1, 1) = -64.0 * xi1 * (4.0 * xi1 + 12.0 * et1 - 5.0);
1091  //
1092  return rResult;
1093  }
1094 
1102  ShapeFunctionsThirdDerivativesType& rResult, const CoordinatesArrayType& rPoint) const override
1103  {
1104  if (rResult.size() != this->PointsNumber())
1105  {
1107  rResult.swap(temp);
1108  }
1109  for (IndexType i = 0; i < rResult.size(); ++i)
1110  {
1112  rResult[i].swap(temp);
1113  }
1114  //
1115  rResult[0][0].resize(2, 2, false);
1116  rResult[0][1].resize(2, 2, false);
1117  rResult[1][0].resize(2, 2, false);
1118  rResult[1][1].resize(2, 2, false);
1119  rResult[2][0].resize(2, 2, false);
1120  rResult[2][1].resize(2, 2, false);
1121  rResult[3][0].resize(2, 2, false);
1122  rResult[3][1].resize(2, 2, false);
1123  rResult[4][0].resize(2, 2, false);
1124  rResult[4][1].resize(2, 2, false);
1125  rResult[5][0].resize(2, 2, false);
1126  rResult[5][1].resize(2, 2, false);
1127  rResult[6][0].resize(2, 2, false);
1128  rResult[6][1].resize(2, 2, false);
1129  rResult[7][0].resize(2, 2, false);
1130  rResult[7][1].resize(2, 2, false);
1131  rResult[8][0].resize(2, 2, false);
1132  rResult[8][1].resize(2, 2, false);
1133  rResult[9][0].resize(2, 2, false);
1134  rResult[9][1].resize(2, 2, false);
1135  rResult[10][0].resize(2, 2, false);
1136  rResult[10][1].resize(2, 2, false);
1137  rResult[11][0].resize(2, 2, false);
1138  rResult[11][1].resize(2, 2, false);
1139  rResult[12][0].resize(2, 2, false);
1140  rResult[12][1].resize(2, 2, false);
1141  rResult[13][0].resize(2, 2, false);
1142  rResult[13][1].resize(2, 2, false);
1143  rResult[14][0].resize(2, 2, false);
1144  rResult[14][1].resize(2, 2, false);
1145  //
1146  double fx3 = 256.0 * rPoint[0] + 256.0 * rPoint[1] - 160.0;
1147  double fx2y = fx3;
1148  double fxy2 = fx3;
1149  double fy3 = fx3;
1150  rResult[0][0](0, 0) = fx3; // dx3
1151  rResult[0][0](0, 1) = fx2y; // dx2 dy
1152  rResult[0][0](1, 0) = fx2y; // dx2 dy
1153  rResult[0][0](1, 1) = fxy2; // dx dy2
1154  rResult[0][1](0, 0) = fx2y; // dx2 dy
1155  rResult[0][1](0, 1) = fxy2; // dx dy2
1156  rResult[0][1](1, 0) = fxy2; // dx dy2
1157  rResult[0][1](1, 1) = fy3; // dy3
1158  //
1159  fx3 = 256.0 * rPoint[0] - 96.0;
1160  rResult[1][0](0, 0) = fx3; // dx3
1161  rResult[1][0](0, 1) = 0.0; // dx2 dy
1162  rResult[1][0](1, 0) = 0.0; // dx2 dy
1163  rResult[1][0](1, 1) = 0.0; // dx dy2
1164  rResult[1][1](0, 0) = 0.0; // dx2 dy
1165  rResult[1][1](0, 1) = 0.0; // dx dy2
1166  rResult[1][1](1, 0) = 0.0; // dx dy2
1167  rResult[1][1](1, 1) = 0.0; // dy3
1168  //
1169  fy3 = 256.0 * rPoint[1] - 96.0;
1170  rResult[2][0](0, 0) = 0.0; // dx3
1171  rResult[2][0](0, 1) = 0.0; // dx2 dy
1172  rResult[2][0](1, 0) = 0.0; // dx2 dy
1173  rResult[2][0](1, 1) = 0.0; // dx dy2
1174  rResult[2][1](0, 0) = 0.0; // dx2 dy
1175  rResult[2][1](0, 1) = 0.0; // dx dy2
1176  rResult[2][1](1, 0) = 0.0; // dx dy2
1177  rResult[2][1](1, 1) = fy3; // dy3
1178  //
1179  fx3 = -768.0 * rPoint[1] - 1024.0 * rPoint[0] + 576.0;
1180  fx2y = -512.0 * rPoint[1] - 768.0 * rPoint[0] + 384.0;
1181  fxy2 = -256.0 * rPoint[1] - 512.0 * rPoint[0] + 192.0;
1182  fy3 = -256.0 * rPoint[0];
1183  rResult[3][0](0, 0) = fx3; // dx3
1184  rResult[3][0](0, 1) = fx2y; // dx2 dy
1185  rResult[3][0](1, 0) = fx2y; // dx2 dy
1186  rResult[3][0](1, 1) = fxy2; // dx dy2
1187  rResult[3][1](0, 0) = fx2y; // dx2 dy
1188  rResult[3][1](0, 1) = fxy2; // dx dy2
1189  rResult[3][1](1, 0) = fxy2; // dx dy2
1190  rResult[3][1](1, 1) = fy3; // dy3
1191  //
1192  fx3 = 1536.0 * rPoint[0] + 768.0 * rPoint[1] - 768.0;
1193  fx2y = 768.0 * rPoint[0] + 256.0 * rPoint[1] - 288.0;
1194  fxy2 = 256.0 * rPoint[0] - 32.0;
1195  rResult[4][0](0, 0) = fx3; // dx3
1196  rResult[4][0](0, 1) = fx2y; // dx2 dy
1197  rResult[4][0](1, 0) = fx2y; // dx2 dy
1198  rResult[4][0](1, 1) = fxy2; // dx dy2
1199  rResult[4][1](0, 0) = fx2y; // dx2 dy
1200  rResult[4][1](0, 1) = fxy2; // dx dy2
1201  rResult[4][1](1, 0) = fxy2; // dx dy2
1202  rResult[4][1](1, 1) = 0.0; // dy3
1203  //
1204  fx3 = 448.0 - 1024.0 * rPoint[0] - 256.0 * rPoint[1];
1205  fx2y = -256.0 * rPoint[0] + 64.0;
1206  rResult[5][0](0, 0) = fx3; // dx3
1207  rResult[5][0](0, 1) = fx2y; // dx2 dy
1208  rResult[5][0](1, 0) = fx2y; // dx2 dy
1209  rResult[5][0](1, 1) = 0.0; // dx dy2
1210  rResult[5][1](0, 0) = fx2y; // dx2 dy
1211  rResult[5][1](0, 1) = 0.0; // dx dy2
1212  rResult[5][1](1, 0) = 0.0; // dx dy2
1213  rResult[5][1](1, 1) = 0.0; // dy3
1214  //
1215  fx3 = 256.0 * rPoint[1];
1216  fx2y = 256.0 * rPoint[0] - 64.0;
1217  rResult[6][0](0, 0) = fx3; // dx3
1218  rResult[6][0](0, 1) = fx2y; // dx2 dy
1219  rResult[6][0](1, 0) = fx2y; // dx2 dy
1220  rResult[6][0](1, 1) = 0.0; // dx dy2
1221  rResult[6][1](0, 0) = fx2y; // dx2 dy
1222  rResult[6][1](0, 1) = 0.0; // dx dy2
1223  rResult[6][1](1, 0) = 0.0; // dx dy2
1224  rResult[6][1](1, 1) = 0.0; // dy3
1225  //
1226  fx2y = 256.0 * rPoint[1] - 32.0;
1227  fxy2 = 256.0 * rPoint[0] - 32.0;
1228  rResult[7][0](0, 0) = 0.0; // dx3
1229  rResult[7][0](0, 1) = fx2y; // dx2 dy
1230  rResult[7][0](1, 0) = fx2y; // dx2 dy
1231  rResult[7][0](1, 1) = fxy2; // dx dy2
1232  rResult[7][1](0, 0) = fx2y; // dx2 dy
1233  rResult[7][1](0, 1) = fxy2; // dx dy2
1234  rResult[7][1](1, 0) = fxy2; // dx dy2
1235  rResult[7][1](1, 1) = 0.0; // dy3
1236  //
1237  fxy2 = 256.0 * rPoint[1] - 64.0;
1238  fy3 = 256.0 * rPoint[0];
1239  rResult[8][0](0, 0) = 0.0; // dx3
1240  rResult[8][0](0, 1) = 0.0; // dx2 dy
1241  rResult[8][0](1, 0) = 0.0; // dx2 dy
1242  rResult[8][0](1, 1) = fxy2; // dx dy2
1243  rResult[8][1](0, 0) = 0.0; // dx2 dy
1244  rResult[8][1](0, 1) = fxy2; // dx dy2
1245  rResult[8][1](1, 0) = fxy2; // dx dy2
1246  rResult[8][1](1, 1) = fy3; // dy3
1247  //
1248  fxy2 = -256.0 * rPoint[1] + 64.0;
1249  fy3 = 448.0 - 256.0 * rPoint[0] - 1024.0 * rPoint[1];
1250  rResult[9][0](0, 0) = 0.0; // dx3
1251  rResult[9][0](0, 1) = 0.0; // dx2 dy
1252  rResult[9][0](1, 0) = 0.0; // dx2 dy
1253  rResult[9][0](1, 1) = fxy2; // dx dy2
1254  rResult[9][1](0, 0) = 0.0; // dx2 dy
1255  rResult[9][1](0, 1) = fxy2; // dx dy2
1256  rResult[9][1](1, 0) = fxy2; // dx dy2
1257  rResult[9][1](1, 1) = fy3; // dy3
1258  //
1259  fx2y = 256.0 * rPoint[1] - 32.0;
1260  fxy2 = 768.0 * rPoint[1] + 256.0 * rPoint[0] - 288.0;
1261  fy3 = 768.0 * rPoint[0] + 1536.0 * rPoint[1] - 768.0;
1262  rResult[10][0](0, 0) = 0.0; // dx3
1263  rResult[10][0](0, 1) = fx2y; // dx2 dy
1264  rResult[10][0](1, 0) = fx2y; // dx2 dy
1265  rResult[10][0](1, 1) = fxy2; // dx dy2
1266  rResult[10][1](0, 0) = fx2y; // dx2 dy
1267  rResult[10][1](0, 1) = fxy2; // dx dy2
1268  rResult[10][1](1, 0) = fxy2; // dx dy2
1269  rResult[10][1](1, 1) = fy3; // dy3
1270  //
1271  fx3 = -256.0 * rPoint[1];
1272  fx2y = -256.0 * rPoint[0] - 512.0 * rPoint[1] + 192.0;
1273  fxy2 = -768.0 * rPoint[1] - 512.0 * rPoint[0] + 384.0;
1274  fy3 = -1024.0 * rPoint[1] - 768.0 * rPoint[0] + 576.0;
1275  rResult[11][0](0, 0) = fx3; // dx3
1276  rResult[11][0](0, 1) = fx2y; // dx2 dy
1277  rResult[11][0](1, 0) = fx2y; // dx2 dy
1278  rResult[11][0](1, 1) = fxy2; // dx dy2
1279  rResult[11][1](0, 0) = fx2y; // dx2 dy
1280  rResult[11][1](0, 1) = fxy2; // dx dy2
1281  rResult[11][1](1, 0) = fxy2; // dx dy2
1282  rResult[11][1](1, 1) = fy3; // dy3
1283  //
1284  fx3 = 768.0 * rPoint[1];
1285  fx2y = 768.0 * rPoint[0] + 1024.0 * rPoint[1] - 448.0;
1286  fxy2 = 1024.0 * rPoint[0] + 768.0 * rPoint[1] - 448.0;
1287  fy3 = 768.0 * rPoint[0];
1288  rResult[12][0](0, 0) = fx3; // dx3
1289  rResult[12][0](0, 1) = fx2y; // dx2 dy
1290  rResult[12][0](1, 0) = fx2y; // dx2 dy
1291  rResult[12][0](1, 1) = fxy2; // dx dy2
1292  rResult[12][1](0, 0) = fx2y; // dx2 dy
1293  rResult[12][1](0, 1) = fxy2; // dx dy2
1294  rResult[12][1](1, 0) = fxy2; // dx dy2
1295  rResult[12][1](1, 1) = fy3; // dy3
1296  //
1297  fx3 = -768.0 * rPoint[1];
1298  fx2y = -768.0 * rPoint[0] - 512.0 * rPoint[1] + 320.0;
1299  fxy2 = -512.0 * rPoint[0] + 64.0;
1300  rResult[13][0](0, 0) = fx3; // dx3
1301  rResult[13][0](0, 1) = fx2y; // dx2 dy
1302  rResult[13][0](1, 0) = fx2y; // dx2 dy
1303  rResult[13][0](1, 1) = fxy2; // dx dy2
1304  rResult[13][1](0, 0) = fx2y; // dx2 dy
1305  rResult[13][1](0, 1) = fxy2; // dx dy2
1306  rResult[13][1](1, 0) = fxy2; // dx dy2
1307  rResult[13][1](1, 1) = 0.0; // dy3
1308  //
1309  fx2y = -512.0 * rPoint[1] + 64.0;
1310  fxy2 = -768.0 * rPoint[1] - 512.0 * rPoint[0] + 320.0;
1311  fy3 = -768.0 * rPoint[0];
1312  rResult[14][0](0, 0) = 0.0; // dx3
1313  rResult[14][0](0, 1) = fx2y; // dx2 dy
1314  rResult[14][0](1, 0) = fx2y; // dx2 dy
1315  rResult[14][0](1, 1) = fxy2; // dx dy2
1316  rResult[14][1](0, 0) = fx2y; // dx2 dy
1317  rResult[14][1](0, 1) = fxy2; // dx dy2
1318  rResult[14][1](1, 0) = fxy2; // dx dy2
1319  rResult[14][1](1, 1) = fy3; // dy3
1320  //
1321  return rResult;
1322  }
1323 
1327 
1329 
1330  protected:
1335  private:
1338 
1339  static const GeometryData msGeometryData;
1340 
1341  static const GeometryDimension msGeometryDimension;
1342 
1346 
1347  friend class Serializer;
1348 
1349  void save(Serializer& rSerializer) const override
1350  {
1352  }
1353 
1354  void load(Serializer& rSerializer) override
1355  {
1357  }
1358 
1359  Triangle2D15() : BaseType(PointsArrayType(), &msGeometryData) {}
1360 
1364 
1368 
1372 
1380  static Matrix CalculateShapeFunctionsIntegrationPointsValues(typename BaseType::IntegrationMethod ThisMethod)
1381  {
1382  IntegrationPointsContainerType all_integration_points = AllIntegrationPoints();
1383  IntegrationPointsArrayType integration_points = all_integration_points[static_cast<int>(ThisMethod)];
1384  //number of integration points
1385  const int integration_points_number = integration_points.size();
1386  //number of nodes in current geometry
1387  const int points_number = 15;
1388  //setting up return matrix
1389  Matrix shape_function_values(integration_points_number, points_number);
1390  //
1391  const double cof1 = 128.0 / 3.0;
1392  const double cof2 = 32.0 / 3.0;
1393  //loop over all integration points
1394  for (int pnt = 0; pnt < integration_points_number; ++pnt)
1395  {
1396  const double xi = integration_points[pnt].X();
1397  const double et = integration_points[pnt].Y();
1398  const double zt = 1.0 - xi - et;
1399  //
1400  shape_function_values(pnt, 0) = zt * (zt - 0.25) * (zt - 0.50) * (zt - 0.75) * cof2;
1401  shape_function_values(pnt, 1) = xi * (xi - 0.25) * (xi - 0.50) * (xi - 0.75) * cof2;
1402  shape_function_values(pnt, 2) = et * (et - 0.25) * (et - 0.50) * (et - 0.75) * cof2;
1403  shape_function_values(pnt, 3) = xi * zt * (zt - 0.25) * (zt - 0.50) * cof1;
1404  shape_function_values(pnt, 4) = xi * zt * (zt - 0.25) * (xi - 0.25) * 64.0;
1405  shape_function_values(pnt, 5) = xi * zt * (xi - 0.25) * (xi - 0.50) * cof1;
1406  shape_function_values(pnt, 6) = xi * et * (xi - 0.25) * (xi - 0.50) * cof1;
1407  shape_function_values(pnt, 7) = xi * et * (xi - 0.25) * (et - 0.25) * 64.0;
1408  shape_function_values(pnt, 8) = xi * et * (et - 0.25) * (et - 0.50) * cof1;
1409  shape_function_values(pnt, 9) = et * zt * (et - 0.25) * (et - 0.50) * cof1;
1410  shape_function_values(pnt, 10) = et * zt * (et - 0.25) * (zt - 0.25) * 64.0;
1411  shape_function_values(pnt, 11) = et * zt * (zt - 0.25) * (zt - 0.50) * cof1;
1412  shape_function_values(pnt, 12) = xi * et * zt * (zt - 0.25) * 128.0;
1413  shape_function_values(pnt, 13) = xi * et * zt * (xi - 0.25) * 128.0;
1414  shape_function_values(pnt, 14) = xi * et * zt * (et - 0.25) * 128.0;
1415  }
1416  return shape_function_values;
1417  }
1418 
1429  CalculateShapeFunctionsIntegrationPointsLocalGradients(typename BaseType::IntegrationMethod ThisMethod)
1430  {
1431  IntegrationPointsContainerType all_integration_points = AllIntegrationPoints();
1432  IntegrationPointsArrayType integration_points = all_integration_points[static_cast<int>(ThisMethod)];
1433  //number of integration points
1434  const int integration_points_number = integration_points.size();
1435  ShapeFunctionsGradientsType d_shape_f_values(integration_points_number);
1436  //
1437  //loop over all integration points
1438  for (int pnt = 0; pnt < integration_points_number; ++pnt)
1439  {
1440  Matrix result(15, 2);
1441  const double xi1 = integration_points[pnt].X();
1442  const double xi2 = xi1 * xi1;
1443  const double xi3 = xi2 * xi1;
1444  const double et1 = integration_points[pnt].Y();
1445  const double et2 = et1 * et1;
1446  const double et3 = et2 * et1;
1447  const double zt1 = 1.0 - xi1 - et1;
1448  const double zt2 = zt1 * zt1;
1449  const double zt3 = zt2 * zt1;
1450  //
1451  noalias(result) = ZeroMatrix(15, 2);
1452  //
1453  result(0, 0) = -(128.0 * zt3 - 144.0 * zt2 + 44.0 * zt1 - 3.0) / 3.0;
1454  result(0, 1) = -(128.0 * zt3 - 144.0 * zt2 + 44.0 * zt1 - 3.0) / 3.0;
1455  //
1456  result(1, 0) = (128.0 * xi3 - 144.0 * xi2 + 44.0 * xi1 - 3.0) / 3.0;
1457  result(1, 1) = 0.0;
1458  //
1459  result(2, 0) = 0.0;
1460  result(2, 1) = (128.0 * et3 - 144.0 * et2 + 44.0 * et1 - 3.0) / 3.0;
1461  //
1462  result(3, 0) = -128.0 * (zt2 - 0.5 * zt1 + 1.0 / 24.0) * xi1 + (128.0 * zt3 - 96.0 * zt2 + 16.0 * zt1) / 3.0;
1463  result(3, 1) = -16.0 * xi1 * (24.0 * zt2 - 12.0 * zt1 + 1.0) / 3.0;
1464  //
1465  result(4, 0) = -128.0 * (xi1 - 0.25) * (zt1 - 0.125) * xi1 + 128.0 * (xi1 - 0.125) * (zt1 - 0.25) * zt1;
1466  result(4, 1) = -4.0 * xi1 * (4.0 * xi1 - 1.0) * (8.0 * zt1 - 1.0);
1467  //
1468  result(5, 0) = -(128.0 * xi3 - 96.0 * xi2 + 16.0 * xi1) / 3.0 + 128.0 * (xi2 - 0.5 * xi1 + 1.0 / 24.0) * zt1;
1469  result(5, 1) = -16.0 * xi1 * (8.0 * xi2 - 6.0 * xi1 + 1.0) / 3.0;
1470  //
1471  result(6, 0) = 16.0 * et1 * (24.0 * xi2 - 12.0 * xi1 + 1.0) / 3.0;
1472  result(6, 1) = (128.0 * xi3 - 96.0 * xi2 + 16.0 * xi1) / 3.0;
1473  //
1474  result(7, 0) = 4.0 * (8.0 * xi1 - 1.0) * (4.0 * et1 - 1.0) * et1;
1475  result(7, 1) = 4.0 * (4.0 * xi1 - 1.0) * (8.0 * et1 - 1.0) * xi1;
1476  //
1477  result(8, 0) = (128.0 * et3 - 96.0 * et2 + 16.0 * et1) / 3.0;
1478  result(8, 1) = 16.0 * xi1 * (24.0 * et2 - 12.0 * et1 + 1.0) / 3.0;
1479  //
1480  result(9, 0) = -16.0 * et1 * (8.0 * et2 - 6.0 * et1 + 1.0) / 3.0;
1481  result(9, 1) = -(128.0 * et3 - 96.0 * et2 + 16.0 * et1) / 3.0 + 128.0 * (et2 - 0.5 * et1 + 1.0 / 24.0) * zt1;
1482  //
1483  result(10, 0) = -4.0 * et1 * (4.0 * et1 - 1.0) * (8.0 * zt1 - 1.0);
1484  result(10, 1) = -128.0 * (et1 - 0.25) * (zt1 - 0.125) * et1 + 128.0 * (et1 - 0.125) * zt1 * (zt1 - 0.25);
1485  //
1486  result(11, 0) = -16.0 * et1 * (24.0 * zt2 - 12.0 * zt1 + 1.0) / 3.0;
1487  result(11, 1) = -128.0 * (zt2 - 0.5 * zt1 + 1.0 / 24.0) * et1 + (128.0 * zt3 - 96.0 * zt2 + 16.0 * zt1) / 3.0;
1488  //
1489  result(12, 0) = 256.0 * et1 * (-xi1 * (zt1 - 0.125) + 0.5 * zt2 - 0.125 * zt1);
1490  result(12, 1) = 256.0 * xi1 * (-et1 * (zt1 - 0.125) + 0.5 * zt2 - 0.125 * zt1);
1491  //
1492  result(13, 0) = -32.0 * et1 * (4.0 * xi2 - xi1) + 256.0 * (xi1 - 0.125) * et1 * zt1;
1493  result(13, 1) = 128.0 * (xi1 - 0.25) * (-et1 + zt1) * xi1;
1494  //
1495  result(14, 0) = 128.0 * (et1 - 0.25) * et1 * (-xi1 + zt1);
1496  result(14, 1) = -32.0 * xi1 * (4.0 * et2 - et1) + 256.0 * (et1 - 0.125) * zt1 * xi1;
1497  //
1498  d_shape_f_values[pnt] = result;
1499  }
1500  return d_shape_f_values;
1501  }
1502 
1503  static const IntegrationPointsContainerType AllIntegrationPoints()
1504  {
1505  IntegrationPointsContainerType integration_points =
1506  {
1507  {
1508  Quadrature<TriangleGaussLegendreIntegrationPoints1, 2, IntegrationPoint<3>>::GenerateIntegrationPoints(),
1509  Quadrature<TriangleGaussLegendreIntegrationPoints2, 2, IntegrationPoint<3>>::GenerateIntegrationPoints(),
1510  Quadrature<TriangleGaussLegendreIntegrationPoints3, 2, IntegrationPoint<3>>::GenerateIntegrationPoints(),
1511  Quadrature<TriangleGaussLegendreIntegrationPoints4, 2, IntegrationPoint<3>>::GenerateIntegrationPoints(),
1512  Quadrature<TriangleGaussLegendreIntegrationPoints5, 2, IntegrationPoint<3>>::GenerateIntegrationPoints()
1513  }
1514  };
1515  return integration_points;
1516  }
1517 
1518  static const ShapeFunctionsValuesContainerType AllShapeFunctionsValues()
1519  {
1520  ShapeFunctionsValuesContainerType shape_functions_values =
1521  {
1522  {
1523  Triangle2D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(GeometryData::IntegrationMethod::GI_GAUSS_1),
1524  Triangle2D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(GeometryData::IntegrationMethod::GI_GAUSS_2),
1525  Triangle2D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(GeometryData::IntegrationMethod::GI_GAUSS_3),
1526  Triangle2D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(GeometryData::IntegrationMethod::GI_GAUSS_4),
1527  Triangle2D15<TPointType>::CalculateShapeFunctionsIntegrationPointsValues(GeometryData::IntegrationMethod::GI_GAUSS_5),
1528  }
1529  };
1530  return shape_functions_values;
1531  }
1532 
1533  static const ShapeFunctionsLocalGradientsContainerType AllShapeFunctionsLocalGradients()
1534  {
1535  ShapeFunctionsLocalGradientsContainerType shape_functions_local_gradients =
1536  {
1537  {
1538  Triangle2D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(GeometryData::IntegrationMethod::GI_GAUSS_1),
1539  Triangle2D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(GeometryData::IntegrationMethod::GI_GAUSS_2),
1540  Triangle2D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(GeometryData::IntegrationMethod::GI_GAUSS_3),
1541  Triangle2D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(GeometryData::IntegrationMethod::GI_GAUSS_4),
1542  Triangle2D15<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients(GeometryData::IntegrationMethod::GI_GAUSS_5),
1543  }
1544  };
1545  return shape_functions_local_gradients;
1546  }
1547 
1551 
1555 
1559 
1560  template<class TOtherPointType> friend class Triangle2D15;
1561 
1565 
1567  }; // Class Geometry
1568 
1572 
1576 
1580  template<class TPointType>
1581  inline std::istream& operator >> (std::istream& rIStream, Triangle2D15<TPointType>& rThis);
1582 
1586  template<class TPointType>
1587  inline std::ostream& operator << (std::ostream& rOStream, const Triangle2D15<TPointType>& rThis)
1588  {
1589  rThis.PrintInfo(rOStream);
1590  rOStream << std::endl;
1591  rThis.PrintData(rOStream);
1592  return rOStream;
1593  }
1594 
1596 
1597  template<class TPointType> const
1598  GeometryData Triangle2D15<TPointType>::msGeometryData(
1601  Triangle2D15<TPointType>::AllIntegrationPoints(),
1602  Triangle2D15<TPointType>::AllShapeFunctionsValues(),
1603  AllShapeFunctionsLocalGradients()
1604  );
1605 
1606  template<class TPointType>
1607  const GeometryDimension Triangle2D15<TPointType>::msGeometryDimension(2, 2);
1608 
1609 }// namespace Kratos.
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
IntegrationMethod DefaultIntegrationMethod() const
Definition: geometry_data.h:425
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
const IntegrationPointsArrayType & IntegrationPoints() const
Definition: geometry.h:2284
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
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 five node 2D line geometry with quartic shape functions.
Definition: line_2d_5.h:62
This class defines the node.
Definition: node.h:65
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
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
A fifteen node 2D triangular geometry with quartic shape functions.
Definition: triangle_2d_15.h:70
GeometryData::KratosGeometryFamily GetGeometryFamily() const override
Definition: triangle_2d_15.h:317
BaseType::CoordinatesArrayType CoordinatesArrayType
Definition: triangle_2d_15.h:131
double DomainSize() const override
Definition: triangle_2d_15.h:532
double Area() const override
Definition: triangle_2d_15.h:509
BaseType::SizeType SizeType
Definition: triangle_2d_15.h:120
BaseType::IntegrationPointsArrayType IntegrationPointsArrayType
Definition: triangle_2d_15.h:146
KRATOS_CLASS_POINTER_DEFINITION(Triangle2D15)
BaseType::ShapeFunctionsThirdDerivativesType ShapeFunctionsThirdDerivativesType
Definition: triangle_2d_15.h:195
ShapeFunctionsSecondDerivativesType & ShapeFunctionsSecondDerivatives(ShapeFunctionsSecondDerivativesType &rResult, const CoordinatesArrayType &rPoint) const override
Definition: triangle_2d_15.h:988
Triangle2D15(Triangle2D15< TOtherPointType > const &rOther)
Definition: triangle_2d_15.h:307
Triangle2D15 & operator=(Triangle2D15< TOtherPointType > const &rOther)
Definition: triangle_2d_15.h:360
BaseType::ShapeFunctionsLocalGradientsContainerType ShapeFunctionsLocalGradientsContainerType
Definition: triangle_2d_15.h:165
GeometriesArrayType GenerateEdges() const override
This method gives you all edges of this geometry.
Definition: triangle_2d_15.h:749
TPointType PointType
Definition: triangle_2d_15.h:105
Triangle2D15(Triangle2D15 const &rOther)
Definition: triangle_2d_15.h:291
ShapeFunctionsThirdDerivativesType & ShapeFunctionsThirdDerivatives(ShapeFunctionsThirdDerivativesType &rResult, const CoordinatesArrayType &rPoint) const override
Definition: triangle_2d_15.h:1101
BaseType::IntegrationPointsContainerType IntegrationPointsContainerType
Definition: triangle_2d_15.h:153
Triangle2D15(const std::string &rGeometryName, const PointsArrayType &rThisPoints)
Constructor with Geometry Name.
Definition: triangle_2d_15.h:275
BaseType::IntegrationPointType IntegrationPointType
Definition: triangle_2d_15.h:137
virtual Matrix & ShapeFunctionsGradients(Matrix &rResult, const CoordinatesArrayType &rPoint)
Definition: triangle_2d_15.h:919
GeometryData::IntegrationMethod IntegrationMethod
Definition: triangle_2d_15.h:94
BaseType::Pointer Create(const BaseType &rGeometry) const override
Creates a new geometry pointer.
Definition: triangle_2d_15.h:397
BaseType::Pointer Create(const IndexType NewGeometryId, PointsArrayType const &rThisPoints) const override
Creates a new geometry pointer.
Definition: triangle_2d_15.h:386
friend class Triangle2D15
Definition: triangle_2d_15.h:1560
void NumberNodesInFaces(DenseVector< unsigned int > &NumberNodesInFaces) const override
Definition: triangle_2d_15.h:764
Matrix & PointsLocalCoordinates(Matrix &rResult) const override
Definition: triangle_2d_15.h:422
BaseType::IndexType IndexType
Definition: triangle_2d_15.h:113
Vector & ShapeFunctionsValues(Vector &rResult, const CoordinatesArrayType &rCoordinates) const override
Returns vector of shape function values at local coordinate.
Definition: triangle_2d_15.h:569
BaseType::PointsArrayType PointsArrayType
Definition: triangle_2d_15.h:126
SizeType EdgesNumber() const override
This method gives you number of all edges of this geometry.
Definition: triangle_2d_15.h:736
Triangle2D15(typename PointType::Pointer pPoint01, typename PointType::Pointer pPoint02, typename PointType::Pointer pPoint03, typename PointType::Pointer pPoint04, typename PointType::Pointer pPoint05, typename PointType::Pointer pPoint06, typename PointType::Pointer pPoint07, typename PointType::Pointer pPoint08, typename PointType::Pointer pPoint09, 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: triangle_2d_15.h:231
BaseType::NormalType NormalType
Definition: triangle_2d_15.h:200
std::string Info() const override
Definition: triangle_2d_15.h:681
Matrix & ShapeFunctionsLocalGradients(Matrix &rResult, const CoordinatesArrayType &rPoint) const override
Definition: triangle_2d_15.h:847
virtual ShapeFunctionsGradientsType ShapeFunctionsLocalGradients()
Definition: triangle_2d_15.h:823
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: triangle_2d_15.h:545
BaseType::JacobiansType JacobiansType
Definition: triangle_2d_15.h:172
Triangle2D15 & operator=(const Triangle2D15 &rOther)
Definition: triangle_2d_15.h:342
Triangle2D15(const PointType &Point01, const PointType &Point02, const PointType &Point03, const PointType &Point04, const PointType &Point05, const PointType &Point06, const PointType &Point07, const PointType &Point08, const PointType &Point09, const PointType &Point10, const PointType &Point11, const PointType &Point12, const PointType &Point13, const PointType &Point14, const PointType &Point15)
Definition: triangle_2d_15.h:206
double Length() const override
Definition: triangle_2d_15.h:493
void PrintData(std::ostream &rOStream) const override
Definition: triangle_2d_15.h:707
Geometry< TPointType > BaseType
Definition: triangle_2d_15.h:79
Triangle2D15(const IndexType GeometryId, const PointsArrayType &rThisPoints)
Constructor with Geometry Id.
Definition: triangle_2d_15.h:267
BaseType::Pointer Create(PointsArrayType const &rThisPoints) const override
Creates a new geometry pointer.
Definition: triangle_2d_15.h:375
BaseType::GeometriesArrayType GeometriesArrayType
Definition: triangle_2d_15.h:100
Triangle2D15(const PointsArrayType &ThisPoints)
Definition: triangle_2d_15.h:258
BaseType::Pointer Create(const IndexType NewGeometryId, const BaseType &rGeometry) const override
Creates a new geometry pointer.
Definition: triangle_2d_15.h:410
SizeType FacesNumber() const override
Returns the number of faces of the current geometry.
Definition: triangle_2d_15.h:758
void NodesInFaces(DenseMatrix< unsigned int > &NodesInFaces) const override
Definition: triangle_2d_15.h:772
BaseType::ShapeFunctionsSecondDerivativesType ShapeFunctionsSecondDerivativesType
Definition: triangle_2d_15.h:187
BaseType::ShapeFunctionsGradientsType ShapeFunctionsGradientsType
Definition: triangle_2d_15.h:179
void PrintInfo(std::ostream &rOStream) const override
Definition: triangle_2d_15.h:692
Line2D5< TPointType > EdgeType
Definition: triangle_2d_15.h:84
double ShapeFunctionValue(IndexType ShapeFunctionIndex, const CoordinatesArrayType &rPoint) const override
Definition: triangle_2d_15.h:607
BaseType::ShapeFunctionsValuesContainerType ShapeFunctionsValuesContainerType
Definition: triangle_2d_15.h:159
GeometryData::KratosGeometryType GetGeometryType() const override
Definition: triangle_2d_15.h:322
virtual ShapeFunctionsGradientsType ShapeFunctionsLocalGradients(IntegrationMethod ThisMethod)
Definition: triangle_2d_15.h:805
~Triangle2D15() override
Definition: triangle_2d_15.h:315
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
const GeometryData Triangle2D15< TPointType >::msGeometryData & msGeometryDimension(), Triangle2D15< 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
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
def load(f)
Definition: ode_solve.py:307
float temp
Definition: rotating_cone.py:85
integer i
Definition: TensorModule.f:17