-
Notifications
You must be signed in to change notification settings - Fork 4
Expand file tree
/
Copy pathHMC_Module_Tools_IO.f90
More file actions
1626 lines (1277 loc) · 72.7 KB
/
HMC_Module_Tools_IO.f90
File metadata and controls
1626 lines (1277 loc) · 72.7 KB
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: HMC_Module_Function_IO.f90
! Author(s): Fabio Delogu, Francesco Silvestro, Simone Gabellani
! Created on May 11, 2015, 10:27 AM
!
! Module to perform IO actions
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Module Header
module HMC_Module_Tools_IO
!------------------------------------------------------------------------------------
! External module(s) and implicit none
#ifdef LIB_NC
use netcdf
#endif
use HMC_Module_Tools_Debug
use HMC_Module_Tools_Generic, only: HMC_Tools_Generic_FindFileLength
implicit none
integer(kind = 4), parameter :: iDeflateLevelNC = 9;
!------------------------------------------------------------------------------------
contains
!------------------------------------------------------------------------------------
! Subroutine for checking variable name NC
#ifdef LIB_NC
#ifdef LIB_DYNARRAY
subroutine HMC_Tools_IO_CheckVar_NC(var_name_list, file_id, var_name_select)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: file_id, var_id
character(len = *) :: var_name_list, var_name_select
character(len = 80) :: string_tmp, string_step, var_name_tmp
integer(kind = 4) :: return_code
integer(kind = 4) :: i, n, m, p, string_idx, max_length
character :: delimiter =';'
type :: string
character(:) , allocatable :: str_val ! Value of variable dynamic string
end type string
type(string), pointer :: p_string_array(:)
type(string) :: string_item
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Initializing variable(s)
var_name_select = ''
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Get sub-string number using a delimiter
string_tmp = trim (adjustl(var_name_list) )
n = count(transfer(trim(string_tmp), 'a', len(trim(string_tmp))) == delimiter) + 1
! Allocate and dump each selected string
m = 1
allocate (p_string_array(n))
do i = 1, n
string_idx = index(string_tmp(m:),delimiter)
string_step = adjustl( string_tmp(m: m + string_idx-2) )
string_item % str_val = string_step
p_string_array(i) = string_item
m = m + string_idx
end do
string_step = adjustl(string_tmp(m:) )
string_item % str_val = string_step
p_string_array(n) = string_item
! Check variable name in the file structure
do i = 1, n
var_name_tmp = p_string_array(i)%str_val
return_code = nf90_inq_varid(file_id, trim(var_name_tmp), var_id)
if (return_code .eq. 0) then
var_name_select = var_name_tmp
exit
end if
end do
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Check variable assignment
if (var_name_select .eq. '') then
call mprintf(.true., iERROR, 'Variable name not correctly defined. Check your input datasets')
else
call mprintf(.true., iINFO_Verbose, 'Select variable '//trim(var_name_select)//' in netcdf file')
endif
!------------------------------------------------------------------------------------
end subroutine HMC_Tools_IO_CheckVar_NC
#endif
#endif
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Subroutine for getting 1d variable NC
#ifdef LIB_NC
subroutine HMC_Tools_IO_Get1d_NC(sVarName, iFileID, a1dVar, sVarUnits, iX, &
bFatalError, iErr)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iFileID, iVarID
integer(kind = 4) :: iX
character(len = 256) :: sVarName
character(len = 256), intent(out) :: sVarUnits
real(kind = 4), dimension(iX), intent(out) :: a1dVar
logical, intent(in) :: bFatalError
integer(kind = 4) :: iErr, iRet
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Initializing units attribute
sVarUnits = "";
iVarID = 0;
a1dVar = -9999.0;
! Finding 2d variable ID
iRet = nf90_inq_varid(iFileID, trim(sVarName), iVarID)
if (iRet /= 0) then
if (bFatalError) then
call mprintf(.true., iERROR, ' INQUIRE FAILED :: Get1d_NC ---> nf90_inq_varid (varname: '//trim(sVarName)//')')
else
call mprintf(.true., iWARN, ' INQUIRE FAILED :: Get1d_NC ---> nf90_inq_varid (varname: '//trim(sVarName)//')')
iErr = iRet
return
endif
else
call mprintf(.true., iINFO_Extra, ' INQUIRE OK :: Get1d_NC ---> nf90_inq_varid (varname: '//trim(sVarName)//')')
endif
! Extracting 2d variable
iRet = nf90_get_var(iFileID, iVarID, a1dVar)
if (iRet /= 0) then
if (bFatalError) then
call mprintf(.true., iERROR, ' READ FAILED :: Get1d_NC ---> nf90_get_var (varname: '//trim(sVarName)//')')
else
call mprintf(.true., iWARN, ' READ FAILED :: Get1d_NC ---> nf90_get_var (varname: '//trim(sVarName)//')')
endif
else
call mprintf(.true., iINFO_Extra, ' READ OK :: Get1d_NC ---> nf90_get_var (varname: '//trim(sVarName)//')')
endif
! Algorithm conclusion
iErr = 0;
!------------------------------------------------------------------------------------
end subroutine HMC_Tools_IO_Get1d_NC
#endif
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Subroutine for getting 2d variable NC
#ifdef LIB_NC
subroutine HMC_Tools_IO_Get2d_NC(sVarName, iFileID, a2dVar, sVarUnits, iX, iY, &
bFatalError, iErr)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iFileID, iVarID
integer(kind = 4) :: iX, iY
character(len = 256) :: sVarName
character(len = 256), intent(out) :: sVarUnits
real(kind = 4), dimension(iX, iY), intent(out) :: a2dVar
logical, intent(in) :: bFatalError
integer(kind = 4) :: iErr, iRet
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Initializing units attribute
sVarUnits = "";
iVarID = 0;
a2dVar = -9999.0;
! Finding 2d variable ID
iRet = nf90_inq_varid(iFileID, trim(sVarName), iVarID)
if (iRet /= 0) then
if (bFatalError) then
call mprintf(.true., iERROR, ' INQUIRE FAILED :: Get2d_NC ---> nf90_inq_varid (varname: '//trim(sVarName)//')')
else
call mprintf(.true., iWARN, ' INQUIRE FAILED :: Get2d_NC ---> nf90_inq_varid (varname: '//trim(sVarName)//')')
iErr = iRet
return
endif
else
call mprintf(.true., iINFO_Extra, ' INQUIRE OK :: Get2d_NC ---> nf90_inq_varid (varname: '//trim(sVarName)//')')
endif
! Extracting 2d variable
iRet = nf90_get_var(iFileID, iVarID, a2dVar)
if (iRet /= 0) then
if (bFatalError) then
call mprintf(.true., iERROR, ' READ FAILED :: Get2d_NC ---> nf90_get_var (varname: '//trim(sVarName)//')')
else
call mprintf(.true., iWARN, ' READ FAILED :: Get2d_NC ---> nf90_get_var (varname: '//trim(sVarName)//')')
endif
else
call mprintf(.true., iINFO_Extra, ' READ OK :: Get2d_NC ---> nf90_get_var (varname: '//trim(sVarName)//')')
endif
! Algorithm conclusion
!iErr = 0;
iErr = iRet
!------------------------------------------------------------------------------------
end subroutine HMC_Tools_IO_Get2d_NC
#endif
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Subroutine to get 3d variable
#ifdef LIB_NC
subroutine HMC_Tools_IO_Get3d_NC(sVarName, iFileID, a3dVar, sVarUnits, dScaleFactor, iT, iX, iY, &
bFatalError, iErr)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iFileID, iVarID
integer(kind = 4) :: iX, iY, iT
character(len = 256) :: sVarName
character(len = 256), intent(out) :: sVarUnits
real(kind = 4), dimension(iX, iY, iT), intent(out) :: a3dVar
real(kind = 4), intent(out) :: dScaleFactor
logical, intent(in) :: bFatalError
integer(kind = 4) :: iErr, iRet
integer(kind = 4) :: iValidLength
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Initializing units attribute
sVarUnits = "";
iVarID = 0;
a3dVar = -9999.0;
! Finding 3d variable ID
iRet = nf90_inq_varid(iFileID, trim(sVarName), iVarID)
if (iRet /= 0) then
if (bFatalError) then
call mprintf(.true., iERROR, ' INQUIRE FAILED :: Get3d_NC ---> nf90_inq_varid (varname: '//trim(sVarName)//')')
else
call mprintf(.true., iWARN, ' INQUIRE FAILED :: Get3d_NC ---> nf90_inq_varid (varname: '//trim(sVarName)//')')
iErr = iRet
return
endif
else
call mprintf(.true., iINFO_Extra, ' INQUIRE OK :: Get3d_NC ---> nf90_inq_varid (varname: '//trim(sVarName)//')')
endif
! Extracting attributes
iRet = nf90_inquire_attribute(iFileID, iVarID, "scale_factor", len = iValidLength)
if (iRet /= 0) then
dScaleFactor = 1
call mprintf(.true., iWARN, ' READ ATTRIBUTE FAILED :: Get3d_NC ---> nf90_get_att (attrname: scale_factor)')
call mprintf(.true., iWARN, ' READ ATTRIBUTE FAILED :: Get3d_NC ---> nf90_get_att --> scale_factor = 1')
else
call check( nf90_get_att(iFileID, iVarID, 'scale_factor', dScaleFactor) )
call mprintf(.true., iINFO_Extra, ' READ ATTRIBUTE OK :: Get3d_NC ---> nf90_get_att (varname: '//trim(sVarName)//')')
endif
! Extracting 3d variable
iRet = nf90_get_var(iFileID, iVarID, a3dVar)
if (iRet /= 0) then
if (bFatalError) then
call mprintf(.true., iERROR, ' READ FAILED :: Get3d_NC ---> nf90_get_var (varname: '//trim(sVarName)//')')
else
call mprintf(.true., iWARN, ' READ FAILED :: Get3d_NC ---> nf90_get_var (varname: '//trim(sVarName)//')')
endif
else
call mprintf(.true., iINFO_Extra, ' READ OK :: Get3d_NC ---> nf90_get_var (varname: '//trim(sVarName)//')')
endif
! Algorithm conclusion
iErr = 0;
!------------------------------------------------------------------------------------
end subroutine HMC_Tools_IO_Get3d_NC
#endif
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Subroutine to put time double variable
#ifdef LIB_NC
subroutine HMC_Tools_IO_PutTime_DBL_NC(iFileID, iID_DimTime, &
sVarName, sVarNameLong, sVarDescription, &
sVarUnits, &
sVarCoords, &
iVarTime, dVarValue)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iFileID
integer(kind = 4) :: iID_DimTime
character(len = 256) :: sVarName, sVarNameLong, sVarDescription
character(len = 256) :: sVarUnits
character(len = 256) :: sVarCoords
integer(kind = 4) :: iVarID, iVarTime
real(kind = 8) :: dVarValue
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Initializing variable(s)
iVarID = 0;
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Definition mode ON - Data mode OFF
call check( nf90_redef(iFileID) )
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Allocating variable
call check( nf90_def_var(iFileID, trim(sVarName), nf90_double, (/iID_DimTime/), iVarID, &
deflate_level = iDeflateLevelNC) )
! Writing variable attribute(s)
call check( nf90_put_att(iFileID, iVarID, 'long_name' , trim(sVarNameLong)) )
call check( nf90_put_att(iFileID, iVarID, 'description' , trim(sVarDescription)) )
call check( nf90_put_att(iFileID, iVarID, 'units' , trim(sVarUnits)) )
call check( nf90_put_att(iFileID, iVarID, 'coordinates' , trim(sVarCoords)) )
! Writing variable value(s)
call check( nf90_inq_varid(iFileID, trim(sVarName), iVarID) )
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Definition mode OFF - Data mode ON
call check( nf90_enddef(iFileID))
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Writing variable(s)
call check( nf90_put_var(iFileID, iVarID, dVarValue))
!------------------------------------------------------------------------------------
end subroutine HMC_Tools_IO_PutTime_DBL_NC
#endif
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Subroutine to put time string variable
#ifdef LIB_NC
subroutine HMC_Tools_IO_PutTime_STR_NC(iFileID, iID_DimTime, &
sVarName, sVarNameLong, sVarDescription, &
sVarUnits, &
sVarCoords, &
iVarTime, sVarValue)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iFileID
integer(kind = 4) :: iID_DimTime
character(len = 256) :: sVarName, sVarNameLong, sVarDescription
character(len = 256) :: sVarUnits
character(len = 256) :: sVarCoords
integer(kind = 4) :: iVarID, iVarTime
character(len = 19) :: sVarValue
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Initializing variable(s)
iVarID = 0;
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Definition mode ON - Data mode OFF
call check( nf90_redef(iFileID) )
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Allocating variable
call check( nf90_def_var(iFileID, trim(sVarName), nf90_char, (/iID_DimTime/), iVarID, &
deflate_level = iDeflateLevelNC) )
! Writing variable attribute(s)
call check( nf90_put_att(iFileID, iVarID, 'long_name' , trim(sVarNameLong)) )
call check( nf90_put_att(iFileID, iVarID, 'description' , trim(sVarDescription)) )
call check( nf90_put_att(iFileID, iVarID, 'units' , trim(sVarUnits)) )
call check( nf90_put_att(iFileID, iVarID, 'coordinates' , trim(sVarCoords)) )
! Writing variable value(s)
call check( nf90_inq_varid(iFileID, trim(sVarName), iVarID) )
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Definition mode OFF - Data mode ON
call check( nf90_enddef(iFileID))
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Writing variable(s)
call check( nf90_put_var(iFileID, iVarID, sVarValue))
!------------------------------------------------------------------------------------
end subroutine HMC_Tools_IO_PutTime_STR_NC
#endif
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Subroutine to put 1d variable
#ifdef LIB_NC
subroutine HMC_Tools_IO_Put1d_NC(iFileID, iID_DimX, &
sVarName, sVarNameLong, sVarDescription, &
sVarUnits, &
sVarCoords, sVarGridMap, &
dVarMissingValue, &
iVarX, a1dVar)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iFileID
integer(kind = 4) :: iID_DimX
character(len = 256) :: sVarName, sVarNameLong, sVarDescription
character(len = 256) :: sVarUnits
character(len = 256) :: sVarCoords, sVarGridMap
real(kind = 4) :: dVarMissingValue
integer(kind = 4) :: iVarID, iVarX
real(kind = 4), dimension (iVarX) :: a1dVar
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Initializing variable(s)
iVarID = 0;
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Definition mode ON - Data mode OFF
call check( nf90_redef(iFileID) )
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Allocating variable
call check( nf90_def_var(iFileID, trim(sVarName), nf90_float, (/iID_DimX/), iVarID, &
deflate_level = iDeflateLevelNC) )
! Writing variable attribute(s)
call check( nf90_put_att(iFileID, iVarID, 'long_name' , trim(sVarNameLong)) )
call check( nf90_put_att(iFileID, iVarID, 'description' , trim(sVarDescription)) )
call check( nf90_put_att(iFileID, iVarID, 'units' , trim(sVarUnits)) )
call check( nf90_put_att(iFileID, iVarID, 'coordinates' , trim(sVarCoords)) )
call check( nf90_put_att(iFileID, iVarID, 'grid_mapping' , trim(sVarGridMap)) )
call check( nf90_put_att(iFileID, iVarID, 'missing_value' , real(dVarMissingValue)) )
! Writing variable value(s)
call check( nf90_inq_varid(iFileID, trim(sVarName), iVarID) )
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Definition mode OFF - Data mode ON
call check( nf90_enddef(iFileID))
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Writing variable(s)
call check( nf90_put_var(iFileID, iVarID, a1dVar, &
start = (/ 1 /), count = (/ iVarX /)) )
!------------------------------------------------------------------------------------
end subroutine HMC_Tools_IO_Put1d_NC
#endif
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Subroutine to put 2d variable
#ifdef LIB_NC
subroutine HMC_Tools_IO_Put2d_NC(iFileID, iID_DimX, iID_DimY, &
sVarName, sVarNameLong, sVarDescription, &
sVarUnits, &
sVarCoords, sVarGridMap, &
dVarMissingValue, &
iVarX, iVarY, a2dVar)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iFileID
integer(kind = 4) :: iID_DimX, iID_DimY
character(len = 256) :: sVarName, sVarNameLong, sVarDescription
character(len = 256) :: sVarUnits
character(len = 256) :: sVarCoords, sVarGridMap
real(kind = 4) :: dVarMissingValue
integer(kind = 4) :: iVarID, iVarX, iVarY
real(kind = 4), dimension (iVarX, iVarY) :: a2dVar
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Subroutine input example
! sVarName = 'SoilM'
! sVarNameLong = 'soil_moisture'
! sVarDescription = 'moisture content'
! sVarUnits = 'm^3/m^3'
! sVarCoords = 'x y z'
! sVarGridMap = 'lambert_conformal_conic'
! dVarMissingValue = -9E15
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Initializing variable(s)
iVarID = 0;
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Definition mode ON - Data mode OFF
call check( nf90_redef(iFileID) )
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Allocating variable
call check( nf90_def_var(iFileID, trim(sVarName), nf90_float, (/iID_DimX, iID_DimY/), iVarID, &
deflate_level = iDeflateLevelNC) )
! Writing variable attribute(s)
call check( nf90_put_att(iFileID, iVarID, 'long_name' , trim(sVarNameLong)) )
call check( nf90_put_att(iFileID, iVarID, 'description' , trim(sVarDescription)) )
call check( nf90_put_att(iFileID, iVarID, 'units' , trim(sVarUnits)) )
call check( nf90_put_att(iFileID, iVarID, 'coordinates' , trim(sVarCoords)) )
call check( nf90_put_att(iFileID, iVarID, 'grid_mapping' , trim(sVarGridMap)) )
call check( nf90_put_att(iFileID, iVarID, 'missing_value' , real(dVarMissingValue)) )
! Writing variable value(s)
call check( nf90_inq_varid(iFileID, trim(sVarName), iVarID) )
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Definition mode OFF - Data mode ON
call check( nf90_enddef(iFileID))
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Writing variable(s)
call check( nf90_put_var(iFileID, iVarID, a2dVar, &
start = (/ 1, 1/), count = (/ iVarX, iVarY/)) )
!------------------------------------------------------------------------------------
end subroutine HMC_Tools_IO_Put2d_NC
#endif
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Subroutine to put 3d variable
#ifdef LIB_NC
subroutine HMC_Tools_IO_Put3d_NC(iFileID, iID_DimX, iID_DimY, iID_DimT, &
sVarName, sVarNameLong, sVarDescription, &
sVarUnits, &
sVarCoords, sVarGridMap, &
dVarMissingValue, dScale_Factor, &
iVarX, iVarY, iVarT, a3dVar)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iFileID
integer(kind = 4) :: iID_DimX, iID_DimY, iID_DimT
character(len = 256) :: sVarName, sVarNameLong, sVarDescription
character(len = 256) :: sVarUnits
character(len = 256) :: sVarCoords, sVarGridMap
real(kind = 4) :: dVarMissingValue, dScale_Factor
integer(kind = 4) :: iT, iVarID, iVarX, iVarY, iVarT
real(kind = 4), dimension (iVarX, iVarY, iVarT) :: a3dVar
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Initializing variable(s)
iVarID = 0;
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Definition mode ON - Data mode OFF
call check( nf90_redef(iFileID) )
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Allocating variable
call check( nf90_def_var(iFileID, trim(sVarName), nf90_int, &
(/iID_DimX, iID_DimY, iID_DimT/), iVarID, &
shuffle = .true., deflate_level = iDeflateLevelNC) )
! Writing variable attribute(s)
call check( nf90_put_att(iFileID, iVarID, 'long_name' , trim(sVarNameLong)) )
call check( nf90_put_att(iFileID, iVarID, 'description' , trim(sVarDescription)) )
call check( nf90_put_att(iFileID, iVarID, 'units' , trim(sVarUnits)) )
call check( nf90_put_att(iFileID, iVarID, 'coordinates' , trim(sVarCoords)) )
call check( nf90_put_att(iFileID, iVarID, 'grid_mapping' , trim(sVarGridMap)) )
call check( nf90_put_att(iFileID, iVarID, 'missing_value' , real(dVarMissingValue)) )
call check( nf90_put_att(iFileID, iVarID, 'scale_factor' , real(1.0/dScale_Factor)) )
! Writing variable value(s)
call check( nf90_inq_varid(iFileID, trim(sVarName), iVarID) )
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Definition mode OFF - Data mode ON
call check( nf90_enddef(iFileID))
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Writing variable(s)
call check( nf90_put_var(iFileID, iVarID, a3dVar, &
start = (/ 1, 1, 1/), count = (/ iVarX, iVarY, iVarT/)) )
!------------------------------------------------------------------------------------
end subroutine HMC_Tools_IO_Put3d_NC
#endif
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Subroutine to get 1d variable ASCII (1 columns)
subroutine HMC_Tools_IO_Get1d_ASCII(sFileName, a1dVar, iX, &
bFatalError, iErr, iSkipRows)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iX, iI, iJ
character(len = 256) :: sFileName
real(kind = 4), dimension(iX), intent(out) :: a1dVar
logical, intent(in) :: bFatalError
logical :: bFileExist
integer(kind = 4) :: iErr, iRet
integer(kind = 4), optional, intent(in) :: iSkipRows
integer(kind = 4) :: iSkipRows_Set, iFileLength
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Check subroutine optional argument(s)
iSkipRows_Set = 0
if (present(iSkipRows)) then
iSkipRows_Set = iSkipRows
endif
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Open and read ascii file
inquire (file = trim(sFileName), exist = bFileExist, iostat = iRet)
if (.not. bFileExist ) then
if (bFatalError) then
call mprintf(.true., iERROR, &
' INQUIRE FAILED :: Get1d_ASCII ---> file not found (filename: '//trim(sFileName)//')')
else
call mprintf(.true., iWARN, &
' INQUIRE FAILED :: Get1d_ASCII ---> file not found (filename: '//trim(sFileName)//')')
iErr = -1
return
endif
elseif ( bFileExist ) then
! Open file
open(unit=1, file=trim(sFileName), status='old', iostat=iRet)
! Check length of data
call HMC_Tools_Generic_FindFileLength(1, iFileLength, iSkipRows_Set)
if (iFileLength .lt. iX) then
call mprintf(.true., iERROR, &
' INQUIRE FAILED :: Get1d_ASCII ---> file length < expected length (filename: '//trim(sFileName)//')')
elseif (iFileLength .gt. iX) then
call mprintf(.true., iWARN, &
' READ OK :: Get1d_ASCII ---> file length > expected length (filename: '//trim(sFileName)//')')
endif
if (iRet == 0) then
call mprintf(.true., iINFO_Extra, &
' READ OK :: Get1d_ASCII ---> file readable (filename: '//trim(sFileName)//')')
do iI = 1,iX
read(1,*) a1dVar(iI)
enddo
iErr = iRet
close(1)
elseif (iRet /= 0) then
if (bFatalError) then
call mprintf(.true., iERROR, &
' READ FAILED :: Get1d_ASCII ---> file not readable (filename: '//trim(sFileName)//')')
else
call mprintf(.true., iWARN, &
' READ FAILED :: Get1d_ASCII ---> file not readable (filename: '//trim(sFileName)//')')
iErr = iRet
return
endif
endif
endif
!------------------------------------------------------------------------------------
end subroutine HMC_Tools_IO_Get1d_ASCII
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Subroutine to get Nd variable ASCII (N columns)
subroutine HMC_Tools_IO_GetNd_ASCII(sFileName, a2dVar, iX, iY, bFatalError, iErr)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iX, iY, iI, iJ
character(len = 256) :: sFileName
real(kind = 4), dimension(iX, iY), intent(out) :: a2dVar
logical, intent(in) :: bFatalError
logical :: bFileExist
integer(kind = 4) :: iErr, iRet
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Initialize variable(s)
a2dVar = -9999.0
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Open and read ascii file
inquire (file = trim(sFileName), exist = bFileExist, iostat = iRet)
if (.not. bFileExist ) then
if (bFatalError) then
call mprintf(.true., iERROR, &
' INQUIRE FAILED :: GetNd_ASCII ---> file not found (filename: '//trim(sFileName)//')')
else
call mprintf(.true., iWARN, &
' INQUIRE FAILED :: GetNd_ASCII ---> file not found (filename: '//trim(sFileName)//')')
iErr = -1
return
endif
elseif ( bFileExist ) then
open(unit=1, file=trim(sFileName), status='old', iostat=iRet)
call mprintf(.true., iINFO_Extra, &
' READ OK :: GetNd_ASCII ---> file readable (filename: '//trim(sFileName)//')')
if (iRet == 0) then
do iI = 1, iX
read(1,*) a2dVar(iI, :)
enddo
close(1)
elseif (iRet /= 0) then
if (bFatalError) then
call mprintf(.true., iERROR, &
' READ FAILED :: GetNd_ASCII ---> file not readable (filename: '//trim(sFileName)//')')
else
call mprintf(.true., iWARN, &
' READ FAILED :: GetNd_ASCII ---> file not readable (filename: '//trim(sFileName)//')')
iErr = iRet
return
endif
endif
endif
!------------------------------------------------------------------------------------
end subroutine HMC_Tools_IO_GetNd_ASCII
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Subroutine to get 2d variable ASCII (2 columns)
subroutine HMC_Tools_IO_Get2d_ASCII(sFileName, a1sVar, a1dVar, iX, bFatalError, iErr)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iX, iI, iJ
character(len = 256) :: sFileName
real(kind = 4), dimension(iX), intent(out) :: a1dVar
character(len = 256), dimension(iX) :: a1sVar
logical, intent(in) :: bFatalError
logical :: bFileExist
integer(kind = 4) :: iErr, iRet
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Initialize variable(s)
a1sVar = ''; a1dVar = -9999.0
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Open and read ascii file
inquire (file = trim(sFileName), exist = bFileExist, iostat = iRet)
if (.not. bFileExist ) then
if (bFatalError) then
call mprintf(.true., iERROR, &
' INQUIRE FAILED :: Get2d_ASCII ---> file not found (filename: '//trim(sFileName)//')')
else
call mprintf(.true., iWARN, &
' INQUIRE FAILED :: Get2d_ASCII ---> file not found (filename: '//trim(sFileName)//')')
iErr = -1
return
endif
elseif ( bFileExist ) then
open(unit=1, file=trim(sFileName), status='old', iostat=iRet)
call mprintf(.true., iINFO_Extra, &
' READ OK :: Get2d_ASCII ---> file readable (filename: '//trim(sFileName)//')')
if (iRet == 0) then
do iI = 1, iX
read(1,*) a1sVar(iI), a1dVar(iI)
enddo
close(1)
elseif (iRet /= 0) then
if (bFatalError) then
call mprintf(.true., iERROR, &
' READ FAILED :: Get2d_ASCII ---> file not readable (filename: '//trim(sFileName)//')')
else
call mprintf(.true., iWARN, &
' READ FAILED :: Get2d_ASCII ---> file not readable (filename: '//trim(sFileName)//')')
iErr = iRet
return
endif
endif
endif
!------------------------------------------------------------------------------------
end subroutine HMC_Tools_IO_Get2d_ASCII
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Subroutine to get ArcGrid variable ASCII
subroutine HMC_Tools_IO_GetArcGrid_ASCII(sFileName, a2dVar, iX, iY, &
bFatalError, iErr)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iX, iY, iI, iJ
character(len = 256) :: sFileName
!real(kind = 4), dimension(iX, iY), intent(out) :: a2dVar
real(kind = 4), dimension(iY, iX), intent(out) :: a2dVar
logical, intent(in) :: bFatalError
logical :: bFileExist
integer(kind = 4) :: iErr, iRet
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Init variable(s)
iErr = 0; a2dVar = 0.0
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Open and read ascii file
inquire (file = trim(sFileName), exist = bFileExist, iostat = iRet)
if (.not. bFileExist ) then
if (bFatalError) then
call mprintf(.true., iERROR, &
' INQUIRE FAILED :: GetArcGrid ---> file not found (filename: '//trim(sFileName)//')')
else
call mprintf(.true., iWARN, &
' INQUIRE FAILED :: GetArcGrid ---> file not found (filename: '//trim(sFileName)//')')
iErr = -1
return
endif
elseif ( bFileExist ) then
open(unit=1, file=trim(sFileName), status='old', iostat=iRet)
if (iRet == 0) then
call mprintf(.true., iINFO_Extra, &
' READ OK :: GetArcGrid ---> file readable (filename: '//trim(sFileName)//')')
do iI = 1,6; read(1,*); enddo
!do iI = 1,iX,1
! read(1,*) (a2dVar(iX - iI + 1,iJ), iJ = 1,iY)
!enddo
do iI = 1,iY,1
read(1,*) (a2dVar(iY - iI + 1,iJ), iJ = 1,iX)
enddo
close(1)
iErr = 0
elseif (iRet /= 0) then
if (bFatalError) then
call mprintf(.true., iERROR, &
' READ FAILED :: GetArcGrid ---> file not readable (filename: '//trim(sFileName)//')')
else
call mprintf(.true., iWARN, &
' READ FAILED :: GetArcGrid ---> file not readable (filename: '//trim(sFileName)//')')
iErr = iRet
return
endif
endif
endif
!------------------------------------------------------------------------------------
end subroutine HMC_Tools_IO_GetArcGrid_ASCII
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Subroutine to put 1d variable ASCII
subroutine HMC_Tools_IO_Put1d_ASCII(sFileName, a1dVar, iX, bFatalError, iErr, sFMT)
!------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iX, iI
character(len = 256) :: sFileName
real(kind = 4), dimension(iX) :: a1dVar
logical, intent(in) :: bFatalError
logical :: bFileExist
integer(kind = 4) :: iErr, iRet
character(len = 20) :: sFMT
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Open and write ASCII file
open(unit = 20, file = trim(sFileName), status = 'new', iostat = iRet)
inquire (file = trim(sFileName), exist = bFileExist, iostat = iRet)
if (.not. bFileExist ) then
if (bFatalError) then
call mprintf(.true., iERROR, &
' INQUIRE FAILED :: Put1d_ASCII ---> file not readable (filename: '//trim(sFileName)//')')
else
call mprintf(.true., iWARN, &
' INQUIRE FAILED :: Put1d_ASCII ---> file not readable (filename: '//trim(sFileName)//')')
iErr = 0
return
endif
elseif ( bFileExist ) then
if (iRet == 0) then
call mprintf(.true., iINFO_Extra, &
' WRITE OK :: Put1d_ASCII ---> file writable (filename: '//trim(sFileName)//')')
do iI = 1,iX
write(20 , sFMT) a1dVar(iI)
enddo
close(20)
elseif (iRet /= 0) then
if (bFatalError) then
call mprintf(.true., iERROR, &
' WRITE FAILED :: Put1d_ASCII ---> file not writable (filename: '//trim(sFileName)//')')
else
call mprintf(.true., iWARN, &
' WRITE FAILED :: Put1d_ASCII ---> file not writable (filename: '//trim(sFileName)//')')
iErr = iRet
return
endif
endif
endif
!------------------------------------------------------------------------------------