- Timestamp:
- 2019-11-22T15:29:17+01:00 (4 years ago)
- Location:
- NEMO/branches/2019/dev_r11943_MERGE_2019/src
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
NEMO/branches/2019/dev_r11943_MERGE_2019/src
- Property svn:mergeinfo deleted
-
NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/step.F90
r11536 r11949 31 31 !! - ! 2015-11 (J. Chanut) free surface simplification (remove filtered free surface) 32 32 !! 4.0 ! 2017-05 (G. Madec) introduction of the vertical physics manager (zdfphy) 33 !! 4.1 ! 2019-08 (A. Coward, D. Storkey) rewrite in preparation for new timestepping scheme 33 34 !!---------------------------------------------------------------------- 34 35 … … 44 45 45 46 PUBLIC stp ! called by nemogcm.F90 47 48 !!---------------------------------------------------------------------- 49 !! time level indices 50 !!---------------------------------------------------------------------- 51 INTEGER, PUBLIC :: Nbb, Nnn, Naa, Nrhs !! used by nemo_init 46 52 47 53 !!---------------------------------------------------------------------- … … 82 88 #if defined key_agrif 83 89 kstp = nit000 + Agrif_Nb_Step() 90 Kbb_a = Nbb; Kmm_a = Nnn; Krhs_a = Nrhs ! agrif_oce module copies of time level indices 84 91 IF( lk_agrif_debug ) THEN 85 92 IF( Agrif_Root() .and. lwp) WRITE(*,*) '---' … … 110 117 ! Update external forcing (tides, open boundaries, and surface boundary condition (including sea-ice) 111 118 !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 112 IF( ln_tide ) CALL sbc_tide( kstp ) ! update tide potential113 IF( ln_apr_dyn ) CALL sbc_apr ( kstp ) ! atmospheric pressure (NB: call before bdy_dta which needs ssh_ib)114 IF( ln_bdy ) CALL bdy_dta ( kstp, kt_offset = +1 ) ! update dynamic & tracer data at open boundaries115 CALL sbc ( kstp )! Sea Boundary Condition (including sea-ice)119 IF( ln_tide ) CALL sbc_tide( kstp ) ! update tide potential 120 IF( ln_apr_dyn ) CALL sbc_apr ( kstp ) ! atmospheric pressure (NB: call before bdy_dta which needs ssh_ib) 121 IF( ln_bdy ) CALL bdy_dta ( kstp, Nnn, kt_offset = +1 ) ! update dynamic & tracer data at open boundaries 122 CALL sbc ( kstp, Nbb, Nnn ) ! Sea Boundary Condition (including sea-ice) 116 123 117 124 !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> … … 119 126 !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 120 127 IF( ln_sto_eos ) CALL sto_par( kstp ) ! Stochastic parameters 121 IF( ln_sto_eos ) CALL sto_pts( ts n) ! Random T/S fluctuations128 IF( ln_sto_eos ) CALL sto_pts( ts(:,:,:,:,Nnn) ) ! Random T/S fluctuations 122 129 123 130 !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> … … 125 132 !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 126 133 ! THERMODYNAMICS 127 CALL eos_rab( ts b, rab_b) ! before local thermal/haline expension ratio at T-points128 CALL eos_rab( ts n, rab_n ) ! now local thermal/haline expension ratio at T-points129 CALL bn2 ( ts b, rab_b, rn2b) ! before Brunt-Vaisala frequency130 CALL bn2 ( ts n, rab_n, rn2) ! now Brunt-Vaisala frequency134 CALL eos_rab( ts(:,:,:,:,Nbb), rab_b, Nnn ) ! before local thermal/haline expension ratio at T-points 135 CALL eos_rab( ts(:,:,:,:,Nnn), rab_n, Nnn ) ! now local thermal/haline expension ratio at T-points 136 CALL bn2 ( ts(:,:,:,:,Nbb), rab_b, rn2b, Nnn ) ! before Brunt-Vaisala frequency 137 CALL bn2 ( ts(:,:,:,:,Nnn), rab_n, rn2, Nnn ) ! now Brunt-Vaisala frequency 131 138 132 139 ! VERTICAL PHYSICS 133 CALL zdf_phy( kstp )! vertical physics update (top/bot drag, avt, avs, avm + MLD)140 CALL zdf_phy( kstp, Nbb, Nnn, Nrhs ) ! vertical physics update (top/bot drag, avt, avs, avm + MLD) 134 141 135 142 ! LATERAL PHYSICS 136 143 ! 137 144 IF( l_ldfslp ) THEN ! slope of lateral mixing 138 CALL eos( ts b, rhd, gdept_0(:,:,:) ) ! before in situ density139 140 IF( ln_zps .AND. .NOT. ln_isfcav) &141 & CALL zps_hde ( kstp, jpts, tsb, gtsu, gtsv, & ! Partial steps: before horizontal gradient142 & rhd, gru , grv ) ! of t, s, rd at the last ocean level143 144 IF( ln_zps .AND. ln_isfcav) &145 & CALL zps_hde_isf( kstp, jpts, tsb, gtsu, gtsv, gtui, gtvi, & ! Partial steps for top cell (ISF)146 & rhd, gru , grv , grui, grvi ) ! of t, s, rd at the first ocean level145 CALL eos( ts(:,:,:,:,Nbb), rhd, gdept_0(:,:,:) ) ! before in situ density 146 147 IF( ln_zps .AND. .NOT. ln_isfcav) & 148 & CALL zps_hde ( kstp, Nnn, jpts, ts(:,:,:,:,Nbb), gtsu, gtsv, & ! Partial steps: before horizontal gradient 149 & rhd, gru , grv ) ! of t, s, rd at the last ocean level 150 151 IF( ln_zps .AND. ln_isfcav) & 152 & CALL zps_hde_isf( kstp, Nnn, jpts, ts(:,:,:,:,Nbb), gtsu, gtsv, gtui, gtvi, & ! Partial steps for top cell (ISF) 153 & rhd, gru , grv , grui, grvi ) ! of t, s, rd at the first ocean level 147 154 IF( ln_traldf_triad ) THEN 148 CALL ldf_slp_triad( kstp )! before slope for triad operator155 CALL ldf_slp_triad( kstp, Nbb, Nnn ) ! before slope for triad operator 149 156 ELSE 150 CALL ldf_slp ( kstp, rhd, rn2b )! before slope for standard operator157 CALL ldf_slp ( kstp, rhd, rn2b, Nbb, Nnn ) ! before slope for standard operator 151 158 ENDIF 152 159 ENDIF 153 ! ! eddy diffusivity coeff.154 IF( l_ldftra_time .OR. l_ldfeiv_time ) CALL ldf_tra( kstp )! and/or eiv coeff.155 IF( l_ldfdyn_time ) CALL ldf_dyn( kstp ) ! eddy viscosity coeff.160 ! ! eddy diffusivity coeff. 161 IF( l_ldftra_time .OR. l_ldfeiv_time ) CALL ldf_tra( kstp, Nbb, Nnn ) ! and/or eiv coeff. 162 IF( l_ldfdyn_time ) CALL ldf_dyn( kstp, Nbb ) ! eddy viscosity coeff. 156 163 157 164 !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> … … 159 166 !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 160 167 161 CALL ssh_nxt ( kstp )! after ssh (includes call to div_hor)162 IF( .NOT.ln_linssh ) CALL dom_vvl_sf_nxt( kstp )! after vertical scale factors163 CALL wzv ( kstp )! now cross-level velocity164 IF( ln_zad_Aimp ) CALL wAimp ( kstp ) ! Adaptive-implicit vertical advection partitioning165 CALL eos ( ts n, rhd, rhop, gdept_n(:,:,:) ) ! now in situ density for hpg computation168 CALL ssh_nxt ( kstp, Nbb, Nnn, ssh, Naa ) ! after ssh (includes call to div_hor) 169 IF( .NOT.ln_linssh ) CALL dom_vvl_sf_nxt( kstp, Nbb, Nnn, Naa ) ! after vertical scale factors 170 CALL wzv ( kstp, Nbb, Nnn, ww, Naa ) ! now cross-level velocity 171 IF( ln_zad_Aimp ) CALL wAimp ( kstp, Nnn ) ! Adaptive-implicit vertical advection partitioning 172 CALL eos ( ts(:,:,:,:,Nnn), rhd, rhop, gdept(:,:,:,Nnn) ) ! now in situ density for hpg computation 166 173 167 174 168 u a(:,:,:) = 0._wp ! set dynamics trends to zero169 v a(:,:,:) = 0._wp175 uu(:,:,:,Nrhs) = 0._wp ! set dynamics trends to zero 176 vv(:,:,:,Nrhs) = 0._wp 170 177 171 178 IF( lk_asminc .AND. ln_asmiau .AND. ln_dyninc ) & 172 & CALL dyn_asm_inc ( kstp ) ! apply dynamics assimilation increment173 IF( ln_bdy ) CALL bdy_dyn3d_dmp ( kstp ) ! bdy damping trends179 & CALL dyn_asm_inc ( kstp, Nbb, Nnn, uu, vv, Nrhs ) ! apply dynamics assimilation increment 180 IF( ln_bdy ) CALL bdy_dyn3d_dmp ( kstp, Nbb, uu, vv, Nrhs ) ! bdy damping trends 174 181 #if defined key_agrif 175 182 IF(.NOT. Agrif_Root()) & 176 183 & CALL Agrif_Sponge_dyn ! momentum sponge 177 184 #endif 178 CALL dyn_adv ( kstp ) ! advection (vector or flux form)179 CALL dyn_vor ( kstp ) ! vorticity term including Coriolis180 CALL dyn_ldf ( kstp) ! lateral mixing181 IF( ln_zdfosm ) CALL dyn_osm ( kstp ) ! OSMOSIS non-local velocity fluxes182 CALL dyn_hpg ( kstp) ! horizontal gradient of Hydrostatic pressure183 CALL dyn_spg ( kstp) ! surface pressure gradient184 185 ! With split-explicit free surface, since now transports have been updated and ssh aas well185 CALL dyn_adv( kstp, Nbb, Nnn , uu, vv, Nrhs ) ! advection (VF or FF) ==> RHS 186 CALL dyn_vor( kstp, Nnn , uu, vv, Nrhs ) ! vorticity ==> RHS 187 CALL dyn_ldf( kstp, Nbb, Nnn , uu, vv, Nrhs ) ! lateral mixing 188 IF( ln_zdfosm ) CALL dyn_osm( kstp, Nnn , uu, vv, Nrhs ) ! OSMOSIS non-local velocity fluxes ==> RHS 189 CALL dyn_hpg( kstp, Nnn , uu, vv, Nrhs ) ! horizontal gradient of Hydrostatic pressure 190 CALL dyn_spg( kstp, Nbb, Nnn, Nrhs, uu, vv, ssh, uu_b, vv_b, Naa ) ! surface pressure gradient 191 192 ! With split-explicit free surface, since now transports have been updated and ssh(:,:,Nrhs) as well 186 193 IF( ln_dynspg_ts ) THEN ! vertical scale factors and vertical velocity need to be updated 187 CALL div_hor ( kstp )! Horizontal divergence (2nd call in time-split case)188 IF(.NOT.ln_linssh) CALL dom_vvl_sf_nxt( kstp, kcall=2 ) ! after vertical scale factors (update depth average component)189 ENDIF 190 CALL dyn_zdf ( kstp ) ! vertical diffusion194 CALL div_hor ( kstp, Nbb, Nnn ) ! Horizontal divergence (2nd call in time-split case) 195 IF(.NOT.ln_linssh) CALL dom_vvl_sf_nxt( kstp, Nbb, Nnn, Naa, kcall=2 ) ! after vertical scale factors (update depth average component) 196 ENDIF 197 CALL dyn_zdf( kstp, Nbb, Nnn, Nrhs, uu, vv, Naa ) ! vertical diffusion ==> after 191 198 192 199 IF( ln_dynspg_ts ) THEN 193 CALL wzv ( kstp ) ! now cross-level velocity 194 IF( ln_zad_Aimp ) CALL wAimp ( kstp ) ! Adaptive-implicit vertical advection partitioning 195 ENDIF 200 CALL wzv ( kstp, Nbb, Nnn, ww, Naa ) ! now cross-level velocity 201 IF( ln_zad_Aimp ) CALL wAimp ( kstp, Nnn ) ! Adaptive-implicit vertical advection partitioning 202 ENDIF 203 196 204 197 205 !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 198 206 ! cool skin 199 207 !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 200 IF ( ln_diurnal ) CALL stp_diurnal( kstp )208 IF ( ln_diurnal ) CALL diurnal_layers( kstp ) 201 209 202 210 !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 203 211 ! diagnostics and outputs 204 212 !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 205 IF( ln_floats ) CALL flo_stp ( kstp )! drifting Floats206 IF( ln_diacfl ) CALL dia_cfl ( kstp )! Courant number diagnostics207 IF( lk_diahth ) CALL dia_hth ( kstp )! Thermocline depth (20 degres isotherm depth)208 IF( ln_diadct ) CALL dia_dct ( kstp )! Transports209 CALL dia_ar5 ( kstp )! ar5 diag210 IF( ln_diaharm ) CALL dia_harm( kstp )! Tidal harmonic analysis211 CALL dia_wri ( kstp )! ocean model: outputs212 ! 213 IF( ln_crs ) CALL crs_fld ( kstp )! ocean model: online field coarsening & output213 IF( ln_floats ) CALL flo_stp ( kstp, Nbb, Nnn ) ! drifting Floats 214 IF( ln_diacfl ) CALL dia_cfl ( kstp, Nnn ) ! Courant number diagnostics 215 IF( lk_diahth ) CALL dia_hth ( kstp, Nnn ) ! Thermocline depth (20 degres isotherm depth) 216 IF( ln_diadct ) CALL dia_dct ( kstp, Nnn ) ! Transports 217 CALL dia_ar5 ( kstp, Nnn ) ! ar5 diag 218 IF( ln_diaharm ) CALL dia_harm( kstp, Nnn ) ! Tidal harmonic analysis 219 CALL dia_wri ( kstp, Nnn ) ! ocean model: outputs 220 ! 221 IF( ln_crs ) CALL crs_fld ( kstp, Nnn ) ! ocean model: online field coarsening & output 214 222 215 223 #if defined key_top … … 217 225 ! Passive Tracer Model 218 226 !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 219 CALL trc_stp ( kstp ) ! time-stepping227 CALL trc_stp ( kstp, Nbb, Nnn, Nrhs, Naa ) ! time-stepping 220 228 #endif 221 229 … … 223 231 ! Active tracers 224 232 !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 225 ts a(:,:,:,:) = 0._wp ! set tracer trends to zero233 ts(:,:,:,:,Nrhs) = 0._wp ! set tracer trends to zero 226 234 227 235 IF( lk_asminc .AND. ln_asmiau .AND. & 228 & ln_trainc ) CALL tra_asm_inc ( kstp) ! apply tracer assimilation increment229 CALL tra_sbc ( kstp) ! surface boundary condition230 IF( ln_traqsr ) CALL tra_qsr ( kstp) ! penetrative solar radiation qsr231 IF( ln_trabbc ) CALL tra_bbc ( kstp) ! bottom heat flux232 IF( ln_trabbl ) CALL tra_bbl ( kstp) ! advective (and/or diffusive) bottom boundary layer scheme233 IF( ln_tradmp ) CALL tra_dmp ( kstp) ! internal damping trends234 IF( ln_bdy ) CALL bdy_tra_dmp ( kstp) ! bdy damping trends236 & ln_trainc ) CALL tra_asm_inc( kstp, Nbb, Nnn, ts, Nrhs ) ! apply tracer assimilation increment 237 CALL tra_sbc ( kstp, Nnn, ts, Nrhs ) ! surface boundary condition 238 IF( ln_traqsr ) CALL tra_qsr ( kstp, Nnn, ts, Nrhs ) ! penetrative solar radiation qsr 239 IF( ln_trabbc ) CALL tra_bbc ( kstp, Nnn, ts, Nrhs ) ! bottom heat flux 240 IF( ln_trabbl ) CALL tra_bbl ( kstp, Nbb, Nnn, ts, Nrhs ) ! advective (and/or diffusive) bottom boundary layer scheme 241 IF( ln_tradmp ) CALL tra_dmp ( kstp, Nbb, Nnn, ts, Nrhs ) ! internal damping trends 242 IF( ln_bdy ) CALL bdy_tra_dmp( kstp, Nbb, ts, Nrhs ) ! bdy damping trends 235 243 #if defined key_agrif 236 244 IF(.NOT. Agrif_Root()) & 237 245 & CALL Agrif_Sponge_tra ! tracers sponge 238 246 #endif 239 CALL tra_adv ( kstp ) ! horizontal & vertical advection240 IF( ln_zdfosm ) CALL tra_osm ( kstp ) ! OSMOSIS non-local tracer fluxes247 CALL tra_adv ( kstp, Nbb, Nnn, ts, Nrhs ) ! hor. + vert. advection ==> RHS 248 IF( ln_zdfosm ) CALL tra_osm ( kstp, Nnn, ts, Nrhs ) ! OSMOSIS non-local tracer fluxes ==> RHS 241 249 IF( lrst_oce .AND. ln_zdfosm ) & 242 & CALL osm_rst ( kstp, 'WRITE' )! write OSMOSIS outputs + wn(so must do here) to restarts243 CALL tra_ldf ( kstp) ! lateral mixing250 & CALL osm_rst ( kstp, Nnn, 'WRITE' ) ! write OSMOSIS outputs + ww (so must do here) to restarts 251 CALL tra_ldf ( kstp, Nbb, Nnn, ts, Nrhs ) ! lateral mixing 244 252 245 253 !!gm : why CALL to dia_ptr has been moved here??? (use trends info?) 246 IF( ln_diaptr ) CALL dia_ptr ! Poleward adv/ldf TRansports diagnostics254 IF( ln_diaptr ) CALL dia_ptr( Nnn ) ! Poleward adv/ldf TRansports diagnostics 247 255 !!gm 248 CALL tra_zdf ( kstp ) ! vertical mixing and after tracer fields249 IF( ln_zdfnpc ) CALL tra_npc ( kstp) ! update after fields by non-penetrative convection250 251 !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 252 ! Set boundary conditions and Swap256 CALL tra_zdf ( kstp, Nbb, Nnn, Nrhs, ts, Naa ) ! vert. mixing & after tracer ==> after 257 IF( ln_zdfnpc ) CALL tra_npc ( kstp, Nnn, Nrhs, ts, Naa ) ! update after fields by non-penetrative convection 258 259 !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 260 ! Set boundary conditions, time filter and swap time levels 253 261 !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 254 262 !!jc1: For agrif, it would be much better to finalize tracers/momentum here (e.g. bdy conditions) and move the swap … … 265 273 !! place. 266 274 !! 267 !!jc2: dynnxt must be the latest call. e3t_b are indeed updated in that routine 268 CALL tra_nxt ( kstp ) ! finalize (bcs) tracer fields at next time step and swap 269 CALL dyn_nxt ( kstp ) ! finalize (bcs) velocities at next time step and swap (always called after tra_nxt) 270 CALL ssh_swp ( kstp ) ! swap of sea surface height 271 IF(.NOT.ln_linssh) CALL dom_vvl_sf_swp( kstp ) ! swap of vertical scale factors 272 ! 273 IF( ln_diahsb ) CALL dia_hsb ( kstp ) ! - ML - global conservation diagnostics 275 !!jc2: dynnxt must be the latest call. e3t(:,:,:,Nbb) are indeed updated in that routine 276 CALL tra_atf ( kstp, Nbb, Nnn, Naa, ts ) ! time filtering of "now" tracer arrays 277 CALL dyn_atf ( kstp, Nbb, Nnn, Naa, uu, vv, e3t, e3u, e3v ) ! time filtering of "now" velocities and scale factors 278 CALL ssh_atf ( kstp, Nbb, Nnn, Naa, ssh ) ! time filtering of "now" sea surface height 279 ! 280 ! Swap time levels 281 Nrhs = Nbb 282 Nbb = Nnn 283 Nnn = Naa 284 Naa = Nrhs 285 ! 286 IF(.NOT.ln_linssh) CALL dom_vvl_sf_update( kstp, Nbb, Nnn, Naa ) ! recompute vertical scale factors 287 ! 288 IF( ln_diahsb ) CALL dia_hsb ( kstp, Nbb, Nnn ) ! - ML - global conservation diagnostics 274 289 275 290 !!gm : This does not only concern the dynamics ==>>> add a new title … … 278 293 !!jc: That would be better, but see comment above 279 294 !! 280 IF( lrst_oce ) CALL rst_write ( kstp ) ! write output ocean restart file295 IF( lrst_oce ) CALL rst_write ( kstp, Nbb, Nnn ) ! write output ocean restart file 281 296 IF( ln_sto_eos ) CALL sto_rst_write( kstp ) ! write restart file for stochastic parameters 282 297 … … 285 300 ! AGRIF 286 301 !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 287 CALL Agrif_Integrate_ChildGrids( stp ) ! allows to finish all the Child Grids before updating 288 289 IF( Agrif_NbStepint() == 0 ) CALL Agrif_update_all( ) ! Update all components 290 #endif 291 IF( ln_diaobs ) CALL dia_obs ( kstp ) ! obs-minus-model (assimilation) diagnostics (call after dynamics update) 302 Kbb_a = Nbb; Kmm_a = Nnn; Krhs_a = Nrhs ! agrif_oce module copies of time level indices 303 CALL Agrif_Integrate_ChildGrids( stp ) ! allows to finish all the Child Grids before updating 304 305 IF( Agrif_NbStepint() == 0 ) THEN 306 CALL Agrif_update_all( ) ! Update all components 307 ENDIF 308 #endif 309 IF( ln_diaobs ) CALL dia_obs ( kstp, Nnn ) ! obs-minus-model (assimilation) diagnostics (call after dynamics update) 292 310 293 311 !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 294 312 ! Control 295 313 !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 296 CALL stp_ctl ( kstp, indic )314 CALL stp_ctl ( kstp, Nbb, Nnn, indic ) 297 315 298 316 IF( kstp == nit000 ) THEN ! 1st time step only … … 306 324 !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 307 325 !!gm why lk_oasis and not lk_cpl ???? 308 IF( lk_oasis ) CALL sbc_cpl_snd( kstp )! coupled mode : field exchanges326 IF( lk_oasis ) CALL sbc_cpl_snd( kstp, Nbb, Nnn ) ! coupled mode : field exchanges 309 327 ! 310 328 #if defined key_iomput … … 319 337 ! 320 338 END SUBROUTINE stp 321 339 ! 322 340 !!====================================================================== 323 341 END MODULE step
Note: See TracChangeset
for help on using the changeset viewer.