-
Notifications
You must be signed in to change notification settings - Fork 12
/
Copy pathice_history_write.F90
1855 lines (1564 loc) · 78.4 KB
/
ice_history_write.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
! SVN:$Id: ice_history_write.F90 567 2013-01-07 02:57:36Z eclare $
!=======================================================================
!
! Writes history in netCDF format
!
! authors Tony Craig and Bruce Briegleb, NCAR
! Elizabeth C. Hunke and William H. Lipscomb, LANL
! C. M. Bitz, UW
!
! 2004 WHL: Block structure added
! 2006 ECH: Accepted some CCSM code into mainstream CICE
! Added ice_present, aicen, vicen; removed aice1...10, vice1...1.
! Added histfreq_n and histfreq='h' options, removed histfreq='w'
! Converted to free source form (F90)
! Added option for binary output instead of netCDF
! 2009 D Bailey and ECH: Generalized for multiple frequency output
! 2010 Alison McLaren and ECH: Added 3D capability
! 2013 ECH split from ice_history.F90
module ice_history_write
use netcdf
use mpi, only: MPI_INFO_NULL, MPI_COMM_WORLD
use ice_kinds_mod
use ice_constants, only: c0, c360, secday, spval, rad_to_deg
use ice_blocks, only: nx_block, ny_block, block, get_block
use ice_exit, only: abort_ice
use ice_domain, only: distrb_info, nblocks, blocks_ice
use ice_domain, only: equal_num_blocks_per_cpu
use ice_communicate, only: my_task, master_task, MPI_COMM_ICE
use ice_broadcast, only: broadcast_scalar
use ice_gather_scatter, only: gather_global
use ice_domain_size, only: nx_global, ny_global, max_nstrm, max_blocks
use ice_grid, only: TLON, TLAT, ULON, ULAT, hm, bm, tarea, uarea, &
dxu, dxt, dyu, dyt, HTN, HTE, ANGLE, ANGLET, &
lont_bounds, latt_bounds, lonu_bounds, latu_bounds
use ice_history_shared
use ice_itd, only: hin_max
use ice_calendar, only: write_ic, histfreq
implicit none
private
public :: ice_write_hist
save
type coord_attributes ! netcdf coordinate attributes
character (len=11) :: short_name
character (len=45) :: long_name
character (len=20) :: units
end type coord_attributes
type req_attributes ! req'd netcdf attributes
type (coord_attributes) :: req
character (len=20) :: coordinates
end type req_attributes
! 4 coordinate variables: TLON, TLAT, ULON, ULAT
INTEGER (kind=int_kind), PARAMETER :: ncoord = 4
! 4 vertices in each grid cell
INTEGER (kind=int_kind), PARAMETER :: nverts = 4
! 4 variables describe T, U grid boundaries:
! lont_bounds, latt_bounds, lonu_bounds, latu_bounds
INTEGER (kind=int_kind), PARAMETER :: nvar_verts = 4
contains
subroutine check(status, msg)
integer, intent (in) :: status
character(len=*), intent (in) :: msg
if(status /= nf90_noerr) then
call abort_ice('ice: NetCDF error '//trim(nf90_strerror(status)//' '//trim(msg)))
end if
end subroutine check
!=======================================================================
!
! write average ice quantities or snapshots
!
! author: Elizabeth C. Hunke, LANL
subroutine ice_write_hist (ns)
use ice_calendar, only: time, sec, idate, idate0, &
dayyr, days_per_year, use_leap_years
use ice_fileunits, only: nu_diag
use ice_restart_shared, only: runid
integer (kind=int_kind), intent(in) :: ns
! local variables
real (kind=dbl_kind), dimension(:,:), allocatable :: work_g1
real (kind=real_kind), dimension(:,:), allocatable :: work_gr
real (kind=real_kind), dimension(:,:,:), allocatable :: work_gr3
real (kind=dbl_kind), dimension(nx_block,ny_block,max_blocks) :: &
work1
integer (kind=int_kind) :: i,k,ic,n,nn, &
ncid,status,imtid,jmtid,kmtidi,kmtids,kmtidb, cmtid,timid,varid, &
nvertexid,ivertex
integer (kind=int_kind), dimension(3) :: dimid
integer (kind=int_kind), dimension(4) :: dimidz
integer (kind=int_kind), dimension(5) :: dimidcz
integer (kind=int_kind), dimension(3) :: dimid_nverts
integer (kind=int_kind), dimension(4) :: dimidex
real (kind=real_kind) :: ltime
character (char_len) :: title
character (char_len_long) :: ncfile(max_nstrm)
integer (kind=int_kind) :: shuffle, deflate, deflate_level
integer (kind=int_kind) :: ind,boundid
character (char_len) :: start_time,current_date,current_time
character (len=8) :: cdate
TYPE(req_attributes), dimension(nvar) :: var
TYPE(coord_attributes), dimension(ncoord) :: coord_var
TYPE(coord_attributes), dimension(nvar_verts) :: var_nverts
TYPE(coord_attributes), dimension(nvarz) :: var_nz
CHARACTER (char_len), dimension(ncoord) :: coord_bounds
! We leave shuffle at 0, this is only useful for integer data.
shuffle = 0
! If history_deflate_level < 0 then don't do deflation,
! otherwise it sets the deflate level
if (history_deflate_level < 0) then
deflate = 0
deflate_level = 0
else
deflate = 1
deflate_level = history_deflate_level
endif
if (my_task == master_task .or. history_parallel_io) then
ltime=time/int(secday)
call construct_filename(ncfile(ns),'nc',ns)
! add local directory path name to ncfile
if (write_ic) then
ncfile(ns) = trim(incond_dir)//ncfile(ns)
else
ncfile(ns) = trim(history_dir)//ncfile(ns)
endif
! create file
if (history_parallel_io) then
call check(nf90_create(ncfile(ns), ior(NF90_NETCDF4, NF90_MPIIO), ncid, &
comm=MPI_COMM_ICE, info=MPI_INFO_NULL), &
'create history ncfile '//ncfile(ns))
if (.not. equal_num_blocks_per_cpu) then
call abort_ice('ice: error history_parallel_io needs equal_num_blocks_per_cpu')
endif
else
call check(nf90_create(ncfile(ns), ior(NF90_CLASSIC_MODEL, NF90_HDF5), ncid), &
'create history ncfile '//ncfile(ns))
endif
!-----------------------------------------------------------------
! define dimensions
!-----------------------------------------------------------------
if (hist_avg) then
call check(nf90_def_dim(ncid,'d2',2,boundid), 'def dim d2')
endif
call check(nf90_def_dim(ncid, 'ni', nx_global, imtid), &
'def dim ni')
call check(nf90_def_dim(ncid, 'nj', ny_global, jmtid), &
'def dim nj')
call check(nf90_def_dim(ncid, 'nc', ncat_hist, cmtid), &
'def dim nc')
call check(nf90_def_dim(ncid, 'nkice', nzilyr, kmtidi), &
'def dim nkice')
call check(nf90_def_dim(ncid, 'nksnow', nzslyr, kmtids), &
'def dim nksnow')
call check(nf90_def_dim(ncid, 'nkbio', nzblyr, kmtidb), &
'def dim nkbio')
call check(nf90_def_dim(ncid, 'time', 1, timid), &
'def dim time')
call check(nf90_def_dim(ncid, 'nvertices', nverts, nvertexid), &
'def dim nverts')
!-----------------------------------------------------------------
! define coordinate variables
!-----------------------------------------------------------------
call check(nf90_def_var(ncid,'time',nf90_float,timid,varid), &
'def var time')
call check(nf90_put_att(ncid,varid,'long_name','model time'), &
'put_att long_name')
write(cdate,'(i8.8)') idate0
write(title,'(a,a,a,a,a,a,a,a)') 'days since ', &
cdate(1:4),'-',cdate(5:6),'-',cdate(7:8),' 00:00:00'
call check(nf90_put_att(ncid,varid,'units',title), &
'put_att time units')
if (days_per_year == 360) then
call check(nf90_put_att(ncid,varid,'calendar','360_day'), &
'att time calendar')
elseif (days_per_year == 365 .and. .not.use_leap_years ) then
call check(nf90_put_att(ncid,varid,'calendar','NoLeap'), &
'att time calendar')
elseif (use_leap_years) then
call check(nf90_put_att(ncid,varid,'calendar','Gregorian'), &
'att time calendar')
else
call abort_ice( 'ice Error: invalid calendar settings')
endif
if (hist_avg) then
call check(nf90_put_att(ncid,varid,'bounds','time_bounds'), &
'att time bounds')
endif
!-----------------------------------------------------------------
! Define attributes for time bounds if hist_avg is true
!-----------------------------------------------------------------
if (hist_avg) then
dimid(1) = boundid
dimid(2) = timid
call check(nf90_def_var(ncid, 'time_bounds', &
nf90_float,dimid(1:2),varid), &
'def var time_bounds')
call check(nf90_put_att(ncid,varid,'long_name', &
'boundaries for time-averaging interval'), &
'att time_bounds long_name')
write(cdate,'(i8.8)') idate0
write(title,'(a,a,a,a,a,a,a,a)') 'days since ', &
cdate(1:4),'-',cdate(5:6),'-',cdate(7:8),' 00:00:00'
call check(nf90_put_att(ncid,varid,'units',title), &
'att time_bounds units')
endif
!-----------------------------------------------------------------
! define information for required time-invariant variables
!-----------------------------------------------------------------
ind = 0
ind = ind + 1
coord_var(ind) = coord_attributes('TLON', &
'T grid center longitude', 'degrees_east')
coord_bounds(ind) = 'lont_bounds'
ind = ind + 1
coord_var(ind) = coord_attributes('TLAT', &
'T grid center latitude', 'degrees_north')
coord_bounds(ind) = 'latt_bounds'
ind = ind + 1
coord_var(ind) = coord_attributes('ULON', &
'U grid center longitude', 'degrees_east')
coord_bounds(ind) = 'lonu_bounds'
ind = ind + 1
coord_var(ind) = coord_attributes('ULAT', &
'U grid center latitude', 'degrees_north')
coord_bounds(ind) = 'latu_bounds'
var_nz(1) = coord_attributes('NCAT', 'category maximum thickness', 'm')
var_nz(2) = coord_attributes('VGRDi', 'vertical ice levels', '1')
var_nz(3) = coord_attributes('VGRDs', 'vertical snow levels', '1')
var_nz(4) = coord_attributes('VGRDb', 'vertical ice-bio levels', '1')
!-----------------------------------------------------------------
! define information for optional time-invariant variables
!-----------------------------------------------------------------
var(n_tarea)%req = coord_attributes('tarea', &
'area of T grid cells', 'm^2')
var(n_tarea)%coordinates = 'TLON TLAT'
var(n_uarea)%req = coord_attributes('uarea', &
'area of U grid cells', 'm^2')
var(n_uarea)%coordinates = 'ULON ULAT'
var(n_dxt)%req = coord_attributes('dxt', &
'T cell width through middle', 'm')
var(n_dxt)%coordinates = 'TLON TLAT'
var(n_dyt)%req = coord_attributes('dyt', &
'T cell height through middle', 'm')
var(n_dyt)%coordinates = 'TLON TLAT'
var(n_dxu)%req = coord_attributes('dxu', &
'U cell width through middle', 'm')
var(n_dxu)%coordinates = 'ULON ULAT'
var(n_dyu)%req = coord_attributes('dyu', &
'U cell height through middle', 'm')
var(n_dyu)%coordinates = 'ULON ULAT'
var(n_HTN)%req = coord_attributes('HTN', &
'T cell width on North side','m')
var(n_HTN)%coordinates = 'TLON TLAT'
var(n_HTE)%req = coord_attributes('HTE', &
'T cell width on East side', 'm')
var(n_HTE)%coordinates = 'TLON TLAT'
var(n_ANGLE)%req = coord_attributes('ANGLE', &
'angle grid makes with latitude line on U grid', &
'radians')
var(n_ANGLE)%coordinates = 'ULON ULAT'
var(n_ANGLET)%req = coord_attributes('ANGLET', &
'angle grid makes with latitude line on T grid', &
'radians')
var(n_ANGLET)%coordinates = 'TLON TLAT'
! These fields are required for CF compliance
! dimensions (nx,ny,nverts)
var_nverts(n_lont_bnds) = coord_attributes('lont_bounds', &
'longitude boundaries of T cells', 'degrees_east')
var_nverts(n_latt_bnds) = coord_attributes('latt_bounds', &
'latitude boundaries of T cells', 'degrees_north')
var_nverts(n_lonu_bnds) = coord_attributes('lonu_bounds', &
'longitude boundaries of U cells', 'degrees_east')
var_nverts(n_latu_bnds) = coord_attributes('latu_bounds', &
'latitude boundaries of U cells', 'degrees_north')
!-----------------------------------------------------------------
! define attributes for time-invariant variables
!-----------------------------------------------------------------
dimid(1) = imtid
dimid(2) = jmtid
dimid(3) = timid
do i = 1, ncoord
call check(nf90_def_var(ncid, coord_var(i)%short_name, nf90_float, &
dimid(1:2), varid), &
'def var '//coord_var(i)%short_name)
if (history_chunksize_x > 0 .and. history_chunksize_y > 0) then
call check(nf90_def_var_chunking(ncid, varid, NF90_CHUNKED, &
(/ history_chunksize_x, history_chunksize_y /)), &
'def var chunking '//coord_var(i)%short_name)
endif
call check(nf90_def_var_deflate(ncid, varid, shuffle, deflate, &
deflate_level), &
'deflate '//coord_var(i)%short_name)
call check(nf90_put_att(ncid,varid,'long_name',coord_var(i)%long_name), &
'put att long_name '//coord_var(i)%short_name)
call check(nf90_put_att(ncid, varid, 'units', coord_var(i)%units), &
'put att units '//coord_var(i)%short_name)
call check(nf90_put_att(ncid,varid,'missing_value',spval), &
'put att missing_value '//coord_var(i)%short_name)
call check(nf90_put_att(ncid, varid, '_FillValue', spval), &
'put att _FillValue '//coord_var(i)%short_name)
if (coord_var(i)%short_name == 'ULAT') then
call check(nf90_put_att(ncid,varid,'comment', &
'Latitude of NE corner of T grid cell'), &
'put att comment for '//coord_var(i)%short_name)
endif
if (f_bounds) then
call check(nf90_put_att(ncid, varid, 'bounds', coord_bounds(i)), &
'put att bounds '//coord_var(i)%short_name)
endif
enddo
! Extra dimensions (NCAT, NZILYR, NZSLYR, NZBLYR)
dimidex(1)=cmtid
dimidex(2)=kmtidi
dimidex(3)=kmtids
dimidex(4)=kmtidb
do i = 1, nvarz
if (igrdz(i)) then
call check(nf90_def_var(ncid, var_nz(i)%short_name, &
nf90_float, dimidex(i), varid), &
'def var '//var_nz(i)%short_name)
call check(nf90_def_var_deflate(ncid, varid, shuffle, deflate, &
deflate_level), &
'deflate '//var_nz(i)%short_name)
call check(nf90_put_att(ncid,varid,'long_name',var_nz(i)%long_name), &
'put att long_name '//var_nz(i)%short_name)
call check(nf90_put_att(ncid, varid, 'units', var_nz(i)%units), &
'for att units '//var_nz(i)%short_name)
endif
enddo
! Attributes for tmask, blkmask defined separately, since they have no units
if (igrd(n_tmask)) then
call check(nf90_def_var(ncid, 'tmask', nf90_float, dimid(1:2), varid), &
'def var tmask')
if (history_chunksize_x > 0 .and. history_chunksize_y > 0) then
call check(nf90_def_var_chunking(ncid, varid, NF90_CHUNKED, &
(/ history_chunksize_x, history_chunksize_y /)), &
'def var chunking tmask')
endif
call check(nf90_def_var_deflate(ncid, varid, shuffle, deflate, &
deflate_level), 'deflating var tmask')
call check(nf90_put_att(ncid,varid, 'long_name', 'ocean grid mask'), &
'put att tmask long_name')
call check(nf90_put_att(ncid, varid, 'coordinates', 'TLON TLAT'), &
'put att tmask units')
call check(nf90_put_att(ncid,varid,'comment', '0 = land, 1 = ocean'), &
'put att tmask comment')
call check(nf90_put_att(ncid,varid,'missing_value',spval), &
'put att missing_value for tmask')
call check(nf90_put_att(ncid,varid,'_FillValue',spval), &
'put att _FillValue for tmask')
endif
if (igrd(n_blkmask)) then
call check(nf90_def_var(ncid, 'blkmask', nf90_float, dimid(1:2), varid), &
'def var blkmask')
if (history_chunksize_x > 0 .and. history_chunksize_y > 0) then
call check(nf90_def_var_chunking(ncid, varid, NF90_CHUNKED, &
(/ history_chunksize_x, history_chunksize_y /)), &
'def var chunking blkmask')
endif
call check(nf90_def_var_deflate(ncid, varid, shuffle, deflate, &
deflate_level), &
'deflating var blkmask')
call check(nf90_put_att(ncid,varid, 'long_name', 'ice grid block mask'), &
'put att blkmask long_name')
call check(nf90_put_att(ncid, varid, 'coordinates', 'TLON TLAT'), &
'put att blkmask coordinates')
call check(nf90_put_att(ncid,varid,'comment', 'mytask + iblk/100'), &
'put att blkmask comment')
call check(nf90_put_att(ncid,varid,'missing_value',spval), &
'put att blkmask missing_value')
call check(nf90_put_att(ncid,varid,'_FillValue',spval), &
'put att blkmask _FillValue')
endif
do i = 3, nvar ! note n_tmask=1, n_blkmask=2
if (igrd(i)) then
call check(nf90_def_var(ncid, var(i)%req%short_name, &
nf90_float, dimid(1:2), varid), &
'def variable '//var(i)%req%short_name)
if (history_chunksize_x > 0 .and. history_chunksize_y > 0) then
call check(nf90_def_var_chunking(ncid, varid, NF90_CHUNKED, &
(/ history_chunksize_x, history_chunksize_y /)), &
'def var chunking '//var(i)%req%short_name)
endif
call check(nf90_def_var_deflate(ncid, varid, shuffle, deflate, &
deflate_level), &
'deflate var '//var(i)%req%short_name)
call check(nf90_put_att(ncid,varid, 'long_name', var(i)%req%long_name), &
'put att long_name '//var(i)%req%short_name)
call check(nf90_put_att(ncid, varid, 'units', var(i)%req%units), &
'put att units '//var(i)%req%short_name)
call check(nf90_put_att(ncid, varid, 'coordinates', var(i)%coordinates), &
'put att coordinates '//var(i)%req%short_name)
call check(nf90_put_att(ncid,varid,'missing_value',spval), &
'put att missing_value '//var(i)%req%short_name)
call check(nf90_put_att(ncid,varid,'_FillValue',spval), &
'put att _FillValue '//var(i)%req%short_name)
endif
enddo
! Fields with dimensions (nverts,nx,ny)
dimid_nverts(1) = nvertexid
dimid_nverts(2) = imtid
dimid_nverts(3) = jmtid
do i = 1, nvar_verts
if (f_bounds) then
call check(nf90_def_var(ncid, var_nverts(i)%short_name, &
nf90_float,dimid_nverts, varid), &
'def var '//var_nverts(i)%short_name)
if (history_chunksize_x > 0 .and. history_chunksize_y > 0) then
call check(nf90_def_var_chunking(ncid, varid, NF90_CHUNKED, &
(/ 1, history_chunksize_x, history_chunksize_y /)), &
'def var chunking '//var_nverts(i)%short_name)
endif
call check(nf90_def_var_deflate(ncid, varid, shuffle, deflate, &
deflate_level), &
'deflate var '//var_nverts(i)%short_name)
call check(nf90_put_att(ncid,varid, 'long_name', &
var_nverts(i)%long_name), &
'put att long_name '//var_nverts(i)%short_name)
call check(nf90_put_att(ncid, varid, 'units', var_nverts(i)%units), &
'put att units '//var_nverts(i)%short_name)
call check(nf90_put_att(ncid,varid,'missing_value',spval), &
'put att missing_value '//var_nverts(i)%short_name)
call check(nf90_put_att(ncid,varid,'_FillValue',spval), &
'put att _FillValue '//var_nverts(i)%short_name)
endif
enddo
do n=1,num_avail_hist_fields_2D
if (avail_hist_fields(n)%vhistfreq == histfreq(ns) .or. write_ic) then
call check(nf90_def_var(ncid, avail_hist_fields(n)%vname, &
nf90_float, dimid, varid), &
'def var '//avail_hist_fields(n)%vname)
if (history_chunksize_x > 0 .and. history_chunksize_y > 0) then
call check(nf90_def_var_chunking(ncid, varid, NF90_CHUNKED, &
(/ history_chunksize_x, history_chunksize_y, 1 /)), &
'def var chunking '//avail_hist_fields(n)%vname)
endif
call check(nf90_def_var_deflate(ncid, varid, shuffle, deflate, &
deflate_level), &
'deflate '//avail_hist_fields(n)%vname)
call check(nf90_put_att(ncid,varid,'units', &
avail_hist_fields(n)%vunit), &
'put att units '//avail_hist_fields(n)%vname)
call check(nf90_put_att(ncid,varid, 'long_name', &
avail_hist_fields(n)%vdesc), &
'put att long_name '//avail_hist_fields(n)%vname)
call check(nf90_put_att(ncid,varid,'coordinates', &
avail_hist_fields(n)%vcoord), &
'put att coordinates '//avail_hist_fields(n)%vname)
call check(nf90_put_att(ncid,varid,'cell_measures', &
avail_hist_fields(n)%vcellmeas), &
'put att cell_measures '//avail_hist_fields(n)%vname)
call check(nf90_put_att(ncid,varid,'missing_value',spval), &
'put att missing_value '//avail_hist_fields(n)%vname)
call check(nf90_put_att(ncid,varid,'_FillValue',spval), &
'put att _FillValue '//avail_hist_fields(n)%vname)
!---------------------------------------------------------------
! Add cell_methods attribute to variables if averaged
!---------------------------------------------------------------
if (hist_avg) then
if (TRIM(avail_hist_fields(n)%vname)/='sig1' .or. &
TRIM(avail_hist_fields(n)%vname)/='sig2') then
call check(nf90_put_att(ncid,varid,'cell_methods','time: mean'), &
'put att cell methods time mean '//avail_hist_fields(n)%vname)
endif
endif
if (histfreq(ns) == '1' .or. .not. hist_avg &
.or. n==n_divu(ns) .or. n==n_shear(ns) & ! snapshots
.or. n==n_sig1(ns) .or. n==n_sig2(ns) &
.or. n==n_trsig(ns) &
.or. n==n_mlt_onset(ns) .or. n==n_frz_onset(ns) &
.or. n==n_hisnap(ns) .or. n==n_aisnap(ns)) then
call check(nf90_put_att(ncid,varid,'time_rep','instantaneous'), &
'put att time_rep instantaneous')
else
call check(nf90_put_att(ncid,varid,'time_rep','averaged'), &
'put att time_rep averaged')
endif
endif
enddo ! num_avail_hist_fields_2D
dimidz(1) = imtid
dimidz(2) = jmtid
dimidz(3) = cmtid
dimidz(4) = timid
do n = n2D + 1, n3Dccum
if (avail_hist_fields(n)%vhistfreq == histfreq(ns) .or. write_ic) then
status = nf90_def_var(ncid, avail_hist_fields(n)%vname, &
nf90_float, dimidz, varid)
if (status /= nf90_noerr) call abort_ice( &
'Error defining variable '//avail_hist_fields(n)%vname)
if (history_chunksize_x > 0 .and. history_chunksize_y > 0) then
call check(nf90_def_var_chunking(ncid, varid, NF90_CHUNKED, &
(/ history_chunksize_x, history_chunksize_y, 1, 1 /)), &
'def var chunking '//avail_hist_fields(n)%vname)
endif
status = nf90_def_var_deflate(ncid, varid, shuffle, deflate, &
deflate_level)
if (status /= nf90_noerr) call abort_ice( &
'Error deflating variable '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'units', &
avail_hist_fields(n)%vunit)
if (status /= nf90_noerr) call abort_ice( &
'Error defining units for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid, 'long_name', &
avail_hist_fields(n)%vdesc)
if (status /= nf90_noerr) call abort_ice( &
'Error defining long_name for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'coordinates', &
avail_hist_fields(n)%vcoord)
if (status /= nf90_noerr) call abort_ice( &
'Error defining coordinates for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'cell_measures', &
avail_hist_fields(n)%vcellmeas)
if (status /= nf90_noerr) call abort_ice( &
'Error defining cell measures for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'missing_value',spval)
if (status /= nf90_noerr) call abort_ice( &
'Error defining missing_value for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'_FillValue',spval)
if (status /= nf90_noerr) call abort_ice( &
'Error defining _FillValue for '//avail_hist_fields(n)%vname)
!---------------------------------------------------------------
! Add cell_methods attribute to variables if averaged
!---------------------------------------------------------------
if (hist_avg) then
status = nf90_put_att(ncid,varid,'cell_methods','time: mean')
if (status /= nf90_noerr) call abort_ice( &
'Error defining cell methods for '//avail_hist_fields(n)%vname)
endif
if (histfreq(ns) == '1' .or. .not. hist_avg) then
status = nf90_put_att(ncid,varid,'time_rep','instantaneous')
else
status = nf90_put_att(ncid,varid,'time_rep','averaged')
endif
endif
enddo ! num_avail_hist_fields_3Dc
dimidz(1) = imtid
dimidz(2) = jmtid
dimidz(3) = kmtidi
dimidz(4) = timid
do n = n3Dccum + 1, n3Dzcum
if (avail_hist_fields(n)%vhistfreq == histfreq(ns) .or. write_ic) then
status = nf90_def_var(ncid, avail_hist_fields(n)%vname, &
nf90_float, dimidz, varid)
if (status /= nf90_noerr) call abort_ice( &
'Error defining variable '//avail_hist_fields(n)%vname)
if (history_chunksize_x > 0 .and. history_chunksize_y > 0) then
call check(nf90_def_var_chunking(ncid, varid, NF90_CHUNKED, &
(/ history_chunksize_x, history_chunksize_y, 1, 1 /)), &
'def var chunking '//avail_hist_fields(n)%vname)
endif
status = nf90_def_var_deflate(ncid, varid, shuffle, deflate, &
deflate_level)
if (status /= nf90_noerr) call abort_ice( &
'Error deflating variable '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'units', &
avail_hist_fields(n)%vunit)
if (status /= nf90_noerr) call abort_ice( &
'Error defining units for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid, 'long_name', &
avail_hist_fields(n)%vdesc)
if (status /= nf90_noerr) call abort_ice( &
'Error defining long_name for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'coordinates', &
avail_hist_fields(n)%vcoord)
if (status /= nf90_noerr) call abort_ice( &
'Error defining coordinates for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'cell_measures', &
avail_hist_fields(n)%vcellmeas)
if (status /= nf90_noerr) call abort_ice( &
'Error defining cell measures for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'missing_value',spval)
if (status /= nf90_noerr) call abort_ice( &
'Error defining missing_value for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'_FillValue',spval)
if (status /= nf90_noerr) call abort_ice( &
'Error defining _FillValue for '//avail_hist_fields(n)%vname)
endif
enddo ! num_avail_hist_fields_3Dz
dimidz(1) = imtid
dimidz(2) = jmtid
dimidz(3) = kmtidb
dimidz(4) = timid
do n = n3Dzcum + 1, n3Dbcum
if (avail_hist_fields(n)%vhistfreq == histfreq(ns) .or. write_ic) then
status = nf90_def_var(ncid, avail_hist_fields(n)%vname, &
nf90_float, dimidz, varid)
if (status /= nf90_noerr) call abort_ice( &
'Error defining variable '//avail_hist_fields(n)%vname)
if (history_chunksize_x > 0 .and. history_chunksize_y > 0) then
call check(nf90_def_var_chunking(ncid, varid, NF90_CHUNKED, &
(/ history_chunksize_x, history_chunksize_y, 1, 1 /)), &
'def var chunking '//avail_hist_fields(n)%vname)
endif
status = nf90_def_var_deflate(ncid, varid, shuffle, deflate, &
deflate_level)
if (status /= nf90_noerr) call abort_ice( &
'Error deflating variable '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'units', &
avail_hist_fields(n)%vunit)
if (status /= nf90_noerr) call abort_ice( &
'Error defining units for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid, 'long_name', &
avail_hist_fields(n)%vdesc)
if (status /= nf90_noerr) call abort_ice( &
'Error defining long_name for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'coordinates', &
avail_hist_fields(n)%vcoord)
if (status /= nf90_noerr) call abort_ice( &
'Error defining coordinates for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'cell_measures', &
avail_hist_fields(n)%vcellmeas)
if (status /= nf90_noerr) call abort_ice( &
'Error defining cell measures for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'missing_value',spval)
if (status /= nf90_noerr) call abort_ice( &
'Error defining missing_value for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'_FillValue',spval)
if (status /= nf90_noerr) call abort_ice( &
'Error defining _FillValue for '//avail_hist_fields(n)%vname)
endif
enddo ! num_avail_hist_fields_3Db
dimidcz(1) = imtid
dimidcz(2) = jmtid
dimidcz(3) = kmtidi
dimidcz(4) = cmtid
dimidcz(5) = timid
do n = n3Dbcum + 1, n4Dicum
if (avail_hist_fields(n)%vhistfreq == histfreq(ns) .or. write_ic) then
status = nf90_def_var(ncid, avail_hist_fields(n)%vname, &
!nf90_float, dimidcz, varid) ! ferret
nf90_float, dimidcz(1:4), varid)
if (status /= nf90_noerr) call abort_ice( &
'Error defining variable '//avail_hist_fields(n)%vname)
if (history_chunksize_x > 0 .and. history_chunksize_y > 0) then
call check(nf90_def_var_chunking(ncid, varid, NF90_CHUNKED, &
(/ history_chunksize_x, history_chunksize_y, 1, 1 /)), &
'def var chunking '//avail_hist_fields(n)%vname)
endif
status = nf90_def_var_deflate(ncid, varid, shuffle, deflate, &
deflate_level)
if (status /= nf90_noerr) call abort_ice( &
'Error deflating variable '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'units', &
avail_hist_fields(n)%vunit)
if (status /= nf90_noerr) call abort_ice( &
'Error defining units for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid, 'long_name', &
avail_hist_fields(n)%vdesc)
if (status /= nf90_noerr) call abort_ice( &
'Error defining long_name for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'coordinates', &
avail_hist_fields(n)%vcoord)
if (status /= nf90_noerr) call abort_ice( &
'Error defining coordinates for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'cell_measures', &
avail_hist_fields(n)%vcellmeas)
if (status /= nf90_noerr) call abort_ice( &
'Error defining cell measures for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'missing_value',spval)
if (status /= nf90_noerr) call abort_ice( &
'Error defining missing_value for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'_FillValue',spval)
if (status /= nf90_noerr) call abort_ice( &
'Error defining _FillValue for '//avail_hist_fields(n)%vname)
!---------------------------------------------------------------
! Add cell_methods attribute to variables if averaged
!---------------------------------------------------------------
if (hist_avg) then
status = nf90_put_att(ncid,varid,'cell_methods','time: mean')
if (status /= nf90_noerr) call abort_ice( &
'Error defining cell methods for '//avail_hist_fields(n)%vname)
endif
if (histfreq(ns) == '1' .or. .not. hist_avg) then
status = nf90_put_att(ncid,varid,'time_rep','instantaneous')
else
status = nf90_put_att(ncid,varid,'time_rep','averaged')
endif
endif
enddo ! num_avail_hist_fields_4Di
dimidcz(1) = imtid
dimidcz(2) = jmtid
dimidcz(3) = kmtids
dimidcz(4) = cmtid
dimidcz(5) = timid
do n = n4Dicum + 1, n4Dscum
if (avail_hist_fields(n)%vhistfreq == histfreq(ns) .or. write_ic) then
status = nf90_def_var(ncid, avail_hist_fields(n)%vname, &
!nf90_float, dimidcz, varid) ! ferret
nf90_float, dimidcz(1:4), varid)
if (status /= nf90_noerr) call abort_ice( &
'Error defining variable '//avail_hist_fields(n)%vname)
if (history_chunksize_x > 0 .and. history_chunksize_y > 0) then
call check(nf90_def_var_chunking(ncid, varid, NF90_CHUNKED, &
(/ history_chunksize_x, history_chunksize_y, 1, 1 /)), &
'def var chunking '//avail_hist_fields(n)%vname)
endif
status = nf90_def_var_deflate(ncid, varid, shuffle, deflate, &
deflate_level)
if (status /= nf90_noerr) call abort_ice( &
'Error deflating variable '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'units', &
avail_hist_fields(n)%vunit)
if (status /= nf90_noerr) call abort_ice( &
'Error defining units for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid, 'long_name', &
avail_hist_fields(n)%vdesc)
if (status /= nf90_noerr) call abort_ice( &
'Error defining long_name for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'coordinates', &
avail_hist_fields(n)%vcoord)
if (status /= nf90_noerr) call abort_ice( &
'Error defining coordinates for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'cell_measures', &
avail_hist_fields(n)%vcellmeas)
if (status /= nf90_noerr) call abort_ice( &
'Error defining cell measures for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'missing_value',spval)
if (status /= nf90_noerr) call abort_ice( &
'Error defining missing_value for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'_FillValue',spval)
if (status /= nf90_noerr) call abort_ice( &
'Error defining _FillValue for '//avail_hist_fields(n)%vname)
!---------------------------------------------------------------
! Add cell_methods attribute to variables if averaged
!---------------------------------------------------------------
if (hist_avg) then
status = nf90_put_att(ncid,varid,'cell_methods','time: mean')
if (status /= nf90_noerr) call abort_ice( &
'Error defining cell methods for '//avail_hist_fields(n)%vname)
endif
if (histfreq(ns) == '1' .or. .not. hist_avg) then
status = nf90_put_att(ncid,varid,'time_rep','instantaneous')
else
status = nf90_put_att(ncid,varid,'time_rep','averaged')
endif
endif
enddo ! num_avail_hist_fields_4Ds
dimidcz(1) = imtid
dimidcz(2) = jmtid
dimidcz(3) = kmtidb
dimidcz(4) = cmtid
dimidcz(5) = timid
do n = n4Dscum + 1, n4Dbcum
if (avail_hist_fields(n)%vhistfreq == histfreq(ns) .or. write_ic) then
status = nf90_def_var(ncid, avail_hist_fields(n)%vname, &
!nf90_float, dimidcz, varid) ! ferret
nf90_float, dimidcz(1:4), varid)
if (status /= nf90_noerr) call abort_ice( &
'Error defining variable '//avail_hist_fields(n)%vname)
if (history_chunksize_x > 0 .and. history_chunksize_y > 0) then
call check(nf90_def_var_chunking(ncid, varid, NF90_CHUNKED, &
(/ history_chunksize_x, history_chunksize_y, 1, 1 /)), &
'def var chunking '//avail_hist_fields(n)%vname)
endif
status = nf90_def_var_deflate(ncid, varid, shuffle, deflate, &
deflate_level)
if (status /= nf90_noerr) call abort_ice( &
'Error deflating variable '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'units', &
avail_hist_fields(n)%vunit)
if (status /= nf90_noerr) call abort_ice( &
'Error defining units for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid, 'long_name', &
avail_hist_fields(n)%vdesc)
if (status /= nf90_noerr) call abort_ice( &
'Error defining long_name for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'coordinates', &
avail_hist_fields(n)%vcoord)
if (status /= nf90_noerr) call abort_ice( &
'Error defining coordinates for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'cell_measures', &
avail_hist_fields(n)%vcellmeas)
if (status /= nf90_noerr) call abort_ice( &
'Error defining cell measures for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'missing_value',spval)
if (status /= nf90_noerr) call abort_ice( &
'Error defining missing_value for '//avail_hist_fields(n)%vname)
status = nf90_put_att(ncid,varid,'_FillValue',spval)
if (status /= nf90_noerr) call abort_ice( &
'Error defining _FillValue for '//avail_hist_fields(n)%vname)
!---------------------------------------------------------------
! Add cell_methods attribute to variables if averaged
!---------------------------------------------------------------
if (hist_avg) then
status = nf90_put_att(ncid,varid,'cell_methods','time: mean')
if (status /= nf90_noerr) call abort_ice( &
'Error defining cell methods for '//avail_hist_fields(n)%vname)
endif
if (histfreq(ns) == '1' .or. .not. hist_avg) then
call check(nf90_put_att(ncid,varid,'time_rep','instantaneous'), &
'put att time_rep instantaneous')
else
call check(nf90_put_att(ncid,varid,'time_rep','averaged'), &
'put att time_rep averaged')
endif
endif
enddo ! num_avail_hist_fields_4Db
!-----------------------------------------------------------------
! global attributes
!-----------------------------------------------------------------
! ... the user should change these to something useful ...
!-----------------------------------------------------------------
#ifdef CCSMCOUPLED
call check(nf90_put_att(ncid,nf90_global,'title',runid), &
'in global attribute title')
#else
title = 'sea ice model output for CICE'
call check(nf90_put_att(ncid,nf90_global,'title',title), &
'global attribute title')
#endif
title = 'Diagnostic and Prognostic Variables'
call check(nf90_put_att(ncid,nf90_global,'contents',title), &
'global attribute contents')
title = 'Los Alamos Sea Ice Model (CICE) Version 5'
call check(nf90_put_att(ncid,nf90_global,'source',title), &
'global attribute source')
#ifdef AusCOM
write(title,'(a,i3,a)') 'This Year Has ',int(dayyr),' days'
#else
if (use_leap_years) then
write(title,'(a,i3,a)') 'This year has ',int(dayyr),' days'
else
write(title,'(a,i3,a)') 'All years have exactly ',int(dayyr),' days'
endif
#endif
call check(nf90_put_att(ncid,nf90_global,'comment',title), &
'global attribute comment')
write(title,'(a,i8.8)') 'File written on model date ',idate
call check(nf90_put_att(ncid,nf90_global,'comment2',title), &
'global attribute date1')
write(title,'(a,i6)') 'seconds elapsed into model date: ',sec
call check(nf90_put_att(ncid,nf90_global,'comment3',title), &
'global attribute date2')
title = 'CF-1.0'
call check(nf90_put_att(ncid,nf90_global,'conventions',title), &
'global attribute conventions')
call date_and_time(date=current_date, time=current_time)
write(start_time,1000) current_date(1:4), current_date(5:6), &
current_date(7:8), current_time(1:2), &
current_time(3:4), current_time(5:8)
1000 format('This dataset was created on ', &
a,'-',a,'-',a,' at ',a,':',a,':',a)
call check(nf90_put_att(ncid,nf90_global,'history',start_time), &
'global attribute history')
call check(nf90_put_att(ncid,nf90_global,'io_flavor','io_netcdf'), &
'global attribute io_flavor')
!-----------------------------------------------------------------
! end define mode
!-----------------------------------------------------------------
call check(nf90_enddef(ncid), 'enddef')
!-----------------------------------------------------------------
! write time variable
!-----------------------------------------------------------------
call check(nf90_inq_varid(ncid,'time',varid), &
'inq varid time')
call check(nf90_put_var(ncid,varid,ltime), &
'put var time')
!-----------------------------------------------------------------
! write time_bounds info
!-----------------------------------------------------------------
if (hist_avg) then
call check(nf90_inq_varid(ncid,'time_bounds',varid), &