-
Notifications
You must be signed in to change notification settings - Fork 248
/
Copy pathquadrilateral_3d_9.h
1778 lines (1489 loc) · 61.4 KB
/
quadrilateral_3d_9.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
// | / |
// ' / __| _` | __| _ \ __|
// . \ | ( | | ( |\__ `
// _|\_\_| \__,_|\__|\___/ ____/
// Multi-Physics
//
// License: BSD License
// Kratos default license: kratos/license.txt
//
// Main authors: Riccardo Rossi
// Janosch Stascheit
// Felix Nagel
// contributors: Hoang Giang Bui
// Josep Maria Carbonell
//
#pragma once
// System includes
// External includes
// Project includes
#include "geometries/line_3d_3.h"
#include "utilities/integration_utilities.h"
#include "integration/quadrilateral_gauss_legendre_integration_points.h"
namespace Kratos
{
/**
* @class Quadrilateral3D9
* @ingroup KratosCore
* @brief A nine node 3D quadrilateral geometry with quadratic shape functions
* @details While the shape functions are only defined in 2D it is possible to define an arbitrary orientation in space. Thus it can be used for defining surfaces on 3D elements.
* The node ordering corresponds with:
* 3-----6-----2
* | |
* | |
* 7 8 5
* | |
* | |
* 0-----4-----1
* @author Riccardo Rossi
* @author Janosch Stascheit
* @author Felix Nagel
*/
template<class TPointType> class Quadrilateral3D9 : public Geometry<TPointType>
{
public:
/**
* Type Definitions
*/
/**
* Base Type: Geometry
*/
typedef Geometry<TPointType> BaseType;
/**
* Type of edge geometry
*/
typedef Line3D3<TPointType> EdgeType;
/**
* Pointer definition of Quadrilateral3D9
*/
KRATOS_CLASS_POINTER_DEFINITION( Quadrilateral3D9 );
/**
* Integration methods implemented in geometry.
*/
typedef GeometryData::IntegrationMethod IntegrationMethod;
/**
* A Vector of counted pointers to Geometries.
* Used for returning edges of the geometry.
*/
typedef typename BaseType::GeometriesArrayType GeometriesArrayType;
/**
* Redefinition of template parameter TPointType.
*/
typedef TPointType PointType;
/**
* Array of coordinates. Can be Nodes, Points or IntegrationPoints
*/
typedef typename BaseType::CoordinatesArrayType CoordinatesArrayType;
/**
* Type used for indexing in geometry class.
* std::size_t used for indexing
* point or integration point access methods and also all other
* methods which need point or integration point index.
*/
typedef typename BaseType::IndexType IndexType;
/**
* This type is used to return size or dimension in
* geometry. Dimension, WorkingDimension, PointsNumber and
* ... return this type as their results.
*/
typedef typename BaseType::SizeType SizeType;
/**
* Array of counted pointers to point.
* This type is used to hold geometry's points.
*/
typedef typename BaseType::PointsArrayType PointsArrayType;
/**
* This type used for representing an integration point in
* geometry. This integration point is a point with an
* additional weight component.
*/
typedef typename BaseType::IntegrationPointType IntegrationPointType;
/**
* A Vector of IntegrationPointType which used to hold
* integration points related to an integration
* method. IntegrationPoints functions used this type to return
* their results.
*/
typedef typename BaseType::IntegrationPointsArrayType IntegrationPointsArrayType;
/**
* A Vector of IntegrationPointsArrayType which used to hold
* integration points related to different integration method
* implemented in geometry.
*/
typedef typename BaseType::IntegrationPointsContainerType IntegrationPointsContainerType;
/**
* A third order tensor used as shape functions' values
* container.
*/
typedef typename BaseType::ShapeFunctionsValuesContainerType
ShapeFunctionsValuesContainerType;
/**
* A fourth order tensor used as shape functions' local
* gradients container in geometry.
*/
typedef typename BaseType::ShapeFunctionsLocalGradientsContainerType
ShapeFunctionsLocalGradientsContainerType;
/**
* A third order tensor to hold jacobian matrices evaluated at
* integration points. Jacobian and InverseOfJacobian functions
* return this type as their result.
*/
typedef typename BaseType::JacobiansType JacobiansType;
/**
* A third order tensor to hold shape functions' local
* gradients. ShapefunctionsLocalGradients function return this
* type as its result.
*/
typedef typename BaseType::ShapeFunctionsGradientsType ShapeFunctionsGradientsType;
/**
* A third order tensor to hold shape functions' local second derivatives.
* ShapefunctionsLocalGradients function return this
* type as its result.
*/
typedef typename BaseType::ShapeFunctionsSecondDerivativesType
ShapeFunctionsSecondDerivativesType;
/**
* A fourth order tensor to hold shape functions' local third derivatives.
* ShapefunctionsLocalGradients function return this
* type as its result.
*/
typedef typename BaseType::ShapeFunctionsThirdDerivativesType
ShapeFunctionsThirdDerivativesType;
/**
* Type of the normal vector used for normal to edges in geomety.
*/
typedef typename BaseType::NormalType NormalType;
/**
* Life Cycle
*/
// Quadrilateral3D9( const PointType& Point1, const PointType& Point2,
// const PointType& Point3, const PointType& Point4,
// const PointType& Point5, const PointType& Point6,
// const PointType& Point7, const PointType& Point8,
// const PointType& Point9 )
// : BaseType( PointsArrayType(), &msGeometryData )
// {
// this->Points().push_back( typename PointType::Pointer( new PointType( Point1 ) ) );
// this->Points().push_back( typename PointType::Pointer( new PointType( Point2 ) ) );
// this->Points().push_back( typename PointType::Pointer( new PointType( Point3 ) ) );
// this->Points().push_back( typename PointType::Pointer( new PointType( Point4 ) ) );
// this->Points().push_back( typename PointType::Pointer( new PointType( Point5 ) ) );
// this->Points().push_back( typename PointType::Pointer( new PointType( Point6 ) ) );
// this->Points().push_back( typename PointType::Pointer( new PointType( Point7 ) ) );
// this->Points().push_back( typename PointType::Pointer( new PointType( Point8 ) ) );
// this->Points().push_back( typename PointType::Pointer( new PointType( Point9 ) ) );
// }
Quadrilateral3D9( 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 )
: BaseType( PointsArrayType(), &msGeometryData )
{
this->Points().push_back( pPoint1 );
this->Points().push_back( pPoint2 );
this->Points().push_back( pPoint3 );
this->Points().push_back( pPoint4 );
this->Points().push_back( pPoint5 );
this->Points().push_back( pPoint6 );
this->Points().push_back( pPoint7 );
this->Points().push_back( pPoint8 );
this->Points().push_back( pPoint9 );
}
Quadrilateral3D9( const PointsArrayType& ThisPoints )
: BaseType( ThisPoints, &msGeometryData )
{
KRATOS_ERROR_IF( this->PointsNumber() != 9 ) << "Invalid points number. Expected 9, given " << this->PointsNumber() << std::endl;
}
/// Constructor with Geometry Id
explicit Quadrilateral3D9(
const IndexType GeometryId,
const PointsArrayType& rThisPoints
) : BaseType(GeometryId, rThisPoints, &msGeometryData)
{
KRATOS_ERROR_IF( this->PointsNumber() != 9 ) << "Invalid points number. Expected 9, given " << this->PointsNumber() << std::endl;
}
/// Constructor with Geometry Name
explicit Quadrilateral3D9(
const std::string& rGeometryName,
const PointsArrayType& rThisPoints
) : BaseType(rGeometryName, rThisPoints, &msGeometryData)
{
KRATOS_ERROR_IF(this->PointsNumber() != 9) << "Invalid points number. Expected 9, given " << this->PointsNumber() << std::endl;
}
/**
* Copy constructor.
* Constructs this geometry as a copy of given geometry.
*
* @note This copy constructor don't copy the points and new
* geometry shares points with given source geometry. It's
* obvious that any change to this new geometry's point affect
* source geometry's points too.
*/
Quadrilateral3D9( Quadrilateral3D9 const& rOther )
: BaseType( rOther )
{
}
/**
* Copy constructor from a geometry with other point type.
* Construct this geometry as a copy of given geometry which
* has different type of points. The given goemetry's
* TOtherPointType* must be implicity convertible to this
* geometry PointType.
*
* @note This copy constructor don't copy the points and new
* geometry shares points with given source geometry. It's
* obvious that any change to this new geometry's point affect
* source geometry's points too.
*/
template<class TOtherPointType> Quadrilateral3D9(
Quadrilateral3D9<TOtherPointType> const& rOther )
: BaseType( rOther )
{
}
/**
* Destructor. Does nothing
*/
~Quadrilateral3D9() override {}
/**
* @brief Gets the geometry family.
* @details This function returns the family type of the geometry. The geometry family categorizes the geometry into a broader classification, aiding in its identification and processing.
* @return GeometryData::KratosGeometryFamily The geometry family.
*/
GeometryData::KratosGeometryFamily GetGeometryFamily() const override
{
return GeometryData::KratosGeometryFamily::Kratos_Quadrilateral;
}
/**
* @brief Gets the geometry type.
* @details This function returns the specific type of the geometry. The geometry type provides a more detailed classification of the geometry.
* @return GeometryData::KratosGeometryType The specific geometry type.
*/
GeometryData::KratosGeometryType GetGeometryType() const override
{
return GeometryData::KratosGeometryType::Kratos_Quadrilateral3D9;
}
/**
* @brief Gets the geometry order type.
* @details This function returns the order type of the geometry. The order type relates to the polynomial degree of the geometry.
* @return GeometryData::KratosGeometryOrderType The geometry order type.
*/
GeometryData::KratosGeometryOrderType GetGeometryOrderType() const override
{
return GeometryData::KratosGeometryOrderType::Kratos_Quadratic_Order;
}
/**
* Operators
*/
/**
* Assignment operator.
*
* @note This operator don't copy the points and this
* geometry shares points with given source geometry. It's
* obvious that any change to this geometry's point affect
* source geometry's points too.
*
* @see Clone
* @see ClonePoints
*/
Quadrilateral3D9& operator=( const Quadrilateral3D9& rOther )
{
BaseType::operator=( rOther );
return *this;
}
/**
* Assignment operator for geometries with different point type.
*
* @note This operator don't copy the points and this
* geometry shares points with given source geometry. It's
* obvious that any change to this geometry's point affect
* source geometry's points too.
*
* @see Clone
* @see ClonePoints
*/
template<class TOtherPointType>
Quadrilateral3D9& operator=( Quadrilateral3D9<TOtherPointType> const & rOther )
{
BaseType::operator=( rOther );
return *this;
}
///@}
///@name Operations
///@{
/**
* @brief Creates a new geometry pointer
* @param NewGeometryId the ID of the new geometry
* @param rThisPoints the nodes of the new geometry
* @return Pointer to the new geometry
*/
typename BaseType::Pointer Create(
const IndexType NewGeometryId,
PointsArrayType const& rThisPoints
) const override
{
return typename BaseType::Pointer( new Quadrilateral3D9(NewGeometryId, rThisPoints ) );
}
/**
* @brief Creates a new geometry pointer
* @param NewGeometryId the ID of the new geometry
* @param rGeometry reference to an existing geometry
* @return Pointer to the new geometry
*/
typename BaseType::Pointer Create(
const IndexType NewGeometryId,
const BaseType& rGeometry
) const override
{
auto p_geometry = typename BaseType::Pointer( new Quadrilateral3D9( NewGeometryId, rGeometry.Points() ) );
p_geometry->SetData(rGeometry.GetData());
return p_geometry;
}
/**
* Informations
*/
/// Returns number of points per direction.
SizeType PointsNumberInDirection(IndexType LocalDirectionIndex) const override
{
if ((LocalDirectionIndex == 0) || (LocalDirectionIndex == 1)) {
return 3;
}
KRATOS_ERROR << "Possible direction index reaches from 0-1. Given direction index: "
<< LocalDirectionIndex << std::endl;
}
/**
* :TODO: the charactereistic sizes have to be reviewed
* by the one who is willing to use them!
*/
/** This method calculates and returns Length or charactereistic
* length of this geometry depending on it's dimension. For one
* dimensional geometry for example Line it returns length of it
* and for the other geometries it gives Characteristic length
* otherwise.
*
* @return double value contains length or Characteristic
* length
* @see Area()
* @see Volume()
* @see DomainSize()
*/
double Length() const override
{
Vector d = this->Points()[2] - this->Points()[0];
return( std::sqrt( d[0]*d[0] + d[1]*d[1] + d[2]*d[2] ) );
}
/** This method calculates and returns area or surface area of
* this geometry depending to it's dimension. For one dimensional
* geometry it returns zero, for two dimensional it gives area
* and for three dimensional geometries it gives surface area.
*
* @return double value contains area or surface
* area.N
* @see Length()
* @see Volume()
* @see DomainSize()
* @todo could be replaced by something more suitable (comment by janosch)
*/
double Area() const override
{
const IntegrationMethod integration_method = msGeometryData.DefaultIntegrationMethod();
return IntegrationUtilities::ComputeDomainSize(*this, integration_method);
}
/**
* @brief This method calculates and returns the volume of this geometry.
* @return Error, the volume of a 2D geometry is not defined (In June 2023)
* @see Length()
* @see Area()
* @see Volume()
*/
double Volume() const override
{
KRATOS_WARNING("Quadrilateral3D9") << "Method not well defined. Replace with DomainSize() instead. This method preserves current behaviour but will be changed in June 2023 (returning error instead)" << std::endl;
return Area();
// TODO: Replace in June 2023
// KRATOS_ERROR << "Quadrilateral3D9:: Method not well defined. Replace with DomainSize() instead." << std::endl;
// return 0.0;
}
/** This method calculates and returns length, area or volume of
* this geometry depending to it's dimension. For one dimensional
* geometry it returns its length, for two dimensional it gives area
* and for three dimensional geometries it gives its volume.
*
* @return double value contains length, area or volume.
* @see Length()
* @see Area()
* @see Volume()
* @todo could be replaced by something more suitable (comment by janosch)
*/
double DomainSize() const override
{
return Area();
}
/**
* @brief Returns whether given arbitrary point is inside the Geometry and the respective
* local point for the given global point
* @param rPoint The point to be checked if is inside o note in global coordinates
* @param rResult The local coordinates of the point
* @param Tolerance The tolerance that will be considered to check if the point is inside or not
* @return True if the point is inside, false otherwise
*/
bool IsInside(
const CoordinatesArrayType& rPoint,
CoordinatesArrayType& rResult,
const double Tolerance = std::numeric_limits<double>::epsilon()
) const override
{
PointLocalCoordinates( rResult, rPoint );
if ( std::abs( rResult[0] ) <= (1.0 + Tolerance) )
{
if ( std::abs( rResult[1] ) <= (1.0 + Tolerance) )
{
return true;
}
}
return false;
}
/**
* @brief Returns the local coordinates of a given arbitrary point
* @param rResult The vector containing the local coordinates of the point
* @param rPoint The point in global coordinates
* @return The vector containing the local coordinates of the point
*/
CoordinatesArrayType& PointLocalCoordinates(
CoordinatesArrayType& rResult,
const CoordinatesArrayType& rPoint
) const override
{
const double tol = 1.0e-8;
const int maxiter = 1000;
//check orientation of surface
std::vector< unsigned int> orientation( 3 );
double dummy = this->GetPoint( 0 ).X();
if ( std::abs( this->GetPoint( 1 ).X() - dummy ) <= tol && std::abs( this->GetPoint( 2 ).X() - dummy ) <= tol && std::abs( this->GetPoint( 3 ).X() - dummy ) <= tol )
orientation[0] = 0;
dummy = this->GetPoint( 0 ).Y();
if ( std::abs( this->GetPoint( 1 ).Y() - dummy ) <= tol && std::abs( this->GetPoint( 2 ).Y() - dummy ) <= tol && std::abs( this->GetPoint( 3 ).Y() - dummy ) <= tol )
orientation[0] = 1;
dummy = this->GetPoint( 0 ).Z();
if ( std::abs( this->GetPoint( 1 ).Z() - dummy ) <= tol && std::abs( this->GetPoint( 2 ).Z() - dummy ) <= tol && std::abs( this->GetPoint( 3 ).Z() - dummy ) <= tol )
orientation[0] = 2;
switch ( orientation[0] )
{
case 0:
orientation[0] = 1;
orientation[1] = 2;
orientation[2] = 0;
break;
case 1:
orientation[0] = 0;
orientation[1] = 2;
orientation[2] = 1;
break;
case 2:
orientation[0] = 0;
orientation[1] = 1;
orientation[2] = 2;
break;
default:
orientation[0] = 0;
orientation[1] = 1;
orientation[2] = 2;
}
Matrix J = ZeroMatrix( 2, 2 );
Matrix invJ = ZeroMatrix( 2, 2 );
if ( rResult.size() != 2 )
rResult.resize( 2, false );
//starting with xi = 0
rResult = ZeroVector( 3 );
Vector DeltaXi = ZeroVector( 3 );
CoordinatesArrayType CurrentGlobalCoords( ZeroVector( 3 ) );
//Newton iteration:
for ( int k = 0; k < maxiter; k++ )
{
CurrentGlobalCoords = ZeroVector( 3 );
this->GlobalCoordinates( CurrentGlobalCoords, rResult );
noalias( CurrentGlobalCoords ) = rPoint - CurrentGlobalCoords;
//Caluclate Inverse of Jacobian
noalias(J) = ZeroMatrix( 2, 2 );
//derivatives of shape functions
Matrix shape_functions_gradients;
shape_functions_gradients = ShapeFunctionsLocalGradients(
shape_functions_gradients, rResult );
//Elements of jacobian matrix (e.g. J(1,1) = dX1/dXi1)
//loop over all nodes
for ( unsigned int i = 0; i < this->PointsNumber(); i++ )
{
Point dummyPoint = this->GetPoint( i );
J( 0, 0 ) += ( dummyPoint[orientation[0]] ) * ( shape_functions_gradients( i, 0 ) );
J( 0, 1 ) += ( dummyPoint[orientation[0]] ) * ( shape_functions_gradients( i, 1 ) );
J( 1, 0 ) += ( dummyPoint[orientation[1]] ) * ( shape_functions_gradients( i, 0 ) );
J( 1, 1 ) += ( dummyPoint[orientation[1]] ) * ( shape_functions_gradients( i, 1 ) );
}
//deteminant of Jacobian
double det_j = J( 0, 0 ) * J( 1, 1 ) - J( 0, 1 ) * J( 1, 0 );
//filling matrix
invJ( 0, 0 ) = ( J( 1, 1 ) ) / ( det_j );
invJ( 1, 0 ) = -( J( 1, 0 ) ) / ( det_j );
invJ( 0, 1 ) = -( J( 0, 1 ) ) / ( det_j );
invJ( 1, 1 ) = ( J( 0, 0 ) ) / ( det_j );
DeltaXi( 0 ) = invJ( 0, 0 ) * CurrentGlobalCoords( orientation[0] ) + invJ( 0, 1 ) * CurrentGlobalCoords( orientation[1] );
DeltaXi( 1 ) = invJ( 1, 0 ) * CurrentGlobalCoords( orientation[0] ) + invJ( 1, 1 ) * CurrentGlobalCoords( orientation[1] );
noalias( rResult ) += DeltaXi;
if ( MathUtils<double>::Norm3( DeltaXi ) > 30 )
{
break;
}
if ( MathUtils<double>::Norm3( DeltaXi ) < tol )
{
if ( !( std::abs( CurrentGlobalCoords( orientation[2] ) ) <= tol ) )
rResult( 0 ) = 2.0;
break;
}
}
return( rResult );
}
/**
* Jacobian
*/
/**
* :TODO: implemented but not yet tested
*/
/**
* Jacobians for given method.
* This method calculates the jacobians matrices in all
* integrations points of given integration method.
*
* @param ThisMethod integration method which jacobians has to
*
* @return JacobiansType a Vector of jacobian
* matrices \f$ J_i \f$ where \f$ i=1,2,...,n \f$ is the integration
* point index of given integration method.
*
* @see DeterminantOfJacobian
* @see InverseOfJacobian
*/
JacobiansType& Jacobian( JacobiansType& rResult,
IntegrationMethod ThisMethod ) const override
{
//getting derivatives of shape functions
ShapeFunctionsGradientsType shape_functions_gradients =
CalculateShapeFunctionsIntegrationPointsLocalGradients( ThisMethod );
//getting values of shape functions
Matrix shape_functions_values =
CalculateShapeFunctionsIntegrationPointsValues( ThisMethod );
//workaround by riccardo...
if ( rResult.size() != this->IntegrationPointsNumber( ThisMethod ) )
{
// KLUDGE: While there is a bug in ublas
// vector resize, I have to put this beside resizing!!
JacobiansType temp( this->IntegrationPointsNumber( ThisMethod ) );
rResult.swap( temp );
}
//loop over all integration points
for ( unsigned int pnt = 0; pnt < this->IntegrationPointsNumber( ThisMethod ); pnt++ )
{
//defining single jacobian matrix
Matrix jacobian = ZeroMatrix( 3, 2 );
//loop over all nodes
for ( unsigned int i = 0; i < this->PointsNumber(); i++ )
{
jacobian( 0, 0 ) += ( this->GetPoint( i ).X() ) * ( shape_functions_gradients[pnt]( i, 0 ) );
jacobian( 0, 1 ) += ( this->GetPoint( i ).X() ) * ( shape_functions_gradients[pnt]( i, 1 ) );
jacobian( 1, 0 ) += ( this->GetPoint( i ).Y() ) * ( shape_functions_gradients[pnt]( i, 0 ) );
jacobian( 1, 1 ) += ( this->GetPoint( i ).Y() ) * ( shape_functions_gradients[pnt]( i, 1 ) );
jacobian( 2, 0 ) += ( this->GetPoint( i ).Z() ) * ( shape_functions_gradients[pnt]( i, 0 ) );
jacobian( 2, 1 ) += ( this->GetPoint( i ).Z() ) * ( shape_functions_gradients[pnt]( i, 1 ) );
}
rResult[pnt] = jacobian;
}//end of loop over all integration points
return rResult;
}
/**
* Jacobians for given method.
* This method calculates the jacobians matrices in all
* integrations points of given integration method.
*
* @param ThisMethod integration method which jacobians has to
*
* @return JacobiansType a Vector of jacobian
* matrices \f$ J_i \f$ where \f$ i=1,2,...,n \f$ is the integration
* point index of given integration method.
*
* @param DeltaPosition Matrix with the nodes position increment which describes
* the configuration where the jacobian has to be calculated.
*
* @see DeterminantOfJacobian
* @see InverseOfJacobian
*/
JacobiansType& Jacobian( JacobiansType& rResult,
IntegrationMethod ThisMethod,
Matrix & DeltaPosition ) const override
{
//getting derivatives of shape functions
ShapeFunctionsGradientsType shape_functions_gradients =
CalculateShapeFunctionsIntegrationPointsLocalGradients( ThisMethod );
//getting values of shape functions
Matrix shape_functions_values =
CalculateShapeFunctionsIntegrationPointsValues( ThisMethod );
//workaround by riccardo...
if ( rResult.size() != this->IntegrationPointsNumber( ThisMethod ) )
{
// KLUDGE: While there is a bug in ublas
// vector resize, I have to put this beside resizing!!
JacobiansType temp( this->IntegrationPointsNumber( ThisMethod ) );
rResult.swap( temp );
}
//loop over all integration points
for ( unsigned int pnt = 0; pnt < this->IntegrationPointsNumber( ThisMethod ); pnt++ )
{
//defining single jacobian matrix
Matrix jacobian = ZeroMatrix( 3, 2 );
//loop over all nodes
for ( unsigned int i = 0; i < this->PointsNumber(); i++ )
{
jacobian( 0, 0 ) += ( this->GetPoint( i ).X() - DeltaPosition(i,0) ) * ( shape_functions_gradients[pnt]( i, 0 ) );
jacobian( 0, 1 ) += ( this->GetPoint( i ).X() - DeltaPosition(i,0) ) * ( shape_functions_gradients[pnt]( i, 1 ) );
jacobian( 1, 0 ) += ( this->GetPoint( i ).Y() - DeltaPosition(i,1) ) * ( shape_functions_gradients[pnt]( i, 0 ) );
jacobian( 1, 1 ) += ( this->GetPoint( i ).Y() - DeltaPosition(i,1) ) * ( shape_functions_gradients[pnt]( i, 1 ) );
jacobian( 2, 0 ) += ( this->GetPoint( i ).Z() - DeltaPosition(i,2) ) * ( shape_functions_gradients[pnt]( i, 0 ) );
jacobian( 2, 1 ) += ( this->GetPoint( i ).Z() - DeltaPosition(i,2) ) * ( shape_functions_gradients[pnt]( i, 1 ) );
}
rResult[pnt] = jacobian;
}//end of loop over all integration points
return rResult;
}
/**
* :TODO: implemented but not yet tested
*/
/**
* Jacobian in specific integration point of given integration
* method. This method calculate jacobian matrix in given
* integration point of given integration method.
*
* @param IntegrationPointIndex index of integration point which jacobians has to
* be calculated in it.
*
* @param ThisMethod integration method which jacobians has to
* be calculated in its integration points.
*
* @return Matrix(double) Jacobian matrix \f$ J_i \f$ where \f$
* i \f$ is the given integration point index of given
* integration method.
*
* @see DeterminantOfJacobian
* @see InverseOfJacobian
*/
Matrix& Jacobian( Matrix& rResult,
IndexType IntegrationPointIndex,
IntegrationMethod ThisMethod ) const override
{
// Setting up size of jacobian matrix
if (rResult.size1() != 3 || rResult.size2() != 2)
rResult.resize( 3, 2 , false);
rResult.clear();
//derivatives of shape functions
ShapeFunctionsGradientsType shape_functions_gradients =
CalculateShapeFunctionsIntegrationPointsLocalGradients( ThisMethod );
Matrix ShapeFunctionsGradientInIntegrationPoint = shape_functions_gradients( IntegrationPointIndex );
//values of shape functions in integration points
DenseVector<double> ShapeFunctionsValuesInIntegrationPoint = ZeroVector( 9 );
/*vector<double>*/
ShapeFunctionsValuesInIntegrationPoint = row(
CalculateShapeFunctionsIntegrationPointsValues( ThisMethod ), IntegrationPointIndex );
//Elements of jacobian matrix (e.g. J(1,1) = dX1/dXi1)
//loop over all nodes
for ( unsigned int i = 0; i < this->PointsNumber(); i++ )
{
rResult( 0, 0 ) += ( this->GetPoint( i ).X() ) * ( ShapeFunctionsGradientInIntegrationPoint( i, 0 ) );
rResult( 0, 1 ) += ( this->GetPoint( i ).X() ) * ( ShapeFunctionsGradientInIntegrationPoint( i, 1 ) );
rResult( 1, 0 ) += ( this->GetPoint( i ).Y() ) * ( ShapeFunctionsGradientInIntegrationPoint( i, 0 ) );
rResult( 1, 1 ) += ( this->GetPoint( i ).Y() ) * ( ShapeFunctionsGradientInIntegrationPoint( i, 1 ) );
rResult( 2, 0 ) += ( this->GetPoint( i ).Z() ) * ( ShapeFunctionsGradientInIntegrationPoint( i, 0 ) );
rResult( 2, 1 ) += ( this->GetPoint( i ).Z() ) * ( ShapeFunctionsGradientInIntegrationPoint( i, 1 ) );
}
return rResult;
}
/**
* :TODO: implemented but not yet tested
*/
/**
* Jacobian in given point. This method calculate jacobian
* matrix in given point.
*
* @param rPoint point which jacobians has to
* be calculated in it.
* @return Matrix of double which is jacobian matrix \f$ J \f$ in given point.
*
* @see DeterminantOfJacobian
* @see InverseOfJacobian
*/
Matrix& Jacobian( Matrix& rResult, const CoordinatesArrayType& rPoint ) const override
{
// Setting up size of jacobian matrix
if (rResult.size1() != 3 || rResult.size2() != 2)
rResult.resize( 3, 2 , false);
rResult.clear();
//derivatives of shape functions
Matrix shape_functions_gradients;
shape_functions_gradients = ShapeFunctionsLocalGradients(
shape_functions_gradients, rPoint );
//Elements of jacobian matrix (e.g. J(1,1) = dX1/dXi1)
//loop over all nodes
for ( unsigned int i = 0; i < this->PointsNumber(); i++ )
{
const auto& coordinates = this->GetPoint(i).Coordinates();
rResult( 0, 0 ) += ( coordinates[0] ) * ( shape_functions_gradients( i, 0 ) );
rResult( 0, 1 ) += ( coordinates[0] ) * ( shape_functions_gradients( i, 1 ) );
rResult( 1, 0 ) += ( coordinates[1] ) * ( shape_functions_gradients( i, 0 ) );
rResult( 1, 1 ) += ( coordinates[1] ) * ( shape_functions_gradients( i, 1 ) );
rResult( 2, 0 ) += ( coordinates[2] ) * ( shape_functions_gradients( i, 0 ) );
rResult( 2, 1 ) += ( coordinates[2] ) * ( shape_functions_gradients( i, 1 ) );
}
return rResult;
}
///@}
///@name Edge
///@{
/**
* @brief This method gives you number of all edges of this geometry.
* @details For example, for a hexahedron, this would be 12
* @return SizeType containes number of this geometry edges.
* @see EdgesNumber()
* @see Edges()
* @see GenerateEdges()
* @see FacesNumber()
* @see Faces()
* @see GenerateFaces()
*/
SizeType EdgesNumber() const override
{
return 4;
}
/**
* @brief This method gives you all edges of this geometry.
* @details This method will gives you all the edges with one dimension less than this geometry.
* For example a triangle would return three lines as its edges or a tetrahedral would return four triangle as its edges but won't return its six edge lines by this method.
* @return GeometriesArrayType containes this geometry edges.
* @see EdgesNumber()
* @see Edge()
*/
GeometriesArrayType GenerateEdges() const override
{
GeometriesArrayType edges = GeometriesArrayType();
edges.push_back( Kratos::make_shared<EdgeType>( this->pGetPoint( 0 ), this->pGetPoint( 1 ), this->pGetPoint( 4 ) ) );
edges.push_back( Kratos::make_shared<EdgeType>( this->pGetPoint( 1 ), this->pGetPoint( 2 ), this->pGetPoint( 5 ) ) );
edges.push_back( Kratos::make_shared<EdgeType>( this->pGetPoint( 2 ), this->pGetPoint( 3 ), this->pGetPoint( 6 ) ) );
edges.push_back( Kratos::make_shared<EdgeType>( this->pGetPoint( 3 ), this->pGetPoint( 0 ), this->pGetPoint( 7 ) ) );
return edges;
}
/**
* Shape Function
*/
/**
* :TODO: implemented but not yet tested
*/
/**
* Calculates the value of a given shape function at a given point.
*
* @param ShapeFunctionIndex The number of the desired shape function
* @param rPoint the given point in local coordinates at which the
* value of the shape function is calculated
*
* @return the value of the shape function at the given point
*/
double ShapeFunctionValue( IndexType ShapeFunctionIndex,
const CoordinatesArrayType& rPoint ) const override
{
double fx1 = 0.5 * ( rPoint[0] - 1 ) * rPoint[0];
double fx2 = 0.5 * ( rPoint[0] + 1 ) * rPoint[0];
double fx3 = 1 - rPoint[0] * rPoint[0];
double fy1 = 0.5 * ( rPoint[1] - 1 ) * rPoint[1];
double fy2 = 0.5 * ( rPoint[1] + 1 ) * rPoint[1];
double fy3 = 1 - rPoint[1] * rPoint[1];
switch ( ShapeFunctionIndex )
{
case 0:
return( fx1*fy1 );
case 1:
return( fx2*fy1 );
case 2:
return( fx2*fy2 );
case 3:
return( fx1*fy2 );
case 4:
return( fx3*fy1 );
case 5:
return( fx2*fy3 );
case 6:
return( fx3*fy2 );
case 7:
return( fx1*fy3 );
case 8:
return( fx3*fy3 );
default:
KRATOS_ERROR << "Wrong index of shape function!" << *this << std::endl;
}
return 0;
}
/** This method gives all non-zero shape functions values
evaluated at the rCoordinates provided
@return Vector of values of shape functions \f$ F_{i} \f$
where i is the shape function index (for NURBS it is the index
of the local enumeration in the element).
@see ShapeFunctionValue
@see ShapeFunctionsLocalGradients
@see ShapeFunctionLocalGradient
*/
Vector& ShapeFunctionsValues(Vector &rResult, const CoordinatesArrayType& rCoordinates) const override
{
if(rResult.size() != 9) rResult.resize(9,false);
double fx1 = 0.5 * ( rCoordinates[0] - 1 ) * rCoordinates[0];
double fx2 = 0.5 * ( rCoordinates[0] + 1 ) * rCoordinates[0];
double fx3 = 1 - rCoordinates[0] * rCoordinates[0];
double fy1 = 0.5 * ( rCoordinates[1] - 1 ) * rCoordinates[1];
double fy2 = 0.5 * ( rCoordinates[1] + 1 ) * rCoordinates[1];
double fy3 = 1 - rCoordinates[1] * rCoordinates[1];
rResult[0] = fx1*fy1 ;
rResult[1] = fx2*fy1 ;
rResult[2] = fx2*fy2 ;
rResult[3] = fx1*fy2 ;
rResult[4] = fx3*fy1 ;
rResult[5] = fx2*fy3 ;
rResult[6] = fx3*fy2 ;
rResult[7] = fx1*fy3 ;
rResult[8] = fx3*fy3 ;
return rResult;
}
///@}
/**
* Input and Output
*/
/**
* Turn back information as a string.
* @return String contains information about this geometry.
* @see PrintData()
* @see PrintInfo()
*/
std::string Info() const override
{
return "2 dimensional quadrilateral with nine nodes in 3D space";
}
/**