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

source: branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/DTA/dtasal.F90 @ 2436

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

v3.3beta: suppress obsolete key_orca_lev10

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