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.
Changeset 2715 for trunk/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 – NEMO

Ignore:
Timestamp:
2011-03-30T17:58:35+02:00 (13 years ago)
Author:
rblod
Message:

First attempt to put dynamic allocation on the trunk

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90

    r2605 r2715  
    4141   LOGICAL , PUBLIC, PARAMETER ::   lk_ldfslp = .TRUE.     !: slopes flag 
    4242   !                                                                             !! Madec operator 
    43    REAL(wp), PUBLIC, DIMENSION(:,:,:)    , ALLOCATABLE ::   uslp, wslpi          !: i_slope at U- and W-points 
    44    REAL(wp), PUBLIC, DIMENSION(:,:,:)    , ALLOCATABLE ::   vslp, wslpj          !: j-slope at V- and W-points 
    45    !                                                                             !! Griffies operator 
    46    REAL(wp), PUBLIC, DIMENSION(:,:,:)    , ALLOCATABLE ::   wslp2                !: wslp**2 from Griffies quarter cells 
    47    REAL(wp), PUBLIC, DIMENSION(:,:,:,:,:), ALLOCATABLE ::   triadi_g, triadj_g   !: skew flux  slopes relative to geopotentials  
    48    REAL(wp), PUBLIC, DIMENSION(:,:,:,:,:), ALLOCATABLE ::   triadi  , triadj     !: isoneutral slopes relative to model-coordinate 
     43   !  Arrays allocated in ldf_slp_init() routine once we know whether we're using the Griffies or Madec operator 
     44   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)     ::   uslp, wslpi          !: i_slope at U- and W-points 
     45   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)     ::   vslp, wslpj          !: j-slope at V- and W-points 
     46   !                                                                !! Griffies operator 
     47   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)     ::   wslp2                !: wslp**2 from Griffies quarter cells 
     48   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:,:) ::   triadi_g, triadj_g   !: skew flux  slopes relative to geopotentials  
     49   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:,:) ::   triadi  , triadj     !: isoneutral slopes relative to model-coordinate 
    4950 
    5051   !                                                              !! Madec operator 
    51    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   omlmask           ! mask of the surface mixed layer at T-pt 
    52    REAL(wp), DIMENSION(:,:)  , ALLOCATABLE ::   uslpml, wslpiml   ! i_slope at U- and W-points just below the mixed layer 
    53    REAL(wp), DIMENSION(:,:)  , ALLOCATABLE ::   vslpml, wslpjml   ! j_slope at V- and W-points just below the mixed layer 
     52   !  Arrays allocated in ldf_slp_init() routine once we know whether we're using the Griffies or Madec operator 
     53   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   omlmask           ! mask of the surface mixed layer at T-pt 
     54   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   uslpml, wslpiml   ! i_slope at U- and W-points just below the mixed layer 
     55   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   vslpml, wslpjml   ! j_slope at V- and W-points just below the mixed layer 
    5456 
    5557   REAL(wp) ::   repsln = 1.e-25_wp       ! tiny value used as minium of di(rho), dj(rho) and dk(rho) 
     58 
     59   ! Workspace arrays for ldf_slp_grif. These could be replaced by several 3D and 2D workspace 
     60   ! arrays from the wrk_nemo module with a bit of code re-writing. The 4D workspace  
     61   ! arrays can't be used here because of the zero-indexing of some of the ranks. ARPDBG. 
     62   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   zdzrho , zdyrho, zdxrho     ! Horizontal and vertical density gradients 
     63   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   zti_mlb, ztj_mlb            ! for Griffies operator only 
    5664 
    5765   !! * Substitutions 
     
    6169#  include "vectopt_loop_substitute.h90" 
    6270   !!---------------------------------------------------------------------- 
    63    !! NEMO/OPA 3.3 , NEMO Consortium (2010) 
     71   !! NEMO/OPA 4.0 , NEMO Consortium (2011) 
    6472   !! $Id$ 
    6573   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt) 
    6674   !!---------------------------------------------------------------------- 
    6775CONTAINS 
     76 
     77   INTEGER FUNCTION ldf_slp_alloc() 
     78      !!---------------------------------------------------------------------- 
     79      !!              ***  FUNCTION ldf_slp_alloc  *** 
     80      !!---------------------------------------------------------------------- 
     81      ! 
     82      ALLOCATE( zdxrho (jpi,jpj,jpk,0:1) , zti_mlb(jpi,jpj,0:1,0:1) ,     & 
     83         &      zdyrho (jpi,jpj,jpk,0:1) , ztj_mlb(jpi,jpj,0:1,0:1) ,     & 
     84         &      zdzrho (jpi,jpj,jpk,0:1)                            , STAT=ldf_slp_alloc ) 
     85         ! 
     86      IF( lk_mpp             )   CALL mpp_sum ( ldf_slp_alloc ) 
     87      IF( ldf_slp_alloc /= 0 )   CALL ctl_warn('ldf_slp_alloc : failed to allocate arrays.') 
     88      ! 
     89   END FUNCTION ldf_slp_alloc 
     90 
    6891 
    6992   SUBROUTINE ldf_slp( kt, prd, pn2 ) 
     
    92115      !!               of now neutral surfaces at u-, w- and v- w-points, resp. 
    93116      !!---------------------------------------------------------------------- 
    94       USE oce , zgru  => ua   ! use ua as workspace 
    95       USE oce , zgrv  => va   ! use va as workspace 
    96       USE oce , zww   => ta   ! use ta as workspace 
    97       USE oce , zwz   => sa   ! use sa as workspace 
    98       !! 
    99       INTEGER , INTENT(in)                         ::   kt    ! ocean time-step index 
    100       REAL(wp), INTENT(in), DIMENSION(jpi,jpj,jpk) ::   prd   ! in situ density 
    101       REAL(wp), INTENT(in), DIMENSION(jpi,jpj,jpk) ::   pn2   ! Brunt-Vaisala frequency (locally ref.) 
     117      USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released 
     118      USE oce     , ONLY:   zgru => ua       , zww => va   ! (ua,va) used as workspace 
     119      USE oce     , ONLY:   zgrv => ta       , zwz => sa   ! (ta,sa) used as workspace 
     120      USE wrk_nemo, ONLY:   zdzr => wrk_3d_1               ! 3D workspace 
     121      !! 
     122      INTEGER , INTENT(in)                   ::   kt    ! ocean time-step index 
     123      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   prd   ! in situ density 
     124      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   pn2   ! Brunt-Vaisala frequency (locally ref.) 
    102125      !! 
    103126      INTEGER  ::   ji , jj , jk    ! dummy loop indices 
     
    108131      REAL(wp) ::   zcj, zfj, zav, zbv, zaj, zbj   !   -      - 
    109132      REAL(wp) ::   zck, zfk,      zbw             !   -      - 
    110       REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zdzr   ! 3D workspace 
    111       !!---------------------------------------------------------------------- 
    112        
     133      !!---------------------------------------------------------------------- 
     134 
     135      IF( wrk_in_use(3, 1) ) THEN 
     136         CALL ctl_stop('ldf_slp: requested workspace arrays are unavailable')   ;   RETURN 
     137      ENDIF 
     138 
    113139      zeps   =  1.e-20_wp        !==   Local constant initialization   ==! 
    114140      z1_16  =  1.0_wp / 16._wp 
     
    342368      ENDIF 
    343369 
    344        
    345370      ! IV. Lateral boundary conditions  
    346371      ! =============================== 
     
    354379      ENDIF 
    355380      ! 
     381      IF( wrk_not_released(3, 1) )   CALL ctl_stop('ldf_slp: failed to release workspace arrays') 
     382      ! 
    356383   END SUBROUTINE ldf_slp 
    357384    
     
    371398      !!             - wslp2              squared slope of neutral surfaces at w-points. 
    372399      !!---------------------------------------------------------------------- 
    373       USE oce,   zdit  => ua   ! use ua as workspace 
    374       USE oce,   zdis  => va   ! use va as workspace 
    375       USE oce,   zdjt  => ta   ! use ta as workspace 
    376       USE oce,   zdjs  => sa   ! use sa as workspace 
    377       !! 
    378       INTEGER, INTENT( in ) ::   kt         ! ocean time-step index 
    379       !! 
     400      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released 
     401      USE oce     , ONLY:   zdit    => ua       , zdis   => va         ! (ua,va) used as workspace 
     402      USE oce     , ONLY:   zdjt    => ta       , zdjs   => sa         ! (ta,sa) used as workspace 
     403      USE wrk_nemo, ONLY:   zdkt    => wrk_3d_2 , zdks   => wrk_3d_3   ! 3D workspace 
     404      USE wrk_nemo, ONLY:   zalpha  => wrk_3d_4 , zbeta => wrk_3d_5    ! alpha, beta at T points, at depth fsgdept 
     405      USE wrk_nemo, ONLY:   z1_mlbw => wrk_2d_1 
     406      ! 
     407      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index 
     408      ! 
    380409      INTEGER  ::   ji, jj, jk, jl, ip, jp, kp  ! dummy loop indices 
    381       INTEGER  ::   iku, ikv                ! temporary integer 
     410      INTEGER  ::   iku, ikv                                  ! local integer 
    382411      REAL(wp) ::   zfacti, zfactj, zatempw,zatempu,zatempv   ! local scalars 
    383       REAL(wp) ::   zbu, zbv, zbti, zbtj 
     412      REAL(wp) ::   zbu, zbv, zbti, zbtj                      !   -      - 
    384413      REAL(wp) ::   zdxrho_raw, zti_coord, zti_raw, zti_lim, zti_lim2, zti_g_raw, zti_g_lim 
    385414      REAL(wp) ::   zdyrho_raw, ztj_coord, ztj_raw, ztj_lim, ztj_lim2, ztj_g_raw, ztj_g_lim 
    386415      REAL(wp) ::   zdzrho_raw 
    387       REAL(wp), DIMENSION(jpi,jpj,jpk,0:1) ::   zdzrho, zdyrho, zdxrho     ! Horizontal and vertical density gradients 
    388       REAL(wp), DIMENSION(jpi,jpj,0:1,0:1) ::   zti_mlb, ztj_mlb 
    389       REAL(wp), DIMENSION(jpi,jpj,jpk)     ::   zdkt, zdks 
    390       REAL(wp), DIMENSION(jpi,jpj,jpk)     ::   zalpha, zbeta       ! alpha, beta at T points, at depth fsgdept 
    391       REAL(wp), DIMENSION(jpi,jpj)         ::   z1_mlbw 
    392       !!---------------------------------------------------------------------- 
     416      !!---------------------------------------------------------------------- 
     417 
     418      IF( wrk_in_use(3, 2,3,4,5) .OR. wrk_in_use(2, 1) )THEN 
     419         CALL ctl_stop('ldf_slp_grif: requested workspace arrays are unavailable')   ;   RETURN 
     420      ENDIF 
    393421 
    394422      !--------------------------------! 
     
    572600      CALL lbc_lnk( wslp2, 'W', 1. )      ! lateral boundary confition on wslp2 only   ==>>> gm : necessary ? to be checked 
    573601      ! 
     602      IF( wrk_not_released(3, 2,3,4,5) .OR.   & 
     603          wrk_not_released(2, 1)       )   CALL ctl_stop('ldf_slp_grif: failed to release workspace arrays') 
     604      ! 
    574605   END SUBROUTINE ldf_slp_grif 
    575606 
     
    591622      !!                omlmask         :  mixed layer mask 
    592623      !!---------------------------------------------------------------------- 
    593       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in) ::   prd            ! in situ density 
    594       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in) ::   pn2            ! Brunt-Vaisala frequency (locally ref.) 
    595       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in) ::   p_gru, p_grv   ! i- & j-gradient of density (u- & v-pts) 
    596       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in) ::   p_dzr          ! z-gradient of density      (T-point) 
     624      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   prd            ! in situ density 
     625      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   pn2            ! Brunt-Vaisala frequency (locally ref.) 
     626      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   p_gru, p_grv   ! i- & j-gradient of density (u- & v-pts) 
     627      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   p_dzr          ! z-gradient of density      (T-point) 
    597628      !! 
    598629      INTEGER  ::   ji , jj , jk         ! dummy loop indices 
     
    704735      !! 
    705736      !! ** Method  :   read the nammbf namelist and check the parameter  
    706       !!      values called by tra_dmp at the first timestep (nit000) 
     737      !!              values called by tra_dmp at the first timestep (nit000) 
    707738      !!---------------------------------------------------------------------- 
    708739      INTEGER ::   ji, jj, jk   ! dummy loop indices 
     
    719750         ALLOCATE( triadi_g(jpi,jpj,jpk,0:1,0:1) , triadj_g(jpi,jpj,jpk,0:1,0:1) , wslp2(jpi,jpj,jpk) , STAT=ierr ) 
    720751         ALLOCATE( triadi  (jpi,jpj,jpk,0:1,0:1) , triadj  (jpi,jpj,jpk,0:1,0:1)                      , STAT=ierr ) 
    721          IF( ierr > 0 ) THEN 
    722             CALL ctl_stop( 'ldf_slp_init : unable to allocate Griffies operator slope ' )   ;   RETURN 
    723          ENDIF 
     752         IF( ierr > 0             )   CALL ctl_stop( 'STOP', 'ldf_slp_init : unable to allocate Griffies operator slope' ) 
     753         IF( ldf_slp_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'ldf_slp_init : unable to allocate workspace arrays' ) 
    724754         ! 
    725755         IF( ln_dynldf_iso )   CALL ctl_stop( 'ldf_slp_init: Griffies operator on momentum not supported' ) 
    726756         ! 
    727          IF( ( ln_traldf_hor .AND. ln_dynldf_hor ) .AND. ln_sco )   & 
    728             &     CALL ctl_stop( 'ldf_slp_init: horizontal Griffies operator ',   & 
    729             &                                              'in s-coordinate not supported' ) 
     757         IF( ( ln_traldf_hor .OR. ln_dynldf_hor ) .AND. ln_sco )   & 
     758            CALL ctl_stop( 'ldf_slp_init: horizontal Griffies operator in s-coordinate not supported' ) 
    730759         ! 
    731760      ELSE                             ! Madec operator : slopes at u-, v-, and w-points 
    732761         ALLOCATE( uslp(jpi,jpj,jpk) , vslp(jpi,jpj,jpk) , wslpi(jpi,jpj,jpk) , wslpj(jpi,jpj,jpk) ,                & 
    733762            &   omlmask(jpi,jpj,jpk) , uslpml(jpi,jpj)   , vslpml(jpi,jpj)    , wslpiml(jpi,jpj)   , wslpjml(jpi,jpj) , STAT=ierr ) 
    734          IF( ierr > 0 ) THEN 
    735             CALL ctl_stop( 'ldf_slp_init : unable to allocate Madec operator slope ' )   ;   RETURN 
    736          ENDIF 
     763         IF( ierr > 0 )   CALL ctl_stop( 'STOP', 'ldf_slp_init : unable to allocate Madec operator slope ' ) 
    737764 
    738765         ! Direction of lateral diffusion (tracers and/or momentum) 
     
    745772!!gm I no longer understand this..... 
    746773         IF( (ln_traldf_hor .OR. ln_dynldf_hor) .AND. .NOT. (lk_vvl .AND. ln_rstart) ) THEN 
    747             IF(lwp) THEN 
    748                WRITE(numout,*) '          Horizontal mixing in s-coordinate: slope = slope of s-surfaces' 
    749             ENDIF 
     774            IF(lwp)   WRITE(numout,*) '          Horizontal mixing in s-coordinate: slope = slope of s-surfaces' 
    750775 
    751776            ! geopotential diffusion in s-coordinates on tracers and/or momentum 
     
    765790               END DO 
    766791            END DO 
    767             ! Lateral boundary conditions on the slopes 
    768             CALL lbc_lnk( uslp , 'U', -1. )      ;      CALL lbc_lnk( vslp , 'V', -1. ) 
    769             CALL lbc_lnk( wslpi, 'W', -1. )      ;      CALL lbc_lnk( wslpj, 'W', -1. ) 
     792            CALL lbc_lnk( uslp , 'U', -1. )   ;   CALL lbc_lnk( vslp , 'V', -1. )      ! Lateral boundary conditions 
     793            CALL lbc_lnk( wslpi, 'W', -1. )   ;   CALL lbc_lnk( wslpj, 'W', -1. ) 
    770794         ENDIF 
    771       ENDIF      ! 
     795      ENDIF 
     796      ! 
    772797   END SUBROUTINE ldf_slp_init 
    773798 
Note: See TracChangeset for help on using the changeset viewer.