-
Notifications
You must be signed in to change notification settings - Fork 93
/
Copy pathinterpolations.jl
1776 lines (1531 loc) · 73.2 KB
/
interpolations.jl
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
"""
Interpolation{ref_shape, order}()
Abstract type for interpolations defined on `ref_shape`
(see [`AbstractRefShape`](@ref)).
`order` corresponds to the order of the interpolation.
The interpolation is used to define shape functions to interpolate
a function between nodes.
The following interpolations are implemented:
* `Lagrange{RefLine,1}`
* `Lagrange{RefLine,2}`
* `Lagrange{RefQuadrilateral,1}`
* `Lagrange{RefQuadrilateral,2}`
* `Lagrange{RefQuadrilateral,3}`
* `Lagrange{RefTriangle,1}`
* `Lagrange{RefTriangle,2}`
* `Lagrange{RefTriangle,3}`
* `Lagrange{RefTriangle,4}`
* `Lagrange{RefTriangle,5}`
* `BubbleEnrichedLagrange{RefTriangle,1}`
* `CrouzeixRaviart{RefTriangle, 1}`
* `CrouzeixRaviart{RefTetrahedron, 1}`
* `RannacherTurek{RefQuadrilateral, 1}`
* `RannacherTurek{RefHexahedron, 1}`
* `Lagrange{RefHexahedron,1}`
* `Lagrange{RefHexahedron,2}`
* `Lagrange{RefTetrahedron,1}`
* `Lagrange{RefTetrahedron,2}`
* `Lagrange{RefPrism,1}`
* `Lagrange{RefPrism,2}`
* `Lagrange{RefPyramid,1}`
* `Lagrange{RefPyramid,2}`
* `Serendipity{RefQuadrilateral,2}`
* `Serendipity{RefHexahedron,2}`
# Examples
```jldoctest
julia> ip = Lagrange{RefTriangle, 2}()
Lagrange{RefTriangle, 2}()
julia> getnbasefunctions(ip)
6
```
"""
abstract type Interpolation{shape #=<: AbstractRefShape=#, order} end
const InterpolationByDim{dim} = Interpolation{<:AbstractRefShape{dim}}
abstract type ScalarInterpolation{refshape, order} <: Interpolation{refshape, order} end
abstract type VectorInterpolation{vdim, refshape, order} <: Interpolation{refshape, order} end
# Number of components for the interpolation.
n_components(::ScalarInterpolation) = 1
n_components(::VectorInterpolation{vdim}) where {vdim} = vdim
# Number of components that are allowed to prescribe in e.g. Dirichlet BC
n_dbc_components(ip::Interpolation) = n_components(ip)
"""
shape_value_type(ip::Interpolation, ::Type{T}) where T<:Number
Return the type of `shape_value(ip::Interpolation, ξ::Vec, ib::Int)`.
"""
shape_value_type(::Interpolation, ::Type{T}) where {T <: Number}
shape_value_type(::ScalarInterpolation, ::Type{T}) where {T <: Number} = T
shape_value_type(::VectorInterpolation{vdim}, ::Type{T}) where {vdim, T <: Number} = Vec{vdim, T}
#shape_value_type(::MatrixInterpolation, T::Type) = Tensor #958
# TODO: Add a fallback that errors if there are multiple dofs per edge/face instead to force
# interpolations to opt-out instead of silently do nothing.
"""
adjust_dofs_during_distribution(::Interpolation)
This function must return `true` if the dofs should be adjusted (i.e. permuted) during dof
distribution. This is in contrast to i) adjusting the dofs during [`reinit!`](@ref) in the
assembly loop, or ii) not adjusting at all (which is not needed for low order
interpolations, generally).
"""
adjust_dofs_during_distribution(::Interpolation)
"""
InterpolationInfo
Gathers all the information needed to distribute dofs for a given interpolation. Note that
this cache is of the same type no matter the interpolation: the purpose is to make
dof-distribution type-stable.
"""
struct InterpolationInfo
nvertexdofs::Vector{Int}
nedgedofs::Vector{Int}
nfacedofs::Vector{Int}
nvolumedofs::Int
reference_dim::Int
adjust_during_distribution::Bool
n_copies::Int
is_discontinuous::Bool
end
function InterpolationInfo(interpolation::Interpolation{shape}, n_copies) where {rdim, shape <: AbstractRefShape{rdim}}
info = InterpolationInfo(
[length(i) for i in vertexdof_indices(interpolation)],
[length(i) for i in edgedof_interior_indices(interpolation)],
[length(i) for i in facedof_interior_indices(interpolation)],
length(volumedof_interior_indices(interpolation)),
rdim,
adjust_dofs_during_distribution(interpolation),
n_copies,
is_discontinuous(interpolation)
)
return info
end
InterpolationInfo(interpolation::Interpolation) = InterpolationInfo(interpolation, 1)
nvertices(::Interpolation{RefShape}) where {RefShape} = nvertices(RefShape)
nedges(::Interpolation{RefShape}) where {RefShape} = nedges(RefShape)
nfaces(::Interpolation{RefShape}) where {RefShape} = nfaces(RefShape)
Base.copy(ip::Interpolation) = ip
"""
Ferrite.getrefdim(::Interpolation)
Return the dimension of the reference element for a given interpolation.
"""
getrefdim(::Interpolation) # To make doc-filtering work
@inline getrefdim(::Interpolation{RefShape}) where {RefShape} = getrefdim(RefShape)
"""
Ferrite.getrefshape(::Interpolation)::AbstractRefShape
Return the reference element shape of the interpolation.
"""
@inline getrefshape(::Interpolation{shape}) where {shape} = shape
"""
Ferrite.getorder(::Interpolation)
Return order of the interpolation.
"""
@inline getorder(::Interpolation{shape, order}) where {shape, order} = order
#####################
# Utility functions #
#####################
"""
Ferrite.getnbasefunctions(ip::Interpolation)
Return the number of base functions for the interpolation `ip`.
"""
getnbasefunctions(::Interpolation)
# The following functions are used to distribute the dofs. Definitions:
# vertexdof: dof on a "corner" of the reference shape
# facedof: dof in the dim-1 dimension (line in 2D, surface in 3D)
# edgedof: dof on a line between 2 vertices (i.e. "corners") (3D only)
# celldof: dof that is local to the element
"""
reference_shape_values!(values::AbstractArray{T}, ip::Interpolation, ξ::Vec)
Evaluate all shape functions of `ip` at once at the reference point `ξ` and store them in
`values`.
"""
@propagate_inbounds function reference_shape_values!(values::AT, ip::IP, ξ::Vec) where {IP <: Interpolation, AT <: AbstractArray}
@boundscheck checkbounds(values, 1:getnbasefunctions(ip))
@inbounds for i in 1:getnbasefunctions(ip)
values[i] = reference_shape_value(ip, ξ, i)
end
return
end
"""
reference_shape_gradients!(gradients::AbstractArray, ip::Interpolation, ξ::Vec)
Evaluate all shape function gradients of `ip` at once at the reference point `ξ` and store
them in `gradients`.
"""
function reference_shape_gradients!(gradients::AT, ip::IP, ξ::Vec) where {IP <: Interpolation, AT <: AbstractArray}
@boundscheck checkbounds(gradients, 1:getnbasefunctions(ip))
@inbounds for i in 1:getnbasefunctions(ip)
gradients[i] = reference_shape_gradient(ip, ξ, i)
end
return
end
"""
reference_shape_gradients_and_values!(gradients::AbstractArray, values::AbstractArray, ip::Interpolation, ξ::Vec)
Evaluate all shape function gradients and values of `ip` at once at the reference point `ξ`
and store them in `values`.
"""
function reference_shape_gradients_and_values!(gradients::GAT, values::SAT, ip::IP, ξ::Vec) where {IP <: Interpolation, SAT <: AbstractArray, GAT <: AbstractArray}
@boundscheck checkbounds(gradients, 1:getnbasefunctions(ip))
@boundscheck checkbounds(values, 1:getnbasefunctions(ip))
@inbounds for i in 1:getnbasefunctions(ip)
gradients[i], values[i] = reference_shape_gradient_and_value(ip, ξ, i)
end
return
end
"""
reference_shape_hessians_gradients_and_values!(hessians::AbstractVector, gradients::AbstractVector, values::AbstractVector, ip::Interpolation, ξ::Vec)
Evaluate all shape function hessians, gradients and values of `ip` at once at the reference point `ξ`
and store them in `hessians`, `gradients`, and `values`.
"""
@propagate_inbounds function reference_shape_hessians_gradients_and_values!(hessians::AbstractVector, gradients::AbstractVector, values::AbstractVector, ip::Interpolation, ξ::Vec)
@boundscheck checkbounds(hessians, 1:getnbasefunctions(ip))
@boundscheck checkbounds(gradients, 1:getnbasefunctions(ip))
@boundscheck checkbounds(values, 1:getnbasefunctions(ip))
@inbounds for i in 1:getnbasefunctions(ip)
hessians[i], gradients[i], values[i] = reference_shape_hessian_gradient_and_value(ip, ξ, i)
end
return
end
"""
reference_shape_value(ip::Interpolation, ξ::Vec, i::Int)
Evaluate the value of the `i`th shape function of the interpolation `ip`
at a point `ξ` on the reference element. The index `i` must
match the index in [`vertices(::Interpolation)`](@ref), [`faces(::Interpolation)`](@ref) and
[`edges(::Interpolation)`](@ref).
For nodal interpolations the indices also must match the
indices of [`reference_coordinates(::Interpolation)`](@ref).
"""
reference_shape_value(ip::Interpolation, ξ::Vec, i::Int)
"""
reference_shape_gradient(ip::Interpolation, ξ::Vec, i::Int)
Evaluate the gradient of the `i`th shape function of the interpolation `ip` in
reference coordinate `ξ`.
"""
function reference_shape_gradient(ip::Interpolation, ξ::Vec, i::Int)
return Tensors.gradient(x -> reference_shape_value(ip, x, i), ξ)
end
"""
reference_shape_gradient_and_value(ip::Interpolation, ξ::Vec, i::Int)
Optimized version combining the evaluation [`Ferrite.reference_shape_value(::Interpolation)`](@ref)
and [`Ferrite.reference_shape_gradient(::Interpolation)`](@ref).
"""
function reference_shape_gradient_and_value(ip::Interpolation, ξ::Vec, i::Int)
return gradient(x -> reference_shape_value(ip, x, i), ξ, :all)
end
"""
reference_shape_hessian_gradient_and_value(ip::Interpolation, ξ::Vec, i::Int)
Optimized version combining the evaluation [`Ferrite.reference_shape_value(::Interpolation)`](@ref),
[`Ferrite.reference_shape_gradient(::Interpolation)`](@ref), and the gradient of the latter.
"""
function reference_shape_hessian_gradient_and_value(ip::Interpolation, ξ::Vec, i::Int)
return hessian(x -> reference_shape_value(ip, x, i), ξ, :all)
end
"""
reference_coordinates(ip::Interpolation)
Returns a vector of coordinates with length [`getnbasefunctions(::Interpolation)`](@ref)
and indices corresponding to the indices of a dof in [`vertices`](@ref), [`faces`](@ref) and
[`edges`](@ref).
Only required for nodal interpolations.
TODO: Separate nodal and non-nodal interpolations.
"""
reference_coordinates(::Interpolation)
"""
vertexdof_indices(ip::Interpolation)
A tuple containing tuples of local dof indices for the respective vertex in local
enumeration on a cell defined by [`vertices(::Cell)`](@ref). The vertex enumeration must
match the vertex enumeration of the corresponding geometrical cell.
!!! note
The dofs appearing in the tuple must be continuous and increasing! The first dof must be
the 1, as vertex dofs are enumerated first.
"""
vertexdof_indices(ip::Interpolation) = ntuple(_ -> (), nvertices(ip))
"""
dirichlet_vertexdof_indices(ip::Interpolation)
A tuple containing tuples of local dof indices for the respective vertex in local
enumeration on a cell defined by [`vertices(::Cell)`](@ref). The vertex enumeration must
match the vertex enumeration of the corresponding geometrical cell.
Used internally in [`ConstraintHandler`](@ref) and defaults to [`vertexdof_indices(ip::Interpolation)`](@ref) for continuous interpolation.
!!! note
The dofs appearing in the tuple must be continuous and increasing! The first dof must be
the 1, as vertex dofs are enumerated first.
"""
dirichlet_vertexdof_indices(ip::Interpolation) = vertexdof_indices(ip)
"""
edgedof_indices(ip::Interpolation)
A tuple containing tuples of local dof indices for the respective edge in local enumeration
on a cell defined by [`edges(::Cell)`](@ref). The edge enumeration must match the edge
enumeration of the corresponding geometrical cell.
The dofs are guaranteed to be aligned with the local ordering of the entities on the oriented edge.
Here the first entries are the vertex dofs, followed by the edge interior dofs.
"""
edgedof_indices(::Interpolation)
"""
dirichlet_edgedof_indices(ip::Interpolation)
A tuple containing tuples of local dof indices for the respective edge in local enumeration
on a cell defined by [`edges(::Cell)`](@ref). The edge enumeration must match the edge
enumeration of the corresponding geometrical cell.
Used internally in [`ConstraintHandler`](@ref) and defaults to [`edgedof_indices(ip::Interpolation)`](@ref) for continuous interpolation.
The dofs are guaranteed to be aligned with the local ordering of the entities on the oriented edge.
Here the first entries are the vertex dofs, followed by the edge interior dofs.
"""
dirichlet_edgedof_indices(ip::Interpolation) = edgedof_indices(ip)
"""
edgedof_interior_indices(ip::Interpolation)
A tuple containing tuples of the local dof indices on the interior of the respective edge in
local enumeration on a cell defined by [`edges(::Cell)`](@ref). The edge enumeration must
match the edge enumeration of the corresponding geometrical cell. Note that the vertex dofs
are included here.
!!! note
The dofs appearing in the tuple must be continuous and increasing! The first dof must be
computed via "last vertex dof index + 1", if edge dofs exist.
"""
edgedof_interior_indices(::Interpolation)
"""
facedof_indices(ip::Interpolation)
A tuple containing tuples of all local dof indices for the respective face in local
enumeration on a cell defined by [`faces(::Cell)`](@ref). The face enumeration must match
the face enumeration of the corresponding geometrical cell.
"""
facedof_indices(::Interpolation)
"""
dirichlet_facedof_indices(ip::Interpolation)
A tuple containing tuples of all local dof indices for the respective face in local
enumeration on a cell defined by [`faces(::Cell)`](@ref). The face enumeration must match
the face enumeration of the corresponding geometrical cell.
Used internally in [`ConstraintHandler`](@ref) and defaults to [`facedof_indices(ip::Interpolation)`](@ref) for continuous interpolation.
"""
dirichlet_facedof_indices(ip::Interpolation) = facedof_indices(ip)
"""
facedof_interior_indices(ip::Interpolation)
A tuple containing tuples of the local dof indices on the interior of the respective face in
local enumeration on a cell defined by [`faces(::Cell)`](@ref). The face enumeration must
match the face enumeration of the corresponding geometrical cell. Note that the vertex and
edge dofs are included here.
!!! note
The dofs appearing in the tuple must be continuous and increasing! The first dof must be
the computed via "last edge interior dof index + 1", if face dofs exist.
"""
facedof_interior_indices(::Interpolation)
"""
volumedof_interior_indices(ip::Interpolation)
Tuple containing the dof indices associated with the interior of a volume.
!!! note
The dofs appearing in the tuple must be continuous and increasing, volumedofs are
enumerated last.
"""
volumedof_interior_indices(::Interpolation) = ()
# Some helpers to skip boilerplate
edgedof_indices(ip::Interpolation) = ntuple(_ -> (), nedges(ip))
edgedof_interior_indices(ip::Interpolation) = ntuple(_ -> (), nedges(ip))
facedof_indices(ip::Interpolation) = ntuple(_ -> (), nfaces(ip))
facedof_interior_indices(ip::Interpolation) = ntuple(_ -> (), nfaces(ip))
"""
boundarydof_indices(::Type{<:BoundaryIndex})
Helper function to generically dispatch on the correct dof sets of a boundary entity.
"""
boundarydof_indices(::Type{<:BoundaryIndex})
boundarydof_indices(::Type{FaceIndex}) = facedof_indices
boundarydof_indices(::Type{EdgeIndex}) = edgedof_indices
boundarydof_indices(::Type{VertexIndex}) = vertexdof_indices
facetdof_indices(ip::InterpolationByDim{3}) = facedof_indices(ip)
facetdof_indices(ip::InterpolationByDim{2}) = edgedof_indices(ip)
facetdof_indices(ip::InterpolationByDim{1}) = vertexdof_indices(ip)
facetdof_interior_indices(ip::InterpolationByDim{3}) = facedof_interior_indices(ip)
facetdof_interior_indices(ip::InterpolationByDim{2}) = edgedof_interior_indices(ip)
facetdof_interior_indices(ip::InterpolationByDim{1}) = ntuple(_ -> (), nvertices(ip))
dirichlet_facetdof_indices(ip::InterpolationByDim{3}) = dirichlet_facedof_indices(ip)
dirichlet_facetdof_indices(ip::InterpolationByDim{2}) = dirichlet_edgedof_indices(ip)
dirichlet_facetdof_indices(ip::InterpolationByDim{1}) = dirichlet_vertexdof_indices(ip)
nfacets(ip::InterpolationByDim{3}) = nfaces(ip)
nfacets(ip::InterpolationByDim{2}) = nedges(ip)
nfacets(ip::InterpolationByDim{1}) = nvertices(ip)
"""
is_discontinuous(::Interpolation)
is_discontinuous(::Type{<:Interpolation})
Checks whether the interpolation is discontinuous (i.e. `DiscontinuousLagrange`)
"""
is_discontinuous(ip::Interpolation) = is_discontinuous(typeof(ip))
is_discontinuous(::Type{<:Interpolation}) = false
"""
dirichlet_boundarydof_indices(::Type{<:BoundaryIndex})
Helper function to generically dispatch on the correct dof sets of a boundary entity.
Used internally in [`ConstraintHandler`](@ref) and defaults to [`boundarydof_indices(ip::Interpolation)`](@ref) for continuous interpolation.
"""
dirichlet_boundarydof_indices(::Type{<:BoundaryIndex})
dirichlet_boundarydof_indices(::Type{FaceIndex}) = dirichlet_facedof_indices
dirichlet_boundarydof_indices(::Type{EdgeIndex}) = dirichlet_edgedof_indices
dirichlet_boundarydof_indices(::Type{VertexIndex}) = dirichlet_vertexdof_indices
dirichlet_boundarydof_indices(::Type{FacetIndex}) = dirichlet_facetdof_indices
#########################
# DiscontinuousLagrange #
#########################
# TODO generalize to arbitrary basis positionings.
"""
Piecewise discontinuous Lagrange basis via Gauss-Lobatto points.
"""
struct DiscontinuousLagrange{shape, order} <: ScalarInterpolation{shape, order}
function DiscontinuousLagrange{shape, order}() where {shape <: AbstractRefShape, order}
return new{shape, order}()
end
end
adjust_dofs_during_distribution(::DiscontinuousLagrange) = false
getlowerorder(::DiscontinuousLagrange{shape, order}) where {shape, order} = DiscontinuousLagrange{shape, order - 1}()
getnbasefunctions(::DiscontinuousLagrange{shape, order}) where {shape, order} = getnbasefunctions(Lagrange{shape, order}())
getnbasefunctions(::DiscontinuousLagrange{shape, 0}) where {shape} = 1
# This just moves all dofs into the interior of the element.
volumedof_interior_indices(ip::DiscontinuousLagrange) = ntuple(i -> i, getnbasefunctions(ip))
# Mirror the Lagrange element for now to avoid repeating.
dirichlet_facedof_indices(ip::DiscontinuousLagrange{shape, order}) where {shape, order} = dirichlet_facedof_indices(Lagrange{shape, order}())
dirichlet_edgedof_indices(ip::DiscontinuousLagrange{shape, order}) where {shape, order} = dirichlet_edgedof_indices(Lagrange{shape, order}())
dirichlet_vertexdof_indices(ip::DiscontinuousLagrange{shape, order}) where {shape, order} = dirichlet_vertexdof_indices(Lagrange{shape, order}())
# Mirror the Lagrange element for now.
function reference_coordinates(ip::DiscontinuousLagrange{shape, order}) where {shape, order}
return reference_coordinates(Lagrange{shape, order}())
end
function reference_shape_value(::DiscontinuousLagrange{shape, order}, ξ::Vec{dim}, i::Int) where {dim, shape <: AbstractRefShape{dim}, order}
return reference_shape_value(Lagrange{shape, order}(), ξ, i)
end
# Excepting the L0 element.
function reference_coordinates(ip::DiscontinuousLagrange{RefHypercube{dim}, 0}) where {dim}
return [Vec{dim, Float64}(ntuple(x -> 0.0, dim))]
end
function reference_coordinates(ip::DiscontinuousLagrange{RefTriangle, 0})
return [Vec{2, Float64}((1 / 3, 1 / 3))]
end
function reference_coordinates(ip::DiscontinuousLagrange{RefTetrahedron, 0})
return [Vec{3, Float64}((1 / 4, 1 / 4, 1 / 4))]
end
function reference_shape_value(ip::DiscontinuousLagrange{shape, 0}, ::Vec{dim, T}, i::Int) where {dim, shape <: AbstractRefShape{dim}, T}
i > 1 && throw(ArgumentError("no shape function $i for interpolation $ip"))
return one(T)
end
is_discontinuous(::Type{<:DiscontinuousLagrange}) = true
############
# Lagrange #
############
"""
Lagrange{refshape, order} <: ScalarInterpolation
Standard continuous Lagrange polynomials with equidistant node placement.
"""
struct Lagrange{shape, order} <: ScalarInterpolation{shape, order}
function Lagrange{shape, order}() where {shape <: AbstractRefShape, order}
return new{shape, order}()
end
end
adjust_dofs_during_distribution(::Lagrange) = true
adjust_dofs_during_distribution(::Lagrange{<:Any, 2}) = false
adjust_dofs_during_distribution(::Lagrange{<:Any, 1}) = false
# Vertices for all Lagrange interpolations are the same
vertexdof_indices(::Lagrange{RefLine}) = ((1,), (2,))
vertexdof_indices(::Lagrange{RefQuadrilateral}) = ((1,), (2,), (3,), (4,))
vertexdof_indices(::Lagrange{RefHexahedron}) = ((1,), (2,), (3,), (4,), (5,), (6,), (7,), (8,))
vertexdof_indices(::Lagrange{RefTriangle}) = ((1,), (2,), (3,))
vertexdof_indices(::Lagrange{RefTetrahedron}) = ((1,), (2,), (3,), (4,))
vertexdof_indices(::Lagrange{RefPrism}) = ((1,), (2,), (3,), (4,), (5,), (6,))
vertexdof_indices(::Lagrange{RefPyramid}) = ((1,), (2,), (3,), (4,), (5,))
getlowerorder(::Lagrange{shape, order}) where {shape, order} = Lagrange{shape, order - 1}()
getlowerorder(::Lagrange{shape, 1}) where {shape} = DiscontinuousLagrange{shape, 0}()
############################
# Lagrange RefLine order 1 #
############################
getnbasefunctions(::Lagrange{RefLine, 1}) = 2
edgedof_indices(::Lagrange{RefLine, 1}) = ((1, 2),)
function reference_coordinates(::Lagrange{RefLine, 1})
return [
Vec{1, Float64}((-1.0,)),
Vec{1, Float64}((1.0,)),
]
end
function reference_shape_value(ip::Lagrange{RefLine, 1}, ξ::Vec{1}, i::Int)
ξ_x = ξ[1]
i == 1 && return (1 - ξ_x) / 2
i == 2 && return (1 + ξ_x) / 2
throw(ArgumentError("no shape function $i for interpolation $ip"))
end
############################
# Lagrange RefLine order 2 #
############################
getnbasefunctions(::Lagrange{RefLine, 2}) = 3
edgedof_indices(::Lagrange{RefLine, 2}) = ((1, 2, 3),)
edgedof_interior_indices(::Lagrange{RefLine, 2}) = (3,)
function reference_coordinates(::Lagrange{RefLine, 2})
return [
Vec{1, Float64}((-1.0,)),
Vec{1, Float64}((1.0,)),
Vec{1, Float64}((0.0,)),
]
end
function reference_shape_value(ip::Lagrange{RefLine, 2}, ξ::Vec{1}, i::Int)
ξ_x = ξ[1]
i == 1 && return ξ_x * (ξ_x - 1) / 2
i == 2 && return ξ_x * (ξ_x + 1) / 2
i == 3 && return 1 - ξ_x^2
throw(ArgumentError("no shape function $i for interpolation $ip"))
end
#####################################
# Lagrange RefQuadrilateral order 1 #
#####################################
getnbasefunctions(::Lagrange{RefQuadrilateral, 1}) = 4
edgedof_indices(::Lagrange{RefQuadrilateral, 1}) = ((1, 2), (2, 3), (3, 4), (4, 1))
facedof_indices(ip::Lagrange{RefQuadrilateral, 1}) = (ntuple(i -> i, getnbasefunctions(ip)),)
function reference_coordinates(::Lagrange{RefQuadrilateral, 1})
return [
Vec{2, Float64}((-1.0, -1.0)),
Vec{2, Float64}((1.0, -1.0)),
Vec{2, Float64}((1.0, 1.0)),
Vec{2, Float64}((-1.0, 1.0)),
]
end
function reference_shape_value(ip::Lagrange{RefQuadrilateral, 1}, ξ::Vec{2}, i::Int)
ξ_x = ξ[1]
ξ_y = ξ[2]
i == 1 && return (1 - ξ_x) * (1 - ξ_y) / 4
i == 2 && return (1 + ξ_x) * (1 - ξ_y) / 4
i == 3 && return (1 + ξ_x) * (1 + ξ_y) / 4
i == 4 && return (1 - ξ_x) * (1 + ξ_y) / 4
throw(ArgumentError("no shape function $i for interpolation $ip"))
end
#####################################
# Lagrange RefQuadrilateral order 2 #
#####################################
getnbasefunctions(::Lagrange{RefQuadrilateral, 2}) = 9
edgedof_indices(::Lagrange{RefQuadrilateral, 2}) = ((1, 2, 5), (2, 3, 6), (3, 4, 7), (4, 1, 8))
edgedof_interior_indices(::Lagrange{RefQuadrilateral, 2}) = ((5,), (6,), (7,), (8,))
facedof_indices(ip::Lagrange{RefQuadrilateral, 2}) = (ntuple(i -> i, getnbasefunctions(ip)),)
facedof_interior_indices(::Lagrange{RefQuadrilateral, 2}) = ((9,))
function reference_coordinates(::Lagrange{RefQuadrilateral, 2})
return [
Vec{2, Float64}((-1.0, -1.0)),
Vec{2, Float64}((1.0, -1.0)),
Vec{2, Float64}((1.0, 1.0)),
Vec{2, Float64}((-1.0, 1.0)),
Vec{2, Float64}((0.0, -1.0)),
Vec{2, Float64}((1.0, 0.0)),
Vec{2, Float64}((0.0, 1.0)),
Vec{2, Float64}((-1.0, 0.0)),
Vec{2, Float64}((0.0, 0.0)),
]
end
function reference_shape_value(ip::Lagrange{RefQuadrilateral, 2}, ξ::Vec{2}, i::Int)
ξ_x = ξ[1]
ξ_y = ξ[2]
i == 1 && return (ξ_x^2 - ξ_x) * (ξ_y^2 - ξ_y) / 4
i == 2 && return (ξ_x^2 + ξ_x) * (ξ_y^2 - ξ_y) / 4
i == 3 && return (ξ_x^2 + ξ_x) * (ξ_y^2 + ξ_y) / 4
i == 4 && return (ξ_x^2 - ξ_x) * (ξ_y^2 + ξ_y) / 4
i == 5 && return (1 - ξ_x^2) * (ξ_y^2 - ξ_y) / 2
i == 6 && return (ξ_x^2 + ξ_x) * (1 - ξ_y^2) / 2
i == 7 && return (1 - ξ_x^2) * (ξ_y^2 + ξ_y) / 2
i == 8 && return (ξ_x^2 - ξ_x) * (1 - ξ_y^2) / 2
i == 9 && return (1 - ξ_x^2) * (1 - ξ_y^2)
throw(ArgumentError("no shape function $i for interpolation $ip"))
end
#####################################
# Lagrange RefQuadrilateral order 3 #
#####################################
getnbasefunctions(::Lagrange{RefQuadrilateral, 3}) = 16
edgedof_indices(::Lagrange{RefQuadrilateral, 3}) = ((1, 2, 5, 6), (2, 3, 7, 8), (3, 4, 9, 10), (4, 1, 11, 12))
edgedof_interior_indices(::Lagrange{RefQuadrilateral, 3}) = ((5, 6), (7, 8), (9, 10), (11, 12))
facedof_indices(ip::Lagrange{RefQuadrilateral, 3}) = (ntuple(i -> i, getnbasefunctions(ip)),)
facedof_interior_indices(::Lagrange{RefQuadrilateral, 3}) = ((13, 14, 15, 16),)
function reference_coordinates(::Lagrange{RefQuadrilateral, 3})
return [
Vec{2, Float64}((-1.0, -1.0)),
Vec{2, Float64}((1.0, -1.0)),
Vec{2, Float64}((1.0, 1.0)),
Vec{2, Float64}((-1.0, 1.0)),
Vec{2, Float64}((-1 / 3, -1.0)),
Vec{2, Float64}((1 / 3, -1.0)),
Vec{2, Float64}((1.0, -1 / 3)),
Vec{2, Float64}((1.0, 1 / 3)),
Vec{2, Float64}((1 / 3, 1.0)),
Vec{2, Float64}((-1 / 3, 1.0)),
Vec{2, Float64}((-1.0, 1 / 3)),
Vec{2, Float64}((-1.0, -1 / 3)),
Vec{2, Float64}((-1 / 3, -1 / 3)),
Vec{2, Float64}((1 / 3, -1 / 3)),
Vec{2, Float64}((-1 / 3, 1 / 3)),
Vec{2, Float64}((1 / 3, 1 / 3)),
]
end
function reference_shape_value(ip::Lagrange{RefQuadrilateral, 3}, ξ::Vec{2}, i::Int)
# See https://defelement.com/elements/examples/quadrilateral-Q-3.html
# Transform domain from [-1, 1] × [-1, 1] to [0, 1] × [0, 1]
ξ_x = (ξ[1] + 1) / 2
ξ_y = (ξ[2] + 1) / 2
i == 1 && return (81 * ξ_x^3 * ξ_y^3) / 4 - (81 * ξ_x^3 * ξ_y^2) / 2 + (99 * ξ_x^3 * ξ_y) / 4 - (9 * ξ_x^3) / 2 - (81 * ξ_x^2 * ξ_y^3) / 2 + (81 * ξ_x^2 * ξ_y^2) - (99 * ξ_x^2 * ξ_y) / 2 + (9 * ξ_x^2) + (99 * ξ_x * ξ_y^3) / 4 - (99 * ξ_x * ξ_y^2) / 2 + (121 * ξ_x * ξ_y) / 4 - (11 * ξ_x) / 2 - (9 * ξ_y^3) / 2 + 9 * ξ_y^2 - (11 * ξ_y) / 2 + 1
i == 2 && return (ξ_x * (- 81 * ξ_x^2 * ξ_y^3 + 162 * ξ_x^2 * ξ_y^2 - 99 * ξ_x^2 * ξ_y + 18 * ξ_x^2 + 81 * ξ_x * ξ_y^3 - 162 * ξ_x * ξ_y^2 + 99 * ξ_x * ξ_y - 18 * ξ_x - 18 * ξ_y^3 + 36 * ξ_y^2 - 22 * ξ_y + 4)) / 4
i == 4 && return (ξ_y * (- 81 * ξ_x^3 * ξ_y^2 + 81 * ξ_x^3 * ξ_y - 18 * ξ_x^3 + 162 * ξ_x^2 * ξ_y^2 - 162 * ξ_x^2 * ξ_y + 36 * ξ_x^2 - 99 * ξ_x * ξ_y^2 + 99 * ξ_x * ξ_y - 22 * ξ_x + 18 * ξ_y^2 - 18 * ξ_y + 4)) / 4
i == 3 && return (ξ_x * ξ_y * (81 * ξ_x^2 * ξ_y^2 - 81 * ξ_x^2 * ξ_y + 18 * ξ_x^2 - 81 * ξ_x * ξ_y^2 + 81 * ξ_x * ξ_y - 18 * ξ_x + 18 * ξ_y^2 - 18 * ξ_y + 4)) / 4
i == 5 && return (9 * ξ_x * (- 27 * ξ_x^2 * ξ_y^3 + 54 * ξ_x^2 * ξ_y^2 - 33 * ξ_x^2 * ξ_y + 6 * ξ_x^2 + 45 * ξ_x * ξ_y^3 - 90 * ξ_x * ξ_y^2 + 55 * ξ_x * ξ_y - 10 * ξ_x - 18 * ξ_y^3 + 36 * ξ_y^2 - 22 * ξ_y + 4)) / 4
i == 6 && return (9 * ξ_x * (27 * ξ_x^2 * ξ_y^3 - 54 * ξ_x^2 * ξ_y^2 + 33 * ξ_x^2 * ξ_y - 6 * ξ_x^2 - 36 * ξ_x * ξ_y^3 + 72 * ξ_x * ξ_y^2 - 44 * ξ_x * ξ_y + 8 * ξ_x + 9 * ξ_y^3 - 18 * ξ_y^2 + 11 * ξ_y - 2)) / 4
i == 12 && return (9 * ξ_y * (- 27 * ξ_x^3 * ξ_y^2 + 45 * ξ_x^3 * ξ_y - 18 * ξ_x^3 + 54 * ξ_x^2 * ξ_y^2 - 90 * ξ_x^2 * ξ_y + 36 * ξ_x^2 - 33 * ξ_x * ξ_y^2 + 55 * ξ_x * ξ_y - 22 * ξ_x + 6 * ξ_y^2 - 10 * ξ_y + 4)) / 4
i == 11 && return (9 * ξ_y * (27 * ξ_x^3 * ξ_y^2 - 36 * ξ_x^3 * ξ_y + 9 * ξ_x^3 - 54 * ξ_x^2 * ξ_y^2 + 72 * ξ_x^2 * ξ_y - 18 * ξ_x^2 + 33 * ξ_x * ξ_y^2 - 44 * ξ_x * ξ_y + 11 * ξ_x - 6 * ξ_y^2 + 8 * ξ_y - 2)) / 4
i == 7 && return (9 * ξ_x * ξ_y * (27 * ξ_x^2 * ξ_y^2 - 45 * ξ_x^2 * ξ_y + 18 * ξ_x^2 - 27 * ξ_x * ξ_y^2 + 45 * ξ_x * ξ_y - 18 * ξ_x + 6 * ξ_y^2 - 10 * ξ_y + 4)) / 4
i == 8 && return (9 * ξ_x * ξ_y * (- 27 * ξ_x^2 * ξ_y^2 + 36 * ξ_x^2 * ξ_y - 9 * ξ_x^2 + 27 * ξ_x * ξ_y^2 - 36 * ξ_x * ξ_y + 9 * ξ_x - 6 * ξ_y^2 + 8 * ξ_y - 2)) / 4
i == 10 && return (9 * ξ_x * ξ_y * (27 * ξ_x^2 * ξ_y^2 - 27 * ξ_x^2 * ξ_y + 6 * ξ_x^2 - 45 * ξ_x * ξ_y^2 + 45 * ξ_x * ξ_y - 10 * ξ_x + 18 * ξ_y^2 - 18 * ξ_y + 4)) / 4
i == 9 && return (9 * ξ_x * ξ_y * (- 27 * ξ_x^2 * ξ_y^2 + 27 * ξ_x^2 * ξ_y - 6 * ξ_x^2 + 36 * ξ_x * ξ_y^2 - 36 * ξ_x * ξ_y + 8 * ξ_x - 9 * ξ_y^2 + 9 * ξ_y - 2)) / 4
i == 13 && return (81 * ξ_x * ξ_y * (9 * ξ_x^2 * ξ_y^2 - 15 * ξ_x^2 * ξ_y + 6 * ξ_x^2 - 15 * ξ_x * ξ_y^2 + 25 * ξ_x * ξ_y - 10 * ξ_x + 6 * ξ_y^2 - 10 * ξ_y + 4)) / 4
i == 14 && return (81 * ξ_x * ξ_y * (- 9 * ξ_x^2 * ξ_y^2 + 15 * ξ_x^2 * ξ_y - 6 * ξ_x^2 + 12 * ξ_x * ξ_y^2 - 20 * ξ_x * ξ_y + 8 * ξ_x - 3 * ξ_y^2 + 5 * ξ_y - 2)) / 4
i == 15 && return (81 * ξ_x * ξ_y * (- 9 * ξ_x^2 * ξ_y^2 + 12 * ξ_x^2 * ξ_y - 3 * ξ_x^2 + 15 * ξ_x * ξ_y^2 - 20 * ξ_x * ξ_y + 5 * ξ_x - 6 * ξ_y^2 + 8 * ξ_y - 2)) / 4
i == 16 && return (81 * ξ_x * ξ_y * (9 * ξ_x^2 * ξ_y^2 - 12 * ξ_x^2 * ξ_y + 3 * ξ_x^2 - 12 * ξ_x * ξ_y^2 + 16 * ξ_x * ξ_y - 4 * ξ_x + 3 * ξ_y^2 - 4 * ξ_y + 1)) / 4
throw(ArgumentError("no shape function $i for interpolation $ip"))
end
################################
# Lagrange RefTriangle order 1 #
################################
getnbasefunctions(::Lagrange{RefTriangle, 1}) = 3
edgedof_indices(::Lagrange{RefTriangle, 1}) = ((1, 2), (2, 3), (3, 1))
facedof_indices(ip::Lagrange{RefTriangle, 1}) = (ntuple(i -> i, getnbasefunctions(ip)),)
function reference_coordinates(::Lagrange{RefTriangle, 1})
return [
Vec{2, Float64}((1.0, 0.0)),
Vec{2, Float64}((0.0, 1.0)),
Vec{2, Float64}((0.0, 0.0)),
]
end
function reference_shape_value(ip::Lagrange{RefTriangle, 1}, ξ::Vec{2}, i::Int)
ξ_x = ξ[1]
ξ_y = ξ[2]
i == 1 && return ξ_x
i == 2 && return ξ_y
i == 3 && return 1 - ξ_x - ξ_y
throw(ArgumentError("no shape function $i for interpolation $ip"))
end
################################
# Lagrange RefTriangle order 2 #
################################
getnbasefunctions(::Lagrange{RefTriangle, 2}) = 6
edgedof_indices(::Lagrange{RefTriangle, 2}) = ((1, 2, 4), (2, 3, 5), (3, 1, 6))
edgedof_interior_indices(::Lagrange{RefTriangle, 2}) = ((4,), (5,), (6,))
facedof_indices(ip::Lagrange{RefTriangle, 2}) = (ntuple(i -> i, getnbasefunctions(ip)),)
function reference_coordinates(::Lagrange{RefTriangle, 2})
return [
Vec{2, Float64}((1.0, 0.0)),
Vec{2, Float64}((0.0, 1.0)),
Vec{2, Float64}((0.0, 0.0)),
Vec{2, Float64}((0.5, 0.5)),
Vec{2, Float64}((0.0, 0.5)),
Vec{2, Float64}((0.5, 0.0)),
]
end
function reference_shape_value(ip::Lagrange{RefTriangle, 2}, ξ::Vec{2}, i::Int)
ξ_x = ξ[1]
ξ_y = ξ[2]
γ = 1 - ξ_x - ξ_y
i == 1 && return ξ_x * (2ξ_x - 1)
i == 2 && return ξ_y * (2ξ_y - 1)
i == 3 && return γ * (2γ - 1)
i == 4 && return 4ξ_x * ξ_y
i == 5 && return 4ξ_y * γ
i == 6 && return 4ξ_x * γ
throw(ArgumentError("no shape function $i for interpolation $ip"))
end
######################################
# Lagrange RefTriangle order 3, 4, 5 #
######################################
# see https://getfem.readthedocs.io/en/latest/userdoc/appendixA.html
const Lagrange2Tri345 = Union{
Lagrange{RefTriangle, 3},
Lagrange{RefTriangle, 4},
Lagrange{RefTriangle, 5},
}
function getnbasefunctions(ip::Lagrange2Tri345)
order = getorder(ip)
return (order + 1) * (order + 2) ÷ 2
end
# Permutation to switch numbering to Ferrite ordering
const permdof2DLagrange2Tri345 = Dict{Int, Vector{Int}}(
1 => [1, 2, 3],
2 => [3, 6, 1, 5, 4, 2],
3 => [4, 10, 1, 7, 9, 8, 5, 2, 3, 6],
4 => [5, 15, 1, 9, 12, 14, 13, 10, 6, 2, 3, 4, 7, 8, 11],
5 => [6, 21, 1, 11, 15, 18, 20, 19, 16, 12, 7, 2, 3, 4, 5, 8, 9, 10, 13, 14, 17],
)
function edgedof_indices(ip::Lagrange2Tri345)
order = getorder(ip)
order == 1 && return ((1, 2), (2, 3), (3, 1))
order == 2 && return ((1, 2, 4), (2, 3, 5), (3, 1, 6))
order == 3 && return ((1, 2, 4, 5), (2, 3, 6, 7), (3, 1, 8, 9))
order == 4 && return ((1, 2, 4, 5, 6), (2, 3, 7, 8, 9), (3, 1, 10, 11, 12))
order == 5 && return ((1, 2, 4, 5, 6, 7), (2, 3, 8, 9, 10, 11), (3, 1, 12, 13, 14, 15))
throw(ArgumentError("Unsupported order $order for Lagrange on triangles."))
end
function edgedof_interior_indices(ip::Lagrange2Tri345)
order = getorder(ip)
order == 1 && return ((), (), ())
order == 2 && return ((4,), (5,), (6,))
order == 3 && return ((4, 5), (6, 7), (8, 9))
order == 4 && return ((4, 5, 6), (7, 8, 9), (10, 11, 12))
order == 5 && return ((4, 5, 6, 7), (8, 9, 10, 11), (12, 13, 14, 15))
throw(ArgumentError("Unsupported order $order for Lagrange on triangles."))
end
facedof_indices(ip::Lagrange2Tri345) = (ntuple(i -> i, getnbasefunctions(ip)),)
function facedof_interior_indices(ip::Lagrange2Tri345)
order = getorder(ip)
ncellintdofs = (order + 1) * (order + 2) ÷ 2 - 3 * order
totaldofs = getnbasefunctions(ip)
return (ntuple(i -> totaldofs - ncellintdofs + i, ncellintdofs),)
end
function reference_coordinates(ip::Lagrange2Tri345)
order = getorder(ip)
coordpts = Vector{Vec{2, Float64}}()
for k in 0:order
for l in 0:(order - k)
push!(coordpts, Vec{2, Float64}((l / order, k / order)))
end
end
return permute!(coordpts, permdof2DLagrange2Tri345[order])
end
function reference_shape_value(ip::Lagrange2Tri345, ξ::Vec{2}, i::Int)
if !(0 < i <= getnbasefunctions(ip))
throw(ArgumentError("no shape function $i for interpolation $ip"))
end
order = getorder(ip)
i = permdof2DLagrange2Tri345[order][i]
ξ_x = ξ[1]
ξ_y = ξ[2]
i1, i2, i3 = _numlin_basis2D(i, order)
val = one(ξ_y)
i1 ≥ 1 && (val *= prod((order - order * (ξ_x + ξ_y) - j) / (j + 1) for j in 0:(i1 - 1)))
i2 ≥ 1 && (val *= prod((order * ξ_x - j) / (j + 1) for j in 0:(i2 - 1)))
i3 ≥ 1 && (val *= prod((order * ξ_y - j) / (j + 1) for j in 0:(i3 - 1)))
return val
end
function _numlin_basis2D(i, order)
c, j1, j2, j3 = 0, 0, 0, 0
for k in 0:order
if i <= c + (order + 1 - k)
j2 = i - c - 1
break
else
j3 += 1
c += order + 1 - k
end
end
j1 = order - j2 - j3
return j1, j2, j3
end
###################################
# Lagrange RefTetrahedron order 1 #
###################################
getnbasefunctions(::Lagrange{RefTetrahedron, 1}) = 4
facedof_indices(::Lagrange{RefTetrahedron, 1}) = ((1, 3, 2), (1, 2, 4), (2, 3, 4), (1, 4, 3))
edgedof_indices(::Lagrange{RefTetrahedron, 1}) = ((1, 2), (2, 3), (3, 1), (1, 4), (2, 4), (3, 4))
function reference_coordinates(::Lagrange{RefTetrahedron, 1})
return [
Vec{3, Float64}((0.0, 0.0, 0.0)),
Vec{3, Float64}((1.0, 0.0, 0.0)),
Vec{3, Float64}((0.0, 1.0, 0.0)),
Vec{3, Float64}((0.0, 0.0, 1.0)),
]
end
function reference_shape_value(ip::Lagrange{RefTetrahedron, 1}, ξ::Vec{3}, i::Int)
ξ_x = ξ[1]
ξ_y = ξ[2]
ξ_z = ξ[3]
i == 1 && return 1 - ξ_x - ξ_y - ξ_z
i == 2 && return ξ_x
i == 3 && return ξ_y
i == 4 && return ξ_z
throw(ArgumentError("no shape function $i for interpolation $ip"))
end
###################################
# Lagrange RefTetrahedron order 2 #
###################################
getnbasefunctions(::Lagrange{RefTetrahedron, 2}) = 10
facedof_indices(::Lagrange{RefTetrahedron, 2}) = ((1, 3, 2, 7, 6, 5), (1, 2, 4, 5, 9, 8), (2, 3, 4, 6, 10, 9), (1, 4, 3, 8, 10, 7))
edgedof_indices(::Lagrange{RefTetrahedron, 2}) = ((1, 2, 5), (2, 3, 6), (3, 1, 7), (1, 4, 8), (2, 4, 9), (3, 4, 10))
edgedof_interior_indices(::Lagrange{RefTetrahedron, 2}) = ((5,), (6,), (7,), (8,), (9,), (10,))
function reference_coordinates(::Lagrange{RefTetrahedron, 2})
return [
Vec{3, Float64}((0.0, 0.0, 0.0)),
Vec{3, Float64}((1.0, 0.0, 0.0)),
Vec{3, Float64}((0.0, 1.0, 0.0)),
Vec{3, Float64}((0.0, 0.0, 1.0)),
Vec{3, Float64}((0.5, 0.0, 0.0)),
Vec{3, Float64}((0.5, 0.5, 0.0)),
Vec{3, Float64}((0.0, 0.5, 0.0)),
Vec{3, Float64}((0.0, 0.0, 0.5)),
Vec{3, Float64}((0.5, 0.0, 0.5)),
Vec{3, Float64}((0.0, 0.5, 0.5)),
]
end
# http://www.colorado.edu/engineering/CAS/courses.d/AFEM.d/AFEM.Ch09.d/AFEM.Ch09.pdf
# http://www.colorado.edu/engineering/CAS/courses.d/AFEM.d/AFEM.Ch10.d/AFEM.Ch10.pdf
function reference_shape_value(ip::Lagrange{RefTetrahedron, 2}, ξ::Vec{3}, i::Int)
ξ_x = ξ[1]
ξ_y = ξ[2]
ξ_z = ξ[3]
i == 1 && return (-2 * ξ_x - 2 * ξ_y - 2 * ξ_z + 1) * (-ξ_x - ξ_y - ξ_z + 1)
i == 2 && return ξ_x * (2 * ξ_x - 1)
i == 3 && return ξ_y * (2 * ξ_y - 1)
i == 4 && return ξ_z * (2 * ξ_z - 1)
i == 5 && return ξ_x * (-4 * ξ_x - 4 * ξ_y - 4 * ξ_z + 4)
i == 6 && return 4 * ξ_x * ξ_y
i == 7 && return 4 * ξ_y * (-ξ_x - ξ_y - ξ_z + 1)
i == 8 && return ξ_z * (-4 * ξ_x - 4 * ξ_y - 4 * ξ_z + 4)
i == 9 && return 4 * ξ_x * ξ_z
i == 10 && return 4 * ξ_y * ξ_z
throw(ArgumentError("no shape function $i for interpolation $ip"))
end
##################################
# Lagrange RefHexahedron order 1 #
##################################
getnbasefunctions(::Lagrange{RefHexahedron, 1}) = 8
facedof_indices(::Lagrange{RefHexahedron, 1}) = ((1, 4, 3, 2), (1, 2, 6, 5), (2, 3, 7, 6), (3, 4, 8, 7), (1, 5, 8, 4), (5, 6, 7, 8))
edgedof_indices(::Lagrange{RefHexahedron, 1}) = ((1, 2), (2, 3), (3, 4), (4, 1), (5, 6), (6, 7), (7, 8), (8, 5), (1, 5), (2, 6), (3, 7), (4, 8))
function reference_coordinates(::Lagrange{RefHexahedron, 1})
return [
Vec{3, Float64}((-1.0, -1.0, -1.0)),
Vec{3, Float64}((1.0, -1.0, -1.0)),
Vec{3, Float64}((1.0, 1.0, -1.0)),
Vec{3, Float64}((-1.0, 1.0, -1.0)),
Vec{3, Float64}((-1.0, -1.0, 1.0)),
Vec{3, Float64}((1.0, -1.0, 1.0)),
Vec{3, Float64}((1.0, 1.0, 1.0)),
Vec{3, Float64}((-1.0, 1.0, 1.0)),
]
end
function reference_shape_value(ip::Lagrange{RefHexahedron, 1}, ξ::Vec{3}, i::Int)
ξ_x = ξ[1]
ξ_y = ξ[2]
ξ_z = ξ[3]
i == 1 && return (1 - ξ_x) * (1 - ξ_y) * (1 - ξ_z) / 8
i == 2 && return (1 + ξ_x) * (1 - ξ_y) * (1 - ξ_z) / 8
i == 3 && return (1 + ξ_x) * (1 + ξ_y) * (1 - ξ_z) / 8
i == 4 && return (1 - ξ_x) * (1 + ξ_y) * (1 - ξ_z) / 8
i == 5 && return (1 - ξ_x) * (1 - ξ_y) * (1 + ξ_z) / 8
i == 6 && return (1 + ξ_x) * (1 - ξ_y) * (1 + ξ_z) / 8
i == 7 && return (1 + ξ_x) * (1 + ξ_y) * (1 + ξ_z) / 8
i == 8 && return (1 - ξ_x) * (1 + ξ_y) * (1 + ξ_z) / 8
throw(ArgumentError("no shape function $i for interpolation $ip"))
end
##################################
# Lagrange RefHexahedron order 2 #
##################################
# Based on vtkTriQuadraticHexahedron (see https://kitware.github.io/vtk-examples/site/Cxx/GeometricObjects/IsoparametricCellsDemo/)
getnbasefunctions(::Lagrange{RefHexahedron, 2}) = 27
facedof_indices(::Lagrange{RefHexahedron, 2}) = (
(1, 4, 3, 2, 12, 11, 10, 9, 21),
(1, 2, 6, 5, 9, 18, 13, 17, 22),
(2, 3, 7, 6, 10, 19, 14, 18, 23),
(3, 4, 8, 7, 11, 20, 15, 19, 24),
(1, 5, 8, 4, 17, 16, 20, 12, 25),
(5, 6, 7, 8, 13, 14, 15, 16, 26),
)
facedof_interior_indices(::Lagrange{RefHexahedron, 2}) = (
(21,), (22,), (23,), (24,), (25,), (26,),
)
edgedof_indices(::Lagrange{RefHexahedron, 2}) = (
(1, 2, 9),
(2, 3, 10),
(3, 4, 11),
(4, 1, 12),
(5, 6, 13),
(6, 7, 14),
(7, 8, 15),
(8, 5, 16),
(1, 5, 17),
(2, 6, 18),
(3, 7, 19),
(4, 8, 20),
)
edgedof_interior_indices(::Lagrange{RefHexahedron, 2}) = (
(9,), (10,), (11,), (12,), (13,), (14,), (15,), (16,), (17), (18,), (19,), (20,),
)
volumedof_interior_indices(::Lagrange{RefHexahedron, 2}) = (27,)
function reference_coordinates(::Lagrange{RefHexahedron, 2})
return [
# vertex
Vec{3, Float64}((-1.0, -1.0, -1.0)), # 1
Vec{3, Float64}((1.0, -1.0, -1.0)), # 2
Vec{3, Float64}((1.0, 1.0, -1.0)), # 3
Vec{3, Float64}((-1.0, 1.0, -1.0)), # 4