forked from NOAA-PSL/stochastic_physics
-
Notifications
You must be signed in to change notification settings - Fork 1
/
setlats_lag_stochy.f
134 lines (129 loc) · 3.83 KB
/
setlats_lag_stochy.f
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
!>@brief The module 'setlats_lag_stochy_mod' contains the subroutine setlats_lag_stochy
module setlats_lag_stochy_mod
implicit none
contains
!>@brief The subroutine 'setlats_a_stochy' selects the latitude points on this task
! and halos
!>@details This code is taken from the legacy spectral GFS
subroutine setlats_lag_stochy(lats_nodes_a, global_lats_a,
& lats_nodes_h, global_lats_h, yhalo)
!
use spectral_layout_mod, only : me,nodes,latg
implicit none
!
integer yhalo
!
integer lats_nodes_a(nodes), lats_nodes_h(nodes)
&, global_lats_a(latg)
&, global_lats_h(latg+2*yhalo*nodes)
!
integer jj,jpt_a,jpt_h,lat_val,nn,nodes_lats
&, j1, j2, iprint
!
lats_nodes_h = 0
!
nodes_lats = 0
do nn=1,nodes
if (lats_nodes_a(nn) > 0) then
lats_nodes_h(nn) = lats_nodes_a(nn) + yhalo + yhalo
nodes_lats = nodes_lats + 1
endif
enddo
!
global_lats_h = 0
!
! set non-yhalo latitudes
!
jpt_a = 0
jpt_h = yhalo
do nn=1,nodes
if (lats_nodes_a(nn) > 0) then
do jj=1,lats_nodes_a(nn)
jpt_a = jpt_a + 1
jpt_h = jpt_h + 1
global_lats_h(jpt_h) = global_lats_a(jpt_a)
enddo
jpt_h = jpt_h + yhalo + yhalo
endif
enddo
!
j1 = latg + (yhalo+yhalo) * nodes_lats
do jj=1,yhalo
j2 = yhalo - jj
global_lats_h(jj) = global_lats_a(1) + j2 ! set north pole yhalo
global_lats_h(j1-j2) = global_lats_a(latg) + 1 - jj ! set south pole yhalo
enddo
!
if (lats_nodes_a(1) /= latg) then
!
! set non-polar south yhalos
jpt_h = 0
do nn=1,nodes-1
if (lats_nodes_h(nn).GT.0) then
jpt_h = jpt_h + lats_nodes_h(nn)
lat_val = global_lats_h(jpt_h-yhalo)
do jj=1,yhalo
global_lats_h(jpt_h-yhalo+jj) = min(lat_val+jj,latg)
enddo
endif
enddo
!
! set non-polar north yhalos
jpt_h = 0
do nn=1,nodes-1
if (lats_nodes_h(nn).GT.0) then
jpt_h = jpt_h + lats_nodes_h(nn)
lat_val = global_lats_h(jpt_h+yhalo+1)
do jj=1,yhalo
global_lats_h(jpt_h+yhalo-(jj-1)) = max(lat_val-jj,1)
enddo
endif
enddo
!
endif
!
iprint = 0
! iprint = 1
if (iprint == 1 .and. me == 0) then
!
write(me+6000,'("setlats_h yhalo=",i3," nodes=",i3/)')
& yhalo,nodes
!
do nn=1,nodes
write(me+6000,'("lats_nodes_a(",i4,")=",i4," ",
& " lats_nodes_h(",i4,")=",i4)')
& nn, lats_nodes_a(nn),
& nn, lats_nodes_h(nn)
enddo
!
jpt_a = 0
do nn=1,nodes
if (lats_nodes_a(nn) > 0) then
write(me+6000,'(" ")')
do jj=1,lats_nodes_a(nn)
jpt_a=jpt_a+1
write(me+6000,'(2i4," global_lats_a(",i4,")=",i4)')
& nn, jj, jpt_a, global_lats_a(jpt_a)
enddo
endif
enddo
!
jpt_h=0
do nn=1,nodes
if (lats_nodes_h(nn).gt.0) then
write(me+6000,'(" ")')
do jj=1,lats_nodes_h(nn)
jpt_h=jpt_h+1
write(me+6000,'(2i4," global_lats_h(",i4,")=",i4)')
& nn, jj, jpt_h, global_lats_h(jpt_h)
enddo
endif
enddo
!
close(6000+me)
endif
! close(6000+me)
!
return
end
end module setlats_lag_stochy_mod