source: trunk/LMDZ.MARS/libf/phymars/planetwide_mod.F90 @ 2151

Last change on this file since 2151 was 1507, checked in by emillour, 9 years ago

Generic model and Mars model:

  • Fix/improvement on planetwide_min/max/sum for the 3D fields which assumed the vertical dimension to be klev. Now works for any (klon,...) field.

EM

File size: 8.0 KB
Line 
1!
2! $Id: $
3!
4module planetwide_mod
5! module which contains functions to obtain max/min/... values over the
6! entire globe (trivial in serial mode, but not in parallel)
7
8use mod_phys_lmdz_para, only : is_master, gather, bcast
9                               
10implicit none
11
12interface planetwide_maxval ! maxval() , over the entire planet
13  module procedure planetwide_maxval_i1, planetwide_maxval_i2, &
14                   planetwide_maxval_r1, planetwide_maxval_r2
15end interface
16
17interface planetwide_minval ! minval() , over the entire planet
18  module procedure planetwide_minval_i1, planetwide_minval_i2, &
19                   planetwide_minval_r1, planetwide_minval_r2
20end interface
21
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
27contains
28
29  subroutine planetwide_maxval_i1(values,values_max)
30  use dimphy, only: klon
31  use mod_grid_phy_lmdz, only : klon_glo
32  implicit none
33  integer,intent(in) :: values(:) ! local grid (klon)
34  integer,intent(out) :: values_max
35#ifdef CPP_PARA
36  integer :: values_glo(klon_glo) ! global grid
37 
38  ! gather field on master:
39  call gather(values,values_glo)
40  ! extract maximum value
41  if (is_master) then
42    values_max=maxval(values_glo)
43  endif
44  ! broadcast information to all cores
45  call bcast(values_max)
46#else
47  values_max=maxval(values)
48#endif
49  end subroutine planetwide_maxval_i1
50
51  subroutine planetwide_maxval_i2(values,values_max)
52  use dimphy, only: klon
53  use mod_grid_phy_lmdz, only : klon_glo
54  implicit none
55  integer,intent(in) :: values(:,:) ! local grid (klon,...)
56  integer,intent(out) :: values_max
57#ifdef CPP_PARA
58  integer :: values_glo(klon_glo,size(values,2)) ! global grid
59 
60  ! gather field on master:
61  call gather(values,values_glo)
62  ! extract maximum value
63  if (is_master) then
64    values_max=maxval(values_glo)
65  endif
66  ! broadcast information to all cores
67  call bcast(values_max)
68#else
69  values_max=maxval(values)
70#endif
71  end subroutine planetwide_maxval_i2
72
73  subroutine planetwide_maxval_r1(values,values_max)
74  use dimphy, only: klon
75  use mod_grid_phy_lmdz, only : klon_glo
76  implicit none
77  real,intent(in) :: values(:) ! local grid (klon)
78  real,intent(out) :: values_max
79#ifdef CPP_PARA
80  real :: values_glo(klon_glo) ! global grid
81 
82  ! gather field on master:
83  call gather(values,values_glo)
84  ! extract maximum value
85  if (is_master) then
86    values_max=maxval(values_glo)
87  endif
88  ! broadcast information to all cores
89  call bcast(values_max)
90#else
91  values_max=maxval(values)
92#endif
93  end subroutine planetwide_maxval_r1
94
95  subroutine planetwide_maxval_r2(values,values_max)
96  use dimphy, only: klon
97  use mod_grid_phy_lmdz, only : klon_glo
98  implicit none
99  real,intent(in) :: values(:,:) ! local grid (klon,...)
100  real,intent(out) :: values_max
101#ifdef CPP_PARA
102  real :: values_glo(klon_glo,size(values,2)) ! global grid
103 
104  ! gather field on master:
105  call gather(values,values_glo)
106  ! extract maximum value
107  if (is_master) then
108    values_max=maxval(values_glo)
109  endif
110  ! broadcast information to all cores
111  call bcast(values_max)
112#else
113  values_max=maxval(values)
114#endif
115  end subroutine planetwide_maxval_r2
116
117!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
118
119  subroutine planetwide_minval_i1(values,values_max)
120  use dimphy, only: klon
121  use mod_grid_phy_lmdz, only : klon_glo
122  implicit none
123  integer,intent(in) :: values(:) ! local grid (klon)
124  integer,intent(out) :: values_max
125#ifdef CPP_PARA
126  integer :: values_glo(klon_glo) ! global grid
127 
128  ! gather field on master:
129  call gather(values,values_glo)
130  ! extract maximum value
131  if (is_master) then
132    values_max=minval(values_glo)
133  endif
134  ! broadcast information to all cores
135  call bcast(values_max)
136#else
137  values_max=minval(values)
138#endif
139  end subroutine planetwide_minval_i1
140
141  subroutine planetwide_minval_i2(values,values_max)
142  use dimphy, only: klon
143  use mod_grid_phy_lmdz, only : klon_glo
144  implicit none
145  integer,intent(in) :: values(:,:) ! local grid (klon,...)
146  integer,intent(out) :: values_max
147#ifdef CPP_PARA
148  integer :: values_glo(klon_glo,size(values,2)) ! global grid
149 
150  ! gather field on master:
151  call gather(values,values_glo)
152  ! extract maximum value
153  if (is_master) then
154    values_max=minval(values_glo)
155  endif
156  ! broadcast information to all cores
157  call bcast(values_max)
158#else
159  values_max=minval(values)
160#endif
161  end subroutine planetwide_minval_i2
162
163  subroutine planetwide_minval_r1(values,values_max)
164  use dimphy, only: klon
165  use mod_grid_phy_lmdz, only : klon_glo
166  implicit none
167  real,intent(in) :: values(:) ! local grid (klon)
168  real,intent(out) :: values_max
169#ifdef CPP_PARA
170  real :: values_glo(klon_glo) ! global grid
171 
172  ! gather field on master:
173  call gather(values,values_glo)
174  ! extract maximum value
175  if (is_master) then
176    values_max=minval(values_glo)
177  endif
178  ! broadcast information to all cores
179  call bcast(values_max)
180#else
181  values_max=minval(values)
182#endif
183  end subroutine planetwide_minval_r1
184
185  subroutine planetwide_minval_r2(values,values_max)
186  use dimphy, only: klon
187  use mod_grid_phy_lmdz, only : klon_glo
188  implicit none
189  real,intent(in) :: values(:,:) ! local grid (klon,...)
190  real,intent(out) :: values_max
191#ifdef CPP_PARA
192  real :: values_glo(klon_glo,size(values,2)) ! global grid
193 
194  ! gather field on master:
195  call gather(values,values_glo)
196  ! extract maximum value
197  if (is_master) then
198    values_max=minval(values_glo)
199  endif
200  ! broadcast information to all cores
201  call bcast(values_max)
202#else
203  values_max=minval(values)
204#endif
205  end subroutine planetwide_minval_r2
206
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
233  use mod_grid_phy_lmdz, only : klon_glo
234  implicit none
235  integer,intent(in) :: values(:,:) ! local grid (klon,...)
236  integer,intent(out) :: values_sum
237#ifdef CPP_PARA
238  integer :: values_glo(klon_glo,size(values,2)) ! 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
277  use mod_grid_phy_lmdz, only : klon_glo
278  implicit none
279  real,intent(in) :: values(:,:) ! local grid (klon,...)
280  real,intent(out) :: values_sum
281#ifdef CPP_PARA
282  real :: values_glo(klon_glo,size(values,2)) ! 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 
298end module planetwide_mod
Note: See TracBrowser for help on using the repository browser.