MOM6
mom_spatial_means Module Reference

Functions/Subroutines

real function, public global_area_mean (var, G)
 
real function, public global_area_integral (var, G)
 
real function, dimension(gv %ke), public global_layer_mean (var, h, G, GV)
 
real function, public global_volume_mean (var, h, G, GV)
 
subroutine, public global_i_mean (array, i_mean, G, mask)
 
subroutine, public global_j_mean (array, j_mean, G, mask)
 
subroutine, public adjust_area_mean_to_zero (array, G, scaling)
 Adjust 2d array such that area mean is zero without moving the zero contour. More...
 

Function/Subroutine Documentation

◆ adjust_area_mean_to_zero()

subroutine, public mom_spatial_means::adjust_area_mean_to_zero ( real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(inout)  array,
type(ocean_grid_type), intent(inout)  G,
real, intent(out), optional  scaling 
)

Adjust 2d array such that area mean is zero without moving the zero contour.

Parameters
[in,out]gGrid structure
[in,out]array2D array to be adjusted
[out]scalingThe scaling factor used

Definition at line 290 of file MOM_spatial_means.F90.

290  type(ocean_grid_type), intent(inout) :: g !< Grid structure
291  real, dimension(SZI_(G),SZJ_(G)), intent(inout) :: array !< 2D array to be adjusted
292  real, optional, intent(out) :: scaling !< The scaling factor used
293  ! Local variables
294  real, dimension(SZI_(G), SZJ_(G)) :: posvals, negvals, areaxposvals, areaxnegvals
295  integer :: i,j
296  real :: areaintposvals, areaintnegvals, posscale, negscale
297 
298  areaxposvals(:,:) = 0.
299  areaxnegvals(:,:) = 0.
300 
301  do j=g%jsc,g%jec ; do i=g%isc,g%iec
302  posvals(i,j) = max(0., array(i,j))
303  areaxposvals(i,j) = g%areaT(i,j) * posvals(i,j)
304  negvals(i,j) = min(0., array(i,j))
305  areaxnegvals(i,j) = g%areaT(i,j) * negvals(i,j)
306  enddo ; enddo
307 
308  areaintposvals = reproducing_sum( areaxposvals )
309  areaintnegvals = reproducing_sum( areaxnegvals )
310 
311  posscale = 0.0 ; negscale = 0.0
312  if ((areaintposvals>0.).and.(areaintnegvals<0.)) then ! Only adjust if possible
313  if (areaintposvals>-areaintnegvals) then ! Scale down positive values
314  posscale = - areaintnegvals / areaintposvals
315  do j=g%jsc,g%jec ; do i=g%isc,g%iec
316  array(i,j) = (posscale * posvals(i,j)) + negvals(i,j)
317  enddo ; enddo
318  elseif (areaintposvals<-areaintnegvals) then ! Scale down negative values
319  negscale = - areaintposvals / areaintnegvals
320  do j=g%jsc,g%jec ; do i=g%isc,g%iec
321  array(i,j) = posvals(i,j) + (negscale * negvals(i,j))
322  enddo ; enddo
323  endif
324  endif
325  if (present(scaling)) scaling = posscale - negscale
326 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19

◆ global_area_integral()

real function, public mom_spatial_means::global_area_integral ( real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(in)  var,
type(ocean_grid_type), intent(in)  G 
)
Parameters
[in]gThe ocean's grid structure

Definition at line 62 of file MOM_spatial_means.F90.

Referenced by mom_forcing_type::forcing_diagnostics().

62  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
63  real, dimension(SZI_(G), SZJ_(G)), intent(in) :: var
64  real, dimension(SZI_(G), SZJ_(G)) :: tmpforsumming
65  real :: global_area_integral
66 
67  integer :: i, j, is, ie, js, je
68  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
69 
70  tmpforsumming(:,:) = 0.
71  do j=js,je ; do i=is, ie
72  tmpforsumming(i,j) = ( var(i,j) * (g%areaT(i,j) * g%mask2dT(i,j)) )
73  enddo ; enddo
74  global_area_integral = reproducing_sum( tmpforsumming )
75 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the caller graph for this function:

