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.
oce.F90 in branches/2017/dev_r8624_ENHANCE4_FREESURFACE/NEMOGCM/NEMO/OPA_SRC – NEMO

source: branches/2017/dev_r8624_ENHANCE4_FREESURFACE/NEMOGCM/NEMO/OPA_SRC/oce.F90 @ 8805

Last change on this file since 8805 was 8805, checked in by jchanut, 6 years ago

Enhance4-freesurface. step 1: recover tracer conservation with split explicit free surface - #1959

  • Property svn:keywords set to Id
File size: 9.0 KB
RevLine 
[3]1MODULE oce
2   !!======================================================================
[15]3   !!                      ***  MODULE  oce  ***
[3]4   !! Ocean        :  dynamics and active tracers defined in memory
5   !!======================================================================
[2528]6   !! History :  1.0  !  2002-11  (G. Madec)  F90: Free form and module
[1438]7   !!            3.1  !  2009-02  (G. Madec, M. Leclair)  pure z* coordinate
[2528]8   !!            3.3  !  2010-09  (C. Ethe) TRA-TRC merge: add ts, gtsu, gtsv 4D arrays
[5836]9   !!            3.7  !  2014-01  (G. Madec) suppression of curl and before hdiv from in-core memory
[3]10   !!----------------------------------------------------------------------
[2715]11   USE par_oce        ! ocean parameters
12   USE lib_mpp        ! MPP library
[3]13
14   IMPLICIT NONE
[15]15   PRIVATE
[3]16
[2715]17   PUBLIC oce_alloc ! routine called by nemo_init in nemogcm.F90
[3]18
[2715]19   !! dynamics and tracer fields                            ! before ! now    ! after  ! the after trends becomes the fields
20   !! --------------------------                            ! fields ! fields ! trends ! only after tra_zdf and dyn_spg
[5836]21   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   ub   ,  un    , ua     !: i-horizontal velocity        [m/s]
22   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   vb   ,  vn    , va     !: j-horizontal velocity        [m/s]
23   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::           wn             !: vertical velocity            [m/s]
24   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::           hdivn          !: horizontal divergence        [s-1]
[7646]25   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   tsb  ,  tsn   , tsa    !: 4D T-S fields                  [Celsius,psu]
26   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   rab_b,  rab_n          !: thermal/haline expansion coef. [Celsius-1,psu-1]
[4990]27   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   rn2b ,  rn2            !: brunt-vaisala frequency**2     [s-2]
[1438]28   !
[2715]29   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   rhd    !: in situ density anomalie rhd=(rho-rau0)/rau0  [no units]
30   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   rhop   !: potential volumic mass                           [kg/m3]
[3]31
[2715]32   !! free surface                                      !  before  ! now    ! after  !
[4354]33   !! ------------                                      !  fields  ! fields ! fields !
34   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ub_b   ,  un_b  ,  ua_b  !: Barotropic velocities at u-point [m/s]
35   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   vb_b   ,  vn_b  ,  va_b  !: Barotropic velocities at v-point [m/s]
36   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshb   ,  sshn  ,  ssha  !: sea surface height at t-point [m]
[359]37
[6140]38   !! Arrays at barotropic time step:                   ! befbefore! before !  now   ! after  !
39   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ubb_e  ,  ub_e  ,  un_e  , ua_e   !: u-external velocity
40   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   vbb_e  ,  vb_e  ,  vn_e  , va_e   !: v-external velocity
41   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshbb_e,  sshb_e,  sshn_e, ssha_e !: external ssh
42   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::                              hu_e   !: external u-depth
43   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::                              hv_e   !: external v-depth
44   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::                              hur_e  !: inverse of u-depth
45   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::                              hvr_e  !: inverse of v-depth
46   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ub2_b  , vb2_b           !: Half step fluxes (ln_bt_fw=T)
[8805]47   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   un_bf  , vn_bf           !: Asselin filter corrective fluxes (ln_bt_fw=T)
[5930]48#if defined key_agrif
[6140]49   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ub2_i_b, vb2_i_b         !: Half step time integrated fluxes
[5930]50#endif
[6140]51   !
52   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   spgu, spgv               !: horizontal surface pressure gradient
[5930]53
[2528]54   !! interpolated gradient (only used in zps case)
55   !! ---------------------
[2715]56   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   gtsu, gtsv   !: horizontal gradient of T, S bottom u-point
57   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   gru , grv    !: horizontal gradient of rd at bottom u-point
[2528]58
[4990]59   !! (ISF) interpolated gradient (only used for ice shelf case)
60   !! ---------------------
61   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   gtui, gtvi   !: horizontal gradient of T, S and rd at top u-point
62   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   grui, grvi   !: horizontal gradient of T, S and rd at top v-point 
63   !! (ISF) ice load
64   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   riceload
65
[4152]66   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   rke          !: kinetic energy
67
[3625]68   !! arrays relating to embedding ice in the ocean. These arrays need to be declared
69   !! even if no ice model is required. In the no ice model or traditional levitating
70   !! ice cases they contain only zeros
71   !! ---------------------
72   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   snwice_mass        !: mass of snow and ice at current  ice time step   [Kg/m2]
73   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   snwice_mass_b      !: mass of snow and ice at previous ice time step   [Kg/m2]
74   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   snwice_fmass       !: time evolution of mass of snow+ice               [Kg/m2/s]
75
[4990]76   !! Energy budget of the leads (open water embedded in sea ice)
77   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   fraqsr_1lev        !: fraction of solar net radiation absorbed in the first ocean level [-]
[4205]78
[3]79   !!----------------------------------------------------------------------
[5836]80   !! NEMO/OPA 3.7 , NEMO Consortium (2014)
[1438]81   !! $Id$
[2715]82   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
83   !!----------------------------------------------------------------------
84CONTAINS
85
86   INTEGER FUNCTION oce_alloc()
87      !!----------------------------------------------------------------------
88      !!                   ***  FUNCTION oce_alloc  ***
89      !!----------------------------------------------------------------------
[6140]90      INTEGER :: ierr(7)
[2715]91      !!----------------------------------------------------------------------
92      !
[6140]93      ierr(:) = 0 
[2715]94      ALLOCATE( ub   (jpi,jpj,jpk)      , un   (jpi,jpj,jpk)      , ua(jpi,jpj,jpk)       ,     &
[5930]95         &      vb   (jpi,jpj,jpk)      , vn   (jpi,jpj,jpk)      , va(jpi,jpj,jpk)       ,     &         
[5836]96         &      wn   (jpi,jpj,jpk)      , hdivn(jpi,jpj,jpk)      ,                             &
[2715]97         &      tsb  (jpi,jpj,jpk,jpts) , tsn  (jpi,jpj,jpk,jpts) , tsa(jpi,jpj,jpk,jpts) ,     &
[4990]98         &      rab_b(jpi,jpj,jpk,jpts) , rab_n(jpi,jpj,jpk,jpts) ,                             &
[5836]99         &      rn2b (jpi,jpj,jpk)      , rn2  (jpi,jpj,jpk)      ,                             &
100         &      rhd  (jpi,jpj,jpk)      , rhop (jpi,jpj,jpk)                              , STAT=ierr(1) )
[2715]101         !
[5836]102      ALLOCATE(rke(jpi,jpj,jpk)  ,                                         &
[4354]103         &     sshb(jpi,jpj)     , sshn(jpi,jpj)   , ssha(jpi,jpj)   ,     &
104         &     ub_b(jpi,jpj)     , un_b(jpi,jpj)   , ua_b(jpi,jpj)   ,     &
105         &     vb_b(jpi,jpj)     , vn_b(jpi,jpj)   , va_b(jpi,jpj)   ,     &
[4292]106         &     spgu  (jpi,jpj)   , spgv(jpi,jpj)   ,                       &
107         &     gtsu(jpi,jpj,jpts), gtsv(jpi,jpj,jpts),                     &
[4990]108         &     gru(jpi,jpj)      , grv(jpi,jpj)      ,                     &
109         &     gtui(jpi,jpj,jpts), gtvi(jpi,jpj,jpts),                     &
110         &     grui(jpi,jpj)     , grvi(jpi,jpj)     ,                     &
111         &     riceload(jpi,jpj),                             STAT=ierr(2) )
[2715]112         !
[4205]113      ALLOCATE( snwice_mass(jpi,jpj) , snwice_mass_b(jpi,jpj), snwice_fmass(jpi,jpj) , STAT=ierr(3) )
[3625]114         !
[4990]115      ALLOCATE( fraqsr_1lev(jpi,jpj) , STAT=ierr(4) )
[4205]116         !
[6140]117      ALLOCATE( ssha_e(jpi,jpj),  sshn_e(jpi,jpj), sshb_e(jpi,jpj), sshbb_e(jpi,jpj), &
118         &        ua_e(jpi,jpj),    un_e(jpi,jpj),   ub_e(jpi,jpj),   ubb_e(jpi,jpj), &
119         &        va_e(jpi,jpj),    vn_e(jpi,jpj),   vb_e(jpi,jpj),   vbb_e(jpi,jpj), &
120         &        hu_e(jpi,jpj),   hur_e(jpi,jpj),   hv_e(jpi,jpj),   hvr_e(jpi,jpj), STAT=ierr(5) )
121         !
[8805]122      ALLOCATE( ub2_b(jpi,jpj), vb2_b(jpi,jpj), un_bf(jpi,jpj), vn_bf(jpi,jpj)      , STAT=ierr(6) )
[6140]123#if defined key_agrif
124      ALLOCATE( ub2_i_b(jpi,jpj), vb2_i_b(jpi,jpj)                                  , STAT=ierr(7) )
125#endif
126         !
[2715]127      oce_alloc = MAXVAL( ierr )
128      IF( oce_alloc /= 0 )   CALL ctl_warn('oce_alloc: failed to allocate arrays')
129      !
130   END FUNCTION oce_alloc
131
[1438]132   !!======================================================================
[3]133END MODULE oce
Note: See TracBrowser for help on using the repository browser.