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

source: branches/2015/dev_r5151_UKMO_ISF/NEMOGCM/NEMO/OPA_SRC/oce.F90 @ 5944

Last change on this file since 5944 was 5189, checked in by mathiot, 9 years ago

ISF cleaning branch: simplification and bug correction in hpg_isf, zps_hde_isf, mixed layer definition, slope, diffusion, vertical advection and top friction

  • Property svn:keywords set to Id
File size: 7.7 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   !!----------------------------------------------------------------------
10   USE par_oce        ! ocean parameters
11   USE lib_mpp        ! MPP library
12
13   IMPLICIT NONE
14   PRIVATE
15
16   PUBLIC oce_alloc ! routine called by nemo_init in nemogcm.F90
17
18   LOGICAL, PUBLIC ::   l_traldf_rot = .FALSE.  !: rotated laplacian operator for lateral diffusion
19
20   !! dynamics and tracer fields                            ! before ! now    ! after  ! the after trends becomes the fields
21   !! --------------------------                            ! fields ! fields ! trends ! only after tra_zdf and dyn_spg
22   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   ub   ,  un    , ua     !: i-horizontal velocity          [m/s]
23   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   vb   ,  vn    , va     !: j-horizontal velocity          [m/s]
24   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   ua_sv,  va_sv          !: Saved trends (time spliting)   [m/s2]
25   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::           wn             !: vertical velocity              [m/s]
26   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   rotb ,  rotn           !: relative vorticity             [s-1]
27   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   hdivb,  hdivn          !: horizontal divergence          [s-1]
28   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   tsb  ,  tsn   , tsa    !: 4D T-S fields                  [Celcius,psu]
29   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   rab_b,  rab_n          !: thermal/haline expansion coef. [Celcius-1,psu-1]
30   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   rn2b ,  rn2            !: brunt-vaisala frequency**2     [s-2]
31   !
32   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   rhd    !: in situ density anomalie rhd=(rho-rau0)/rau0  [no units]
33   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   rhop   !: potential volumic mass                           [kg/m3]
34
35   !! free surface                                      !  before  ! now    ! after  !
36   !! ------------                                      !  fields  ! fields ! fields !
37   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ub_b   ,  un_b  ,  ua_b  !: Barotropic velocities at u-point [m/s]
38   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   vb_b   ,  vn_b  ,  va_b  !: Barotropic velocities at v-point [m/s]
39   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshb   ,  sshn  ,  ssha  !: sea surface height at t-point [m]
40   !
41   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   spgu, spgv               !: horizontal surface pressure gradient
42
43   !! interpolated gradient (only used in zps case)
44   !! ---------------------
45   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   gtsu, gtsv   !: horizontal gradient of T, S bottom u-point
46   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   gru , grv    !: horizontal gradient of rd at bottom u-point
47
48   !! (ISF) interpolated gradient (only used for ice shelf case)
49   !! ---------------------
50   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   gtui, gtvi   !: horizontal gradient of T, S and rd at top u-point
51   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   grui, grvi   !: horizontal gradient of T, S and rd at top v-point 
52   !! (ISF) ice load
53   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   riceload
54
55   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   rke          !: kinetic energy
56
57   !! arrays relating to embedding ice in the ocean. These arrays need to be declared
58   !! even if no ice model is required. In the no ice model or traditional levitating
59   !! ice cases they contain only zeros
60   !! ---------------------
61   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   snwice_mass        !: mass of snow and ice at current  ice time step   [Kg/m2]
62   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   snwice_mass_b      !: mass of snow and ice at previous ice time step   [Kg/m2]
63   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   snwice_fmass       !: time evolution of mass of snow+ice               [Kg/m2/s]
64
65   !! Energy budget of the leads (open water embedded in sea ice)
66   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   fraqsr_1lev        !: fraction of solar net radiation absorbed in the first ocean level [-]
67
68   !!----------------------------------------------------------------------
69   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
70   !! $Id$
71   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
72   !!----------------------------------------------------------------------
73CONTAINS
74
75   INTEGER FUNCTION oce_alloc()
76      !!----------------------------------------------------------------------
77      !!                   ***  FUNCTION oce_alloc  ***
78      !!----------------------------------------------------------------------
79      INTEGER :: ierr(4)
80      !!----------------------------------------------------------------------
81      !
82      ALLOCATE( ub   (jpi,jpj,jpk)      , un   (jpi,jpj,jpk)      , ua(jpi,jpj,jpk)       ,     &
83         &      vb   (jpi,jpj,jpk)      , vn   (jpi,jpj,jpk)      , va(jpi,jpj,jpk)       ,     &     
84         &      ua_sv(jpi,jpj,jpk)      , va_sv(jpi,jpj,jpk)      ,                             &     
85         &      wn   (jpi,jpj,jpk)      ,                                                       &
86         &      rotb (jpi,jpj,jpk)      , rotn (jpi,jpj,jpk)      ,                             &   
87         &      hdivb(jpi,jpj,jpk)      , hdivn(jpi,jpj,jpk)      ,                             &
88         &      tsb  (jpi,jpj,jpk,jpts) , tsn  (jpi,jpj,jpk,jpts) , tsa(jpi,jpj,jpk,jpts) ,     &
89         &      rab_b(jpi,jpj,jpk,jpts) , rab_n(jpi,jpj,jpk,jpts) ,                             &
90         &      rn2b (jpi,jpj,jpk)      , rn2  (jpi,jpj,jpk)                              , STAT=ierr(1) )
91         !
92      ALLOCATE(rhd (jpi,jpj,jpk) ,                                         &
93         &     rhop(jpi,jpj,jpk) ,                                         &
94         &     rke(jpi,jpj,jpk)  ,                                         &
95         &     sshb(jpi,jpj)     , sshn(jpi,jpj)   , ssha(jpi,jpj)   ,     &
96         &     ub_b(jpi,jpj)     , un_b(jpi,jpj)   , ua_b(jpi,jpj)   ,     &
97         &     vb_b(jpi,jpj)     , vn_b(jpi,jpj)   , va_b(jpi,jpj)   ,     &
98         &     spgu  (jpi,jpj)   , spgv(jpi,jpj)   ,                       &
99         &     gtsu(jpi,jpj,jpts), gtsv(jpi,jpj,jpts),                     &
100         &     gru(jpi,jpj)      , grv(jpi,jpj)      ,                     &
101         &     gtui(jpi,jpj,jpts), gtvi(jpi,jpj,jpts),                     &
102         &     grui(jpi,jpj)     , grvi(jpi,jpj)     ,                     &
103         &     riceload(jpi,jpj),                             STAT=ierr(2) )
104         !
105      ALLOCATE( snwice_mass(jpi,jpj) , snwice_mass_b(jpi,jpj), snwice_fmass(jpi,jpj) , STAT=ierr(3) )
106         !
107      ALLOCATE( fraqsr_1lev(jpi,jpj) , STAT=ierr(4) )
108         !
109      oce_alloc = MAXVAL( ierr )
110      IF( oce_alloc /= 0 )   CALL ctl_warn('oce_alloc: failed to allocate arrays')
111      !
112   END FUNCTION oce_alloc
113
114   !!======================================================================
115END MODULE oce
Note: See TracBrowser for help on using the repository browser.