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

source: trunk/NEMOGCM/NEMO/OPA_SRC/DTA/dtasal.F90 @ 2715

Last change on this file since 2715 was 2715, checked in by rblod, 13 years ago

First attempt to put dynamic allocation on the trunk

  • Property svn:keywords set to Id
File size: 10.2 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, ierr, 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
92                                   ! Allocate salinity data array
93                                ALLOCATE( s_dta(jpi,jpj,jpk)           , STAT=ierr  )
94         IF( ierr > 0              )   CALL ctl_stop( 'STOP', 'dta_sal: unable to allocate s_dta array' )
95                                   ! Allocate sf_tem structure
96                                ierr2 = 0
97                                ALLOCATE( sf_sal(1)                    , STAT=ierr0 )
98                                ALLOCATE( sf_sal(1)%fnow(jpi,jpj,jpk)  , STAT=ierr1 )
99         IF( sn_sal%ln_tint )   ALLOCATE( sf_sal(1)%fdta(jpi,jpj,jpk,2), STAT=ierr2 )
100         IF( ierr0+ierr1+ierr2 > 0 )   CALL ctl_stop( 'STOP', 'dta_sal: unable to allocate sf_sal structure' )
101         !                         ! fill sf_sal with sn_sal and control print
102         CALL fld_fill( sf_sal, (/ sn_sal /), cn_dir, 'dta_sal', 'Salinity data', 'namdta_sal' )
103         linit_sal = .TRUE.       
104      ENDIF
105     
106      ! 2. Read monthly file
107      ! -------------------
108     
109      CALL fld_read( kt, 1, sf_sal )
110
111      IF( lwp .AND. kt == nit000 ) THEN
112         WRITE(numout,*)
113         WRITE(numout,*) ' read Levitus salinity ok'
114         WRITE(numout,*)
115      ENDIF
116
117#if defined key_tradmp
118      IF( cp_cfg == "orca"  .AND. jp_cfg == 2 ) THEN     !  ORCA_R2 configuration
119         !
120         ij0 = 101   ;   ij1 = 109
121         ii0 = 141   ;   ii1 = 155   
122         DO jj = mj0(ij0), mj1(ij1)                  ! Reduced salinity in the Alboran Sea
123            DO ji = mi0(ii0), mi1(ii1)
124               sf_sal(1)%fnow(ji,jj,13:13) = sf_sal(1)%fnow(ji,jj,13:13) - 0.15
125               sf_sal(1)%fnow(ji,jj,14:15) = sf_sal(1)%fnow(ji,jj,14:15) - 0.25
126               sf_sal(1)%fnow(ji,jj,16:17) = sf_sal(1)%fnow(ji,jj,16:17) - 0.30
127               sf_sal(1)%fnow(ji,jj,18:25) = sf_sal(1)%fnow(ji,jj,18:25) - 0.35
128            END DO
129         END DO
130         !
131         IF( nn_cla == 1 ) THEN 
132            !                                         ! New salinity profile at Gibraltar
133            il0 = 138   ;   il1 = 138   
134            ij0 = 101   ;   ij1 = 102
135            ii0 = 139   ;   ii1 = 139   
136            DO jl = mi0(il0), mi1(il1)
137               DO jj = mj0(ij0), mj1(ij1)
138                  DO ji = mi0(ii0), mi1(ii1)
139                        sf_sal(1)%fnow(ji,jj,:) = sf_sal(1)%fnow(jl,jj,:)
140                  END DO
141               END DO
142            END DO
143            !                                         ! New salinity profile at Bab el Mandeb
144            il0 = 164   ;   il1 = 164   
145            ij0 =  87   ;   ij1 =  88
146            ii0 = 161   ;   ii1 = 163   
147            DO jl = mi0(il0), mi1(il1)
148               DO jj = mj0(ij0), mj1(ij1)
149                  DO ji = mi0(ii0), mi1(ii1)
150                     sf_sal(1)%fnow(ji,jj,:) = sf_sal(1)%fnow(jl,jj,:)
151                  END DO
152               END DO
153            END DO
154            !
155         ENDIF
156            !
157      ENDIF
158#endif   
159       
160      s_dta(:,:,:)=sf_sal(1)%fnow(:,:,:)
161       
162      IF( ln_sco ) THEN
163         DO jj = 1, jpj                  ! interpolation of salinites
164            DO ji = 1, jpi
165               DO jk = 1, jpk
166                  zl=fsdept_0(ji,jj,jk)
167                  IF(zl < gdept_0(1)  ) zsaldta(jk) =  s_dta(ji,jj,1    ) 
168                  IF(zl > gdept_0(jpk)) zsaldta(jk) =  s_dta(ji,jj,jpkm1) 
169                  DO jkk = 1, jpkm1
170                     IF((zl-gdept_0(jkk))*(zl-gdept_0(jkk+1)).le.0.0) THEN
171                          zsaldta(jk) = s_dta(ji,jj,jkk)                                 &
172                                     &           + (zl-gdept_0(jkk))/(gdept_0(jkk+1)-gdept_0(jkk))      &
173                                     &                              *(s_dta(ji,jj,jkk+1) - s_dta(ji,jj,jkk))
174                     ENDIF
175                  END DO
176               END DO
177               DO jk = 1, jpkm1
178                  s_dta(ji,jj,jk) = zsaldta(jk) 
179               END DO
180               s_dta(ji,jj,jpk) = 0.0 
181            END DO
182         END DO
183           
184         IF( lwp .AND. kt==nn_it000 ) THEN
185            WRITE(numout,*)
186            WRITE(numout,*) ' Levitus salinity data interpolated to s-coordinate'
187            WRITE(numout,*)
188         ENDIF
189
190      ELSE
191         !                                  ! Mask
192         s_dta(:,:,:) = s_dta(:,:,:) * tmask(:,:,:)
193         s_dta(:,:,jpk) = 0.e0
194         IF( ln_zps ) THEN               ! z-coord. partial steps
195            DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
196               DO ji = 1, jpi
197                  ik = mbkt(ji,jj)
198                  IF( ik > 1 ) THEN
199                     zl = ( gdept_0(ik) - fsdept_0(ji,jj,ik) ) / ( gdept_0(ik) - gdept_0(ik-1) )
200                     s_dta(ji,jj,ik) = (1.-zl) * s_dta(ji,jj,ik) + zl * s_dta(ji,jj,ik-1)
201                  ENDIF
202               END DO
203            END DO
204         ENDIF
205      ENDIF
206       
207      IF( lwp .AND. kt == nit000 ) THEN
208         WRITE(numout,*)' salinity Levitus '
209         WRITE(numout,*)
210         WRITE(numout,*)'  level = 1'
211         CALL prihre(s_dta(:,:,1),    jpi,jpj,1,jpi,20,1,jpj,20,1.,numout)
212         WRITE(numout,*)'  level = ',jpk/2
213         CALL prihre(s_dta(:,:,jpk/2),jpi,jpj,1,jpi,20,1,jpj,20,1.,numout)           
214         WRITE(numout,*) '  level = ',jpkm1
215         CALL prihre(s_dta(:,:,jpkm1),jpi,jpj,1,jpi,20,1,jpj,20,1.,numout)
216      ENDIF
217      !
218   END SUBROUTINE dta_sal
219
220#else
221   !!----------------------------------------------------------------------
222   !!   Default option:                                    NO salinity data
223   !!----------------------------------------------------------------------
224   LOGICAL , PUBLIC, PARAMETER ::   lk_dtasal = .FALSE.   !: salinity data flag
225CONTAINS
226   SUBROUTINE dta_sal( kt )        ! Empty routine
227      WRITE(*,*) 'dta_sal: You should not have seen this print! error?', kt
228   END SUBROUTINE dta_sal
229#endif
230   !!======================================================================
231END MODULE dtasal
Note: See TracBrowser for help on using the repository browser.