-
Notifications
You must be signed in to change notification settings - Fork 139
/
Copy pathice_grid.F90
2562 lines (2065 loc) · 92 KB
/
ice_grid.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
#ifdef ncdf
#define USE_NETCDF
#endif
!=======================================================================
! Spatial grids, masks, and boundary conditions
!
! authors: Elizabeth C. Hunke and William H. Lipscomb, LANL
! Tony Craig, NCAR
!
! 2004: Block structure added by William Lipscomb
! init_grid split into two parts as in POP 2.0
! Boundary update routines replaced by POP versions
! 2006: Converted to free source form (F90) by Elizabeth Hunke
! 2007: Option to read from netcdf files (A. Keen, Met Office)
! Grid reading routines reworked by E. Hunke for boundary values
module ice_grid
use ice_kinds_mod
use ice_broadcast, only: broadcast_scalar, broadcast_array
use ice_boundary, only: ice_HaloUpdate, ice_HaloExtrapolate
use ice_communicate, only: my_task, master_task
use ice_blocks, only: block, get_block, nx_block, ny_block, nghost
use ice_domain_size, only: nx_global, ny_global, max_blocks
use ice_domain, only: blocks_ice, nblocks, halo_info, distrb_info, &
ew_boundary_type, ns_boundary_type, init_domain_distribution
use ice_fileunits, only: nu_diag, nu_grid, nu_kmt, &
get_fileunit, release_fileunit
use ice_gather_scatter, only: gather_global, scatter_global
use ice_read_write, only: ice_read, ice_read_nc, ice_read_global, &
ice_read_global_nc, ice_open, ice_open_nc, ice_close_nc
use ice_timers, only: timer_bound, ice_timer_start, ice_timer_stop
use ice_exit, only: abort_ice
use ice_global_reductions, only: global_minval, global_maxval
use icepack_intfc, only: icepack_warnings_flush, icepack_warnings_aborted
use icepack_intfc, only: icepack_query_parameters
implicit none
private
public :: init_grid1, init_grid2, &
t2ugrid_vector, u2tgrid_vector, &
to_ugrid, to_tgrid, alloc_grid
character (len=char_len_long), public :: &
grid_format , & ! file format ('bin'=binary or 'nc'=netcdf)
gridcpl_file , & ! input file for POP coupling grid info
grid_file , & ! input file for POP grid info
kmt_file , & ! input file for POP grid info
bathymetry_file, & ! input bathymetry for basalstress
bathymetry_format, & ! bathymetry file format (default or pop)
grid_spacing , & ! default of 30.e3m or set by user in namelist
grid_type ! current options are rectangular (default),
! displaced_pole, tripole, regional
real (kind=dbl_kind), dimension (:,:,:), allocatable, public :: &
dxt , & ! width of T-cell through the middle (m)
dyt , & ! height of T-cell through the middle (m)
dxu , & ! width of U-cell through the middle (m)
dyu , & ! height of U-cell through the middle (m)
HTE , & ! length of eastern edge of T-cell (m)
HTN , & ! length of northern edge of T-cell (m)
tarea , & ! area of T-cell (m^2)
uarea , & ! area of U-cell (m^2)
tarear , & ! 1/tarea
uarear , & ! 1/uarea
tinyarea,& ! puny*tarea
tarean , & ! area of NH T-cells
tareas , & ! area of SH T-cells
ULON , & ! longitude of velocity pts (radians)
ULAT , & ! latitude of velocity pts (radians)
TLON , & ! longitude of temp pts (radians)
TLAT , & ! latitude of temp pts (radians)
ANGLE , & ! for conversions between POP grid and lat/lon
ANGLET , & ! ANGLE converted to T-cells
bathymetry , & ! ocean depth, for grounding keels and bergs (m)
ocn_gridcell_frac ! only relevant for lat-lon grids
! gridcell value of [1 - (land fraction)] (T-cell)
real (kind=dbl_kind), dimension (:,:,:), allocatable, public :: &
cyp , & ! 1.5*HTE - 0.5*HTE
cxp , & ! 1.5*HTN - 0.5*HTN
cym , & ! 0.5*HTE - 1.5*HTE
cxm , & ! 0.5*HTN - 1.5*HTN
dxhy , & ! 0.5*(HTE - HTE)
dyhx ! 0.5*(HTN - HTN)
! grid dimensions for rectangular grid
real (kind=dbl_kind), public :: &
dxrect, & ! user_specified spacing (cm) in x-direction (uniform HTN)
dyrect ! user_specified spacing (cm) in y-direction (uniform HTE)
! Corners of grid boxes for history output
real (kind=dbl_kind), dimension (:,:,:,:), allocatable, public :: &
lont_bounds, & ! longitude of gridbox corners for T point
latt_bounds, & ! latitude of gridbox corners for T point
lonu_bounds, & ! longitude of gridbox corners for U point
latu_bounds ! latitude of gridbox corners for U point
! geometric quantities used for remapping transport
real (kind=dbl_kind), dimension (:,:,:), allocatable, public :: &
xav , & ! mean T-cell value of x
yav , & ! mean T-cell value of y
xxav , & ! mean T-cell value of xx
! xyav , & ! mean T-cell value of xy
! yyav , & ! mean T-cell value of yy
yyav ! mean T-cell value of yy
! xxxav, & ! mean T-cell value of xxx
! xxyav, & ! mean T-cell value of xxy
! xyyav, & ! mean T-cell value of xyy
! yyyav ! mean T-cell value of yyy
real (kind=dbl_kind), &
dimension (:,:,:,:,:), allocatable, public :: &
mne, & ! matrices used for coordinate transformations in remapping
mnw, & ! ne = northeast corner, nw = northwest, etc.
mse, &
msw
! masks
real (kind=dbl_kind), dimension (:,:,:), allocatable, public :: &
hm , & ! land/boundary mask, thickness (T-cell)
bm , & ! task/block id
uvm , & ! land/boundary mask, velocity (U-cell)
kmt ! ocean topography mask for bathymetry (T-cell)
logical (kind=log_kind), public :: &
use_bathymetry ! flag for reading in bathymetry_file
logical (kind=log_kind), &
dimension (:,:,:), allocatable, public :: &
tmask , & ! land/boundary mask, thickness (T-cell)
umask , & ! land/boundary mask, velocity (U-cell)
lmask_n, & ! northern hemisphere mask
lmask_s ! southern hemisphere mask
real (kind=dbl_kind), dimension (:,:,:), allocatable, public :: &
rndex_global ! global index for local subdomain (dbl)
logical (kind=log_kind), private :: &
l_readCenter ! If anglet exist in grid file read it otherwise calculate it
!=======================================================================
contains
!=======================================================================
!
! Allocate space for all variables
!
subroutine alloc_grid
integer (int_kind) :: ierr
allocate( &
dxt (nx_block,ny_block,max_blocks), & ! width of T-cell through the middle (m)
dyt (nx_block,ny_block,max_blocks), & ! height of T-cell through the middle (m)
dxu (nx_block,ny_block,max_blocks), & ! width of U-cell through the middle (m)
dyu (nx_block,ny_block,max_blocks), & ! height of U-cell through the middle (m)
HTE (nx_block,ny_block,max_blocks), & ! length of eastern edge of T-cell (m)
HTN (nx_block,ny_block,max_blocks), & ! length of northern edge of T-cell (m)
tarea (nx_block,ny_block,max_blocks), & ! area of T-cell (m^2)
uarea (nx_block,ny_block,max_blocks), & ! area of U-cell (m^2)
tarear (nx_block,ny_block,max_blocks), & ! 1/tarea
uarear (nx_block,ny_block,max_blocks), & ! 1/uarea
tinyarea (nx_block,ny_block,max_blocks), & ! puny*tarea
tarean (nx_block,ny_block,max_blocks), & ! area of NH T-cells
tareas (nx_block,ny_block,max_blocks), & ! area of SH T-cells
ULON (nx_block,ny_block,max_blocks), & ! longitude of velocity pts (radians)
ULAT (nx_block,ny_block,max_blocks), & ! latitude of velocity pts (radians)
TLON (nx_block,ny_block,max_blocks), & ! longitude of temp pts (radians)
TLAT (nx_block,ny_block,max_blocks), & ! latitude of temp pts (radians)
ANGLE (nx_block,ny_block,max_blocks), & ! for conversions between POP grid and lat/lon
ANGLET (nx_block,ny_block,max_blocks), & ! ANGLE converted to T-cells
bathymetry(nx_block,ny_block,max_blocks),& ! ocean depth, for grounding keels and bergs (m)
ocn_gridcell_frac(nx_block,ny_block,max_blocks),& ! only relevant for lat-lon grids
cyp (nx_block,ny_block,max_blocks), & ! 1.5*HTE - 0.5*HTE
cxp (nx_block,ny_block,max_blocks), & ! 1.5*HTN - 0.5*HTN
cym (nx_block,ny_block,max_blocks), & ! 0.5*HTE - 1.5*HTE
cxm (nx_block,ny_block,max_blocks), & ! 0.5*HTN - 1.5*HTN
dxhy (nx_block,ny_block,max_blocks), & ! 0.5*(HTE - HTE)
dyhx (nx_block,ny_block,max_blocks), & ! 0.5*(HTN - HTN)
xav (nx_block,ny_block,max_blocks), & ! mean T-cell value of x
yav (nx_block,ny_block,max_blocks), & ! mean T-cell value of y
xxav (nx_block,ny_block,max_blocks), & ! mean T-cell value of xx
yyav (nx_block,ny_block,max_blocks), & ! mean T-cell value of yy
hm (nx_block,ny_block,max_blocks), & ! land/boundary mask, thickness (T-cell)
bm (nx_block,ny_block,max_blocks), & ! task/block id
uvm (nx_block,ny_block,max_blocks), & ! land/boundary mask, velocity (U-cell)
kmt (nx_block,ny_block,max_blocks), & ! ocean topography mask for bathymetry (T-cell)
tmask (nx_block,ny_block,max_blocks), & ! land/boundary mask, thickness (T-cell)
umask (nx_block,ny_block,max_blocks), & ! land/boundary mask, velocity (U-cell)
lmask_n (nx_block,ny_block,max_blocks), & ! northern hemisphere mask
lmask_s (nx_block,ny_block,max_blocks), & ! southern hemisphere mask
rndex_global(nx_block,ny_block,max_blocks), & ! global index for local subdomain (dbl)
lont_bounds(4,nx_block,ny_block,max_blocks), & ! longitude of gridbox corners for T point
latt_bounds(4,nx_block,ny_block,max_blocks), & ! latitude of gridbox corners for T point
lonu_bounds(4,nx_block,ny_block,max_blocks), & ! longitude of gridbox corners for U point
latu_bounds(4,nx_block,ny_block,max_blocks), & ! latitude of gridbox corners for U point
mne (2,2,nx_block,ny_block,max_blocks), & ! matrices used for coordinate transformations in remapping
mnw (2,2,nx_block,ny_block,max_blocks), & ! ne = northeast corner, nw = northwest, etc.
mse (2,2,nx_block,ny_block,max_blocks), &
msw (2,2,nx_block,ny_block,max_blocks), &
stat=ierr)
if (ierr/=0) call abort_ice('(alloc_grid): Out of memory')
end subroutine alloc_grid
!=======================================================================
! Distribute blocks across processors. The distribution is optimized
! based on latitude and topography, contained in the ULAT and KMT arrays.
!
! authors: William Lipscomb and Phil Jones, LANL
subroutine init_grid1
use ice_blocks, only: nx_block, ny_block
use ice_broadcast, only: broadcast_array
use ice_constants, only: c1
integer (kind=int_kind) :: &
fid_grid, & ! file id for netCDF grid file
fid_kmt ! file id for netCDF kmt file
character (char_len) :: &
fieldname ! field name in netCDF file
real (kind=dbl_kind), dimension(:,:), allocatable :: &
work_g1, work_g2
real (kind=dbl_kind) :: &
rad_to_deg
character(len=*), parameter :: subname = '(init_grid1)'
!-----------------------------------------------------------------
! Get global ULAT and KMT arrays used for block decomposition.
!-----------------------------------------------------------------
call icepack_query_parameters(rad_to_deg_out=rad_to_deg)
call icepack_warnings_flush(nu_diag)
if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
file=__FILE__, line=__LINE__)
allocate(work_g1(nx_global,ny_global))
allocate(work_g2(nx_global,ny_global))
if (trim(grid_type) == 'displaced_pole' .or. &
trim(grid_type) == 'tripole' .or. &
trim(grid_type) == 'regional' ) then
if (trim(grid_format) == 'nc') then
call ice_open_nc(grid_file,fid_grid)
call ice_open_nc(kmt_file,fid_kmt)
fieldname='ulat'
call ice_read_global_nc(fid_grid,1,fieldname,work_g1,.true.)
fieldname='kmt'
call ice_read_global_nc(fid_kmt,1,fieldname,work_g2,.true.)
if (my_task == master_task) then
call ice_close_nc(fid_grid)
call ice_close_nc(fid_kmt)
endif
else
call ice_open(nu_grid,grid_file,64) ! ULAT
call ice_open(nu_kmt, kmt_file, 32) ! KMT
call ice_read_global(nu_grid,1,work_g1,'rda8',.true.) ! ULAT
call ice_read_global(nu_kmt, 1,work_g2,'ida4',.true.) ! KMT
if (my_task == master_task) then
close (nu_grid)
close (nu_kmt)
endif
endif
else ! rectangular grid
work_g1(:,:) = 75._dbl_kind/rad_to_deg ! arbitrary polar latitude
work_g2(:,:) = c1
endif
call broadcast_array(work_g1, master_task) ! ULAT
call broadcast_array(work_g2, master_task) ! KMT
!-----------------------------------------------------------------
! distribute blocks among processors
!-----------------------------------------------------------------
call init_domain_distribution(work_g2, work_g1) ! KMT, ULAT
deallocate(work_g1)
deallocate(work_g2)
!-----------------------------------------------------------------
! write additional domain information
!-----------------------------------------------------------------
if (my_task == master_task) then
write(nu_diag,'(a26,i6)') ' Block size: nx_block = ',nx_block
write(nu_diag,'(a26,i6)') ' ny_block = ',ny_block
endif
end subroutine init_grid1
!=======================================================================
! Horizontal grid initialization:
!
! U{LAT,LONG} = true {latitude,longitude} of U points
! HT{N,E} = cell widths on {N,E} sides of T cell
! ANGLE = angle between local x direction and true east
! hm = land mask (c1 for ocean points, c0 for land points)
! D{X,Y}{T,U} = {x,y} spacing centered at {T,U} points
! T-grid and ghost cell values
! Various grid quantities needed for dynamics and transport
!
! author: Elizabeth C. Hunke, LANL
subroutine init_grid2
use ice_blocks, only: get_block, block, nx_block, ny_block
use ice_constants, only: c0, c1, c2, p5, p25, c1p5, &
field_loc_center, field_loc_NEcorner, &
field_type_scalar, field_type_vector, field_type_angle
use ice_domain_size, only: max_blocks
integer (kind=int_kind) :: &
i, j, iblk, &
ilo,ihi,jlo,jhi ! beginning and end of physical domain
real (kind=dbl_kind) :: &
angle_0, angle_w, angle_s, angle_sw, &
pi, pi2, puny
logical (kind=log_kind), dimension(nx_block,ny_block,max_blocks):: &
out_of_range
real (kind=dbl_kind), dimension(:,:), allocatable :: &
work_g1
type (block) :: &
this_block ! block information for current block
character(len=*), parameter :: subname = '(init_grid2)'
!-----------------------------------------------------------------
! lat, lon, cell widths, angle, land mask
!-----------------------------------------------------------------
call icepack_query_parameters(pi_out=pi, pi2_out=pi2, puny_out=puny)
call icepack_warnings_flush(nu_diag)
if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
file=__FILE__, line=__LINE__)
if (trim(grid_type) == 'displaced_pole' .or. &
trim(grid_type) == 'tripole' .or. &
trim(grid_type) == 'regional' ) then
if (trim(grid_format) == 'nc') then
call popgrid_nc ! read POP grid lengths from nc file
else
call popgrid ! read POP grid lengths directly
endif
#ifdef CESMCOUPLED
elseif (trim(grid_type) == 'latlon') then
call latlongrid ! lat lon grid for sequential CESM (CAM mode)
return
#endif
elseif (trim(grid_type) == 'cpom_grid') then
call cpomgrid ! cpom model orca1 type grid
else
call rectgrid ! regular rectangular grid
endif
!-----------------------------------------------------------------
! T-grid cell and U-grid cell quantities
!-----------------------------------------------------------------
tarea(:,:,:) = c0
!$OMP PARALLEL DO PRIVATE(iblk,i,j,ilo,ihi,jlo,jhi,this_block)
do iblk = 1, nblocks
this_block = get_block(blocks_ice(iblk),iblk)
ilo = this_block%ilo
ihi = this_block%ihi
jlo = this_block%jlo
jhi = this_block%jhi
do j = jlo, jhi
do i = ilo, ihi
tarea(i,j,iblk) = dxt(i,j,iblk)*dyt(i,j,iblk)
uarea(i,j,iblk) = dxu(i,j,iblk)*dyu(i,j,iblk)
if (tarea(i,j,iblk) > c0) then
tarear(i,j,iblk) = c1/tarea(i,j,iblk)
else
tarear(i,j,iblk) = c0 ! possible on boundaries
endif
if (uarea(i,j,iblk) > c0) then
uarear(i,j,iblk) = c1/uarea(i,j,iblk)
else
uarear(i,j,iblk) = c0 ! possible on boundaries
endif
tinyarea(i,j,iblk) = puny*tarea(i,j,iblk)
dxhy(i,j,iblk) = p5*(HTE(i,j,iblk) - HTE(i-1,j,iblk))
dyhx(i,j,iblk) = p5*(HTN(i,j,iblk) - HTN(i,j-1,iblk))
enddo
enddo
do j = jlo, jhi+1
do i = ilo, ihi+1
cyp(i,j,iblk) = (c1p5*HTE(i,j,iblk) - p5*HTE(i-1,j,iblk))
cxp(i,j,iblk) = (c1p5*HTN(i,j,iblk) - p5*HTN(i,j-1,iblk))
! match order of operations in cyp, cxp for tripole grids
cym(i,j,iblk) = -(c1p5*HTE(i-1,j,iblk) - p5*HTE(i,j,iblk))
cxm(i,j,iblk) = -(c1p5*HTN(i,j-1,iblk) - p5*HTN(i,j,iblk))
enddo
enddo
enddo ! iblk
!$OMP END PARALLEL DO
!-----------------------------------------------------------------
! Ghost cell updates
! On the tripole grid, one must be careful with updates of
! quantities that involve a difference of cell lengths.
! For example, dyhx and dxhy are cell-centered vector components.
! Also note that on the tripole grid, cxp and cxm would swap places,
! as would cyp and cym. These quantities are computed only
! in north and east ghost cells (above), not south and west.
!-----------------------------------------------------------------
call ice_timer_start(timer_bound)
call ice_HaloUpdate (tarea, halo_info, &
field_loc_center, field_type_scalar, &
fillValue=c1)
call ice_HaloUpdate (uarea, halo_info, &
field_loc_NEcorner, field_type_scalar, &
fillValue=c1)
call ice_HaloUpdate (tarear, halo_info, &
field_loc_center, field_type_scalar, &
fillValue=c1)
call ice_HaloUpdate (uarear, halo_info, &
field_loc_NEcorner, field_type_scalar, &
fillValue=c1)
call ice_HaloUpdate (tinyarea, halo_info, &
field_loc_center, field_type_scalar, &
fillValue=c1)
call ice_HaloUpdate (dxhy, halo_info, &
field_loc_center, field_type_vector, &
fillValue=c1)
call ice_HaloUpdate (dyhx, halo_info, &
field_loc_center, field_type_vector, &
fillValue=c1)
call ice_timer_stop(timer_bound)
!-----------------------------------------------------------------
! Calculate ANGLET to be compatible with POP ocean model
! First, ensure that -pi <= ANGLE <= pi
!-----------------------------------------------------------------
out_of_range = .false.
where (ANGLE < -pi .or. ANGLE > pi) out_of_range = .true.
if (count(out_of_range) > 0) then
write(nu_diag,*) subname,' angle = ',minval(ANGLE),maxval(ANGLE),count(out_of_range)
call abort_ice (subname//' ANGLE out of expected range', &
file=__FILE__, line=__LINE__)
endif
!-----------------------------------------------------------------
! Compute ANGLE on T-grid
!-----------------------------------------------------------------
if (trim(grid_type) == 'cpom_grid') then
ANGLET(:,:,:) = ANGLE(:,:,:)
else if (.not. (l_readCenter)) then
ANGLET = c0
!$OMP PARALLEL DO PRIVATE(iblk,i,j,ilo,ihi,jlo,jhi,this_block, &
!$OMP angle_0,angle_w,angle_s,angle_sw)
do iblk = 1, nblocks
this_block = get_block(blocks_ice(iblk),iblk)
ilo = this_block%ilo
ihi = this_block%ihi
jlo = this_block%jlo
jhi = this_block%jhi
do j = jlo, jhi
do i = ilo, ihi
angle_0 = ANGLE(i ,j ,iblk) ! w----0
angle_w = ANGLE(i-1,j ,iblk) ! | |
angle_s = ANGLE(i, j-1,iblk) ! | |
angle_sw = ANGLE(i-1,j-1,iblk) ! sw---s
ANGLET(i,j,iblk) = atan2(p25*(sin(angle_0)+ &
sin(angle_w)+ &
sin(angle_s)+ &
sin(angle_sw)),&
p25*(cos(angle_0)+ &
cos(angle_w)+ &
cos(angle_s)+ &
cos(angle_sw)))
enddo
enddo
enddo
!$OMP END PARALLEL DO
endif ! cpom_grid
if (trim(grid_type) == 'regional' .and. &
(.not. (l_readCenter))) then
! for W boundary extrapolate from interior
!$OMP PARALLEL DO PRIVATE(iblk,i,j,ilo,ihi,jlo,jhi,this_block)
do iblk = 1, nblocks
this_block = get_block(blocks_ice(iblk),iblk)
ilo = this_block%ilo
ihi = this_block%ihi
jlo = this_block%jlo
jhi = this_block%jhi
i = ilo
if (this_block%i_glob(i) == 1) then
do j = jlo, jhi
ANGLET(i,j,iblk) = c2*ANGLET(i+1,j,iblk)-ANGLET(i+2,j,iblk)
enddo
endif
enddo
!$OMP END PARALLEL DO
endif ! regional
call ice_timer_start(timer_bound)
call ice_HaloUpdate (ANGLET, halo_info, &
field_loc_center, field_type_angle, &
fillValue=c1)
call ice_timer_stop(timer_bound)
call makemask ! velocity mask, hemisphere masks
if (.not. (l_readCenter)) then
call Tlatlon ! get lat, lon on the T grid
endif
!-----------------------------------------------------------------
! bathymetry
!-----------------------------------------------------------------
if (trim(bathymetry_format) == 'default') then
call get_bathymetry
elseif (trim(bathymetry_format) == 'pop') then
call get_bathymetry_popfile
else
call abort_ice(subname//'ERROR: bathymetry_format value must be default or pop', &
file=__FILE__, line=__LINE__)
endif
!----------------------------------------------------------------
! Corner coordinates for CF compliant history files
!----------------------------------------------------------------
call gridbox_corners
!-----------------------------------------------------------------
! Compute global index (used for unpacking messages from coupler)
!-----------------------------------------------------------------
if (my_task==master_task) then
allocate(work_g1(nx_global,ny_global))
do j=1,ny_global
do i=1,nx_global
work_g1(i,j) = real((j-1)*nx_global + i,kind=dbl_kind)
enddo
enddo
else
allocate(work_g1(1,1)) ! to save memory
endif
call scatter_global(rndex_global, work_g1, &
master_task, distrb_info, &
field_loc_center, field_type_scalar)
deallocate(work_g1)
end subroutine init_grid2
!=======================================================================
! POP displaced pole grid and land mask (or tripole).
! Grid record number, field and units are: \\
! (1) ULAT (radians) \\
! (2) ULON (radians) \\
! (3) HTN (cm) \\
! (4) HTE (cm) \\
! (5) HUS (cm) \\
! (6) HUW (cm) \\
! (7) ANGLE (radians)
!
! Land mask record number and field is (1) KMT.
!
! author: Elizabeth C. Hunke, LANL
subroutine popgrid
use ice_blocks, only: nx_block, ny_block
use ice_constants, only: c0, c1, &
field_loc_center, field_loc_NEcorner, &
field_type_scalar, field_type_angle
use ice_domain_size, only: max_blocks
integer (kind=int_kind) :: &
i, j, iblk, &
ilo,ihi,jlo,jhi ! beginning and end of physical domain
logical (kind=log_kind) :: diag
real (kind=dbl_kind), dimension(:,:), allocatable :: &
work_g1
real (kind=dbl_kind), dimension (nx_block,ny_block,max_blocks) :: &
work1
type (block) :: &
this_block ! block information for current block
character(len=*), parameter :: subname = '(popgrid)'
call ice_open(nu_grid,grid_file,64)
call ice_open(nu_kmt,kmt_file,32)
diag = .true. ! write diagnostic info
!-----------------------------------------------------------------
! topography
!-----------------------------------------------------------------
call ice_read(nu_kmt,1,work1,'ida4',diag, &
field_loc=field_loc_center, &
field_type=field_type_scalar)
hm (:,:,:) = c0
kmt(:,:,:) = c0
!$OMP PARALLEL DO PRIVATE(iblk,i,j,ilo,ihi,jlo,jhi,this_block)
do iblk = 1, nblocks
this_block = get_block(blocks_ice(iblk),iblk)
ilo = this_block%ilo
ihi = this_block%ihi
jlo = this_block%jlo
jhi = this_block%jhi
do j = jlo, jhi
do i = ilo, ihi
kmt(i,j,iblk) = work1(i,j,iblk)
if (kmt(i,j,iblk) >= c1) hm(i,j,iblk) = c1
enddo
enddo
enddo
!$OMP END PARALLEL DO
!-----------------------------------------------------------------
! lat, lon, angle
!-----------------------------------------------------------------
allocate(work_g1(nx_global,ny_global))
call ice_read_global(nu_grid,1,work_g1,'rda8',.true.) ! ULAT
call gridbox_verts(work_g1,latt_bounds)
call scatter_global(ULAT, work_g1, master_task, distrb_info, &
field_loc_NEcorner, field_type_scalar)
call ice_HaloExtrapolate(ULAT, distrb_info, &
ew_boundary_type, ns_boundary_type)
call ice_read_global(nu_grid,2,work_g1,'rda8',.true.) ! ULON
call gridbox_verts(work_g1,lont_bounds)
call scatter_global(ULON, work_g1, master_task, distrb_info, &
field_loc_NEcorner, field_type_scalar)
call ice_HaloExtrapolate(ULON, distrb_info, &
ew_boundary_type, ns_boundary_type)
call ice_read_global(nu_grid,7,work_g1,'rda8',.true.) ! ANGLE
call scatter_global(ANGLE, work_g1, master_task, distrb_info, &
field_loc_NEcorner, field_type_angle)
!-----------------------------------------------------------------
! cell dimensions
! calculate derived quantities from global arrays to preserve
! information on boundaries
!-----------------------------------------------------------------
call ice_read_global(nu_grid,3,work_g1,'rda8',.true.) ! HTN
call primary_grid_lengths_HTN(work_g1) ! dxu, dxt
call ice_read_global(nu_grid,4,work_g1,'rda8',.true.) ! HTE
call primary_grid_lengths_HTE(work_g1) ! dyu, dyt
deallocate(work_g1)
if (my_task == master_task) then
close (nu_grid)
close (nu_kmt)
endif
end subroutine popgrid
!=======================================================================
! POP displaced pole grid and land mask.
! Grid record number, field and units are: \\
! (1) ULAT (radians) \\
! (2) ULON (radians) \\
! (3) HTN (cm) \\
! (4) HTE (cm) \\
! (5) HUS (cm) \\
! (6) HUW (cm) \\
! (7) ANGLE (radians)
!
! Land mask record number and field is (1) KMT.
!
! author: Elizabeth C. Hunke, LANL
! Revised for netcdf input: Ann Keen, Met Office, May 2007
subroutine popgrid_nc
use ice_blocks, only: nx_block, ny_block
use ice_constants, only: c0, c1, &
field_loc_center, field_loc_NEcorner, &
field_type_scalar, field_type_angle
use ice_domain_size, only: max_blocks
#ifdef USE_NETCDF
use netcdf
#endif
integer (kind=int_kind) :: &
i, j, iblk, &
ilo,ihi,jlo,jhi, & ! beginning and end of physical domain
fid_grid, & ! file id for netCDF grid file
fid_kmt ! file id for netCDF kmt file
logical (kind=log_kind) :: diag
character (char_len) :: &
fieldname ! field name in netCDF file
real (kind=dbl_kind) :: &
pi
real (kind=dbl_kind), dimension(:,:), allocatable :: &
work_g1
real (kind=dbl_kind), dimension (nx_block,ny_block,max_blocks) :: &
work1
type (block) :: &
this_block ! block information for current block
integer(kind=int_kind) :: &
varid
integer (kind=int_kind) :: &
status ! status flag
character(len=*), parameter :: subname = '(popgrid_nc)'
#ifdef USE_NETCDF
call icepack_query_parameters(pi_out=pi)
call icepack_warnings_flush(nu_diag)
if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
file=__FILE__, line=__LINE__)
call ice_open_nc(grid_file,fid_grid)
call ice_open_nc(kmt_file,fid_kmt)
diag = .true. ! write diagnostic info
l_readCenter = .false.
!-----------------------------------------------------------------
! topography
!-----------------------------------------------------------------
fieldname='kmt'
call ice_read_nc(fid_kmt,1,fieldname,work1,diag, &
field_loc=field_loc_center, &
field_type=field_type_scalar)
hm (:,:,:) = c0
kmt(:,:,:) = c0
!$OMP PARALLEL DO PRIVATE(iblk,i,j,ilo,ihi,jlo,jhi,this_block)
do iblk = 1, nblocks
this_block = get_block(blocks_ice(iblk),iblk)
ilo = this_block%ilo
ihi = this_block%ihi
jlo = this_block%jlo
jhi = this_block%jhi
do j = jlo, jhi
do i = ilo, ihi
kmt(i,j,iblk) = work1(i,j,iblk)
if (kmt(i,j,iblk) >= c1) hm(i,j,iblk) = c1
enddo
enddo
enddo
!$OMP END PARALLEL DO
!-----------------------------------------------------------------
! lat, lon, angle
!-----------------------------------------------------------------
allocate(work_g1(nx_global,ny_global))
fieldname='ulat'
call ice_read_global_nc(fid_grid,1,fieldname,work_g1,diag) ! ULAT
call gridbox_verts(work_g1,latt_bounds)
call scatter_global(ULAT, work_g1, master_task, distrb_info, &
field_loc_NEcorner, field_type_scalar)
call ice_HaloExtrapolate(ULAT, distrb_info, &
ew_boundary_type, ns_boundary_type)
fieldname='ulon'
call ice_read_global_nc(fid_grid,1,fieldname,work_g1,diag) ! ULON
call gridbox_verts(work_g1,lont_bounds)
call scatter_global(ULON, work_g1, master_task, distrb_info, &
field_loc_NEcorner, field_type_scalar)
call ice_HaloExtrapolate(ULON, distrb_info, &
ew_boundary_type, ns_boundary_type)
fieldname='angle'
call ice_read_global_nc(fid_grid,1,fieldname,work_g1,diag) ! ANGLE
call scatter_global(ANGLE, work_g1, master_task, distrb_info, &
field_loc_NEcorner, field_type_angle)
! fix ANGLE: roundoff error due to single precision
where (ANGLE > pi) ANGLE = pi
where (ANGLE < -pi) ANGLE = -pi
! if grid file includes anglet then read instead
fieldname='anglet'
if (my_task == master_task) then
status = nf90_inq_varid(fid_grid, trim(fieldname) , varid)
if (status /= nf90_noerr) then
write(nu_diag,*) subname//' CICE will calculate angleT, TLON and TLAT'
else
write(nu_diag,*) subname//' angleT, TLON and TLAT is read from grid file'
l_readCenter = .true.
endif
endif
call broadcast_scalar(l_readCenter,master_task)
if (l_readCenter) then
call ice_read_global_nc(fid_grid,1,fieldname,work_g1,diag)
call scatter_global(ANGLET, work_g1, master_task, distrb_info, &
field_loc_center, field_type_angle)
where (ANGLET > pi) ANGLET = pi
where (ANGLET < -pi) ANGLET = -pi
fieldname="tlon"
call ice_read_global_nc(fid_grid,1,fieldname,work_g1,diag)
call scatter_global(TLON, work_g1, master_task, distrb_info, &
field_loc_center, field_type_scalar)
fieldname="tlat"
call ice_read_global_nc(fid_grid,1,fieldname,work_g1,diag)
call scatter_global(TLAT, work_g1, master_task, distrb_info, &
field_loc_center, field_type_scalar)
endif
!-----------------------------------------------------------------
! cell dimensions
! calculate derived quantities from global arrays to preserve
! information on boundaries
!-----------------------------------------------------------------
fieldname='htn'
call ice_read_global_nc(fid_grid,1,fieldname,work_g1,diag) ! HTN
call primary_grid_lengths_HTN(work_g1) ! dxu, dxt
fieldname='hte'
call ice_read_global_nc(fid_grid,1,fieldname,work_g1,diag) ! HTE
call primary_grid_lengths_HTE(work_g1) ! dyu, dyt
deallocate(work_g1)
if (my_task == master_task) then
call ice_close_nc(fid_grid)
call ice_close_nc(fid_kmt)
endif
#else
call abort_ice(subname//'ERROR: USE_NETCDF cpp not defined', &
file=__FILE__, line=__LINE__)
#endif
end subroutine popgrid_nc
#ifdef CESMCOUPLED
!=======================================================================
! Read in kmt file that matches CAM lat-lon grid and has single column
! functionality
! author: Mariana Vertenstein
! 2007: Elizabeth Hunke upgraded to netcdf90 and cice ncdf calls
subroutine latlongrid
! use ice_boundary
use ice_domain_size
use ice_scam, only : scmlat, scmlon, single_column
use ice_constants, only: c0, c1, p5, p25, &
field_loc_center, field_type_scalar, radius
#ifdef USE_NETCDF
use netcdf
#endif
integer (kind=int_kind) :: &
i, j, iblk
integer (kind=int_kind) :: &
ni, nj, ncid, dimid, varid, ier
type (block) :: &
this_block ! block information for current block
integer (kind=int_kind) :: &
ilo,ihi,jlo,jhi ! beginning and end of physical domain
real (kind=dbl_kind) :: &
closelat, & ! Single-column latitude value
closelon, & ! Single-column longitude value
closelatidx, & ! Single-column latitude index to retrieve
closelonidx ! Single-column longitude index to retrieve
integer (kind=int_kind) :: &
start(2), & ! Start index to read in
count(2) ! Number of points to read in
integer (kind=int_kind) :: &
start3(3), & ! Start index to read in
count3(3) ! Number of points to read in
integer (kind=int_kind) :: &
status ! status flag
real (kind=dbl_kind), allocatable :: &
lats(:),lons(:),pos_lons(:), glob_grid(:,:) ! temporaries
real (kind=dbl_kind) :: &
pos_scmlon,& ! temporary
pi, &
puny, &
scamdata ! temporary
character(len=*), parameter :: subname = '(lonlatgrid)'
#ifdef USE_NETCDF
!-----------------------------------------------------------------
! - kmt file is actually clm fractional land file
! - Determine consistency of dimensions
! - Read in lon/lat centers in degrees from kmt file
! - Read in ocean from "kmt" file (1 for ocean, 0 for land)
!-----------------------------------------------------------------
call icepack_query_parameters(pi_out=pi, puny_out=puny)
call icepack_warnings_flush(nu_diag)
if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
file=__FILE__, line=__LINE__)
! Determine dimension of domain file and check for consistency
if (my_task == master_task) then
call ice_open_nc(kmt_file, ncid)
status = nf90_inq_dimid (ncid, 'ni', dimid)
status = nf90_inquire_dimension(ncid, dimid, len=ni)
status = nf90_inq_dimid (ncid, 'nj', dimid)
status = nf90_inquire_dimension(ncid, dimid, len=nj)
end if
! Determine start/count to read in for either single column or global lat-lon grid
! If single_column, then assume that only master_task is used since there is only one task
if (single_column) then
! Check for consistency
if (my_task == master_task) then
if ((nx_global /= 1).or. (ny_global /= 1)) then
write(nu_diag,*) 'Because you have selected the column model flag'
write(nu_diag,*) 'Please set nx_global=ny_global=1 in file'
write(nu_diag,*) 'ice_domain_size.F and recompile'
call abort_ice (subname//'ERROR: check nx_global, ny_global')
endif
end if
! Read in domain file for single column
allocate(lats(nj))
allocate(lons(ni))
allocate(pos_lons(ni))
allocate(glob_grid(ni,nj))
start3=(/1,1,1/)
count3=(/ni,nj,1/)
status = nf90_inq_varid(ncid, 'xc' , varid)
if (status /= nf90_noerr) call abort_ice (subname//' inq_varid xc')
status = nf90_get_var(ncid, varid, glob_grid, start3, count3)
if (status /= nf90_noerr) call abort_ice (subname//' get_var xc')
do i = 1,ni
lons(i) = glob_grid(i,1)
end do
status = nf90_inq_varid(ncid, 'yc' , varid)
if (status /= nf90_noerr) call abort_ice (subname//' inq_varid yc')