Line data Source code
1 : module mean_module
2 :
3 : use shr_kind_mod, only : r8 => shr_kind_r8
4 : use shr_reprosum_mod, only : shr_reprosum_calc, shr_reprosum_tolExceeded, &
5 : shr_reprosum_recompute
6 : use perf_mod
7 : use cam_logfile, only : iulog
8 :
9 : public gmean, gmeanxy
10 :
11 : private
12 : real(r8), parameter :: D0_0 = 0.0_r8
13 :
14 : contains
15 : !-----------------------------------------------------------------------
16 : !BOP
17 : ! !ROUTINE: gmean --- Calculate the mean of a 2D field
18 : !
19 : ! !INTERFACE:
20 :
21 0 : subroutine gmean(grid, q, qmean)
22 :
23 : ! !USES:
24 : use commap, only : w
25 : use dynamics_vars, only : T_FVDYCORE_GRID
26 :
27 : #if defined( SPMD )
28 : use parutilitiesmodule, only : parcollective, sumop
29 : #endif
30 :
31 : implicit none
32 :
33 : ! !INPUT PARAMETERS:
34 :
35 : type (T_FVDYCORE_GRID), intent(in) :: grid ! Grid information
36 : real(r8), intent(in) :: q(grid%im,grid%jfirst:grid%jlast) ! 2D field
37 :
38 : real(r8) qmean
39 :
40 : ! !DESCRIPTION:
41 : ! Calculate the mean of a 2D field
42 : !
43 : ! !REVISION HISTORY:
44 : ! 00.08.01 Lin Creation
45 : ! 01.01.10 Lin Revised
46 : ! 01.06.27 Mirin Use y communicator
47 : ! 05.07.12 Sawyer Simplified interface with grid argument
48 : !
49 : !EOP
50 : !-----------------------------------------------------------------------
51 : !BOC
52 :
53 0 : real(r8) :: xsum(grid%jm)
54 : integer :: i, j, im, jm, jfirst, jlast
55 :
56 0 : im = grid%im
57 0 : jm = grid%jm
58 0 : jfirst = grid%jfirst
59 0 : jlast = grid%jlast
60 :
61 0 : do j=1,jm
62 0 : xsum(j) = D0_0
63 : enddo
64 0 : do j=jfirst,jlast
65 0 : do i=1,im
66 0 : xsum(j) = xsum(j) + q(i,j)
67 : enddo
68 0 : xsum(j) = xsum(j)*w(j)
69 : enddo
70 :
71 : #if defined( SPMD )
72 0 : if (grid%npr_y .ne. 1) then
73 0 : call parcollective( grid%comm_y, sumop, jm, xsum )
74 : endif
75 : #endif
76 :
77 0 : qmean = D0_0
78 0 : do j=1,jm
79 0 : qmean = qmean + xsum(j)
80 : enddo
81 0 : qmean = qmean / (2*im)
82 :
83 0 : return
84 : !EOC
85 : end subroutine gmean
86 : !-----------------------------------------------------------------------
87 :
88 : !-----------------------------------------------------------------------
89 : !BOP
90 : ! !ROUTINE: gmeanxy --- Calculate the mean of a 2D field (XY decomp)
91 : !
92 : ! !INTERFACE:
93 :
94 5376 : subroutine gmeanxy(grid, q, qmean)
95 :
96 : ! !USES:
97 : use commap, only : w
98 : use dynamics_vars, only : T_FVDYCORE_GRID
99 :
100 : #if defined( SPMD )
101 : use parutilitiesmodule, only : parcollective, sumop
102 : #endif
103 :
104 : implicit none
105 :
106 : ! !INPUT PARAMETERS:
107 :
108 : type (T_FVDYCORE_GRID), intent(in) :: grid ! Grid information
109 : real(r8), intent(in), target :: q(grid%ifirstxy:grid%ilastxy, &
110 : grid%jfirstxy:grid%jlastxy) ! 2D field
111 :
112 : real(r8) qmean
113 :
114 : ! !DESCRIPTION:
115 : ! Calculate the mean of a 2D field on an XY decomposition
116 : ! This is inefficiently programmed (global collective operation),
117 : ! and is therefore only intended for initialization phase.
118 : !
119 : ! PW: gmeanxy is called in fv_prints, so replaced inefficient algorithm
120 : ! with shr_reprosum_calc.
121 : !
122 : ! !REVISION HISTORY:
123 : ! 00.08.01 Lin Creation
124 : ! 01.01.10 Lin Revised
125 : ! 01.06.27 Mirin Use y communicator
126 : ! 05.07.12 Sawyer Simplified interface with grid argument
127 : ! 05.08.26 Sawyer Modified for XY decomposition
128 : ! 08.07.03 Worley Introduced repro_sum logic
129 : ! 12.10.29 Santos repro_sum is now shr_reprosum
130 : !
131 : !EOP
132 : !-----------------------------------------------------------------------
133 : !BOC
134 :
135 : real(r8) :: q_tmp(grid%ifirstxy:grid%ilastxy, &
136 10752 : grid%jfirstxy:grid%jlastxy)
137 : real(r8) :: rel_diff(2), qmean_tmp(1), xsum
138 5376 : real(r8), allocatable :: q_global(:,:)
139 :
140 : integer :: i, j, im, jm, ifirstxy, ilastxy, jfirstxy, jlastxy
141 : integer :: lim, ljm, lijm
142 :
143 : logical :: write_warning
144 :
145 5376 : im = grid%im
146 5376 : jm = grid%jm
147 5376 : ifirstxy = grid%ifirstxy
148 5376 : ilastxy = grid%ilastxy
149 5376 : jfirstxy = grid%jfirstxy
150 5376 : jlastxy = grid%jlastxy
151 :
152 5376 : lim = ilastxy - ifirstxy + 1
153 5376 : ljm = jlastxy - jfirstxy + 1
154 5376 : lijm = lim*ljm
155 :
156 21504 : do j=jfirstxy,jlastxy
157 408576 : do i=ifirstxy,ilastxy
158 403200 : q_tmp(i,j) = q(i,j)*w(j)
159 : enddo
160 : enddo
161 :
162 5376 : call t_startf("gmeanxy_reprosum")
163 : call shr_reprosum_calc(q_tmp, qmean_tmp, lijm, lijm, 1, gbl_count=im*jm, &
164 5376 : commid=grid%commxy, rel_diff=rel_diff)
165 5376 : qmean = qmean_tmp(1)
166 5376 : call t_stopf("gmeanxy_reprosum")
167 :
168 : ! check that "fast" reproducible sum is accurate enough. If not, calculate
169 : ! using old method
170 5376 : write_warning = .false.
171 5376 : if (grid%iam == 0) write_warning = .true.
172 5376 : if ( shr_reprosum_tolExceeded('gmeanxy', 1, write_warning, &
173 : iulog, rel_diff) ) then
174 0 : if ( shr_reprosum_recompute ) then
175 0 : call t_startf("gmeanxy_sumfix")
176 0 : allocate( q_global(im,jm) )
177 0 : q_global = D0_0
178 0 : do j=jfirstxy,jlastxy
179 0 : do i=ifirstxy,ilastxy
180 0 : q_global(i,j) = q_tmp(i,j)
181 : enddo
182 : enddo
183 :
184 : #if defined( SPMD )
185 0 : call parcollective( grid%commxy, sumop, im, jm, q_global )
186 : #endif
187 0 : qmean = D0_0
188 0 : do j=1,jm
189 0 : xsum = D0_0
190 0 : do i=1,im
191 0 : xsum = xsum + q_global(i,j)
192 : enddo
193 0 : qmean = qmean + xsum
194 : enddo
195 :
196 0 : deallocate( q_global )
197 0 : call t_stopf("gmeanxy_sumfix")
198 : endif
199 : endif
200 :
201 5376 : qmean = qmean / (2*im)
202 :
203 5376 : return
204 : !EOC
205 5376 : end subroutine gmeanxy
206 : !-----------------------------------------------------------------------
207 :
208 : end module mean_module
|