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/OBC/obcdyn_bt.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/OBC/obcdyn_bt.F90

    r2528 r2715  
    11MODULE obcdyn_bt 
     2   !!====================================================================== 
     3   !!                       ***  MODULE  obcdyn_bt  *** 
     4   !! Ocean dynamics:   Radiation/prescription of sea surface heights on each open boundary 
     5   !!====================================================================== 
     6   !! History :  1.0  ! 2005-12  (V. Garnier) original code 
     7   !!---------------------------------------------------------------------- 
    28#if ( defined key_dynspg_ts || defined key_dynspg_exp ) && defined key_obc 
    3    !!================================================================================= 
    4    !!                       ***  MODULE  obcdyn_bt  *** 
    5    !! Ocean dynamics:   Radiation/prescription of sea surface heights 
    6    !!                   on each open boundary 
    7    !!================================================================================= 
    8  
    9    !!--------------------------------------------------------------------------------- 
     9   !!---------------------------------------------------------------------- 
     10   !!   'key_dynspg_ts'     OR                   time spliting free surface 
     11   !!   'key_dynspg_exp'    AND                       explicit free surface 
     12   !!   'key_obc'                                   Open Boundary Condition 
     13   !!---------------------------------------------------------------------- 
    1014   !!   obc_dyn_bt        : call the subroutine for each open boundary 
    1115   !!   obc_dyn_bt_east   : Flather's algorithm at the east open boundary 
     
    1317   !!   obc_dyn_bt_north  : Flather's algorithm at the north open boundary 
    1418   !!   obc_dyn_bt_south  : Flather's algorithm at the south open boundary 
    15    !!---------------------------------------------------------------------------------- 
    16  
    17    !!---------------------------------------------------------------------------------- 
    18    !! * Modules used 
     19   !!---------------------------------------------------------------------- 
    1920   USE oce             ! ocean dynamics and tracers  
    2021   USE dom_oce         ! ocean space and time domain 
     
    3031   PRIVATE 
    3132 
    32    !! * Accessibility 
    33    PUBLIC obc_dyn_bt  ! routine called in dynnxt (explicit free surface case) 
    34  
    35    !!--------------------------------------------------------------------------------- 
     33   PUBLIC   obc_dyn_bt  ! routine called in dynnxt (explicit free surface case) 
     34 
     35   !!---------------------------------------------------------------------- 
    3636   !! NEMO/OPA 3.3 , NEMO Consortium (2010) 
    3737   !! $Id$  
    38    !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) 
    39    !!---------------------------------------------------------------------- 
    40  
     38   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt) 
     39   !!---------------------------------------------------------------------- 
    4140CONTAINS 
    4241 
    43    SUBROUTINE obc_dyn_bt ( kt ) 
    44       !!------------------------------------------------------------------------------ 
    45       !!                      SUBROUTINE obc_dyn_bt 
    46       !!                     *********************** 
    47       !! ** Purpose : 
    48       !!      Apply Flather's algorithm at open boundaries for the explicit 
    49       !!      free surface case and free surface case with time-splitting 
     42   SUBROUTINE obc_dyn_bt( kt ) 
     43      !!---------------------------------------------------------------------- 
     44      !!                 ***  SUBROUTINE obc_dyn_bt  *** 
     45      !! 
     46      !! ** Purpose :   Apply Flather's algorithm at open boundaries for the explicit 
     47      !!              free surface case and free surface case with time-splitting 
    5048      !! 
    5149      !!      This routine is called in dynnxt.F routine and updates ua, va and sshn.  
     
    5553      !!      open one (must be done in the param_obc.h90 file). 
    5654      !! 
    57       !! ** Reference :  
    58       !!         Flather, R. A., 1976, Mem. Soc. R. Sci. Liege, Ser. 6, 10, 141-164 
    59       !! 
    60       !! History : 
    61       !!   9.0  !  05-12  (V. Garnier) original  
    62       !!---------------------------------------------------------------------- 
    63       !! * Arguments 
    64       INTEGER, INTENT( in ) ::   kt 
    65  
     55      !! Reference :   Flather, R. A., 1976, Mem. Soc. R. Sci. Liege, Ser. 6, 10, 141-164 
     56      !!---------------------------------------------------------------------- 
     57      INTEGER, INTENT(in) ::   kt 
    6658      !!---------------------------------------------------------------------- 
    6759 
     
    8577 
    8678# if defined key_dynspg_exp 
     79 
    8780   SUBROUTINE obc_dyn_bt_east  
    88       !!------------------------------------------------------------------------------ 
     81      !!---------------------------------------------------------------------- 
    8982      !!                  ***  SUBROUTINE obc_dyn_bt_east  *** 
    9083      !!               
     
    9386      !!      Fix sea surface height (sshn) on east open boundary 
    9487      !!      The logical lfbceast must be .TRUE. 
    95       !! 
    96       !!  History : 
    97       !!   9.0  !  05-12  (V. Garnier) original 
    98       !!------------------------------------------------------------------------------ 
    99       !! * Local declaration 
    100       INTEGER ::   ji, jj, jk ! dummy loop indices 
    101       !!------------------------------------------------------------------------------ 
     88      !!---------------------------------------------------------------------- 
     89      INTEGER, INTENT(in) ::   kt 
     90      !!---------------------------------------------------------------------- 
     91      INTEGER ::   ji, jj, jk   ! dummy loop indices 
     92      !!---------------------------------------------------------------------- 
    10293 
    10394      DO ji = nie0, nie1 
     
    120111 
    121112   SUBROUTINE obc_dyn_bt_west  
    122       !!------------------------------------------------------------------------------ 
     113      !!---------------------------------------------------------------------- 
    123114      !!                  ***  SUBROUTINE obc_dyn_bt_west  *** 
    124115      !!                   
     
    127118      !!      Fix sea surface height (sshn) on west open boundary 
    128119      !!      The logical lfbcwest must be .TRUE. 
    129       !! 
    130       !!  History : 
    131       !!   9.0  !  05-12  (V. Garnier) original 
    132       !!------------------------------------------------------------------------------ 
    133       !! * Local declaration 
    134       INTEGER ::   ji, jj, jk ! dummy loop indices 
    135       !!------------------------------------------------------------------------------ 
    136  
     120      !!---------------------------------------------------------------------- 
     121      INTEGER ::   ji, jj, jk   ! dummy loop indices 
     122      !!---------------------------------------------------------------------- 
     123      ! 
    137124      DO ji = niw0, niw1 
    138125         DO jk = 1, jpkm1 
     
    147134         END DO 
    148135      END DO 
    149  
     136      ! 
    150137   END SUBROUTINE obc_dyn_bt_west 
     138 
    151139 
    152140   SUBROUTINE obc_dyn_bt_north  
     
    158146      !!      Fix sea surface height (sshn) on north open boundary 
    159147      !!      The logical lfbcnorth must be .TRUE. 
    160       !! 
    161       !!  History : 
    162       !!   9.0  !  05-12  (V. Garnier) original 
    163       !!------------------------------------------------------------------------------ 
    164       !! * Local declaration 
    165       INTEGER ::   ji, jj, jk ! dummy loop indices 
    166       !!------------------------------------------------------------------------------ 
    167  
     148      !!---------------------------------------------------------------------- 
     149      INTEGER ::   ji, jj, jk   ! dummy loop indices 
     150      !!---------------------------------------------------------------------- 
     151      ! 
    168152      DO jj = njn0, njn1 
    169153         DO jk = 1, jpkm1 
     
    180164         END DO 
    181165      END DO 
    182  
     166      ! 
    183167   END SUBROUTINE obc_dyn_bt_north 
    184168 
     169 
    185170   SUBROUTINE obc_dyn_bt_south  
    186       !!------------------------------------------------------------------------------ 
     171      !!---------------------------------------------------------------------- 
    187172      !!                ***  SUBROUTINE obc_dyn_bt_south  *** 
    188173      !!                     
     
    191176      !!      Fix sea surface height (sshn) on south open boundary 
    192177      !!      The logical lfbcsouth must be .TRUE. 
    193       !! 
    194       !!  History : 
    195       !!   9.0  !  05-12  (V. Garnier) original 
    196       !!------------------------------------------------------------------------------ 
    197       !! * Local declaration 
    198       INTEGER ::   ji, jj, jk ! dummy loop indices 
    199  
    200       !!------------------------------------------------------------------------------ 
    201  
     178      !!---------------------------------------------------------------------- 
     179      INTEGER ::   ji, jj, jk   ! dummy loop indices 
     180      !!---------------------------------------------------------------------- 
     181      ! 
    202182      DO jj = njs0, njs1 
    203183         DO jk = 1, jpkm1 
     
    212192         END DO 
    213193      END DO 
    214  
     194      ! 
    215195   END SUBROUTINE obc_dyn_bt_south 
    216196 
     
    225205      !!      Fix sea surface height (sshn) on east open boundary 
    226206      !!      The logical lfbceast must be .TRUE. 
    227       !! 
    228       !!  History : 
    229       !!   9.0  !  05-12  (V. Garnier) original 
    230       !!------------------------------------------------------------------------------ 
    231       !! * Local declaration 
    232       INTEGER ::   ji, jj, jk ! dummy loop indices 
    233       !!------------------------------------------------------------------------------ 
    234  
     207      !!---------------------------------------------------------------------- 
     208      INTEGER ::   ji, jj, jk   ! dummy loop indices 
     209      !!---------------------------------------------------------------------- 
     210      ! 
    235211      DO ji = nie0, nie1 
    236212         DO jk = 1, jpkm1 
     
    245221         END DO 
    246222      END DO 
    247  
     223      ! 
    248224   END SUBROUTINE obc_dyn_bt_east 
    249225 
     226 
    250227   SUBROUTINE obc_dyn_bt_west  
    251       !!------------------------------------------------------------------------------ 
     228      !!--------------------------------------------------------------------- 
    252229      !!                  ***  SUBROUTINE obc_dyn_bt_west  *** 
    253230      !! 
    254       !! ** Purpose : 
    255       !! ** Purpose : 
    256       !!      Apply Flather algorithm on west OBC velocities ua, va 
     231      !! ** Purpose :   Apply Flather algorithm on west OBC velocities ua, va 
    257232      !!      Fix sea surface height (sshn) on west open boundary 
    258233      !!      The logical lfbcwest must be .TRUE. 
    259       !! 
    260       !!  History : 
    261       !!   9.0  !  05-12  (V. Garnier) original 
    262       !!------------------------------------------------------------------------------ 
    263       !! * Local declaration 
    264       INTEGER ::   ji, jj, jk ! dummy loop indices 
    265       !!------------------------------------------------------------------------------ 
    266  
     234      !!---------------------------------------------------------------------- 
     235      INTEGER ::   ji, jj, jk   ! dummy loop indices 
     236      !!---------------------------------------------------------------------- 
     237      ! 
    267238      DO ji = niw0, niw1 
    268239         DO jk = 1, jpkm1 
     
    275246         END DO 
    276247      END DO 
    277  
     248      ! 
    278249   END SUBROUTINE obc_dyn_bt_west 
     250 
    279251 
    280252   SUBROUTINE obc_dyn_bt_north  
    281253      !!------------------------------------------------------------------------------ 
    282       !!                     SUBROUTINE obc_dyn_bt_north 
    283       !!                    ************************* 
     254      !!                ***  SUBROUTINE obc_dyn_bt_north  *** 
     255      !!                 
    284256      !! ** Purpose : 
    285257      !!      Apply Flather algorithm on north OBC velocities ua, va 
    286258      !!      Fix sea surface height (sshn) on north open boundary 
    287259      !!      The logical lfbcnorth must be .TRUE. 
    288       !! 
    289       !!  History : 
    290       !!   9.0  !  05-12  (V. Garnier) original 
    291       !!------------------------------------------------------------------------------ 
    292       !! * Local declaration 
    293       INTEGER ::   ji, jj, jk ! dummy loop indices 
    294       !!------------------------------------------------------------------------------ 
    295  
     260      !!---------------------------------------------------------------------- 
     261      INTEGER ::   ji, jj, jk   ! dummy loop indices 
     262      !!---------------------------------------------------------------------- 
     263      ! 
    296264      DO jj = njn0, njn1 
    297265         DO jk = 1, jpkm1 
     
    306274         END DO 
    307275      END DO 
    308  
     276      ! 
    309277   END SUBROUTINE obc_dyn_bt_north 
     278 
    310279 
    311280   SUBROUTINE obc_dyn_bt_south  
    312281      !!------------------------------------------------------------------------------ 
    313       !!                     SUBROUTINE obc_dyn_bt_south 
    314       !!                    ************************* 
     282      !!                ***  SUBROUTINE obc_dyn_bt_south  *** 
     283      !!                   
    315284      !! ** Purpose : 
    316285      !!      Apply Flather algorithm on south OBC velocities ua, va 
    317286      !!      Fix sea surface height (sshn) on south open boundary 
    318287      !!      The logical lfbcsouth must be .TRUE. 
    319       !! 
    320       !!  History : 
    321       !!   9.0  !  05-12  (V. Garnier) original 
    322       !!------------------------------------------------------------------------------ 
    323       !! * Local declaration 
    324       INTEGER ::   ji, jj, jk ! dummy loop indices 
    325  
    326       !!------------------------------------------------------------------------------ 
    327  
     288      !!---------------------------------------------------------------------- 
     289      INTEGER ::   ji, jj, jk   ! dummy loop indices 
     290      !!---------------------------------------------------------------------- 
     291      ! 
    328292      DO jj = njs0, njs1 
    329293         DO jk = 1, jpkm1 
     
    336300         END DO 
    337301      END DO 
    338  
     302      ! 
    339303   END SUBROUTINE obc_dyn_bt_south 
    340304 
    341305# endif 
     306 
    342307#else 
    343    !!================================================================================= 
    344    !!                       ***  MODULE  obcdyn_bt  *** 
    345    !! Ocean dynamics:   Radiation of velocities on each open boundary 
    346    !!================================================================================= 
     308   !!---------------------------------------------------------------------- 
     309   !!   Default option       No Open Boundaries or not explicit fre surface 
     310   !!---------------------------------------------------------------------- 
    347311CONTAINS 
    348  
    349    SUBROUTINE obc_dyn_bt 
    350                               ! No open boundaries ==> empty routine 
     312   SUBROUTINE obc_dyn_bt      ! Dummy routine 
    351313   END SUBROUTINE obc_dyn_bt 
    352314#endif 
    353315 
     316   !!====================================================================== 
    354317END MODULE obcdyn_bt 
Note: See TracChangeset for help on using the changeset viewer.