6 use mpp_domains_mod
, only : center, corner, north, east
7 use data_override_mod
, only : data_override_init, data_override
14 use mom_io, only : read_data
17 use astronomy_mod
, only : orbital_time, diurnal_solar, daily_mean_solar
38 #include "MOM_memory.h" 39 #include "version_variable.h" 48 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(in) :: uhtr
49 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(in) :: vhtr
50 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) ,
intent(in) :: h_pre
51 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) ,
intent(inout) :: h_new
54 integer :: i, j, k, m, is, ie, js, je, nz
56 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
59 do i=is-1,ie+1 ;
do j=js-1,je+1
61 h_new(i,j,k) = max(0.0, g%areaT(i,j)*h_pre(i,j,k) + &
62 ((uhtr(i-1,j,k) - uhtr(i,j,k)) + (vhtr(i,j-1,k) - vhtr(i,j,k))))
67 h_new(i,j,k) = h_new(i,j,k) + &
68 max(gv%Angstrom, 1.0e-13*h_new(i,j,k) - g%areaT(i,j)*h_pre(i,j,k))
71 h_new(i,j,k) = h_new(i,j,k)/g%areaT(i,j)
82 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) ,
intent(in) :: ea
83 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) ,
intent(in) :: eb
84 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) ,
intent(in) :: h_pre
85 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) ,
intent(inout) :: h_new
88 integer :: i, j, k, m, is, ie, js, je, nz
90 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
97 h_new(i,j,1) = max(0.0, h_pre(i,j,1) + (eb(i,j,1) - ea(i,j,2) + ea(i,j,1) ))
98 h_new(i,j,1) = h_new(i,j,1) + &
99 max(0.0, 1.0e-13*h_new(i,j,1) - h_pre(i,j,1))
103 h_new(i,j,nz) = max(0.0, h_pre(i,j,nz) + (ea(i,j,nz) - eb(i,j,nz-1)+eb(i,j,nz)))
104 h_new(i,j,nz) = h_new(i,j,nz) + &
105 max(0.0, 1.0e-13*h_new(i,j,nz) - h_pre(i,j,nz))
110 do k=2,nz-1 ;
do i=is-1,ie+1
112 h_new(i,j,k) = max(0.0, h_pre(i,j,k) + ((ea(i,j,k) - eb(i,j,k-1)) + &
113 (eb(i,j,k) - ea(i,j,k+1))))
114 h_new(i,j,k) = h_new(i,j,k) + &
115 max(0.0, 1.0e-13*h_new(i,j,k) - h_pre(i,j,k))
128 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: uh
129 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: vh
130 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) ,
intent(inout) :: ea
131 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) ,
intent(inout) :: eb
132 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) ,
intent(in) :: h_pre
135 integer :: i, j, k, m, is, ie, js, je, nz
136 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: top_flux, bottom_flux
137 real :: pos_flux, hvol, h_neglect, scale_factor, max_off_cfl
146 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
151 do j=js-1,je+1 ;
do i=is-1,ie+1
152 top_flux(i,j,k) = -ea(i,j,k)
153 bottom_flux(i,j,k) = -(eb(i,j,k)-ea(i,j,k+1))
156 do k=2, nz-1 ;
do j=js-1,je+1 ;
do i=is-1,ie+1
157 top_flux(i,j,k) = -(ea(i,j,k)-eb(i,j,k-1))
158 bottom_flux(i,j,k) = -(eb(i,j,k)-ea(i,j,k+1))
159 enddo ;
enddo ;
enddo 162 do j=js-1,je+1 ;
do i=is-1,ie+1
163 top_flux(i,j,k) = -(ea(i,j,k)-eb(i,j,k-1))
164 bottom_flux(i,j,k) = -eb(i,j,k)
170 do k = 1, nz ;
do j=js-1,je+1 ;
do i=is-1,ie+1
172 hvol = h_pre(i,j,k)*g%areaT(i,j)
173 pos_flux = max(0.0,-uh(i-1,j,k)) + max(0.0, -vh(i,j-1,k)) + &
174 max(0.0, uh(i,j,k)) + max(0.0, vh(i,j,k)) + &
175 max(0.0, top_flux(i,j,k)*g%areaT(i,j)) + max(0.0, bottom_flux(i,j,k)*g%areaT(i,j))
177 if (pos_flux>hvol .and. pos_flux>0.0)
then 178 scale_factor = ( hvol )/pos_flux*max_off_cfl
184 if (-uh(i-1,j,k)>0) uh(i-1,j,k) = uh(i-1,j,k)*scale_factor
185 if (uh(i,j,k)>0) uh(i,j,k) = uh(i,j,k)*scale_factor
186 if (-vh(i,j-1,k)>0) vh(i,j-1,k) = vh(i,j-1,k)*scale_factor
187 if (vh(i,j,k)>0) vh(i,j,k) = vh(i,j,k)*scale_factor
189 if (k>1 .and. k<nz)
then 191 if (top_flux(i,j,k)>0.0)
then 192 ea(i,j,k) = ea(i,j,k)*scale_factor
193 eb(i,j,k-1) = eb(i,j,k-1)*scale_factor
195 if (bottom_flux(i,j,k)>0.0)
then 196 eb(i,j,k) = eb(i,j,k)*scale_factor
197 ea(i,j,k+1) = ea(i,j,k+1)*scale_factor
201 if (top_flux(i,j,k)>0.0) ea(i,j,k) = ea(i,j,k)*scale_factor
202 if (bottom_flux(i,j,k)>0.0)
then 203 eb(i,j,k) = eb(i,j,k)*scale_factor
204 ea(i,j,k+1) = ea(i,j,k+1)*scale_factor
208 if (top_flux(i,j,k)>0.0)
then 209 ea(i,j,k) = ea(i,j,k)*scale_factor
210 eb(i,j,k-1) = eb(i,j,k-1)*scale_factor
212 if (bottom_flux(i,j,k)>0.0) eb(i,j,k)=eb(i,j,k)*scale_factor
214 enddo ;
enddo ;
enddo 223 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in ) :: hvol
224 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: uh
226 real,
dimension(SZIB_(G),SZK_(G)) :: uh2d
227 real,
dimension(SZIB_(G)) :: uh2d_sum
228 real,
dimension(SZI_(G),SZK_(G)) :: h2d
229 real,
dimension(SZI_(G)) :: h2d_sum
231 integer :: i, j, k, m, is, ie, js, je, nz
235 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
240 do k=1,nz ;
do i=is-1,ie
241 uh2d(i,k) = uh(i,j,k)
242 uh2d_sum(i) = uh2d_sum(i) + uh2d(i,k)
247 do k=1,nz ;
do i=is-1,ie+1
248 h2d(i,k) = hvol(i,j,k)
249 if (hvol(i,j,k)>0.)
then 250 h2d_sum(i) = h2d_sum(i) + h2d(i,k)
252 h2d(i,k) = gv%H_subroundoff
259 if ( uh2d_sum(i)>0.0 )
then 261 uh2d(i,k) = uh2d_sum(i)*(h2d(i,k)/h2d_sum(i))
263 elseif (uh2d_sum(i)<0.0)
then 265 uh2d(i,k) = uh2d_sum(i)*(h2d(i+1,k)/h2d_sum(i+1))
274 uh_neglect = gv%Angstrom*min(g%areaT(i,j),g%areaT(i+1,j))
275 if ( abs(sum(uh2d(i,:))-uh2d_sum(i)) > uh_neglect) &
276 call mom_error(warning,
"Column integral of uh does not match after "//&
277 "barotropic redistribution")
280 do k=1,nz ;
do i=is-1,ie
281 uh(i,j,k) = uh2d(i,k)
289 type(ocean_grid_type),
pointer :: G
290 type(verticalgrid_type),
pointer :: GV
291 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in ) :: hvol
292 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: vh
294 real,
dimension(SZJB_(G),SZK_(G)) :: vh2d
295 real,
dimension(SZJB_(G)) :: vh2d_sum
296 real,
dimension(SZJ_(G),SZK_(G)) :: h2d
297 real,
dimension(SZJ_(G)) :: h2d_sum
299 integer :: i, j, k, m, is, ie, js, je, nz
303 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
308 do k=1,nz ;
do j=js-1,je
309 vh2d(j,k) = vh(i,j,k)
310 vh2d_sum(j) = vh2d_sum(j) + vh2d(j,k)
315 do k=1,nz ;
do j=js-1,je+1
316 h2d(j,k) = hvol(i,j,k)
317 if (hvol(i,j,k)>0.)
then 318 h2d_sum(j) = h2d_sum(j) + h2d(j,k)
320 h2d(j,k) = gv%H_subroundoff
326 if ( vh2d_sum(j)>0.0 )
then 328 vh2d(j,k) = vh2d_sum(j)*(h2d(j,k)/h2d_sum(j))
330 elseif (vh2d_sum(j)<0.0)
then 332 vh2d(j,k) = vh2d_sum(j)*(h2d(j+1,k)/h2d_sum(j+1))
341 vh_neglect = gv%Angstrom*min(g%areaT(i,j),g%areaT(i,j+1))
342 if ( abs(sum(vh2d(j,:))-vh2d_sum(j)) > vh_neglect)
then 343 call mom_error(warning,
"Column integral of vh does not match after "//&
344 "barotropic redistribution")
349 do k=1,nz ;
do j=js-1,je
350 vh(i,j,k) = vh2d(j,k)
359 type(ocean_grid_type),
pointer :: G
360 type(verticalgrid_type),
pointer :: GV
361 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in ) :: hvol
362 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: uh
364 real,
dimension(SZIB_(G),SZK_(G)) :: uh2d
365 real,
dimension(SZI_(G),SZK_(G)) :: h2d
367 real :: uh_neglect, uh_remain, uh_add, uh_sum, uh_col, uh_max
368 real :: hup, hdown, hlos, min_h
369 integer :: i, j, k, m, is, ie, js, je, nz, k_rev
372 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
374 min_h = gv%Angstrom*0.1
378 do k=1,nz ;
do i=is-2,ie+1
379 uh2d(i,k) = uh(i,j,k)
381 do k=1,nz ;
do i=is-1,ie+1
383 h2d(i,k) = hvol(i,j,k)-min_h*g%areaT(i,j)
387 uh_col = sum(uh2d(i,:))
389 uh_remain = uh2d(i,k)
391 if (abs(uh_remain)>0.0)
then 393 uh_sum = uh_remain + uh2d(i,k_rev)
396 hlos = max(0.0,uh2d(i+1,k_rev))
397 if ((((hup - hlos) + uh_sum) < 0.0) .and. &
398 ((0.5*hup + uh_sum) < 0.0))
then 399 uh2d(i,k_rev) = min(-0.5*hup,-hup+hlos,0.0)
400 uh_remain = uh_sum - uh2d(i,k_rev)
402 uh2d(i,k_rev) = uh_sum
408 hlos = max(0.0,-uh2d(i-1,k_rev))
409 if ((((hup - hlos) - uh_sum) < 0.0) .and. &
410 ((0.5*hup - uh_sum) < 0.0))
then 411 uh2d(i,k_rev) = max(0.5*hup,hup-hlos,0.0)
412 uh_remain = uh_sum - uh2d(i,k_rev)
414 uh2d(i,k_rev) = uh_sum
422 if (abs(uh_remain)>0.0)
then 424 uh2d(i,k+1) = uh2d(i,k+1) + uh_remain
426 uh2d(i,k) = uh2d(i,k) + uh_remain
427 call mom_error(warning,
"Water column cannot accommodate UH redistribution. Tracer may not be conserved")
434 uh_neglect = gv%Angstrom*min(g%areaT(i,j),g%areaT(i+1,j))
435 if (abs(uh_col - sum(uh2d(i,:)))>uh_neglect)
then 436 call mom_error(warning,
"Column integral of uh does not match after "//&
437 "upwards redistribution")
442 do k=1,nz ;
do i=is-1,ie
443 uh(i,j,k) = uh2d(i,k)
452 type(ocean_grid_type),
pointer :: G
453 type(verticalgrid_type),
pointer :: GV
454 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in ) :: hvol
455 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: vh
457 real,
dimension(SZJB_(G),SZK_(G)) :: vh2d
458 real,
dimension(SZJB_(G)) :: vh2d_sum
459 real,
dimension(SZJ_(G),SZK_(G)) :: h2d
460 real,
dimension(SZJ_(G)) :: h2d_sum
462 real :: vh_neglect, vh_remain, vh_col, vh_sum
463 real :: hup, hlos, min_h
464 integer :: i, j, k, m, is, ie, js, je, nz, k_rev
467 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
469 min_h = 0.1*gv%Angstrom
473 do k=1,nz ;
do j=js-2,je+1
474 vh2d(j,k) = vh(i,j,k)
476 do k=1,nz ;
do j=js-1,je+1
477 h2d(j,k) = hvol(i,j,k)-min_h*g%areaT(i,j)
481 vh_col = sum(vh2d(j,:))
483 vh_remain = vh2d(j,k)
485 if (abs(vh_remain)>0.0)
then 487 vh_sum = vh_remain + vh2d(j,k_rev)
490 hlos = max(0.0,vh2d(j+1,k_rev))
491 if ((((hup - hlos) + vh_sum) < 0.0) .and. &
492 ((0.5*hup + vh_sum) < 0.0))
then 493 vh2d(j,k_rev) = min(-0.5*hup,-hup+hlos,0.0)
494 vh_remain = vh_sum - vh2d(j,k_rev)
496 vh2d(j,k_rev) = vh_sum
502 hlos = max(0.0,-vh2d(j-1,k_rev))
503 if ((((hup - hlos) - vh_sum) < 0.0) .and. &
504 ((0.5*hup - vh_sum) < 0.0))
then 505 vh2d(j,k_rev) = max(0.5*hup,hup-hlos,0.0)
506 vh_remain = vh_sum - vh2d(j,k_rev)
508 vh2d(j,k_rev) = vh_sum
517 if (abs(vh_remain)>0.0)
then 519 vh2d(j,k+1) = vh2d(j,k+1) + vh_remain
521 vh2d(j,k) = vh2d(j,k) + vh_remain
522 call mom_error(warning,
"Water column cannot accommodate VH redistribution. Tracer will not be conserved")
529 vh_neglect = gv%Angstrom*min(g%areaT(i,j),g%areaT(i,j+1))
530 if ( abs(vh_col-sum(vh2d(j,:))) > vh_neglect)
then 531 call mom_error(warning,
"Column integral of vh does not match after "//&
532 "upwards redistribution")
536 do k=1,nz ;
do j=js-1,je
537 vh(i,j,k) = vh2d(j,k)
546 type(forcing),
intent(inout) :: fluxes
547 type(ocean_grid_type),
intent(in) :: G
548 type(time_type),
intent(in) :: Time_start
549 type(time_type),
intent(in) :: Time_end
551 real :: diurnal_factor, time_since_ae, rad
552 real :: fracday_dt, fracday_day
553 real :: cosz_day, cosz_dt, rrsun_day, rrsun_dt
554 type(time_type) :: dt_here
556 integer :: i, j, k, i2, j2, isc, iec, jsc, jec, i_off, j_off
558 isc = g%isc ; iec = g%iec ; jsc = g%jsc ; jec = g%jec
559 i_off = lbound(fluxes%sens,1) - g%isc ; j_off = lbound(fluxes%sens,2) - g%jsc
563 time_since_ae = orbital_time(time_start)
564 dt_here = time_end - time_start
571 do j=jsc,jec ;
do i=isc,iec
578 call diurnal_solar(g%geoLatT(i,j)*rad, g%geoLonT(i,j)*rad, time_start, cosz=cosz_dt, &
579 fracday=fracday_dt, rrsun=rrsun_dt, dt_time=dt_here)
580 call daily_mean_solar (g%geoLatT(i,j)*rad, time_since_ae, cosz_day, fracday_day, rrsun_day)
581 diurnal_factor = cosz_dt*fracday_dt*rrsun_dt / &
582 max(1e-30, cosz_day*fracday_day*rrsun_day)
584 i2 = i+i_off ; j2 = j+j_off
585 fluxes%sw(i2,j2) = fluxes%sw(i2,j2) * diurnal_factor
586 fluxes%sw_vis_dir(i2,j2) = fluxes%sw_vis_dir(i2,j2) * diurnal_factor
587 fluxes%sw_vis_dif (i2,j2) = fluxes%sw_vis_dif (i2,j2) * diurnal_factor
588 fluxes%sw_nir_dir(i2,j2) = fluxes%sw_nir_dir(i2,j2) * diurnal_factor
589 fluxes%sw_nir_dif (i2,j2) = fluxes%sw_nir_dif (i2,j2) * diurnal_factor
597 uhtr, vhtr, temp_mean, salt_mean, mld, Kd, fluxes, ridx_sum, ridx_snap, read_mld, read_sw, &
598 read_ts_uvh, do_ale_in)
600 type(ocean_grid_type),
pointer,
intent(inout) :: G
601 type(verticalgrid_type),
pointer,
intent(in ) :: GV
602 integer,
intent(in ) :: nk_input
603 character(len=*),
intent(in ) :: mean_file
604 character(len=*),
intent(in ) :: sum_file
605 character(len=*),
intent(in ) :: snap_file
606 character(len=*),
intent(in ) :: surf_file
607 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: uhtr
608 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: vhtr
609 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: h_end
610 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: temp_mean
611 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: salt_mean
612 real,
dimension(SZI_(G),SZJ_(G)),
intent(inout) :: mld
613 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(inout) :: Kd
614 type(forcing),
intent(inout) :: fluxes
615 integer,
intent(in ) :: ridx_sum
616 integer,
intent(in ) :: ridx_snap
617 logical,
intent(in ) :: read_mld
618 logical,
intent(in ) :: read_sw
619 logical,
intent(in ) :: read_ts_uvh
620 logical,
optional,
intent(in ) :: do_ale_in
623 integer :: i, j, k, is, ie, js, je, nz
627 if (
present(do_ale_in) ) do_ale = do_ale_in
629 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
632 if (read_ts_uvh)
then 634 temp_mean(:,:,:) = 0.0
635 salt_mean(:,:,:) = 0.0
639 call read_data(sum_file,
'uhtr_sum', uhtr(:,:,1:nk_input),domain=g%Domain%mpp_domain, &
640 timelevel=ridx_sum, position=east)
641 call read_data(sum_file,
'vhtr_sum', vhtr(:,:,1:nk_input), domain=g%Domain%mpp_domain, &
642 timelevel=ridx_sum, position=north)
643 call read_data(snap_file,
'h_end', h_end(:,:,1:nk_input), domain=g%Domain%mpp_domain, &
644 timelevel=ridx_snap,position=center)
645 call read_data(mean_file,
'temp', temp_mean(:,:,1:nk_input), domain=g%Domain%mpp_domain, &
646 timelevel=ridx_sum,position=center)
647 call read_data(mean_file,
'salt', salt_mean(:,:,1:nk_input), domain=g%Domain%mpp_domain, &
648 timelevel=ridx_sum,position=center)
651 do j=js,je ;
do i=is,ie
652 if (g%mask2dT(i,j)>0.)
then 653 temp_mean(:,:,nk_input:nz) = temp_mean(i,j,nk_input)
654 salt_mean(:,:,nk_input:nz) = salt_mean(i,j,nk_input)
659 call read_data( mean_file,
'Kd_interface', kd(:,:,1:nk_input+1), domain=g%Domain%mpp_domain, &
660 timelevel=ridx_sum,position=center)
665 if (.not.
ASSOCIATED(fluxes%netMassOut))
then 666 allocate(fluxes%netMassOut(g%isd:g%ied,g%jsd:g%jed))
667 fluxes%netMassOut(:,:) = 0.0
669 if (.not.
ASSOCIATED(fluxes%netMassIn))
then 670 allocate(fluxes%netMassIn(g%isd:g%ied,g%jsd:g%jed))
671 fluxes%netMassIn(:,:) = 0.0
674 fluxes%netMassOut(:,:) = 0.0
675 fluxes%netMassIn(:,:) = 0.0
676 call read_data(surf_file,
'massout_flux_sum',fluxes%netMassOut, domain=g%Domain%mpp_domain, &
678 call read_data(surf_file,
'massin_flux_sum', fluxes%netMassIn, domain=g%Domain%mpp_domain, &
681 do j=js,je ;
do i=is,ie
682 if (g%mask2dT(i,j)<1.0)
then 683 fluxes%netMassOut(i,j) = 0.0
684 fluxes%netMassIn(i,j) = 0.0
691 call read_data(surf_file,
'ePBL_h_ML', mld, domain=g%Domain%mpp_domain, timelevel=ridx_sum)
699 call read_data(mean_file,
'sw_vis',fluxes%sw_vis_dir, domain=g%Domain%mpp_domain, &
701 call read_data(mean_file,
'sw_nir',fluxes%sw_nir_dir, domain=g%Domain%mpp_domain, &
703 fluxes%sw_vis_dir(:,:) = fluxes%sw_vis_dir(:,:)*0.5
704 fluxes%sw_vis_dif (:,:) = fluxes%sw_vis_dir
705 fluxes%sw_nir_dir(:,:) = fluxes%sw_nir_dir(:,:)*0.5
706 fluxes%sw_nir_dif (:,:) = fluxes%sw_nir_dir
707 fluxes%sw = fluxes%sw_vis_dir + fluxes%sw_vis_dif + fluxes%sw_nir_dir + fluxes%sw_nir_dif
708 do j=js,je ;
do i=is,ie
709 if (g%mask2dT(i,j)<1.0)
then 711 fluxes%sw_vis_dir(i,j) = 0.0
712 fluxes%sw_nir_dir(i,j) = 0.0
713 fluxes%sw_vis_dif (i,j) = 0.0
714 fluxes%sw_nir_dif (i,j) = 0.0
717 call pass_var(fluxes%sw,g%Domain)
718 call pass_var(fluxes%sw_vis_dir,g%Domain)
719 call pass_var(fluxes%sw_vis_dif,g%Domain)
720 call pass_var(fluxes%sw_nir_dir,g%Domain)
721 call pass_var(fluxes%sw_nir_dif,g%Domain)
728 hend, uhtr_all, vhtr_all, hend_all, temp, salt, temp_all, salt_all )
729 type(ocean_grid_type),
intent(inout) :: G
730 type(verticalgrid_type),
intent(in ) :: GV
731 integer,
intent(in ) :: nk_input
732 integer,
intent(in ) :: ridx_sum
733 character(len=200),
intent(in ) :: mean_file
734 character(len=200),
intent(in ) :: sum_file
735 character(len=200),
intent(in ) :: snap_file
736 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: uhtr
737 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: vhtr
738 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: hend
739 real,
dimension(:,:,:,:),
allocatable,
intent(inout) :: uhtr_all
740 real,
dimension(:,:,:,:),
allocatable,
intent(inout) :: vhtr_all
741 real,
dimension(:,:,:,:),
allocatable,
intent(inout) :: hend_all
742 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: temp
743 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: salt
744 real,
dimension(:,:,:,:),
allocatable,
intent(inout) :: temp_all
745 real,
dimension(:,:,:,:),
allocatable,
intent(inout) :: salt_all
747 integer :: i, j, k, is, ie, js, je, nz
748 real,
parameter :: fill_value = 0.
749 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
752 if (.not.
allocated(uhtr_all)) &
753 call mom_error(fatal,
"uhtr_all not allocated before call to update_transport_from_arrays")
754 if (.not.
allocated(vhtr_all)) &
755 call mom_error(fatal,
"vhtr_all not allocated before call to update_transport_from_arrays")
756 if (.not.
allocated(hend_all)) &
757 call mom_error(fatal,
"hend_all not allocated before call to update_transport_from_arrays")
758 if (.not.
allocated(temp_all)) &
759 call mom_error(fatal,
"temp_all not allocated before call to update_transport_from_arrays")
760 if (.not.
allocated(salt_all)) &
761 call mom_error(fatal,
"salt_all not allocated before call to update_transport_from_arrays")
764 do k=1,nk_input ;
do j=js,je ;
do i=is,ie
765 uhtr(i,j,k) = uhtr_all(i,j,k,ridx_sum)
766 vhtr(i,j,k) = vhtr_all(i,j,k,ridx_sum)
767 hend(i,j,k) = hend_all(i,j,k,ridx_sum)
768 temp(i,j,k) = temp_all(i,j,k,ridx_sum)
769 salt(i,j,k) = salt_all(i,j,k,ridx_sum)
770 enddo ;
enddo ;
enddo 773 do k=nk_input+1,nz ;
do j=js,je ;
do i=is,ie
774 uhtr(i,j,k) = fill_value
775 vhtr(i,j,k) = fill_value
776 hend(i,j,k) = fill_value
777 temp(i,j,k) = fill_value
778 salt(i,j,k) = fill_value
779 enddo ;
enddo ;
enddo 790 integer :: read_index
793 integer :: next_modulo_time
795 read_index = mod(inidx+1,numtime)
796 if (read_index < 0) read_index = inidx-read_index
797 if (read_index == 0) read_index = numtime
799 next_modulo_time = read_index
Provides kernels for single-column interpolation, re-integration (re-mapping of integrated quantities...
logical function, public check_column_integrals(nk_1, field_1, nk_2, field_2, missing_value)
Returns false if the column integrals of two given quantities are within roundoff of each other...
This module implements boundary forcing for MOM6.
subroutine, public distribute_residual_vh_upwards(G, GV, hvol, vh)
In the case where offline advection has failed to converge, redistribute the u-flux into layers above...
subroutine, public update_offline_from_files(G, GV, nk_input, mean_file, sum_file, snap_file, surf_file, h_end, uhtr, vhtr, temp_mean, salt_mean, mld, Kd, fluxes, ridx_sum, ridx_snap, read_mld, read_sw, read_ts_uvh, do_ale_in)
Controls the reading in 3d mass fluxes, diffusive fluxes, and other fields stored in a previous integ...
integer, parameter, public to_all
Ocean grid type. See mom_grid for details.
Contains routines related to offline transport of tracers. These routines are likely to be called fro...
Provides the ocean grid type.
The vertvisc_type structure contains vertical viscosities, drag coefficients, and related fields...
subroutine, public calltree_leave(mesg)
Writes a message about leaving a subroutine if call tree reporting is active.
subroutine, public distribute_residual_vh_barotropic(G, GV, hvol, vh)
Redistribute the v-flux as a barotropic equivalent.
This module contains I/O framework code.
subroutine, public offline_add_diurnal_sw(fluxes, G, Time_start, Time_end)
add_diurnal_SW adjusts the shortwave fluxes in an forcying_type variable to add a synthetic diurnal c...
subroutine, public reintegrate_column(nsrc, h_src, uh_src, ndest, h_dest, missing_value, uh_dest)
Conservatively calculate integrated data, uh_dest, on grid h_dest, from layer-integrated data...
logical function, public is_root_pe()
integer function, public next_modulo_time(inidx, numtime)
Calculates the next timelevel to read from the input fields. This allows the 'looping' of the fields...
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
subroutine, public distribute_residual_uh_barotropic(G, GV, hvol, uh)
In the case where offline advection has failed to converge, redistribute the u-flux into remainder of...
subroutine, public update_h_horizontal_flux(G, GV, uhtr, vhtr, h_pre, h_new)
This updates thickness based on the convergence of horizontal mass fluxes NOTE: Only used in non-ALE ...
subroutine, public update_h_vertical_flux(G, GV, ea, eb, h_pre, h_new)
Updates layer thicknesses due to vertical mass transports NOTE: Only used in non-ALE configuration...
subroutine, public distribute_residual_uh_upwards(G, GV, hvol, uh)
In the case where offline advection has failed to converge, redistribute the u-flux into layers above...
subroutine, public update_offline_from_arrays(G, GV, nk_input, ridx_sum, mean_file, sum_file, snap_file, uhtr, vhtr, hend, uhtr_all, vhtr_all, hend_all, temp, salt, temp_all, salt_all)
Fields for offline transport are copied from the stored arrays read during initialization.
subroutine, public limit_mass_flux_3d(G, GV, uh, vh, ea, eb, h_pre)
This routine limits the mass fluxes so that the a layer cannot be completely depleted. NOTE: Only used in non-ALE mode.
subroutine, public mom_error(level, message, all_print)
subroutine, public calltree_enter(mesg, n)
Writes a message about entering a subroutine if call tree reporting is active.