File: libmeshb7_mod.f90

package info (click to toggle)
libmeshb 7.80-4
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 2,848 kB
  • sloc: ansic: 12,810; f90: 1,146; fortran: 406; makefile: 218
file content (846 lines) | stat: -rw-r--r-- 34,659 bytes parent folder | download
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
!----------------------------------------------------------
!
!                       LIBMESH V 7.79
! 
!----------------------------------------------------------
!
!   Description:         handles .meshb file format I/O
!   Author:              Loic MARECHAL
!   Creation date:       dec 08 2015
!   Last modification:   feb 12 2024
!
!----------------------------------------------------------

module libmeshb7
  
  use iso_fortran_env
  use, intrinsic :: iso_c_binding, only: c_int,c_long,c_loc,c_ptr,c_null_ptr
  
  implicit none
  
  ! Procedures definition
  integer(int64) , external :: gmfopenmeshf77
  integer(int32) , external :: gmfclosemeshf77
  integer(int32) , external :: GmfStatKwdf77
  integer(int32) , external :: gmfsetkwdf77
  integer(int32) , external :: gmfgotokwdf77
  integer(int32) , external :: gmfsethonodesorderingf77
  integer(int32) , external :: gmfgetlinef77
  integer(int32) , external :: gmfsetlinef77
  integer(int32) , external :: gmfgetblockf77
  integer(int32) , external :: gmfsetblockf77
  
  
  ! Parameters definition
  integer(int32), parameter :: gmfmaxtyp=1000
  integer(int32), parameter :: gmfmaxkwd=227
  integer(int32), parameter :: gmfread=1
  integer(int32), parameter :: gmfwrite=2
  integer(int32), parameter :: gmfsca=1
  integer(int32), parameter :: gmfvec=2
  integer(int32), parameter :: gmfsymmat=3
  integer(int32), parameter :: gmfmat=4
  integer(int32), parameter :: gmffloat=8
  integer(int32), parameter :: gmfdouble=9
  integer(int32), parameter :: gmfint=10
  integer(int32), parameter :: gmflong=11
  integer(int32), parameter :: gmfinttab=14
  integer(int32), parameter :: gmflongtab=15
  integer(int32), parameter :: gmffloatvec=12
  integer(int32), parameter :: gmfdoublevec=13
  integer(int32), parameter :: gmfintvec=14
  integer(int32), parameter :: gmflongvec=15
  integer(int32), parameter :: gmfargtab=100
  integer(int32), parameter :: gmfarglst=101
  
  ! Keywords list  
  integer(int32), parameter :: gmfdimension=3
  integer(int32), parameter :: gmfmeshversionformatted=1
  integer(int32), parameter :: gmfvertices=4
  integer(int32), parameter :: gmfedges=5
  integer(int32), parameter :: gmftriangles=6
  integer(int32), parameter :: gmfquadrilaterals=7
  integer(int32), parameter :: gmftetrahedra=8
  integer(int32), parameter :: gmfprisms=9
  integer(int32), parameter :: gmfhexahedra=10
  integer(int32), parameter :: gmfiterationsall=11
  integer(int32), parameter :: gmftimesall=12
  integer(int32), parameter :: gmfcorners=13
  integer(int32), parameter :: gmfridges=14
  integer(int32), parameter :: gmfrequiredvertices=15
  integer(int32), parameter :: gmfrequirededges=16
  integer(int32), parameter :: gmfrequiredtriangles=17
  integer(int32), parameter :: gmfrequiredquadrilaterals=18
  integer(int32), parameter :: gmftangentatedgevertices=19
  integer(int32), parameter :: gmfnormalatvertices=20
  integer(int32), parameter :: gmfnormalattrianglevertices=21
  integer(int32), parameter :: gmfnormalatquadrilateralvertices=22
  integer(int32), parameter :: gmfangleofcornerbound=23
  integer(int32), parameter :: gmftrianglesp2=24
  integer(int32), parameter :: gmfedgesp2=25
  integer(int32), parameter :: gmfsolatpyramids=26
  integer(int32), parameter :: gmfquadrilateralsq2=27
  integer(int32), parameter :: gmfisolatpyramids=28
  integer(int32), parameter :: gmfsubdomainfromgeom=29
  integer(int32), parameter :: gmftetrahedrap2=30
  integer(int32), parameter :: gmffault_neartri=31
  integer(int32), parameter :: gmffault_inter=32
  integer(int32), parameter :: gmfhexahedraq2=33
  integer(int32), parameter :: gmfextraverticesatedges=34
  integer(int32), parameter :: gmfextraverticesattriangles=35
  integer(int32), parameter :: gmfextraverticesatquadrilaterals=36
  integer(int32), parameter :: gmfextraverticesattetrahedra=37
  integer(int32), parameter :: gmfextraverticesatprisms=38
  integer(int32), parameter :: gmfextraverticesathexahedra=39
  integer(int32), parameter :: gmfverticesongeometricvertices=40
  integer(int32), parameter :: gmfverticesongeometricedges=41
  integer(int32), parameter :: gmfverticesongeometrictriangles=42
  integer(int32), parameter :: gmfverticesongeometricquadrilaterals=43
  integer(int32), parameter :: gmfedgesongeometricedges=44
  integer(int32), parameter :: gmffault_freeedge=45
  integer(int32), parameter :: gmfpolyhedra=46
  integer(int32), parameter :: gmfpolygons=47
  integer(int32), parameter :: gmffault_overlap=48
  integer(int32), parameter :: gmfpyramids=49
  integer(int32), parameter :: gmfboundingbox=50
  integer(int32), parameter :: gmfbody=51
  integer(int32), parameter :: gmfprivatetable=52
  integer(int32), parameter :: gmffault_badshape=53
  integer(int32), parameter :: gmfend=54
  integer(int32), parameter :: gmftrianglesongeometrictriangles=55
  integer(int32), parameter :: gmftrianglesongeometricquadrilaterals=56
  integer(int32), parameter :: gmfquadrilateralsongeometrictriangles=57
  integer(int32), parameter :: gmfquadrilateralsongeometricquadrilaterals=58
  integer(int32), parameter :: gmftangents=59
  integer(int32), parameter :: gmfnormals=60
  integer(int32), parameter :: gmftangentatvertices=61
  integer(int32), parameter :: gmfsolatvertices=62
  integer(int32), parameter :: gmfsolatedges=63
  integer(int32), parameter :: gmfsolattriangles=64
  integer(int32), parameter :: gmfsolatquadrilaterals=65
  integer(int32), parameter :: gmfsolattetrahedra=66
  integer(int32), parameter :: gmfsolatprisms=67
  integer(int32), parameter :: gmfsolathexahedra=68
  integer(int32), parameter :: gmfdsolatvertices=69
  integer(int32), parameter :: gmfisolatvertices=70
  integer(int32), parameter :: gmfisolatedges=71
  integer(int32), parameter :: gmfisolattriangles=72
  integer(int32), parameter :: gmfisolatquadrilaterals=73
  integer(int32), parameter :: gmfisolattetrahedra=74
  integer(int32), parameter :: gmfisolatprisms=75
  integer(int32), parameter :: gmfisolathexahedra=76
  integer(int32), parameter :: gmfiterations=77
  integer(int32), parameter :: gmftime=78
  integer(int32), parameter :: gmffault_smalltri=79
  integer(int32), parameter :: gmfcoarsehexahedra=80
  integer(int32), parameter :: gmfcomments=81
  integer(int32), parameter :: gmfperiodicvertices=82
  integer(int32), parameter :: gmfperiodicedges=83
  integer(int32), parameter :: gmfperiodictriangles=84
  integer(int32), parameter :: gmfperiodicquadrilaterals=85
  integer(int32), parameter :: gmfprismsp2=86
  integer(int32), parameter :: gmfpyramidsp2=87
  integer(int32), parameter :: gmfquadrilateralsq3=88
  integer(int32), parameter :: gmfquadrilateralsq4=89
  integer(int32), parameter :: gmftrianglesp3=90
  integer(int32), parameter :: gmftrianglesp4=91
  integer(int32), parameter :: gmfedgesp3=92
  integer(int32), parameter :: gmfedgesp4=93
  integer(int32), parameter :: gmfirefgroups=94
  integer(int32), parameter :: gmfdrefgroups=95
  integer(int32), parameter :: gmftetrahedrap3=96
  integer(int32), parameter :: gmftetrahedrap4=97
  integer(int32), parameter :: gmfhexahedraq3=98
  integer(int32), parameter :: gmfhexahedraq4=99
  integer(int32), parameter :: gmfpyramidsp3=100
  integer(int32), parameter :: gmfpyramidsp4=101
  integer(int32), parameter :: gmfprismsp3=102
  integer(int32), parameter :: gmfprismsp4=103
  integer(int32), parameter :: gmfhosolatedgesp1=104
  integer(int32), parameter :: gmfhosolatedgesp2=105
  integer(int32), parameter :: gmfhosolatedgesp3=106
  integer(int32), parameter :: gmfhosolattrianglesp1=107
  integer(int32), parameter :: gmfhosolattrianglesp2=108
  integer(int32), parameter :: gmfhosolattrianglesp3=109
  integer(int32), parameter :: gmfhosolatquadrilateralsq1=110
  integer(int32), parameter :: gmfhosolatquadrilateralsq2=111
  integer(int32), parameter :: gmfhosolatquadrilateralsq3=112
  integer(int32), parameter :: gmfhosolattetrahedrap1=113
  integer(int32), parameter :: gmfhosolattetrahedrap2=114
  integer(int32), parameter :: gmfhosolattetrahedrap3=115
  integer(int32), parameter :: gmfhosolatpyramidsp1=116
  integer(int32), parameter :: gmfhosolatpyramidsp2=117
  integer(int32), parameter :: gmfhosolatpyramidsp3=118
  integer(int32), parameter :: gmfhosolatprismsp1=119
  integer(int32), parameter :: gmfhosolatprismsp2=120
  integer(int32), parameter :: gmfhosolatprismsp3=121
  integer(int32), parameter :: gmfhosolathexahedraq1=122
  integer(int32), parameter :: gmfhosolathexahedraq2=123
  integer(int32), parameter :: gmfhosolathexahedraq3=124
  integer(int32), parameter :: gmfbezierbasis=125
  integer(int32), parameter :: gmfbyteflow=126
  integer(int32), parameter :: gmfedgesp2ordering=127
  integer(int32), parameter :: gmfedgesp3ordering=128
  integer(int32), parameter :: gmftrianglesp2ordering=129
  integer(int32), parameter :: gmftrianglesp3ordering=130
  integer(int32), parameter :: gmfquadrilateralsq2ordering=131
  integer(int32), parameter :: gmfquadrilateralsq3ordering=132
  integer(int32), parameter :: gmftetrahedrap2ordering=133
  integer(int32), parameter :: gmftetrahedrap3ordering=134
  integer(int32), parameter :: gmfpyramidsp2ordering=135
  integer(int32), parameter :: gmfpyramidsp3ordering=136
  integer(int32), parameter :: gmfprismsp2ordering=137
  integer(int32), parameter :: gmfprismsp3ordering=138
  integer(int32), parameter :: gmfhexahedraq2ordering=139
  integer(int32), parameter :: gmfhexahedraq3ordering=140
  integer(int32), parameter :: gmfedgesp1ordering=141
  integer(int32), parameter :: gmfedgesp4ordering=142
  integer(int32), parameter :: gmftrianglesp1ordering=143
  integer(int32), parameter :: gmftrianglesp4ordering=144
  integer(int32), parameter :: gmfquadrilateralsq1ordering=145
  integer(int32), parameter :: gmfquadrilateralsq4ordering=146
  integer(int32), parameter :: gmftetrahedrap1ordering=147
  integer(int32), parameter :: gmftetrahedrap4ordering=148
  integer(int32), parameter :: gmfpyramidsp1ordering=149
  integer(int32), parameter :: gmfpyramidsp4ordering=150
  integer(int32), parameter :: gmfprismsp1ordering=151
  integer(int32), parameter :: gmfprismsp4ordering=152
  integer(int32), parameter :: gmfhexahedraq1ordering=153
  integer(int32), parameter :: gmfhexahedraq4ordering=154
  integer(int32), parameter :: gmffloatingpointprecision=155
  integer(int32), parameter :: gmfhosolatedgesp4=156
  integer(int32), parameter :: gmfhosolattrianglesp4=157
  integer(int32), parameter :: gmfhosolatquadrilateralsq4=158
  integer(int32), parameter :: gmfhosolattetrahedrap4=159
  integer(int32), parameter :: gmfhosolatpyramidsp4=160
  integer(int32), parameter :: gmfhosolatprismsp4=161
  integer(int32), parameter :: gmfhosolathexahedraq4=162
  integer(int32), parameter :: gmfhosolatedgesp1nodespositions=163
  integer(int32), parameter :: gmfhosolatedgesp2nodespositions=164
  integer(int32), parameter :: gmfhosolatedgesp3nodespositions=165
  integer(int32), parameter :: gmfhosolatedgesp4nodespositions=166
  integer(int32), parameter :: gmfhosolattrianglesp1nodespositions=167
  integer(int32), parameter :: gmfhosolattrianglesp2nodespositions=168
  integer(int32), parameter :: gmfhosolattrianglesp3nodespositions=169
  integer(int32), parameter :: gmfhosolattrianglesp4nodespositions=170
  integer(int32), parameter :: gmfhosolatquadrilateralsq1nodespositions=171
  integer(int32), parameter :: gmfhosolatquadrilateralsq2nodespositions=172
  integer(int32), parameter :: gmfhosolatquadrilateralsq3nodespositions=173
  integer(int32), parameter :: gmfhosolatquadrilateralsq4nodespositions=174
  integer(int32), parameter :: gmfhosolattetrahedrap1nodespositions=175
  integer(int32), parameter :: gmfhosolattetrahedrap2nodespositions=176
  integer(int32), parameter :: gmfhosolattetrahedrap3nodespositions=177
  integer(int32), parameter :: gmfhosolattetrahedrap4nodespositions=178
  integer(int32), parameter :: gmfhosolatpyramidsp1nodespositions=179
  integer(int32), parameter :: gmfhosolatpyramidsp2nodespositions=180
  integer(int32), parameter :: gmfhosolatpyramidsp3nodespositions=181
  integer(int32), parameter :: gmfhosolatpyramidsp4nodespositions=182
  integer(int32), parameter :: gmfhosolatprismsp1nodespositions=183
  integer(int32), parameter :: gmfhosolatprismsp2nodespositions=184
  integer(int32), parameter :: gmfhosolatprismsp3nodespositions=185
  integer(int32), parameter :: gmfhosolatprismsp4nodespositions=186
  integer(int32), parameter :: gmfhosolathexahedraq1nodespositions=187
  integer(int32), parameter :: gmfhosolathexahedraq2nodespositions=188
  integer(int32), parameter :: gmfhosolathexahedraq3nodespositions=189
  integer(int32), parameter :: gmfhosolathexahedraq4nodespositions=190
  integer(int32), parameter :: gmfedgesreferenceelement=191
  integer(int32), parameter :: gmftrianglereferenceelement=192
  integer(int32), parameter :: gmfquadrilateralreferenceelement=193
  integer(int32), parameter :: gmftetrahedronreferenceelement=194
  integer(int32), parameter :: gmfpyramidreferenceelement=195
  integer(int32), parameter :: gmfprismreferenceelement=196
  integer(int32), parameter :: gmfhexahedronreferenceelement=197
  integer(int32), parameter :: gmfboundarylayers=198
  integer(int32), parameter :: gmfreferencestrings=199
  integer(int32), parameter :: gmfprisms9=200
  integer(int32), parameter :: gmfhexahedra12=201
  integer(int32), parameter :: gmfquadrilaterals6=202
  integer(int32), parameter :: gmfboundarypolygonheaders=203
  integer(int32), parameter :: gmfboundarypolygonvertices=204
  integer(int32), parameter :: gmfinnerpolygonheaders=205
  integer(int32), parameter :: gmfinnerpolygonvertices=206
  integer(int32), parameter :: gmfpolyhedraheaders=207
  integer(int32), parameter :: gmfpolyhedrafaces=208
  integer(int32), parameter :: gmfdomains=209
  integer(int32), parameter :: gmfverticesgid=210
  integer(int32), parameter :: gmfedgesgid=211
  integer(int32), parameter :: gmftrianglesgid=212
  integer(int32), parameter :: gmfquadrilateralsgid=213
  integer(int32), parameter :: gmftetrahedragid=214
  integer(int32), parameter :: gmfpyramidsgid=215
  integer(int32), parameter :: gmfprismsgid=216
  integer(int32), parameter :: gmfhexahedragid=217
  integer(int32), parameter :: gmfsolatboundarypolygons=218
  integer(int32), parameter :: gmfsolatpolyhedra=219
  integer(int32), parameter :: gmfverticesongeometrynodes=220
  integer(int32), parameter :: gmfverticesongeometryedges=221
  integer(int32), parameter :: gmfedgesongeometryedges=222
  integer(int32), parameter :: gmfverticesongeometryfaces=223
  integer(int32), parameter :: gmfedgesongeometryfaces=224
  integer(int32), parameter :: gmftrianglesongeometryfaces=225
  integer(int32), parameter :: gmfquadrialteralsongeometryfaces=226
  integer(int32), parameter :: gmfmeshongeometry=227
  
  interface     GmfStatKwdF90
    module procedure GmfStatKwdF90_0 !> vertices & nodes
    module procedure GmfStatKwdF90_1 !> solutions       
  end interface GmfStatKwdF90
  
  interface     GmfSetKwdF90
    module procedure GmfSetKwdF90_0 !> vertices & nodes
    module procedure GmfSetKwdF90_1 !> solutions       
  end interface GmfSetKwdF90
  
  interface     GmfGetLineF90
    module procedure GmfGetLineF90_i
    module procedure GmfGetLineF90_d
  end interface GmfGetLineF90
  
  interface     GmfSetLineF90
    module procedure GmfSetLineF90_i
    module procedure GmfSetLineF90_d
    module procedure GmfSetLineF90_sol_i
    module procedure GmfSetLineF90_sol_d
  end interface GmfSetLineF90
  
  interface     GmfGetBlockF90
   !module procedure GmfGetBlockF90_00
    module procedure GmfGetBlockF90_01    !> nodes    + ref
    module procedure GmfGetBlockF90_01Bis !> nodes   
    module procedure GmfGetBlockF90_02    !> vertices + ref
    module procedure GmfGetBlockF90_02Bis !> solutions     
  end interface GmfGetBlockF90
  
  interface     GmfSetBlockF90
   !module procedure GmfGetBlockF90_00
    module procedure GmfSetBlockF90_01    !> nodes    + ref
    module procedure GmfSetBlockF90_02    !> vertices + ref
    module procedure GmfSetBlockF90_02Bis !> solutions     
  end interface GmfSetBlockF90
  
