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/2015/dev_MetOffice_merge_2015/NEMOGCM/NEMO/OPA_SRC – NEMO

source: branches/2015/dev_MetOffice_merge_2015/NEMOGCM/NEMO/OPA_SRC/oce.F90 @ 6012

Last change on this file since 6012 was 6012, checked in by mathiot, 8 years ago

merge MetO branch with dev_r5151_UKMO_ISF

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