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.
dtasal.F90 in branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OPA_SRC/DTA – NEMO

source: branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OPA_SRC/DTA/dtasal.F90 @ 2636

Last change on this file since 2636 was 2636, checked in by gm, 13 years ago

dynamic mem: #785 ; move ctl_stop & warn in lib_mpp to avoid a circular dependency + ctl_stop improvment

  • Property svn:keywords set to Id
File size: 9.9 KB
Line 
1MODULE dtasal
2   !!======================================================================
3   !!                     ***  MODULE  dtasal  ***
4   !! Ocean data  :  read ocean salinity data from monthly atlas data
5   !!=====================================================================
6   !! History :  OPA  ! 1991-03  ()  Original code
7   !!             -   ! 1992-07  (M. Imbard)
8   !!            8.0  ! 1999-10  (M.A. Foujols, M. Imbard)  NetCDF FORMAT
9   !!   NEMO     1.0  ! 2002-06  (G. Madec)  F90: Free form and module
10   !!            3.3  ! 2010-10  (C. Bricaud, S. Masson)  use of fldread
11   !!----------------------------------------------------------------------
12#if defined key_dtasal   ||   defined key_esopa
13   !!----------------------------------------------------------------------
14   !!   'key_dtasal'                                          salinity data
15   !!----------------------------------------------------------------------
16   !!   dta_sal      : read ocean salinity data
17   !!----------------------------------------------------------------------
18   USE oce             ! ocean dynamics and tracers
19   USE dom_oce         ! ocean space and time domain
20   USE phycst          ! physical constants
21   USE fldread         ! read input fields
22   USE in_out_manager  ! I/O manager
23   USE lib_mpp         ! MPP library
24
25   IMPLICIT NONE
26   PRIVATE
27
28   PUBLIC   dta_sal    ! called by step.F90 and inidta.F90
29
30   LOGICAL , PUBLIC, PARAMETER                     ::   lk_dtasal = .TRUE. !: salinity data flag
31   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:,:) ::   s_dta              !: salinity data at given time-step
32
33   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_sal   ! structure of input SST (file informations, fields read)
34
35   !! * Substitutions
36#  include "domzgr_substitute.h90"
37   !!----------------------------------------------------------------------
38   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
39   !! $Id$
40   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
41   !!----------------------------------------------------------------------
42CONTAINS
43
44   SUBROUTINE dta_sal( kt )
45      !!----------------------------------------------------------------------
46      !!                   ***  ROUTINE dta_sal  ***
47      !!       
48      !! ** Purpose :   Reads monthly salinity data
49      !!             
50      !! ** Method  : - Read on unit numsdt the monthly salinity data interpo-
51      !!                lated onto the model grid.
52      !!              - At each time step, a linear interpolation is applied
53      !!                between two monthly values.
54      !!----------------------------------------------------------------------
55      INTEGER, INTENT(in) ::   kt   ! ocean time step
56      !
57      INTEGER ::   ji, jj, jk, jl, jkk       ! local loop indicies
58      INTEGER ::   ik, ierr0, ierr1, ierr2   ! local integers
59#if defined key_tradmp
60      INTEGER ::   il0, il1, ii0, ii1, ij0, ij1   ! local integers
61#endif
62      REAL(wp)::   zl
63      REAL(wp), DIMENSION(jpk) :: zsaldta         ! auxiliary array for interpolation
64      CHARACTER(len=100)       :: cn_dir          ! Root directory for location of ssr files
65      TYPE(FLD_N)              :: sn_sal
66      LOGICAL , SAVE           :: linit_sal = .FALSE.
67      !!
68      NAMELIST/namdta_sal/   cn_dir, sn_sal
69      !!----------------------------------------------------------------------
70     
71      ! 1. Initialization
72      ! -----------------------
73     
74      IF( kt == nit000 .AND. ( .NOT. linit_sal ) ) THEN
75       
76         !                         ! set file information
77         cn_dir = './'             ! directory in which the model is executed
78         ! ... default values (NB: frequency positive => hours, negative => months)
79         !            !   file    ! frequency ! variable ! time intep !  clim   ! 'yearly' or ! weights  ! rotation   !
80         !            !   name    !  (hours)  !  name    !   (T/F)    !  (T/F)  !  'monthly'  ! filename ! pairs      !
81         sn_sal = FLD_N( 'salinity',  -1.     ,'vosaline',  .false.   , .true.  ,  'monthly'  , ''       , ''         )
82
83         REWIND ( numnam )         ! read in namlist namdta_sal
84         READ( numnam, namdta_sal ) 
85
86         IF(lwp) THEN              ! control print
87            WRITE(numout,*)
88            WRITE(numout,*) 'dta_sal : Salinity Climatology '
89            WRITE(numout,*) '~~~~~~~ '
90         ENDIF
91                                ALLOCATE( sf_sal(1)                    , STAT=ierr0 )
92                                ALLOCATE( sf_sal(1)%fnow(jpi,jpj,jpk)  , STAT=ierr1 )
93         IF( sn_sal%ln_tint )   ALLOCATE( sf_sal(1)%fdta(jpi,jpj,jpk,2), STAT=ierr2 )
94         IF( ierr0+ierr1+ierr2 > 0 )   CALL ctl_stop( 'STOP', 'dta_sal: unable to allocate sf_sal structure' )
95         !                         ! fill sf_sal with sn_sal and control print
96         CALL fld_fill( sf_sal, (/ sn_sal /), cn_dir, 'dta_sal', 'Salinity data', 'namdta_sal' )
97         linit_sal = .TRUE.       
98      ENDIF
99     
100      ! 2. Read monthly file
101      ! -------------------
102     
103      CALL fld_read( kt, 1, sf_sal )
104
105      IF( lwp .AND. kt == nit000 ) THEN
106         WRITE(numout,*)
107         WRITE(numout,*) ' read Levitus salinity ok'
108         WRITE(numout,*)
109      ENDIF
110
111#if defined key_tradmp
112      IF( cp_cfg == "orca"  .AND. jp_cfg == 2 ) THEN     !  ORCA_R2 configuration
113         !
114         ij0 = 101   ;   ij1 = 109
115         ii0 = 141   ;   ii1 = 155   
116         DO jj = mj0(ij0), mj1(ij1)                  ! Reduced salinity in the Alboran Sea
117            DO ji = mi0(ii0), mi1(ii1)
118               sf_sal(1)%fnow(ji,jj,13:13) = sf_sal(1)%fnow(ji,jj,13:13) - 0.15
119               sf_sal(1)%fnow(ji,jj,14:15) = sf_sal(1)%fnow(ji,jj,14:15) - 0.25
120               sf_sal(1)%fnow(ji,jj,16:17) = sf_sal(1)%fnow(ji,jj,16:17) - 0.30
121               sf_sal(1)%fnow(ji,jj,18:25) = sf_sal(1)%fnow(ji,jj,18:25) - 0.35
122            END DO
123         END DO
124         !
125         IF( nn_cla == 1 ) THEN 
126            !                                         ! New salinity profile at Gibraltar
127            il0 = 138   ;   il1 = 138   
128            ij0 = 101   ;   ij1 = 102
129            ii0 = 139   ;   ii1 = 139   
130            DO jl = mi0(il0), mi1(il1)
131               DO jj = mj0(ij0), mj1(ij1)
132                  DO ji = mi0(ii0), mi1(ii1)
133                        sf_sal(1)%fnow(ji,jj,:) = sf_sal(1)%fnow(jl,jj,:)
134                  END DO
135               END DO
136            END DO
137            !                                         ! New salinity profile at Bab el Mandeb
138            il0 = 164   ;   il1 = 164   
139            ij0 =  87   ;   ij1 =  88
140            ii0 = 161   ;   ii1 = 163   
141            DO jl = mi0(il0), mi1(il1)
142               DO jj = mj0(ij0), mj1(ij1)
143                  DO ji = mi0(ii0), mi1(ii1)
144                     sf_sal(1)%fnow(ji,jj,:) = sf_sal(1)%fnow(jl,jj,:)
145                  END DO
146               END DO
147            END DO
148            !
149         ENDIF
150            !
151      ENDIF
152#endif   
153       
154      s_dta(:,:,:)=sf_sal(1)%fnow(:,:,:)
155       
156      IF( ln_sco ) THEN
157         DO jj = 1, jpj                  ! interpolation of salinites
158            DO ji = 1, jpi
159               DO jk = 1, jpk
160                  zl=fsdept_0(ji,jj,jk)
161                  IF(zl < gdept_0(1)  ) zsaldta(jk) =  s_dta(ji,jj,1    ) 
162                  IF(zl > gdept_0(jpk)) zsaldta(jk) =  s_dta(ji,jj,jpkm1) 
163                  DO jkk = 1, jpkm1
164                     IF((zl-gdept_0(jkk))*(zl-gdept_0(jkk+1)).le.0.0) THEN
165                          zsaldta(jk) = s_dta(ji,jj,jkk)                                 &
166                                     &           + (zl-gdept_0(jkk))/(gdept_0(jkk+1)-gdept_0(jkk))      &
167                                     &                              *(s_dta(ji,jj,jkk+1) - s_dta(ji,jj,jkk))
168                     ENDIF
169                  END DO
170               END DO
171               DO jk = 1, jpkm1
172                  s_dta(ji,jj,jk) = zsaldta(jk) 
173               END DO
174               s_dta(ji,jj,jpk) = 0.0 
175            END DO
176         END DO
177           
178         IF( lwp .AND. kt==nn_it000 ) THEN
179            WRITE(numout,*)
180            WRITE(numout,*) ' Levitus salinity data interpolated to s-coordinate'
181            WRITE(numout,*)
182         ENDIF
183
184      ELSE
185         !                                  ! Mask
186         s_dta(:,:,:) = s_dta(:,:,:) * tmask(:,:,:)
187         s_dta(:,:,jpk) = 0.e0
188         IF( ln_zps ) THEN               ! z-coord. partial steps
189            DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
190               DO ji = 1, jpi
191                  ik = mbkt(ji,jj)
192                  IF( ik > 1 ) THEN
193                     zl = ( gdept_0(ik) - fsdept_0(ji,jj,ik) ) / ( gdept_0(ik) - gdept_0(ik-1) )
194                     s_dta(ji,jj,ik) = (1.-zl) * s_dta(ji,jj,ik) + zl * s_dta(ji,jj,ik-1)
195                  ENDIF
196               END DO
197            END DO
198         ENDIF
199      ENDIF
200       
201      IF( lwp .AND. kt == nit000 ) THEN
202         WRITE(numout,*)' salinity Levitus '
203         WRITE(numout,*)
204         WRITE(numout,*)'  level = 1'
205         CALL prihre(s_dta(:,:,1),    jpi,jpj,1,jpi,20,1,jpj,20,1.,numout)
206         WRITE(numout,*)'  level = ',jpk/2
207         CALL prihre(s_dta(:,:,jpk/2),jpi,jpj,1,jpi,20,1,jpj,20,1.,numout)           
208         WRITE(numout,*) '  level = ',jpkm1
209         CALL prihre(s_dta(:,:,jpkm1),jpi,jpj,1,jpi,20,1,jpj,20,1.,numout)
210      ENDIF
211      !
212   END SUBROUTINE dta_sal
213
214#else
215   !!----------------------------------------------------------------------
216   !!   Default option:                                    NO salinity data
217   !!----------------------------------------------------------------------
218   LOGICAL , PUBLIC, PARAMETER ::   lk_dtasal = .FALSE.   !: salinity data flag
219CONTAINS
220   SUBROUTINE dta_sal( kt )        ! Empty routine
221      WRITE(*,*) 'dta_sal: You should not have seen this print! error?', kt
222   END SUBROUTINE dta_sal
223#endif
224   !!======================================================================
225END MODULE dtasal
Note: See TracBrowser for help on using the repository browser.