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

source: branches/UKMO/dev_r5518_amm15_test/NEMOGCM/NEMO/OPA_SRC/oce.F90 @ 5886

Last change on this file since 5886 was 5886, checked in by jgraham, 8 years ago

Removed keywords

File size: 8.8 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   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   aru , arv   
48   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   gzu , gzv   
49   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   ge3ru, ge3rv   !: horizontal gradient of T, S and rd at top v-point 
50
51   !! (ISF) interpolated gradient (only used for ice shelf case)
52   !! ---------------------
53   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   gtui, gtvi   !: horizontal gradient of T, S and rd at top u-point
54   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   grui, grvi   !: horizontal gradient of T, S and rd at top v-point 
55   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   arui, arvi   !: horizontal average  of rd          at top v-point 
56   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   gzui, gzvi   !: horizontal gradient of z           at top v-point 
57   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   ge3rui, ge3rvi   !: horizontal gradient of T, S and rd at top v-point 
58   !! (ISF) ice load
59   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   riceload
60
61   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   rke          !: kinetic energy
62
63   !! arrays relating to embedding ice in the ocean. These arrays need to be declared
64   !! even if no ice model is required. In the no ice model or traditional levitating
65   !! ice cases they contain only zeros
66   !! ---------------------
67   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   snwice_mass        !: mass of snow and ice at current  ice time step   [Kg/m2]
68   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   snwice_mass_b      !: mass of snow and ice at previous ice time step   [Kg/m2]
69   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   snwice_fmass       !: time evolution of mass of snow+ice               [Kg/m2/s]
70
71   !! Energy budget of the leads (open water embedded in sea ice)
72   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   fraqsr_1lev        !: fraction of solar net radiation absorbed in the first ocean level [-]
73
74   !!----------------------------------------------------------------------
75   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
76   !! $Id$
77   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
78   !!----------------------------------------------------------------------
79CONTAINS
80
81   INTEGER FUNCTION oce_alloc()
82      !!----------------------------------------------------------------------
83      !!                   ***  FUNCTION oce_alloc  ***
84      !!----------------------------------------------------------------------
85      INTEGER :: ierr(4)
86      !!----------------------------------------------------------------------
87      !
88      ALLOCATE( ub   (jpi,jpj,jpk)      , un   (jpi,jpj,jpk)      , ua(jpi,jpj,jpk)       ,     &
89         &      vb   (jpi,jpj,jpk)      , vn   (jpi,jpj,jpk)      , va(jpi,jpj,jpk)       ,     &     
90         &      ua_sv(jpi,jpj,jpk)      , va_sv(jpi,jpj,jpk)      ,                             &     
91         &      wn   (jpi,jpj,jpk)      ,                                                       &
92         &      rotb (jpi,jpj,jpk)      , rotn (jpi,jpj,jpk)      ,                             &   
93         &      hdivb(jpi,jpj,jpk)      , hdivn(jpi,jpj,jpk)      ,                             &
94         &      tsb  (jpi,jpj,jpk,jpts) , tsn  (jpi,jpj,jpk,jpts) , tsa(jpi,jpj,jpk,jpts) ,     &
95         &      rab_b(jpi,jpj,jpk,jpts) , rab_n(jpi,jpj,jpk,jpts) ,                             &
96         &      rn2b (jpi,jpj,jpk)      , rn2  (jpi,jpj,jpk)                              , STAT=ierr(1) )
97         !
98      ALLOCATE(rhd (jpi,jpj,jpk) ,                                         &
99         &     rhop(jpi,jpj,jpk) ,                                         &
100         &     rke(jpi,jpj,jpk)  ,                                         &
101         &     sshb(jpi,jpj)     , sshn(jpi,jpj)   , ssha(jpi,jpj)   ,     &
102         &     ub_b(jpi,jpj)     , un_b(jpi,jpj)   , ua_b(jpi,jpj)   ,     &
103         &     vb_b(jpi,jpj)     , vn_b(jpi,jpj)   , va_b(jpi,jpj)   ,     &
104         &     spgu  (jpi,jpj)   , spgv(jpi,jpj)   ,                       &
105         &     gtsu(jpi,jpj,jpts), gtsv(jpi,jpj,jpts),                     &
106         &     aru(jpi,jpj)      , arv(jpi,jpj)      ,                     &
107         &     gzu(jpi,jpj)      , gzv(jpi,jpj)      ,                     &
108         &     gru(jpi,jpj)      , grv(jpi,jpj)      ,                     &
109         &     ge3ru(jpi,jpj)    , ge3rv(jpi,jpj)    ,                     &
110         &     gtui(jpi,jpj,jpts), gtvi(jpi,jpj,jpts),                     &
111         &     arui(jpi,jpj)     , arvi(jpi,jpj)     ,                     &
112         &     gzui(jpi,jpj)     , gzvi(jpi,jpj)     ,                     &
113         &     ge3rui(jpi,jpj)   , ge3rvi(jpi,jpj)   ,                     &
114         &     grui(jpi,jpj)     , grvi(jpi,jpj)     ,                     &
115         &     riceload(jpi,jpj),                             STAT=ierr(2) )
116         !
117      ALLOCATE( snwice_mass(jpi,jpj) , snwice_mass_b(jpi,jpj), snwice_fmass(jpi,jpj) , STAT=ierr(3) )
118         !
119      ALLOCATE( fraqsr_1lev(jpi,jpj) , STAT=ierr(4) )
120         !
121      oce_alloc = MAXVAL( ierr )
122      IF( oce_alloc /= 0 )   CALL ctl_warn('oce_alloc: failed to allocate arrays')
123      !
124   END FUNCTION oce_alloc
125
126   !!======================================================================
127END MODULE oce
Note: See TracBrowser for help on using the repository browser.