-
Notifications
You must be signed in to change notification settings - Fork 144
/
state_structure_mod.f90
2493 lines (1729 loc) · 81 KB
/
state_structure_mod.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
! DART software - Copyright UCAR. This open source software is provided
! by UCAR, "as is", without charge, subject to all terms of use at
! http://www.image.ucar.edu/DAReS/DART/DART_download
!-------------------------------------------------------------------------------
module state_structure_mod
!> \defgroup state_structure_mod state_structure_mod
!> @{ @brief Contains the stucture of the state vector and routines to access this structure
!>
!> Usage:
!> call add_domain() in static_init_model()
!> - must be called.
!> - may be called more than once in static_init_model if a model has
!> more than one domain (e.g. WRF).
!>
!> The add_domain() call adds a 'domain' to the state. This may be a component in
!> the case of XCESM or another coupled model.
!>
!> There are three ways to add a domain (these are overloaded as add_domain):
!> * add_domain_blank. This takes model size as an argument.
!>
!> * add_domain_from_file. This takes a netcdf file and a list of variables
!>
!> * add_domain_from_spec. This makes a skeleton structure for a domain. Dimensions
!> for each variable must be added using add_dimension_to_variable(). This is intended
!> to be used to create netcdf output for models like bgrid_solo that are spun up.
!> Usage: call add_domain_from_spec
!> for each variable call add_dimension_to_variable() for each dimension in the variable
!> call finished_adding_domain_from_spec
!>
!> There are optional arguments to add_domain_from_spec and add_domain_from_file. These
!> are:
!> * kind_list - dart kind of each variable
!> * clamp_vals - upper and lower bounds each variable (missing_r8 for no clamping)
!> * update_list - true/false. Write out the variable. Default is .true.
!>
!> The state_structure hierarchy is:
!> state -> domain -> variable
!> Each time add_domain is called another domain is added to the state.
!> Outside code can query with the state structure using domain_id and variable_id as
!> arguments to accessor funcitons.
!> Domain_id and variable_id are integers currently. It may be better to make domain_id
!> a private type so other information can be stored in the domain_id.
!>
!> Interaction with the state_stucture depends on which module is accessing the structure.
!> For example, model_mod querying the number of dimensions a variable is not concerned
!> with the unlimited dimension, but state_vector_io_mod is concerned with the unlimtied
!> dimension.
!> The variable type has an io_type inside it to hold io information
!> The io accessor functions have get_io_* in their name.
!>
!> get_dart_vector_index() and its inverse get_model_variable_indices() link
!> a model x,y,z to dart index. The order of the state vector is no longer under model_mod
!> control. Beware when converting model_mods such as CAM that transform the order of
!> variables after reading from a netcdf file. There may be calculations in model_mod
!> that are assuming a transformed order which no longer exists.
use utilities_mod, only : E_ERR, error_handler, do_output
use obs_kind_mod, only : get_name_for_quantity, get_index_for_quantity
use types_mod, only : r8, r4, i8, digits12, MISSING_R8, MISSING_R4, MISSING_I, &
obstypelength, MAX_NUM_DOMS
use netcdf_utilities_mod, only : nc_check
use sort_mod, only : index_sort
use netcdf
implicit none
private
public :: add_domain, &
get_domain_size, &
get_num_domains, &
get_variable_size, &
get_variable_name, &
get_kind_string, &
get_kind_index, &
get_varid_from_varname, &
get_varid_from_kind, &
get_varids_from_kind, &
get_num_variables, &
get_num_dims, &
get_dim_lengths, &
get_dim_length, &
get_dim_name, &
get_io_num_dims, &
get_io_dim_ids, &
get_io_dim_lengths, &
get_io_num_unique_dims, &
get_io_unique_dim_name, &
get_io_unique_dim_length, &
get_io_dim_length, &
add_time_unlimited, &
get_unlimited_dimid, &
set_var_id, &
get_io_clamping_maxval, &
get_io_clamping_minval, &
do_io_clamping, &
do_io_update, &
get_index_start, &
get_index_end, &
get_sum_variables, &
get_sum_variables_below, &
get_model_variable_indices, &
get_dart_vector_index, &
get_num_varids_from_kind, &
get_xtype, &
get_units, &
get_long_name, &
get_short_name, &
get_has_missing_value, &
get_FillValue, &
get_missing_value, &
get_has_FillValue, &
get_add_offset, &
get_scale_factor, &
set_dart_kinds, &
set_clamping, &
set_update_list, &
add_dimension_to_variable, &
finished_adding_domain, &
state_structure_info, &
hyperslice_domain, &
has_unlimited_dim
! diagnostic files
!>@todo FIXME these routines are deprecated because we are no supporting 'diagnostic'
!> files, but they will likely be useful for the single file (multiple member) input.
public :: create_diagnostic_structure, &
end_diagnostic_structure
character(len=*), parameter :: source = 'state_structure_mod.f90'
character(len=512) :: string1, string2, string3
!-------------------------------------------------------------------------------
! global variables
!-------------------------------------------------------------------------------
integer, parameter :: diagnostic_domain = MAX_NUM_DOMS + 1 ! Need to separate this from state
!-------------------------------------------------------------------------------
! Describes information pertaining to the IO portion of the variable information
!-------------------------------------------------------------------------------
type io_information
private
! netcdf variable id
integer :: varid ! HK Do you want one for reading, one for writing?
integer :: xtype ! netCDF variable type (NF90_double, etc.)
! clamping variables
logical :: clamping = .false. ! does variable need to be range-restricted before
real(r8) :: minvalue = missing_r8 ! min value for clamping
real(r8) :: maxvalue = missing_r8 ! max value for clamping
! logical :: out_of_range_fail = .false. ! is out of range fatal if range-checking?
! dimension information, including unlimited dimensions
integer :: io_numdims = 0
integer, dimension(NF90_MAX_VAR_DIMS) :: io_dimIds
integer, dimension(NF90_MAX_VAR_DIMS) :: dimlens
! update information
logical :: update = .true. ! default to update variables
! CF-Conventions
character(len=NF90_MAX_NAME) :: units = ' '
character(len=NF90_MAX_NAME) :: short_name = ' '
character(len=NF90_MAX_NAME) :: long_name = ' '
logical :: has_missing_value = .false.
logical :: has_FillValue = .false.
integer :: missingINT = MISSING_I ! missing values
real(r4) :: missingR4 = MISSING_R4
real(r8) :: missingR8 = MISSING_R8
integer :: spvalINT = MISSING_I ! fill values
real(r4) :: spvalR4 = MISSING_R4
real(r8) :: spvalR8 = MISSING_R8
real(r8) :: scale_factor = MISSING_R8
real(r8) :: add_offset = MISSING_R8
end type io_information
!-------------------------------------------------------------------------------
! state_type -> domain_type -> variable_type
! #domains #variables info
!-------------------------------------------------------------------------------
! Describes a variable, enough for netcdf io purposes to read and write
!-------------------------------------------------------------------------------
type variable_type
private
! variable information
character(len=NF90_MAX_NAME) :: varname
integer :: var_size = 0
!>@todo FIXME : should we have a missingI8
integer(i8) :: index_start = -1_i8 ! first occurance of variable in state
integer(i8) :: index_end = -1_i8 ! last occurance of variable in state
! dimension information
integer :: numdims = 0 ! number of dims - excluding TIME?
logical :: var_has_unlim = .false.
character(len=NF90_MAX_NAME), dimension(NF90_MAX_VAR_DIMS) :: dimname
integer, dimension(NF90_MAX_VAR_DIMS) :: dimlens
! dart information
integer :: dart_kind = -1
character(len=obstypelength) :: kind_string = 'unknown'
type(io_information) :: io_info
end type variable_type
!-------------------------------------------------------------------------------
! Describes a domain
!-------------------------------------------------------------------------------
type domain_type
private
! variables in domain
integer :: num_variables
type(variable_type), allocatable :: variable(:)
! dimension information for domain
integer :: num_unique_dims
character(len=NF90_MAX_NAME), allocatable :: unique_dim_names(:)
integer, allocatable :: unique_dim_length(:)
integer, allocatable :: original_dim_IDs(:)
integer :: unlimDimId = -1 ! initialize to no unlimited dimension
logical :: has_unlimited = .false.
! number of elements in the domain
integer(i8) :: dom_size = 0_i8
! netcdf file describing the shape of the variable
character(len=256) :: info_file = 'NULL'
! string identifying the manner in which the domain was created
! 'blank', 'file', or 'spec'
character(len=6) :: method = 'none'
end type domain_type
!-------------------------------------------------------------------------------
! Describes a state vector
!-------------------------------------------------------------------------------
type state_type
private
! domains or separate files
integer :: num_domains = 0
! number of domains + 1 for diagnostic domain
type(domain_type) :: domain(MAX_NUM_DOMS + 1)
! size of the state vector
integer(i8) :: model_size = 0_i8
!>@todo is this different to character(len=128) :: name? Redundant?
!character(len=1024) :: domain_name(MAX_NUM_DOMS)
end type state_type
!-------------------------------------------------------------------------------
! module handle
!-------------------------------------------------------------------------------
type(state_type) :: state
! This is used when dart writes the state to diagnostic files
! It is a transformation of the state_type
! One diagnostic file contains all domains.
! It is not part of the state type
logical :: diagnostic_initialized = .false.
! debug flag for get_state_indices
logical :: debug = .false.
interface add_domain
module procedure add_domain_blank
module procedure add_domain_from_file
module procedure add_domain_from_spec
end interface
interface get_index_start
module procedure get_index_start_from_varname
module procedure get_index_start_from_varid
end interface
interface get_index_end
module procedure get_index_end_from_varname
module procedure get_index_end_from_varid
end interface
interface get_missing_value
module procedure get_missing_value_r8
module procedure get_missing_value_r4
module procedure get_missing_value_int
end interface
interface get_FillValue
module procedure get_spval_r8
module procedure get_spval_r4
module procedure get_spval_int
end interface
contains
!-------------------------------------------------------------------------------
!> Given an info_file, reads in a list var_names(num_vars)
!> into the state_strucutre.
!>
!> Returns a dom_id that can be used to harvest information of a particular domain
!>
!> Does this need to be a function or a subroutine?
function add_domain_from_file(info_file, num_vars, var_names, kind_list, clamp_vals, update_list) result(dom_id)
character(len=*), intent(in) :: info_file
integer, intent(in) :: num_vars
character(len=*), intent(in) :: var_names(num_vars)
integer, intent(in), optional :: kind_list(num_vars)
real(r8), intent(in), optional :: clamp_vals(num_vars, 2)
logical, intent(in), optional :: update_list(num_vars)
integer :: dom_id
integer :: ivar
! add to domains
call assert_below_max_num_domains('add_domain_from_file')
state%num_domains = state%num_domains + 1
!>@todo dom_id should be a handle.
dom_id = state%num_domains
! save information about the information file
state%domain(dom_id)%info_file = info_file
state%domain(dom_id)%method = 'file'
! set number of variables in this domain
state%domain(dom_id)%num_variables = num_vars
! load up the variable names
allocate(state%domain(dom_id)%variable(num_vars))
do ivar = 1, num_vars
state%domain(dom_id)%variable(ivar)%varname = var_names(ivar)
enddo
! load up variable id's and sizes
call load_state_variable_info(state%domain(dom_id),dom_id)
! load up the domain unique dimension info
call load_unique_dim_info(dom_id)
! load up any cf-conventions if they exist
call load_common_cf_conventions(state%domain(dom_id))
if ( present(kind_list) ) call set_dart_kinds (dom_id, num_vars, kind_list)
if ( present(clamp_vals) ) call set_clamping (dom_id, num_vars, clamp_vals)
if ( present(update_list) ) call set_update_list(dom_id, num_vars, update_list)
end function add_domain_from_file
!-------------------------------------------------------------------------------
!> Defines a skeleton structure for the state structure. Dimension can be
!> added to variables with add_dimension_to_variable.
!>
!> Returns a dom_id that can be used to harvest information of a particular domain
function add_domain_from_spec(num_vars, var_names, kind_list, clamp_vals, update_list) result(dom_id)
integer, intent(in) :: num_vars
character(len=*), intent(in) :: var_names(num_vars)
integer, intent(in), optional :: kind_list(num_vars)
real(r8), intent(in), optional :: clamp_vals(num_vars, 2)
logical, intent(in), optional :: update_list(num_vars)
integer :: dom_id
integer :: ivar
! add to domains
call assert_below_max_num_domains('add_domain_from_spec')
state%num_domains = state%num_domains + 1
dom_id = state%num_domains
state%domain(dom_id)%method = 'spec'
! set number of variables in this domain
state%domain(dom_id)%num_variables = num_vars
state%domain(dom_id)%num_unique_dims = 0
! load up the variable names
allocate(state%domain(dom_id)%variable(num_vars))
do ivar = 1, num_vars
state%domain(dom_id)%variable(ivar)%varname = var_names(ivar)
state%domain(dom_id)%variable(ivar)%io_info%xtype = NF90_DOUBLE
enddo
if ( present(kind_list) ) call set_dart_kinds (dom_id, num_vars, kind_list)
if ( present(clamp_vals) ) call set_clamping (dom_id, num_vars, clamp_vals)
if ( present(update_list) ) call set_update_list(dom_id, num_vars, update_list)
end function add_domain_from_spec
!-------------------------------------------------------------------------------
!> Add a blank domain - one variable called state, length = domain_size
! HK the above comment is not true, there are three dimensions created in this function.
! HK this should set has_unlimited = .true.
function add_domain_blank(domain_size) result(dom_id)
integer(i8), intent(in) :: domain_size
integer :: dom_id
integer(i8) :: domain_offset
! add to domains
call assert_below_max_num_domains('add_domain_blank')
state%num_domains = state%num_domains + 1
dom_id = state%num_domains
if (state%num_domains > 1 ) then
domain_offset = get_index_end(dom_id-1,get_num_variables(dom_id-1))
else
domain_offset = 0
endif
! domain
state%domain(dom_id)%method = 'blank'
state%domain(dom_id)%num_variables = 1
state%domain(dom_id)%dom_size = domain_size
state%model_size = state%model_size + domain_size
state%domain(dom_id)%num_unique_dims = 3
allocate(state%domain(dom_id)%original_dim_IDs(3))
allocate(state%domain(dom_id)%unique_dim_names(3))
allocate(state%domain(dom_id)%unique_dim_length(3))
state%domain(dom_id)%unique_dim_names(1) = 'location'
state%domain(dom_id)%unique_dim_names(2) = 'member'
state%domain(dom_id)%unique_dim_names(3) = 'time'
state%domain(dom_id)%unique_dim_length(1) = domain_size
state%domain(dom_id)%unique_dim_length(2) = 1
state%domain(dom_id)%unique_dim_length(3) = 1
state%domain(dom_id)%original_dim_IDs(1) = 1
state%domain(dom_id)%original_dim_IDs(2) = 2
state%domain(dom_id)%original_dim_IDs(3) = NF90_UNLIMITED
! variable
allocate(state%domain(dom_id)%variable(1))
state%domain(dom_id)%variable(1)%varname = 'state'
state%domain(dom_id)%variable(1)%numdims = 1
state%domain(dom_id)%variable(1)%var_size = domain_size
state%domain(dom_id)%variable(1)%index_start = domain_offset + 1
state%domain(dom_id)%variable(1)%index_end = domain_offset + domain_size
!>@todo FIXME : what is a good default for kind_string
state%domain(dom_id)%variable(1)%kind_string = 'QTY_STATE_VARIABLE'
state%domain(dom_id)%variable(1)%dart_kind = &
get_index_for_quantity(state%domain(dom_id)%variable(1)%kind_string)
state%domain(dom_id)%variable(1)%dimname(1) = 'location'
state%domain(dom_id)%variable(1)%dimname(2) = 'member'
state%domain(dom_id)%variable(1)%dimname(3) = 'time'
state%domain(dom_id)%variable(1)%dimlens(1) = domain_size
state%domain(dom_id)%variable(1)%dimlens(2) = 1
state%domain(dom_id)%variable(1)%dimlens(3) = 1
state%domain(dom_id)%variable(1)%io_info%dimlens(1) = domain_size
state%domain(dom_id)%variable(1)%io_info%dimlens(2) = 1
state%domain(dom_id)%variable(1)%io_info%dimlens(3) = 1
state%domain(dom_id)%variable(1)%io_info%xtype = NF90_DOUBLE
state%domain(dom_id)%variable(1)%io_info%units = 'none'
state%domain(dom_id)%variable(1)%io_info%io_numdims = 3
state%domain(dom_id)%variable(1)%io_info%io_dimids(1) = 1
state%domain(dom_id)%variable(1)%io_info%io_dimids(2) = 2
state%domain(dom_id)%variable(1)%io_info%io_dimids(3) = NF90_UNLIMITED
end function add_domain_blank
!-------------------------------------------------------------------------------
!> Load metadata from netcdf file info state_strucutre
subroutine load_state_variable_info(domain, domain_index)
type(domain_type), intent(inout) :: domain
integer, intent(in) :: domain_index
! netcdf variables
integer :: ret, ncfile
character(len=256) :: nc_filename
nc_filename = domain%info_file
! open netcdf file - all restart files in a domain have the same info?
ret = nf90_open(nc_filename, NF90_NOWRITE, ncfile)
call nc_check(ret, 'load_state_variable_info nf90_open', trim(nc_filename))
! get the dimension id of the unlimited dimension if it exists
ret = nf90_inquire(ncfile, unlimitedDimId=domain%unlimDimId)
call nc_check(ret, 'load_state_variable_info, nf90_inquire')
if ( domain%unlimDimID /= -1 ) domain%has_unlimited = .true.
! get variable ids
call load_variable_ids(ncfile, domain, domain_index)
! get all variable sizes, only readers store dimensions?
call load_variable_sizes(ncfile, domain)
! close netcdf file
ret = nf90_close(ncfile)
call nc_check(ret, 'load_state_variable_info nf90_close', trim(nc_filename))
end subroutine load_state_variable_info
!-------------------------------------------------------------------------------
!> Load netcdf variable ids
subroutine load_variable_ids(ncfile, domain, domain_index)
integer, intent(in) :: ncfile ! netdcf file id - should this be part of the domain handle?
type(domain_type), intent(inout) :: domain
integer, intent(in) :: domain_index
integer :: ret ! netcdf return value
integer :: ivar, num_vars
num_vars = domain%num_variables
do ivar = 1, num_vars
! load netcdf id from variable name
ret = nf90_inq_varid(ncfile, domain%variable(ivar)%varname, &
domain%variable(ivar)%io_info%varid)
write(string1,*)'domain ',domain_index,', variable #',ivar,' "', &
trim(domain%variable(ivar)%varname)//'" from file "'//trim(domain%info_file)//'"'
call nc_check(ret, 'load_variable_ids, nf90_inq_var_id', string1)
enddo
end subroutine load_variable_ids
!-------------------------------------------------------------------------------
!> load dimension information and calculate variable and domain sizes
!-------------------------------------------------------------------------------
subroutine load_variable_sizes(ncfile, domain)
integer, intent(in) :: ncfile ! netdcf file id - should this be part of the domain handle?
type(domain_type), intent(inout) :: domain
integer :: ivar, jdim, num_vars, num_dims !< loop variables
integer :: variable_size
integer :: ret ! netcdf return value
integer(i8) :: index_start, domain_size
domain_size = 0
num_vars = domain%num_variables
index_start = state%model_size + 1
do ivar = 1, num_vars
! from netcdf id load variable dimension and ids
ret = nf90_inquire_variable(ncfile, domain%variable(ivar)%io_info%varid, &
ndims = domain%variable(ivar)%io_info%io_numdims, &
dimids = domain%variable(ivar)%io_info%io_dimIds, &
xtype = domain%variable(ivar)%io_info%xtype)
call nc_check(ret, 'load_variable_sizes, inq_variable', &
trim(domain%variable(ivar)%varname))
variable_size = 1
num_dims = domain%variable(ivar)%io_info%io_numdims
do jdim = 1, num_dims
! load dimension names and lengths
ret = nf90_inquire_dimension(ncfile, domain%variable(ivar)%io_info%io_dimIds(jdim), &
name = domain%variable(ivar)%dimname(jdim), &
len = domain%variable(ivar)%io_info%dimlens(jdim))
domain%variable(ivar)%dimlens(jdim) = domain%variable(ivar)%io_info%dimlens(jdim)
call nc_check(ret, 'load_variable_sizes, inq_dimension', &
trim(domain%variable(ivar)%dimname(jdim)))
!>@todo FIXME we'll have to document that no user can have a 'member' dimension
!>in their own netcdf files. it's more reasonable to indicate 'time' is special
!>but this code doesn't know if its a user file or one we wrote. but without
!>skipping these dimensions you get the wrong variable size because it includes
!>all times and all members in the size. this needs to be revisited later.
if ((domain%variable(ivar)%dimname(jdim) == 'time') .or. &
(domain%variable(ivar)%dimname(jdim) == 'member')) cycle
variable_size = variable_size * domain%variable(ivar)%dimlens(jdim)
enddo
! to be consistent this needs to ignore both 'time' and 'member' for
! files we write, and newer netcdf libs support multiple unlimited dims
! so this does need to change. but it works, somehow, as-is so leave it
! for now.
! subtract the unlimited dimension if it exists
!>@todo : how to handle models with multiple unlimited dimensions?
!> nf90_inquire returns the first unlimited dimension id which
!> is the slowest varying dimension. For now am assuming that
!> there can only be one unlimited dimension. Just subtract
!> to get 'spatial' dimensions.
if ( any(domain%variable(ivar)%io_info%io_dimIds(1:num_dims) == domain%unlimDimId) ) then
domain%variable(ivar)%numdims = num_dims - 1
domain%variable(ivar)%var_has_unlim = .TRUE.
else
domain%variable(ivar)%numdims = num_dims
endif
! member is not a spatial domain but could be included in a single file
do jdim = 1, num_dims
if ( domain%variable(ivar)%dimname(jdim) == 'member') then
domain%variable(ivar)%numdims = domain%variable(ivar)%numdims - 1
endif
enddo
domain%variable(ivar)%var_size = variable_size
! first and last location of variable in the state index
domain%variable(ivar)%index_start = index_start
domain%variable(ivar)%index_end = index_start + variable_size - 1
! update counters
domain_size = domain_size + variable_size
index_start = domain%variable(ivar)%index_end + 1
enddo
domain%dom_size = domain_size
state%model_size = state%model_size + domain_size
end subroutine load_variable_sizes
!-------------------------------------------------------------------------------
!> Identify the unique dimensions within a domain
!> I think you can sort dimension ids and remove duplicates
!> What if there are 0 dimensions?
!-------------------------------------------------------------------------------
subroutine load_unique_dim_info(dom_id)
integer, intent(in) :: dom_id ! domain identifier
integer, allocatable :: array_of_dimids(:)
integer, allocatable :: array_of_lengths(:)
integer, allocatable :: array_of_indices(:)
character(len=NF90_MAX_NAME), allocatable :: array_of_names(:)
logical, allocatable :: unique(:)
integer :: ivar, jdim, jdim_dom
integer :: count_dims, ndims, nvars, ndims_dom
integer :: num_unique
num_unique = 0
ndims_dom = get_domain_num_dims(dom_id)
nvars = get_num_variables(dom_id)
allocate(array_of_dimids(ndims_dom))
allocate(array_of_names(ndims_dom))
allocate(array_of_lengths(ndims_dom))
allocate(array_of_indices(ndims_dom))
allocate(unique(ndims_dom))
count_dims = 1
do ivar = 1, nvars
ndims = get_io_num_dims(dom_id, ivar)
do jdim = 1, ndims
array_of_dimids(count_dims) = state%domain(dom_id)%variable(ivar)%io_info%io_dimIds(jdim)
array_of_names(count_dims) = state%domain(dom_id)%variable(ivar)%dimname(jdim)
array_of_lengths(count_dims) = state%domain(dom_id)%variable(ivar)%dimlens(jdim)
count_dims = count_dims + 1
enddo
enddo
call index_sort(array_of_dimids, array_of_indices, ndims_dom)
count_dims = 1
unique(:) = .false.
unique(1) = .true.
do jdim_dom = 2, ndims_dom
if( array_of_dimids(array_of_indices(jdim_dom)) /= &
array_of_dimids(array_of_indices(jdim_dom-1)) ) then
count_dims = count_dims + 1
unique(jdim_dom) = .true.
endif
enddo
state%domain(dom_id)%num_unique_dims = count_dims
allocate(state%domain(dom_id)%unique_dim_names(count_dims))
allocate(state%domain(dom_id)%unique_dim_length(count_dims))
allocate(state%domain(dom_id)%original_dim_IDs(count_dims))
count_dims = 1
do jdim_dom = 1, ndims_dom
if(unique(jdim_dom)) then
state%domain(dom_id)%unique_dim_names( count_dims) = array_of_names( array_of_indices(jdim_dom))
state%domain(dom_id)%unique_dim_length(count_dims) = array_of_lengths(array_of_indices(jdim_dom))
state%domain(dom_id)%original_dim_IDs( count_dims) = array_of_dimids( array_of_indices(jdim_dom))
count_dims = count_dims + 1
endif
enddo
deallocate(array_of_dimids, array_of_names, array_of_lengths, array_of_indices, unique)
end subroutine load_unique_dim_info
!-------------------------------------------------------------------------------
!> Check to see if the template file has any of the common cf-conventions
!>
!> * units
!> * short_name
!> * long_name
!> * short_name
!> * _FillValue
!> * missing_value
!> * add_offset
!> * scale_factor
!>
!> If they exist, load them up into the state structure.
!-------------------------------------------------------------------------------
subroutine load_common_cf_conventions(domain)
type(domain_type), intent(inout) :: domain
integer :: ivar
integer :: nvars
! netcdf variables
integer :: ret, ncid, VarID
integer :: var_xtype
integer :: cf_spvalINT
real(r4) :: cf_spvalR4
real(digits12) :: cf_spvalR8
real(digits12) :: cf_scale_factor, cf_add_offset
character(len=512) :: ncFilename
character(len=NF90_MAX_NAME) :: var_name
character(len=NF90_MAX_NAME) :: cf_long_name, cf_short_name, cf_units
character(len=*), parameter :: routine = 'load_common_cf_conventions'
ncFilename = domain%info_file
ret = nf90_open(ncFilename, NF90_NOWRITE, ncid)
call nc_check(ret, routine,'nf90_open '//trim(ncFilename))
! determine attributes of each variable in turn
nvars = domain%num_variables
do ivar = 1, nvars
var_name = domain%variable(ivar)%varname
ret = nf90_inq_varid(ncid, trim(var_name), VarID)
call nc_check(ret, routine, 'inq_varid '//trim(var_name))
! If the short_name, long_name and/or units attributes are set, get them.
! They are not REQUIRED by DART but are nice to keep around if they are present.
if( nf90_inquire_attribute( ncid, VarID, 'long_name') == NF90_NOERR ) then
ret = nf90_get_att(ncid, VarID, 'long_name' , cf_long_name)
call nc_check(ret, routine, 'get_att long_name '//trim(var_name))
domain%variable(ivar)%io_info%long_name = cf_long_name
endif
if( nf90_inquire_attribute( ncid, VarID, 'short_name') == NF90_NOERR ) then
ret = nf90_get_att(ncid, VarID, 'short_name' , cf_short_name)
call nc_check(ret, routine, 'get_att short_name '//trim(var_name))
domain%variable(ivar)%io_info%short_name = cf_short_name
endif
if( nf90_inquire_attribute( ncid, VarID, 'units') == NF90_NOERR ) then
ret = nf90_get_att(ncid, VarID, 'units' , cf_units)
call nc_check(ret, routine, 'get_att units '//trim(var_name))
domain%variable(ivar)%io_info%units = cf_units
endif
! Saving any FillValue, missing_value attributes ...
! Also stuff them into the 'r8' slots to facilitate simpler clamp_variable
! implementation. (Since we coerce the DART state to 'r8')
var_xtype = domain%variable(ivar)%io_info%xtype
select case (var_xtype)
case ( NF90_INT )
! Sometimes the attributes are specified as GLOBAL attributes
ret = nf90_get_att(ncid, NF90_GLOBAL, '_FillValue', cf_spvalINT)
if (ret == NF90_NOERR) then
domain%variable(ivar)%io_info%spvalINT = cf_spvalINT
domain%variable(ivar)%io_info%spvalR8 = real(cf_spvalINT,r8)
domain%variable(ivar)%io_info%has_FillValue = .true.
endif
ret = nf90_get_att(ncid, NF90_GLOBAL, 'missing_value', cf_spvalINT)
if (ret == NF90_NOERR) then
domain%variable(ivar)%io_info%missingINT = cf_spvalINT
domain%variable(ivar)%io_info%spvalR8 = real(cf_spvalINT,r8)
domain%variable(ivar)%io_info%has_missing_value = .true.
endif
! Usually the attributes are specified as variable attributes
ret = nf90_get_att(ncid, VarID, '_FillValue', cf_spvalINT)
if (ret == NF90_NOERR) then
domain%variable(ivar)%io_info%spvalINT = cf_spvalINT
domain%variable(ivar)%io_info%spvalR8 = real(cf_spvalINT,r8)
domain%variable(ivar)%io_info%has_FillValue = .true.
endif
ret = nf90_get_att(ncid, VarID, 'missing_value', cf_spvalINT)
if (ret == NF90_NOERR) then
domain%variable(ivar)%io_info%missingINT = cf_spvalINT
domain%variable(ivar)%io_info%missingR8 = real(cf_spvalINT,r8)
domain%variable(ivar)%io_info%has_missing_value = .true.
endif
case ( NF90_FLOAT )
ret = nf90_get_att(ncid, NF90_GLOBAL, '_FillValue', cf_spvalR4)
if (ret == NF90_NOERR) then
domain%variable(ivar)%io_info%spvalR4 = cf_spvalR4
domain%variable(ivar)%io_info%spvalR8 = real(cf_spvalR4,r8)
domain%variable(ivar)%io_info%has_FillValue = .true.
endif
ret = nf90_get_att(ncid, NF90_GLOBAL, 'missing_value', cf_spvalR4)
if (ret == NF90_NOERR) then
domain%variable(ivar)%io_info%missingR4 = cf_spvalR4
domain%variable(ivar)%io_info%spvalR8 = real(cf_spvalR4,r8)
domain%variable(ivar)%io_info%has_missing_value = .true.
endif
ret = nf90_get_att(ncid, VarID, '_FillValue', cf_spvalR4)
if (ret == NF90_NOERR) then
domain%variable(ivar)%io_info%spvalR4 = cf_spvalR4
domain%variable(ivar)%io_info%spvalR8 = real(cf_spvalR4,r8)
domain%variable(ivar)%io_info%has_FillValue = .true.
endif
ret = nf90_get_att(ncid, VarID, 'missing_value', cf_spvalR4)
if (ret == NF90_NOERR) then
domain%variable(ivar)%io_info%missingR4 = cf_spvalR4
domain%variable(ivar)%io_info%missingR8 = real(cf_spvalR4,r8)
domain%variable(ivar)%io_info%has_missing_value = .true.
endif
case ( NF90_DOUBLE )
! If r8 = r4,
! the missing_value must be present in both missingR4 and missingR8
! ditto for _FillValue.
! This satisfies the overloaded operator 'get_missing_value, get_FillValue'
ret = nf90_get_att(ncid, NF90_GLOBAL, '_FillValue', cf_spvalR8)
if (ret == NF90_NOERR) then
domain%variable(ivar)%io_info%spvalR4 = cf_spvalR8
domain%variable(ivar)%io_info%spvalR8 = cf_spvalR8
domain%variable(ivar)%io_info%has_FillValue = .true.
endif
ret = nf90_get_att(ncid, NF90_GLOBAL, 'missing_value', cf_spvalR8)
if (ret == NF90_NOERR) then
domain%variable(ivar)%io_info%missingR4 = cf_spvalR8
domain%variable(ivar)%io_info%missingR8 = cf_spvalR8
domain%variable(ivar)%io_info%has_missing_value = .true.
endif
ret = nf90_get_att(ncid, VarID, '_FillValue', cf_spvalR8)
if (ret == NF90_NOERR) then
domain%variable(ivar)%io_info%spvalR4 = cf_spvalR8
domain%variable(ivar)%io_info%spvalR8 = cf_spvalR8
domain%variable(ivar)%io_info%has_FillValue = .true.
endif
ret = nf90_get_att(ncid, VarID, 'missing_value', cf_spvalR8)
if (ret == NF90_NOERR) then
domain%variable(ivar)%io_info%missingR4 = cf_spvalR8
domain%variable(ivar)%io_info%missingR8 = cf_spvalR8
domain%variable(ivar)%io_info%has_missing_value = .true.
endif
case DEFAULT
write(string1,*) ' unsupported netcdf variable type : ', var_xtype
call error_handler(E_ERR,routine,string1,source)
end select
! If the variable has one or the other, no problem.
! If the variable has both _FillValue and missing_value attributes, the
! values must be the same or we are lost. DART only supports one missing
! value code. When we go to write, we have no way of knowing which value
! to use as a replacement for the DART missing code.
if ( domain%variable(ivar)%io_info%has_missing_value .and. &
domain%variable(ivar)%io_info%has_FillValue ) then
if ( domain%variable(ivar)%io_info%missingR8 /= &
domain%variable(ivar)%io_info%spvalR8 ) then
write(string1, *) trim(var_name)//' missing_value /= _FillValue '
write(string2,*) 'missing_value is ', domain%variable(ivar)%io_info%missingR8
write(string3,*) '_FillValue is ', domain%variable(ivar)%io_info%spvalR8
call error_handler(E_ERR,'set_dart_missing_value:',string1, &
source, text2=string2, text3=string3)
endif
endif
!>@todo FIXME : Not supporting scale factor or offset at the moment, so just error.
!> To fully support netCDF I/O we need to
!> pack and unpack the variable if these attributes exist.
if (nf90_get_att(ncid, VarID, 'scale_factor', cf_scale_factor) == NF90_NOERR) then
domain%variable(ivar)%io_info%scale_factor = cf_scale_factor
write(string1,*) 'scale_factor not supported at the moment'
write(string2,*) 'contact DART if you would like to get this to work'
call error_handler(E_ERR,routine,string1,source,text2=string2)
endif
if (nf90_get_att(ncid, VarID, 'add_offset', cf_add_offset) == NF90_NOERR) then
domain%variable(ivar)%io_info%add_offset = cf_add_offset
write(string1,*) 'add_offset not supported at the moment'
write(string2,*) 'contact DART if you would like to get this to work'
call error_handler(E_ERR,routine,string1,source,text2=string2)
endif
enddo
! close netcdf file
ret = nf90_close(ncid)
call nc_check(ret, routine, 'nf90_close', trim(ncFilename))
end subroutine load_common_cf_conventions
!-------------------------------------------------------------------------------
!> TIEGCM, the top level of each variable is the boundary condition.
!> So the top level should not be part of the state.
!> How general does this need to be? At the moment it is just a hack to
!> be able to proceed with TIEGCM.
!>
!> This slicing needs to be on the the DART state structure, not the io netcdf
!> structure.
subroutine hyperslice_domain(dom_id, dim_name, new_length)
integer, intent(in) :: dom_id
character(len=*), intent(in) :: dim_name
integer, intent(in) :: new_length
integer :: ivar, jdim
integer(i8) :: orig_domain_size, orig_model_size, index_start
integer(i8) :: domain_size, variable_size
orig_domain_size = get_domain_size(dom_id)
orig_model_size = state%model_size
! Start of the domain remains the same
index_start = state%domain(dom_id)%variable(1)%index_start
! Reduce the dimension size
do ivar = 1, get_num_variables(dom_id)
do jdim = 1, get_num_dims(dom_id, ivar)
if (get_dim_name(dom_id, ivar, jdim) == dim_name) then
state%domain(dom_id)%variable(ivar)%dimlens(jdim) = new_length
endif
enddo
enddo
! need to modify state index accessing for each variable:
! variable_size
! index_start, index_end
! model size
! domain size
domain_size = 0
! variable size
do ivar = 1, get_num_variables(dom_id)
variable_size = 1
do jdim = 1, get_num_dims(dom_id, ivar)
!HK there is a todo about this, what is going on with time and member?
if ((state%domain(dom_id)%variable(ivar)%dimname(jdim) == 'time') .or. &
(state%domain(dom_id)%variable(ivar)%dimname(jdim) == 'member')) cycle
variable_size = variable_size *state%domain(dom_id)%variable(ivar)%dimlens(jdim)