contains
  
  function  GmfOpenMeshF90(name, GmfKey, ver, dim) result(unit)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    character(*)  , intent(in)    :: name
    integer(int32), intent(in)    :: GmfKey
    integer(int32), intent(inout) :: ver
    integer(int32), intent(inout) :: dim
    integer(int64)                :: unit
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    unit = GmfOpenMeshf77(trim(name), GmfKey, ver, dim)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfOpenMeshF90
  
  function     GmfCloseMeshF90(unit) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64) :: unit
    integer(int32) :: res
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    res=GmfCloseMeshF77(unit)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfCloseMeshF90
  
  function     GmfStatKwdF90_0(unit, GmfKey) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64) :: unit
    integer(int32) :: GmfKey
    integer(int32) :: Nmb
    integer(int32) :: res
    !>
    integer(int32) :: t(1),d,ho,s
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    res = GmfStatKwdf77(unit, GmfKey, 0, 0, t(1), 0, 0)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfStatKwdF90_0
  
  function     GmfStatKwdF90_1(unit, GmfKey, r, s, t, d, ho) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64) :: unit
    integer(int32) :: GmfKey
    integer(int32) :: Nmb
    integer(int32) :: res
    !>
    integer(int32) :: r,s,t(:),d,ho
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    res = GmfStatKwdf77(unit, GmfKey, r, s, t(1), d, ho)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfStatKwdF90_1
  
  function     GmfGotoKwdF90(unit, GmfKey) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64) :: unit
    integer(int32) :: GmfKey
    integer(int32) :: res
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    res=GmfgotokwdF77(unit, GmfKey)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfGotoKwdF90
  
  function    GmfSetKwdF90_0(unit, GmfKey, Nmb) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64) :: unit
    integer(int32) :: GmfKey
    integer(int32) :: Nmb
    !>
    integer(int32) :: t(1),d,ho,s
    integer(int32) :: res
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    res = GmfSetKwdF77(unit, GmfKey, Nmb, 0, t(1), 0, ho)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfSetKwdF90_0
  
  function     GmfSetKwdF90_1(unit, GmfKey, Nmb, d, t, s, ho) result(res)
   !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
   integer(int64) :: unit
   integer(int32) :: GmfKey
   integer(int32) :: Nmb
   integer(int32) :: t(:)
   integer(int32) :: d,ho,s
   integer(int32) :: res
   !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
   !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
   res = GmfSetKwdF77(unit, GmfKey, Nmb, d, t, s, ho)
   !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
   return
  end function GmfSetKwdF90_1
  
  function     GmfSetHONodesOrderingF90(unit, GmfKey, BasTab, OrdTab) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64) :: unit
    integer(int32) :: GmfKey
    integer(int32) :: BasTab(:,:)
    integer(int32) :: OrdTab(:,:)
    integer(int32) :: res
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    res=GmfSetHONodesOrderingF77(unit,GmfKey,BasTab,OrdTab)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfSetHONodesOrderingF90
  
  function     GmfGetLineF90_i(unit, GmfKey, Tab, Ref) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    !> Reading Nodes and Ref
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64) :: unit
    integer(int32) :: GmfKey
    integer(int32) :: Tab(:)
    integer(int32) :: Ref
    integer(int32) :: res
    !>
    real(real64)   :: dTab(1)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    res=GmfGetLineF77(unit, GmfKey, Tab(1), dTab(1), Ref)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfGetLineF90_i
  
  function     GmfGetLineF90_d(unit, GmfKey, Tab, Ref) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    !> Reading Vertices and Ref
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64) :: unit
    integer(int32) :: GmfKey
    real(real64)   :: Tab(:)
    integer(int32) :: Ref
    integer(int32) :: res
    !>
    integer(int32) :: iTab(1)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    res=GmfGetLineF77(unit, GmfKey, iTab(1), Tab(1), Ref)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
  end function GmfGetLineF90_d
  
  function     GmfSetLineF90_i(unit, GmfKey, Tab, Ref) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    !> Writting Nodes and Ref
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64) :: unit
    integer(int32) :: GmfKey
    integer(int32) :: Tab(:)
    integer(int32) :: Ref
    integer(int32) :: res
    !>
    real(real64)   :: dTab(1)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    res=GmfSetLineF77(unit, GmfKey, Tab(1), dTab(1), Ref)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfSetLineF90_i
  
  function     GmfSetLineF90_d(unit, GmfKey, Tab, Ref) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    !> Writting Vertices and Ref
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64) :: unit
    integer(int32) :: GmfKey
    real(real64)   :: Tab(:)
    integer(int32) :: Ref
    integer(int32) :: res
    !>
    integer(int32) :: iTab(1)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    res=GmfSetLineF77(unit, GmfKey, iTab(1), Tab(1), Ref)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfSetLineF90_d
  
  function     GmfSetLineF90_sol_i(unit, GmfKey, iTab) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    !> Writting Nodes and Ref
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64) :: unit
    integer(int32) :: GmfKey
    integer(int32) :: iTab(:)
    integer(int32) :: iRef
    integer(int32) :: res
    !>
    real(real64)   :: dTab(1)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    res=GmfSetLineF77(unit, GmfKey, iTab(1), dTab(1), iRef)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfSetLineF90_sol_i
  
  function     GmfSetLineF90_sol_d(unit, GmfKey, dTab) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    !> Writting Vertices and Ref
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64) :: unit
    integer(int32) :: GmfKey
    real(real64)   :: dTab(:)
    integer(int32) :: res
    !>
    integer(int32) :: iRef
    integer(int32) :: iTab(1)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    res=GmfSetLineF77(unit, GmfKey, iTab(1), dTab(1), iRef)
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfSetLineF90_sol_d
    
  function     GmfGetBlockF90_00(unit, GmfKey, ad0, ad1, iTab, dTab, Ref) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64), intent(in)    :: unit
    integer(int32), intent(in)    :: GmfKey
    integer(int32), intent(in)    :: ad0
    integer(int32), intent(in)    :: ad1
    integer(int32), intent(inout) :: iTab(:,:)
    real(real64)  , intent(inout) :: dTab(:,:)
    integer(int32), intent(inout) :: Ref(  :)
    integer(int32)                :: res
    !>
    integer(int32)                :: Nmb
    integer(int32), pointer       :: map(:)=>null()
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    Nmb=ad1-ad0+1
    
    res=GmfGetBlockF77(unit       ,&
    &                  GmfKey     ,&
    &                  ad0        ,&
    &                  ad1        ,&
    &                  int32      ,&
    &                  map        ,&
    &                  iTab(1,  1),&
    &                  iTab(1,Nmb),&
    &                  dTab(1,  1),&
    &                  dTab(1,Nmb),&
    &                   Ref(    1),&
    &                   Ref(  Nmb) )
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfGetBlockF90_00
  
  function     GmfGetBlockF90_01(unit, GmfKey, ad0, ad1, Tab, Ref) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64), intent(in)    :: unit
    integer(int32), intent(in)    :: GmfKey
    integer(int32), intent(in)    :: ad0
    integer(int32), intent(in)    :: ad1
    integer(int32), intent(inout) :: Tab(:,:)
    integer(int32), intent(inout) :: Ref(  :)
    integer(int32)                :: res
    !>
    integer(int32)                :: Nmb
    real(real64)                  :: dTab(1)
    integer(int32), pointer       :: map(:)=>null()

    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    
    Nmb=ad1-ad0+1
    
    print '("GmfGetBlockF90_01 (ad0,ad1)=(",i0,",",i0,") Nmb=",i0)',ad0,ad1,Nmb
    print '("GmfGetBlockF90_01 size(Tab)=",i0,"x",i0)',size(Tab,1),size(Tab,2)
    print '("GmfGetBlockF90_01 size(Ref)=  ",i0)',size(Ref)
    
    res=GmfGetBlockF77(unit       ,&
    &                  GmfKey     ,&
    &                  ad0        ,&
    &                  ad1        ,&
    &                  int32      ,&
    &                  map        ,&
    &                  Tab(1,  1) ,&
    &                  Tab(1,Nmb) ,&
    &                  dTab(1)    ,&
    &                  dTab(1)    ,&
    &                  Ref(    1) ,&
    &                  Ref(  Nmb)  )
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfGetBlockF90_01

  function     GmfGetBlockF90_01Bis(unit, GmfKey, ad0, ad1, Tab) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64), intent(in)    :: unit
    integer(int32), intent(in)    :: GmfKey
    integer(int32), intent(in)    :: ad0
    integer(int32), intent(in)    :: ad1
    integer(int32), intent(inout) :: Tab(:,:)
    integer(int32)                :: res
    !>
    integer(int32)                :: Nmb
    real(real64)                  :: dTab(1)
    integer(int32)                :: Ref(1)
    integer(int32), pointer       :: map(:)=>null()

    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    
    Nmb=ad1-ad0+1
    
    print '("GmfGetBlockF90_01Bis (ad0,ad1)=(",i0,",",i0,") Nmb=",i0)',ad0,ad1,Nmb
    print '("GmfGetBlockF90_01Bis size(Tab)=",i0,"x",i0)',size(Tab,1),size(Tab,2)
    print '("GmfGetBlockF90_01Bis size(Ref)=  ",i0)',size(Ref)
    
    res=GmfGetBlockF77(unit       ,&
    &                  GmfKey     ,&
    &                  ad0        ,&
    &                  ad1        ,&
    &                  int32      ,&
    &                  map        ,&
    &                  Tab(1,  1) ,&
    &                  Tab(1,Nmb) ,&
    &                  dTab(1)    ,&
    &                  dTab(1)    ,&
    &                  Ref(    1) ,&
    &                  Ref(    1)  )
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfGetBlockF90_01Bis  
  
  function     GmfGetBlockF90_02(unit, GmfKey, ad0, ad1, Tab, Ref) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64), intent(in)    :: unit
    integer(int32), intent(in)    :: GmfKey
    integer(int32), intent(in)    :: ad0
    integer(int32), intent(in)    :: ad1
    real(real64)  , intent(inout) :: Tab(:,:)
    integer(int32), intent(inout) :: Ref(  :)
    integer(int32)                :: res
    !>
    integer(int32) :: iTab(1)
    integer(int32) :: Nmb
    integer(int32), pointer       :: map(:)=>null()
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    Nmb=ad1-ad0+1
    
    print '("GmfGetBlockF90_02 (ad0,ad1)=(",i0,",",i0,") Nmb=",i0)',ad0,ad1,Nmb
    print '("GmfGetBlockF90_02 size(Tab)=",i0,"x",i0)',size(Tab,1),size(Tab,2)
    print '("GmfGetBlockF90_02 size(Ref)=  ",i0)',size(Ref)
    
    res=GmfGetBlockF77(unit       ,&
    &                  GmfKey     ,&
    &                  ad0        ,&
    &                  ad1        ,&
    &                  int32      ,&
    &                  map        ,&
    &                  iTab(1)    ,&
    &                  iTab(1)    ,&
    &                  Tab(1,  1) ,&
    &                  Tab(1,Nmb) ,&
    &                  Ref(    1) ,&
    &                  Ref(  Nmb)  )
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfGetBlockF90_02
  
  function     GmfGetBlockF90_02Bis(unit, GmfKey, ad0, ad1, Tab) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64), intent(in)    :: unit
    integer(int32), intent(in)    :: GmfKey
    integer(int32), intent(in)    :: ad0
    integer(int32), intent(in)    :: ad1
    real(real64)  , intent(inout) :: Tab(:,:)
    integer(int32)                :: res
    !>
    integer(int32) :: Ref (1)
    integer(int32) :: iTab(1)
    integer(int32) :: Nmb
    integer(int32), pointer       :: map(:)=>null()
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    Nmb=ad1-ad0+1
    
    print '("GmfGetBlockF90_02Bis (ad0,ad1)=(",i0,",",i0,") Nmb=",i0)',ad0,ad1,Nmb
    print '("GmfGetBlockF90_02Bis size(Tab)=",i0,"x",i0)',size(Tab,1),size(Tab,2)
    print '("GmfGetBlockF90_02Bis size(Ref)=  ",i0)',size(Ref)
    
    res=GmfGetBlockF77(unit       ,&
    &                  GmfKey     ,&
    &                  ad0        ,&
    &                  ad1        ,&
    &                  int32      ,&
    &                  map        ,&
    &                  iTab(1)    ,&
    &                  iTab(1)    ,&
    &                  Tab(1,  1) ,&
    &                  Tab(1,Nmb) ,&
    &                  Ref(    1) ,&
    &                  Ref(    1)  )
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfGetBlockF90_02Bis
  
  function     GmfSetBlockF90_01(unit, GmfKey, ad0, ad1, Tab, Ref) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64), intent(in)    :: unit
    integer(int32), intent(in)    :: GmfKey
    integer(int32), intent(in)    :: ad0
    integer(int32), intent(in)    :: ad1
    integer(int32), intent(inout) :: Tab(:,:)
    integer(int32), intent(inout) :: Ref(  :)
    integer(int32)                :: res
    !>
    integer(int32)                :: Nmb
    real(real64)                  :: dTab(1)
    integer(int32), pointer       :: map(:)=>null()

    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    Nmb=ad1-ad0+1
    
    print '("GmfSetBlockF90_01 (ad0,ad1)=(",i0,",",i0,") Nmb=",i0)',ad0,ad1,Nmb
    print '("GmfSetBlockF90_01 size(Tab)=",i0,"x",i0)',size(Tab,1),size(Tab,2)
    print '("GmfSetBlockF90_01 size(Ref)=  ",i0)',size(Ref)
    
    res=GmfSetBlockF77(unit       ,&
    &                  GmfKey     ,&
    &                  ad0        ,&
    &                  ad1        ,&
    &                  int32      ,&
    &                  map        ,&
    &                  Tab(1,  1) ,&
    &                  Tab(1,Nmb) ,&
    &                  dTab(1)    ,&
    &                  dTab(1)    ,&
    &                  Ref(    1) ,&
    &                  Ref(  Nmb)  )
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfSetBlockF90_01  
  
  function     GmfSetBlockF90_02(unit, GmfKey, ad0, ad1, Tab, Ref) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64), intent(in)    :: unit
    integer(int32), intent(in)    :: GmfKey
    integer(int32), intent(in)    :: ad0
    integer(int32), intent(in)    :: ad1
    real(real64)  , intent(inout) :: Tab(:,:)
    integer(int32), intent(inout) :: Ref(  :)
    integer(int32)                :: res
    !>
    integer(int32)                :: iTab(1)
    integer(int32)                :: Nmb
    integer(int32), pointer       :: map(:)=>null()
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    Nmb=ad1-ad0+1
    
    print '("GmfSetBlockF90_02 (ad0,ad1)=(",i0,",",i0,") Nmb=",i0)',ad0,ad1,Nmb
    print '("GmfSetBlockF90_02 size(Tab)=",i0,"x",i0)',size(Tab,1),size(Tab,2)
    print '("GmfSetBlockF90_02 size(Ref)=  ",i0)',size(Ref)
    
    res=GmfSetBlockF77(unit       ,&
    &                  GmfKey     ,&
    &                  ad0        ,&
    &                  ad1        ,&
    &                  int32      ,&
    &                  map        ,&
    &                  iTab(1)    ,&
    &                  iTab(1)    ,&
    &                  Tab(1,  1) ,&
    &                  Tab(1,Nmb) ,&
    &                  Ref(    1) ,&
    &                  Ref(  Nmb)  )
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfSetBlockF90_02
  
  function     GmfSetBlockF90_02Bis(unit, GmfKey, ad0, ad1, Tab) result(res)
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    integer(int64), intent(in)    :: unit
    integer(int32), intent(in)    :: GmfKey
    integer(int32), intent(in)    :: ad0
    integer(int32), intent(in)    :: ad1
    real(real64)  , intent(inout) :: Tab(:,:)
    integer(int32)                :: res
    !>
    integer(int32)                :: Ref(1)
    integer(int32)                :: iTab(1)
    integer(int32)                :: Nmb
    integer(int32), pointer       :: map(:)=>null()
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
    Nmb=ad1-ad0+1
    
    print '("GmfSetBlockF90_02Bis (ad0,ad1)=(",i0,",",i0,") Nmb=",i0)',ad0,ad1,Nmb
    print '("GmfSetBlockF90_02Bis size(Tab)=",i0,"x",i0)',size(Tab,1),size(Tab,2)
    print '("GmfSetBlockF90_02Bis size(Ref)=  ",i0)',size(Ref)
    
    res=GmfSetBlockF77(unit       ,&
    &                  GmfKey     ,&
    &                  ad0        ,&
    &                  ad1        ,&
    &                  int32      ,&
    &                  map        ,&
    &                  iTab(   1) ,&
    &                  iTab(   1) ,&
    &                  Tab(1,  1) ,&
    &                  Tab(1,Nmb) ,&
    &                  Ref(    1) ,&
    &                  Ref(    1)  )
    !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    return
  end function GmfSetBlockF90_02Bis

end module libmeshb7