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/UKMO/dev_r5518_GO6_package_OMP/NEMOGCM/NEMO/OPA_SRC – NEMO

source: branches/UKMO/dev_r5518_GO6_package_OMP/NEMOGCM/NEMO/OPA_SRC/lib_fortran.F90 @ 9176

Last change on this file since 9176 was 9176, checked in by andmirek, 6 years ago

#2001: OMP directives

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