◆ global_area_mean()

real function, public mom_spatial_means::global_area_mean ( real, dimension(szi_(g), szj_(g)), intent(in)  var,
type(ocean_grid_type), intent(in)  G 
)
Parameters
[in]gThe ocean's grid structure

Definition at line 45 of file MOM_spatial_means.F90.

Referenced by mom_diagnostics::calculate_diagnostic_fields(), and mom_forcing_type::forcing_diagnostics().

45  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
46  real, dimension(SZI_(G), SZJ_(G)), intent(in) :: var
47  real, dimension(SZI_(G), SZJ_(G)) :: tmpforsumming
48  real :: global_area_mean
49 
50  integer :: i, j, is, ie, js, je
51  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
52 
53  tmpforsumming(:,:) = 0.
54  do j=js,je ; do i=is, ie
55  tmpforsumming(i,j) = ( var(i,j) * (g%areaT(i,j) * g%mask2dT(i,j)) )
56  enddo ; enddo
57  global_area_mean = reproducing_sum( tmpforsumming ) * g%IareaT_global
58 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the caller graph for this function:

◆ global_i_mean()

subroutine, public mom_spatial_means::global_i_mean ( real, dimension(szi_(g),szj_(g)), intent(in)  array,
real, dimension(szj_(g)), intent(out)  i_mean,
type(ocean_grid_type), intent(inout)  G,
real, dimension(szi_(g),szj_(g)), intent(in), optional  mask 
)
Parameters
[in,out]gThe ocean's grid structure

Definition at line 133 of file MOM_spatial_means.F90.

References mom_error_handler::mom_error(), mom_coms::query_efp_overflow_error(), and mom_coms::reset_efp_overflow_error().

Referenced by mom_sponge::apply_sponge().

133  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure
134  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: array
135  real, dimension(SZJ_(G)), intent(out) :: i_mean
136  real, dimension(SZI_(G),SZJ_(G)), optional, intent(in) :: mask
137 
138 ! This subroutine determines the global mean of a field along rows of
139 ! constant i, returning it in a 1-d array using the local indexing.
140 
141 ! Arguments: array - The 2-d array whose i-mean is to be taken.
142 ! (out) i_mean - Global mean of array along its i-axis.
143 ! (in) G - The ocean's grid structure.
144 ! (in) mask - An array used for weighting the i-mean.
145 
146  type(efp_type), allocatable, dimension(:) :: asum, mask_sum
147  real :: mask_sum_r
148  integer :: is, ie, js, je, idg_off, jdg_off
149  integer :: i, j
150 
151  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
152  idg_off = g%idg_offset ; jdg_off = g%jdg_offset
153 
154  call reset_efp_overflow_error()
155 
156  allocate(asum(g%jsg:g%jeg))
157  if (present(mask)) then
158  allocate(mask_sum(g%jsg:g%jeg))
159 
160  do j=g%jsg,g%jeg
161  asum(j) = real_to_efp(0.0) ; mask_sum(j) = real_to_efp(0.0)
162  enddo
163 
164  do i=is,ie ; do j=js,je
165  asum(j+jdg_off) = asum(j+jdg_off) + real_to_efp(array(i,j)*mask(i,j))
166  mask_sum(j+jdg_off) = mask_sum(j+jdg_off) + real_to_efp(mask(i,j))
167  enddo ; enddo
168 
169  if (query_efp_overflow_error()) call mom_error(fatal, &
170  "global_i_mean overflow error occurred before sums across PEs.")
171 
172  call efp_list_sum_across_pes(asum(g%jsg:g%jeg), g%jeg-g%jsg+1)
173  call efp_list_sum_across_pes(mask_sum(g%jsg:g%jeg), g%jeg-g%jsg+1)
174 
175  if (query_efp_overflow_error()) call mom_error(fatal, &
176  "global_i_mean overflow error occurred during sums across PEs.")
177 
178  do j=js,je
179  mask_sum_r = efp_to_real(mask_sum(j+jdg_off))
180  if (mask_sum_r == 0.0 ) then ; i_mean(j) = 0.0 ; else
181  i_mean(j) = efp_to_real(asum(j+jdg_off)) / mask_sum_r
182  endif
183  enddo
184 
185  deallocate(mask_sum)
186  else
187  do j=g%jsg,g%jeg ; asum(j) = real_to_efp(0.0) ; enddo
188 
189  do i=is,ie ; do j=js,je
190  asum(j+jdg_off) = asum(j+jdg_off) + real_to_efp(array(i,j))
191  enddo ; enddo
192 
193  if (query_efp_overflow_error()) call mom_error(fatal, &
194  "global_i_mean overflow error occurred before sum across PEs.")
195 
196  call efp_list_sum_across_pes(asum(g%jsg:g%jeg), g%jeg-g%jsg+1)
197 
198  if (query_efp_overflow_error()) call mom_error(fatal, &
199  "global_i_mean overflow error occurred during sum across PEs.")
200 
201  do j=js,je
202  i_mean(j) = efp_to_real(asum(j+jdg_off)) / real(g%ieg-g%isg+1)
203  enddo
204  endif
205 
206  deallocate(asum)
207 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:
Here is the caller graph for this function:

