Ignore:
Timestamp:
Jun 13, 2014, 5:28:30 PM (11 years ago)
Author:
milmd
Message:

LMDZ.GENERIC. Diagnostics are done globally. They slightly differ from previous version because of some shorcuts. Correct a bug when using meanOLR option.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phystd/planetwide_mod.F90

    r1216 r1295  
    2020end interface
    2121
     22interface planetwide_sumval ! sum() , over the entire planet
     23  module procedure planetwide_sumval_i1, planetwide_sumval_i2, &
     24                   planetwide_sumval_r1, planetwide_sumval_r2
     25end interface
     26
    2227contains
    2328
     
    200205  end subroutine planetwide_minval_r2
    201206
    202 
     207!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     208
     209  subroutine planetwide_sumval_i1(values,values_sum)
     210  use dimphy, only: klon
     211  use mod_grid_phy_lmdz, only : klon_glo
     212  implicit none
     213  integer,intent(in) :: values(:) ! local grid (klon)
     214  integer,intent(out) :: values_sum
     215#ifdef CPP_PARA
     216  integer :: values_glo(klon_glo) ! global grid
     217 
     218  ! gather field on master:
     219  call gather(values,values_glo)
     220  ! calculate sum value
     221  if (is_master) then
     222    values_sum=SUM(values_glo(:))
     223  endif
     224  ! broadcast information to all cores
     225  call bcast(values_sum)
     226#else
     227  values_sum=SUM(values(:))
     228#endif
     229  end subroutine planetwide_sumval_i1
     230
     231  subroutine planetwide_sumval_i2(values,values_sum)
     232  use dimphy, only: klon, klev
     233  use mod_grid_phy_lmdz, only : klon_glo
     234  implicit none
     235  integer,intent(in) :: values(:,:) ! local grid (klon,klev)
     236  integer,intent(out) :: values_sum
     237#ifdef CPP_PARA
     238  integer :: values_glo(klon_glo,klev) ! global grid
     239 
     240  ! gather field on master:
     241  call gather(values,values_glo)
     242  ! calculate sum value
     243  if (is_master) then
     244    values_sum=SUM(values_glo)
     245  endif
     246  ! broadcast information to all cores
     247  call bcast(values_sum)
     248#else
     249  values_sum=SUM(values)
     250#endif
     251  end subroutine planetwide_sumval_i2
     252
     253  subroutine planetwide_sumval_r1(values,values_sum)
     254  use dimphy, only: klon
     255  use mod_grid_phy_lmdz, only : klon_glo
     256  implicit none
     257  real,intent(in) :: values(:) ! local grid (klon)
     258  real,intent(out) :: values_sum
     259#ifdef CPP_PARA
     260  real :: values_glo(klon_glo) ! global grid
     261 
     262  ! gather field on master:
     263  call gather(values,values_glo)
     264  ! calculate sum value
     265  if (is_master) then
     266    values_sum=SUM(values_glo)
     267  endif
     268  ! broadcast information to all cores
     269  call bcast(values_sum)
     270#else
     271  values_sum=SUM(values)
     272#endif
     273  end subroutine planetwide_sumval_r1
     274
     275  subroutine planetwide_sumval_r2(values,values_sum)
     276  use dimphy, only: klon, klev
     277  use mod_grid_phy_lmdz, only : klon_glo
     278  implicit none
     279  real,intent(in) :: values(:,:) ! local grid (klon,klev)
     280  real,intent(out) :: values_sum
     281#ifdef CPP_PARA
     282  real :: values_glo(klon_glo,klev) ! global grid
     283 
     284  ! gather field on master:
     285  call gather(values,values_glo)
     286  ! calculate sum value
     287  if (is_master) then
     288    values_sum=SUM(values_glo)
     289  endif
     290  ! broadcast information to all cores
     291  call bcast(values_sum)
     292#else
     293  values_sum=SUM(values)
     294#endif
     295  end subroutine planetwide_sumval_r2
     296 
     297 
    203298end module planetwide_mod
Note: See TracChangeset for help on using the changeset viewer.