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 2590 for branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OPA_SRC/OBC – NEMO

Ignore:
Timestamp:
2011-02-18T13:49:27+01:00 (13 years ago)
Author:
trackstand2
Message:

Merge branch 'dynamic_memory' into master-svn-dyn

Location:
branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OPA_SRC/OBC
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OPA_SRC/OBC/obc_oce.F90

    r2528 r2590  
    7171   REAL(wp), PUBLIC ::   obcsurftot       !: Total lateral surface of open boundaries 
    7272    
    73    REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: &  !: 
     73   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: &  !: 
    7474      obctmsk,            &  !: mask array identical to tmask, execpt along OBC where it is set to 0 
    7575      !                      !  it used to calculate the cumulate flux E-P in the obcvol.F90 routine 
     
    8787   INTEGER ::   nje1m2, nje0m1    !: do loop index in mpp case for jpjefm1-1,jpjed 
    8888 
    89    REAL(wp), DIMENSION(jpj) ::   &  !: 
     89   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:) ::   &  !: 
    9090      sshfoe,           & !: now climatology of the east boundary sea surface height 
    9191      ubtfoe,vbtfoe       !: now climatology of the east boundary barotropic transport 
    9292      
    93    REAL(wp), DIMENSION(jpj,jpk) ::   &  !: 
     93   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   &  !: 
    9494      ufoe, vfoe,       & !: now climatology of the east boundary velocities  
    9595      tfoe, sfoe,       & !: now climatology of the east boundary temperature and salinity 
     
    9797      !                   ! in the obcdyn.F90 routine 
    9898 
    99    REAL(wp), DIMENSION(jpi,jpj) ::   sshfoe_b      !: east boundary ssh correction averaged over the barotropic loop 
     99   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshfoe_b      !: east boundary ssh correction averaged over the barotropic loop 
    100100      !                                            !  (if Flather's algoritm applied at open boundary) 
    101101 
     
    103103   !! Arrays for radiative East OBC:  
    104104   !!------------------------------- 
    105    REAL(wp), DIMENSION(jpj,jpk,3,3) ::   uebnd, vebnd      !: baroclinic u & v component of the velocity over 3 rows  
     105   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   uebnd, vebnd      !: baroclinic u & v component of the velocity over 3 rows  
    106106      !                                                    !  and 3 time step (now, before, and before before) 
    107    REAL(wp), DIMENSION(jpj,jpk,2,2) ::   tebnd, sebnd      !: East boundary temperature and salinity over 2 rows  
     107   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   tebnd, sebnd      !: East boundary temperature and salinity over 2 rows  
    108108      !                                                    !  and 2 time step (now and before) 
    109    REAL(wp), DIMENSION(jpj,jpk) ::   u_cxebnd, v_cxebnd    !: Zonal component of the phase speed ratio computed with  
     109   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   u_cxebnd, v_cxebnd    !: Zonal component of the phase speed ratio computed with  
    110110      !                                                    !  radiation of u and v velocity (respectively) at the  
    111111      !                                                    !  east open boundary (u_cxebnd = cx rdt ) 
    112    REAL(wp), DIMENSION(jpj,jpk) ::   uemsk, vemsk, temsk   !: 2D mask for the East OB 
     112   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   uemsk, vemsk, temsk   !: 2D mask for the East OB 
    113113 
    114114   ! Note that those arrays are optimized for mpp case  
     
    124124   INTEGER ::   njw1m2, njw0m1     !: do loop index in mpp case for jpjwfm2,jpjwd 
    125125 
    126    REAL(wp), DIMENSION(jpj) ::   &  !: 
     126   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:) ::   &  !: 
    127127      sshfow,           & !: now climatology of the west boundary sea surface height 
    128128      ubtfow,vbtfow       !: now climatology of the west boundary barotropic transport 
    129129 
    130    REAL(wp), DIMENSION(jpj,jpk) ::   &  !: 
     130   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   &  !: 
    131131      ufow, vfow,       & !: now climatology of the west velocities  
    132132      tfow, sfow,       & !: now climatology of the west temperature and salinity 
     
    134134      !                   !  in the obcdyn.F90 routine 
    135135 
    136    REAL(wp), DIMENSION(jpi,jpj) ::   sshfow_b    !: west boundary ssh correction averaged over the barotropic loop 
     136   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshfow_b    !: west boundary ssh correction averaged over the barotropic loop 
    137137      !                                          !  (if Flather's algoritm applied at open boundary) 
    138138 
     
    140140   !! Arrays for radiative West OBC 
    141141   !!------------------------------- 
    142    REAL(wp), DIMENSION(jpj,jpk,3,3) ::   uwbnd, vwbnd     !: baroclinic u & v components of the velocity over 3 rows  
     142   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   uwbnd, vwbnd     !: baroclinic u & v components of the velocity over 3 rows  
    143143      !                                                   !  and 3 time step (now, before, and before before) 
    144    REAL(wp), DIMENSION(jpj,jpk,2,2) ::   twbnd, swbnd     !: west boundary temperature and salinity over 2 rows and  
     144   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   twbnd, swbnd     !: west boundary temperature and salinity over 2 rows and  
    145145      !                                                   !  2 time step (now and before) 
    146    REAL(wp), DIMENSION(jpj,jpk) ::   u_cxwbnd, v_cxwbnd   !: Zonal component of the phase speed ratio computed with  
     146   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   u_cxwbnd, v_cxwbnd   !: Zonal component of the phase speed ratio computed with  
    147147      !                                                   !  radiation of zonal and meridional velocity (respectively)  
    148148      !                                                   !  at the west open boundary (u_cxwbnd = cx rdt ) 
    149    REAL(wp), DIMENSION(jpj,jpk) ::   uwmsk, vwmsk, twmsk  !: 2D mask for the West OB 
     149   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   uwmsk, vwmsk, twmsk  !: 2D mask for the West OB 
    150150 
    151151   ! Note that those arrays are optimized for mpp case  
     
    162162   INTEGER ::   njn0m1, njn1m1     !: do loop index in mpp case for jpnob-1 
    163163 
    164    REAL(wp), DIMENSION(jpi) ::   &  !: 
     164   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:) ::   &  !: 
    165165      sshfon,           & !: now climatology of the north boundary sea surface height 
    166166      ubtfon,vbtfon       !: now climatology of the north boundary barotropic transport 
    167167 
    168    REAL(wp), DIMENSION(jpi,jpk) ::   &    !: 
     168   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   &    !: 
    169169      ufon, vfon,       & !: now climatology of the north boundary velocities 
    170170      tfon, sfon,       & !: now climatology of the north boundary temperature and salinity 
     
    172172      !                   !  in yhe obcdyn.F90 routine 
    173173 
    174    REAL(wp), DIMENSION(jpi,jpj) ::   sshfon_b      !: north boundary ssh correction averaged over the barotropic loop 
     174   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshfon_b      !: north boundary ssh correction averaged over the barotropic loop 
    175175      !                                            !  (if Flather's algoritm applied at open boundary) 
    176176 
     
    178178   !! Arrays for radiative North OBC 
    179179   !!-------------------------------- 
    180    REAL(wp), DIMENSION(jpi,jpk,3,3) ::   unbnd, vnbnd      !: baroclinic u & v components of the velocity over 3 
     180   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   unbnd, vnbnd      !: baroclinic u & v components of the velocity over 3 
    181181      !                                                    !  rows and 3 time step (now, before, and before before) 
    182    REAL(wp), DIMENSION(jpi,jpk,2,2) ::   tnbnd, snbnd      !: north boundary temperature and salinity over 
     182   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   tnbnd, snbnd      !: north boundary temperature and salinity over 
    183183      !                                                    !  2 rows and 2 time step (now and before) 
    184    REAL(wp), DIMENSION(jpi,jpk) ::   u_cynbnd, v_cynbnd    !: Meridional component of the phase speed ratio compu- 
     184   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   u_cynbnd, v_cynbnd    !: Meridional component of the phase speed ratio compu- 
    185185      !                                                    !  ted with radiation of zonal and meridional velocity  
    186186      !                                                    !  (respectively) at the north OB (u_cynbnd = cx rdt ) 
    187    REAL(wp), DIMENSION(jpi,jpk) ::   unmsk, vnmsk, tnmsk   !: 2D mask for the North OB 
     187   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   unmsk, vnmsk, tnmsk   !: 2D mask for the North OB 
    188188 
    189189   ! Note that those arrays are optimized for mpp case  
     
    199199   INTEGER ::   njs0p1, njs1p1     !: do loop index in mpp case for jpsob+1 
    200200 
    201    REAL(wp), DIMENSION(jpi) ::    &   !: 
     201   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:) ::    &   !: 
    202202      sshfos,           & !: now climatology of the south boundary sea surface height 
    203203      ubtfos,vbtfos       !: now climatology of the south boundary barotropic transport 
    204204 
    205    REAL(wp), DIMENSION(jpi,jpk) ::    &   !: 
     205   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::    &   !: 
    206206      ufos, vfos,       & !: now climatology of the south boundary velocities  
    207207      tfos, sfos,       & !: now climatology of the south boundary temperature and salinity 
     
    209209      !                   !  in the obcdyn.F90 routine 
    210210 
    211    REAL(wp), DIMENSION(jpi,jpj) ::   sshfos_b     !: south boundary ssh correction averaged over the barotropic loop 
     211   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshfos_b     !: south boundary ssh correction averaged over the barotropic loop 
    212212      !                                           !  (if Flather's algoritm applied at open boundary) 
    213213 
     
    215215   !! Arrays for radiative South OBC   (computed by the forward time step in dynspg) 
    216216   !!-------------------------------- 
    217    REAL(wp), DIMENSION(jpi,jpk,3,3) ::   usbnd, vsbnd     !: baroclinic u & v components of the velocity over 3  
     217   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   usbnd, vsbnd     !: baroclinic u & v components of the velocity over 3  
    218218      !                                                   !  rows and 3 time step (now, before, and before before) 
    219    REAL(wp), DIMENSION(jpi,jpk,2,2) ::   tsbnd, ssbnd     !: south boundary temperature and salinity over 
     219   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   tsbnd, ssbnd     !: south boundary temperature and salinity over 
    220220      !                                                   !  2 rows and 2 time step (now and before) 
    221    REAL(wp), DIMENSION(jpi,jpk) ::   u_cysbnd, v_cysbnd   !: Meridional component of the phase speed ratio 
     221   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   u_cysbnd, v_cysbnd   !: Meridional component of the phase speed ratio 
    222222      !                                                   !  computed with radiation of zonal and meridional velocity  
    223223      !                                                   !  (repsectively) at the south OB (u_cynbnd = cx rdt ) 
    224    REAL(wp), DIMENSION(jpi,jpk) ::   usmsk, vsmsk, tsmsk  !: 2D mask for the South OB 
     224   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   usmsk, vsmsk, tsmsk  !: 2D mask for the South OB 
    225225 
    226226#else 
     
    235235   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) 
    236236   !!====================================================================== 
     237#if defined key_obc 
     238CONTAINS 
     239 
     240   FUNCTION obc_oce_alloc() 
     241     IMPLICIT none 
     242 
     243     ALLOCATE(                                                               & 
     244              !! East open boundary 
     245              obctmsk(jpi,jpj), obcumask(jpi,jpj), obcvmask(jpi,jpj),        & 
     246              sshfoe(jpjed:jpjef), ubtfoe(jpjed:jpjef), vbtfoe(jpjed:jpjef), & 
     247              ufoe(jpj,jpk), vfoe(jpj,jpk), tfoe(jpj,jpk), sfoe(jpj,jpk),    & 
     248              uclie(jpj,jpk), sshfoe_b(jpjed:jpjef,jpj),                     & 
     249              !! Arrays for radiative East OBC 
     250              uebnd(jpj,jpk,3,3), vebnd(jpj,jpk,3,3) ,                       & 
     251              tebnd(jpj,jpk,2,2), sebnd(jpj,jpk,2,2),                        & 
     252              u_cxebnd(jpj,jpk), v_cxebnd(jpj,jpk),                          & 
     253              uemsk(jpj,jpk), vemsk(jpj,jpk), temsk(jpj,jpk),                & 
     254              !! West open boundary 
     255              sshfow(jpjwd:jpjwf), ubtfow(jpjwd:jpjwf), vbtfow(jpjwd:jpjwf), & 
     256              ufow(jpj,jpk), vfow(jpj,jpk), tfow(jpj,jpk),                   & 
     257              sfow(jpj,jpk), ucliw(jpj,jpk), sshfow_b(jpjwd:jpjwf,jpj),      & 
     258              !! Arrays for radiative West OBC 
     259              uwbnd(jpj,jpk,3,3), vwbnd(jpj,jpk,3,3),                        & 
     260              twbnd(jpj,jpk,2,2), swbnd(jpj,jpk,2,2),                        & 
     261              u_cxwbnd(jpj,jpk), v_cxwbnd(jpj,jpk),                          & 
     262              uwmsk(jpj,jpk), vwmsk(jpj,jpk), twmsk(jpj,jpk),                & 
     263              !! North open boundary 
     264              sshfon(jpind:jpinf), ubtfon(jpind:jpinf), vbtfon(jpind:jpinf), & 
     265              ufon(jpi,jpk), vfon(jpi,jpk), tfon(jpi,jpk),                   & 
     266              sfon(jpi,jpk), vclin(jpi,jpk), sshfon_b(jpind:jpinf,jpj),      & 
     267              !! Arrays for radiative North OBC 
     268              unbnd(jpi,jpk,3,3), vnbnd(jpi,jpk,3,3),                        & 
     269              tnbnd(jpi,jpk,2,2), snbnd(jpi,jpk,2,2),                        & 
     270              u_cynbnd(jpi,jpk), v_cynbnd(jpi,jpk),                          & 
     271              unmsk(jpi,jpk), vnmsk(jpi,jpk), tnmsk (jpi,jpk),               & 
     272              !! South open boundary 
     273              sshfos(jpisd:jpisf), ubtfos(jpisd:jpisf), vbtfos(jpisd:jpisf), & 
     274              ufos(jpi,jpk), vfos(jpi,jpk), tfos(jpi,jpk),                   & 
     275              sfos(jpi,jpk), vclis(jpi,jpk),                                 & 
     276              sshfos_b(jpisd:jpisf,jpj),                                     & 
     277              !! Arrays for radiative South OBC  
     278              usbnd(jpi,jpk,3,3), vsbnd(jpi,jpk,3,3),                        & 
     279              tsbnd(jpi,jpk,2,2), ssbnd(jpi,jpk,2,2),                        & 
     280              u_cysbnd(jpi,jpk), v_cysbnd(jpi,jpk),                          & 
     281              usmsk(jpi,jpk), vsmsk(jpi,jpk), tsmsk(jpi,jpk),                & 
     282              !! 
     283              Stat=obc_oce_alloc ) 
     284 
     285   END FUNCTION obc_oce_alloc 
     286#endif ! Defined key_obc 
     287 
    237288END MODULE obc_oce 
  • branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OPA_SRC/OBC/obcdta.F90

    r2528 r2590  
    5353   ! bt arrays for interpolating time dependent data on the boundaries 
    5454   INTEGER :: nt_m=0, ntobc_m 
    55    REAL(wp), DIMENSION(jpj,0:jptobc) :: ubtedta, vbtedta, sshedta  ! East 
    56    REAL(wp), DIMENSION(jpj,0:jptobc) :: ubtwdta, vbtwdta, sshwdta ! West 
    57    REAL(wp), DIMENSION(jpi,0:jptobc) :: ubtndta, vbtndta, sshndta ! North 
    58    REAL(wp), DIMENSION(jpi,0:jptobc) :: ubtsdta, vbtsdta, sshsdta ! South 
     55   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubtedta, vbtedta, sshedta       ! East 
     56   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubtwdta, vbtwdta, sshwdta   ! West 
     57   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubtndta, vbtndta, sshndta   ! North 
     58   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubtsdta, vbtsdta, sshsdta   ! South 
    5959   ! arrays used for interpolating time dependent data on the boundaries 
    60    REAL(wp), DIMENSION(jpj,jpk,0:jptobc) :: uedta, vedta, tedta, sedta    ! East 
    61    REAL(wp), DIMENSION(jpj,jpk,0:jptobc) :: uwdta, vwdta, twdta, swdta    ! West 
    62    REAL(wp), DIMENSION(jpi,jpk,0:jptobc) :: undta, vndta, tndta, sndta    ! North 
    63    REAL(wp), DIMENSION(jpi,jpk,0:jptobc) :: usdta, vsdta, tsdta, ssdta    ! South 
     60   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: uedta, vedta, tedta, sedta    ! East 
     61   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: uwdta, vwdta, twdta, swdta    ! West 
     62   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: undta, vndta, tndta, sndta    ! North 
     63   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: usdta, vsdta, tsdta, ssdta    ! South 
    6464# else 
    6565   ! bt arrays for interpolating time dependent data on the boundaries 
    66    REAL(wp), DIMENSION(jpj,jptobc) :: ubtedta, vbtedta, sshedta  ! East 
    67    REAL(wp), DIMENSION(jpj,jptobc) :: ubtwdta, vbtwdta, sshwdta        ! West 
    68    REAL(wp), DIMENSION(jpi,jptobc) :: ubtndta, vbtndta, sshndta        ! North 
    69    REAL(wp), DIMENSION(jpi,jptobc) :: ubtsdta, vbtsdta, sshsdta        ! South 
     66   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubtedta, vbtedta, sshedta        ! East 
     67   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubtwdta, vbtwdta, sshwdta        ! West 
     68   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubtndta, vbtndta, sshndta        ! North 
     69   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubtsdta, vbtsdta, sshsdta        ! South 
    7070   ! arrays used for interpolating time dependent data on the boundaries 
    71    REAL(wp), DIMENSION(jpj,jpk,jptobc) :: uedta, vedta, tedta, sedta    ! East 
    72    REAL(wp), DIMENSION(jpj,jpk,jptobc) :: uwdta, vwdta, twdta, swdta    ! West 
    73    REAL(wp), DIMENSION(jpi,jpk,jptobc) :: undta, vndta, tndta, sndta    ! North 
    74    REAL(wp), DIMENSION(jpi,jpk,jptobc) :: usdta, vsdta, tsdta, ssdta    ! South 
     71   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: uedta, vedta, tedta, sedta    ! East 
     72   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: uwdta, vwdta, twdta, swdta    ! West 
     73   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: undta, vndta, tndta, sndta    ! North 
     74   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: usdta, vsdta, tsdta, ssdta    ! South 
    7575# endif 
    76    LOGICAL, DIMENSION (jpj,jpk ) :: ltemsk=.TRUE., luemsk=.TRUE., lvemsk=.TRUE.  ! boolean msks 
    77    LOGICAL, DIMENSION (jpj,jpk ) :: ltwmsk=.TRUE., luwmsk=.TRUE., lvwmsk=.TRUE.  ! used for outliers 
    78    LOGICAL, DIMENSION (jpi,jpk ) :: ltnmsk=.TRUE., lunmsk=.TRUE., lvnmsk=.TRUE.  ! checks 
    79    LOGICAL, DIMENSION (jpi,jpk ) :: ltsmsk=.TRUE., lusmsk=.TRUE., lvsmsk=.TRUE. 
     76   ! Masks set to .TRUE. after successful allocation below 
     77   LOGICAL, ALLOCATABLE, SAVE, DIMENSION (:,: ) :: ltemsk, luemsk, lvemsk  ! boolean msks 
     78   LOGICAL, ALLOCATABLE, SAVE, DIMENSION (:,: ) :: ltwmsk, luwmsk, lvwmsk  ! used for outliers 
     79   LOGICAL, ALLOCATABLE, SAVE, DIMENSION (:,: ) :: ltnmsk, lunmsk, lvnmsk  ! checks 
     80   LOGICAL, ALLOCATABLE, SAVE, DIMENSION (:,: ) :: ltsmsk, lusmsk, lvsmsk 
    8081 
    8182   !! * Substitutions 
     
    8990 
    9091CONTAINS 
     92 
     93   FUNCTION obc_dta_alloc() 
     94      !!------------------------------------------------------------------- 
     95      !!                     ***  ROUTINE obc_dta_alloc  *** 
     96      !!                     
     97      !!------------------------------------------------------------------- 
     98      IMPLICIT none 
     99      INTEGER :: obc_dta_alloc 
     100      INTEGER :: ierr(2) 
     101      !!------------------------------------------------------------------- 
     102 
     103# if defined key_dynspg_ts 
     104      ALLOCATE(ubtedta(jpj,0:jptobc), vbtedta(jpj,0:jptobc), & 
     105               sshedta(jpj,0:jptobc), ubtwdta(jpj,0:jptobc), & 
     106               vbtwdta(jpj,0:jptobc), sshwdta(jpj,0:jptobc), & 
     107               ubtndta(jpi,0:jptobc), vbtndta(jpi,0:jptobc), & 
     108               sshndta(jpi,0:jptobc), ubtsdta(jpi,0:jptobc), & 
     109               vbtsdta(jpi,0:jptobc), sshsdta(jpi,0:jptobc), & 
     110               ! arrays used for interpolating time dependent data on the boundaries 
     111               uedta(jpj,jpk,0:jptobc), vedta(jpj,jpk,0:jptobc), & 
     112               tedta(jpj,jpk,0:jptobc), sedta(jpj,jpk,0:jptobc), & 
     113               uwdta(jpj,jpk,0:jptobc), vwdta(jpj,jpk,0:jptobc), & 
     114               twdta(jpj,jpk,0:jptobc), swdta(jpj,jpk,0:jptobc), & 
     115               undta(jpi,jpk,0:jptobc), vndta(jpi,jpk,0:jptobc), & 
     116               tndta(jpi,jpk,0:jptobc), sndta(jpi,jpk,0:jptobc), & 
     117               usdta(jpi,jpk,0:jptobc), vsdta(jpi,jpk,0:jptobc), & 
     118               tsdta(jpi,jpk,0:jptobc), ssdta(jpi,jpk,0:jptobc), Stat=ierr(1) ) 
     119# else 
     120               ! bt arrays for interpolating time dependent data on the boundaries 
     121      ALLOCATE(ubtedta(jpj,jptobc), vbtedta(jpj,jptobc), sshedta(jpj,jptobc), & 
     122               ubtwdta(jpj,jptobc), vbtwdta(jpj,jptobc), sshwdta(jpj,jptobc), & 
     123               ubtndta(jpi,jptobc), vbtndta(jpi,jptobc), sshndta(jpi,jptobc), & 
     124               ubtsdta(jpi,jptobc), vbtsdta(jpi,jptobc), sshsdta(jpi,jptobc), & 
     125               ! arrays used for interpolating time dependent data on the boundaries 
     126               uedta(jpj,jpk,jptobc), vedta(jpj,jpk,jptobc),  & 
     127               tedta(jpj,jpk,jptobc), sedta(jpj,jpk,jptobc),  & 
     128               uwdta(jpj,jpk,jptobc), vwdta(jpj,jpk,jptobc),  & 
     129               twdta(jpj,jpk,jptobc), swdta(jpj,jpk,jptobc),  & 
     130               undta(jpi,jpk,jptobc), vndta(jpi,jpk,jptobc),  & 
     131               tndta(jpi,jpk,jptobc), sndta(jpi,jpk,jptobc),  & 
     132               usdta(jpi,jpk,jptobc), vsdta(jpi,jpk,jptobc),  & 
     133               tsdta(jpi,jpk,jptobc), ssdta(jpi,jpk,jptobc), Stat=ierr(1) ) 
     134# endif 
     135 
     136      ALLOCATE(uedta(jpj,jpk,jptobc), vedta(jpj,jpk,jptobc), & 
     137               tedta(jpj,jpk,jptobc), sedta(jpj,jpk,jptobc), & 
     138               uwdta(jpj,jpk,jptobc), vwdta(jpj,jpk,jptobc), & 
     139               twdta(jpj,jpk,jptobc), swdta(jpj,jpk,jptobc), & 
     140               undta(jpj,jpk,jptobc), vndta(jpj,jpk,jptobc), & 
     141               tndta(jpj,jpk,jptobc), sndta(jpj,jpk,jptobc), & 
     142               usdta(jpj,jpk,jptobc), vsdta(jpj,jpk,jptobc), & 
     143               tsdta(jpj,jpk,jptobc), ssdta(jpj,jpk,jptobc), & 
     144               ltemsk(jpj,jpk), luemsk(jpj,jpk), lvemsk(jpj,jpk), & 
     145               ltwmsk(jpj,jpk), luwmsk(jpj,jpk), lvwmsk(jpj,jpk), & 
     146               ltnmsk(jpj,jpk), lunmsk(jpj,jpk), lvnmsk(jpj,jpk), & 
     147               ltsmsk(jpj,jpk), lusmsk(jpj,jpk), lvsmsk(jpj,jpk), & 
     148               Stat=ierr(2)) 
     149 
     150      obc_dta_alloc = MAXVAL(ierr) 
     151 
     152      IF(obc_dta_alloc == 0)THEN 
     153         ! Initialise mask values following successful allocation 
     154         ltemsk(:)=.TRUE. 
     155         luemsk(:)=.TRUE. 
     156         lvemsk(:)=.TRUE. 
     157         ltwmsk(:)=.TRUE. 
     158         luwmsk(:)=.TRUE. 
     159         lvwmsk(:)=.TRUE. 
     160         ltnmsk(:)=.TRUE. 
     161         lunmsk(:)=.TRUE. 
     162         lvnmsk(:)=.TRUE. 
     163         ltsmsk(:)=.TRUE. 
     164         lusmsk(:)=.TRUE. 
     165         lvsmsk(:)=.TRUE. 
     166      END IF 
     167 
     168   END FUNCTION obc_dta_alloc 
     169 
    91170 
    92171   SUBROUTINE obc_dta( kt ) 
Note: See TracChangeset for help on using the changeset viewer.