◆ global_j_mean()

subroutine, public mom_spatial_means::global_j_mean ( real, dimension(szi_(g),szj_(g)), intent(in)  array,
real, dimension(szi_(g)), intent(out)  j_mean,
type(ocean_grid_type), intent(inout)  G,
real, dimension(szi_(g),szj_(g)), intent(in), optional  mask 
)
Parameters
[in,out]gThe ocean's grid structure

Definition at line 211 of file MOM_spatial_means.F90.

References mom_error_handler::mom_error(), mom_coms::query_efp_overflow_error(), and mom_coms::reset_efp_overflow_error().

211  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure
212  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: array
213  real, dimension(SZI_(G)), intent(out) :: j_mean
214  real, dimension(SZI_(G),SZJ_(G)), optional, intent(in) :: mask
215 
216 ! This subroutine determines the global mean of a field along rows of
217 ! constant j, returning it in a 1-d array using the local indexing.
218 
219 ! Arguments: array - The 2-d array whose j-mean is to be taken.
220 ! (out) j_mean - Global mean of array along its j-axis.
221 ! (in) G - The ocean's grid structure.
222 ! (in) mask - An array used for weighting the j-mean.
223 
224  type(efp_type), allocatable, dimension(:) :: asum, mask_sum
225  real :: mask_sum_r
226  integer :: is, ie, js, je, idg_off, jdg_off
227  integer :: i, j
228 
229  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
230  idg_off = g%idg_offset ; jdg_off = g%jdg_offset
231 
232  call reset_efp_overflow_error()
233 
234  allocate(asum(g%isg:g%ieg))
235  if (present(mask)) then
236  allocate (mask_sum(g%isg:g%ieg))
237 
238  do i=g%isg,g%ieg
239  asum(i) = real_to_efp(0.0) ; mask_sum(i) = real_to_efp(0.0)
240  enddo
241 
242  do i=is,ie ; do j=js,je
243  asum(i+idg_off) = asum(i+idg_off) + real_to_efp(array(i,j)*mask(i,j))
244  mask_sum(i+idg_off) = mask_sum(i+idg_off) + real_to_efp(mask(i,j))
245  enddo ; enddo
246 
247  if (query_efp_overflow_error()) call mom_error(fatal, &
248  "global_j_mean overflow error occurred before sums across PEs.")
249 
250  call efp_list_sum_across_pes(asum(g%isg:g%ieg), g%ieg-g%isg+1)
251  call efp_list_sum_across_pes(mask_sum(g%isg:g%ieg), g%ieg-g%isg+1)
252 
253  if (query_efp_overflow_error()) call mom_error(fatal, &
254  "global_j_mean overflow error occurred during sums across PEs.")
255 
256  do i=is,ie
257  mask_sum_r = efp_to_real(mask_sum(i+idg_off))
258  if (mask_sum_r == 0.0 ) then ; j_mean(i) = 0.0 ; else
259  j_mean(i) = efp_to_real(asum(i+idg_off)) / mask_sum_r
260  endif
261  enddo
262 
263  deallocate(mask_sum)
264  else
265  do i=g%isg,g%ieg ; asum(i) = real_to_efp(0.0) ; enddo
266 
267  do i=is,ie ; do j=js,je
268  asum(i+idg_off) = asum(i+idg_off) + real_to_efp(array(i,j))
269  enddo ; enddo
270 
271  if (query_efp_overflow_error()) call mom_error(fatal, &
272  "global_j_mean overflow error occurred before sum across PEs.")
273 
274  call efp_list_sum_across_pes(asum(g%isg:g%ieg), g%ieg-g%isg+1)
275 
276  if (query_efp_overflow_error()) call mom_error(fatal, &
277  "global_j_mean overflow error occurred during sum across PEs.")
278 
279  do i=is,ie
280  j_mean(i) = efp_to_real(asum(i+idg_off)) / real(g%jeg-g%jsg+1)
281  enddo
282  endif
283 
284  deallocate(asum)
285 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:

