-
Notifications
You must be signed in to change notification settings - Fork 62
/
basis_routines.f90
8320 lines (7711 loc) · 401 KB
/
basis_routines.f90
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
!> \file
!> \author Chris Bradley
!> \brief This module contains all basis function routines.
!>
!> \section LICENSE
!>
!> Version: MPL 1.1/GPL 2.0/LGPL 2.1
!>
!> The contents of this file are subject to the Mozilla Public License
!> Version 1.1 (the "License"); you may not use this file except in
!> compliance with the License. You may obtain a copy of the License at
!> http://www.mozilla.org/MPL/
!>
!> Software distributed under the License is distributed on an "AS IS"
!> basis, WITHOUT WARRANTY OF ANY KIND, either express or implied. See the
!> License for the specific language governing rights and limitations
!> under the License.
!>
!> The Original Code is OpenCMISS
!>
!> The Initial Developer of the Original Code is University of Auckland,
!> Auckland, New Zealand, the University of Oxford, Oxford, United
!> Kingdom and King's College, London, United Kingdom. Portions created
!> by the University of Auckland, the University of Oxford and King's
!> College, London are Copyright (C) 2007-2010 by the University of
!> Auckland, the University of Oxford and King's College, London.
!> All Rights Reserved.
!>
!> Contributor(s):
!>
!> Alternatively, the contents of this file may be used under the terms of
!> either the GNU General Public License Version 2 or later (the "GPL"), or
!> the GNU Lesser General Public License Version 2.1 or later (the "LGPL"),
!> in which case the provisions of the GPL or the LGPL are applicable instead
!> of those above. If you wish to allow use of your version of this file only
!> under the terms of either the GPL or the LGPL, and not to allow others to
!> use your version of this file under the terms of the MPL, indicate your
!> decision by deleting the provisions above and replace them with the notice
!> and other provisions required by the GPL or the LGPL. If you do not delete
!> the provisions above, a recipient may use your version of this file under
!> the terms of any one of the MPL, the GPL or the LGPL.
!>
!> This module contains all basis function routines.
MODULE BasisRoutines
USE BaseRoutines
USE BasisAccessRoutines
USE CONSTANTS
USE INPUT_OUTPUT
USE ISO_VARYING_STRING
USE KINDS
USE STRINGS
USE TYPES
#include "macros.h"
IMPLICIT NONE
PRIVATE
!Module parameters
!> \addtogroup BASIS_ROUTINES_BasisTypes BASIS_ROUTINES::BasisTypes
!> \brief Basis definition type parameters
!> \todo Combine simplex and serendipity elements???
!> \see BASIS_ROUTINES,OPENCMISS_BasisTypes
!>@{
INTEGER(INTG), PARAMETER :: BASIS_LAGRANGE_HERMITE_TP_TYPE=1 !<Lagrange-Hermite tensor product basis type \see BASIS_ROUTINES_BasisTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_SIMPLEX_TYPE=2 !<Simplex basis type \see BASIS_ROUTINES_BasisTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_SERENDIPITY_TYPE=3 !<Serendipity basis type \see BASIS_ROUTINES_BasisTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_AUXILLIARY_TYPE=4 !<Auxillary basis type \see BASIS_ROUTINES_BasisTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_B_SPLINE_TP_TYPE=5 !<B-spline basis type \see BASIS_ROUTINES_BasisTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_FOURIER_LAGRANGE_HERMITE_TP_TYPE=6 !<Fourier-Lagrange tensor product basis type \see BASIS_ROUTINES_BasisTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_EXTENDED_LAGRANGE_TP_TYPE=7 !< Extendend Lagrange tensor product basis type \see BASIS_ROUTINES_BasisTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_RADIAL_TYPE=7 !< Radial basis typee \see BASIS_ROUTINES_BasisTypes,BASIS_ROUTINES
!>@}
!> \addtogroup BASIS_ROUTINES_InterpolationSpecifications BASIS_ROUTINES::InterpolationSpecifications
!> \brief Interpolation specification parameters
!> \see BASIS_ROUTINES,OPENCMISS_InterpolationSpecifications
!>@{
INTEGER(INTG), PARAMETER :: BASIS_LINEAR_LAGRANGE_INTERPOLATION=1 !<Linear Lagrange interpolation specification \see BASIS_ROUTINES_InterpolationSpecifications,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_QUADRATIC_LAGRANGE_INTERPOLATION=2 !<Quadratic Lagrange interpolation specification \see BASIS_ROUTINES_InterpolationSpecifications,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_CUBIC_LAGRANGE_INTERPOLATION=3 !<Cubic Lagrange interpolation specification \see BASIS_ROUTINES_InterpolationSpecifications,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_CUBIC_HERMITE_INTERPOLATION=4 !<Cubic Hermite interpolation specification \see BASIS_ROUTINES_InterpolationSpecifications,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_QUADRATIC1_HERMITE_INTERPOLATION=5 !<Quadratic Hermite (no derivative at xi=0) interpolation specification \see BASIS_ROUTINES_InterpolationSpecifications,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_QUADRATIC2_HERMITE_INTERPOLATION=6 !<Quadratic Hermite (no derivative at xi=1) interpolation specification \see BASIS_ROUTINES_InterpolationSpecifications,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_LINEAR_SIMPLEX_INTERPOLATION=7 !<Linear Simplex interpolation specification \see BASIS_ROUTINES_InterpolationSpecifications,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_QUADRATIC_SIMPLEX_INTERPOLATION=8 !<Quadratic Simplex interpolation specification \see BASIS_ROUTINES_InterpolationSpecifications,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_CUBIC_SIMPLEX_INTERPOLATION=9 !<Cubic Simplex interpolation specification \see BASIS_ROUTINES_InterpolationSpecifications,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_GAUSSIAN_RADIAL_INTERPOLATION=10 !<Gaussian Radial interpolation specification \see BASIS_ROUTINES_InterpolationSpecifications,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_MULTIQUARTIC_RADIAL_INTERPOLATION=11 !<Multiquartic Radial interpolation specification \see BASIS_ROUTINES_InterpolationSpecifications,BASIS_ROUTINES
!>@}
!> \addtogroup BASIS_ROUTINES_InterpolationTypes BASIS_ROUTINES::InterpolationTypes
!> \brief Interpolation type parameters for a Xi direction
!> \see BASIS_ROUTINES
!>@{
INTEGER(INTG), PARAMETER :: BASIS_LAGRANGE_INTERPOLATION=1 !<Lagrange interpolation \see BASIS_ROUTINES_InterpolationTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_HERMITE_INTERPOLATION=2 !<Hermite interpolation \see BASIS_ROUTINES_InterpolationTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_SIMPLEX_INTERPOLATION=3 !<Simplex interpolation \see BASIS_ROUTINES_InterpolationTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_SERENDIPITY_INTERPOLATION=4 !<Serendipity interpolation \see BASIS_ROUTINES_InterpolationTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_TRANSITION_INTERPOLATION=5 !<Transition interpolation \see BASIS_ROUTINES_InterpolationTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_SINGULAR_INTERPOLATION=6 !<Singular interpolation \see BASIS_ROUTINES_InterpolationTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_FOURIER_INTERPOLATION=7 !<Fourier interpolation \see BASIS_ROUTINES_InterpolationTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_RADIAL_INTERPOLATION=8 !<Radial interpolation \see BASIS_ROUTINES_InterpolationTypes,BASIS_ROUTINES
!>@}
!> \addtogroup BASIS_ROUTINES_InterpolationOrder BASIS_ROUTINES::InterpolationOrder
!> \brief Interpolation order for a Xi direction
!> \see BASIS_ROUTINES
!>@{
INTEGER(INTG), PARAMETER :: BASIS_LINEAR_INTERPOLATION_ORDER=1 !<Linear interpolation order \see BASIS_ROUTINES_InterpolationOrder,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_QUADRATIC_INTERPOLATION_ORDER=2 !<Quadratic interpolation order \see BASIS_ROUTINES_InterpolationOrder,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_CUBIC_INTERPOLATION_ORDER=3 !<Cubic interpolation order \see BASIS_ROUTINES_InterpolationOrder,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_QUADRATIC1_INTERPOLATION_ORDER=4 !<Quadratic (no derivative at xi=0) interpolation order \see BASIS_ROUTINES_InterpolationOrder,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_QUADRATIC2_INTERPOLATION_ORDER=5 !<Quadratic (no derivative at xi=1) interpolation order \see BASIS_ROUTINES_InterpolationOrder,BASIS_ROUTINES
!>@}
!> \addtogroup BASIS_ROUTINES_QuadratureTypes BASIS_ROUTINES::QuadratureTypes
!> \brief Quadrature type parameters
!> \see BASIS_ROUTINES,OPENCMISS_QuadratureTypes
!>@{
INTEGER(INTG), PARAMETER :: BASIS_GAUSS_LEGENDRE_QUADRATURE=1 !<Gauss-Legendre quadrature \see BASIS_ROUTINES_QuadratureTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_GAUSS_LAGUERRE_QUADRATURE=2 !<Gauss-Laguerre quadrature \see BASIS_ROUTINES_QuadratureTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_GUASS_HERMITE_QUADRATURE=3 !<Gauss-Hermite quadrature \see BASIS_ROUTINES_QuadratureTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_ADAPTIVE_GAUSS_LEGENDRE_QUADRATURE=4 !<Adaptive Gauss-Legendre quadrature \see BASIS_ROUTINES_QuadratureTypes,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_GAUSS_SIMPLEX_QUADRATURE=5 !<Gauss-Legendre for Simplex elements quadrature \see BASIS_ROUTINES_QuadratureTypes,BASIS_ROUTINES
!>@}
!> \addtogroup BASIS_ROUTINES_XiCollapse BASIS_ROUTINES::XiCollapse
!> \brief Xi collapse parameters
!> \see BASIS_ROUTINES
!>@{
INTEGER(INTG), PARAMETER :: BASIS_XI_COLLAPSED=1 !<The Xi direction is collapsed \see BASIS_ROUTINES_XiCollapse,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_COLLAPSED_AT_XI0=2 !<The Xi direction at the xi=0 end of this Xi direction is collapsed \see BASIS_ROUTINES_XiCollapse,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_COLLAPSED_AT_XI1=3 !<The Xi direction at the xi=1 end of this Xi direction is collapsed \see BASIS_ROUTINES_XiCollapse,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_NOT_COLLAPSED=4 !<The Xi direction is not collapsed \see BASIS_ROUTINES_XiCollapse,BASIS_ROUTINES
!>@}
!> \addtogroup BASIS_ROUTINES_BoundaryXiType BASIS_ROUTINES::BoundaryXiType
!> \brief The boundary xi types
!> \see BASIS_ROUTINES
!>@{
INTEGER(INTG), PARAMETER :: BASIS_NO_BOUNDARY_XI=0 !<The xi is not on a boundary \see BASIS_ROUTINES_XiBoundaryType,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_LINE_BOUNDARY_XI=1 !<The xi is on a boundary line \see BASIS_ROUTINES_XiBoundaryType,BASIS_ROUTINES
INTEGER(INTG), PARAMETER :: BASIS_FACE_BOUNDARY_XI=2 !<The xi is on a boundary face \see BASIS_ROUTINES_XiBoundaryType,BASIS_ROUTINES
!>@}
!!Module types
!
!!>Contains information on the defined basis functions
!TYPE BASIS_FUNCTIONS_TYPE
! PRIVATE
! INTEGER(INTG) :: numberOfBasisFunctions !<The number of basis functions definegd
! TYPE(BASIS_PTR_TYPE), POINTER :: BASES(:) !<The array of pointers to the defined basis functions
!END TYPE BASIS_FUNCTIONS_TYPE
!Module variables
!Interfaces
INTERFACE Basis_CollapsedXiGet
MODULE PROCEDURE BASIS_COLLAPSED_XI_GET
END INTERFACE Basis_CollapsedXiGet
!>Sets/changes the collapsed Xi flags for a basis.
INTERFACE BASIS_COLLAPSED_XI_SET
MODULE PROCEDURE BASIS_COLLAPSED_XI_SET_NUMBER
MODULE PROCEDURE BASIS_COLLAPSED_XI_SET_PTR
END INTERFACE BASIS_COLLAPSED_XI_SET
!>Sets/changes the collapsed Xi flags for a basis.
INTERFACE Basis_CollapsedXiSet
MODULE PROCEDURE BASIS_COLLAPSED_XI_SET_NUMBER
MODULE PROCEDURE BASIS_COLLAPSED_XI_SET_PTR
END INTERFACE Basis_CollapsedXiSet
INTERFACE BASIS_CREATE_START
MODULE PROCEDURE Basis_CreateStart
END INTERFACE BASIS_CREATE_START
INTERFACE BASIS_CREATE_FINISH
MODULE PROCEDURE Basis_CreateFinish
END INTERFACE BASIS_CREATE_FINISH
!>Evaluates the appropriate partial derivative index for the specificied basis function at a Xi location \see BASIS_ROUTINES
INTERFACE BASIS_EVALUATE_XI
MODULE PROCEDURE BASIS_EVALUATE_XI_DP
END INTERFACE BASIS_EVALUATE_XI
!>Evaluates the appropriate partial derivative index for the specificied basis function at a Xi location \see BASIS_ROUTINES
INTERFACE Basis_EvaluateXi
MODULE PROCEDURE BASIS_EVALUATE_XI_DP
END INTERFACE Basis_EvaluateXi
!>Interpolates the appropriate partial derivative index of the elements parameters for basis function at a Gauss point \see BASIS_ROUTINES
INTERFACE BASIS_INTERPOLATE_GAUSS
MODULE PROCEDURE BASIS_INTERPOLATE_GAUSS_DP
END INTERFACE BASIS_INTERPOLATE_GAUSS
!>Interpolates the appropriate partial derivative index of the elements parameters for basis function at a Gauss point \see BASIS_ROUTINES
INTERFACE Basis_InterpolateGauss
MODULE PROCEDURE BASIS_INTERPOLATE_GAUSS_DP
END INTERFACE Basis_InterpolateGauss
!>Interpolates the appropriate partial derivative index of the elements parameters for basis function at Xi location \see BASIS_ROUTINES
INTERFACE BASIS_INTERPOLATE_XI
MODULE PROCEDURE BASIS_INTERPOLATE_XI_DP
END INTERFACE BASIS_INTERPOLATE_XI
!>Interpolates the appropriate partial derivative index of the elements parameters for basis function at Xi location \see BASIS_ROUTINES
INTERFACE Basis_InterpolateXi
MODULE PROCEDURE BASIS_INTERPOLATE_XI_DP
END INTERFACE Basis_InterpolateXi
!>Interpolates the requested partial derivative index(ices) of the element parameters for basis function at a face Gauss point \see BASIS_ROUTINES
INTERFACE BASIS_INTERPOLATE_LOCAL_FACE_GAUSS
MODULE PROCEDURE BASIS_INTERPOLATE_LOCAL_FACE_GAUSS_DP
END INTERFACE BASIS_INTERPOLATE_LOCAL_FACE_GAUSS
!>Interpolates the requested partial derivative index(ices) of the element parameters for basis function at a face Gauss point \see BASIS_ROUTINES
INTERFACE Basis_InterpolateLocalFaceGauss
MODULE PROCEDURE BASIS_INTERPOLATE_LOCAL_FACE_GAUSS_DP
END INTERFACE Basis_InterpolateLocalFaceGauss
INTERFACE Basis_LocalNodeXiCalculate
MODULE PROCEDURE BASIS_LOCAL_NODE_XI_CALCULATE
END INTERFACE Basis_LocalNodeXiCalculate
INTERFACE Basis_InterpolationXiGet
MODULE PROCEDURE BASIS_INTERPOLATION_XI_GET
END INTERFACE Basis_InterpolationXiGet
!>Sets/changes the interpolation type in each Xi direction for a basis
INTERFACE BASIS_INTERPOLATION_XI_SET
MODULE PROCEDURE BASIS_INTERPOLATION_XI_SET_NUMBER
MODULE PROCEDURE BASIS_INTERPOLATION_XI_SET_PTR
END INTERFACE BASIS_INTERPOLATION_XI_SET
!>Sets/changes the interpolation type in each Xi direction for a basis
INTERFACE Basis_InterpolationXiSet
MODULE PROCEDURE BASIS_INTERPOLATION_XI_SET_NUMBER
MODULE PROCEDURE BASIS_INTERPOLATION_XI_SET_PTR
END INTERFACE Basis_InterpolationXiSet
INTERFACE Basis_NumberOfLocalNodesGet
MODULE PROCEDURE BASIS_NUMBER_OF_LOCAL_NODES_GET
END INTERFACE Basis_NumberOfLocalNodesGet
INTERFACE BASIS_NUMBER_OF_XI_GET
MODULE PROCEDURE Basis_NumberOfXiGet
END INTERFACE BASIS_NUMBER_OF_XI_GET
INTERFACE BASIS_NUMBER_OF_XI_SET
MODULE PROCEDURE Basis_NumberOfXiSet
END INTERFACE BASIS_NUMBER_OF_XI_SET
INTERFACE Basis_QuadratureDestroy
MODULE PROCEDURE BASIS_QUADRATURE_DESTROY
END INTERFACE Basis_QuadratureDestroy
INTERFACE Basis_QuadratureMultipleGaussXiGet
MODULE PROCEDURE BASIS_QUADRATURE_MULTIPLE_GAUSS_XI_GET
END INTERFACE Basis_QuadratureMultipleGaussXiGet
INTERFACE Basis_QuadratureNumberOfGaussXiGet
MODULE PROCEDURE BASIS_QUADRATURE_NUMBER_OF_GAUSS_XI_GET
END INTERFACE Basis_QuadratureNUmberOfGaussXiGet
INTERFACE Basis_QuadratureNumberOfGaussXiSet
MODULE PROCEDURE BASIS_QUADRATURE_NUMBER_OF_GAUSS_XI_SET
END INTERFACE Basis_QuadratureNUmberOfGaussXiSet
INTERFACE Basis_QuadratureOrderGet
MODULE PROCEDURE BASIS_QUADRATURE_ORDER_GET
END INTERFACE Basis_QuadratureOrderGet
!>Sets/changes the order of a quadrature for a basis quadrature.
INTERFACE BASIS_QUADRATURE_ORDER_SET
MODULE PROCEDURE BASIS_QUADRATURE_ORDER_SET_NUMBER
MODULE PROCEDURE BASIS_QUADRATURE_ORDER_SET_PTR
END INTERFACE BASIS_QUADRATURE_ORDER_SET
!>Sets/changes the order of a quadrature for a basis quadrature.
INTERFACE Basis_QuadratureOrderSet
MODULE PROCEDURE BASIS_QUADRATURE_ORDER_SET_NUMBER
MODULE PROCEDURE BASIS_QUADRATURE_ORDER_SET_PTR
END INTERFACE Basis_QuadratureOrderSet
INTERFACE Basis_QuadratureSingleGaussXiGet
MODULE PROCEDURE BASIS_QUADRATURE_SINGLE_GAUSS_XI_GET
END INTERFACE Basis_QuadratureSingleGaussXiGet
INTERFACE Basis_QuadratureTypeGet
MODULE PROCEDURE BASIS_QUADRATURE_TYPE_GET
END INTERFACE Basis_QuadratureTypeGet
!>Sets/changes the quadrature type for a basis
INTERFACE BASIS_QUADRATURE_TYPE_SET
MODULE PROCEDURE BASIS_QUADRATURE_TYPE_SET_NUMBER
MODULE PROCEDURE BASIS_QUADRATURE_TYPE_SET_PTR
END INTERFACE BASIS_QUADRATURE_TYPE_SET
!>Sets/changes the quadrature type for a basis
INTERFACE Basis_QuadratureTypeSet
MODULE PROCEDURE BASIS_QUADRATURE_TYPE_SET_NUMBER
MODULE PROCEDURE BASIS_QUADRATURE_TYPE_SET_PTR
END INTERFACE Basis_QuadratureTypeSet
INTERFACE Basis_TypeGet
MODULE PROCEDURE BASIS_TYPE_GET
END INTERFACE Basis_TypeGet
!>Sets/changes the type for a basis.
INTERFACE BASIS_TYPE_SET
MODULE PROCEDURE BASIS_TYPE_SET_NUMBER
MODULE PROCEDURE BASIS_TYPE_SET_PTR
END INTERFACE BASIS_TYPE_SET
!>Sets/changes the type for a basis.
INTERFACE Basis_TypeSet
MODULE PROCEDURE BASIS_TYPE_SET_NUMBER
MODULE PROCEDURE BASIS_TYPE_SET_PTR
END INTERFACE Basis_TypeSet
!>Evaluates a linear Simplex basis function
INTERFACE SIMPLEX_LINEAR_EVALUATE
MODULE PROCEDURE SIMPLEX_LINEAR_EVALUATE_DP
END INTERFACE SIMPLEX_LINEAR_EVALUATE
!>Evaluates a quadratic Simplex basis function
INTERFACE SIMPLEX_QUADRATIC_EVALUATE
MODULE PROCEDURE SIMPLEX_QUADRATIC_EVALUATE_DP
END INTERFACE SIMPLEX_QUADRATIC_EVALUATE
!>Evaluates a cubic Simplex basis function
INTERFACE SIMPLEX_CUBIC_EVALUATE
MODULE PROCEDURE SIMPLEX_CUBIC_EVALUATE_DP
END INTERFACE SIMPLEX_CUBIC_EVALUATE
!>Evaluates the Lagrange/Hermite/Fourier tensor product basis function for the given basis
INTERFACE BASIS_LHTP_BASIS_EVALUATE
MODULE PROCEDURE BASIS_LHTP_BASIS_EVALUATE_DP
END INTERFACE BASIS_LHTP_BASIS_EVALUATE
!>Evaluates the Lagrange/Hermite/Fourier tensor product basis function for the given basis
INTERFACE Basis_LHTPBasisEvaluate
MODULE PROCEDURE BASIS_LHTP_BASIS_EVALUATE_DP
END INTERFACE Basis_LHTPBasisEvaluate
PUBLIC BASIS_LAGRANGE_HERMITE_TP_TYPE,BASIS_SIMPLEX_TYPE,BASIS_SERENDIPITY_TYPE,BASIS_AUXILLIARY_TYPE,BASIS_B_SPLINE_TP_TYPE, &
& BASIS_FOURIER_LAGRANGE_HERMITE_TP_TYPE,BASIS_EXTENDED_LAGRANGE_TP_TYPE, BASIS_RADIAL_TYPE
PUBLIC BASIS_LINEAR_LAGRANGE_INTERPOLATION,BASIS_QUADRATIC_LAGRANGE_INTERPOLATION,BASIS_CUBIC_LAGRANGE_INTERPOLATION, &
& BASIS_CUBIC_HERMITE_INTERPOLATION,BASIS_QUADRATIC1_HERMITE_INTERPOLATION,BASIS_QUADRATIC2_HERMITE_INTERPOLATION, &
& BASIS_LINEAR_SIMPLEX_INTERPOLATION,BASIS_QUADRATIC_SIMPLEX_INTERPOLATION,BASIS_CUBIC_SIMPLEX_INTERPOLATION, &
& BASIS_GAUSSIAN_RADIAL_INTERPOLATION, BASIS_MULTIQUARTIC_RADIAL_INTERPOLATION
PUBLIC BASIS_LINEAR_INTERPOLATION_ORDER,BASIS_QUADRATIC_INTERPOLATION_ORDER,BASIS_CUBIC_INTERPOLATION_ORDER, &
& BASIS_QUADRATIC1_INTERPOLATION_ORDER,BASIS_QUADRATIC2_INTERPOLATION_ORDER
PUBLIC BASIS_GAUSS_LEGENDRE_QUADRATURE,BASIS_GAUSS_LAGUERRE_QUADRATURE,BASIS_GUASS_HERMITE_QUADRATURE,&
& BASIS_ADAPTIVE_GAUSS_LEGENDRE_QUADRATURE,BASIS_GAUSS_SIMPLEX_QUADRATURE
PUBLIC BASIS_XI_COLLAPSED,BASIS_COLLAPSED_AT_XI0,BASIS_COLLAPSED_AT_XI1,BASIS_NOT_COLLAPSED
PUBLIC BASIS_NO_BOUNDARY_XI,BASIS_LINE_BOUNDARY_XI,BASIS_FACE_BOUNDARY_XI
PUBLIC Basis_AreaToXiCoordinates
PUBLIC Basis_BoundaryXiToXi,Basis_XiToBoundaryXi
PUBLIC BASIS_COLLAPSED_XI_SET
PUBLIC Basis_CollapsedXiSet
PUBLIC BASIS_EVALUATE_XI
PUBLIC Basis_EvaluateXi
PUBLIC Basis_GaussPointsCalculate
PUBLIC BASIS_INTERPOLATE_GAUSS,BASIS_INTERPOLATE_XI,BASIS_INTERPOLATE_LOCAL_FACE_GAUSS
PUBLIC Basis_InterpolateGauss,Basis_InterpolateXi,Basis_InterpolateLocalFaceGauss
PUBLIC BASIS_LOCAL_NODE_XI_CALCULATE
PUBLIC Basis_LocalNodeXiCalculate
PUBLIC BASIS_NUMBER_OF_LOCAL_NODES_GET
PUBLIC Basis_NumberOfLocalNodesGet
PUBLIC BASIS_INTERPOLATION_XI_SET
PUBLIC Basis_InterpolationXiSet
PUBLIC BASIS_NUMBER_OF_XI_SET
PUBLIC Basis_NumberOfXiSet
PUBLIC BASIS_QUADRATURE_NUMBER_OF_GAUSS_XI_SET
PUBLIC Basis_QuadratureNumberOfGaussXiSet
PUBLIC BASIS_QUADRATURE_DESTROY,BASIS_QUADRATURE_ORDER_SET,BASIS_QUADRATURE_TYPE_SET
PUBLIC Basis_QuadratureDestroy,Basis_QuadratureOrderSet,Basis_QuadratureTypeSet
PUBLIC BASIS_TYPE_SET
PUBLIC Basis_TypeSet
PUBLIC Basis_QuadratureLocalFaceGaussEvaluateSet
PUBLIC BASIS_CREATE_START,BASIS_CREATE_FINISH
PUBLIC Basis_CreateStart,Basis_CreateFinish
PUBLIC Basis_Destroy
PUBLIC Bases_Finalise,Bases_Initialise
PUBLIC BASIS_COLLAPSED_XI_GET
PUBLIC Basis_CollapsedXiGet
PUBLIC BASIS_INTERPOLATION_XI_GET
PUBLIC Basis_InterpolationXiGet
PUBLIC BASIS_NUMBER_OF_XI_GET
PUBLIC Basis_NumberOfXiGet
PUBLIC BASIS_QUADRATURE_NUMBER_OF_GAUSS_XI_GET
PUBLIC Basis_QuadratureNumberOfGaussXiGet
!!TODO: These should be changed to operate on an array rather than have separate single/multiple routines
PUBLIC BASIS_QUADRATURE_SINGLE_GAUSS_XI_GET, BASIS_QUADRATURE_MULTIPLE_GAUSS_XI_GET
PUBLIC Basis_QuadratureSingleGaussXiGet,Basis_QuadratureMultipleGaussXiGet
PUBLIC BASIS_QUADRATURE_ORDER_GET
PUBLIC Basis_QuadratureOrderGet
PUBLIC BASIS_QUADRATURE_TYPE_GET
PUBLIC Basis_QuadratureTypeGet
PUBLIC BASIS_TYPE_GET
PUBLIC Basis_TypeGet
PUBLIC Basis_XiToAreaCoordinates
CONTAINS
!
!================================================================================================================================
!
!>Finalises the bases and deallocates all memory
SUBROUTINE Bases_Finalise(err,error,*)
!Argument variables
INTEGER(INTG), INTENT(OUT) :: err !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: error !<The error string
!Local Variables
ENTERS("Bases_Finalise",err,error,*999)
!Destroy any created basis functions
DO WHILE(basisFunctions%numberOfBasisFunctions>0)
CALL BASIS_DESTROY(basisFunctions%bases(1)%ptr,err,error,*999)
ENDDO !nb
!Destroy basis functions and deallocated any memory allocated
basisFunctions%numberOfBasisFunctions=0
IF(ALLOCATED(basisFunctions%bases)) DEALLOCATE(basisFunctions%bases)
EXITS("Bases_Finalise")
RETURN
999 ERRORSEXITS("Bases_Finalise",err,error)
RETURN 1
END SUBROUTINE Bases_Finalise
!
!================================================================================================================================
!
!>Initialises the bases.
SUBROUTINE Bases_Initialise(err,error,*)
!Argument variables
INTEGER(INTG), INTENT(OUT) :: err !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: error !<The error string
!Local Variables
ENTERS("Bases_Initialise",err,error,*999)
basisFunctions%numberOfBasisFunctions=0
EXITS("Bases_Initialise")
RETURN
999 ERRORSEXITS("Bases_Initialise",err,error)
RETURN 1
END SUBROUTINE Bases_Initialise
!
!================================================================================================================================
!
!>Converts area coordinates to xi coordinates. \see BASIS_ROUTINES::Basis_XiToAreaCoordinates
SUBROUTINE Basis_AreaToXiCoordinates(areaCoordinates,xiCoordinates,err,error,*)
!Argument variables
REAL(DP), INTENT(IN) :: areaCoordinates(:) !<The area coordinates to convert
REAL(DP), INTENT(OUT) :: xiCoordinates(:) !<On return, the xi coordinates corresponding to the area coordinates
INTEGER(INTG), INTENT(OUT) :: err !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: error !<The error string
!Local Variables
TYPE(VARYING_STRING) :: localError
ENTERS("Basis_AreaToXiCoordinates",err,error,*999)
IF(SIZE(areaCoordinates,1)/=(SIZE(xiCoordinates,1)+1)) THEN
localError="Invalid number of coordinates. The number of area coordinates of "// &
& TRIM(NumberToVString(SIZE(areaCoordinates,1),"*",err,error))// &
& " should be equal to the number of xi coordinates of "// &
& TRIM(NumberToVString(SIZE(xiCoordinates,1),"*",err,error))//" plus one."
CALL FlagError(localError,err,error,*999)
ENDIF
SELECT CASE(SIZE(xiCoordinates,1))
CASE(1)
xiCoordinates(1)=1.0_DP-areaCoordinates(1)
CASE(2)
xiCoordinates(1)=1.0_DP-areaCoordinates(1)
xiCoordinates(2)=1.0_DP-areaCoordinates(2)
CASE(3)
xiCoordinates(1)=1.0_DP-areaCoordinates(1)
xiCoordinates(2)=1.0_DP-areaCoordinates(2)
xiCoordinates(3)=1.0_DP-areaCoordinates(3)
CASE DEFAULT
localError="The number of xi coordinates of "//TRIM(NumberToVString(SIZE(xiCoordinates,1),"*",err,error))// &
& " is invalid. The number must be >= 1 and <= 3."
CALL FlagError(localError,err,error,*999)
END SELECT
EXITS("Basis_AreaToXiCoordinates")
RETURN
999 ERRORSEXITS("Basis_AreaToXiCoordinates",err,error)
RETURN 1
END SUBROUTINE Basis_AreaToXiCoordinates
!
!================================================================================================================================
!
!>Finishes the creation of a new basis \see BASIS_ROUTINES::Basis_CreateStart,OpenCMISS::BasisCreateFinish
SUBROUTINE Basis_CreateFinish(basis,err,error,*)
!Argument variables
TYPE(BASIS_TYPE), POINTER :: basis !<A pointer to the basis to finish the creation of
INTEGER(INTG), INTENT(OUT) :: err !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: error !<The error string
!Local Variables
INTEGER(INTG) :: xiIdx,xiCoordIdx,localNodeIdx,localNodeIdx1,localNodeIdx2,localNodeIdx3,localNodeIdx4,elemParamIdx, &
& localLineIdx,localFaceIdx,columnStart,columnStart2,columnStop,columnStop2
TYPE(VARYING_STRING) :: localError
ENTERS("Basis_CreateFinish",err,error,*999)
IF(.NOT.ASSOCIATED(basis)) CALL FlagError("Basis is not associated.",err,error,*999)
IF(basis%BASIS_FINISHED) CALL FlagError("Basis has already been finished.",err,error,*999)
SELECT CASE(basis%type)
CASE(BASIS_LAGRANGE_HERMITE_TP_TYPE)
CALL Basis_LHTPFamilyCreate(basis,err,error,*999)
CASE(BASIS_SIMPLEX_TYPE)
CALL Basis_SimplexFamilyCreate(basis,err,error,*999)
CASE(BASIS_RADIAL_TYPE)
CALL BASIS_RADIAL_FAMILY_CREATE(basis,err,error,*999)
CASE DEFAULT
localError="Basis type "//TRIM(NumberToVString(BASIS%TYPE,"*",err,error))//" is invalid or not implemented"
CALL FlagError(localError,err,error,*999)
END SELECT
basis%BASIS_FINISHED=.TRUE.
IF(diagnostics1) THEN
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE,"Basis user number = ",basis%USER_NUMBER,err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Basis family number = ",basis%FAMILY_NUMBER,err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Basis global number = ",basis%GLOBAL_NUMBER,err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Basis type = ",basis%type,err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Basis degenerate = ",basis%degenerate,err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Number of Xi directions = ",basis%NUMBER_OF_XI,err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Number of Xi coordinates = ",basis%NUMBER_OF_XI_COORDINATES,err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_XI_COORDINATES,4,4,basis%INTERPOLATION_TYPE, &
& '(" Interpolation type(:) :",4(X,I2))','(25X,4(X,I2))',err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_XI_COORDINATES,4,4,basis%INTERPOLATION_ORDER, &
& '(" Interpolation order(:) :",4(X,I2))','(26X,4(X,I2))',err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_XI,3,3,basis%COLLAPSED_XI, &
& '(" Collapsed xi(:) :",3(X,I2))','(26X,3(X,I2))',err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Number of partial derivatives = ",basis%NUMBER_OF_PARTIAL_DERIVATIVES, &
& err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Total number of nodes = ",basis%NUMBER_OF_NODES,err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_XI_COORDINATES,4,4,basis%NUMBER_OF_NODES_XIC, &
& '(" Number of nodes(:) :",4(X,I2))','(22X,4(X,I2))',err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_NODES,8,8,basis%NUMBER_OF_DERIVATIVES, &
& '(" Number of derivatives(:) :",8(X,I2))','(28X,8(X,I2))',err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Number of element parameters = ",basis%NUMBER_OF_ELEMENT_PARAMETERS, &
& err,error,*999)
! CPB 23/07/07 Doxygen may or may not like this line for some reason????
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_NODES,8,8,basis%NODE_AT_COLLAPSE, &
& '(" Node at collapse(:) :",8(X,L1))','(23X,8(X,L1))',err,error,*999)
CALL WriteString(DIAGNOSTIC_OUTPUT_TYPE," Node position index:",err,error,*999)
DO xiCoordIdx=1,basis%NUMBER_OF_XI_COORDINATES
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Xi coordinate = ",xiCoordIdx,err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_NODES,16,16,basis%NODE_POSITION_INDEX(:,xiCoordIdx), &
& '(" Index(:) :",16(X,I2))','(18X,16(X,I2))',err,error,*999)
ENDDO !xiCoordIdx
CALL WriteString(DIAGNOSTIC_OUTPUT_TYPE," Inverse node position index:",err,error,*999)
SELECT CASE(basis%NUMBER_OF_XI_COORDINATES)
CASE(1)
DO localNodeIdx1=1,basis%NUMBER_OF_NODES_XIC(1)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Xi coordinate = 1, Node position index = ",localNodeIdx1, &
& err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Index = ",basis%NODE_POSITION_INDEX_INV(localNodeIdx1,1,1,1), &
& err,error,*999)
ENDDO !localNodeIdx1
CASE(2)
DO localNodeIdx2=1,basis%NUMBER_OF_NODES_XIC(2)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Xi coordinate = 2, Node position index = ",localNodeIdx2, &
& err,error,*999)
DO localNodeIdx1=1,basis%NUMBER_OF_NODES_XIC(1)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Xi coordinate = 1, Node position index = ",localNodeIdx1, &
& err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Index = ",basis%NODE_POSITION_INDEX_INV(localNodeIdx1, &
& localNodeIdx2,1,1),err,error,*999)
ENDDO !localNodeIdx1
ENDDO !localNodeIdx2
CASE(3)
DO localNodeIdx3=1,basis%NUMBER_OF_NODES_XIC(3)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Xi coordinate = 3, Node position index = ",localNodeIdx3, &
& err,error,*999)
DO localNodeIdx2=1,basis%NUMBER_OF_NODES_XIC(2)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Xi coordinate = 2, Node position index = ",localNodeIdx2, &
& err,error,*999)
DO localNodeIdx1=1,basis%NUMBER_OF_NODES_XIC(1)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Xi coordinate = 1, Node position index = ",localNodeIdx1, &
& err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Index = ",basis%NODE_POSITION_INDEX_INV(localNodeIdx1, &
& localNodeIdx2,localNodeIdx3,1),err,error,*999)
ENDDO !localNodeIdx1
ENDDO !localNodeIdx2
ENDDO !localNodeIdx3
CASE(4)
DO localNodeIdx4=1,basis%NUMBER_OF_NODES_XIC(4)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Xi coordinate = 4, Node position index = ",localNodeIdx4, &
& err,error,*999)
DO localNodeIdx3=1,basis%NUMBER_OF_NODES_XIC(3)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Xi coordinate = 3, Node position index = ",localNodeIdx3, &
& err,error,*999)
DO localNodeIdx2=1,basis%NUMBER_OF_NODES_XIC(2)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Xi coordinate = 2, Node position index = ",localNodeIdx2, &
& err,error,*999)
DO localNodeIdx1=1,basis%NUMBER_OF_NODES_XIC(1)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Xi coordinate = 1, Node position index = ",localNodeIdx1, &
& err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Index = " &
& ,basis%NODE_POSITION_INDEX_INV(localNodeIdx1,localNodeIdx2,localNodeIdx3,localNodeIdx4),err,error,*999)
ENDDO !localNodeIdx1
ENDDO !localNodeIdx2
ENDDO !localNodeIdx3
ENDDO !localNodeIdx4
CASE DEFAULT
CALL FlagError("Invalid number of xi coordinates",err,error,*999)
END SELECT
CALL WriteString(DIAGNOSTIC_OUTPUT_TYPE," Derivative order index:",err,error,*999)
DO xiIdx=1,basis%NUMBER_OF_XI
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Xi = ",xiIdx,err,error,*999)
DO localNodeIdx=1,basis%NUMBER_OF_NODES
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Node = ",localNodeIdx,err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_DERIVATIVES(localNodeIdx),8,8, &
& basis%DERIVATIVE_ORDER_INDEX(:,localNodeIdx,xiIdx),'(" Index(:) :",8(X,I2))','(18X,8(X,I2))',err,error,*999)
ENDDO !localNodeIdx
ENDDO !xiIdx
CALL WriteString(DIAGNOSTIC_OUTPUT_TYPE," Inverse derivative order index:",err,error,*999)
CALL WriteString(DIAGNOSTIC_OUTPUT_TYPE," Element parameter index:",err,error,*999)
DO localNodeIdx=1,basis%NUMBER_OF_NODES
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Node = ",localNodeIdx,err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_DERIVATIVES(localNodeIdx),8,8, &
& basis%ELEMENT_PARAMETER_INDEX(:,localNodeIdx),'(" Index(:) :",8(X,I2))','(18X,8(X,I2))',err,error,*999)
ENDDO !localNodeIdx
CALL WriteString(DIAGNOSTIC_OUTPUT_TYPE," Inverse element parameter index:",err,error,*999)
DO elemParamIdx=1,basis%NUMBER_OF_ELEMENT_PARAMETERS
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Element parameter = ",elemParamIdx,err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,2,2,2, &
& basis%ELEMENT_PARAMETER_INDEX_INV(:,elemParamIdx),'(" Index(:) :",2(X,I2))','(18X,2(X,I2))',err,error,*999)
ENDDO !elemParamIdx
CALL WriteString(DIAGNOSTIC_OUTPUT_TYPE," Partial derivative index:",err,error,*999)
DO localNodeIdx=1,basis%NUMBER_OF_NODES
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Node = ",localNodeIdx,err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_DERIVATIVES(localNodeIdx),8,8, &
& basis%PARTIAL_DERIVATIVE_INDEX(:,localNodeIdx),'(" Index(:) :",8(X,I2))','(18X,8(X,I2))',err,error,*999)
ENDDO !localNodeIdx
IF(basis%NUMBER_OF_XI==3) THEN
CALL WriteString(DIAGNOSTIC_OUTPUT_TYPE," Local faces:",err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Number of local faces = ",basis%NUMBER_OF_LOCAL_FACES,err,error,*999)
DO localFaceIdx=1,basis%NUMBER_OF_LOCAL_FACES
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Local face = ",localFaceIdx,err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Number of nodes in local face = ", &
& basis%NUMBER_OF_NODES_IN_LOCAL_FACE(localFaceIdx),err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_NODES_IN_LOCAL_FACE(localFaceIdx),4,4, &
& basis%NODE_NUMBERS_IN_LOCAL_FACE(:,localFaceIdx),'(" Nodes in local face :",4(X,I2))','(33X,4(X,I2))', &
& err,error,*999)
DO localNodeIdx=1,basis%NUMBER_OF_NODES_IN_LOCAL_FACE(localFaceIdx)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Local face node: ",localNodeIdx,err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%DERIVATIVE_NUMBERS_IN_LOCAL_FACE(0,localNodeIdx, &
& localFaceIdx),4,4,basis%DERIVATIVE_NUMBERS_IN_LOCAL_FACE(1:,localNodeIdx,localFaceIdx), &
& '(" Derivatives in local face :",4(X,I2))','(33X,4(X,I2))',err,error,*999)
ENDDO !localNodeIdx
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_XI_COORDINATES-1,3,3, &
& basis%localFaceXiDirections(:,localFaceIdx),'(" Local face xi directions :",3(X,I2))','(33X,3(X,I2))', &
& err,error,*999)
CALL WriteStringFmtValue(DIAGNOSTIC_OUTPUT_TYPE," Local face xi normal : ", &
& basis%localFaceXiNormal(localFaceIdx),"I2",err,error,*999)
ENDDO !localFaceIdx
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,2*basis%NUMBER_OF_XI_COORDINATES+1,9,9, &
& basis%xiNormalLocalFace(-basis%NUMBER_OF_XI_COORDINATES:basis%NUMBER_OF_XI_COORDINATES), &
& '(" Xi normal local face :",9(X,I2))','(26X,9(X,I2))',err,error,*999)
ENDIF
CALL WriteString(DIAGNOSTIC_OUTPUT_TYPE," Local lines:",err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Number of local lines = ",basis%NUMBER_OF_LOCAL_LINES,err,error,*999)
DO localLineIdx=1,basis%NUMBER_OF_LOCAL_LINES
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Local line = ",localLineIdx,err,error,*999)
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Number of nodes in local line = ", &
& basis%NUMBER_OF_NODES_IN_LOCAL_LINE(localLineIdx),err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_NODES_IN_LOCAL_LINE(localLineIdx),4,4, &
& basis%NODE_NUMBERS_IN_LOCAL_LINE(:,localLineIdx),'(" Nodes in local line :",4(X,I2))','(33X,4(X,I2))', &
& err,error,*999)
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_NODES_IN_LOCAL_LINE(localLineIdx),4,4, &
& basis%DERIVATIVE_NUMBERS_IN_LOCAL_LINE(:,localLineIdx),'(" Derivatives in local line :",4(X,I2))', &
& '(33X,4(X,I2))',err,error,*999)
CALL WriteStringFmtValue(DIAGNOSTIC_OUTPUT_TYPE," Local line xi direction : ", &
& basis%localLineXiDirection(localLineIdx),"I2",err,error,*999)
IF(basis%NUMBER_OF_XI>1) THEN
CALL WriteStringVector(DIAGNOSTIC_OUTPUT_TYPE,1,1,basis%NUMBER_OF_XI-1,2,2,basis%localLineXiNormals(:,localLineIdx), &
& '(" Local line xi normals :",2(X,I2))','(31X,2(X,I2))',err,error,*999)
ENDIF
ENDDO !localLineIdx
IF(basis%NUMBER_OF_XI>=2) THEN
IF(basis%NUMBER_OF_XI==3) THEN
columnStart=1
columnStart2=-basis%NUMBER_OF_XI_COORDINATES
columnStop=2*basis%NUMBER_OF_XI_COORDINATES+1
columnStop2=basis%NUMBER_OF_XI_COORDINATES
ELSE
columnStart=1
columnStart2=1
columnStop=1
columnStop2=1
ENDIF
CALL WriteStringMatrix(DIAGNOSTIC_OUTPUT_TYPE,1,1,2*basis%NUMBER_OF_XI_COORDINATES+1, &
& columnStart,1,columnStop,9,9,basis%xiNormalsLocalLine(-basis%NUMBER_OF_XI_COORDINATES:basis%NUMBER_OF_XI_COORDINATES, &
& columnStart2:columnStop2),WRITE_STRING_MATRIX_NAME_AND_INDICES, &
& '(" Xi normal local line','(",I2,",:)',':",9(X,I2))','(31X,9(X,I2))',err,error,*999)
ENDIF
CALL WriteStringValue(DIAGNOSTIC_OUTPUT_TYPE," Number of sub-bases = ",basis%numberOfSubBases,err,error,*999)
ENDIF
EXITS("Basis_CreateFinish")
RETURN
999 ERRORSEXITS("Basis_CreateFinish",err,error)
RETURN 1
END SUBROUTINE Basis_CreateFinish
!
!================================================================================================================================
!
!>Starts the creation of a new basis
!>The default values of the BASIS attributes are:
!>- TYPE: 1 (BASIS_LAGRANGE_HERMITE_TP_TYPE)
!>- NUMBER_OF_XI: 3
!>- INTERPOLATION_XI: (1,1,1) (BASIS_LINEAR_LAGRANGE_INTERPOLATIONs)
!>- INTERPOLATION_TYPE: (1,1,1) (BASIS_LAGRANGE_INTERPOLATIONs)
!>- INTERPOLATION_ORDER: (1,1,1) (BASIS_LINEAR_INTERPOLATION_ORDERs)
!>- DEGENERATE: false
!>- COLLAPSED_XI: (4,4,4) (BASIS_NOT_COLLAPSEDs)
!>- QUADRATURE:
!> - TYPE: 1 (BASIS_LAGRANGE_HERMITE_TP_TYPE)
!> - NUMBER_OF_GAUSS_XI: (2,2,2)
!> - GAUSS_ORDER: 0
SUBROUTINE Basis_CreateStart(userNumber,basis,err,error,*)
!Argument variables
INTEGER(INTG), INTENT(IN) :: userNumber !<The user number of the basis to start the creation of
TYPE(BASIS_TYPE), POINTER :: basis !<On return, A pointer to the created basis. Must not be associated on entry.
INTEGER(INTG), INTENT(OUT) :: err !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: error !<The error string
!Local Variables
INTEGER(INTG) :: basisIdx
TYPE(BASIS_TYPE), POINTER :: newBasis
TYPE(BASIS_PTR_TYPE), ALLOCATABLE :: newBases(:)
TYPE(VARYING_STRING) :: localError
ENTERS("Basis_CreateStart",err,error,*999)
IF(ASSOCIATED(basis)) CALL FlagError("Basis is already associated",err,error,*999)
!See if basis number has already been created
CALL Basis_UserNumberFind(userNumber,basis,err,error,*999)
IF(ASSOCIATED(basis)) THEN
localError="A basis with a user number of "//TRIM(NumberToVString(userNumber,"*",err,error))// &
& " already exists."
CALL FlagError(localError,err,error,*999)
ENDIF
!Allocate new basis function and add it to the basis functions
ALLOCATE(newBases(basisFunctions%numberOfBasisFunctions+1),STAT=err)
IF(err/=0) CALL FlagError("Could not allocate new bases.",err,error,*999)
NULLIFY(newBasis)
CALL Basis_Initialise(newBasis,err,error,*999)
DO basisIdx=1,basisFunctions%numberOfBasisFunctions
newBases(basisIdx)%ptr=>basisFunctions%bases(basisIdx)%ptr
ENDDO !basisIdx
basisFunctions%numberOfBasisFunctions=basisFunctions%numberOfBasisFunctions+1
newBases(basisFunctions%numberOfBasisFunctions)%ptr=>newBasis
CALL MOVE_ALLOC(newBases,basisFunctions%bases)
!Set the basis parameters
newBasis%USER_NUMBER=userNumber
newBasis%FAMILY_NUMBER=0
newBasis%GLOBAL_NUMBER=basisFunctions%numberOfBasisFunctions
!Set the default basis parameters
newBasis%type=BASIS_LAGRANGE_HERMITE_TP_TYPE
newBasis%NUMBER_OF_XI=3
ALLOCATE(newBasis%INTERPOLATION_XI(3),STAT=err)
IF(err/=0) CALL FlagError("Could not allocate basis interpolation xi.",err,error,*999)
newBasis%INTERPOLATION_XI(1:3)=[BASIS_LINEAR_LAGRANGE_INTERPOLATION,BASIS_LINEAR_LAGRANGE_INTERPOLATION, &
& BASIS_LINEAR_LAGRANGE_INTERPOLATION]
ALLOCATE(newBasis%COLLAPSED_XI(3),STAT=err)
IF(err/=0) CALL FlagError("Could not allocate basis collapsed xi.",err,error,*999)
newBasis%COLLAPSED_XI(1:3)=BASIS_NOT_COLLAPSED
!Initialise the basis quadrature
NULLIFY(newBasis%QUADRATURE%basis)
CALL BASIS_QUADRATURE_INITIALISE(newBasis,err,error,*999)
basis=>newBasis
EXITS("Basis_CreateStart")
RETURN
999 IF(ASSOCIATED(newBasis)) CALL BASIS_DESTROY(newBasis,err,error,*998)
998 IF(ALLOCATED(newBases)) DEALLOCATE(newBases)
ERRORSEXITS("Basis_CreateStart",err,error)
RETURN 1
END SUBROUTINE Basis_CreateStart
!
!================================================================================================================================
!
!>Destroys a basis identified by its basis user number \see BASIS_ROUTINES::BASIS_DESTROY_FAMILY,OpenCMISS::Iron::cmfe_BasisDestroy
RECURSIVE SUBROUTINE BASIS_DESTROY_NUMBER(USER_NUMBER,err,error,*)
!Argument variables
INTEGER(INTG), INTENT(IN) :: USER_NUMBER !<The user number of the basis to destroy
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
!Local Variables
ENTERS("BASIS_DESTROY_NUMBER",err,error,*999)
CALL Basis_FamilyDestroy(USER_NUMBER,0,err,error,*999)
EXITS("BASIS_DESTROY_NUMBER")
RETURN
999 ERRORSEXITS("BASIS_DESTROY_NUMBER",err,error)
RETURN 1
END SUBROUTINE BASIS_DESTROY_NUMBER
!
!================================================================================================================================
!
!>Destroys a basis. \see BASIS_ROUTINES::BASIS_DESTROY_FAMILY,OpenCMISS::Iron::cmfe_BasisDestroy
RECURSIVE SUBROUTINE BASIS_DESTROY(BASIS,err,error,*)
!Argument variables
TYPE(BASIS_TYPE), POINTER :: BASIS !<A pointer to the basis to destroy
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
!Local Variables
INTEGER(INTG) :: USER_NUMBER
ENTERS("BASIS_DESTROY",err,error,*999)
IF(ASSOCIATED(BASIS)) THEN
USER_NUMBER=BASIS%USER_NUMBER
CALL Basis_FamilyDestroy(USER_NUMBER,0,err,error,*999)
!NULLIFY(BASIS)
ELSE
CALL FlagError("Basis is not associated.",err,error,*999)
ENDIF
EXITS("BASIS_DESTROY")
RETURN
999 ERRORSEXITS("BASIS_DESTROY",err,error)
RETURN 1
END SUBROUTINE BASIS_DESTROY
!
!================================================================================================================================
!
!>Evaluates the appropriate partial derivative index at position XI for the basis for double precision arguments.
!>Note for simplex basis functions the XI coordinates should exclude the last area coordinate.
FUNCTION BASIS_EVALUATE_XI_DP(BASIS,ELEMENT_PARAMETER_INDEX,PARTIAL_DERIV_INDEX,XI,err,error)
!Argument variables
TYPE(BASIS_TYPE), POINTER :: BASIS !<A pointer to the basis
INTEGER(INTG), INTENT(IN) :: ELEMENT_PARAMETER_INDEX !<The element parameter index to evaluate i.e., the local basis index within the element basis.
INTEGER(INTG), INTENT(IN) :: PARTIAL_DERIV_INDEX !<The partial derivative index to evaluate \see CONSTANTS_PartialDerivativeConstants
REAL(DP), INTENT(IN) :: XI(:) !<The Xi position to evaluate the basis function at
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
!Function variable
REAL(DP) :: BASIS_EVALUATE_XI_DP
!Local Variables
INTEGER(INTG) :: nn,nk
REAL(DP) :: XIL(SIZE(XI,1)+1)
TYPE(VARYING_STRING) :: LOCAL_ERROR
ENTERS("BASIS_EVALUATE_XI_DP",err,error,*999)
BASIS_EVALUATE_XI_DP=0.0_DP
IF(ASSOCIATED(BASIS)) THEN
IF(ELEMENT_PARAMETER_INDEX>0.AND.ELEMENT_PARAMETER_INDEX<=BASIS%NUMBER_OF_ELEMENT_PARAMETERS) THEN
SELECT CASE(BASIS%TYPE)
CASE(BASIS_LAGRANGE_HERMITE_TP_TYPE)
nn=BASIS%ELEMENT_PARAMETER_INDEX_INV(1,ELEMENT_PARAMETER_INDEX)
nk=BASIS%ELEMENT_PARAMETER_INDEX_INV(2,ELEMENT_PARAMETER_INDEX)
BASIS_EVALUATE_XI_DP=BASIS_LHTP_BASIS_EVALUATE(BASIS,nn,nk,PARTIAL_DERIV_INDEX,XI,err,error)
IF(ERR/=0) GOTO 999
CASE(BASIS_SIMPLEX_TYPE)
!Create the area coordinates from the xi coordinates
CALL Basis_XiToAreaCoordinates(XI(1:SIZE(XI,1)),XIL(1:SIZE(XI,1)+1),err,error,*999)
nn=BASIS%ELEMENT_PARAMETER_INDEX_INV(1,ELEMENT_PARAMETER_INDEX)
BASIS_EVALUATE_XI_DP=BASIS_SIMPLEX_BASIS_EVALUATE(BASIS,nn,PARTIAL_DERIV_INDEX,XIL,err,error)
IF(ERR/=0) GOTO 999
CASE(BASIS_RADIAL_TYPE)
IF(ERR/=0) GOTO 999
CASE DEFAULT
LOCAL_ERROR="Basis type "//TRIM(NumberToVString(BASIS%TYPE,"*",err,error))//" is invalid or not implemented."
CALL FlagError(LOCAL_ERROR,err,error,*999)
END SELECT
ELSE
LOCAL_ERROR="The specified element parameter index of "// &
& TRIM(NumberToVString(ELEMENT_PARAMETER_INDEX,"*",err,error))// &
& " is invalid. The index must be > 0 and <= "// &
& TRIM(NumberToVString(BASIS%NUMBER_OF_ELEMENT_PARAMETERS,"*",err,error))//"."
CALL FlagError(LOCAL_ERROR,err,error,*999)
ENDIF
ELSE
CALL FlagError("Basis is not associated.",err,error,*999)
ENDIF
EXITS("BASIS_EVALUATE_XI_DP")
RETURN
999 ERRORSEXITS("BASIS_EVALUATE_XI_DP",err,error)
RETURN
END FUNCTION BASIS_EVALUATE_XI_DP
!
!================================================================================================================================
!
!>Converts a xi location on a boundary face or line to a full xi location.
SUBROUTINE Basis_BoundaryXiToXi(basis,localLineFaceNumber,boundaryXi,fullXi,err,error,*)
!Argument variables
TYPE(BASIS_TYPE), POINTER :: basis !<A pointer to the basis to convert the boundary xi for
INTEGER(INTG), INTENT(IN) :: localLineFaceNumber !<The local line/face number containing the boundary xi
REAL(DP), INTENT(IN) :: boundaryXi(:) !<The boundary xi location to convert to the full xi location. Note the size of the boundary xi array will determine if we are dealing with a line xi or a face xi.