-
Notifications
You must be signed in to change notification settings - Fork 118
/
fv_grid_tools.F90
2408 lines (2045 loc) · 86 KB
/
fv_grid_tools.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
!***********************************************************************
!* GNU Lesser General Public License
!*
!* This file is part of the FV3 dynamical core.
!*
!* The FV3 dynamical core is free software: you can redistribute it
!* and/or modify it under the terms of the
!* GNU Lesser General Public License as published by the
!* Free Software Foundation, either version 3 of the License, or
!* (at your option) any later version.
!*
!* The FV3 dynamical core is distributed in the hope that it will be
!* useful, but WITHOUT ANYWARRANTY; without even the implied warranty
!* of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
!* See the GNU General Public License for more details.
!*
!* You should have received a copy of the GNU Lesser General Public
!* License along with the FV3 dynamical core.
!* If not, see <http://www.gnu.org/licenses/>.
!***********************************************************************
module fv_grid_tools_mod
use constants_mod, only: grav, omega, pi=>pi_8, cnst_radius=>radius
use fv_arrays_mod, only: fv_atmos_type, fv_grid_type, fv_grid_bounds_type, R_GRID
use fv_grid_utils_mod, only: gnomonic_grids, great_circle_dist, &
mid_pt_sphere, spherical_angle, &
cell_center2, get_area, inner_prod, fill_ghost, &
direct_transform, cube_transform, dist2side_latlon, &
spherical_linear_interpolation, big_number
use fv_timing_mod, only: timing_on, timing_off
use fv_mp_mod, only: is_master, fill_corners, XDir, YDir
use fv_mp_mod, only: mp_gather, mp_bcst, mp_reduce_max, mp_stop, grids_master_procs
use sorted_index_mod, only: sorted_inta, sorted_intb
use mpp_mod, only: mpp_error, FATAL, get_unit, mpp_chksum, mpp_pe, stdout, &
mpp_send, mpp_recv, mpp_sync_self, EVENT_RECV, mpp_npes, &
mpp_sum, mpp_max, mpp_min, mpp_root_pe, mpp_broadcast
use mpp_domains_mod, only: mpp_update_domains, mpp_get_boundary, &
mpp_get_ntile_count, mpp_get_pelist, &
mpp_get_compute_domains, mpp_global_field, &
mpp_get_data_domain, mpp_get_compute_domain, &
mpp_get_global_domain, mpp_global_sum, mpp_global_max, mpp_global_min
use mpp_domains_mod, only: domain2d
use mpp_io_mod, only: mpp_get_att_value
use mpp_parameter_mod, only: AGRID_PARAM=>AGRID, &
DGRID_NE_PARAM=>DGRID_NE, &
CGRID_NE_PARAM=>CGRID_NE, &
CGRID_SW_PARAM=>CGRID_SW, &
BGRID_NE_PARAM=>BGRID_NE, &
BGRID_SW_PARAM=>BGRID_SW, &
SCALAR_PAIR, &
CORNER, CENTER, XUPDATE
use fms_mod, only: get_mosaic_tile_grid
use fms_io_mod, only: file_exist, field_exist, read_data, &
get_global_att_value, get_var_att_value
use mosaic_mod, only : get_mosaic_ntiles
use mpp_mod, only: mpp_transmit, mpp_recv
implicit none
private
#include <netcdf.inc>
real(kind=R_GRID), parameter:: radius = cnst_radius
real(kind=R_GRID) , parameter:: todeg = 180.0d0/pi ! convert to degrees
real(kind=R_GRID) , parameter:: torad = pi/180.0d0 ! convert to radians
real(kind=R_GRID) , parameter:: missing = 1.d25
real(kind=R_GRID) :: csFac
logical, parameter :: debug_message_size = .false.
logical :: write_grid_char_file = .false.
public :: todeg, missing, init_grid, spherical_to_cartesian
contains
subroutine read_grid(Atm, grid_file, ndims, nregions, ng)
! read_grid :: read grid from mosaic grid file.
type(fv_atmos_type), intent(inout), target :: Atm
character(len=*), intent(IN) :: grid_file
integer, intent(IN) :: ndims
integer, intent(IN) :: nregions
integer, intent(IN) :: ng
real, allocatable, dimension(:,:) :: tmpx, tmpy
real(kind=R_GRID), pointer, dimension(:,:,:) :: grid
character(len=128) :: units = ""
character(len=256) :: atm_mosaic, atm_hgrid, grid_form
character(len=1024) :: attvalue
integer :: ntiles, i, j, stdunit
integer :: isc2, iec2, jsc2, jec2
integer :: start(4), nread(4)
integer :: is, ie, js, je
integer :: isd, ied, jsd, jed
integer,save :: halo=3 ! for regional domain external tools
is = Atm%bd%is
ie = Atm%bd%ie
js = Atm%bd%js
je = Atm%bd%je
isd = Atm%bd%isd
ied = Atm%bd%ied
jsd = Atm%bd%jsd
jed = Atm%bd%jed
grid => Atm%gridstruct%grid_64
if(.not. file_exist(grid_file)) call mpp_error(FATAL, 'fv_grid_tools(read_grid): file '// &
trim(grid_file)//' does not exist')
!--- make sure the grid file is mosaic file.
if( field_exist(grid_file, 'atm_mosaic_file') .OR. field_exist(grid_file, 'gridfiles') ) then
stdunit = stdout()
write(stdunit,*) '==>Note from fv_grid_tools_mod(read_grid): read atmosphere grid from mosaic version grid'
else
call mpp_error(FATAL, 'fv_grid_tools(read_grid): neither atm_mosaic_file nor gridfiles exists in file ' &
//trim(grid_file))
endif
if(field_exist(grid_file, 'atm_mosaic_file')) then
call read_data(grid_file, "atm_mosaic_file", atm_mosaic)
atm_mosaic = "INPUT/"//trim(atm_mosaic)
else
atm_mosaic = trim(grid_file)
endif
call get_mosaic_tile_grid(atm_hgrid, atm_mosaic, Atm%domain)
grid_form = "none"
if( get_global_att_value(atm_hgrid, "history", attvalue) ) then
if( index(attvalue, "gnomonic_ed") > 0) grid_form = "gnomonic_ed"
endif
if(grid_form .NE. "gnomonic_ed") call mpp_error(FATAL, &
"fv_grid_tools(read_grid): the grid should be 'gnomonic_ed' when reading from grid file, contact developer")
!FIXME: Doesn't work for a nested grid
ntiles = get_mosaic_ntiles(atm_mosaic)
if( .not. Atm%gridstruct%bounded_domain) then !<-- The regional setup has only 1 tile so do not shutdown in that case.
if(ntiles .NE. 6) call mpp_error(FATAL, &
'fv_grid_tools(read_grid): ntiles should be 6 in mosaic file '//trim(atm_mosaic) )
if(nregions .NE. 6) call mpp_error(FATAL, &
'fv_grid_tools(read_grid): nregions should be 6 when reading from mosaic file '//trim(grid_file) )
endif
call get_var_att_value(atm_hgrid, 'x', 'units', units)
!--- get the geographical coordinates of super-grid.
isc2 = 2*is-1; iec2 = 2*ie+1
jsc2 = 2*js-1; jec2 = 2*je+1
if( Atm%gridstruct%bounded_domain ) then
isc2 = 2*(isd+halo)-1; iec2 = 2*(ied+1+halo)-1 ! For the regional domain the cell corner locations must be transferred
jsc2 = 2*(jsd+halo)-1; jec2 = 2*(jed+1+halo)-1 ! from the entire supergrid to the compute grid, including the halo region.
endif
allocate(tmpx(isc2:iec2, jsc2:jec2) )
allocate(tmpy(isc2:iec2, jsc2:jec2) )
start = 1; nread = 1
start(1) = isc2; nread(1) = iec2 - isc2 + 1
start(2) = jsc2; nread(2) = jec2 - jsc2 + 1
call read_data(atm_hgrid, 'x', tmpx, start, nread, no_domain=.TRUE.) !<-- tmpx (lon, deg east) is on the supergrid
call read_data(atm_hgrid, 'y', tmpy, start, nread, no_domain=.TRUE.) !<-- tmpy (lat, deg) is on the supergrid
!--- geographic grid at cell corner
grid(isd: is-1, jsd:js-1,1:ndims)=0.
grid(isd: is-1, je+2:jed+1,1:ndims)=0.
grid(ie+2:ied+1,jsd:js-1,1:ndims)=0.
grid(ie+2:ied+1,je+2:jed+1,1:ndims)=0.
if(len_trim(units) < 6) call mpp_error(FATAL, &
"fv_grid_tools_mod(read_grid): the length of units must be no less than 6")
if(units(1:6) == 'degree') then
if( .not. Atm%gridstruct%bounded_domain) then
do j = js, je+1
do i = is, ie+1
grid(i,j,1) = tmpx(2*i-1,2*j-1)*pi/180.
grid(i,j,2) = tmpy(2*i-1,2*j-1)*pi/180.
enddo
enddo
else
!
!*** In the regional case the halo surrounding the domain was included in the read.
!*** Transfer the compute and halo regions to the compute grid.
!
do j = jsd, jed+1
do i = isd, ied+1
grid(i,j,1) = tmpx(2*i+halo+2,2*j+halo+2)*pi/180.
grid(i,j,2) = tmpy(2*i+halo+2,2*j+halo+2)*pi/180.
enddo
enddo
endif
else if(units(1:6) == 'radian') then
do j = js, je+1
do i = is, ie+1
grid(i,j,1) = tmpx(2*i-1,2*j-1)
grid(i,j,2) = tmpy(2*i-1,2*j-1)
enddo
enddo
else
print*, 'units is ' , trim(units), len_trim(units), mpp_pe()
call mpp_error(FATAL, 'fv_grid_tools_mod(read_grid): units must start with degree or radian')
endif
deallocate(tmpx, tmpy)
nullify(grid)
end subroutine read_grid
!#################################################################################
subroutine get_symmetry(data_in, data_out, ishift, jshift, npes_x, npes_y, domain, tile, npx_g, bd)
type(fv_grid_bounds_type), intent(IN) :: bd
integer, intent(in) :: ishift, jshift, npes_x, npes_y
real(kind=R_GRID), dimension(bd%is:bd%ie+ishift, bd%js:bd%je+jshift ), intent(in) :: data_in
real(kind=R_GRID), dimension(bd%is:bd%ie+jshift, bd%js:bd%je+ishift ), intent(out) :: data_out
real(kind=R_GRID), dimension(:), allocatable :: send_buffer
real(kind=R_GRID), dimension(:), allocatable :: recv_buffer
integer, dimension(:), allocatable :: is_recv, ie_recv, js_recv, je_recv, pe_recv
integer, dimension(:), allocatable :: is_send, ie_send, js_send, je_send, pe_send
integer, dimension(:), allocatable :: isl, iel, jsl, jel, pelist, msg1, msg2
integer :: msgsize, pos, ntiles, npes_per_tile, npes
integer :: send_buf_size, recv_buf_size, buffer_pos
integer :: is0, ie0, js0, je0
integer :: is1, ie1, js1, je1
integer :: is2, ie2, js2, je2
integer :: i, j, p, nrecv, nsend, tile_you, is3, ie3, nlist
integer :: start_pe, ipos, jpos, from_pe, to_pe
type(domain2d) :: domain
integer :: tile, npx_g
integer :: is, ie, js, je
integer :: isd, ied, jsd, jed
is = bd%is
ie = bd%ie
js = bd%js
je = bd%je
isd = bd%isd
ied = bd%ied
jsd = bd%jsd
jed = bd%jed
!--- This routine will be called only for cubic sphere grid. so 6 tiles will be assumed
!--- also number of processors on each tile will be the same.
ntiles = mpp_get_ntile_count(domain)
npes = mpp_npes()
if(ntiles .NE. 6 ) call mpp_error(FATAL, 'fv_grid_tools(get_symmetry): ntiles should be 6 ')
if(mod(npes,ntiles) /= 0) call mpp_error(FATAL, 'fv_grid_tools(get_symmetry): npes should be divided by ntiles')
npes_per_tile = npes/ntiles
! if(npes_x == npes_y) then ! even, simple communication
if(npes_x == npes_y .AND. mod(npx_g-1,npes_x) == 0 ) then ! even,
msgsize = (ie-is+1+jshift)*(je-js+1+ishift)
pos = mod((mpp_pe()-mpp_root_pe()), npes_x*npes_y)
start_pe = mpp_pe() - pos
ipos = mod(pos, npes_x)
jpos = pos/npes_x
from_pe = start_pe + ipos*npes_x + jpos
to_pe = from_pe
allocate(recv_buffer(msgsize))
call mpp_recv(recv_buffer(1), glen=msgsize, from_pe=from_pe, block=.FALSE. )
pos = 0
allocate(send_buffer(msgsize))
do j = js, je+jshift
do i = is, ie+ishift
pos = pos + 1
send_buffer(pos) = data_in(i,j)
enddo
enddo
call mpp_send(send_buffer(1), plen=msgsize, to_pe=to_pe)
call mpp_sync_self(check=EVENT_RECV) ! To ensure recv is completed.
!--unpack buffer
pos = 0
do i = is, ie+jshift
do j = js, je+ishift
pos = pos + 1
data_out(i,j) = recv_buffer(pos)
enddo
enddo
call mpp_sync_self()
deallocate(send_buffer, recv_buffer)
else
allocate(is_recv(0:npes_per_tile-1), ie_recv(0:npes_per_tile-1))
allocate(js_recv(0:npes_per_tile-1), je_recv(0:npes_per_tile-1))
allocate(is_send(0:npes_per_tile-1), ie_send(0:npes_per_tile-1))
allocate(js_send(0:npes_per_tile-1), je_send(0:npes_per_tile-1))
allocate(pe_send(0:npes_per_tile-1), pe_recv(0:npes_per_tile-1))
if(debug_message_size) then
allocate(msg1 (0:npes_per_tile-1), msg2 (0:npes_per_tile-1))
msg1 = 0
msg2 = 0
endif
allocate(pelist(0:npes-1))
call mpp_get_pelist(domain, pelist)
allocate(isl(0:npes-1), iel(0:npes-1), jsl(0:npes-1), jel(0:npes-1) )
call mpp_get_compute_domains(domain, xbegin=isl, xend=iel, ybegin=jsl, yend=jel)
!--- pre-post receiving
buffer_pos = 0
nrecv = 0
nsend = 0
recv_buf_size = 0
!--- first set up the receiving index
nlist = 0
do p = 0, npes-1
tile_you = p/(npes_x*npes_y) + 1
if(tile_you .NE. tile) cycle
!--- my index for data_out after rotation
is1 = js; ie1 = je + ishift;
js1 = is; je1 = ie + jshift;
!--- your index for data_out
is2 = isl(p); ie2 = iel(p) + ishift;
js2 = jsl(p); je2 = jel(p) + jshift;
is0 = max(is1,is2); ie0 = min(ie1,ie2)
js0 = max(js1,js2); je0 = min(je1,je2)
msgsize = 0
if(ie0 .GE. is0 .AND. je0 .GE. js0) then
msgsize = (ie0-is0+1)*(je0-js0+1)
recv_buf_size = recv_buf_size + msgsize
pe_recv(nrecv) = pelist(p)
!--- need to rotate back the index
is_recv(nrecv) = js0; ie_recv(nrecv) = je0
js_recv(nrecv) = is0; je_recv(nrecv) = ie0
nrecv = nrecv+1
endif
if(debug_message_size) then
msg1(nlist) = msgsize
call mpp_recv(msg2(nlist), glen=1, from_pe=pelist(p), block=.FALSE. )
nlist = nlist + 1
endif
enddo
!--- Then setup the sending index.
send_buf_size = 0
do p = 0, npes-1
tile_you = p/(npes_x*npes_y) + 1
if(tile_you .NE. tile) cycle
!--- my index on data_in
is1 = is; ie1 = ie + ishift;
js1 = js; je1 = je + jshift;
!--- your index on data_out after rotate
is2 = jsl(p); ie2 = jel(p) + ishift;
js2 = isl(p); je2 = iel(p) + jshift;
is0 = max(is1,is2); ie0 = min(ie1,ie2)
js0 = max(js1,js2); je0 = min(je1,je2)
msgsize = 0
if(ie0 .GE. is0 .AND. je0 .GE. js0 )then
msgsize = (ie0-is0+1)*(je0-js0+1)
send_buf_size = send_buf_size + msgsize
pe_send(nsend) = pelist(p)
is_send(nsend) = is0; ie_send(nsend) = ie0
js_send(nsend) = js0; je_send(nsend) = je0
nsend = nsend+1
endif
IF(debug_message_size) call mpp_send(msgsize, plen=1, to_pe=pelist(p) )
enddo
!--- check to make sure send and recv size match.
if(debug_message_size) then
call mpp_sync_self(check=EVENT_RECV) ! To ensure recv is completed.
do p = 0, nlist-1
if(msg1(p) .NE. msg2(p)) then
call mpp_error(FATAL, "fv_grid_tools_mod(get_symmetry): mismatch on send and recv size")
endif
enddo
call mpp_sync_self()
deallocate(msg1, msg2)
endif
!--- pre-post data
allocate(recv_buffer(recv_buf_size))
buffer_pos = 0
do p = 0, nrecv-1
is0 = is_recv(p); ie0 = ie_recv(p)
js0 = js_recv(p); je0 = je_recv(p)
msgsize = (ie0-is0+1)*(je0-js0+1)
call mpp_recv(recv_buffer(buffer_pos+1), glen=msgsize, from_pe=pe_recv(p), block=.FALSE. )
buffer_pos = buffer_pos + msgsize
enddo
!--- send the data
buffer_pos = 0
allocate(send_buffer(send_buf_size))
do p = 0, nsend-1
is0 = is_send(p); ie0 = ie_send(p)
js0 = js_send(p); je0 = je_send(p)
msgsize = (ie0-is0+1)*(je0-js0+1)
pos = buffer_pos
do j = js0, je0
do i = is0, ie0
pos = pos+1
send_buffer(pos) = data_in(i,j)
enddo
enddo
call mpp_send(send_buffer(buffer_pos+1), plen=msgsize, to_pe=pe_send(p) )
buffer_pos = buffer_pos + msgsize
enddo
call mpp_sync_self(check=EVENT_RECV) ! To ensure recv is completed.
!--- unpack buffer
pos = 0
do p = 0, nrecv-1
is0 = is_recv(p); ie0 = ie_recv(p)
js0 = js_recv(p); je0 = je_recv(p)
do i = is0, ie0
do j = js0, je0
pos = pos + 1
data_out(i,j) = recv_buffer(pos)
enddo
enddo
enddo
call mpp_sync_self()
deallocate(isl, iel, jsl, jel, pelist)
deallocate(is_recv, ie_recv, js_recv, je_recv, pe_recv)
deallocate(is_send, ie_send, js_send, je_send, pe_send)
deallocate(recv_buffer, send_buffer)
endif
end subroutine get_symmetry
subroutine init_grid(Atm, grid_name, grid_file, npx, npy, npz, ndims, nregions, ng, tile_coarse)
! init_grid :: read grid from input file and setup grid descriptors
!--------------------------------------------------------
type(fv_atmos_type), intent(inout), target :: Atm
character(len=80), intent(IN) :: grid_name
character(len=120),intent(IN) :: grid_file
integer, intent(IN) :: npx, npy, npz
integer, intent(IN) :: ndims
integer, intent(IN) :: nregions
integer, intent(IN) :: ng
integer, intent(IN) :: tile_coarse(:)
!--------------------------------------------------------
real(kind=R_GRID) :: xs(npx,npy)
real(kind=R_GRID) :: ys(npx,npy)
real(kind=R_GRID) :: dp, dl
real(kind=R_GRID) :: x1,x2,y1,y2,z1,z2
integer :: i,j,k,n,nreg
integer :: fileLun
real(kind=R_GRID) :: p1(3), p2(3), p3(3), p4(3)
real(kind=R_GRID) :: dist,dist1,dist2, pa(2), pa1(2), pa2(2), pb(2)
real(kind=R_GRID) :: pt(3), pt1(3), pt2(3), pt3(3)
real(kind=R_GRID) :: ee1(3), ee2(3)
real(kind=R_GRID) :: angN,angM,angAV,ang
real(kind=R_GRID) :: aspN,aspM,aspAV,asp
real(kind=R_GRID) :: dxN, dxM, dxAV
real(kind=R_GRID) :: dx_local, dy_local
real(kind=R_GRID) :: vec1(3), vec2(3), vec3(3), vec4(3)
real(kind=R_GRID) :: vecAvg(3), vec3a(3), vec3b(3), vec4a(3), vec4b(3)
real(kind=R_GRID) :: xyz1(3), xyz2(3)
! real(kind=R_GRID) :: grid_global(1-ng:npx +ng,1-ng:npy +ng,ndims,1:nregions)
integer :: ios, ip, jp
integer :: igrid
integer :: tmplun
character(len=80) :: tmpFile
real(kind=R_GRID), dimension(Atm%bd%is:Atm%bd%ie) :: sbuffer, nbuffer
real(kind=R_GRID), dimension(Atm%bd%js:Atm%bd%je) :: wbuffer, ebuffer
real(kind=R_GRID), pointer, dimension(:,:,:) :: agrid, grid
real(kind=R_GRID), pointer, dimension(:,:) :: area, area_c
real(kind=R_GRID), pointer, dimension(:,:) :: sina, cosa, dx, dy, dxc, dyc, dxa, dya
real, pointer, dimension(:,:,:) :: e1, e2
real, pointer, dimension(:,:) :: rarea, rarea_c
real, pointer, dimension(:,:) :: rdx, rdy, rdxc, rdyc, rdxa, rdya
integer, pointer, dimension(:,:,:) :: iinta, jinta, iintb, jintb
real(kind=R_GRID), pointer, dimension(:,:,:,:) :: grid_global
integer, pointer :: npx_g, npy_g, ntiles_g, tile
logical, pointer :: sw_corner, se_corner, ne_corner, nw_corner
logical, pointer :: latlon, cubed_sphere, have_south_pole, have_north_pole, stretched_grid
type(domain2d), pointer :: domain
integer :: is, ie, js, je
integer :: isd, ied, jsd, jed
integer :: istart, iend, jstart, jend
is = Atm%bd%is
ie = Atm%bd%ie
js = Atm%bd%js
je = Atm%bd%je
isd = Atm%bd%isd
ied = Atm%bd%ied
jsd = Atm%bd%jsd
jed = Atm%bd%jed
!!! Associate pointers
agrid => Atm%gridstruct%agrid_64
grid => Atm%gridstruct%grid_64
area => Atm%gridstruct%area_64
area_c => Atm%gridstruct%area_c_64
rarea => Atm%gridstruct%rarea
rarea_c => Atm%gridstruct%rarea_c
sina => Atm%gridstruct%sina_64
cosa => Atm%gridstruct%cosa_64
dx => Atm%gridstruct%dx_64
dy => Atm%gridstruct%dy_64
dxc => Atm%gridstruct%dxc_64
dyc => Atm%gridstruct%dyc_64
dxa => Atm%gridstruct%dxa_64
dya => Atm%gridstruct%dya_64
rdx => Atm%gridstruct%rdx
rdy => Atm%gridstruct%rdy
rdxc => Atm%gridstruct%rdxc
rdyc => Atm%gridstruct%rdyc
rdxa => Atm%gridstruct%rdxa
rdya => Atm%gridstruct%rdya
e1 => Atm%gridstruct%e1
e2 => Atm%gridstruct%e2
if (Atm%neststruct%nested .or. ANY(Atm%neststruct%child_grids)) then
grid_global => Atm%grid_global
else if( trim(grid_file) .NE. 'INPUT/grid_spec.nc') then
allocate(grid_global(1-ng:npx +ng,1-ng:npy +ng,ndims,1:nregions))
endif
iinta => Atm%gridstruct%iinta
jinta => Atm%gridstruct%jinta
iintb => Atm%gridstruct%iintb
jintb => Atm%gridstruct%jintb
npx_g => Atm%gridstruct%npx_g
npy_g => Atm%gridstruct%npy_g
ntiles_g => Atm%gridstruct%ntiles_g
sw_corner => Atm%gridstruct%sw_corner
se_corner => Atm%gridstruct%se_corner
ne_corner => Atm%gridstruct%ne_corner
nw_corner => Atm%gridstruct%nw_corner
latlon => Atm%gridstruct%latlon
cubed_sphere => Atm%gridstruct%cubed_sphere
have_south_pole => Atm%gridstruct%have_south_pole
have_north_pole => Atm%gridstruct%have_north_pole
stretched_grid => Atm%gridstruct%stretched_grid
tile => Atm%tile_of_mosaic
domain => Atm%domain
npx_g = npx
npy_g = npy
ntiles_g = nregions
latlon = .false.
cubed_sphere = .false.
if ( (Atm%flagstruct%do_schmidt .or. Atm%flagstruct%do_cube_transform) .and. abs(atm%flagstruct%stretch_fac-1.) > 1.E-5 ) then
stretched_grid = .true.
if (Atm%flagstruct%do_schmidt .and. Atm%flagstruct%do_cube_transform) then
call mpp_error(FATAL, ' Cannot set both do_schmidt and do_cube_transform to .true.')
endif
endif
if (Atm%flagstruct%grid_type>3) then
if (Atm%flagstruct%grid_type == 4) then
call setup_cartesian(npx, npy, Atm%flagstruct%dx_const, Atm%flagstruct%dy_const, &
Atm%flagstruct%deglat, Atm%bd)
else
call mpp_error(FATAL, 'init_grid: unsupported grid type')
endif
else
cubed_sphere = .true.
if (Atm%neststruct%nested) then
!Read grid if it exists
! still need to set up
call setup_aligned_nest(Atm)
else
if(trim(grid_file) == 'INPUT/grid_spec.nc') then
call read_grid(Atm, grid_file, ndims, nregions, ng)
else
if (Atm%flagstruct%grid_type>=0) call gnomonic_grids(Atm%flagstruct%grid_type, npx-1, xs, ys)
if (is_master()) then
if (Atm%flagstruct%grid_type>=0) then
do j=1,npy
do i=1,npx
grid_global(i,j,1,1) = xs(i,j)
grid_global(i,j,2,1) = ys(i,j)
enddo
enddo
! mirror_grid assumes that the tile=1 is centered on equator and greenwich meridian Lon[-pi,pi]
call mirror_grid(grid_global, ng, npx, npy, 2, 6)
do n=1,nregions
do j=1,npy
do i=1,npx
!---------------------------------
! Shift the corner away from Japan
!---------------------------------
!--------------------- This will result in the corner close to east coast of China ------------------
if ( .not. ( Atm%flagstruct%do_schmidt .or. Atm%flagstruct%do_cube_transform) .and. (Atm%flagstruct%shift_fac)>1.E-4 ) &
grid_global(i,j,1,n) = grid_global(i,j,1,n) - pi/Atm%flagstruct%shift_fac
!----------------------------------------------------------------------------------------------------
if ( grid_global(i,j,1,n) < 0. ) &
grid_global(i,j,1,n) = grid_global(i,j,1,n) + 2.*pi
if (ABS(grid_global(i,j,1,1)) < 1.d-10) grid_global(i,j,1,1) = 0.0
if (ABS(grid_global(i,j,2,1)) < 1.d-10) grid_global(i,j,2,1) = 0.0
enddo
enddo
enddo
else
call mpp_error(FATAL, "fv_grid_tools: reading of ASCII grid files no longer supported")
endif
grid_global( 1,1:npy,:,2)=grid_global(npx,1:npy,:,1)
grid_global( 1,1:npy,:,3)=grid_global(npx:1:-1,npy,:,1)
grid_global(1:npx,npy,:,5)=grid_global(1,npy:1:-1,:,1)
grid_global(1:npx,npy,:,6)=grid_global(1:npx,1,:,1)
grid_global(1:npx, 1,:,3)=grid_global(1:npx,npy,:,2)
grid_global(1:npx, 1,:,4)=grid_global(npx,npy:1:-1,:,2)
grid_global(npx,1:npy,:,6)=grid_global(npx:1:-1,1,:,2)
grid_global( 1,1:npy,:,4)=grid_global(npx,1:npy,:,3)
grid_global( 1,1:npy,:,5)=grid_global(npx:1:-1,npy,:,3)
grid_global(npx,1:npy,:,3)=grid_global(1,1:npy,:,4)
grid_global(1:npx, 1,:,5)=grid_global(1:npx,npy,:,4)
grid_global(1:npx, 1,:,6)=grid_global(npx,npy:1:-1,:,4)
grid_global( 1,1:npy,:,6)=grid_global(npx,1:npy,:,5)
!------------------------
! Schmidt transformation:
!------------------------
if ( Atm%flagstruct%do_schmidt ) then
do n=1,nregions
call direct_transform(Atm%flagstruct%stretch_fac, 1, npx, 1, npy, &
Atm%flagstruct%target_lon, Atm%flagstruct%target_lat, &
n, grid_global(1:npx,1:npy,1,n), grid_global(1:npx,1:npy,2,n))
enddo
elseif (Atm%flagstruct%do_cube_transform) then
do n=1,nregions
call cube_transform(Atm%flagstruct%stretch_fac, 1, npx, 1, npy, &
Atm%flagstruct%target_lon, Atm%flagstruct%target_lat, &
n, grid_global(1:npx,1:npy,1,n), grid_global(1:npx,1:npy,2,n))
enddo
endif
endif !is master
call mpp_broadcast(grid_global, size(grid_global), mpp_root_pe())
!--- copy grid to compute domain
do n=1,ndims
do j=js,je+1
do i=is,ie+1
grid(i,j,n) = grid_global(i,j,n,tile)
enddo
enddo
enddo
endif !(trim(grid_file) == 'INPUT/grid_spec.nc')
!
! SJL: For phys/exchange grid, etc
!
call mpp_update_domains( grid, Atm%domain, position=CORNER)
if (.not. (Atm%gridstruct%bounded_domain)) then
call fill_corners(grid(:,:,1), npx, npy, FILL=XDir, BGRID=.true.)
call fill_corners(grid(:,:,2), npx, npy, FILL=XDir, BGRID=.true.)
endif
!--- dx and dy
if( .not. Atm%gridstruct%bounded_domain) then
istart=is
iend=ie
jstart=js
jend=je
else
istart=isd
iend=ied
jstart=jsd
jend=jed
endif
do j = jstart, jend+1
do i = istart, iend
p1(1) = grid(i ,j,1)
p1(2) = grid(i ,j,2)
p2(1) = grid(i+1,j,1)
p2(2) = grid(i+1,j,2)
dx(i,j) = great_circle_dist( p2, p1, radius )
enddo
enddo
if( stretched_grid .or. Atm%gridstruct%bounded_domain ) then
do j = jstart, jend
do i = istart, iend+1
p1(1) = grid(i,j, 1)
p1(2) = grid(i,j, 2)
p2(1) = grid(i,j+1,1)
p2(2) = grid(i,j+1,2)
dy(i,j) = great_circle_dist( p2, p1, radius )
enddo
enddo
else
call get_symmetry(dx(is:ie,js:je+1), dy(is:ie+1,js:je), 0, 1, Atm%layout(1), Atm%layout(2), &
Atm%domain, Atm%tile_of_mosaic, Atm%gridstruct%npx_g, Atm%bd)
endif
call mpp_get_boundary( dy, dx, Atm%domain, ebufferx=ebuffer, wbufferx=wbuffer, sbuffery=sbuffer, nbuffery=nbuffer,&
flags=SCALAR_PAIR+XUPDATE, gridtype=CGRID_NE_PARAM)
if( .not. Atm%gridstruct%bounded_domain ) then
if(is == 1 .AND. mod(tile,2) .NE. 0) then ! on the west boundary
dy(is, js:je) = wbuffer(js:je)
endif
if(ie == npx-1) then ! on the east boundary
dy(ie+1, js:je) = ebuffer(js:je)
endif
endif
call mpp_update_domains( dy, dx, Atm%domain, flags=SCALAR_PAIR, &
gridtype=CGRID_NE_PARAM, complete=.true.)
if (cubed_sphere .and. (.not. (Atm%gridstruct%bounded_domain))) then
call fill_corners(dx, dy, npx, npy, DGRID=.true.)
endif
if( .not. stretched_grid ) &
call sorted_inta(isd, ied, jsd, jed, cubed_sphere, grid, iinta, jinta)
agrid(:,:,:) = -1.e25
!--- compute agrid (use same indices as for dx/dy above)
do j=jstart,jend
do i=istart,iend
if ( stretched_grid ) then
call cell_center2(grid(i,j, 1:2), grid(i+1,j, 1:2), &
grid(i,j+1,1:2), grid(i+1,j+1,1:2), &
agrid(i,j,1:2) )
else
call cell_center2(grid(iinta(1,i,j),jinta(1,i,j),1:2), &
grid(iinta(2,i,j),jinta(2,i,j),1:2), &
grid(iinta(3,i,j),jinta(3,i,j),1:2), &
grid(iinta(4,i,j),jinta(4,i,j),1:2), &
agrid(i,j,1:2) )
endif
enddo
enddo
call mpp_update_domains( agrid, Atm%domain, position=CENTER, complete=.true. )
if (.not. (Atm%gridstruct%bounded_domain)) then
call fill_corners(agrid(:,:,1), npx, npy, XDir, AGRID=.true.)
call fill_corners(agrid(:,:,2), npx, npy, YDir, AGRID=.true.)
endif
do j=jsd,jed
do i=isd,ied
call mid_pt_sphere(grid(i, j,1:2), grid(i, j+1,1:2), p1)
call mid_pt_sphere(grid(i+1,j,1:2), grid(i+1,j+1,1:2), p2)
dxa(i,j) = great_circle_dist( p2, p1, radius )
!
call mid_pt_sphere(grid(i,j ,1:2), grid(i+1,j ,1:2), p1)
call mid_pt_sphere(grid(i,j+1,1:2), grid(i+1,j+1,1:2), p2)
dya(i,j) = great_circle_dist( p2, p1, radius )
enddo
enddo
! call mpp_update_domains( dxa, dya, Atm%domain, flags=SCALAR_PAIR, gridtype=AGRID_PARAM)
if (cubed_sphere .and. (.not. (Atm%gridstruct%bounded_domain))) then
call fill_corners(dxa, dya, npx, npy, AGRID=.true.)
endif
end if !if nested
! do j=js,je
! do i=is,ie+1
do j=jsd,jed
do i=isd+1,ied
dxc(i,j) = great_circle_dist(agrid(i,j,:), agrid(i-1,j,:), radius)
enddo
!xxxxxx
!Are the following 2 lines appropriate for the regional domain?
!xxxxxx
dxc(isd,j) = dxc(isd+1,j)
dxc(ied+1,j) = dxc(ied,j)
enddo
! do j=js,je+1
! do i=is,ie
do j=jsd+1,jed
do i=isd,ied
dyc(i,j) = great_circle_dist(agrid(i,j,:), agrid(i,j-1,:), radius)
enddo
enddo
!xxxxxx
!Are the following 2 lines appropriate for the regional domain?
!xxxxxx
do i=isd,ied
dyc(i,jsd) = dyc(i,jsd+1)
dyc(i,jed+1) = dyc(i,jed)
end do
if( .not. stretched_grid ) &
call sorted_intb(isd, ied, jsd, jed, is, ie, js, je, npx, npy, &
cubed_sphere, agrid, iintb, jintb)
call grid_area( npx, npy, ndims, nregions, Atm%gridstruct%bounded_domain, Atm%gridstruct, Atm%domain, Atm%bd )
! stretched_grid = .false.
!----------------------------------
! Compute area_c, rarea_c, dxc, dyc
!----------------------------------
if ( .not. stretched_grid .and. (.not. (Atm%gridstruct%bounded_domain))) then
! For symmetrical grids:
if ( is==1 ) then
i = 1
do j=js,je+1
call mid_pt_sphere(grid(i,j-1,1:2), grid(i,j, 1:2), p1)
call mid_pt_sphere(grid(i,j ,1:2), grid(i,j+1,1:2), p4)
p2(1:2) = agrid(i,j-1,1:2)
p3(1:2) = agrid(i,j, 1:2)
area_c(i,j) = 2.*get_area(p1, p4, p2, p3, radius)
enddo
do j=js,je
call mid_pt_sphere(grid(i,j,1:2), grid(i,j+1,1:2), p1)
p2(1:2) = agrid(i,j,1:2)
dxc(i,j) = 2.*great_circle_dist( p1, p2, radius )
enddo
endif
if ( (ie+1)==npx ) then
i = npx
do j=js,je+1
p1(1:2) = agrid(i-1,j-1,1:2)
call mid_pt_sphere(grid(i,j-1,1:2), grid(i,j, 1:2), p2)
call mid_pt_sphere(grid(i,j ,1:2), grid(i,j+1,1:2), p3)
p4(1:2) = agrid(i-1,j,1:2)
area_c(i,j) = 2.*get_area(p1, p4, p2, p3, radius)
enddo
do j=js,je
p1(1:2) = agrid(i-1,j,1:2)
call mid_pt_sphere(grid(i,j,1:2), grid(i,j+1,1:2), p2)
dxc(i,j) = 2.*great_circle_dist( p1, p2, radius )
enddo
endif
if ( js==1 ) then
j = 1
do i=is,ie+1
call mid_pt_sphere(grid(i-1,j,1:2), grid(i, j,1:2), p1)
call mid_pt_sphere(grid(i, j,1:2), grid(i+1,j,1:2), p2)
p3(1:2) = agrid(i, j,1:2)
p4(1:2) = agrid(i-1,j,1:2)
area_c(i,j) = 2.*get_area(p1, p4, p2, p3, radius)
enddo
do i=is,ie
call mid_pt_sphere(grid(i,j,1:2), grid(i+1,j,1:2), p1)
p2(1:2) = agrid(i,j,1:2)
dyc(i,j) = 2.*great_circle_dist( p1, p2, radius )
enddo
endif
if ( (je+1)==npy ) then
j = npy
do i=is,ie+1
p1(1:2) = agrid(i-1,j-1,1:2)
p2(1:2) = agrid(i ,j-1,1:2)
call mid_pt_sphere(grid(i, j,1:2), grid(i+1,j,1:2), p3)
call mid_pt_sphere(grid(i-1,j,1:2), grid(i, j,1:2), p4)
area_c(i,j) = 2.*get_area(p1, p4, p2, p3, radius)
enddo
do i=is,ie
p1(1:2) = agrid(i,j-1,1:2)
call mid_pt_sphere(grid(i,j,1:2), grid(i+1,j,1:2), p2)
dyc(i,j) = 2.*great_circle_dist( p1, p2, radius )
enddo
endif
if ( sw_corner ) then
i=1; j=1
p1(1:2) = grid(i,j,1:2)
call mid_pt_sphere(grid(i,j,1:2), grid(i+1,j,1:2), p2)
p3(1:2) = agrid(i,j,1:2)
call mid_pt_sphere(grid(i,j,1:2), grid(i,j+1,1:2), p4)
area_c(i,j) = 3.*get_area(p1, p4, p2, p3, radius)
endif
if ( se_corner ) then
i=npx; j=1
call mid_pt_sphere(grid(i-1,j,1:2), grid(i,j,1:2), p1)
p2(1:2) = grid(i,j,1:2)
call mid_pt_sphere(grid(i,j,1:2), grid(i,j+1,1:2), p3)
p4(1:2) = agrid(i,j,1:2)
area_c(i,j) = 3.*get_area(p1, p4, p2, p3, radius)
endif
if ( ne_corner ) then
i=npx; j=npy
p1(1:2) = agrid(i-1,j-1,1:2)
call mid_pt_sphere(grid(i,j-1,1:2), grid(i,j,1:2), p2)
p3(1:2) = grid(i,j,1:2)
call mid_pt_sphere(grid(i-1,j,1:2), grid(i,j,1:2), p4)
area_c(i,j) = 3.*get_area(p1, p4, p2, p3, radius)
endif
if ( nw_corner ) then
i=1; j=npy
call mid_pt_sphere(grid(i,j-1,1:2), grid(i,j,1:2), p1)
p2(1:2) = agrid(i,j-1,1:2)
call mid_pt_sphere(grid(i,j,1:2), grid(i+1,j,1:2), p3)
p4(1:2) = grid(i,j,1:2)
area_c(i,j) = 3.*get_area(p1, p4, p2, p3, radius)
endif
endif
!-----------------
call mpp_update_domains( dxc, dyc, Atm%domain, flags=SCALAR_PAIR, &
gridtype=CGRID_NE_PARAM, complete=.true.)
if (cubed_sphere .and. (.not. (Atm%gridstruct%bounded_domain))) then
call fill_corners(dxc, dyc, npx, npy, CGRID=.true.)
endif
call mpp_update_domains( area, Atm%domain, complete=.true. )
!Handling outermost ends for area_c
if (Atm%gridstruct%bounded_domain) then
if (is == 1) then
do j=jsd,jed
area_c(isd,j) = area_c(isd+1,j)
end do
if (js == 1) area_c(isd,jsd) = area_c(isd+1,jsd+1)
if (js == npy-1) area_c(isd,jed+1) = area_c(isd+1,jed)
end if
if (ie == npx-1) then
do j=jsd,jed
area_c(ied+1,j) = area_c(ied,j)
end do
if (js == 1) area_c(ied+1,jsd) = area_c(ied,jsd+1)
if (js == npy-1) area_c(ied+1,jed+1) = area_c(ied,jed)
end if
if (js == 1) then
do i=isd,ied
area_c(i,jsd) = area_c(i,jsd+1)
end do
end if
if (je == npy-1) then
do i=isd,ied
area_c(i,jed+1) = area_c(i,jed)
end do
end if
end if
call mpp_update_domains( area_c, Atm%domain, position=CORNER, complete=.true.)
! Handle corner Area ghosting
if (cubed_sphere .and. (.not. (Atm%gridstruct%bounded_domain))) then
call fill_ghost(area, npx, npy, -big_number, Atm%bd) ! fill in garbage values
call fill_corners(area_c, npx, npy, FILL=XDir, BGRID=.true.)
endif
do j=jsd,jed+1
do i=isd,ied
rdx(i,j) = 1.0/dx(i,j)
enddo
enddo
do j=jsd,jed
do i=isd,ied+1
rdy(i,j) = 1.0/dy(i,j)
enddo
enddo
do j=jsd,jed
do i=isd,ied+1
rdxc(i,j) = 1.0/dxc(i,j)
enddo
enddo
do j=jsd,jed+1
do i=isd,ied
rdyc(i,j) = 1.0/dyc(i,j)
enddo
enddo
do j=jsd,jed
do i=isd,ied
rarea(i,j) = 1.0/area(i,j)
rdxa(i,j) = 1./dxa(i,j)
rdya(i,j) = 1./dya(i,j)
enddo
enddo
do j=jsd,jed+1
do i=isd,ied+1
rarea_c(i,j) = 1.0/area_c(i,j)
enddo
enddo