◆ global_layer_mean()

real function, dimension( gv %ke), public mom_spatial_means::global_layer_mean ( real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  var,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  h,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV 
)
Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]hLayer thicknesses, in H (usually m or kg m-2)

Definition at line 79 of file MOM_spatial_means.F90.

Referenced by mom_diagnostics::calculate_diagnostic_fields().

79  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
80  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
81  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: var
82  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
83  real, dimension(SZK_(GV)) :: global_layer_mean
84 
85  real, dimension(SZI_(G), SZJ_(G), SZK_(GV)) :: tmpforsumming, weight
86  real, dimension(SZK_(GV)) :: scalarij, weightij
87  real, dimension(SZK_(GV)) :: global_temp_scalar, global_weight_scalar
88  integer :: i, j, k, is, ie, js, je, nz
89  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
90 
91  tmpforsumming(:,:,:) = 0. ; weight(:,:,:) = 0.
92 
93  do k=1,nz ; do j=js,je ; do i=is,ie
94  weight(i,j,k) = h(i,j,k) * (g%areaT(i,j) * g%mask2dT(i,j))
95  tmpforsumming(i,j,k) = var(i,j,k) * weight(i,j,k)
96  enddo ; enddo ; enddo
97 
98  global_temp_scalar = reproducing_sum(tmpforsumming,sums=scalarij)
99  global_weight_scalar = reproducing_sum(weight,sums=weightij)
100 
101  do k=1, nz
102  global_layer_mean(k) = scalarij(k) / weightij(k)
103  enddo
104 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the caller graph for this function:

◆ global_volume_mean()

real function, public mom_spatial_means::global_volume_mean ( real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  var,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  h,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV 
)
Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]hLayer thicknesses, in H (usually m or kg m-2)

Definition at line 108 of file MOM_spatial_means.F90.

Referenced by mom_diagnostics::calculate_diagnostic_fields().

108  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
109  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
110  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: var
111  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
112  real :: global_volume_mean
113 
114  real :: weight_here
115  real, dimension(SZI_(G), SZJ_(G)) :: tmpforsumming, sum_weight
116  integer :: i, j, k, is, ie, js, je, nz
117  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
118 
119  tmpforsumming(:,:) = 0. ; sum_weight(:,:) = 0.
120 
121  do k=1,nz ; do j=js,je ; do i=is,ie
122  weight_here = h(i,j,k) * (g%areaT(i,j) * g%mask2dT(i,j))
123  tmpforsumming(i,j) = tmpforsumming(i,j) + var(i,j,k) * weight_here
124  sum_weight(i,j) = sum_weight(i,j) + weight_here
125  enddo ; enddo ; enddo
126  global_volume_mean = (reproducing_sum(tmpforsumming)) / &
127  (reproducing_sum(sum_weight))
128 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the caller graph for this function: