Line data Source code
1 : module polar_avg
2 : !-----------------------------------------------------------------------
3 : !
4 : ! Purpose:
5 : ! These routines are used by the fv dycore to set the collocated
6 : ! pole points at the limits of the latitude dimension to the same
7 : ! value.
8 : !
9 : ! Methods:
10 : ! The reprosum reproducible distributed sum is used for these
11 : ! calculations.
12 : !
13 : ! Author: A. Mirin
14 : !
15 : !-----------------------------------------------------------------------
16 :
17 : !-----------------------------------------------------------------------
18 : !- use statements ------------------------------------------------------
19 : !-----------------------------------------------------------------------
20 : use shr_kind_mod, only: r8 => shr_kind_r8
21 : use dycore, only: dycore_is
22 : use dyn_grid, only: get_dyn_grid_parm
23 : use phys_grid, only: get_ncols_p, get_lat_all_p
24 : use ppgrid, only: begchunk, endchunk, pcols
25 : use shr_reprosum_mod, only: shr_reprosum_calc
26 : #if ( defined SPMD )
27 : use mpishorthand, only: mpicom
28 : #endif
29 :
30 : !-----------------------------------------------------------------------
31 : !- module boilerplate --------------------------------------------------
32 : !-----------------------------------------------------------------------
33 : implicit none
34 : private
35 : save
36 :
37 : !-----------------------------------------------------------------------
38 : ! Public interfaces ----------------------------------------------------
39 : !-----------------------------------------------------------------------
40 : public :: &
41 : polar_average ! support for LR dycore polar averaging
42 :
43 : interface polar_average
44 : module procedure polar_average2d, polar_average3d
45 : end interface
46 :
47 : CONTAINS
48 : !
49 : !========================================================================
50 : !
51 297984 : subroutine polar_average2d(field)
52 : !-----------------------------------------------------------------------
53 : ! Purpose: Set the collocated pole points at the limits of the latitude
54 : ! dimension to the same value.
55 : ! Author: J. Edwards
56 : !-----------------------------------------------------------------------
57 : !
58 : ! Arguments
59 : !
60 : real(r8), intent(inout) :: field(pcols,begchunk:endchunk)
61 : !
62 : ! Local workspace
63 : !
64 : integer :: i, c, ln, ls, ncols
65 : integer :: plat, plon
66 297984 : integer, allocatable :: lats(:)
67 : #if (! defined SPMD)
68 : integer :: mpicom = 0
69 : #endif
70 :
71 : real(r8) :: sum(2)
72 297984 : real(r8), allocatable :: n_pole(:), s_pole(:)
73 : !
74 : !-----------------------------------------------------------------------
75 : !
76 297984 : if(.not. dycore_is('LR')) return
77 :
78 297984 : plon = get_dyn_grid_parm('plon')
79 297984 : plat = get_dyn_grid_parm('plat')
80 1191936 : allocate(lats(pcols), n_pole(plon), s_pole(plon))
81 297984 : ln=0
82 297984 : ls=0
83 86117376 : n_pole = 0._r8
84 86117376 : s_pole = 0._r8
85 :
86 1787904 : do c=begchunk,endchunk
87 1489920 : call get_lat_all_p(c,pcols,lats)
88 1489920 : ncols = get_ncols_p(c)
89 :
90 23242752 : do i=1,ncols
91 22944768 : if(lats(i).eq.1) then
92 111744 : ln=ln+1
93 111744 : n_pole(ln) = field(i,c)
94 21343104 : else if(lats(i).eq.plat) then
95 111744 : ls=ls+1
96 111744 : s_pole(ls) = field(i,c)
97 : end if
98 : enddo
99 :
100 : end do
101 :
102 : call shr_reprosum_calc(n_pole, sum(1:1), ln, plon, 1, &
103 297984 : gbl_count=plon, commid=mpicom)
104 :
105 : call shr_reprosum_calc(s_pole, sum(2:2), ls, plon, 1, &
106 297984 : gbl_count=plon, commid=mpicom)
107 :
108 297984 : ln=0
109 297984 : ls=0
110 1787904 : do c=begchunk,endchunk
111 1489920 : call get_lat_all_p(c,pcols,lats)
112 1489920 : ncols = get_ncols_p(c)
113 :
114 23242752 : do i=1,ncols
115 22944768 : if(lats(i).eq.1) then
116 111744 : ln=ln+1
117 111744 : field(i,c) = sum(1)/plon
118 21343104 : else if(lats(i).eq.plat) then
119 111744 : ls=ls+1
120 111744 : field(i,c) = sum(2)/plon
121 : end if
122 : enddo
123 :
124 : end do
125 :
126 297984 : deallocate(lats, n_pole, s_pole)
127 :
128 297984 : end subroutine polar_average2d
129 :
130 : !
131 : !========================================================================
132 : !
133 :
134 52224 : subroutine polar_average3d(nlev, field)
135 : !-----------------------------------------------------------------------
136 : ! Purpose: Set the collocated pole points at the limits of the latitude
137 : ! dimension to the same value.
138 : ! Author: J. Edwards
139 : !-----------------------------------------------------------------------
140 : !
141 : ! Arguments
142 : !
143 : integer, intent(in) :: nlev
144 : real(r8), intent(inout) :: field(pcols,nlev,begchunk:endchunk)
145 : !
146 : ! Local workspace
147 : !
148 : integer :: i, c, ln, ls, ncols, k
149 : integer :: plat, plon
150 52224 : integer, allocatable :: lats(:)
151 : #if (! defined SPMD)
152 : integer :: mpicom = 0
153 : #endif
154 :
155 104448 : real(r8) :: sum(nlev,2)
156 52224 : real(r8), allocatable :: n_pole(:,:), s_pole(:,:)
157 : !
158 : !-----------------------------------------------------------------------
159 : !
160 52224 : if(.not. dycore_is('LR')) return
161 :
162 52224 : plon = get_dyn_grid_parm('plon')
163 52224 : plat = get_dyn_grid_parm('plat')
164 313344 : allocate(lats(pcols), n_pole(plon,nlev), s_pole(plon,nlev))
165 52224 : ln=0
166 52224 : ls=0
167 307233792 : n_pole = 0._r8
168 307233792 : s_pole = 0._r8
169 :
170 313344 : do c=begchunk,endchunk
171 261120 : call get_lat_all_p(c,pcols,lats)
172 261120 : ncols = get_ncols_p(c)
173 :
174 4073472 : do i=1,ncols
175 4021248 : if(lats(i).eq.1) then
176 19584 : ln=ln+1
177 418176 : do k=1,nlev
178 418176 : n_pole(ln,k) = field(i,k,c)
179 : end do
180 3740544 : else if(lats(i).eq.plat) then
181 19584 : ls=ls+1
182 418176 : do k=1,nlev
183 418176 : s_pole(ls,k) = field(i,k,c)
184 : end do
185 : end if
186 : enddo
187 : end do
188 :
189 : call shr_reprosum_calc(n_pole, sum(:,1), ln, plon, nlev, &
190 52224 : gbl_count=plon, commid=mpicom)
191 :
192 : call shr_reprosum_calc(s_pole, sum(:,2), ls, plon, nlev, &
193 52224 : gbl_count=plon, commid=mpicom)
194 :
195 52224 : ln=0
196 52224 : ls=0
197 313344 : do c=begchunk,endchunk
198 261120 : call get_lat_all_p(c,pcols,lats)
199 261120 : ncols = get_ncols_p(c)
200 :
201 4073472 : do i=1,ncols
202 4021248 : if(lats(i).eq.1) then
203 19584 : ln=ln+1
204 418176 : do k=1,nlev
205 418176 : field(i,k,c) = sum(k,1)/plon
206 : end do
207 3740544 : else if(lats(i).eq.plat) then
208 19584 : ls=ls+1
209 418176 : do k=1,nlev
210 418176 : field(i,k,c) = sum(k,2)/plon
211 : end do
212 : end if
213 : enddo
214 :
215 : end do
216 :
217 52224 : deallocate(lats, n_pole, s_pole)
218 :
219 52224 : end subroutine polar_average3d
220 :
221 : end module polar_avg
|