Changeset 353
- Timestamp:
- 2005-12-12T15:20:26+01:00 (19 years ago)
- Location:
- trunk/NEMO/OPA_SRC/OBC
- Files:
-
- 6 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/NEMO/OPA_SRC/OBC/obc_oce.F90
r247 r353 11 11 !! 8.0 01/91 (CLIPPER) Original code 12 12 !! 8.5 06/02 (C. Talandier) modules 13 !! 06/04 (F. Durand) ORCA2_ZIND config 14 !! 06/04 (F. Durand) Dimensions of arrays vsdta, tsdta, ssdta, 15 !! vndta, tndta, sndta, uwdta, twdta, swdta, uedta, tedta, sedta 16 !! are defined to the actual dimensions of the OBs i.e. 17 !! (jpisd:jpisf,jpk,jptobc) for the South OB 18 !! (jpind:jpinf,jpk,jptobc) for the North OB 19 !! (jpjwd:jpjwf,jpk,jptobc) for the West OB 20 !! (jpjed:jpjef,jpk,jptobc) for the East OB 21 !! 13 22 !!---------------------------------------------------------------------- 14 23 !! * Modules used … … 34 43 nobc_dta = 0 , & !: = 0 use the initial state as obc data 35 44 ! ! = 1 read obc data in obcxxx.dta files 36 nmoisold , & !: number of the last read month on the OBC 37 nbef, naft !: index of the aftera and before fields on the OBC45 46 LOGICAL :: ln_obc_clim = .true. !: obc data files are climatological 38 47 39 48 REAL(wp) :: & !!: open boundary namelist (namobc) … … 116 125 ! ! in the obcdyn.F90 routine 117 126 118 REAL(wp), DIMENSION(jpj glo,jpk,1) :: & !:127 REAL(wp), DIMENSION(jpjed:jpjef,jpk,jptobc) :: & !: 119 128 uedta, tedta, sedta !: array used for interpolating monthly data on the east boundary 120 129 … … 167 176 ! ! in the obcdyn.F90 routine 168 177 169 REAL(wp), DIMENSION(jpj glo,jpk,1) :: & !:178 REAL(wp), DIMENSION(jpjwd:jpjwf,jpk,jptobc) :: & !: 170 179 uwdta, twdta, swdta !: array used for interpolating monthly data on the west boundary 171 180 … … 219 228 ! ! in yhe obcdyn.F90 routine 220 229 221 REAL(wp), DIMENSION(jpi glo,jpk,1) :: & !:230 REAL(wp), DIMENSION(jpind:jpinf,jpk,jptobc) :: & !: 222 231 vndta, tndta, sndta !: array used for interpolating monthly data on the north boundary 223 232 … … 270 279 ! ! in the obcdyn.F90 routine 271 280 272 REAL(wp), DIMENSION(jpi glo,jpk,1) :: & !:281 REAL(wp), DIMENSION(jpisd:jpisf,jpk,jptobc) :: & !: 273 282 vsdta, tsdta, ssdta !: array used for interpolating monthly data on the south boundary 274 283 -
trunk/NEMO/OPA_SRC/OBC/obc_par.F90
r247 r353 11 11 !! 8.0 01/91 (CLIPPER) Original code 12 12 !! 9.0 06/02 (C. Talandier) modules 13 !! 06/04 (F. Durand) ORCA_R2_ZIND config 14 !! 06/04 (F. Durand) jptobc is defined as a parameter, 15 !! in order to allow time-dependent OBCs fields on input 13 16 !!---------------------------------------------------------------------- 14 17 !! * Modules used … … 34 37 !! open boundary parameter 35 38 !!--------------------------------------------------------------------- 39 INTEGER, PARAMETER :: & !: time dimension of the BCS fields on input 40 jptobc = 2 36 41 !! * EAST open boundary 37 42 LOGICAL, PARAMETER :: & !: -
trunk/NEMO/OPA_SRC/OBC/obc_par_EEL_R5.h90
r247 r353 9 9 !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt 10 10 !!---------------------------------------------------------------------- 11 INTEGER, PARAMETER :: & !: time dimension for the BCS fields on input 12 jptobc = 2 11 13 12 14 !! * EAST open boundary -
trunk/NEMO/OPA_SRC/OBC/obc_vectopt_loop_substitute.h90
r247 r353 10 10 !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt 11 11 !!---------------------------------------------------------------------- 12 #if defined key_vectopt_loop && defined key_obc 12 #if defined key_vectopt_loop && defined key_obc && ! defined key_mpp_mpi && ! defined key_mpp_shmem 13 13 # define fs_niw0 jpiwob 14 14 # define fs_niw1 jpiwob -
trunk/NEMO/OPA_SRC/OBC/obcdta.F90
r247 r353 19 19 USE daymod ! calendar 20 20 USE in_out_manager ! I/O logical units 21 USE lib_mpp ! distribu ed memory computing21 USE lib_mpp ! distributed memory computing 22 22 USE dynspg_rl ! 23 USE ioipsl 23 24 24 25 … … 30 31 PRIVATE 31 32 32 !! * Accessibility 33 PUBLIC obc_dta ! routines called by step.F90 33 !! * Accessibility 34 PUBLIC obc_dta ! routines called by step.F90 35 36 !! * Shared module variables 37 INTEGER :: & 38 nlecto = 0, & ! switch for the first read 39 ntobc1 , & ! first record used 40 ntobc2 ! second record used 34 41 35 42 !! * Substitutions 36 43 # include "obc_vectopt_loop_substitute.h90" 37 44 !!--------------------------------------------------------------------------------- 38 !! OPA 9.0 , LOCEAN-IPSL (2005) 39 !! $Header$ 40 !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt 45 !! OPA 9.0 , LODYC-IPSL (2003) 41 46 !!--------------------------------------------------------------------------------- 42 47 43 48 CONTAINS 44 49 45 SUBROUTINE obc_dta( kt ) 46 !!--------------------------------------------------------------------------- 47 !! *** SUBROUTINE obc_dta *** 48 !! 49 !! ** Purpose : Find the climatological boundary arrays for the specified date, 50 !! Originally this routine interpolated between monthly fields 51 !! of a climatology. 52 !! When using hydrological section data, we have only on snapshot 53 !! and do not need to interpolate. 54 !! 55 !! ** Method : Determine the current month from kdat, and interpolate for 56 !! the current day. 50 SUBROUTINE obc_dta (kt) 51 !!-------------------------------------------------------------------- 52 !! *** SUBROUTINE obc_dta *** 53 !! 54 !! ** Purpose : 55 !! Find the climatological boundary arrays for the specified date, 56 !! The boundary arrays are netcdf files. Three possible cases: 57 !! - one time frame only in the file (time dimension = 1). 58 !! in that case the boundary data does not change in time. 59 !! - many time frames. In that case, if we have 12 frames 60 !! we assume monthly fields. 61 !! Else, we assume that time_counter is in seconds 62 !! since the beginning of either the current year or a reference 63 !! year given in the namelist. 64 !! (no check is done so far but one would have to check the "unit" 65 !! attribute of variable time_counter). 57 66 !! 58 67 !! History : 59 68 !! ! 98-05 (J.M. Molines) Original code 60 69 !! 8.5 ! 02-10 (C. Talandier, A-M. Treguier) Free surface, F90 61 !!--------------------------------------------------------------------------- 70 !! 9.0 ! 04-06 (F. Durand, A-M. Treguier) Netcdf BC files on input 71 !!-------------------------------------------------------------------- 62 72 !! * Arguments 63 73 INTEGER, INTENT( in ) :: kt ! ocean time-step index 64 74 65 75 !! * Local declarations 66 INTEGER :: ji, jj, jk, ii, ij ! dummy loop indices 67 INTEGER :: inum = 11 ! temporary logical unit 68 INTEGER :: imois, iman, imoisp1 76 INTEGER :: ji, jj, jk, ii, ij ! dummy loop indices 77 INTEGER :: itimo, iman, imois 69 78 INTEGER :: i15 70 INTEGER :: irecl, irec, ios71 INTEGER, PARAMETER :: &72 jpmois = 1, &73 jpf = 374 79 REAL(wp) :: zxy 75 CHARACTER (len=4 ) :: clversion 76 CHARACTER (len=80) :: clcom 77 !!--------------------------------------------------------------------- 78 79 80 IF( lk_dynspg_rl ) CALL obc_dta_psi( kt ) ! update bsf data at open boundaries 81 82 83 ! 0. Initialization of date 84 ! imois is the index (1 to 12) of the first month to be used in the 85 ! interpolation. ndastp is the date (integer format yymmdd) calculated 86 ! in routine day.F starting from ndate0 given in namelist. 87 ! ----------------------------------------------------------------------- 88 89 iman = jpmois 90 i15 = nday/16 91 imois = nmonth + i15 - 1 92 IF( imois == 0 ) imois = iman 93 imoisp1 = MOD( imois + 1, iman ) 94 IF( imoisp1 == 0 ) imoisp1 = iman 95 ! ... zxy is the weight of the after field 96 zxy = FLOAT( nday + 15 - 30 * i15 ) / 30. 97 IF( zxy > 1.01 .OR. zxy < 0. ) THEN 98 IF(lwp) WRITE(numout,*)' ' 99 IF(lwp) WRITE(numout,*)'obc_dta: Pbm with the the weight of the after field zxy ' 100 IF(lwp) WRITE(numout,*)'~~~~~~~~~~~' 101 nstop = nstop + 1 102 END IF 103 104 ! 1. First call: 105 ! ---------------- 106 107 IF( kt == nit000 ) THEN 108 IF(lwp) WRITE(numout,*)' ' 109 IF(lwp) WRITE(numout,*)'obcdta: initial step in obc_dta' 110 IF(lwp) WRITE(numout,*)'~~~~~~ months ',imois,' and', imoisp1,' read' 111 112 ! 1.1 Tangential velocities set to zero 113 ! -------------------------------------- 114 IF( lp_obc_east ) vfoe(:,1:jpkm1) = 0.e0 115 IF( lp_obc_west ) vfow(:,1:jpkm1) = 0.e0 116 IF( lp_obc_south ) ufos(:,1:jpkm1) = 0.e0 117 IF( lp_obc_north ) ufon(:,1:jpkm1) = 0.e0 118 119 ! 1.2 Set the Normal velocity and tracers data for the EAST OBC 120 ! ------------------------------------------------------------- 121 122 IF( lp_obc_east ) THEN 123 80 !! * Ajouts FD 81 INTEGER :: iyrel ! number of years since 1992 82 INTEGER :: isrel ! number of seconds since 1/1/1992 83 INTEGER, SAVE :: itobce, itobcw, & ! number of time steps in OBC files 84 itobcs, itobcn, & ! " " " " 85 itobc 86 INTEGER :: ikprint ! frequency for printouts. 87 INTEGER :: fid_e, fid_w, fid_n, fid_s, fid ! file identifiers 88 LOGICAL :: l_exv 89 INTEGER, DIMENSION(flio_max_dims) :: & 90 f_d , & ! dimensions lenght 91 start, & ! starting index read 92 count ! number of indices to be read 93 ! time_counter variable of BCs 94 REAL(wp),DIMENSION(:),ALLOCATABLE :: ztcobc 95 96 CHARACTER(LEN=25) :: f_name,v_name 97 !!-------------------------------------------------------------------- 98 99 IF( lk_dynspg_rl ) THEN 100 CALL obc_dta_psi (kt) ! update bsf data at open boundaries 101 IF( nobc_dta == 1 .AND. kt == nit000 ) THEN 102 IF(lwp) WRITE(numout,*) ' time-variable psi boundary data not allowed yet' 103 STOP 104 ENDIF 105 ENDIF 106 107 CALL ipslnlf (new_number=numout) 108 109 ! 1. First call: check time frames available in files. 110 ! ------------------------------------------------------- 111 112 IF( kt == nit000 ) THEN 113 114 IF(lwp) WRITE(numout,*) 115 IF(lwp) WRITE(numout,*) 'obc_dta : find boundary data' 116 IF(lwp) WRITE(numout,*) '~~~~~~~' 117 118 IF( nobc_dta == 0 ) THEN 119 IF(lwp) WRITE(numout,*) ' OBC data taken from initial conditions.' 120 ntobc1 = 1 121 ntobc2 = 1 122 ELSE 123 IF(lwp) WRITE(numout,*) ' OBC data taken from netcdf files.' 124 IF(lwp) WRITE(numout,*) ' climatology (T/F):',ln_obc_clim 125 ! check the number of time steps in the files. 126 itobce =0 ; itobcw = 0; itobcn = 0; itobcs = 0 127 v_name = 'time_counter' 128 IF( lp_obc_east ) THEN 129 CALL flioopfd ('obceast_TS.nc',fid_e) 130 CALL flioinqv (fid_e,TRIM(v_name),l_exv,len_dims=f_d) 131 IF( l_exv ) THEN 132 itobce = f_d(1) 133 ELSE 134 WRITE(numout,*) ' Variable ',TRIM(v_name),' not found in file ','obceast_TS.nc' 135 ENDIF 136 ENDIF 137 IF( lp_obc_west ) THEN 138 CALL flioopfd ('obcwest_TS.nc',fid_w) 139 CALL flioinqv (fid_w,TRIM(v_name),l_exv,len_dims=f_d) 140 IF( l_exv ) THEN 141 itobcw = f_d(1) 142 ELSE 143 WRITE(numout,*) ' Variable ',TRIM(v_name),' not found in file ','obcwest_TS.nc' 144 ENDIF 145 ENDIF 146 IF( lp_obc_north ) THEN 147 CALL flioopfd ('obcnorth_TS.nc',fid_n) 148 CALL flioinqv (fid_n,TRIM(v_name),l_exv,len_dims=f_d) 149 IF( l_exv ) THEN 150 itobcn = f_d(1) 151 ELSE 152 WRITE(numout,*) ' Variable ',TRIM(v_name),' not found in file ','obcnorth_TS.nc' 153 ENDIF 154 ENDIF 155 IF( lp_obc_south ) THEN 156 CALL flioopfd ('obcsouth_TS.nc',fid_s) 157 CALL flioinqv (fid_s,TRIM(v_name),l_exv,len_dims=f_d) 158 IF( l_exv ) THEN 159 itobcs = f_d(1) 160 ELSE 161 WRITE(numout,*) ' Variable ',TRIM(v_name),' not found in file ','obcsouth_TS.nc' 162 ENDIF 163 ENDIF 164 165 itobc = MAX(itobce,itobcw,itobcn,itobcs) 166 nstop = 0 167 IF( lp_obc_east .AND. itobce /= itobc ) nstop = nstop+1 168 IF( lp_obc_west .AND. itobcw /= itobc ) nstop = nstop+1 169 IF( lp_obc_north .AND. itobcn /= itobc ) nstop = nstop+1 170 IF( lp_obc_south .AND. itobcs /= itobc ) nstop = nstop+1 171 IF( nstop /= 0 ) THEN 172 IF( lwp ) THEN 173 WRITE(numout,*) ' obcdta : all files must have the same number of time steps' 174 WRITE(numout,*) ' east, west, north, south: ', itobce, itobcw, itobcn, itobcs 175 ENDIF 176 STOP 177 ENDIF 178 IF( itobc == 1 ) THEN 179 IF( lwp ) WRITE(numout,*) ' obcdta found one time step only in the OBC files' 180 ELSE 181 ALLOCATE (ztcobc(itobc)) 182 l_exv = .TRUE. 183 IF( lp_obc_east ) THEN 184 IF( l_exv ) THEN 185 CALL fliogetv (fid_e,TRIM(v_name),ztcobc) 186 l_exv = .FALSE. 187 ENDIF 188 CALL flioclo (fid_e) 189 ENDIF 190 IF( lp_obc_west ) THEN 191 IF( l_exv ) THEN 192 CALL fliogetv (fid_w,TRIM(v_name),ztcobc) 193 l_exv = .FALSE. 194 ENDIF 195 CALL flioclo (fid_w) 196 ENDIF 197 IF( lp_obc_north ) THEN 198 IF( l_exv ) THEN 199 CALL fliogetv (fid_n,TRIM(v_name),ztcobc) 200 l_exv = .FALSE. 201 ENDIF 202 CALL flioclo (fid_n) 203 ENDIF 204 IF( lp_obc_south ) THEN 205 IF( l_exv ) THEN 206 CALL fliogetv (fid_s,TRIM(v_name),ztcobc) 207 l_exv = .FALSE. 208 ENDIF 209 CALL flioclo (fid_s) 210 ENDIF 211 IF( lwp ) WRITE(numout,*) ' obcdta found', itobc,' time steps in the OBC files' 212 IF( .NOT. ln_obc_clim .AND. itobc == 12 ) THEN 213 IF ( lwp ) WRITE(numout,*) ' WARNING: With monthly data we assume climatology' 214 ln_obc_clim = .true. 215 ENDIF 216 ENDIF 217 ENDIF 218 219 ! 1.1 Tangential velocities set to zero 220 ! -------------------------------------- 221 IF( lp_obc_east ) vfoe = 0.0 222 IF( lp_obc_west ) vfow = 0.0 223 IF( lp_obc_south ) ufos = 0.0 224 IF( lp_obc_north ) ufon = 0.0 225 226 ! 1.2 Data temperature, salinity, normal velocities set to zero 227 ! or initial conditions if nobc_dta == 0 228 ! -------------------------------------------------------------- 229 230 IF( lp_obc_east ) THEN 124 231 ! initialisation to zero 125 sedta(:,:,1) = 0.e0 126 tedta(:,:,1) = 0.e0 127 uedta(:,:,1) = 0.e0 128 ! ! ================== ! 129 IF( nobc_dta == 0 ) THEN ! initial state used 130 ! ! ================== ! 131 DO ji = fs_nie0, fs_nie1 ! Vector opt. 232 sedta(:,:,:) = 0.e0 233 tedta(:,:,:) = 0.e0 234 uedta(:,:,:) = 0.e0 235 ! ! ================== ! 236 IF( nobc_dta == 0 ) THEN ! initial state used 237 ! ! ================== ! 238 ! Fills sedta, tedta, uedta (global arrays) 239 ! Remark: this works for njzoom = 1. 240 ! Should the definition of ij include njzoom? 241 DO ji = nie0, nie1 132 242 DO jk = 1, jpkm1 133 DO jj = 1, jpj243 DO jj = nje0p1, nje1m1 134 244 ij = jj -1 + njmpp 135 245 sedta(ij,jk,1) = sn(ji,jj,jk)*tmask(ji,jj,jk) … … 139 249 END DO 140 250 END DO 141 142 DO jk = 1, jpkm1 143 DO jj = 1, jpj 144 ij = jj -1 + njmpp 145 sfoe(jj,jk) = sedta(ij,jk,1)*temsk(jj,jk) 146 tfoe(jj,jk) = tedta(ij,jk,1)*temsk(jj,jk) 147 ufoe(jj,jk) = uedta(ij,jk,1)*uemsk(jj,jk) 148 END DO 149 END DO 150 ! ! =================== ! 151 ELSE ! read in obceast.dta 152 ! ! =================== ! 153 OPEN(UNIT = inum, & 154 IOSTAT = ios, & 155 FILE ='obceast.dta', & 156 FORM ='UNFORMATTED', & 157 ACCESS ='DIRECT', & 158 RECL = 4096 ) 159 IF( ios > 0 ) THEN 160 IF(lwp) WRITE(numout,*) 'obc_dta: Pbm to OPEN the obceast.dta file ' 161 IF(lwp) WRITE(numout,*) '~~~~~~~' 162 nstop = nstop + 1 163 END IF 164 READ(inum,REC=1) clversion,clcom,irecl 165 CLOSE(inum) 166 IF(lwp) WRITE(numout,*)' ' 167 IF(lwp) WRITE(numout,*)' opening obceast.dta with irecl=',irecl 168 OPEN(UNIT = inum, & 169 IOSTAT = ios, & 170 FILE ='obceast.dta', & 171 FORM ='UNFORMATTED', & 172 ACCESS ='DIRECT', & 173 RECL = irecl ) 174 IF( ios > 0 ) THEN 175 IF(lwp) WRITE(numout,*) 'obc_dta: Pbm to OPEN the obceast.dta file ' 176 IF(lwp) WRITE(numout,*) '~~~~~~~' 177 nstop = nstop + 1 178 END IF 179 180 ! ... Read datafile and set temperature, salinity and normal velocity 181 ! ... initialise the sedta, tedta arrays 182 ! ... index 1 refer to before, 2 to after 183 184 DO jk = 1, jpkm1 185 irec = 2 + (jk -1)* jpf 186 READ(inum,REC=irec )((sedta(jj,jk,1),ji=1,1),jj=jpjed, jpjef) 187 READ(inum,REC=irec+1)((tedta(jj,jk,1),ji=1,1),jj=jpjed, jpjef) 188 READ(inum,REC=irec+2)((uedta(jj,jk,1),ji=1,1),jj=jpjed, jpjef) 189 ! ... set climatological temperature and salinity at the boundary 190 ! ... do not interpolate between months : impose values of climatology 191 DO jj = 1, jpj 192 ij = jj -1 + njmpp 193 sfoe(jj,jk) = sedta(ij,jk,1)*temsk(jj,jk) 194 tfoe(jj,jk) = tedta(ij,jk,1)*temsk(jj,jk) 195 END DO 196 END DO 197 CLOSE(inum) 198 199 # if ! defined key_dynspg_fsc 200 ! ... Rigid lid case: make sure uedta is baroclinic velocity 201 ! ... In rigid lid case uedta needs to be the baroclinic component. 202 203 CALL obc_cli( uedta, uclie, fs_nie0, fs_nie1, 0, jpj, njmpp ) 204 205 # endif 206 ! ... Set normal velocity (on nie0, nie1 <=> jpieob) 207 DO jk = 1, jpkm1 208 DO jj = 1, jpj 209 ij = jj -1 + njmpp 210 ufoe(jj,jk) = uedta(ij,jk,1)*uemsk(jj,jk) 211 END DO 212 END DO 213 ENDIF 214 ENDIF 215 216 ! 1.3 Set the Normal velocity and tracers data for the WEST OBC 217 ! ------------------------------------------------------------- 218 219 IF( lp_obc_west ) THEN 220 251 ENDIF 252 ENDIF 253 254 IF( lp_obc_west ) THEN 221 255 ! initialisation to zero 222 swdta(:,:,1) = 0.e0 223 twdta(:,:,1) = 0.e0 224 uwdta(:,:,1) = 0.e0 225 ! ! ================== ! 226 IF( nobc_dta == 0 ) THEN ! initial state used 227 ! ! ================== ! 228 DO ji = fs_niw0, fs_niw1 ! Vector opt. 256 swdta(:,:,:) = 0.e0 257 twdta(:,:,:) = 0.e0 258 uwdta(:,:,:) = 0.e0 259 ! ! ================== ! 260 IF( nobc_dta == 0 ) THEN ! initial state used ! 261 ! ! ================== ! 262 ! Fills swdta, twdta, uwdta (global arrays) 263 ! Remark: this works for njzoom = 1. 264 ! Should the definition of ij include njzoom? 265 DO ji = niw0, niw1 229 266 DO jk = 1, jpkm1 230 DO jj = 1, jpj267 DO jj = njw0p1, njw1m1 231 268 ij = jj -1 + njmpp 232 269 swdta(ij,jk,1) = sn(ji,jj,jk)*tmask(ji,jj,jk) … … 236 273 END DO 237 274 END DO 238 239 DO jk = 1, jpkm1 240 DO jj = 1, jpj 241 ij = jj -1 + njmpp 242 sfow(jj,jk) = swdta(ij,jk,1)*twmsk(jj,jk) 243 tfow(jj,jk) = twdta(ij,jk,1)*twmsk(jj,jk) 244 ufow(jj,jk) = uwdta(ij,jk,1)*uwmsk(jj,jk) 275 ENDIF 276 ENDIF 277 278 IF( lp_obc_north) THEN 279 ! initialisation to zero 280 sndta(:,:,:) = 0.e0 281 tndta(:,:,:) = 0.e0 282 vndta(:,:,:) = 0.e0 283 ! ! ================== ! 284 IF( nobc_dta == 0 ) THEN ! initial state used 285 ! ! ================== ! 286 ! Fills sndta, tndta, vndta (global arrays) 287 ! Remark: this works for njzoom = 1. 288 ! Should the definition of ij include njzoom? 289 DO jj = njn0, njn1 290 DO jk = 1, jpkm1 291 DO ji = nin0p1, nin1m1 292 ii = ji -1 + nimpp 293 sndta(ii,jk,1) = sn(ji,jj,jk)*tmask(ji,jj,jk) 294 tndta(ii,jk,1) = tn(ji,jj,jk)*tmask(ji,jj,jk) 295 vndta(ii,jk,1) = vn(ji,jj,jk)*vmask(ji,jj,jk) 296 END DO 245 297 END DO 246 298 END DO 247 ! ! =================== ! 248 ELSE ! read in obceast.dta 249 ! ! =================== ! 250 OPEN(UNIT = inum, & 251 IOSTAT = ios, & 252 FILE ='obcwest.dta', & 253 FORM ='UNFORMATTED', & 254 ACCESS ='DIRECT', & 255 RECL = 4096 ) 256 IF( ios > 0 ) THEN 257 IF(lwp) WRITE(numout,*) 'obc_dta: Pbm to OPEN the obcwest.dta file ' 258 IF(lwp) WRITE(numout,*) '~~~~~~~' 259 nstop = nstop + 1 260 END IF 261 READ(inum,REC=1) clversion, clcom,irecl 262 CLOSE(inum) 263 IF(lwp) WRITE(numout,*)' ' 264 IF(lwp) WRITE(numout,*)' opening obcwest.dta with irecl=',irecl 265 OPEN(UNIT = inum, & 266 IOSTAT = ios, & 267 FILE ='obcwest.dta', & 268 FORM ='UNFORMATTED', & 269 ACCESS ='DIRECT', & 270 RECL = irecl ) 271 IF( ios > 0 ) THEN 272 IF(lwp) WRITE(numout,*) 'obc_dta: Pbm to OPEN the obcwest.dta file ' 273 IF(lwp) WRITE(numout,*) '~~~~~~~' 274 nstop = nstop + 1 275 END IF 276 277 ! ... Read datafile and set temperature, salinity and normal velocity 278 ! ... initialise the swdta, twdta arrays 279 ! ... index 1 refer to before, 2 to after 280 DO jk = 1, jpkm1 281 irec = 2 + (jk -1)* jpf 282 READ(inum,REC=irec )((swdta(jj,jk,1),ji=1,1),jj=jpjwd, jpjwf) 283 READ(inum,REC=irec+1)((twdta(jj,jk,1),ji=1,1),jj=jpjwd, jpjwf) 284 READ(inum,REC=irec+2)((uwdta(jj,jk,1),ji=1,1),jj=jpjwd, jpjwf) 285 DO jj = 1, jpj 286 ij = jj -1 + njmpp 287 sfow(jj,jk) = swdta(ij,jk,1)*twmsk(jj,jk) 288 tfow(jj,jk) = twdta(ij,jk,1)*twmsk(jj,jk) 299 ENDIF 300 ENDIF 301 302 IF( lp_obc_south ) THEN 303 ! initialisation to zero 304 ssdta(:,:,:) = 0.e0 305 tsdta(:,:,:) = 0.e0 306 vsdta(:,:,:) = 0.e0 307 ! ! ================== ! 308 IF( nobc_dta == 0 ) THEN ! initial state used 309 ! ! ================== ! 310 ! Fills ssdta, tsdta, vsdta (global arrays) 311 ! Remark: this works for njzoom = 1. 312 ! Should the definition of ij include njzoom? 313 DO jj = njs0, njs1 314 DO jk = 1, jpkm1 315 DO ji = nis0p1, nis1m1 316 ii = ji -1 + nimpp 317 ssdta(ii,jk,1) = sn(ji,jj,jk)*tmask(ji,jj,jk) 318 tsdta(ii,jk,1) = tn(ji,jj,jk)*tmask(ji,jj,jk) 319 vsdta(ii,jk,1) = vn(ji,jj,jk)*vmask(ji,jj,jk) 320 END DO 289 321 END DO 290 322 END DO 291 CLOSE(inum) 292 293 #if ! defined key_dynspg_fsc 294 ! ... Rigid lid case: make sure uwdta is baroclinic velocity 295 ! ... In rigid lid case uwdta needs to be the baroclinic component. 296 297 CALL obc_cli( uwdta, ucliw, fs_niw0, fs_niw1, 0, jpj, njmpp ) 298 299 # endif 300 ! ... Set normal velocity (on niw0, niw1 <=> jpiwob) 301 DO jk = 1, jpkm1 302 DO jj = 1, jpj 303 ij = jj -1 + njmpp 304 ufow(jj,jk) = uwdta(ij,jk,1)*uwmsk(jj,jk) 305 END DO 306 END DO 307 ENDIF 308 ENDIF 309 310 ! 1.4 Set the Normal velocity and tracers data for the NORTH OBC 311 ! --------------------------------------------------------------- 312 313 IF( lp_obc_north ) THEN 314 315 ! initialisation to zero 316 sndta(:,:,1) = 0.e0 317 tndta(:,:,1) = 0.e0 318 vndta(:,:,1) = 0.e0 319 320 OPEN(UNIT = inum, & 321 IOSTAT = ios, & 322 FILE ='obcnorth.dta', & 323 FORM ='UNFORMATTED', & 324 ACCESS ='DIRECT', & 325 RECL = 4096 ) 326 IF( ios > 0 ) THEN 327 IF(lwp) WRITE(numout,*) 'obc_dta: Pbm to OPEN the obcnorth.dta file ' 328 IF(lwp) WRITE(numout,*) '~~~~~~~' 329 nstop = nstop + 1 330 END IF 331 READ(inum,REC=1) clversion,clcom,irecl 332 CLOSE(inum) 333 IF(lwp) WRITE(numout,*)' ' 334 IF(lwp) WRITE(numout,*)' opening obcnorth.dta with irecl=',irecl 335 OPEN(UNIT = inum, & 336 IOSTAT = ios, & 337 FILE ='obcnorth.dta', & 338 FORM ='UNFORMATTED', & 339 ACCESS ='DIRECT', & 340 RECL = irecl ) 341 IF( ios > 0 ) THEN 342 IF(lwp) WRITE(numout,*) 'obc_dta: Pbm to OPEN the obcnorth.dta file ' 343 IF(lwp) WRITE(numout,*) '~~~~~~~' 344 nstop = nstop + 1 345 END IF 346 347 ! ... Read datafile and set temperature and salinity 348 ! ... initialise the sndta, tndta arrays 349 ! ... index 1 refer to before, 2 to after 350 DO jk = 1, jpkm1 351 irec = 2 + (jk -1)* jpf 352 READ(inum,REC=irec )((sndta(jj,jk,1),ji=1,1),jj=jpind, jpinf) 353 READ(inum,REC=irec+1)((tndta(jj,jk,1),ji=1,1),jj=jpind, jpinf) 354 READ(inum,REC=irec+2)((vndta(jj,jk,1),ji=1,1),jj=jpind, jpinf) 355 ! ... do not interpolate: impose values of climatology 356 DO ji = 1, jpi 357 ii = ji -1 + nimpp 358 sfon(ji,jk) = sndta(ii,jk,1)*tnmsk(ji,jk) 359 tfon(ji,jk) = tndta(ii,jk,1)*tnmsk(ji,jk) 360 END DO 361 END DO 362 CLOSE(inum) 363 364 # if ! defined key_dynspg_fsc 365 ! ... Rigid lid case: make sure vndta is baroclinic velocity 366 ! In rigid lid case vndta needs to be the baroclinic component. 367 ! substract the barotropic velocity component (zvnbtpe) 368 ! from the total one (vndta). 369 370 CALL obc_cli( vndta, vclin, fs_njn0, fs_njn1, 1, jpi, nimpp ) 371 # endif 372 ! ... Set normal velocity (on njn0, njn1 <=> jpjnob) 373 DO jk = 1, jpkm1 374 DO ji = 1, jpi 375 ii = ji -1 + nimpp 376 vfon(ji,jk) = vndta(ii,jk,1)*vnmsk(ji,jk) 377 END DO 378 END DO 379 380 END IF 381 382 ! 1.5 Set the Normal velocity and tracers data for the SOUTH OBC 383 ! --------------------------------------------------------------- 384 385 IF( lp_obc_south ) THEN 386 387 ! initialisation to zero 388 ssdta(:,:,1) = 0.e0 389 tsdta(:,:,1) = 0.e0 390 vsdta(:,:,1) = 0.e0 391 392 OPEN(UNIT = inum, & 393 IOSTAT = ios, & 394 FILE ='obcsouth.dta', & 395 FORM ='UNFORMATTED', & 396 ACCESS ='DIRECT', & 397 RECL = 4096 ) 398 IF( ios > 0 ) THEN 399 IF(lwp) WRITE(numout,*) 'obc_dta: Pbm to OPEN the obcsouth.dta file ' 400 IF(lwp) WRITE(numout,*) '~~~~~~~' 401 nstop = nstop + 1 402 END IF 403 READ(inum,REC=1) clversion,clcom,irecl 404 CLOSE(inum) 405 IF(lwp) WRITE(numout,*)' ' 406 IF(lwp) WRITE(numout,*)' opening obcsouth.dta with irecl=',irecl 407 OPEN(UNIT = inum, & 408 IOSTAT = ios, & 409 FILE ='obcsouth.dta', & 410 FORM ='UNFORMATTED', & 411 ACCESS ='DIRECT', & 412 RECL = irecl ) 413 IF( ios > 0 ) THEN 414 IF(lwp) WRITE(numout,*) 'obc_dta: Pbm to OPEN the obcsouth.dta file ' 415 IF(lwp) WRITE(numout,*) '~~~~~~~' 416 nstop = nstop + 1 417 END IF 418 323 ENDIF 324 ENDIF 325 326 ENDIF ! end if kt == nit000 327 328 ! 2. Initialize the time we are at. 329 ! Does this every time the routine is called, 330 ! excepted when nobc_dta = 0 331 !--------------------------------------------------------------------- 332 IF( nobc_dta == 0 ) THEN 333 itimo = 1 334 zxy = 0 335 ELSE 336 IF( itobc == 1 ) THEN 337 itimo = 1 338 ELSE IF( itobc == 12 ) THEN ! BC are monthly 339 ! we assume we have climatology in that case 340 iman = 12 341 i15 = nday / 16 342 imois = nmonth + i15 - 1 343 IF( imois == 0 ) imois = iman 344 itimo = imois 345 ELSE 346 IF(lwp) WRITE(numout,*) 'data other than constant or monthly not written yet' 347 STOP 348 ENDIF 349 ENDIF 350 351 ! 2.1 Read two records in the file if necessary 352 ! --------------------------------------------- 353 IF( ( nobc_dta == 1 ) .AND. ( ( kt == nit000 .AND. nlecto == 0 ) .OR. itimo /= ntobc1 ) ) THEN 354 nlecto = 1 355 356 ! Calendar computation 357 IF( itobc == 1 ) THEN ! BC are constant in time 358 ntobc1 = 1 359 ntobc2 = 1 360 ELSE IF( itobc == 12 ) THEN ! BC are monthly 361 ntobc1 = itimo ! first file record used 362 ntobc2 = ntobc1 + 1 ! last file record used 363 ntobc1 = MOD( ntobc1, iman ) 364 IF( ntobc1 == 0 ) ntobc1 = iman 365 ntobc2 = MOD( ntobc2, iman ) 366 IF( ntobc2 == 0 ) ntobc2 = iman 367 IF( lwp ) THEN 368 WRITE(numout,*) ' read monthly obc first record file used ntobc1 ', ntobc1 369 WRITE(numout,*) ' read monthly obc last record file used ntobc2 ', ntobc2 370 ENDIF 371 ELSE 372 !!!!!!!!!!!!!ATTENTION el: A verifier en fction de la convention choisie pour 373 !!!!!!!!!!!!! le codage de nyear, pour les runs interannuels !!!!!!!!!!!!!! 374 !!! attention if ln_obc_clim is true, go back to jan 1st after december 31st 375 iyrel=nyear-1991 376 IF( ( iyrel < 1 ) .OR. ( iyrel > 13 ) ) THEN 377 IF( lwp ) WRITE(numout,*) 'Pb OBCDTA : iyrel' 378 STOP 379 ENDIF 380 ! Compute nb of seconds from 1/1/1992 00:00 : 381 isrel=(365*(iyrel-1)+nday_year)*86400 382 IF( lwp ) THEN 383 WRITE(numout,*)'Nbre de secondes ecoulees depuis le 1/1/1992:' 384 WRITE(numout,*) isrel 385 ENDIF 386 387 ! need to calculate here ntobc1 and ntobc2, the two time steps to be read 388 389 ENDIF 390 ! ======================= ! 391 ! BCs read ! 392 ! ! ======================= ! 393 IF( lp_obc_east ) THEN 419 394 ! ... Read datafile and set temperature, salinity and normal velocity 420 ! ... initialise the ssdta, tsdta arrays 421 ! ... index 1 refer to before, 2 to after 422 DO jk = 1, jpkm1 423 irec = 2 + (jk -1)* jpf 424 READ(inum,REC=irec )((ssdta(jj,jk,1),ji=1,1),jj=jpisd, jpisf) 425 READ(inum,REC=irec+1)((tsdta(jj,jk,1),ji=1,1),jj=jpisd, jpisf) 426 READ(inum,REC=irec+2)((vsdta(jj,jk,1),ji=1,1),jj=jpisd, jpisf) 427 ! ... do not interpolate: impose values of climatology 428 DO ji = 1, jpi 429 ii = ji -1 + nimpp 430 sfos(ji,jk) = ssdta(ii,jk,1)*tsmsk(ji,jk) 431 tfos(ji,jk) = tsdta(ii,jk,1)*tsmsk(ji,jk) 432 END DO 433 END DO 434 CLOSE(inum) 435 436 # if ! defined key_dynspg_fsc 437 ! ... Rigid lid case: make sure vsdta is baroclinic velocity 438 ! ... In rigid lid case vsdta needs to be the baroclinic component. 439 ! ... substract the barotropic velocity component (zvsbtpe) 440 ! ... from the total one (vsdta). 441 442 CALL obc_cli( vsdta, vclis, fs_njs0, fs_njs1, 1, jpi, nimpp ) 443 # endif 444 ! ... Set normal velocity (on njs0, njs1 <=> jpjsob) 445 DO jk = 1, jpkm1 446 DO ji = 1, jpi 447 ii = ji -1 + nimpp 448 vfos(ji,jk) = vsdta(ii,jk,1)*vsmsk(ji,jk) 449 END DO 450 END DO 451 452 END IF 453 454 ELSE 455 456 ! 2. CALL not the first step ... 457 ! do not read but interpolate between months. 458 ! here no interpolation is done but the code is kept as a reminder 459 ! ---------------------------------------------------------------------- 460 461 IF( lp_obc_east ) THEN 395 ! ... initialise the sedta, tedta, uedta arrays 396 CALL flioopfd ('obceast_TS.nc',fid_e) 397 CALL obc_dta_gv (fid_e,'y','vosaline',sedta(:,:,1),jpjef-jpjed+1,jpk,ntobc1) 398 CALL obc_dta_gv (fid_e,'y','vosaline',sedta(:,:,2),jpjef-jpjed+1,jpk,ntobc2) 399 CALL obc_dta_gv (fid_e,'y','votemper',tedta(:,:,1),jpjef-jpjed+1,jpk,ntobc1) 400 CALL obc_dta_gv (fid_e,'y','votemper',tedta(:,:,2),jpjef-jpjed+1,jpk,ntobc2) 401 CALL flioclo (fid_e) 402 403 CALL flioopfd ('obceast_U.nc',fid_e) 404 CALL obc_dta_gv (fid_e,'y','vozocrtx',uedta(:,:,1),jpjef-jpjed+1,jpk,ntobc1) 405 CALL obc_dta_gv (fid_e,'y','vozocrtx',uedta(:,:,2),jpjef-jpjed+1,jpk,ntobc2) 406 CALL flioclo (fid_e) 407 ! Usually printout is done only once at kt = nit000, 408 ! unless nprint (namelist) > 1 409 IF( lwp .AND. ( kt == nit000 .OR. nprint /= 0 ) ) THEN 410 WRITE(numout,*) 411 WRITE(numout,*) ' Read East OBC data records ', ntobc1, ntobc2 412 ikprint = (jpjef-jpjed+1)/20 +1 413 WRITE(numout,*) ' Temperature record 1 - printout every 3 level' 414 CALL prihre( tedta(:,:,1),jpjef-jpjed+1,jpk,1,jpjef-jpjed+1,ikprint, & 415 & jpk, 1, -3, 1., numout ) 416 WRITE(numout,*) 417 WRITE(numout,*) ' Salinity record 1 - printout every 3 level' 418 CALL prihre( sedta(:,:,1), jpjef-jpjed+1, jpk, 1, jpjef-jpjed+1, ikprint, & 419 & jpk, 1, -3, 1., numout ) 420 WRITE(numout,*) 421 WRITE(numout,*) ' Normal velocity U record 1 - printout every 3 level' 422 CALL prihre( uedta(:,:,1), jpjef-jpjed+1, jpk, 1, jpjef-jpjed+1, ikprint, & 423 & jpk, 1, -3, 1., numout ) 424 ENDIF 425 ENDIF 426 427 IF( lp_obc_west ) THEN 428 ! ... Read datafile and set temperature, salinity and normal velocity 429 ! ... initialise the swdta, twdta, uwdta arrays 430 CALL flioopfd ('obcwest_TS.nc',fid_w) 431 CALL obc_dta_gv (fid_w,'y','vosaline',swdta(:,:,1),jpjwf-jpjwd+1,jpk,ntobc1) 432 CALL obc_dta_gv (fid_w,'y','vosaline',swdta(:,:,2),jpjwf-jpjwd+1,jpk,ntobc2) 433 CALL obc_dta_gv (fid_w,'y','votemper',twdta(:,:,1),jpjwf-jpjwd+1,jpk,ntobc1) 434 CALL obc_dta_gv (fid_w,'y','votemper',twdta(:,:,2),jpjwf-jpjwd+1,jpk,ntobc2) 435 CALL flioclo (fid_w) 436 437 CALL flioopfd ('obcwest_U.nc',fid_w) 438 CALL obc_dta_gv (fid_w,'y','vozocrtx',uwdta(:,:,1),jpjwf-jpjwd+1,jpk,ntobc1) 439 CALL obc_dta_gv (fid_w,'y','vozocrtx',uwdta(:,:,2),jpjwf-jpjwd+1,jpk,ntobc2) 440 CALL flioclo (fid_w) 441 442 IF( lwp .AND. ( kt == nit000 .OR. nprint /= 0 ) ) THEN 443 WRITE(numout,*) 444 WRITE(numout,*) ' Read West OBC data records ', ntobc1, ntobc2 445 ikprint = (jpjwf-jpjwd+1)/20 +1 446 WRITE(numout,*) ' Temperature record 1 - printout every 3 level' 447 CALL prihre( twdta(:,:,1), jpjwf-jpjwd+1, jpk, 1, jpjwf-jpjwd+1, ikprint, & 448 & jpk, 1, -3, 1., numout ) 449 WRITE(numout,*) 450 WRITE(numout,*) ' Salinity record 1 - printout every 3 level' 451 CALL prihre( swdta(:,:,1), jpjwf-jpjwd+1, jpk, 1, jpjwf-jpjwd+1, ikprint, & 452 & jpk, 1, -3, 1., numout ) 453 WRITE(numout,*) 454 WRITE(numout,*) ' Normal velocity U record 1 - printout every 3 level' 455 CALL prihre( uwdta(:,:,1), jpjwf-jpjwd+1, jpk, 1, jpjwf-jpjwd+1, ikprint, & 456 & jpk, 1, -3, 1., numout ) 457 ENDIF 458 ENDIF 459 460 IF( lp_obc_north ) THEN 461 CALL flioopfd ('obcnorth_TS.nc',fid_n) 462 CALL obc_dta_gv (fid_n,'x','vosaline',sndta(:,:,1),jpinf-jpind+1,jpk,ntobc1) 463 CALL obc_dta_gv (fid_n,'x','vosaline',sndta(:,:,2),jpinf-jpind+1,jpk,ntobc2) 464 CALL obc_dta_gv (fid_n,'x','votemper',tndta(:,:,1),jpinf-jpind+1,jpk,ntobc1) 465 CALL obc_dta_gv (fid_n,'x','votemper',tndta(:,:,2),jpinf-jpind+1,jpk,ntobc2) 466 CALL flioclo (fid_n) 467 468 CALL flioopfd ('obcnorth_V.nc',fid_n) 469 CALL obc_dta_gv (fid_n,'x','vomecrty',vndta(:,:,1),jpinf-jpind+1,jpk,ntobc1) 470 CALL obc_dta_gv (fid_n,'x','vomecrty',vndta(:,:,2),jpinf-jpind+1,jpk,ntobc2) 471 CALL flioclo (fid_n) 472 473 IF( lwp .AND. ( kt == nit000 .OR. nprint /= 0 ) ) THEN 474 WRITE(numout,*) 475 WRITE(numout,*) ' Read North OBC data records ', ntobc1, ntobc2 476 ikprint = (jpinf-jpind+1)/20 +1 477 WRITE(numout,*) ' Temperature record 1 - printout every 3 level' 478 CALL prihre( tndta(:,:,1), jpinf-jpind+1, jpk, 1, jpinf-jpind+1, ikprint, & 479 & jpk, 1, -3, 1., numout ) 480 WRITE(numout,*) 481 WRITE(numout,*) ' Salinity record 1 - printout every 3 level' 482 CALL prihre( sndta(:,:,1), jpinf-jpind+1, jpk, 1, jpinf-jpind+1, ikprint, & 483 & jpk, 1, -3, 1., numout ) 484 WRITE(numout,*) 485 WRITE(numout,*) ' Normal velocity V record 1 - printout every 3 level' 486 CALL prihre( vndta(:,:,1), jpinf-jpind+1, jpk, 1, jpinf-jpind+1, ikprint, & 487 & jpk, 1, -3, 1., numout ) 488 ENDIF 489 ENDIF 490 491 IF( lp_obc_south ) THEN 492 CALL flioopfd ('obcsouth_TS.nc',fid_s) 493 CALL obc_dta_gv (fid_s,'x','vosaline',ssdta(:,:,1),jpisf-jpisd+1,jpk,ntobc1) 494 CALL obc_dta_gv (fid_s,'x','vosaline',ssdta(:,:,2),jpisf-jpisd+1,jpk,ntobc2) 495 CALL obc_dta_gv (fid_s,'x','votemper',tsdta(:,:,1),jpisf-jpisd+1,jpk,ntobc1) 496 CALL obc_dta_gv (fid_s,'x','votemper',tsdta(:,:,2),jpisf-jpisd+1,jpk,ntobc2) 497 CALL flioclo (fid_s) 498 499 CALL flioopfd ('obcsouth_V.nc',fid_s) 500 CALL obc_dta_gv (fid_s,'x','vomecrty',vsdta(:,:,1),jpisf-jpisd+1,jpk,ntobc1) 501 CALL obc_dta_gv (fid_s,'x','vomecrty',vsdta(:,:,2),jpisf-jpisd+1,jpk,ntobc2) 502 CALL flioclo (fid_s) 503 504 IF( lwp .AND. ( kt == nit000 .OR. nprint /= 0 ) ) THEN 505 WRITE(numout,*) 506 WRITE(numout,*) ' Read South OBC data records ', ntobc1, ntobc2 507 ikprint = (jpisf-jpisd+1)/20 +1 508 WRITE(numout,*) ' Temperature record 1 - printout every 3 level' 509 CALL prihre( tsdta(:,:,1), jpisf-jpisd+1, jpk, 1, jpisf-jpisd+1, ikprint, & 510 & jpk, 1, -3, 1., numout ) 511 WRITE(numout,*) 512 WRITE(numout,*) ' Salinity record 1 - printout every 3 level' 513 CALL prihre( ssdta(:,:,1), jpisf-jpisd+1, jpk, 1, jpisf-jpisd+1, ikprint, & 514 & jpk, 1, -3, 1., numout ) 515 WRITE(numout,*) 516 WRITE(numout,*) ' Normal velocity V record 1 - printout every 3 level' 517 CALL prihre( vsdta(:,:,1), jpisf-jpisd+1, jpk, 1, jpisf-jpisd+1, ikprint, & 518 & jpk, 1, -3, 1., numout ) 519 ENDIF 520 ENDIF 521 522 ENDIF ! end of the test on the condition to read or not the files 523 524 ! 3. Call at every time step : 525 ! Linear interpolation of BCs to current time step 526 ! ---------------------------------------------------- 527 528 IF( itobc == 1 .OR. nobc_dta == 0 ) THEN 529 zxy = 0. 530 ELSE IF( itobc == 12 ) THEN 531 zxy = FLOAT( nday + 15 - 30 * i15 ) / 30. 532 ELSE 533 zxy = (ztcobc(ntobc1)-FLOAT(isrel))/(ztcobc(ntobc1)-ztcobc(ntobc2)) 534 ENDIF 535 536 IF( lp_obc_east ) THEN 537 ! fills sfoe, tfoe, ufoe (local to each processor) 462 538 DO jk = 1, jpkm1 463 DO jj = 1, jpj539 DO jj = nje0p1, nje1m1 464 540 ij = jj -1 + njmpp 465 sfoe(jj,jk) = sedta(ij,jk,1)*temsk(jj,jk) 466 tfoe(jj,jk) = tedta(ij,jk,1)*temsk(jj,jk) 467 ufoe(jj,jk) = uedta(ij,jk,1)*uemsk(jj,jk) 541 sfoe(jj,jk) = ( zxy * sedta(ij,jk,2) + & 542 & (1.-zxy) * sedta(ij,jk,1) ) * temsk(jj,jk) 543 tfoe(jj,jk) = ( zxy * tedta(ij,jk,2) + & 544 & (1.-zxy) * tedta(ij,jk,1) ) * temsk(jj,jk) 545 ufoe(jj,jk) = ( zxy * uedta(ij,jk,2) + & 546 & (1.-zxy) * uedta(ij,jk,1) ) * uemsk(jj,jk) 468 547 END DO 469 548 END DO 470 END IF 471 472 IF( lp_obc_west ) THEN 549 ENDIF 550 551 IF( lp_obc_west ) THEN 552 ! fills sfow, tfow, ufow (local to each processor) 473 553 DO jk = 1, jpkm1 474 DO jj = 1, jpj554 DO jj = njw0p1, njw1m1 475 555 ij = jj -1 + njmpp 476 sfow(jj,jk) = swdta(ij,jk,1)*twmsk(jj,jk) 477 tfow(jj,jk) = twdta(ij,jk,1)*twmsk(jj,jk) 478 ufow(jj,jk) = uwdta(ij,jk,1)*uwmsk(jj,jk) 556 sfow(jj,jk) = ( zxy * swdta(ij,jk,2) + & 557 & (1.-zxy) * swdta(ij,jk,1) ) * twmsk(jj,jk) 558 tfow(jj,jk) = ( zxy * twdta(ij,jk,2) + & 559 & (1.-zxy) * twdta(ij,jk,1) ) * twmsk(jj,jk) 560 ufow(jj,jk) = ( zxy * uwdta(ij,jk,2) + & 561 & (1.-zxy) * uwdta(ij,jk,1) ) * uwmsk(jj,jk) 479 562 END DO 480 563 END DO 481 END IF 482 483 IF( lp_obc_north ) THEN 564 ENDIF 565 566 IF( lp_obc_north ) THEN 567 ! fills sfon, tfon, vfon (local to each processor) 484 568 DO jk = 1, jpkm1 485 DO ji = 1, jpi569 DO ji = nin0p1, nin1m1 486 570 ii = ji -1 + nimpp 487 sfon(ji,jk) = sndta(ii,jk,1)*tnmsk(ji,jk) 488 tfon(ji,jk) = tndta(ii,jk,1)*tnmsk(ji,jk) 489 vfon(ji,jk) = vndta(ii,jk,1)*vnmsk(ji,jk) 571 sfon(ji,jk) = ( zxy * sndta(ii,jk,2) + & 572 & (1.-zxy) * sndta(ii,jk,1) ) * tnmsk(ji,jk) 573 tfon(ji,jk) = ( zxy * tndta(ii,jk,2) + & 574 & (1.-zxy) * tndta(ii,jk,1) ) * tnmsk(ji,jk) 575 vfon(ji,jk) = ( zxy * vndta(ii,jk,2) + & 576 & (1.-zxy) * vndta(ii,jk,1) ) * vnmsk(ji,jk) 490 577 END DO 491 578 END DO 492 END IF 493 494 IF( lp_obc_south ) THEN 579 ENDIF 580 581 IF( lp_obc_south ) THEN 582 ! fills sfos, tfos, vfos (local to each processor) 495 583 DO jk = 1, jpkm1 496 DO ji = 1, jpi 497 ii = ji -1 + nimpp 498 sfos(ji,jk) = ssdta(ii,jk,1)*tsmsk(ji,jk) 499 tfos(ji,jk) = tsdta(ii,jk,1)*tsmsk(ji,jk) 500 vfos(ji,jk) = vsdta(ii,jk,1)*vsmsk(ji,jk) 501 END DO 502 END DO 503 END IF 504 505 END IF 506 584 DO ji = nis0p1, nis1m1 585 ii = ji -1 + nimpp 586 sfos(ji,jk) = ( zxy * ssdta(ii,jk,2) + & 587 & (1.-zxy) * ssdta(ii,jk,1) ) * tsmsk(ji,jk) 588 tfos(ji,jk) = ( zxy * tsdta(ii,jk,2) + & 589 & (1.-zxy) * tsdta(ii,jk,1) ) * tsmsk(ji,jk) 590 vfos(ji,jk) = ( zxy * vsdta(ii,jk,2) + & 591 & (1.-zxy) * vsdta(ii,jk,1) ) * vsmsk(ji,jk) 592 END DO 593 END DO 594 ENDIF 595 507 596 END SUBROUTINE obc_dta 508 597 509 598 # if defined key_dynspg_fsc 510 599 !!----------------------------------------------------------------------------- … … 565 654 zsver1 = bsfic0(1) 566 655 zsver2 = zsver1 567 IF( kt <= kbsfstart ) THEN656 IF( kt <= kbsfstart ) THEN 568 657 zcoef = float(kt)/float(kbsfstart) 569 658 ELSE 570 zcoef = 1. e0659 zcoef = 1. 571 660 END IF 572 661 bsfic(1) = zsver1*zcoef 573 IF( lwp .AND. ( kt <= kbsfstart ) ) THEN662 IF( lwp .AND. ( kt <= kbsfstart ) ) THEN 574 663 IF(lwp) WRITE(numout,*)' ' 575 664 IF(lwp) WRITE(numout,*)'obcdta: spinup phase in obc_dta_psi routine' … … 584 673 ! ---------------------------------------------------------------------------- 585 674 586 IF( ( neuler == 0 ) .AND. ( kt == nit000 ) ) THEN675 IF( ( neuler == 0 ) .AND. ( kt == nit000 ) ) THEN 587 676 z2dtr = 1./rdt 588 677 ELSE … … 593 682 ! ... However, bsfb(ii,ij) is constant along the same coastlines 594 683 ! ... ---> can be improved using an extra array for storing bsficb (before) 595 IF( nbobc > 1 ) THEN684 IF( nbobc > 1 ) THEN 596 685 DO jnic = 1,nbobc - 1 597 gcbic(jnic) = 0. e0686 gcbic(jnic) = 0. 598 687 ip=mnic(0,jnic) 599 688 DO jip = 1,ip … … 601 690 ij = mjic(jip,0,jnic) 602 691 IF( ii >= nldi+ nimpp - 1 .AND. ii <= nlei+ nimpp - 1 .AND. & 603 ij >= nldj+ njmpp - 1 .AND. ij <= nlej+ njmpp - 1 ) THEN692 ij >= nldj+ njmpp - 1 .AND. ij <= nlej+ njmpp - 1 ) THEN 604 693 iii=ii-nimpp+1 605 694 ijj=ij-njmpp+1 … … 615 704 ! ---------------------------------------------------------------- 616 705 617 IF( lp _obc_east )THEN618 619 IF( kt == nit000 .OR. kt <= kbsfstart ) THEN706 IF( lpeastobc ) THEN 707 708 IF( kt == nit000 .OR. kt <= kbsfstart ) THEN 620 709 OPEN(inum,file='obceastbsf.dta') 621 710 READ(inum,*) … … 633 722 END IF 634 723 635 IF( lp _obc_west )THEN724 IF( lpwestobc) THEN 636 725 637 726 IF( kt == nit000 .OR. kt <= kbsfstart ) then … … 646 735 END IF 647 736 DO jj=jpjwd, jpjwfm1 648 bfow(jj) = bfow(jj) *zcoef737 bfow(jj)=bfow(jj)*zcoef 649 738 END DO 650 739 651 740 END IF 652 741 653 IF( lp _obc_south )THEN654 655 IF( kt == nit000 .OR. kt <= kbsfstart )THEN742 IF( lpsouthobc) THEN 743 744 IF( kt == nit000.OR.kt <= kbsfstart ) THEN 656 745 OPEN(inum,file='obcsouthbsf.dta') 657 746 READ(inum,*) … … 669 758 END IF 670 759 671 IF( lp _obc_north )THEN672 IF( kt == nit000 .OR. kt <= kbsfstart )THEN760 IF( lpnorthobc) THEN 761 IF( kt == nit000.OR.kt <= kbsfstart ) THEN 673 762 OPEN(inum,file='obcnorthbsf.dta') 674 763 READ(inum,*) … … 690 779 # endif 691 780 781 SUBROUTINE obc_dta_gv (ifid,cldim,clobc,pdta,kobcij,kobck,ktobc) 782 !!----------------------------------------------------------------------------- 783 !! *** SUBROUTINE obc_dta_gv *** 784 !! 785 !! ** Purpose : Read a OBC forcing field from netcdf file 786 !! Input file are supposed to be 3D e.g. 787 !! - for a South or North OB : longitude x depth x time 788 !! - for a West or East OB : latitude x depth x time 789 !! 790 !! History : 791 !! ! 04-06 (A.-M. Treguier, F. Durand) Original code 792 !! ! 05-02 (J. Bellier, C. Talandier) use fliocom CALL 793 !!---------------------------------------------------------------------------- 794 !! * Arguments 795 INTEGER, INTENT(IN) :: & 796 ifid, & ! netcdf file name identifier 797 kobcij, & ! Horizontal (i or j) dimension of the array 798 kobck, & ! vertical dimension 799 ktobc ! starting time index read 800 CHARACTER(LEN=*), INTENT(IN) :: & 801 cldim, & ! dimension along which is the open boundary ('x' or 'y') 802 clobc ! name of the netcdf variable read 803 REAL, DIMENSION(kobcij,kobck,1), INTENT(OUT) :: & 804 pdta ! 3D array of OBC forcing field 805 806 !! * Local declarations 807 INTEGER :: indim 808 LOGICAL :: l_exv 809 INTEGER,DIMENSION(4) :: f_d, istart, icount 810 REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE :: v_tmp_4 811 !---------------------------------------------------------------------- 812 813 CALL flioinqv (ifid,TRIM(clobc),l_exv,nb_dims=indim,len_dims=f_d) 814 IF( l_exv ) THEN 815 ! checks the number of dimensions 816 IF( indim == 3 ) THEN 817 istart(1:3) = (/ 1, 1, ktobc /) 818 icount(1:3) = (/ kobcij, kobck, 1 /) 819 CALL fliogetv (ifid,TRIM(clobc),pdta,start=istart(1:3),count=icount(1:3)) 820 ELSE IF( indim == 4 ) THEN 821 istart(1:4) = (/ 1, 1, 1, ktobc /) 822 IF( TRIM(cldim) == 'y' ) THEN 823 icount(1:4) = (/ 1, kobcij, kobck, 1 /) 824 ELSE 825 icount(1:4) = (/ kobcij, 1, kobck, 1 /) 826 ENDIF 827 ALLOCATE (v_tmp_4(icount(1),icount(2),icount(3),icount(4))) 828 CALL fliogetv (ifid,TRIM(clobc),v_tmp_4,start=istart(1:4),count=icount(1:4)) 829 IF( TRIM(cldim) == 'y' ) THEN 830 pdta(1:kobcij,1:kobck,1:1) = v_tmp_4(1,1:kobcij,1:kobck,1:1) 831 ELSE 832 pdta(1:kobcij,1:kobck,1:1) = v_tmp_4(1:kobcij,1,1:kobck,1:1) 833 ENDIF 834 DEALLOCATE (v_tmp_4) 835 ELSE 836 IF( lwp ) THEN 837 WRITE(numout,*) ' Problem in OBC file for ',TRIM(clobc),' :' 838 WRITE(numout,*) ' number of dimensions (not 3 or 4) =',indim 839 ENDIF 840 STOP 841 ENDIF 842 ELSE 843 WRITE(numout,*) ' Variable ',TRIM(clobc),' not found' 844 ENDIF 845 846 END SUBROUTINE obc_dta_gv 847 692 848 #else 693 !!-------------------------------------------------------------------- ----------694 !! default option: Dummy moduleNO Open Boundary Conditions695 !!-------------------------------------------------------------------- ----------849 !!-------------------------------------------------------------------- 850 !! default option : Dummy module NO Open Boundary Conditions 851 !!-------------------------------------------------------------------- 696 852 CONTAINS 697 853 SUBROUTINE obc_dta( kt ) ! Dummy routine 698 699 854 INTEGER, INTENT (in) :: kt 855 WRITE(*,*) 'obc_dta: You should not have seen this print! error?', kt 700 856 END SUBROUTINE obc_dta 701 857 #endif 702 858 703 !!===================================================================== =========859 !!===================================================================== 704 860 END MODULE obcdta -
trunk/NEMO/OPA_SRC/OBC/obcini.F90
r292 r353 61 61 62 62 !! * Local declarations 63 INTEGER :: ji, jj, istop 63 INTEGER :: ji, jj, istop , inumfbc 64 64 INTEGER, DIMENSION(4) :: icorner 65 65 REAL(wp) :: zbsic1, zbsic2, zbsic3 … … 69 69 & rdpeob, rdpwob, rdpnob, rdpsob, & 70 70 & zbsic1, zbsic2, zbsic3, & 71 & nbic, volemp, nobc_dta 71 & nbic, volemp, nobc_dta, ln_obc_clim 72 72 !!---------------------------------------------------------------------- 73 73 … … 134 134 IF(lwp) WRITE(numout,*) ' data in file (=1) or nobc_dta = ', nobc_dta 135 135 IF(lwp) WRITE(numout,*) ' initial state used (=0) ' 136 IF(lwp) WRITE(numout,*) ' climatology (true) or not:', ln_obc_clim 136 137 IF( lwp.AND.lp_obc_east ) THEN 137 138 WRITE(numout,*) ' East open boundary :' … … 179 180 lfbcnorth = .FALSE. 180 181 lfbcsouth = .FALSE. 182 inumfbc = 0 181 183 ! ... look for Fixed Boundaries (rdp = 0 ) 182 184 ! ... When specified, lbcxxx flags are set to TRUE and rdpxxx are set to … … 187 189 rdpein = 1e-3 188 190 rdpeob = 1e-3 191 inumfbc = inumfbc+1 189 192 END IF 190 193 IF( lp_obc_west .AND. ( rdpwin == 0 .AND. rdpwob == 0 ) ) THEN … … 192 195 rdpwin = 1e-3 193 196 rdpwob = 1e-3 197 inumfbc = inumfbc+1 194 198 END IF 195 199 IF( lp_obc_north .AND. ( rdpnin == 0 .AND. rdpnob == 0 ) ) THEN … … 197 201 rdpnin = 1e-3 198 202 rdpnob = 1e-3 203 inumfbc = inumfbc+1 199 204 END IF 200 205 IF( lp_obc_south .AND. ( rdpsin == 0 .AND. rdpsob == 0 ) ) THEN … … 202 207 rdpsin = 1e-3 203 208 rdpsob = 1e-3 209 inumfbc = inumfbc+1 204 210 END IF 205 211 … … 682 688 ! 6. Initialization of open boundary variables (u, v, bsf, t, s) 683 689 ! -------------------------------------------------------------- 690 ! only if at least one boundary is radiative 684 691 685 692 ! ... Restart from restart.obc 686 IF (ln_rstart ) THEN693 IF ( inumfbc < nbobc .AND. ln_rstart ) THEN 687 694 CALL obc_rst_lec 688 695 ELSE … … 730 737 ENDIF 731 738 ELSE 732 IF( tmask(jpieob+1,jpjed ,1) /= 0. .OR. &733 tmask(jpieob+1,jpjed+1,1) /= 1. ) THEN734 IF(lwp) WRITE(numout,cform_err)735 IF(lwp) WRITE(numout,*) ' starting point is not a land T-point.'736 IF(lwp) WRITE(numout,*) ' or starting point + 1 is not a ocean T-point.'737 istop = istop + 1738 END IF739 IF( tmask(jpieob+1,jpjef ,1) /= 0. .OR. &740 tmask(jpieob+1,jpjef-1,1) /= 1. ) THEN741 IF(lwp) WRITE(numout,cform_err)742 IF(lwp) WRITE(numout,*) ' ending point is not a land T-point.'743 IF(lwp) WRITE(numout,*) ' or ending point - 1 is not a ocean T-point.'744 istop = istop + 1745 END IF746 739 747 740 ! ... stop if e r r o r (s) detected … … 773 766 ENDIF 774 767 ELSE 775 IF( tmask(jpiwob,jpjwd ,1) /= 0. .OR. &776 tmask(jpiwob,jpjwd+1,1) /= 1. ) THEN777 IF(lwp) WRITE(numout,cform_err)778 IF(lwp) WRITE(numout,*) ' starting point is not a land T-point.'779 IF(lwp) WRITE(numout,*) ' or starting point + 1 is not a ocean T-point.'780 istop = istop + 1781 END IF782 IF ( tmask(jpieob+1,jpjef ,1) /= 0. .OR. &783 tmask(jpieob+1,jpjef-1,1) /= 1. ) THEN784 IF(lwp) WRITE(numout,cform_err)785 IF(lwp) WRITE(numout,*) ' ending point is not a land T-point.'786 IF(lwp) WRITE(numout,*) ' or ending point - 1 is not a ocean T-point.'787 istop = istop + 1788 END IF789 768 790 769 ! ... stop if e r r o r (s) detected … … 816 795 ENDIF 817 796 ELSE 818 IF( tmask(jpind , jpjnob+1,1) /= 0. .OR. &819 tmask(jpind+1, jpjnob+1,1) /= 1. ) THEN820 IF(lwp) WRITE(numout,cform_err)821 IF(lwp) WRITE(numout,*) ' starting point is not a land T-point.'822 IF(lwp) WRITE(numout,*) ' or starting point + 1 is not a ocean T-point.'823 istop = istop + 1824 END IF825 IF( tmask(jpinf ,jpjnob+1,1) /= 0. .OR. &826 tmask(jpinf-1,jpjnob+1,1) /= 1. ) THEN827 IF(lwp) WRITE(numout,cform_err)828 IF(lwp) WRITE(numout,*) ' ending point is not a land T-point.'829 IF(lwp) WRITE(numout,*) ' or ending point - 1 is not a ocean T-point.'830 istop = istop + 1831 END IF832 797 833 798 ! ... stop if e r r o r (s) detected … … 859 824 ENDIF 860 825 ELSE 861 IF( tmask(jpisd , jpjsob,1) /= 0. .OR. &862 tmask(jpisd+1, jpjsob,1) /= 1. ) THEN863 IF(lwp) WRITE(numout,cform_err)864 IF(lwp) WRITE(numout,*) ' starting point is not a land T-point.'865 IF(lwp) WRITE(numout,*) ' or starting point + 1 is not a ocean T-point.'866 istop = istop + 1867 END IF868 IF( tmask(jpisf ,jpjsob,1) /= 0. .OR. &869 tmask(jpisf-1,jpjsob,1) /= 1. ) THEN870 IF(lwp) WRITE(numout,cform_err)871 IF(lwp) WRITE(numout,*) ' ending point is not a land T-point.'872 IF(lwp) WRITE(numout,*) ' or ending point - 1 is not a ocean T-point.'873 istop = istop + 1874 END IF875 826 876 827 ! ... stop if e r r o r (s) detected
Note: See TracChangeset
for help on using the changeset viewer.