istate.F90 27 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571
  1. MODULE istate
  2. !!======================================================================
  3. !! *** MODULE istate ***
  4. !! Ocean state : initial state setting
  5. !!=====================================================================
  6. !! History : OPA ! 1989-12 (P. Andrich) Original code
  7. !! 5.0 ! 1991-11 (G. Madec) rewritting
  8. !! 6.0 ! 1996-01 (G. Madec) terrain following coordinates
  9. !! 8.0 ! 2001-09 (M. Levy, M. Ben Jelloul) istate_eel
  10. !! 8.0 ! 2001-09 (M. Levy, M. Ben Jelloul) istate_uvg
  11. !! NEMO 1.0 ! 2003-08 (G. Madec, C. Talandier) F90: Free form, modules + EEL R5
  12. !! - ! 2004-05 (A. Koch-Larrouy) istate_gyre
  13. !! 2.0 ! 2006-07 (S. Masson) distributed restart using iom
  14. !! 3.3 ! 2010-10 (C. Ethe) merge TRC-TRA
  15. !! 3.4 ! 2011-04 (G. Madec) Merge of dtatem and dtasal & suppression of tb,tn/sb,sn
  16. !!----------------------------------------------------------------------
  17. !!----------------------------------------------------------------------
  18. !! istate_init : initial state setting
  19. !! istate_tem : analytical profile for initial Temperature
  20. !! istate_sal : analytical profile for initial Salinity
  21. !! istate_eel : initial state setting of EEL R5 configuration
  22. !! istate_gyre : initial state setting of GYRE configuration
  23. !! istate_uvg : initial velocity in geostropic balance
  24. !!----------------------------------------------------------------------
  25. USE oce ! ocean dynamics and active tracers
  26. USE dom_oce ! ocean space and time domain
  27. USE c1d ! 1D vertical configuration
  28. USE daymod ! calendar
  29. USE eosbn2 ! eq. of state, Brunt Vaisala frequency (eos routine)
  30. USE ldftra_oce ! ocean active tracers: lateral physics
  31. USE zdf_oce ! ocean vertical physics
  32. USE phycst ! physical constants
  33. USE dtatsd ! data temperature and salinity (dta_tsd routine)
  34. USE dtauvd ! data: U & V current (dta_uvd routine)
  35. USE zpshde ! partial step: hor. derivative (zps_hde routine)
  36. USE eosbn2 ! equation of state (eos bn2 routine)
  37. USE domvvl ! varying vertical mesh
  38. USE dynspg_oce ! pressure gradient schemes
  39. USE dynspg_flt ! filtered free surface
  40. USE sol_oce ! ocean solver variables
  41. !
  42. USE in_out_manager ! I/O manager
  43. USE iom ! I/O library
  44. USE lib_mpp ! MPP library
  45. USE restart ! restart
  46. USE wrk_nemo ! Memory allocation
  47. USE timing ! Timing
  48. IMPLICIT NONE
  49. PRIVATE
  50. PUBLIC istate_init ! routine called by step.F90
  51. !! * Substitutions
  52. # include "domzgr_substitute.h90"
  53. # include "vectopt_loop_substitute.h90"
  54. !!----------------------------------------------------------------------
  55. !! NEMO/OPA 3.7 , NEMO Consortium (2014)
  56. !! $Id: istate.F90 4990 2014-12-15 16:42:49Z timgraham $
  57. !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
  58. !!----------------------------------------------------------------------
  59. CONTAINS
  60. SUBROUTINE istate_init
  61. !!----------------------------------------------------------------------
  62. !! *** ROUTINE istate_init ***
  63. !!
  64. !! ** Purpose : Initialization of the dynamics and tracer fields.
  65. !!----------------------------------------------------------------------
  66. INTEGER :: ji, jj, jk ! dummy loop indices
  67. REAL(wp), POINTER, DIMENSION(:,:,:,:) :: zuvd ! U & V data workspace
  68. !!----------------------------------------------------------------------
  69. !
  70. IF( nn_timing == 1 ) CALL timing_start('istate_init')
  71. !
  72. IF(lwp) WRITE(numout,*) ' '
  73. IF(lwp) WRITE(numout,*) 'istate_ini : Initialization of the dynamics and tracers'
  74. IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
  75. CALL dta_tsd_init ! Initialisation of T & S input data
  76. IF( lk_c1d ) CALL dta_uvd_init ! Initialization of U & V input data
  77. rhd (:,:,: ) = 0._wp ; rhop (:,:,: ) = 0._wp ! set one for all to 0 at level jpk
  78. rn2b (:,:,: ) = 0._wp ; rn2 (:,:,: ) = 0._wp ! set one for all to 0 at levels 1 and jpk
  79. tsa (:,:,:,:) = 0._wp ! set one for all to 0 at level jpk
  80. rab_b(:,:,:,:) = 0._wp ; rab_n(:,:,:,:) = 0._wp ! set one for all to 0 at level jpk
  81. IF( ln_rstart ) THEN ! Restart from a file
  82. ! ! -------------------
  83. CALL rst_read ! Read the restart file
  84. CALL day_init ! model calendar (using both namelist and restart infos)
  85. ELSE
  86. ! ! Start from rest
  87. ! ! ---------------
  88. numror = 0 ! define numror = 0 -> no restart file to read
  89. neuler = 0 ! Set time-step indicator at nit000 (euler forward)
  90. CALL day_init ! model calendar (using both namelist and restart infos)
  91. ! ! Initialization of ocean to zero
  92. ! before fields ! now fields
  93. sshb (:,:) = 0._wp ; sshn (:,:) = 0._wp
  94. ub (:,:,:) = 0._wp ; un (:,:,:) = 0._wp
  95. vb (:,:,:) = 0._wp ; vn (:,:,:) = 0._wp
  96. rotb (:,:,:) = 0._wp ; rotn (:,:,:) = 0._wp
  97. hdivb(:,:,:) = 0._wp ; hdivn(:,:,:) = 0._wp
  98. !
  99. IF( cp_cfg == 'eel' ) THEN
  100. CALL istate_eel ! EEL configuration : start from pre-defined U,V T-S fields
  101. ELSEIF( cp_cfg == 'gyre' ) THEN
  102. CALL istate_gyre ! GYRE configuration : start from pre-defined T-S fields
  103. ELSE ! Initial T-S, U-V fields read in files
  104. IF ( ln_tsd_init ) THEN ! read 3D T and S data at nit000
  105. CALL dta_tsd( nit000, tsb )
  106. tsn(:,:,:,:) = tsb(:,:,:,:)
  107. !
  108. ELSE ! Initial T-S fields defined analytically
  109. CALL istate_t_s
  110. ENDIF
  111. IF ( ln_uvd_init .AND. lk_c1d ) THEN ! read 3D U and V data at nit000
  112. CALL wrk_alloc( jpi, jpj, jpk, 2, zuvd )
  113. CALL dta_uvd( nit000, zuvd )
  114. ub(:,:,:) = zuvd(:,:,:,1) ; un(:,:,:) = ub(:,:,:)
  115. vb(:,:,:) = zuvd(:,:,:,2) ; vn(:,:,:) = vb(:,:,:)
  116. CALL wrk_dealloc( jpi, jpj, jpk, 2, zuvd )
  117. ENDIF
  118. ENDIF
  119. !
  120. CALL eos( tsb, rhd, rhop, gdept_0(:,:,:) ) ! before potential and in situ densities
  121. #if ! defined key_c1d
  122. IF( ln_zps .AND. .NOT. ln_isfcav) &
  123. & CALL zps_hde ( nit000, jpts, tsb, gtsu, gtsv, & ! Partial steps: before horizontal gradient
  124. & rhd, gru , grv ) ! of t, s, rd at the last ocean level
  125. IF( ln_zps .AND. ln_isfcav) &
  126. & CALL zps_hde_isf( nit000, jpts, tsb, gtsu, gtsv, & ! Partial steps for top cell (ISF)
  127. & rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv , &
  128. & gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi ) ! of t, s, rd at the last ocean level
  129. #endif
  130. !
  131. ! - ML - sshn could be modified by istate_eel, so that initialization of fse3t_b is done here
  132. IF( lk_vvl ) THEN
  133. DO jk = 1, jpk
  134. fse3t_b(:,:,jk) = fse3t_n(:,:,jk)
  135. ENDDO
  136. ENDIF
  137. !
  138. ENDIF
  139. !
  140. IF( lk_agrif ) THEN ! read free surface arrays in restart file
  141. IF( ln_rstart ) THEN
  142. IF( lk_dynspg_flt ) THEN ! read or initialize the following fields
  143. ! ! gcx, gcxb for agrif_opa_init
  144. IF( sol_oce_alloc() > 0 ) CALL ctl_stop('agrif sol_oce_alloc: allocation of arrays failed')
  145. CALL flt_rst( nit000, 'READ' )
  146. ENDIF
  147. ENDIF ! explicit case not coded yet with AGRIF
  148. ENDIF
  149. !
  150. !
  151. ! Initialize "now" and "before" barotropic velocities:
  152. ! Do it whatever the free surface method, these arrays
  153. ! being eventually used
  154. !
  155. !
  156. un_b(:,:) = 0._wp ; vn_b(:,:) = 0._wp
  157. ub_b(:,:) = 0._wp ; vb_b(:,:) = 0._wp
  158. !
  159. DO jk = 1, jpkm1
  160. DO jj = 1, jpj
  161. DO ji = 1, jpi
  162. un_b(ji,jj) = un_b(ji,jj) + fse3u_n(ji,jj,jk) * un(ji,jj,jk) * umask(ji,jj,jk)
  163. vn_b(ji,jj) = vn_b(ji,jj) + fse3v_n(ji,jj,jk) * vn(ji,jj,jk) * vmask(ji,jj,jk)
  164. !
  165. ub_b(ji,jj) = ub_b(ji,jj) + fse3u_b(ji,jj,jk) * ub(ji,jj,jk) * umask(ji,jj,jk)
  166. vb_b(ji,jj) = vb_b(ji,jj) + fse3v_b(ji,jj,jk) * vb(ji,jj,jk) * vmask(ji,jj,jk)
  167. END DO
  168. END DO
  169. END DO
  170. !
  171. un_b(:,:) = un_b(:,:) * hur (:,:)
  172. vn_b(:,:) = vn_b(:,:) * hvr (:,:)
  173. !
  174. ub_b(:,:) = ub_b(:,:) * hur_b(:,:)
  175. vb_b(:,:) = vb_b(:,:) * hvr_b(:,:)
  176. !
  177. !
  178. IF( nn_timing == 1 ) CALL timing_stop('istate_init')
  179. !
  180. END SUBROUTINE istate_init
  181. SUBROUTINE istate_t_s
  182. !!---------------------------------------------------------------------
  183. !! *** ROUTINE istate_t_s ***
  184. !!
  185. !! ** Purpose : Intialization of the temperature field with an
  186. !! analytical profile or a file (i.e. in EEL configuration)
  187. !!
  188. !! ** Method : - temperature: use Philander analytic profile
  189. !! - salinity : use to a constant value 35.5
  190. !!
  191. !! References : Philander ???
  192. !!----------------------------------------------------------------------
  193. INTEGER :: ji, jj, jk
  194. REAL(wp) :: zsal = 35.50
  195. !!----------------------------------------------------------------------
  196. !
  197. IF(lwp) WRITE(numout,*)
  198. IF(lwp) WRITE(numout,*) 'istate_t_s : Philander s initial temperature profile'
  199. IF(lwp) WRITE(numout,*) '~~~~~~~~~~ and constant salinity (',zsal,' psu)'
  200. !
  201. DO jk = 1, jpk
  202. tsn(:,:,jk,jp_tem) = ( ( ( 7.5 - 0. * ABS( gphit(:,:) )/30. ) * ( 1.-TANH((fsdept(:,:,jk)-80.)/30.) ) &
  203. & + 10. * ( 5000. - fsdept(:,:,jk) ) /5000.) ) * tmask(:,:,jk)
  204. tsb(:,:,jk,jp_tem) = tsn(:,:,jk,jp_tem)
  205. END DO
  206. tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:)
  207. tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
  208. !
  209. END SUBROUTINE istate_t_s
  210. SUBROUTINE istate_eel
  211. !!----------------------------------------------------------------------
  212. !! *** ROUTINE istate_eel ***
  213. !!
  214. !! ** Purpose : Initialization of the dynamics and tracers for EEL R5
  215. !! configuration (channel with or without a topographic bump)
  216. !!
  217. !! ** Method : - set temprature field
  218. !! - set salinity field
  219. !! - set velocity field including horizontal divergence
  220. !! and relative vorticity fields
  221. !!----------------------------------------------------------------------
  222. USE divcur ! hor. divergence & rel. vorticity (div_cur routine)
  223. USE iom
  224. !
  225. INTEGER :: inum ! temporary logical unit
  226. INTEGER :: ji, jj, jk ! dummy loop indices
  227. INTEGER :: ijloc
  228. REAL(wp) :: zh1, zh2, zslope, zcst, zfcor ! temporary scalars
  229. REAL(wp) :: zt1 = 15._wp ! surface temperature value (EEL R5)
  230. REAL(wp) :: zt2 = 5._wp ! bottom temperature value (EEL R5)
  231. REAL(wp) :: zsal = 35.0_wp ! constant salinity (EEL R2, R5 and R6)
  232. REAL(wp) :: zueel = 0.1_wp ! constant uniform zonal velocity (EEL R5)
  233. REAL(wp), DIMENSION(jpiglo,jpjglo) :: zssh ! initial ssh over the global domain
  234. !!----------------------------------------------------------------------
  235. !
  236. SELECT CASE ( jp_cfg )
  237. ! ! ====================
  238. CASE ( 5 ) ! EEL R5 configuration
  239. ! ! ====================
  240. !
  241. ! set temperature field with a linear profile
  242. ! -------------------------------------------
  243. IF(lwp) WRITE(numout,*)
  244. IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: linear temperature profile'
  245. IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
  246. !
  247. zh1 = gdept_1d( 1 )
  248. zh2 = gdept_1d(jpkm1)
  249. !
  250. zslope = ( zt1 - zt2 ) / ( zh1 - zh2 )
  251. zcst = ( zt1 * ( zh1 - zh2) - ( zt1 - zt2 ) * zh1 ) / ( zh1 - zh2 )
  252. !
  253. DO jk = 1, jpk
  254. tsn(:,:,jk,jp_tem) = ( zt2 + zt1 * exp( - fsdept(:,:,jk) / 1000 ) ) * tmask(:,:,jk)
  255. tsb(:,:,jk,jp_tem) = tsn(:,:,jk,jp_tem)
  256. END DO
  257. !
  258. IF(lwp) CALL prizre( tsn(:,:,:,jp_tem), jpi , jpj , jpk , jpj/2 , &
  259. & 1 , jpi , 5 , 1 , jpk , &
  260. & 1 , 1. , numout )
  261. !
  262. ! set salinity field to a constant value
  263. ! --------------------------------------
  264. IF(lwp) WRITE(numout,*)
  265. IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
  266. IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
  267. !
  268. tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:)
  269. tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
  270. !
  271. ! set the dynamics: U,V, hdiv, rot (and ssh if necessary)
  272. ! ----------------
  273. ! Start EEL5 configuration with barotropic geostrophic velocities
  274. ! according the sshb and sshn SSH imposed.
  275. ! we assume a uniform grid (hence the use of e1t(1,1) for delta_y)
  276. ! we use the Coriolis frequency at mid-channel.
  277. ub(:,:,:) = zueel * umask(:,:,:)
  278. un(:,:,:) = ub(:,:,:)
  279. ijloc = mj0(INT(jpjglo-1)/2)
  280. zfcor = ff(1,ijloc)
  281. !
  282. DO jj = 1, jpjglo
  283. zssh(:,jj) = - (FLOAT(jj)- FLOAT(jpjglo-1)/2.)*zueel*e1t(1,1)*zfcor/grav
  284. END DO
  285. !
  286. IF(lwp) THEN
  287. WRITE(numout,*) ' Uniform zonal velocity for EEL R5:',zueel
  288. WRITE(numout,*) ' Geostrophic SSH profile as a function of y:'
  289. WRITE(numout,'(12(1x,f6.2))') zssh(1,:)
  290. ENDIF
  291. !
  292. DO jj = 1, nlcj
  293. DO ji = 1, nlci
  294. sshb(ji,jj) = zssh( mig(ji) , mjg(jj) ) * tmask(ji,jj,1)
  295. END DO
  296. END DO
  297. sshb(nlci+1:jpi, : ) = 0.e0 ! set to zero extra mpp columns
  298. sshb( : ,nlcj+1:jpj) = 0.e0 ! set to zero extra mpp rows
  299. !
  300. sshn(:,:) = sshb(:,:) ! set now ssh to the before value
  301. !
  302. IF( nn_rstssh /= 0 ) THEN
  303. nn_rstssh = 0 ! hand-made initilization of ssh
  304. CALL ctl_warn( 'istate_eel: force nn_rstssh = 0' )
  305. ENDIF
  306. !
  307. CALL div_cur( nit000 ) ! horizontal divergence and relative vorticity (curl)
  308. ! N.B. the vertical velocity will be computed from the horizontal divergence field
  309. ! in istate by a call to wzv routine
  310. ! ! ==========================
  311. CASE ( 2 , 6 ) ! EEL R2 or R6 configuration
  312. ! ! ==========================
  313. !
  314. ! set temperature field with a NetCDF file
  315. ! ----------------------------------------
  316. IF(lwp) WRITE(numout,*)
  317. IF(lwp) WRITE(numout,*) 'istate_eel : EEL R2 or R6: read initial temperature in a NetCDF file'
  318. IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
  319. !
  320. CALL iom_open ( 'eel.initemp', inum )
  321. CALL iom_get ( inum, jpdom_data, 'initemp', tsb(:,:,:,jp_tem) ) ! read before temprature (tb)
  322. CALL iom_close( inum )
  323. !
  324. tsn(:,:,:,jp_tem) = tsb(:,:,:,jp_tem) ! set nox temperature to tb
  325. !
  326. IF(lwp) CALL prizre( tsn(:,:,:,jp_tem), jpi , jpj , jpk , jpj/2 , &
  327. & 1 , jpi , 5 , 1 , jpk , &
  328. & 1 , 1. , numout )
  329. !
  330. ! set salinity field to a constant value
  331. ! --------------------------------------
  332. IF(lwp) WRITE(numout,*)
  333. IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
  334. IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
  335. !
  336. tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:)
  337. tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
  338. !
  339. ! ! ===========================
  340. CASE DEFAULT ! NONE existing configuration
  341. ! ! ===========================
  342. WRITE(ctmp1,*) 'EEL with a ', jp_cfg,' km resolution is not coded'
  343. CALL ctl_stop( ctmp1 )
  344. !
  345. END SELECT
  346. !
  347. END SUBROUTINE istate_eel
  348. SUBROUTINE istate_gyre
  349. !!----------------------------------------------------------------------
  350. !! *** ROUTINE istate_gyre ***
  351. !!
  352. !! ** Purpose : Initialization of the dynamics and tracers for GYRE
  353. !! configuration (double gyre with rotated domain)
  354. !!
  355. !! ** Method : - set temprature field
  356. !! - set salinity field
  357. !!----------------------------------------------------------------------
  358. INTEGER :: ji, jj, jk ! dummy loop indices
  359. INTEGER :: inum ! temporary logical unit
  360. INTEGER, PARAMETER :: ntsinit = 0 ! (0/1) (analytical/input data files) T&S initialization
  361. !!----------------------------------------------------------------------
  362. !
  363. SELECT CASE ( ntsinit)
  364. !
  365. CASE ( 0 ) ! analytical T/S profil deduced from LEVITUS
  366. IF(lwp) WRITE(numout,*)
  367. IF(lwp) WRITE(numout,*) 'istate_gyre : initial analytical T and S profil deduced from LEVITUS '
  368. IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
  369. !
  370. DO jk = 1, jpk
  371. DO jj = 1, jpj
  372. DO ji = 1, jpi
  373. tsn(ji,jj,jk,jp_tem) = ( 16. - 12. * TANH( (fsdept(ji,jj,jk) - 400) / 700 ) ) &
  374. & * (-TANH( (500-fsdept(ji,jj,jk)) / 150 ) + 1) / 2 &
  375. & + ( 15. * ( 1. - TANH( (fsdept(ji,jj,jk)-50.) / 1500.) ) &
  376. & - 1.4 * TANH((fsdept(ji,jj,jk)-100.) / 100.) &
  377. & + 7. * (1500. - fsdept(ji,jj,jk)) / 1500. ) &
  378. & * (-TANH( (fsdept(ji,jj,jk) - 500) / 150) + 1) / 2
  379. tsn(ji,jj,jk,jp_tem) = tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk)
  380. tsb(ji,jj,jk,jp_tem) = tsn(ji,jj,jk,jp_tem)
  381. tsn(ji,jj,jk,jp_sal) = ( 36.25 - 1.13 * TANH( (fsdept(ji,jj,jk) - 305) / 460 ) ) &
  382. & * (-TANH((500 - fsdept(ji,jj,jk)) / 150) + 1) / 2 &
  383. & + ( 35.55 + 1.25 * (5000. - fsdept(ji,jj,jk)) / 5000. &
  384. & - 1.62 * TANH( (fsdept(ji,jj,jk) - 60. ) / 650. ) &
  385. & + 0.2 * TANH( (fsdept(ji,jj,jk) - 35. ) / 100. ) &
  386. & + 0.2 * TANH( (fsdept(ji,jj,jk) - 1000.) / 5000.) ) &
  387. & * (-TANH((fsdept(ji,jj,jk) - 500) / 150) + 1) / 2
  388. tsn(ji,jj,jk,jp_sal) = tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk)
  389. tsb(ji,jj,jk,jp_sal) = tsn(ji,jj,jk,jp_sal)
  390. END DO
  391. END DO
  392. END DO
  393. !
  394. CASE ( 1 ) ! T/S data fields read in dta_tem.nc/data_sal.nc files
  395. IF(lwp) WRITE(numout,*)
  396. IF(lwp) WRITE(numout,*) 'istate_gyre : initial T and S read from dta_tem.nc/data_sal.nc files'
  397. IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
  398. IF(lwp) WRITE(numout,*) ' NetCDF FORMAT'
  399. ! Read temperature field
  400. ! ----------------------
  401. CALL iom_open ( 'data_tem', inum )
  402. CALL iom_get ( inum, jpdom_data, 'votemper', tsn(:,:,:,jp_tem) )
  403. CALL iom_close( inum )
  404. tsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) * tmask(:,:,:)
  405. tsb(:,:,:,jp_tem) = tsn(:,:,:,jp_tem)
  406. ! Read salinity field
  407. ! -------------------
  408. CALL iom_open ( 'data_sal', inum )
  409. CALL iom_get ( inum, jpdom_data, 'vosaline', tsn(:,:,:,jp_sal) )
  410. CALL iom_close( inum )
  411. tsn(:,:,:,jp_sal) = tsn(:,:,:,jp_sal) * tmask(:,:,:)
  412. tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
  413. !
  414. END SELECT
  415. !
  416. IF(lwp) THEN
  417. WRITE(numout,*)
  418. WRITE(numout,*) ' Initial temperature and salinity profiles:'
  419. WRITE(numout, "(9x,' level gdept_1d temperature salinity ')" )
  420. WRITE(numout, "(10x, i4, 3f10.2)" ) ( jk, gdept_1d(jk), tsn(2,2,jk,jp_tem), tsn(2,2,jk,jp_sal), jk = 1, jpk )
  421. ENDIF
  422. !
  423. END SUBROUTINE istate_gyre
  424. SUBROUTINE istate_uvg
  425. !!----------------------------------------------------------------------
  426. !! *** ROUTINE istate_uvg ***
  427. !!
  428. !! ** Purpose : Compute the geostrophic velocities from (tn,sn) fields
  429. !!
  430. !! ** Method : Using the hydrostatic hypothesis the now hydrostatic
  431. !! pressure is computed by integrating the in-situ density from the
  432. !! surface to the bottom.
  433. !! p=integral [ rau*g dz ]
  434. !!----------------------------------------------------------------------
  435. USE dynspg ! surface pressure gradient (dyn_spg routine)
  436. USE divcur ! hor. divergence & rel. vorticity (div_cur routine)
  437. USE lbclnk ! ocean lateral boundary condition (or mpp link)
  438. !
  439. INTEGER :: ji, jj, jk ! dummy loop indices
  440. INTEGER :: indic ! ???
  441. REAL(wp) :: zmsv, zphv, zmsu, zphu, zalfg ! temporary scalars
  442. REAL(wp), POINTER, DIMENSION(:,:,:) :: zprn
  443. !!----------------------------------------------------------------------
  444. !
  445. CALL wrk_alloc( jpi, jpj, jpk, zprn)
  446. !
  447. IF(lwp) WRITE(numout,*)
  448. IF(lwp) WRITE(numout,*) 'istate_uvg : Start from Geostrophy'
  449. IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
  450. ! Compute the now hydrostatic pressure
  451. ! ------------------------------------
  452. zalfg = 0.5 * grav * rau0
  453. zprn(:,:,1) = zalfg * fse3w(:,:,1) * ( 1 + rhd(:,:,1) ) ! Surface value
  454. DO jk = 2, jpkm1 ! Vertical integration from the surface
  455. zprn(:,:,jk) = zprn(:,:,jk-1) &
  456. & + zalfg * fse3w(:,:,jk) * ( 2. + rhd(:,:,jk) + rhd(:,:,jk-1) )
  457. END DO
  458. ! Compute geostrophic balance
  459. ! ---------------------------
  460. DO jk = 1, jpkm1
  461. DO jj = 2, jpjm1
  462. DO ji = fs_2, fs_jpim1 ! vertor opt.
  463. zmsv = 1. / MAX( umask(ji-1,jj+1,jk) + umask(ji ,jj+1,jk) &
  464. + umask(ji-1,jj ,jk) + umask(ji ,jj ,jk) , 1. )
  465. zphv = ( zprn(ji ,jj+1,jk) - zprn(ji-1,jj+1,jk) ) * umask(ji-1,jj+1,jk) / e1u(ji-1,jj+1) &
  466. + ( zprn(ji+1,jj+1,jk) - zprn(ji ,jj+1,jk) ) * umask(ji ,jj+1,jk) / e1u(ji ,jj+1) &
  467. + ( zprn(ji ,jj ,jk) - zprn(ji-1,jj ,jk) ) * umask(ji-1,jj ,jk) / e1u(ji-1,jj ) &
  468. + ( zprn(ji+1,jj ,jk) - zprn(ji ,jj ,jk) ) * umask(ji ,jj ,jk) / e1u(ji ,jj )
  469. zphv = 1. / rau0 * zphv * zmsv * vmask(ji,jj,jk)
  470. zmsu = 1. / MAX( vmask(ji+1,jj ,jk) + vmask(ji ,jj ,jk) &
  471. + vmask(ji+1,jj-1,jk) + vmask(ji ,jj-1,jk) , 1. )
  472. zphu = ( zprn(ji+1,jj+1,jk) - zprn(ji+1,jj ,jk) ) * vmask(ji+1,jj ,jk) / e2v(ji+1,jj ) &
  473. + ( zprn(ji ,jj+1,jk) - zprn(ji ,jj ,jk) ) * vmask(ji ,jj ,jk) / e2v(ji ,jj ) &
  474. + ( zprn(ji+1,jj ,jk) - zprn(ji+1,jj-1,jk) ) * vmask(ji+1,jj-1,jk) / e2v(ji+1,jj-1) &
  475. + ( zprn(ji ,jj ,jk) - zprn(ji ,jj-1,jk) ) * vmask(ji ,jj-1,jk) / e2v(ji ,jj-1)
  476. zphu = 1. / rau0 * zphu * zmsu * umask(ji,jj,jk)
  477. ! Compute the geostrophic velocities
  478. un(ji,jj,jk) = -2. * zphu / ( ff(ji,jj) + ff(ji ,jj-1) )
  479. vn(ji,jj,jk) = 2. * zphv / ( ff(ji,jj) + ff(ji-1,jj ) )
  480. END DO
  481. END DO
  482. END DO
  483. IF(lwp) WRITE(numout,*) ' we force to zero bottom velocity'
  484. ! Susbtract the bottom velocity (level jpk-1 for flat bottom case)
  485. ! to have a zero bottom velocity
  486. DO jk = 1, jpkm1
  487. un(:,:,jk) = ( un(:,:,jk) - un(:,:,jpkm1) ) * umask(:,:,jk)
  488. vn(:,:,jk) = ( vn(:,:,jk) - vn(:,:,jpkm1) ) * vmask(:,:,jk)
  489. END DO
  490. CALL lbc_lnk( un, 'U', -1. )
  491. CALL lbc_lnk( vn, 'V', -1. )
  492. ub(:,:,:) = un(:,:,:)
  493. vb(:,:,:) = vn(:,:,:)
  494. ! WARNING !!!!!
  495. ! after initializing u and v, we need to calculate the initial streamfunction bsf.
  496. ! Otherwise, only the trend will be computed and the model will blow up (inconsistency).
  497. ! to do that, we call dyn_spg with a special trick:
  498. ! we fill ua and va with the velocities divided by dt, and the streamfunction will be brought to the
  499. ! right value assuming the velocities have been set up in one time step.
  500. ! we then set bsfd to zero (first guess for next step is d(psi)/dt = 0.)
  501. ! sets up s false trend to calculate the barotropic streamfunction.
  502. ua(:,:,:) = ub(:,:,:) / rdt
  503. va(:,:,:) = vb(:,:,:) / rdt
  504. ! calls dyn_spg. we assume euler time step, starting from rest.
  505. indic = 0
  506. CALL dyn_spg( nit000, indic ) ! surface pressure gradient
  507. ! the new velocity is ua*rdt
  508. CALL lbc_lnk( ua, 'U', -1. )
  509. CALL lbc_lnk( va, 'V', -1. )
  510. ub(:,:,:) = ua(:,:,:) * rdt
  511. vb(:,:,:) = va(:,:,:) * rdt
  512. ua(:,:,:) = 0.e0
  513. va(:,:,:) = 0.e0
  514. un(:,:,:) = ub(:,:,:)
  515. vn(:,:,:) = vb(:,:,:)
  516. ! Compute the divergence and curl
  517. CALL div_cur( nit000 ) ! now horizontal divergence and curl
  518. hdivb(:,:,:) = hdivn(:,:,:) ! set the before to the now value
  519. rotb (:,:,:) = rotn (:,:,:) ! set the before to the now value
  520. !
  521. CALL wrk_dealloc( jpi, jpj, jpk, zprn)
  522. !
  523. END SUBROUTINE istate_uvg
  524. !!=====================================================================
  525. END MODULE istate