New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
lib_fortran.F90 in branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC – NEMO

source: branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/lib_fortran.F90 @ 7508

Last change on this file since 7508 was 7508, checked in by mocavero, 7 years ago

changes on code duplication and workshare construct

  • Property svn:keywords set to Id
File size: 32.8 KB
Line 
1MODULE lib_fortran
2   !!======================================================================
3   !!                       ***  MODULE  lib_fortran  ***
4   !! Fortran utilities:  includes some low levels fortran functionality
5   !!======================================================================
6   !! History :  3.2  !  2010-05  (M. Dunphy, R. Benshila)  Original code
7   !!            3.4  !  2013-06  (C. Rousset)  add glob_min, glob_max
8   !!                                           + 3d dim. of input is fexible (jpk, jpl...)
9   !!----------------------------------------------------------------------
10
11   !!----------------------------------------------------------------------
12   !!   glob_sum    : generic interface for global masked summation over
13   !!                 the interior domain for 1 or 2 2D or 3D arrays
14   !!                 it works only for T points
15   !!   SIGN        : generic interface for SIGN to overwrite f95 behaviour
16   !!                 of intrinsinc sign function
17   !!----------------------------------------------------------------------
18   USE par_oce         ! Ocean parameter
19   USE dom_oce         ! ocean domain
20   USE in_out_manager  ! I/O manager
21   USE lib_mpp         ! distributed memory computing
22
23   IMPLICIT NONE
24   PRIVATE
25
26   PUBLIC   glob_sum      ! used in many places (masked with tmask_i)
27   PUBLIC   glob_sum_full ! used in many places (masked with tmask_h, ie omly over the halos)
28   PUBLIC   DDPDD         ! also used in closea module
29   PUBLIC   glob_min, glob_max
30#if defined key_nosignedzero
31   PUBLIC SIGN
32#endif
33
34   INTERFACE glob_sum
35      MODULE PROCEDURE glob_sum_1d, glob_sum_2d, glob_sum_3d, &
36         &             glob_sum_2d_a, glob_sum_3d_a
37   END INTERFACE
38   INTERFACE glob_sum_full
39      MODULE PROCEDURE glob_sum_full_2d, glob_sum_full_3d
40   END INTERFACE
41   INTERFACE glob_min
42      MODULE PROCEDURE glob_min_2d, glob_min_3d,glob_min_2d_a, glob_min_3d_a 
43   END INTERFACE
44   INTERFACE glob_max
45      MODULE PROCEDURE glob_max_2d, glob_max_3d,glob_max_2d_a, glob_max_3d_a 
46   END INTERFACE
47
48#if defined key_nosignedzero
49   INTERFACE SIGN
50      MODULE PROCEDURE SIGN_SCALAR, SIGN_ARRAY_1D, SIGN_ARRAY_2D, SIGN_ARRAY_3D,   &
51         &             SIGN_ARRAY_1D_A, SIGN_ARRAY_2D_A, SIGN_ARRAY_3D_A,          &
52         &             SIGN_ARRAY_1D_B, SIGN_ARRAY_2D_B, SIGN_ARRAY_3D_B
53   END INTERFACE
54#endif
55
56   !!----------------------------------------------------------------------
57   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
58   !! $Id$
59   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
60   !!----------------------------------------------------------------------
61CONTAINS
62
63#if ! defined key_mpp_rep
64   ! --- SUM ---
65
66   FUNCTION glob_sum_1d( ptab, kdim )
67      !!-----------------------------------------------------------------------
68      !!                  ***  FUNCTION  glob_sum_1D  ***
69      !!
70      !! ** Purpose : perform a masked sum on the inner global domain of a 1D array
71      !!-----------------------------------------------------------------------
72      INTEGER :: kdim
73      REAL(wp), INTENT(in), DIMENSION(kdim) ::   ptab        ! input 1D array
74      REAL(wp)                              ::   glob_sum_1d ! global sum
75      !!-----------------------------------------------------------------------
76      !
77      glob_sum_1d = SUM( ptab(:) )
78      IF( lk_mpp )   CALL mpp_sum( glob_sum_1d )
79      !
80   END FUNCTION glob_sum_1d
81
82   FUNCTION glob_sum_2d( ptab )
83      !!-----------------------------------------------------------------------
84      !!                  ***  FUNCTION  glob_sum_2D  ***
85      !!
86      !! ** Purpose : perform a masked sum on the inner global domain of a 2D array
87      !!-----------------------------------------------------------------------
88      REAL(wp), INTENT(in), DIMENSION(:,:) ::   ptab          ! input 2D array
89      REAL(wp)                             ::   glob_sum_2d   ! global masked sum
90      !!-----------------------------------------------------------------------
91      !
92      glob_sum_2d = SUM( ptab(:,:)*tmask_i(:,:) )
93      IF( lk_mpp )   CALL mpp_sum( glob_sum_2d )
94      !
95   END FUNCTION glob_sum_2d
96
97
98   FUNCTION glob_sum_3d( ptab )
99      !!-----------------------------------------------------------------------
100      !!                  ***  FUNCTION  glob_sum_3D  ***
101      !!
102      !! ** Purpose : perform a masked sum on the inner global domain of a 3D array
103      !!-----------------------------------------------------------------------
104      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   ptab          ! input 3D array
105      REAL(wp)                               ::   glob_sum_3d, a, b   ! global masked sum
106      !!
107      INTEGER :: jk, jj, ji
108      INTEGER :: ijpk ! local variable: size of the 3d dimension of ptab
109      !!-----------------------------------------------------------------------
110      !
111      ijpk = SIZE(ptab,3)
112      !
113      glob_sum_3d = 0.e0
114!$OMP PARALLEL
115!$OMP DO schedule(static) private(jk,jj,ji,b) REDUCTION(+:a,glob_sum_3d)
116      DO jk = 1, ijpk
117      a = 0.e0
118      b = 0.e0
119         DO jj = 1, jpj
120            DO ji = 1, jpi
121               b = ptab(ji,jj,jk)*tmask_i(ji,jj)
122               a = a + b
123           END DO
124         END DO
125         glob_sum_3d = glob_sum_3d + a
126      END DO
127!$OMP END PARALLEL
128
129      IF( lk_mpp )   CALL mpp_sum( glob_sum_3d )
130      !
131   END FUNCTION glob_sum_3d
132
133
134   FUNCTION glob_sum_2d_a( ptab1, ptab2 )
135      !!-----------------------------------------------------------------------
136      !!                  ***  FUNCTION  glob_sum_2D _a ***
137      !!
138      !! ** Purpose : perform a masked sum on the inner global domain of two 2D array
139      !!-----------------------------------------------------------------------
140      REAL(wp), INTENT(in), DIMENSION(:,:) ::   ptab1, ptab2    ! input 2D array
141      REAL(wp)            , DIMENSION(2)   ::   glob_sum_2d_a   ! global masked sum
142      !!-----------------------------------------------------------------------
143      !
144      glob_sum_2d_a(1) = SUM( ptab1(:,:)*tmask_i(:,:) )
145      glob_sum_2d_a(2) = SUM( ptab2(:,:)*tmask_i(:,:) )
146      IF( lk_mpp )   CALL mpp_sum( glob_sum_2d_a, 2 )
147      !
148   END FUNCTION glob_sum_2d_a
149
150
151   FUNCTION glob_sum_3d_a( ptab1, ptab2 )
152      !!-----------------------------------------------------------------------
153      !!                  ***  FUNCTION  glob_sum_3D_a ***
154      !!
155      !! ** Purpose : perform a masked sum on the inner global domain of two 3D array
156      !!-----------------------------------------------------------------------
157      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   ptab1, ptab2    ! input 3D array
158      REAL(wp)            , DIMENSION(2)     ::   glob_sum_3d_a   ! global masked sum
159      !!
160      INTEGER :: jk
161      INTEGER :: ijpk ! local variable: size of the 3d dimension of ptab
162      !!-----------------------------------------------------------------------
163      !
164      ijpk = SIZE(ptab1,3)
165      !
166      glob_sum_3d_a(:) = 0.e0
167      DO jk = 1, ijpk
168         glob_sum_3d_a(1) = glob_sum_3d_a(1) + SUM( ptab1(:,:,jk)*tmask_i(:,:) )
169         glob_sum_3d_a(2) = glob_sum_3d_a(2) + SUM( ptab2(:,:,jk)*tmask_i(:,:) )
170      END DO
171      IF( lk_mpp )   CALL mpp_sum( glob_sum_3d_a, 2 )
172      !
173   END FUNCTION glob_sum_3d_a
174
175   FUNCTION glob_sum_full_2d( ptab )
176      !!----------------------------------------------------------------------
177      !!                  ***  FUNCTION  glob_sum_full_2d ***
178      !!
179      !! ** Purpose : perform a sum in calling DDPDD routine (nomask)
180      !!----------------------------------------------------------------------
181      REAL(wp), INTENT(in), DIMENSION(:,:) ::   ptab
182      REAL(wp)                             ::   glob_sum_full_2d   ! global sum
183      !!
184      !!-----------------------------------------------------------------------
185      !
186      glob_sum_full_2d = SUM( ptab(:,:) * tmask_h(:,:) )
187      IF( lk_mpp )   CALL mpp_sum( glob_sum_full_2d )
188      !
189   END FUNCTION glob_sum_full_2d
190
191   FUNCTION glob_sum_full_3d( ptab )
192      !!----------------------------------------------------------------------
193      !!                  ***  FUNCTION  glob_sum_full_3d ***
194      !!
195      !! ** Purpose : perform a sum on a 3D array in calling DDPDD routine (nomask)
196      !!----------------------------------------------------------------------
197      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   ptab
198      REAL(wp)                               ::   glob_sum_full_3d   ! global sum
199      !!
200      INTEGER    ::   ji, jj, jk   ! dummy loop indices
201      INTEGER    ::   ijpk ! local variables: size of ptab
202      !!-----------------------------------------------------------------------
203      !
204      ijpk = SIZE(ptab,3)
205      !
206      glob_sum_full_3d = 0.e0
207      DO jk = 1, ijpk
208         glob_sum_full_3d = glob_sum_full_3d + SUM( ptab(:,:,jk) * tmask_h(:,:) )
209      END DO
210      IF( lk_mpp )   CALL mpp_sum( glob_sum_full_3d )
211      !
212   END FUNCTION glob_sum_full_3d
213
214
215#else 
216   !!----------------------------------------------------------------------
217   !!   'key_mpp_rep'                                   MPP reproducibility
218   !!----------------------------------------------------------------------
219   
220   ! --- SUM ---
221   FUNCTION glob_sum_1d( ptab, kdim )
222      !!----------------------------------------------------------------------
223      !!                  ***  FUNCTION  glob_sum_1d ***
224      !!
225      !! ** Purpose : perform a sum in calling DDPDD routine
226      !!----------------------------------------------------------------------
227      INTEGER , INTENT(in) :: kdim
228      REAL(wp), INTENT(in), DIMENSION(kdim) ::   ptab
229      REAL(wp)                              ::   glob_sum_1d   ! global sum
230      !!
231      COMPLEX(wp)::   ctmp
232      REAL(wp)   ::   ztmp
233      INTEGER    ::   ji   ! dummy loop indices
234      !!-----------------------------------------------------------------------
235      !
236      ztmp = 0.e0
237      ctmp = CMPLX( 0.e0, 0.e0, wp )
238      DO ji = 1, kdim
239         ztmp =  ptab(ji)
240         CALL DDPDD( CMPLX( ztmp, 0.e0, wp ), ctmp )
241         END DO
242      IF( lk_mpp )   CALL mpp_sum( ctmp )   ! sum over the global domain
243      glob_sum_1d = REAL(ctmp,wp)
244      !
245   END FUNCTION glob_sum_1d
246
247   FUNCTION glob_sum_2d( ptab )
248      !!----------------------------------------------------------------------
249      !!                  ***  FUNCTION  glob_sum_2d ***
250      !!
251      !! ** Purpose : perform a sum in calling DDPDD routine
252      !!----------------------------------------------------------------------
253      REAL(wp), INTENT(in), DIMENSION(:,:) ::   ptab
254      REAL(wp)                             ::   glob_sum_2d   ! global masked sum
255      !!
256      COMPLEX(wp)::   ctmp
257      REAL(wp)   ::   ztmp
258      INTEGER    ::   ji, jj   ! dummy loop indices
259      !!-----------------------------------------------------------------------
260      !
261      ztmp = 0.e0
262      ctmp = CMPLX( 0.e0, 0.e0, wp )
263      DO jj = 1, jpj
264         DO ji =1, jpi
265         ztmp =  ptab(ji,jj) * tmask_i(ji,jj)
266         CALL DDPDD( CMPLX( ztmp, 0.e0, wp ), ctmp )
267         END DO
268      END DO
269      IF( lk_mpp )   CALL mpp_sum( ctmp )   ! sum over the global domain
270      glob_sum_2d = REAL(ctmp,wp)
271      !
272   END FUNCTION glob_sum_2d
273
274
275   FUNCTION glob_sum_3d( ptab )
276      !!----------------------------------------------------------------------
277      !!                  ***  FUNCTION  glob_sum_3d ***
278      !!
279      !! ** Purpose : perform a sum on a 3D array in calling DDPDD routine
280      !!----------------------------------------------------------------------
281      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   ptab
282      REAL(wp)                               ::   glob_sum_3d   ! global masked sum
283      !!
284      COMPLEX(wp)::   ctmp
285      REAL(wp)   ::   ztmp
286      INTEGER    ::   ji, jj, jk   ! dummy loop indices
287      INTEGER    ::   ijpk ! local variables: size of ptab
288      !!-----------------------------------------------------------------------
289      !
290      ijpk = SIZE(ptab,3)
291      !
292      ztmp = 0.e0
293      ctmp = CMPLX( 0.e0, 0.e0, wp )
294      DO jk = 1, ijpk
295         DO jj = 1, jpj
296            DO ji =1, jpi
297            ztmp =  ptab(ji,jj,jk) * tmask_i(ji,jj)
298            CALL DDPDD( CMPLX( ztmp, 0.e0, wp ), ctmp )
299            END DO
300         END DO
301      END DO
302      IF( lk_mpp )   CALL mpp_sum( ctmp )   ! sum over the global domain
303      glob_sum_3d = REAL(ctmp,wp)
304      !
305   END FUNCTION glob_sum_3d
306
307
308   FUNCTION glob_sum_2d_a( ptab1, ptab2 )
309      !!----------------------------------------------------------------------
310      !!                  ***  FUNCTION  glob_sum_2d_a ***
311      !!
312      !! ** Purpose : perform a sum on two 2D arrays in calling DDPDD routine
313      !!----------------------------------------------------------------------
314      REAL(wp), INTENT(in), DIMENSION(:,:) ::   ptab1, ptab2
315      REAL(wp)                             ::   glob_sum_2d_a   ! global masked sum
316      !!
317      COMPLEX(wp)::   ctmp
318      REAL(wp)   ::   ztmp
319      INTEGER    ::   ji, jj   ! dummy loop indices
320      !!-----------------------------------------------------------------------
321      !
322      ztmp = 0.e0
323      ctmp = CMPLX( 0.e0, 0.e0, wp )
324      DO jj = 1, jpj
325         DO ji =1, jpi
326         ztmp =  ptab1(ji,jj) * tmask_i(ji,jj)
327         CALL DDPDD( CMPLX( ztmp, 0.e0, wp ), ctmp )
328         ztmp =  ptab2(ji,jj) * tmask_i(ji,jj)
329         CALL DDPDD( CMPLX( ztmp, 0.e0, wp ), ctmp )
330         END DO
331      END DO
332      IF( lk_mpp )   CALL mpp_sum( ctmp )   ! sum over the global domain
333      glob_sum_2d_a = REAL(ctmp,wp)
334      !
335   END FUNCTION glob_sum_2d_a
336
337
338   FUNCTION glob_sum_3d_a( ptab1, ptab2 )
339      !!----------------------------------------------------------------------
340      !!                  ***  FUNCTION  glob_sum_3d_a ***
341      !!
342      !! ** Purpose : perform a sum on two 3D array in calling DDPDD routine
343      !!----------------------------------------------------------------------
344      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   ptab1, ptab2
345      REAL(wp)                               ::   glob_sum_3d_a   ! global masked sum
346      !!
347      COMPLEX(wp)::   ctmp
348      REAL(wp)   ::   ztmp
349      INTEGER    ::   ji, jj, jk   ! dummy loop indices
350      INTEGER    ::   ijpk ! local variables: size of ptab
351      !!-----------------------------------------------------------------------
352      !
353      ijpk = SIZE(ptab1,3)
354      !
355      ztmp = 0.e0
356      ctmp = CMPLX( 0.e0, 0.e0, wp )
357      DO jk = 1, ijpk
358         DO jj = 1, jpj
359            DO ji = 1, jpi
360               ztmp =  ptab1(ji,jj,jk) * tmask_i(ji,jj)
361               CALL DDPDD( CMPLX( ztmp, 0.e0, wp ), ctmp )
362               ztmp =  ptab2(ji,jj,jk) * tmask_i(ji,jj)
363               CALL DDPDD( CMPLX( ztmp, 0.e0, wp ), ctmp )
364            END DO
365         END DO   
366      END DO
367      IF( lk_mpp )   CALL mpp_sum( ctmp )   ! sum over the global domain
368      glob_sum_3d_a = REAL(ctmp,wp)
369      !
370   END FUNCTION glob_sum_3d_a   
371
372   FUNCTION glob_sum_full_2d( ptab )
373      !!----------------------------------------------------------------------
374      !!                  ***  FUNCTION  glob_sum_full_2d ***
375      !!
376      !! ** Purpose : perform a sum in calling DDPDD routine
377      !!----------------------------------------------------------------------
378      REAL(wp), INTENT(in), DIMENSION(:,:) ::   ptab
379      REAL(wp)                             ::   glob_sum_full_2d   ! global sum (nomask)
380      !!
381      COMPLEX(wp)::   ctmp
382      REAL(wp)   ::   ztmp
383      INTEGER    ::   ji, jj   ! dummy loop indices
384      !!-----------------------------------------------------------------------
385      !
386      ztmp = 0.e0
387      ctmp = CMPLX( 0.e0, 0.e0, wp )
388      DO jj = 1, jpj
389         DO ji =1, jpi
390         ztmp =  ptab(ji,jj) * tmask_h(ji,jj) 
391         CALL DDPDD( CMPLX( ztmp, 0.e0, wp ), ctmp )
392         END DO
393      END DO
394      IF( lk_mpp )   CALL mpp_sum( ctmp )   ! sum over the global domain
395      glob_sum_full_2d = REAL(ctmp,wp)
396      !
397   END FUNCTION glob_sum_full_2d
398
399   FUNCTION glob_sum_full_3d( ptab )
400      !!----------------------------------------------------------------------
401      !!                  ***  FUNCTION  glob_sum_full_3d ***
402      !!
403      !! ** Purpose : perform a sum on a 3D array in calling DDPDD routine
404      !!----------------------------------------------------------------------
405      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   ptab
406      REAL(wp)                               ::   glob_sum_full_3d   ! global sum (nomask)
407      !!
408      COMPLEX(wp)::   ctmp
409      REAL(wp)   ::   ztmp
410      INTEGER    ::   ji, jj, jk   ! dummy loop indices
411      INTEGER    ::   ijpk ! local variables: size of ptab
412      !!-----------------------------------------------------------------------
413      !
414      ijpk = SIZE(ptab,3)
415      !
416      ztmp = 0.e0
417      ctmp = CMPLX( 0.e0, 0.e0, wp )
418      DO jk = 1, ijpk
419         DO jj = 1, jpj
420            DO ji =1, jpi
421            ztmp =  ptab(ji,jj,jk) * tmask_h(ji,jj)
422            CALL DDPDD( CMPLX( ztmp, 0.e0, wp ), ctmp )
423            END DO
424         END DO
425      END DO
426      IF( lk_mpp )   CALL mpp_sum( ctmp )   ! sum over the global domain
427      glob_sum_full_3d = REAL(ctmp,wp)
428      !
429   END FUNCTION glob_sum_full_3d
430
431
432
433#endif
434
435   ! --- MIN ---
436   FUNCTION glob_min_2d( ptab ) 
437      !!-----------------------------------------------------------------------
438      !!                  ***  FUNCTION  glob_min_2D  ***
439      !!
440      !! ** Purpose : perform a masked min on the inner global domain of a 2D array
441      !!-----------------------------------------------------------------------
442      REAL(wp), INTENT(in), DIMENSION(:,:) ::   ptab          ! input 2D array
443      REAL(wp)                             ::   glob_min_2d   ! global masked min
444      !!-----------------------------------------------------------------------
445      !
446      glob_min_2d = MINVAL( ptab(:,:)*tmask_i(:,:) )
447      IF( lk_mpp )   CALL mpp_min( glob_min_2d )
448      !
449   END FUNCTION glob_min_2d
450 
451   FUNCTION glob_min_3d( ptab ) 
452      !!-----------------------------------------------------------------------
453      !!                  ***  FUNCTION  glob_min_3D  ***
454      !!
455      !! ** Purpose : perform a masked min on the inner global domain of a 3D array
456      !!-----------------------------------------------------------------------
457      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   ptab          ! input 3D array
458      REAL(wp)                               ::   glob_min_3d   ! global masked min
459      !!
460      INTEGER :: jk
461      INTEGER :: ijpk ! local variable: size of the 3d dimension of ptab
462      !!-----------------------------------------------------------------------
463      !
464      ijpk = SIZE(ptab,3)
465      !
466      glob_min_3d = MINVAL( ptab(:,:,1)*tmask_i(:,:) )
467      DO jk = 2, ijpk
468         glob_min_3d = MIN( glob_min_3d, MINVAL( ptab(:,:,jk)*tmask_i(:,:) ) )
469      END DO
470      IF( lk_mpp )   CALL mpp_min( glob_min_3d )
471      !
472   END FUNCTION glob_min_3d
473
474
475   FUNCTION glob_min_2d_a( ptab1, ptab2 ) 
476      !!-----------------------------------------------------------------------
477      !!                  ***  FUNCTION  glob_min_2D _a ***
478      !!
479      !! ** Purpose : perform a masked min on the inner global domain of two 2D array
480      !!-----------------------------------------------------------------------
481      REAL(wp), INTENT(in), DIMENSION(:,:) ::   ptab1, ptab2    ! input 2D array
482      REAL(wp)            , DIMENSION(2)   ::   glob_min_2d_a   ! global masked min
483      !!-----------------------------------------------------------------------
484      !             
485      glob_min_2d_a(1) = MINVAL( ptab1(:,:)*tmask_i(:,:) )
486      glob_min_2d_a(2) = MINVAL( ptab2(:,:)*tmask_i(:,:) )
487      IF( lk_mpp )   CALL mpp_min( glob_min_2d_a, 2 )
488      !
489   END FUNCTION glob_min_2d_a
490 
491 
492   FUNCTION glob_min_3d_a( ptab1, ptab2 ) 
493      !!-----------------------------------------------------------------------
494      !!                  ***  FUNCTION  glob_min_3D_a ***
495      !!
496      !! ** Purpose : perform a masked min on the inner global domain of two 3D array
497      !!-----------------------------------------------------------------------
498      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   ptab1, ptab2    ! input 3D array
499      REAL(wp)            , DIMENSION(2)     ::   glob_min_3d_a   ! global masked min
500      !!
501      INTEGER :: jk
502      INTEGER :: ijpk ! local variable: size of the 3d dimension of ptab
503      !!-----------------------------------------------------------------------
504      !
505      ijpk = SIZE(ptab1,3)
506      !
507      glob_min_3d_a(1) = MINVAL( ptab1(:,:,1)*tmask_i(:,:) )
508      glob_min_3d_a(2) = MINVAL( ptab2(:,:,1)*tmask_i(:,:) )
509      DO jk = 2, ijpk
510         glob_min_3d_a(1) = MIN( glob_min_3d_a(1), MINVAL( ptab1(:,:,jk)*tmask_i(:,:) ) )
511         glob_min_3d_a(2) = MIN( glob_min_3d_a(2), MINVAL( ptab2(:,:,jk)*tmask_i(:,:) ) )
512      END DO
513      IF( lk_mpp )   CALL mpp_min( glob_min_3d_a, 2 )
514      !
515   END FUNCTION glob_min_3d_a
516
517   ! --- MAX ---
518   FUNCTION glob_max_2d( ptab ) 
519      !!-----------------------------------------------------------------------
520      !!                  ***  FUNCTION  glob_max_2D  ***
521      !!
522      !! ** Purpose : perform a masked max on the inner global domain of a 2D array
523      !!-----------------------------------------------------------------------
524      REAL(wp), INTENT(in), DIMENSION(:,:) ::   ptab          ! input 2D array
525      REAL(wp)                             ::   glob_max_2d   ! global masked max
526      !!-----------------------------------------------------------------------
527      !
528      glob_max_2d = MAXVAL( ptab(:,:)*tmask_i(:,:) )
529      IF( lk_mpp )   CALL mpp_max( glob_max_2d )
530      !
531   END FUNCTION glob_max_2d
532 
533   FUNCTION glob_max_3d( ptab ) 
534      !!-----------------------------------------------------------------------
535      !!                  ***  FUNCTION  glob_max_3D  ***
536      !!
537      !! ** Purpose : perform a masked max on the inner global domain of a 3D array
538      !!-----------------------------------------------------------------------
539      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   ptab          ! input 3D array
540      REAL(wp)                               ::   glob_max_3d   ! global masked max
541      !!
542      INTEGER :: jk
543      INTEGER :: ijpk ! local variable: size of the 3d dimension of ptab
544      !!-----------------------------------------------------------------------
545      !
546      ijpk = SIZE(ptab,3)
547      !
548      glob_max_3d = MAXVAL( ptab(:,:,1)*tmask_i(:,:) )
549      DO jk = 2, ijpk
550         glob_max_3d = MAX( glob_max_3d, MAXVAL( ptab(:,:,jk)*tmask_i(:,:) ) )
551      END DO
552      IF( lk_mpp )   CALL mpp_max( glob_max_3d )
553      !
554   END FUNCTION glob_max_3d
555
556
557   FUNCTION glob_max_2d_a( ptab1, ptab2 ) 
558      !!-----------------------------------------------------------------------
559      !!                  ***  FUNCTION  glob_max_2D _a ***
560      !!
561      !! ** Purpose : perform a masked max on the inner global domain of two 2D array
562      !!-----------------------------------------------------------------------
563      REAL(wp), INTENT(in), DIMENSION(:,:) ::   ptab1, ptab2    ! input 2D array
564      REAL(wp)            , DIMENSION(2)   ::   glob_max_2d_a   ! global masked max
565      !!-----------------------------------------------------------------------
566      !             
567      glob_max_2d_a(1) = MAXVAL( ptab1(:,:)*tmask_i(:,:) )
568      glob_max_2d_a(2) = MAXVAL( ptab2(:,:)*tmask_i(:,:) )
569      IF( lk_mpp )   CALL mpp_max( glob_max_2d_a, 2 )
570      !
571   END FUNCTION glob_max_2d_a
572 
573 
574   FUNCTION glob_max_3d_a( ptab1, ptab2 ) 
575      !!-----------------------------------------------------------------------
576      !!                  ***  FUNCTION  glob_max_3D_a ***
577      !!
578      !! ** Purpose : perform a masked max on the inner global domain of two 3D array
579      !!-----------------------------------------------------------------------
580      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   ptab1, ptab2    ! input 3D array
581      REAL(wp)            , DIMENSION(2)     ::   glob_max_3d_a   ! global masked max
582      !!
583      INTEGER :: jk
584      INTEGER :: ijpk ! local variable: size of the 3d dimension of ptab
585      !!-----------------------------------------------------------------------
586      !
587      ijpk = SIZE(ptab1,3)
588      !
589      glob_max_3d_a(1) = MAXVAL( ptab1(:,:,1)*tmask_i(:,:) )
590      glob_max_3d_a(2) = MAXVAL( ptab2(:,:,1)*tmask_i(:,:) )
591      DO jk = 2, ijpk
592         glob_max_3d_a(1) = MAX( glob_max_3d_a(1), MAXVAL( ptab1(:,:,jk)*tmask_i(:,:) ) )
593         glob_max_3d_a(2) = MAX( glob_max_3d_a(2), MAXVAL( ptab2(:,:,jk)*tmask_i(:,:) ) )
594      END DO
595      IF( lk_mpp )   CALL mpp_max( glob_max_3d_a, 2 )
596      !
597   END FUNCTION glob_max_3d_a
598
599
600   SUBROUTINE DDPDD( ydda, yddb )
601      !!----------------------------------------------------------------------
602      !!               ***  ROUTINE DDPDD ***
603      !!
604      !! ** Purpose : Add a scalar element to a sum
605      !!
606      !!
607      !! ** Method  : The code uses the compensated summation with doublet
608      !!              (sum,error) emulated useing complex numbers. ydda is the
609      !!               scalar to add to the summ yddb
610      !!
611      !! ** Action  : This does only work for MPI.
612      !!
613      !! References : Using Acurate Arithmetics to Improve Numerical
614      !!              Reproducibility and Sability in Parallel Applications
615      !!              Yun HE and Chris H. Q. DING, Journal of Supercomputing 18, 259-277, 2001
616      !!----------------------------------------------------------------------
617      COMPLEX(wp), INTENT(in   ) ::   ydda
618      COMPLEX(wp), INTENT(inout) ::   yddb
619      !
620      REAL(wp) :: zerr, zt1, zt2  ! local work variables
621      !!-----------------------------------------------------------------------
622      !
623      ! Compute ydda + yddb using Knuth's trick.
624      zt1  = REAL(ydda) + REAL(yddb)
625      zerr = zt1 - REAL(ydda)
626      zt2  = ( (REAL(yddb) - zerr) + (REAL(ydda) - (zt1 - zerr)) )   &
627         &   + AIMAG(ydda)         + AIMAG(yddb)
628      !
629      ! The result is t1 + t2, after normalization.
630      yddb = CMPLX( zt1 + zt2, zt2 - ((zt1 + zt2) - zt1), wp )
631      !
632   END SUBROUTINE DDPDD
633
634#if defined key_nosignedzero
635   !!----------------------------------------------------------------------
636   !!   'key_nosignedzero'                                         F90 SIGN
637   !!----------------------------------------------------------------------
638
639   FUNCTION SIGN_SCALAR( pa, pb )
640      !!-----------------------------------------------------------------------
641      !!                  ***  FUNCTION SIGN_SCALAR  ***
642      !!
643      !! ** Purpose : overwrite f95 behaviour of intrinsinc sign function
644      !!-----------------------------------------------------------------------
645      REAL(wp) :: pa,pb          ! input
646      REAL(wp) :: SIGN_SCALAR    ! result
647      !!-----------------------------------------------------------------------
648      IF ( pb >= 0.e0) THEN   ;   SIGN_SCALAR = ABS(pa)
649      ELSE                    ;   SIGN_SCALAR =-ABS(pa)
650      ENDIF
651   END FUNCTION SIGN_SCALAR
652
653
654   FUNCTION SIGN_ARRAY_1D( pa, pb )
655      !!-----------------------------------------------------------------------
656      !!                  ***  FUNCTION SIGN_ARRAY_1D  ***
657      !!
658      !! ** Purpose : overwrite f95 behaviour of intrinsinc sign function
659      !!-----------------------------------------------------------------------
660      REAL(wp) :: pa,pb(:)                   ! input
661      REAL(wp) :: SIGN_ARRAY_1D(SIZE(pb,1))  ! result
662      !!-----------------------------------------------------------------------
663      WHERE ( pb >= 0.e0 )   ;   SIGN_ARRAY_1D = ABS(pa)
664      ELSEWHERE              ;   SIGN_ARRAY_1D =-ABS(pa)
665      END WHERE
666   END FUNCTION SIGN_ARRAY_1D
667
668
669   FUNCTION SIGN_ARRAY_2D(pa,pb)
670      !!-----------------------------------------------------------------------
671      !!                  ***  FUNCTION SIGN_ARRAY_2D  ***
672      !!
673      !! ** Purpose : overwrite f95 behaviour of intrinsinc sign function
674      !!-----------------------------------------------------------------------
675      REAL(wp) :: pa,pb(:,:)      ! input
676      REAL(wp) :: SIGN_ARRAY_2D(SIZE(pb,1),SIZE(pb,2))  ! result
677      !!-----------------------------------------------------------------------
678      WHERE ( pb >= 0.e0 )   ;   SIGN_ARRAY_2D = ABS(pa)
679      ELSEWHERE              ;   SIGN_ARRAY_2D =-ABS(pa)
680      END WHERE
681   END FUNCTION SIGN_ARRAY_2D
682
683   FUNCTION SIGN_ARRAY_3D(pa,pb)
684      !!-----------------------------------------------------------------------
685      !!                  ***  FUNCTION SIGN_ARRAY_3D  ***
686      !!
687      !! ** Purpose : overwrite f95 behaviour of intrinsinc sign function
688      !!-----------------------------------------------------------------------
689      REAL(wp) :: pa,pb(:,:,:)      ! input
690      REAL(wp) :: SIGN_ARRAY_3D(SIZE(pb,1),SIZE(pb,2),SIZE(pb,3))  ! result
691      !!-----------------------------------------------------------------------
692      WHERE ( pb >= 0.e0 )   ;   SIGN_ARRAY_3D = ABS(pa)
693      ELSEWHERE              ;   SIGN_ARRAY_3D =-ABS(pa)
694      END WHERE
695   END FUNCTION SIGN_ARRAY_3D
696
697
698   FUNCTION SIGN_ARRAY_1D_A(pa,pb)
699      !!-----------------------------------------------------------------------
700      !!                  ***  FUNCTION SIGN_ARRAY_1D_A  ***
701      !!
702      !! ** Purpose : overwrite f95 behaviour of intrinsinc sign function
703      !!-----------------------------------------------------------------------
704      REAL(wp) :: pa(:),pb(:)      ! input
705      REAL(wp) :: SIGN_ARRAY_1D_A(SIZE(pb,1))  ! result
706      !!-----------------------------------------------------------------------
707      WHERE ( pb >= 0.e0 )   ;   SIGN_ARRAY_1D_A = ABS(pa)
708      ELSEWHERE              ;   SIGN_ARRAY_1D_A =-ABS(pa)
709      END WHERE
710   END FUNCTION SIGN_ARRAY_1D_A
711
712
713   FUNCTION SIGN_ARRAY_2D_A(pa,pb)
714      !!-----------------------------------------------------------------------
715      !!                  ***  FUNCTION SIGN_ARRAY_2D_A  ***
716      !!
717      !! ** Purpose : overwrite f95 behaviour of intrinsinc sign function
718      !!-----------------------------------------------------------------------
719      REAL(wp) :: pa(:,:),pb(:,:)      ! input
720      REAL(wp) :: SIGN_ARRAY_2D_A(SIZE(pb,1),SIZE(pb,2))  ! result
721      !!-----------------------------------------------------------------------
722      WHERE ( pb >= 0.e0 )   ;   SIGN_ARRAY_2D_A = ABS(pa)
723      ELSEWHERE              ;   SIGN_ARRAY_2D_A =-ABS(pa)
724      END WHERE
725   END FUNCTION SIGN_ARRAY_2D_A
726
727
728   FUNCTION SIGN_ARRAY_3D_A(pa,pb)
729      !!-----------------------------------------------------------------------
730      !!                  ***  FUNCTION SIGN_ARRAY_3D_A  ***
731      !!
732      !! ** Purpose : overwrite f95 behaviour of intrinsinc sign function
733      !!-----------------------------------------------------------------------
734      REAL(wp) :: pa(:,:,:),pb(:,:,:)  ! input
735      REAL(wp) :: SIGN_ARRAY_3D_A(SIZE(pb,1),SIZE(pb,2),SIZE(pb,3)) ! result
736      !!-----------------------------------------------------------------------
737      WHERE ( pb >= 0.e0 )   ;   SIGN_ARRAY_3D_A = ABS(pa)
738      ELSEWHERE              ;   SIGN_ARRAY_3D_A =-ABS(pa)
739      END WHERE
740   END FUNCTION SIGN_ARRAY_3D_A
741
742
743   FUNCTION SIGN_ARRAY_1D_B(pa,pb)
744      !!-----------------------------------------------------------------------
745      !!                  ***  FUNCTION SIGN_ARRAY_1D_B  ***
746      !!
747      !! ** Purpose : overwrite f95 behaviour of intrinsinc sign function
748      !!-----------------------------------------------------------------------
749      REAL(wp) :: pa(:),pb      ! input
750      REAL(wp) :: SIGN_ARRAY_1D_B(SIZE(pa,1))  ! result
751      !!-----------------------------------------------------------------------
752      IF( pb >= 0.e0 ) THEN   ;   SIGN_ARRAY_1D_B = ABS(pa)
753      ELSE                    ;   SIGN_ARRAY_1D_B =-ABS(pa)
754      ENDIF
755   END FUNCTION SIGN_ARRAY_1D_B
756
757
758   FUNCTION SIGN_ARRAY_2D_B(pa,pb)
759      !!-----------------------------------------------------------------------
760      !!                  ***  FUNCTION SIGN_ARRAY_2D_B  ***
761      !!
762      !! ** Purpose : overwrite f95 behaviour of intrinsinc sign function
763      !!-----------------------------------------------------------------------
764      REAL(wp) :: pa(:,:),pb      ! input
765      REAL(wp) :: SIGN_ARRAY_2D_B(SIZE(pa,1),SIZE(pa,2))  ! result
766      !!-----------------------------------------------------------------------
767      IF( pb >= 0.e0 ) THEN   ;   SIGN_ARRAY_2D_B = ABS(pa)
768      ELSE                    ;   SIGN_ARRAY_2D_B =-ABS(pa)
769      ENDIF
770   END FUNCTION SIGN_ARRAY_2D_B
771
772
773   FUNCTION SIGN_ARRAY_3D_B(pa,pb)
774      !!-----------------------------------------------------------------------
775      !!                  ***  FUNCTION SIGN_ARRAY_3D_B  ***
776      !!
777      !! ** Purpose : overwrite f95 behaviour of intrinsinc sign function
778      !!-----------------------------------------------------------------------
779      REAL(wp) :: pa(:,:,:),pb      ! input
780      REAL(wp) :: SIGN_ARRAY_3D_B(SIZE(pa,1),SIZE(pa,2),SIZE(pa,3))  ! result
781      !!-----------------------------------------------------------------------
782      IF( pb >= 0.e0 ) THEN   ;   SIGN_ARRAY_3D_B = ABS(pa)
783      ELSE                    ;   SIGN_ARRAY_3D_B =-ABS(pa)
784      ENDIF
785   END FUNCTION SIGN_ARRAY_3D_B
786#endif
787
788   !!======================================================================
789END MODULE lib_fortran
Note: See TracBrowser for help on using the repository browser.