13 use mom_cpu_clock,
only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
24 use mom_forcing_type,
only : set_net_mass_forcing, copy_common_forcing_fields
32 use mom_io,
only : east_face, north_face, num_timelevels
34 use mom_restart,
only : restart_init_end, save_restart, restore_state
35 use mom_time_manager,
only : time_type,
operator(+),
operator(/), get_time, time_type_to_real
47 use idealized_hurricane,
only : idealized_hurricane_wind_forcing, scm_idealized_hurricane_wind_forcing
57 use data_override_mod,
only : data_override_init, data_override
59 implicit none ;
private
61 #include <MOM_memory.h>
64 public surface_forcing_init
65 public forcing_save_restart
71 logical :: use_temperature
72 logical :: restorebuoy
74 logical :: variable_winds
75 logical :: variable_buoyforce
84 real :: latent_heat_fusion
85 real :: latent_heat_vapor
90 logical :: read_gust_2d
91 real,
pointer :: gust(:,:) => null()
94 real,
pointer :: t_restore(:,:) => null()
95 real,
pointer :: s_restore(:,:) => null()
96 real,
pointer :: dens_restore(:,:) => null()
98 integer :: buoy_last_lev_read = -1
102 real :: gyres_taux_const
103 real :: gyres_taux_sin_amp
104 real :: gyres_taux_cos_amp
105 real :: gyres_taux_n_pis
106 logical :: answers_2018
110 logical :: fix_ustar_gustless_bug
113 real :: scurves_ydata(20) = 90.
114 real :: scurves_taux(20) = 0.
121 logical :: first_call_set_forcing = .true.
122 logical :: archaic_omip_file = .true.
124 logical :: dataoverrideisinitialized = .false.
127 real :: constantheatforcing
129 character(len=8) :: wind_stagger
138 character(len=200) :: inputdir
139 character(len=200) :: wind_config
140 character(len=200) :: wind_file
141 character(len=200) :: buoy_config
143 character(len=200) :: longwave_file =
''
144 character(len=200) :: shortwave_file =
''
145 character(len=200) :: evaporation_file =
''
146 character(len=200) :: sensibleheat_file =
''
147 character(len=200) :: latentheat_file =
''
149 character(len=200) :: rain_file =
''
150 character(len=200) :: snow_file =
''
151 character(len=200) :: runoff_file =
''
153 character(len=200) :: longwaveup_file =
''
154 character(len=200) :: shortwaveup_file =
''
156 character(len=200) :: sstrestore_file =
''
158 character(len=200) :: salinityrestore_file =
''
161 character(len=80) :: stress_x_var =
''
162 character(len=80) :: stress_y_var =
''
163 character(len=80) :: ustar_var =
''
164 character(len=80) :: lw_var =
''
165 character(len=80) :: sw_var =
''
166 character(len=80) :: latent_var =
''
167 character(len=80) :: sens_var =
''
168 character(len=80) :: evap_var =
''
169 character(len=80) :: rain_var =
''
170 character(len=80) :: snow_var =
''
171 character(len=80) :: lrunoff_var =
''
172 character(len=80) :: frunoff_var =
''
173 character(len=80) :: sst_restore_var =
''
174 character(len=80) :: sss_restore_var =
''
177 integer :: wind_nlev = -1
178 integer :: sw_nlev = -1
179 integer :: lw_nlev = -1
180 integer :: latent_nlev = -1
181 integer :: sens_nlev = -1
182 integer :: evap_nlev = -1
183 integer :: precip_nlev = -1
184 integer :: runoff_nlev = -1
185 integer :: sst_nlev = -1
186 integer :: sss_nlev = -1
189 integer :: wind_last_lev = -1
190 integer :: sw_last_lev = -1
191 integer :: lw_last_lev = -1
192 integer :: latent_last_lev = -1
193 integer :: sens_last_lev = -1
194 integer :: evap_last_lev = -1
195 integer :: precip_last_lev = -1
196 integer :: runoff_last_lev = -1
197 integer :: sst_last_lev = -1
198 integer :: sss_last_lev = -1
214 integer :: id_clock_forcing
222 subroutine set_forcing(sfc_state, forces, fluxes, day_start, day_interval, G, US, CS)
226 type(
forcing),
intent(inout) :: fluxes
227 type(time_type),
intent(in) :: day_start
228 type(time_type),
intent(in) :: day_interval
235 type(time_type) :: day_center
236 integer :: isd, ied, jsd, jed
237 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
239 call cpu_clock_begin(id_clock_forcing)
240 call calltree_enter(
"set_forcing, MOM_surface_forcing.F90")
242 day_center = day_start + day_interval/2
243 dt = time_type_to_real(day_interval)
245 if (cs%first_call_set_forcing)
then
250 if (trim(cs%buoy_config) /=
"NONE")
then
251 if ( cs%use_temperature )
then
253 if (cs%restorebuoy)
then
254 call safe_alloc_ptr(cs%T_Restore,isd, ied, jsd, jed)
255 call safe_alloc_ptr(fluxes%heat_added, isd, ied, jsd, jed)
256 call safe_alloc_ptr(cs%S_Restore, isd, ied, jsd, jed)
259 call safe_alloc_ptr(fluxes%buoy, isd, ied, jsd, jed)
261 if (cs%restorebuoy)
call safe_alloc_ptr(cs%Dens_Restore, isd, ied, jsd, jed)
267 if (cs%variable_winds .or. cs%first_call_set_forcing)
then
269 if (trim(cs%wind_config) ==
"file")
then
270 call wind_forcing_from_file(sfc_state, forces, day_center, g, us, cs)
271 elseif (trim(cs%wind_config) ==
"data_override")
then
272 call wind_forcing_by_data_override(sfc_state, forces, day_center, g, us, cs)
273 elseif (trim(cs%wind_config) ==
"2gyre")
then
274 call wind_forcing_2gyre(sfc_state, forces, day_center, g, us, cs)
275 elseif (trim(cs%wind_config) ==
"1gyre")
then
276 call wind_forcing_1gyre(sfc_state, forces, day_center, g, us, cs)
277 elseif (trim(cs%wind_config) ==
"gyres")
then
278 call wind_forcing_gyres(sfc_state, forces, day_center, g, us, cs)
279 elseif (trim(cs%wind_config) ==
"zero")
then
280 call wind_forcing_const(sfc_state, forces, 0., 0., day_center, g, us, cs)
281 elseif (trim(cs%wind_config) ==
"const")
then
282 call wind_forcing_const(sfc_state, forces, cs%tau_x0, cs%tau_y0, day_center, g, us, cs)
283 elseif (trim(cs%wind_config) ==
"Neverworld" .or. trim(cs%wind_config) ==
"Neverland")
then
284 call neverworld_wind_forcing(sfc_state, forces, day_center, g, us, cs)
285 elseif (trim(cs%wind_config) ==
"scurves")
then
286 call scurve_wind_forcing(sfc_state, forces, day_center, g, us, cs)
287 elseif (trim(cs%wind_config) ==
"ideal_hurr")
then
288 call idealized_hurricane_wind_forcing(sfc_state, forces, day_center, g, us, cs%idealized_hurricane_CSp)
289 elseif (trim(cs%wind_config) ==
"SCM_ideal_hurr")
then
290 call scm_idealized_hurricane_wind_forcing(sfc_state, forces, day_center, g, us, cs%idealized_hurricane_CSp)
291 elseif (trim(cs%wind_config) ==
"SCM_CVmix_tests")
then
292 call scm_cvmix_tests_wind_forcing(sfc_state, forces, day_center, g, us, cs%SCM_CVmix_tests_CSp)
293 elseif (trim(cs%wind_config) ==
"USER")
then
294 call user_wind_forcing(sfc_state, forces, day_center, g, us, cs%user_forcing_CSp)
295 elseif (cs%variable_winds .and. .not.cs%first_call_set_forcing)
then
296 call mom_error(fatal, &
297 "MOM_surface_forcing: Variable winds defined with no wind config")
299 call mom_error(fatal, &
300 "MOM_surface_forcing:Unrecognized wind config "//trim(cs%wind_config))
305 if (cs%restorebuoy .and. .not.cs%variable_buoyforce)
then
306 call mom_error(fatal,
"With RESTOREBUOY = True, VARIABLE_BUOYFORCE = True should be used. "//&
307 "Otherwise, this can lead to diverging solutions when a simulation "//&
308 "is continued using a restart file.")
311 if ((cs%variable_buoyforce .or. cs%first_call_set_forcing) .and. &
312 (.not.cs%adiabatic))
then
313 if (trim(cs%buoy_config) ==
"file")
then
314 call buoyancy_forcing_from_files(sfc_state, fluxes, day_center, dt, g, us, cs)
315 elseif (trim(cs%buoy_config) ==
"data_override")
then
316 call buoyancy_forcing_from_data_override(sfc_state, fluxes, day_center, dt, g, us, cs)
317 elseif (trim(cs%buoy_config) ==
"zero")
then
318 call buoyancy_forcing_zero(sfc_state, fluxes, day_center, dt, g, cs)
319 elseif (trim(cs%buoy_config) ==
"const")
then
320 call buoyancy_forcing_const(sfc_state, fluxes, day_center, dt, g, us, cs)
321 elseif (trim(cs%buoy_config) ==
"linear")
then
322 call buoyancy_forcing_linear(sfc_state, fluxes, day_center, dt, g, us, cs)
323 elseif (trim(cs%buoy_config) ==
"MESO")
then
324 call meso_buoyancy_forcing(sfc_state, fluxes, day_center, dt, g, us, cs%MESO_forcing_CSp)
325 elseif (trim(cs%buoy_config) ==
"SCM_CVmix_tests")
then
326 call scm_cvmix_tests_buoyancy_forcing(sfc_state, fluxes, day_center, g, us, cs%SCM_CVmix_tests_CSp)
327 elseif (trim(cs%buoy_config) ==
"USER")
then
328 call user_buoyancy_forcing(sfc_state, fluxes, day_center, dt, g, us, cs%user_forcing_CSp)
329 elseif (trim(cs%buoy_config) ==
"BFB")
then
330 call bfb_buoyancy_forcing(sfc_state, fluxes, day_center, dt, g, us, cs%BFB_forcing_CSp)
331 elseif (trim(cs%buoy_config) ==
"dumbbell")
then
332 call dumbbell_buoyancy_forcing(sfc_state, fluxes, day_center, dt, g, us, cs%dumbbell_forcing_CSp)
333 elseif (trim(cs%buoy_config) ==
"NONE")
then
334 call mom_mesg(
"MOM_surface_forcing: buoyancy forcing has been set to omitted.")
335 elseif (cs%variable_buoyforce .and. .not.cs%first_call_set_forcing)
then
336 call mom_error(fatal, &
337 "MOM_surface_forcing: Variable buoy defined with no buoy config.")
339 call mom_error(fatal, &
340 "MOM_surface_forcing: Unrecognized buoy config "//trim(cs%buoy_config))
344 if (
associated(cs%tracer_flow_CSp))
then
345 call call_tracer_set_forcing(sfc_state, fluxes, day_start, day_interval, g, cs%tracer_flow_CSp)
349 call user_alter_forcing(sfc_state, fluxes, day_center, g, cs%urf_CS)
352 if (cs%variable_winds .or. cs%first_call_set_forcing)
then
353 call copy_common_forcing_fields(forces, fluxes, g)
354 call set_derived_forcing_fields(forces, fluxes, g, us, cs%Rho0)
357 if ((cs%variable_buoyforce .or. cs%first_call_set_forcing) .and. &
358 (.not.cs%adiabatic))
then
359 call set_net_mass_forcing(fluxes, forces, g, us)
362 cs%first_call_set_forcing = .false.
364 call cpu_clock_end(id_clock_forcing)
365 call calltree_leave(
"set_forcing")
367 end subroutine set_forcing
370 subroutine wind_forcing_const(sfc_state, forces, tau_x0, tau_y0, day, G, US, CS)
374 real,
intent(in) :: tau_x0
375 real,
intent(in) :: tau_y0
376 type(time_type),
intent(in) :: day
382 real :: Pa_conversion
384 integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
386 call calltree_enter(
"wind_forcing_const, MOM_surface_forcing.F90")
387 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
388 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
389 pa_conversion = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
392 mag_tau = pa_conversion * sqrt( tau_x0**2 + tau_y0**2)
394 do j=js,je ;
do i=is-1,ieq
395 forces%taux(i,j) = tau_x0 * pa_conversion
398 do j=js-1,jeq ;
do i=is,ie
399 forces%tauy(i,j) = tau_y0 * pa_conversion
402 if (cs%read_gust_2d)
then
403 if (
associated(forces%ustar))
then ;
do j=js,je ;
do i=is,ie
404 forces%ustar(i,j) = sqrt( us%L_to_Z * ( mag_tau + cs%gust(i,j) ) / cs%Rho0 )
405 enddo ;
enddo ;
endif
407 if (
associated(forces%ustar))
then ;
do j=js,je ;
do i=is,ie
408 forces%ustar(i,j) = sqrt( us%L_to_Z * ( mag_tau + cs%gust_const ) / cs%Rho0 )
409 enddo ;
enddo ;
endif
412 call calltree_leave(
"wind_forcing_const")
413 end subroutine wind_forcing_const
417 subroutine wind_forcing_2gyre(sfc_state, forces, day, G, US, CS)
421 type(time_type),
intent(in) :: day
428 integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
430 call calltree_enter(
"wind_forcing_2gyre, MOM_surface_forcing.F90")
431 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
432 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
437 do j=js,je ;
do i=is-1,ieq
438 forces%taux(i,j) = 0.1*us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z * &
439 (1.0 - cos(2.0*pi*(g%geoLatCu(i,j)-cs%South_lat) / cs%len_lat))
442 do j=js-1,jeq ;
do i=is,ie
443 forces%tauy(i,j) = 0.0
446 call calltree_leave(
"wind_forcing_2gyre")
447 end subroutine wind_forcing_2gyre
451 subroutine wind_forcing_1gyre(sfc_state, forces, day, G, US, CS)
455 type(time_type),
intent(in) :: day
462 integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
464 call calltree_enter(
"wind_forcing_1gyre, MOM_surface_forcing.F90")
465 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
466 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
471 do j=js,je ;
do i=is-1,ieq
472 forces%taux(i,j) = -0.2*us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z * &
473 cos(pi*(g%geoLatCu(i,j)-cs%South_lat)/cs%len_lat)
476 do j=js-1,jeq ;
do i=is,ie
477 forces%tauy(i,j) = 0.0
480 call calltree_leave(
"wind_forcing_1gyre")
481 end subroutine wind_forcing_1gyre
484 subroutine wind_forcing_gyres(sfc_state, forces, day, G, US, CS)
488 type(time_type),
intent(in) :: day
495 integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
497 call calltree_enter(
"wind_forcing_gyres, MOM_surface_forcing.F90")
498 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
499 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
504 do j=js-1,je+1 ;
do i=is-1,ieq
505 y = (g%geoLatCu(i,j)-cs%South_lat) / cs%len_lat
506 forces%taux(i,j) = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z * &
507 (cs%gyres_taux_const + &
508 ( cs%gyres_taux_sin_amp*sin(cs%gyres_taux_n_pis*pi*y) &
509 + cs%gyres_taux_cos_amp*cos(cs%gyres_taux_n_pis*pi*y) ))
512 do j=js-1,jeq ;
do i=is-1,ie+1
513 forces%tauy(i,j) = 0.0
517 if (cs%answers_2018)
then
518 do j=js,je ;
do i=is,ie
519 forces%ustar(i,j) = sqrt(us%L_to_Z * ((cs%gust_const/cs%Rho0) + &
520 sqrt(0.5*(forces%tauy(i,j-1)*forces%tauy(i,j-1) + forces%tauy(i,j)*forces%tauy(i,j) + &
521 forces%taux(i-1,j)*forces%taux(i-1,j) + forces%taux(i,j)*forces%taux(i,j)))/cs%Rho0) )
524 i_rho = us%L_to_Z / cs%Rho0
525 do j=js,je ;
do i=is,ie
526 forces%ustar(i,j) = sqrt( (cs%gust_const + &
527 sqrt(0.5*((forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2) + &
528 (forces%taux(i-1,j)**2 + forces%taux(i,j)**2))) ) * i_rho )
532 call calltree_leave(
"wind_forcing_gyres")
533 end subroutine wind_forcing_gyres
537 subroutine neverworld_wind_forcing(sfc_state, forces, day, G, US, CS)
541 type(time_type),
intent(in) :: day
547 integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
548 integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
553 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
554 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
555 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
556 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
567 forces%taux(:,:) = 0.0
568 tau_max = 0.2 * us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
570 do j=js,je ;
do i=is-1,ieq
571 y = (g%geoLatT(i,j)-g%south_lat)/g%len_lat
574 forces%taux(i,j) = forces%taux(i,j) + tau_max * ( (1/0.29)*y - ( 1/(2*pi) )*sin( (2*pi*y) / 0.29 ) )
576 if ((y > 0.29) .and. (y <= (0.8-off)))
then
577 forces%taux(i,j) = forces%taux(i,j) + tau_max *(0.35+0.65*cos(pi*(y-0.29)/(0.51-off)) )
579 if ((y > (0.8-off)) .and. (y <= (1-off)))
then
580 forces%taux(i,j) = forces%taux(i,j) + tau_max *( 1.5*( (y-1+off) - (0.1/pi)*sin(10.0*pi*(y-0.8+off)) ) )
582 forces%taux(i,j) = g%mask2dCu(i,j) * forces%taux(i,j)
585 do j=js-1,jeq ;
do i=is,ie
586 forces%tauy(i,j) = g%mask2dCv(i,j) * 0.0
590 if (
associated(forces%ustar))
then
591 i_rho = us%L_to_Z / cs%Rho0
592 do j=js,je ;
do i=is,ie
593 forces%ustar(i,j) = sqrt( (cs%gust_const + &
594 sqrt(0.5*((forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2) + &
595 (forces%taux(i-1,j)**2 + forces%taux(i,j)**2))) ) * i_rho )
599 end subroutine neverworld_wind_forcing
602 subroutine scurve_wind_forcing(sfc_state, forces, day, G, US, CS)
606 type(time_type),
intent(in) :: day
612 integer :: i, j, kseg
613 real :: lon, lat, I_rho, y, L
621 do j=g%jsd,g%jed ;
do i=g%IsdB,g%IedB
622 lon = g%geoLonCu(i,j)
623 lat = g%geoLatCu(i,j)
626 do while (lat>=cs%scurves_ydata(kseg+1) .and. kseg<6)
629 do while (lat<cs%scurves_ydata(kseg) .and. kseg>1)
633 y = lat - cs%scurves_ydata(kseg)
634 l = cs%scurves_ydata(kseg+1) - cs%scurves_ydata(kseg)
635 forces%taux(i,j) = cs%scurves_taux(kseg) + &
636 ( cs%scurves_taux(kseg+1) - cs%scurves_taux(kseg) ) * scurve(y, l)
637 forces%taux(i,j) = g%mask2dCu(i,j) * forces%taux(i,j)
640 do j=g%JsdB,g%JedB ;
do i=g%isd,g%ied
641 forces%tauy(i,j) = g%mask2dCv(i,j) * 0.0
645 if (
associated(forces%ustar))
then
646 i_rho = us%L_to_Z / cs%Rho0
647 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
648 forces%ustar(i,j) = sqrt( (cs%gust_const + &
649 sqrt(0.5*((forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2) + &
650 (forces%taux(i-1,j)**2 + forces%taux(i,j)**2))) ) * i_rho )
654 end subroutine scurve_wind_forcing
657 real function scurve(x,L)
658 real ,
intent(in) :: x
659 real ,
intent(in) :: l
663 scurve = (3. - 2.*s) * (s*s)
667 subroutine wind_forcing_from_file(sfc_state, forces, day, G, US, CS)
671 type(time_type),
intent(in) :: day
677 character(len=200) :: filename
678 real :: temp_x(SZI_(G),SZJ_(G))
679 real :: temp_y(SZI_(G),SZJ_(G))
680 real :: Pa_conversion
682 integer :: time_lev_daily
683 integer :: time_lev_monthly
685 integer :: days, seconds
686 integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
687 logical :: read_Ustar
689 call calltree_enter(
"wind_forcing_from_file, MOM_surface_forcing.F90")
690 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
691 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
692 pa_conversion = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
694 call get_time(day, seconds, days)
695 time_lev_daily = days - 365*floor(real(days) / 365.0)
697 if (time_lev_daily < 31)
then ; time_lev_monthly = 0
698 elseif (time_lev_daily < 59)
then ; time_lev_monthly = 1
699 elseif (time_lev_daily < 90)
then ; time_lev_monthly = 2
700 elseif (time_lev_daily < 120)
then ; time_lev_monthly = 3
701 elseif (time_lev_daily < 151)
then ; time_lev_monthly = 4
702 elseif (time_lev_daily < 181)
then ; time_lev_monthly = 5
703 elseif (time_lev_daily < 212)
then ; time_lev_monthly = 6
704 elseif (time_lev_daily < 243)
then ; time_lev_monthly = 7
705 elseif (time_lev_daily < 273)
then ; time_lev_monthly = 8
706 elseif (time_lev_daily < 304)
then ; time_lev_monthly = 9
707 elseif (time_lev_daily < 334)
then ; time_lev_monthly = 10
708 else ; time_lev_monthly = 11
711 time_lev_daily = time_lev_daily+1
712 time_lev_monthly = time_lev_monthly+1
714 select case (cs%wind_nlev)
715 case (12) ; time_lev = time_lev_monthly
716 case (365) ; time_lev = time_lev_daily
717 case default ; time_lev = 1
720 if (time_lev /= cs%wind_last_lev)
then
721 filename = trim(cs%wind_file)
722 read_ustar = (len_trim(cs%ustar_var) > 0)
726 select case ( uppercase(cs%wind_stagger(1:1)) )
728 temp_x(:,:) = 0.0 ; temp_y(:,:) = 0.0
730 temp_x(:,:), temp_y(:,:), g%Domain, stagger=agrid, &
731 timelevel=time_lev, scale=pa_conversion)
733 call pass_vector(temp_x, temp_y, g%Domain, to_all, agrid)
734 do j=js,je ;
do i=is-1,ieq
735 forces%taux(i,j) = 0.5 * cs%wind_scale * (temp_x(i,j) + temp_x(i+1,j))
737 do j=js-1,jeq ;
do i=is,ie
738 forces%tauy(i,j) = 0.5 * cs%wind_scale * (temp_y(i,j) + temp_y(i,j+1))
741 if (.not.read_ustar)
then
742 if (cs%read_gust_2d)
then
743 do j=js,je ;
do i=is,ie
744 forces%ustar(i,j) = sqrt((cs%gust(i,j) + &
745 sqrt(temp_x(i,j)*temp_x(i,j) + temp_y(i,j)*temp_y(i,j))) * us%L_to_Z / cs%Rho0)
748 do j=js,je ;
do i=is,ie
749 forces%ustar(i,j) = sqrt(us%L_to_Z * (cs%gust_const/cs%Rho0 + &
750 sqrt(temp_x(i,j)*temp_x(i,j) + temp_y(i,j)*temp_y(i,j)) / cs%Rho0) )
755 if (g%symmetric)
then
756 if (.not.
associated(g%Domain_aux))
call mom_error(fatal, &
757 " wind_forcing_from_file with C-grid input and symmetric memory "//&
758 " called with a non-associated auxiliary domain in the grid type.")
761 temp_x(:,:) = 0.0 ; temp_y(:,:) = 0.0
763 temp_x(:,:), temp_y(:,:), &
764 g%Domain_aux, stagger=cgrid_ne, timelevel=time_lev, &
766 do j=js,je ;
do i=is,ie
767 forces%taux(i,j) = cs%wind_scale * temp_x(i,j)
768 forces%tauy(i,j) = cs%wind_scale * temp_y(i,j)
773 forces%taux(:,:), forces%tauy(:,:), &
774 g%Domain, stagger=cgrid_ne, timelevel=time_lev, &
777 if (cs%wind_scale /= 1.0)
then
778 do j=js,je ;
do i=isq,ieq
779 forces%taux(i,j) = cs%wind_scale * forces%taux(i,j)
781 do j=jsq,jeq ;
do i=is,ie
782 forces%tauy(i,j) = cs%wind_scale * forces%tauy(i,j)
787 call pass_vector(forces%taux, forces%tauy, g%Domain, to_all)
788 if (.not.read_ustar)
then
789 if (cs%read_gust_2d)
then
790 do j=js, je ;
do i=is, ie
791 forces%ustar(i,j) = sqrt((cs%gust(i,j) + &
792 sqrt(0.5*((forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2) + &
793 (forces%taux(i-1,j)**2 + forces%taux(i,j)**2))) ) * us%L_to_Z / cs%Rho0 )
796 do j=js, je ;
do i=is, ie
797 forces%ustar(i,j) = sqrt(us%L_to_Z * ( (cs%gust_const/cs%Rho0) + &
798 sqrt(0.5*((forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2) + &
799 (forces%taux(i-1,j)**2 + forces%taux(i,j)**2)))/cs%Rho0))
804 call mom_error(fatal,
"wind_forcing_from_file: Unrecognized stagger "//&
805 trim(cs%wind_stagger)//
" is not 'A' or 'C'.")
809 call mom_read_data(filename, cs%Ustar_var, forces%ustar(:,:), &
810 g%Domain, timelevel=time_lev, scale=us%m_to_Z*us%T_to_s)
813 cs%wind_last_lev = time_lev
817 call calltree_leave(
"wind_forcing_from_file")
818 end subroutine wind_forcing_from_file
822 subroutine wind_forcing_by_data_override(sfc_state, forces, day, G, US, CS)
826 type(time_type),
intent(in) :: day
832 real :: temp_x(SZI_(G),SZJ_(G))
833 real :: temp_y(SZI_(G),SZJ_(G))
834 real :: temp_ustar(SZI_(G),SZJ_(G))
835 real :: Pa_conversion
836 integer :: i, j, is_in, ie_in, js_in, je_in
837 logical :: read_uStar
839 call calltree_enter(
"wind_forcing_by_data_override, MOM_surface_forcing.F90")
841 if (.not.cs%dataOverrideIsInitialized)
then
843 call data_override_init(ocean_domain_in=g%Domain%mpp_domain)
844 cs%dataOverrideIsInitialized = .true.
847 is_in = g%isc - g%isd + 1 ; ie_in = g%iec - g%isd + 1
848 js_in = g%jsc - g%jsd + 1 ; je_in = g%jec - g%jsd + 1
849 pa_conversion = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
851 temp_x(:,:) = 0.0 ; temp_y(:,:) = 0.0
852 call data_override(
'OCN',
'taux', temp_x, day, is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
853 call data_override(
'OCN',
'tauy', temp_y, day, is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
854 call pass_vector(temp_x, temp_y, g%Domain, to_all, agrid)
856 do j=g%jsc,g%jec ;
do i=g%isc-1,g%IecB
857 forces%taux(i,j) = pa_conversion * 0.5 * (temp_x(i,j) + temp_x(i+1,j))
859 do j=g%jsc-1,g%JecB ;
do i=g%isc,g%iec
860 forces%tauy(i,j) = pa_conversion * 0.5 * (temp_y(i,j) + temp_y(i,j+1))
863 read_ustar = (len_trim(cs%ustar_var) > 0)
865 do j=g%jsc,g%jec ;
do i=g%isc,g%iec ; temp_ustar(i,j) = us%Z_to_m*us%s_to_T*forces%ustar(i,j) ;
enddo ;
enddo
866 call data_override(
'OCN',
'ustar', temp_ustar, day, is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
867 do j=g%jsc,g%jec ;
do i=g%isc,g%iec ; forces%ustar(i,j) = us%m_to_Z*us%T_to_s*temp_ustar(i,j) ;
enddo ;
enddo
869 if (cs%read_gust_2d)
then
870 call data_override(
'OCN',
'gust', cs%gust, day, is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
871 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
872 forces%ustar(i,j) = sqrt((pa_conversion * sqrt(temp_x(i,j)*temp_x(i,j) + &
873 temp_y(i,j)*temp_y(i,j)) + cs%gust(i,j)) * us%L_to_Z / cs%Rho0)
876 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
877 forces%ustar(i,j) = sqrt(us%L_to_Z * (pa_conversion*sqrt(temp_x(i,j)*temp_x(i,j) + &
878 temp_y(i,j)*temp_y(i,j))/cs%Rho0 + cs%gust_const/cs%Rho0 ))
883 call pass_vector(forces%taux, forces%tauy, g%Domain, to_all)
886 call calltree_leave(
"wind_forcing_by_data_override")
887 end subroutine wind_forcing_by_data_override
891 subroutine buoyancy_forcing_from_files(sfc_state, fluxes, day, dt, G, US, CS)
894 type(
forcing),
intent(inout) :: fluxes
895 type(time_type),
intent(in) :: day
896 real,
intent(in) :: dt
903 real,
dimension(SZI_(G),SZJ_(G)) :: &
913 real :: kg_m2_s_conversion
917 integer :: time_lev_daily
918 integer :: time_lev_monthly
921 integer :: days, seconds
922 integer :: i, j, is, ie, js, je
924 call calltree_enter(
"buoyancy_forcing_from_files, MOM_surface_forcing.F90")
926 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
927 kg_m2_s_conversion = us%kg_m2s_to_RZ_T
929 if (cs%use_temperature) rhoxcp = cs%Rho0 * fluxes%C_p
932 call get_time(day, seconds, days)
934 time_lev_daily = days - 365*floor(real(days) / 365.0)
936 if (time_lev_daily < 31)
then ; time_lev_monthly = 0
937 elseif (time_lev_daily < 59)
then ; time_lev_monthly = 1
938 elseif (time_lev_daily < 90)
then ; time_lev_monthly = 2
939 elseif (time_lev_daily < 120)
then ; time_lev_monthly = 3
940 elseif (time_lev_daily < 151)
then ; time_lev_monthly = 4
941 elseif (time_lev_daily < 181)
then ; time_lev_monthly = 5
942 elseif (time_lev_daily < 212)
then ; time_lev_monthly = 6
943 elseif (time_lev_daily < 243)
then ; time_lev_monthly = 7
944 elseif (time_lev_daily < 273)
then ; time_lev_monthly = 8
945 elseif (time_lev_daily < 304)
then ; time_lev_monthly = 9
946 elseif (time_lev_daily < 334)
then ; time_lev_monthly = 10
947 else ; time_lev_monthly = 11
950 time_lev_daily = time_lev_daily +1
951 time_lev_monthly = time_lev_monthly+1
953 if (time_lev_daily /= cs%buoy_last_lev_read)
then
956 select case (cs%LW_nlev)
957 case (12) ; time_lev = time_lev_monthly
958 case (365) ; time_lev = time_lev_daily
959 case default ; time_lev = 1
961 call mom_read_data(cs%longwave_file, cs%LW_var, fluxes%lw(:,:), &
962 g%Domain, timelevel=time_lev, scale=us%W_m2_to_QRZ_T)
963 if (cs%archaic_OMIP_file)
then
964 call mom_read_data(cs%longwaveup_file,
"lwup_sfc", temp(:,:), g%Domain, &
965 timelevel=time_lev, scale=us%W_m2_to_QRZ_T)
966 do j=js,je ;
do i=is,ie ; fluxes%LW(i,j) = fluxes%LW(i,j) - temp(i,j) ;
enddo ;
enddo
968 cs%LW_last_lev = time_lev
971 select case (cs%evap_nlev)
972 case (12) ; time_lev = time_lev_monthly
973 case (365) ; time_lev = time_lev_daily
974 case default ; time_lev = 1
976 if (cs%archaic_OMIP_file)
then
977 call mom_read_data(cs%evaporation_file, cs%evap_var, fluxes%evap(:,:), &
978 g%Domain, timelevel=time_lev, scale=-kg_m2_s_conversion)
979 do j=js,je ;
do i=is,ie
980 fluxes%latent(i,j) = cs%latent_heat_vapor*fluxes%evap(i,j)
981 fluxes%latent_evap_diag(i,j) = fluxes%latent(i,j)
984 call mom_read_data(cs%evaporation_file, cs%evap_var, fluxes%evap(:,:), &
985 g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
987 cs%evap_last_lev = time_lev
989 select case (cs%latent_nlev)
990 case (12) ; time_lev = time_lev_monthly
991 case (365) ; time_lev = time_lev_daily
992 case default ; time_lev = 1
994 if (.not.cs%archaic_OMIP_file)
then
995 call mom_read_data(cs%latentheat_file, cs%latent_var, fluxes%latent(:,:), &
996 g%Domain, timelevel=time_lev, scale=us%W_m2_to_QRZ_T)
997 do j=js,je ;
do i=is,ie
998 fluxes%latent_evap_diag(i,j) = fluxes%latent(i,j)
1001 cs%latent_last_lev = time_lev
1003 select case (cs%sens_nlev)
1004 case (12) ; time_lev = time_lev_monthly
1005 case (365) ; time_lev = time_lev_daily
1006 case default ; time_lev = 1
1008 if (cs%archaic_OMIP_file)
then
1009 call mom_read_data(cs%sensibleheat_file, cs%sens_var, fluxes%sens(:,:), &
1010 g%Domain, timelevel=time_lev, scale=-us%W_m2_to_QRZ_T)
1012 call mom_read_data(cs%sensibleheat_file, cs%sens_var, fluxes%sens(:,:), &
1013 g%Domain, timelevel=time_lev, scale=us%W_m2_to_QRZ_T)
1015 cs%sens_last_lev = time_lev
1017 select case (cs%SW_nlev)
1018 case (12) ; time_lev = time_lev_monthly
1019 case (365) ; time_lev = time_lev_daily
1020 case default ; time_lev = 1
1022 call mom_read_data(cs%shortwave_file, cs%SW_var, fluxes%sw(:,:), g%Domain, &
1023 timelevel=time_lev, scale=us%W_m2_to_QRZ_T)
1024 if (cs%archaic_OMIP_file)
then
1025 call mom_read_data(cs%shortwaveup_file,
"swup_sfc", temp(:,:), g%Domain, &
1026 timelevel=time_lev, scale=us%W_m2_to_QRZ_T)
1027 do j=js,je ;
do i=is,ie
1028 fluxes%sw(i,j) = fluxes%sw(i,j) - temp(i,j)
1031 cs%SW_last_lev = time_lev
1033 select case (cs%precip_nlev)
1034 case (12) ; time_lev = time_lev_monthly
1035 case (365) ; time_lev = time_lev_daily
1036 case default ; time_lev = 1
1039 fluxes%fprec(:,:), g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
1041 fluxes%lprec(:,:), g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
1042 if (cs%archaic_OMIP_file)
then
1043 do j=js,je ;
do i=is,ie
1044 fluxes%lprec(i,j) = fluxes%lprec(i,j) - fluxes%fprec(i,j)
1047 cs%precip_last_lev = time_lev
1049 select case (cs%runoff_nlev)
1050 case (12) ; time_lev = time_lev_monthly
1051 case (365) ; time_lev = time_lev_daily
1052 case default ; time_lev = 1
1054 if (cs%archaic_OMIP_file)
then
1055 call mom_read_data(cs%runoff_file, cs%lrunoff_var, temp(:,:), &
1056 g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
1057 do j=js,je ;
do i=is,ie
1058 fluxes%lrunoff(i,j) = temp(i,j)*us%m_to_L**2*g%IareaT(i,j)
1060 call mom_read_data(cs%runoff_file, cs%frunoff_var, temp(:,:), &
1061 g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
1062 do j=js,je ;
do i=is,ie
1063 fluxes%frunoff(i,j) = temp(i,j)*us%m_to_L**2*g%IareaT(i,j)
1066 call mom_read_data(cs%runoff_file, cs%lrunoff_var, fluxes%lrunoff(:,:), &
1067 g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
1068 call mom_read_data(cs%runoff_file, cs%frunoff_var, fluxes%frunoff(:,:), &
1069 g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
1071 cs%runoff_last_lev = time_lev
1074 if (cs%restorebuoy)
then
1075 select case (cs%SST_nlev)
1076 case (12) ; time_lev = time_lev_monthly
1077 case (365) ; time_lev = time_lev_daily
1078 case default ; time_lev = 1
1080 call mom_read_data(cs%SSTrestore_file, cs%SST_restore_var, &
1081 cs%T_Restore(:,:), g%Domain, timelevel=time_lev)
1082 cs%SST_last_lev = time_lev
1084 select case (cs%SSS_nlev)
1085 case (12) ; time_lev = time_lev_monthly
1086 case (365) ; time_lev = time_lev_daily
1087 case default ; time_lev = 1
1089 call mom_read_data(cs%salinityrestore_file, cs%SSS_restore_var, &
1090 cs%S_Restore(:,:), g%Domain, timelevel=time_lev)
1091 cs%SSS_last_lev = time_lev
1093 cs%buoy_last_lev_read = time_lev_daily
1101 do j=js,je ;
do i=is,ie
1102 fluxes%evap(i,j) = fluxes%evap(i,j) * g%mask2dT(i,j)
1103 fluxes%lprec(i,j) = fluxes%lprec(i,j) * g%mask2dT(i,j)
1104 fluxes%fprec(i,j) = fluxes%fprec(i,j) * g%mask2dT(i,j)
1105 fluxes%lrunoff(i,j) = fluxes%lrunoff(i,j) * g%mask2dT(i,j)
1106 fluxes%frunoff(i,j) = fluxes%frunoff(i,j) * g%mask2dT(i,j)
1107 fluxes%lw(i,j) = fluxes%lw(i,j) * g%mask2dT(i,j)
1108 fluxes%sens(i,j) = fluxes%sens(i,j) * g%mask2dT(i,j)
1109 fluxes%sw(i,j) = fluxes%sw(i,j) * g%mask2dT(i,j)
1110 fluxes%latent(i,j) = fluxes%latent(i,j) * g%mask2dT(i,j)
1112 fluxes%latent_evap_diag(i,j) = fluxes%latent_evap_diag(i,j) * g%mask2dT(i,j)
1113 fluxes%latent_fprec_diag(i,j) = -fluxes%fprec(i,j)*cs%latent_heat_fusion
1114 fluxes%latent_frunoff_diag(i,j) = -fluxes%frunoff(i,j)*cs%latent_heat_fusion
1121 if (cs%restorebuoy)
then
1123 if (cs%use_temperature)
then
1124 do j=js,je ;
do i=is,ie
1125 if (g%mask2dT(i,j) > 0)
then
1126 fluxes%heat_added(i,j) = g%mask2dT(i,j) * &
1127 ((cs%T_Restore(i,j) - sfc_state%SST(i,j)) * rhoxcp * cs%Flux_const_T)
1128 fluxes%vprec(i,j) = - (cs%Rho0*cs%Flux_const_S) * &
1129 (cs%S_Restore(i,j) - sfc_state%SSS(i,j)) / &
1130 (0.5*(sfc_state%SSS(i,j) + cs%S_Restore(i,j)))
1132 fluxes%heat_added(i,j) = 0.0
1133 fluxes%vprec(i,j) = 0.0
1137 do j=js,je ;
do i=is,ie
1138 if (g%mask2dT(i,j) > 0)
then
1139 fluxes%buoy(i,j) = (cs%Dens_Restore(i,j) - sfc_state%sfc_density(i,j)) * &
1140 (cs%G_Earth * cs%Flux_const / cs%Rho0)
1142 fluxes%buoy(i,j) = 0.0
1148 if (.not.cs%use_temperature)
then
1149 call mom_error(fatal,
"buoyancy_forcing in MOM_surface_forcing: "// &
1150 "The fluxes need to be defined without RESTOREBUOY.")
1165 call calltree_leave(
"buoyancy_forcing_from_files")
1166 end subroutine buoyancy_forcing_from_files
1169 subroutine buoyancy_forcing_from_data_override(sfc_state, fluxes, day, dt, G, US, CS)
1172 type(
forcing),
intent(inout) :: fluxes
1173 type(time_type),
intent(in) :: day
1174 real,
intent(in) :: dt
1181 real,
dimension(SZI_(G),SZJ_(G)) :: &
1190 real :: kg_m2_s_conversion
1194 integer :: time_lev_daily
1195 integer :: time_lev_monthly
1196 integer :: itime_lev
1198 integer :: days, seconds
1199 integer :: i, j, is, ie, js, je, isd, ied, jsd, jed
1200 integer :: is_in, ie_in, js_in, je_in
1202 call calltree_enter(
"buoyancy_forcing_from_data_override, MOM_surface_forcing.F90")
1204 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1205 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1206 kg_m2_s_conversion = us%kg_m2s_to_RZ_T
1208 if (cs%use_temperature) rhoxcp = cs%Rho0 * fluxes%C_p
1210 if (.not.cs%dataOverrideIsInitialized)
then
1211 call data_override_init(ocean_domain_in=g%Domain%mpp_domain)
1212 cs%dataOverrideIsInitialized = .true.
1215 is_in = g%isc - g%isd + 1
1216 ie_in = g%iec - g%isd + 1
1217 js_in = g%jsc - g%jsd + 1
1218 je_in = g%jec - g%jsd + 1
1220 call data_override(
'OCN',
'lw', fluxes%lw(:,:), day, &
1221 is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1222 if (us%QRZ_T_to_W_m2 /= 1.0)
then ;
do j=js,je ;
do i=is,ie
1223 fluxes%lw(i,j) = fluxes%lw(i,j) * us%W_m2_to_QRZ_T
1224 enddo ;
enddo ;
endif
1225 call data_override(
'OCN',
'evap', fluxes%evap(:,:), day, &
1226 is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1229 do j=js,je ;
do i=is,ie
1233 fluxes%evap(i,j) = -kg_m2_s_conversion*fluxes%evap(i,j)
1234 fluxes%latent(i,j) = cs%latent_heat_vapor*fluxes%evap(i,j)
1235 fluxes%latent_evap_diag(i,j) = fluxes%latent(i,j)
1238 call data_override(
'OCN',
'sens', fluxes%sens(:,:), day, &
1239 is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1242 do j=js,je ;
do i=is,ie
1243 fluxes%sens(i,j) = -us%W_m2_to_QRZ_T * fluxes%sens(i,j)
1247 call data_override(
'OCN',
'sw', fluxes%sw(:,:), day, &
1248 is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1249 if (us%QRZ_T_to_W_m2 /= 1.0)
then ;
do j=js,je ;
do i=is,ie
1250 fluxes%sw(i,j) = fluxes%sw(i,j) * us%W_m2_to_QRZ_T
1251 enddo ;
enddo ;
endif
1253 call data_override(
'OCN',
'snow', fluxes%fprec(:,:), day, &
1254 is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1256 call data_override(
'OCN',
'rain', fluxes%lprec(:,:), day, &
1257 is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1259 call data_override(
'OCN',
'runoff', fluxes%lrunoff(:,:), day, &
1260 is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1262 call data_override(
'OCN',
'calving', fluxes%frunoff(:,:), day, &
1263 is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1265 if (kg_m2_s_conversion /= 1.0)
then ;
do j=js,je ;
do i=is,ie
1266 fluxes%lprec(i,j) = fluxes%lprec(i,j) * kg_m2_s_conversion
1267 fluxes%fprec(i,j) = fluxes%fprec(i,j) * kg_m2_s_conversion
1268 fluxes%lrunoff(i,j) = fluxes%lrunoff(i,j) * kg_m2_s_conversion
1269 fluxes%frunoff(i,j) = fluxes%frunoff(i,j) * kg_m2_s_conversion
1270 enddo ;
enddo ;
endif
1273 if (cs%restorebuoy)
then
1274 call data_override(
'OCN',
'SST_restore', cs%T_restore(:,:), day, &
1275 is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1277 call data_override(
'OCN',
'SSS_restore', cs%S_restore(:,:), day, &
1278 is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1283 if (cs%restorebuoy)
then
1284 if (cs%use_temperature)
then
1285 do j=js,je ;
do i=is,ie
1286 if (g%mask2dT(i,j) > 0)
then
1287 fluxes%heat_added(i,j) = g%mask2dT(i,j) * &
1288 ((cs%T_Restore(i,j) - sfc_state%SST(i,j)) * rhoxcp * cs%Flux_const_T)
1289 fluxes%vprec(i,j) = - (cs%Rho0*cs%Flux_const_S) * &
1290 (cs%S_Restore(i,j) - sfc_state%SSS(i,j)) / &
1291 (0.5*(sfc_state%SSS(i,j) + cs%S_Restore(i,j)))
1293 fluxes%heat_added(i,j) = 0.0
1294 fluxes%vprec(i,j) = 0.0
1298 do j=js,je ;
do i=is,ie
1299 if (g%mask2dT(i,j) > 0)
then
1300 fluxes%buoy(i,j) = (cs%Dens_Restore(i,j) - sfc_state%sfc_density(i,j)) * &
1301 (cs%G_Earth * cs%Flux_const / cs%Rho0)
1303 fluxes%buoy(i,j) = 0.0
1308 if (.not.cs%use_temperature)
then
1309 call mom_error(fatal,
"buoyancy_forcing in MOM_surface_forcing: "// &
1310 "The fluxes need to be defined without RESTOREBUOY.")
1321 do j=js,je ;
do i=is,ie
1322 fluxes%evap(i,j) = fluxes%evap(i,j) * g%mask2dT(i,j)
1323 fluxes%lprec(i,j) = fluxes%lprec(i,j) * g%mask2dT(i,j)
1324 fluxes%fprec(i,j) = fluxes%fprec(i,j) * g%mask2dT(i,j)
1325 fluxes%lrunoff(i,j) = fluxes%lrunoff(i,j) * g%mask2dT(i,j)
1326 fluxes%frunoff(i,j) = fluxes%frunoff(i,j) * g%mask2dT(i,j)
1327 fluxes%lw(i,j) = fluxes%lw(i,j) * g%mask2dT(i,j)
1328 fluxes%latent(i,j) = fluxes%latent(i,j) * g%mask2dT(i,j)
1329 fluxes%sens(i,j) = fluxes%sens(i,j) * g%mask2dT(i,j)
1330 fluxes%sw(i,j) = fluxes%sw(i,j) * g%mask2dT(i,j)
1332 fluxes%latent_evap_diag(i,j) = fluxes%latent_evap_diag(i,j) * g%mask2dT(i,j)
1333 fluxes%latent_fprec_diag(i,j) = -fluxes%fprec(i,j)*cs%latent_heat_fusion
1334 fluxes%latent_frunoff_diag(i,j) = -fluxes%frunoff(i,j)*cs%latent_heat_fusion
1348 call calltree_leave(
"buoyancy_forcing_from_data_override")
1349 end subroutine buoyancy_forcing_from_data_override
1352 subroutine buoyancy_forcing_zero(sfc_state, fluxes, day, dt, G, CS)
1355 type(
forcing),
intent(inout) :: fluxes
1356 type(time_type),
intent(in) :: day
1357 real,
intent(in) :: dt
1363 integer :: i, j, is, ie, js, je
1365 call calltree_enter(
"buoyancy_forcing_zero, MOM_surface_forcing.F90")
1366 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1368 if (cs%use_temperature)
then
1369 do j=js,je ;
do i=is,ie
1370 fluxes%evap(i,j) = 0.0
1371 fluxes%lprec(i,j) = 0.0
1372 fluxes%fprec(i,j) = 0.0
1373 fluxes%vprec(i,j) = 0.0
1374 fluxes%lrunoff(i,j) = 0.0
1375 fluxes%frunoff(i,j) = 0.0
1376 fluxes%lw(i,j) = 0.0
1377 fluxes%latent(i,j) = 0.0
1378 fluxes%sens(i,j) = 0.0
1379 fluxes%sw(i,j) = 0.0
1380 fluxes%latent_evap_diag(i,j) = 0.0
1381 fluxes%latent_fprec_diag(i,j) = 0.0
1382 fluxes%latent_frunoff_diag(i,j) = 0.0
1385 do j=js,je ;
do i=is,ie
1386 fluxes%buoy(i,j) = 0.0
1390 call calltree_leave(
"buoyancy_forcing_zero")
1391 end subroutine buoyancy_forcing_zero
1395 subroutine buoyancy_forcing_const(sfc_state, fluxes, day, dt, G, US, CS)
1398 type(
forcing),
intent(inout) :: fluxes
1399 type(time_type),
intent(in) :: day
1400 real,
intent(in) :: dt
1407 integer :: i, j, is, ie, js, je
1408 call calltree_enter(
"buoyancy_forcing_const, MOM_surface_forcing.F90")
1409 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1411 if (cs%use_temperature)
then
1412 do j=js,je ;
do i=is,ie
1413 fluxes%evap(i,j) = 0.0
1414 fluxes%lprec(i,j) = 0.0
1415 fluxes%fprec(i,j) = 0.0
1416 fluxes%vprec(i,j) = 0.0
1417 fluxes%lrunoff(i,j) = 0.0
1418 fluxes%frunoff(i,j) = 0.0
1419 fluxes%lw(i,j) = 0.0
1420 fluxes%latent(i,j) = 0.0
1421 fluxes%sens(i,j) = cs%constantHeatForcing * g%mask2dT(i,j)
1422 fluxes%sw(i,j) = 0.0
1423 fluxes%latent_evap_diag(i,j) = 0.0
1424 fluxes%latent_fprec_diag(i,j) = 0.0
1425 fluxes%latent_frunoff_diag(i,j) = 0.0
1428 do j=js,je ;
do i=is,ie
1429 fluxes%buoy(i,j) = 0.0
1433 call calltree_leave(
"buoyancy_forcing_const")
1434 end subroutine buoyancy_forcing_const
1438 subroutine buoyancy_forcing_linear(sfc_state, fluxes, day, dt, G, US, CS)
1441 type(
forcing),
intent(inout) :: fluxes
1442 type(time_type),
intent(in) :: day
1443 real,
intent(in) :: dt
1450 real :: y, T_restore, S_restore
1451 integer :: i, j, is, ie, js, je
1453 call calltree_enter(
"buoyancy_forcing_linear, MOM_surface_forcing.F90")
1454 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1457 if (cs%use_temperature)
then
1458 do j=js,je ;
do i=is,ie
1459 fluxes%evap(i,j) = 0.0
1460 fluxes%lprec(i,j) = 0.0
1461 fluxes%fprec(i,j) = 0.0
1462 fluxes%vprec(i,j) = 0.0
1463 fluxes%lrunoff(i,j) = 0.0
1464 fluxes%frunoff(i,j) = 0.0
1465 fluxes%lw(i,j) = 0.0
1466 fluxes%latent(i,j) = 0.0
1467 fluxes%sens(i,j) = 0.0
1468 fluxes%sw(i,j) = 0.0
1469 fluxes%latent_evap_diag(i,j) = 0.0
1470 fluxes%latent_fprec_diag(i,j) = 0.0
1471 fluxes%latent_frunoff_diag(i,j) = 0.0
1474 do j=js,je ;
do i=is,ie
1475 fluxes%buoy(i,j) = 0.0
1479 if (cs%restorebuoy)
then
1480 if (cs%use_temperature)
then
1481 do j=js,je ;
do i=is,ie
1482 y = (g%geoLatCu(i,j)-cs%South_lat)/cs%len_lat
1483 t_restore = cs%T_south + (cs%T_north-cs%T_south)*y
1484 s_restore = cs%S_south + (cs%S_north-cs%S_south)*y
1485 if (g%mask2dT(i,j) > 0)
then
1486 fluxes%heat_added(i,j) = g%mask2dT(i,j) * &
1487 ((t_restore - sfc_state%SST(i,j)) * ((cs%Rho0 * fluxes%C_p) * cs%Flux_const))
1488 fluxes%vprec(i,j) = - (cs%Rho0*cs%Flux_const) * &
1489 (s_restore - sfc_state%SSS(i,j)) / &
1490 (0.5*(sfc_state%SSS(i,j) + s_restore))
1492 fluxes%heat_added(i,j) = 0.0
1493 fluxes%vprec(i,j) = 0.0
1497 call mom_error(fatal,
"buoyancy_forcing_linear in MOM_surface_forcing: "// &
1498 "RESTOREBUOY to linear not written yet.")
1509 if (.not.cs%use_temperature)
then
1510 call mom_error(fatal,
"buoyancy_forcing_linear in MOM_surface_forcing: "// &
1511 "The fluxes need to be defined without RESTOREBUOY.")
1515 call calltree_leave(
"buoyancy_forcing_linear")
1516 end subroutine buoyancy_forcing_linear
1519 subroutine forcing_save_restart(CS, G, Time, directory, time_stamped, &
1524 type(time_type),
intent(in) :: time
1525 character(len=*),
intent(in) :: directory
1526 logical,
optional,
intent(in) :: time_stamped
1528 character(len=*),
optional,
intent(in) :: filename_suffix
1531 if (.not.
associated(cs))
return
1532 if (.not.
associated(cs%restart_CSp))
return
1534 call save_restart(directory, time, g, cs%restart_CSp, time_stamped)
1536 end subroutine forcing_save_restart
1539 subroutine surface_forcing_init(Time, G, US, param_file, diag, CS, tracer_flow_CSp)
1540 type(time_type),
intent(in) :: time
1544 type(
diag_ctrl),
target,
intent(inout) :: diag
1552 type(time_type) :: time_frc
1554 # include "version_variable.h"
1555 real :: flux_const_default
1556 logical :: default_2018_answers
1557 character(len=40) :: mdl =
"MOM_surface_forcing"
1558 character(len=200) :: filename, gust_file
1560 if (
associated(cs))
then
1561 call mom_error(warning,
"surface_forcing_init called with an associated "// &
1562 "control structure.")
1567 id_clock_forcing=cpu_clock_id(
'(Ocean surface forcing)', grain=clock_module)
1568 call cpu_clock_begin(id_clock_forcing)
1571 if (
associated(tracer_flow_csp)) cs%tracer_flow_CSp => tracer_flow_csp
1575 call get_param(param_file, mdl,
"ENABLE_THERMODYNAMICS", cs%use_temperature, &
1576 "If true, Temperature and salinity are used as state "//&
1577 "variables.", default=.true.)
1578 call get_param(param_file, mdl,
"INPUTDIR", cs%inputdir, &
1579 "The directory in which all input files are found.", &
1581 cs%inputdir = slasher(cs%inputdir)
1583 call get_param(param_file, mdl,
"ADIABATIC", cs%adiabatic, &
1584 "There are no diapycnal mass fluxes if ADIABATIC is "//&
1585 "true. This assumes that KD = KDML = 0.0 and that "//&
1586 "there is no buoyancy forcing, but makes the model "//&
1587 "faster by eliminating subroutine calls.", default=.false.)
1588 call get_param(param_file, mdl,
"VARIABLE_WINDS", cs%variable_winds, &
1589 "If true, the winds vary in time after the initialization.", &
1591 call get_param(param_file, mdl,
"VARIABLE_BUOYFORCE", cs%variable_buoyforce, &
1592 "If true, the buoyancy forcing varies in time after the "//&
1593 "initialization of the model.", default=.true.)
1595 call get_param(param_file, mdl,
"BUOY_CONFIG", cs%buoy_config, &
1596 "The character string that indicates how buoyancy forcing "//&
1597 "is specified. Valid options include (file), (zero), "//&
1598 "(linear), (USER), (BFB) and (NONE).", default=
"zero")
1599 if (trim(cs%buoy_config) ==
"file")
then
1600 call get_param(param_file, mdl,
"ARCHAIC_OMIP_FORCING_FILE", cs%archaic_OMIP_file, &
1601 "If true, use the forcing variable decomposition from "//&
1602 "the old German OMIP prescription that predated CORE. If "//&
1603 "false, use the variable groupings available from MOM "//&
1604 "output diagnostics of forcing variables.", default=.true.)
1605 if (cs%archaic_OMIP_file)
then
1606 call get_param(param_file, mdl,
"LONGWAVEDOWN_FILE", cs%longwave_file, &
1607 "The file with the downward longwave heat flux, in "//&
1608 "variable lwdn_sfc.", fail_if_missing=.true.)
1609 call get_param(param_file, mdl,
"LONGWAVEUP_FILE", cs%longwaveup_file, &
1610 "The file with the upward longwave heat flux, in "//&
1611 "variable lwup_sfc.", fail_if_missing=.true.)
1612 call get_param(param_file, mdl,
"EVAPORATION_FILE", cs%evaporation_file, &
1613 "The file with the evaporative moisture flux, in "//&
1614 "variable evap.", fail_if_missing=.true.)
1615 call get_param(param_file, mdl,
"SENSIBLEHEAT_FILE", cs%sensibleheat_file, &
1616 "The file with the sensible heat flux, in "//&
1617 "variable shflx.", fail_if_missing=.true.)
1618 call get_param(param_file, mdl,
"SHORTWAVEUP_FILE", cs%shortwaveup_file, &
1619 "The file with the upward shortwave heat flux.", &
1620 fail_if_missing=.true.)
1621 call get_param(param_file, mdl,
"SHORTWAVEDOWN_FILE", cs%shortwave_file, &
1622 "The file with the downward shortwave heat flux.", &
1623 fail_if_missing=.true.)
1624 call get_param(param_file, mdl,
"SNOW_FILE", cs%snow_file, &
1625 "The file with the downward frozen precip flux, in "//&
1626 "variable snow.", fail_if_missing=.true.)
1627 call get_param(param_file, mdl,
"PRECIP_FILE", cs%rain_file, &
1628 "The file with the downward total precip flux, in "//&
1629 "variable precip.", fail_if_missing=.true.)
1630 call get_param(param_file, mdl,
"FRESHDISCHARGE_FILE", cs%runoff_file, &
1631 "The file with the fresh and frozen runoff/calving fluxes, "//&
1632 "invariables disch_w and disch_s.", fail_if_missing=.true.)
1635 cs%latentheat_file = cs%evaporation_file ; cs%latent_var =
"evap"
1636 cs%LW_var =
"lwdn_sfc"; cs%SW_var =
"swdn_sfc"; cs%sens_var =
"shflx"
1637 cs%evap_var =
"evap"; cs%rain_var =
"precip"; cs%snow_var =
"snow"
1638 cs%lrunoff_var =
"disch_w"; cs%frunoff_var =
"disch_s"
1641 call get_param(param_file, mdl,
"LONGWAVE_FILE", cs%longwave_file, &
1642 "The file with the longwave heat flux, in the variable "//&
1643 "given by LONGWAVE_FORCING_VAR.", fail_if_missing=.true.)
1644 call get_param(param_file, mdl,
"LONGWAVE_FORCING_VAR", cs%LW_var, &
1645 "The variable with the longwave forcing field.", default=
"LW")
1647 call get_param(param_file, mdl,
"SHORTWAVE_FILE", cs%shortwave_file, &
1648 "The file with the shortwave heat flux, in the variable "//&
1649 "given by SHORTWAVE_FORCING_VAR.", fail_if_missing=.true.)
1650 call get_param(param_file, mdl,
"SHORTWAVE_FORCING_VAR", cs%SW_var, &
1651 "The variable with the shortwave forcing field.", default=
"SW")
1653 call get_param(param_file, mdl,
"EVAPORATION_FILE", cs%evaporation_file, &
1654 "The file with the evaporative moisture flux, in the "//&
1655 "variable given by EVAP_FORCING_VAR.", fail_if_missing=.true.)
1656 call get_param(param_file, mdl,
"EVAP_FORCING_VAR", cs%evap_var, &
1657 "The variable with the evaporative moisture flux.", &
1660 call get_param(param_file, mdl,
"LATENTHEAT_FILE", cs%latentheat_file, &
1661 "The file with the latent heat flux, in the variable "//&
1662 "given by LATENT_FORCING_VAR.", fail_if_missing=.true.)
1663 call get_param(param_file, mdl,
"LATENT_FORCING_VAR", cs%latent_var, &
1664 "The variable with the latent heat flux.", default=
"latent")
1666 call get_param(param_file, mdl,
"SENSIBLEHEAT_FILE", cs%sensibleheat_file, &
1667 "The file with the sensible heat flux, in the variable "//&
1668 "given by SENSIBLE_FORCING_VAR.", fail_if_missing=.true.)
1669 call get_param(param_file, mdl,
"SENSIBLE_FORCING_VAR", cs%sens_var, &
1670 "The variable with the sensible heat flux.", default=
"sensible")
1672 call get_param(param_file, mdl,
"RAIN_FILE", cs%rain_file, &
1673 "The file with the liquid precipitation flux, in the "//&
1674 "variable given by RAIN_FORCING_VAR.", fail_if_missing=.true.)
1675 call get_param(param_file, mdl,
"RAIN_FORCING_VAR", cs%rain_var, &
1676 "The variable with the liquid precipitation flux.", &
1677 default=
"liq_precip")
1678 call get_param(param_file, mdl,
"SNOW_FILE", cs%snow_file, &
1679 "The file with the frozen precipitation flux, in the "//&
1680 "variable given by SNOW_FORCING_VAR.", fail_if_missing=.true.)
1681 call get_param(param_file, mdl,
"SNOW_FORCING_VAR", cs%snow_var, &
1682 "The variable with the frozen precipitation flux.", &
1683 default=
"froz_precip")
1685 call get_param(param_file, mdl,
"RUNOFF_FILE", cs%runoff_file, &
1686 "The file with the fresh and frozen runoff/calving "//&
1687 "fluxes, in variables given by LIQ_RUNOFF_FORCING_VAR "//&
1688 "and FROZ_RUNOFF_FORCING_VAR.", fail_if_missing=.true.)
1689 call get_param(param_file, mdl,
"LIQ_RUNOFF_FORCING_VAR", cs%lrunoff_var, &
1690 "The variable with the liquid runoff flux.", &
1691 default=
"liq_runoff")
1692 call get_param(param_file, mdl,
"FROZ_RUNOFF_FORCING_VAR", cs%frunoff_var, &
1693 "The variable with the frozen runoff flux.", &
1694 default=
"froz_runoff")
1697 call get_param(param_file, mdl,
"SSTRESTORE_FILE", cs%SSTrestore_file, &
1698 "The file with the SST toward which to restore in the "//&
1699 "variable given by SST_RESTORE_VAR.", fail_if_missing=.true.)
1700 call get_param(param_file, mdl,
"SALINITYRESTORE_FILE", cs%salinityrestore_file, &
1701 "The file with the surface salinity toward which to "//&
1702 "restore in the variable given by SSS_RESTORE_VAR.", &
1703 fail_if_missing=.true.)
1705 if (cs%archaic_OMIP_file)
then
1706 cs%SST_restore_var =
"TEMP" ; cs%SSS_restore_var =
"SALT"
1708 call get_param(param_file, mdl,
"SST_RESTORE_VAR", cs%SST_restore_var, &
1709 "The variable with the SST toward which to restore.", &
1711 call get_param(param_file, mdl,
"SSS_RESTORE_VAR", cs%SSS_restore_var, &
1712 "The variable with the SSS toward which to restore.", &
1717 cs%shortwave_file = trim(cs%inputdir)//trim(cs%shortwave_file)
1718 cs%longwave_file = trim(cs%inputdir)//trim(cs%longwave_file)
1719 cs%sensibleheat_file = trim(cs%inputdir)//trim(cs%sensibleheat_file)
1720 cs%latentheat_file = trim(cs%inputdir)//trim(cs%latentheat_file)
1721 cs%evaporation_file = trim(cs%inputdir)//trim(cs%evaporation_file)
1722 cs%snow_file = trim(cs%inputdir)//trim(cs%snow_file)
1723 cs%rain_file = trim(cs%inputdir)//trim(cs%rain_file)
1724 cs%runoff_file = trim(cs%inputdir)//trim(cs%runoff_file)
1726 cs%shortwaveup_file = trim(cs%inputdir)//trim(cs%shortwaveup_file)
1727 cs%longwaveup_file = trim(cs%inputdir)//trim(cs%longwaveup_file)
1729 cs%SSTrestore_file = trim(cs%inputdir)//trim(cs%SSTrestore_file)
1730 cs%salinityrestore_file = trim(cs%inputdir)//trim(cs%salinityrestore_file)
1731 elseif (trim(cs%buoy_config) ==
"const")
then
1732 call get_param(param_file, mdl,
"SENSIBLE_HEAT_FLUX", cs%constantHeatForcing, &
1733 "A constant heat forcing (positive into ocean) applied "//&
1734 "through the sensible heat flux field. ", &
1735 units=
'W/m2', scale=us%W_m2_to_QRZ_T, fail_if_missing=.true.)
1737 call get_param(param_file, mdl,
"WIND_CONFIG", cs%wind_config, &
1738 "The character string that indicates how wind forcing "//&
1739 "is specified. Valid options include (file), (2gyre), "//&
1740 "(1gyre), (gyres), (zero), and (USER).", default=
"zero")
1741 if (trim(cs%wind_config) ==
"file")
then
1742 call get_param(param_file, mdl,
"WIND_FILE", cs%wind_file, &
1743 "The file in which the wind stresses are found in "//&
1744 "variables STRESS_X and STRESS_Y.", fail_if_missing=.true.)
1745 call get_param(param_file, mdl,
"WINDSTRESS_X_VAR",cs%stress_x_var, &
1746 "The name of the x-wind stress variable in WIND_FILE.", &
1748 call get_param(param_file, mdl,
"WINDSTRESS_Y_VAR", cs%stress_y_var, &
1749 "The name of the y-wind stress variable in WIND_FILE.", &
1751 call get_param(param_file, mdl,
"WIND_STAGGER",cs%wind_stagger, &
1752 "A character indicating how the wind stress components "//&
1753 "are staggered in WIND_FILE. This may be A or C for now.", &
1755 call get_param(param_file, mdl,
"WINDSTRESS_SCALE", cs%wind_scale, &
1756 "A value by which the wind stresses in WIND_FILE are rescaled.", &
1757 default=1.0, units=
"nondim")
1758 call get_param(param_file, mdl,
"USTAR_FORCING_VAR", cs%ustar_var, &
1759 "The name of the friction velocity variable in WIND_FILE "//&
1760 "or blank to get ustar from the wind stresses plus the "//&
1761 "gustiness.", default=
" ", units=
"nondim")
1762 cs%wind_file = trim(cs%inputdir) // trim(cs%wind_file)
1764 if (trim(cs%wind_config) ==
"gyres")
then
1765 call get_param(param_file, mdl,
"TAUX_CONST", cs%gyres_taux_const, &
1766 "With the gyres wind_config, the constant offset in the "//&
1767 "zonal wind stress profile: "//&
1768 " A in taux = A + B*sin(n*pi*y/L) + C*cos(n*pi*y/L).", &
1769 units=
"Pa", default=0.0)
1770 call get_param(param_file, mdl,
"TAUX_SIN_AMP",cs%gyres_taux_sin_amp, &
1771 "With the gyres wind_config, the sine amplitude in the "//&
1772 "zonal wind stress profile: "//&
1773 " B in taux = A + B*sin(n*pi*y/L) + C*cos(n*pi*y/L).", &
1774 units=
"Pa", default=0.0)
1775 call get_param(param_file, mdl,
"TAUX_COS_AMP",cs%gyres_taux_cos_amp, &
1776 "With the gyres wind_config, the cosine amplitude in "//&
1777 "the zonal wind stress profile: "//&
1778 " C in taux = A + B*sin(n*pi*y/L) + C*cos(n*pi*y/L).", &
1779 units=
"Pa", default=0.0)
1780 call get_param(param_file, mdl,
"TAUX_N_PIS",cs%gyres_taux_n_pis, &
1781 "With the gyres wind_config, the number of gyres in "//&
1782 "the zonal wind stress profile: "//&
1783 " n in taux = A + B*sin(n*pi*y/L) + C*cos(n*pi*y/L).", &
1784 units=
"nondim", default=0.0)
1785 call get_param(param_file, mdl,
"DEFAULT_2018_ANSWERS", default_2018_answers, &
1786 "This sets the default value for the various _2018_ANSWERS parameters.", &
1788 call get_param(param_file, mdl,
"WIND_GYRES_2018_ANSWERS", cs%answers_2018, &
1789 "If true, use the order of arithmetic and expressions that recover the answers "//&
1790 "from the end of 2018. Otherwise, use expressions for the gyre friction velocities "//&
1791 "that are rotationally invariant and more likely to be the same between compilers.", &
1792 default=default_2018_answers)
1794 cs%answers_2018 = .false.
1796 if (trim(cs%wind_config) ==
"scurves")
then
1797 call get_param(param_file, mdl,
"WIND_SCURVES_LATS", cs%scurves_ydata, &
1798 "A list of latitudes defining a piecewise scurve profile "//&
1799 "for zonal wind stress.", &
1800 units=
"degrees N", fail_if_missing=.true.)
1801 call get_param(param_file, mdl,
"WIND_SCURVES_TAUX", cs%scurves_taux, &
1802 "A list of zonal wind stress values at latitudes "//&
1803 "WIND_SCURVES_LATS defining a piecewise scurve profile.", &
1804 units=
"Pa", fail_if_missing=.true.)
1806 if ((trim(cs%wind_config) ==
"2gyre") .or. &
1807 (trim(cs%wind_config) ==
"1gyre") .or. &
1808 (trim(cs%wind_config) ==
"gyres") .or. &
1809 (trim(cs%buoy_config) ==
"linear"))
then
1810 cs%south_lat = g%south_lat
1811 cs%len_lat = g%len_lat
1814 call get_param(param_file, mdl,
"RHO_0", cs%Rho0, &
1815 "The mean ocean density used with BOUSSINESQ true to "//&
1816 "calculate accelerations and the mass for conservation "//&
1817 "properties, or with BOUSSINSEQ false to convert some "//&
1818 "parameters from vertical units of m to kg m-2.", &
1819 units=
"kg m-3", default=1035.0, scale=us%kg_m3_to_R)
1820 call get_param(param_file, mdl,
"RESTOREBUOY", cs%restorebuoy, &
1821 "If true, the buoyancy fluxes drive the model back "//&
1822 "toward some specified surface state with a rate "//&
1823 "given by FLUXCONST.", default= .false.)
1824 call get_param(param_file, mdl,
"LATENT_HEAT_FUSION", cs%latent_heat_fusion, &
1825 "The latent heat of fusion.", default=hlf, &
1826 units=
"J/kg", scale=us%J_kg_to_Q)
1827 call get_param(param_file, mdl,
"LATENT_HEAT_VAPORIZATION", cs%latent_heat_vapor, &
1828 "The latent heat of fusion.", default=hlv, units=
"J/kg", scale=us%J_kg_to_Q)
1829 if (cs%restorebuoy)
then
1831 call get_param(param_file, mdl,
"FLUXCONST", cs%Flux_const, &
1832 "The constant that relates the restoring surface fluxes to the relative "//&
1833 "surface anomalies (akin to a piston velocity). Note the non-MKS units.", &
1834 default=0.0, units=
"m day-1", scale=us%m_to_Z*us%T_to_s/86400.0, &
1835 unscaled=flux_const_default)
1837 if (cs%use_temperature)
then
1838 call get_param(param_file, mdl,
"FLUXCONST_T", cs%Flux_const_T, &
1839 "The constant that relates the restoring surface temperature "//&
1840 "flux to the relative surface anomaly (akin to a piston "//&
1841 "velocity). Note the non-MKS units.", &
1842 units=
"m day-1", scale=us%m_to_Z*us%T_to_s/86400.0, &
1843 default=flux_const_default)
1844 call get_param(param_file, mdl,
"FLUXCONST_S", cs%Flux_const_S, &
1845 "The constant that relates the restoring surface salinity "//&
1846 "flux to the relative surface anomaly (akin to a piston "//&
1847 "velocity). Note the non-MKS units.", &
1848 units=
"m day-1", scale=us%m_to_Z*us%T_to_s/86400.0, &
1849 default=flux_const_default)
1852 if (trim(cs%buoy_config) ==
"linear")
then
1853 call get_param(param_file, mdl,
"SST_NORTH", cs%T_north, &
1854 "With buoy_config linear, the sea surface temperature "//&
1855 "at the northern end of the domain toward which to "//&
1856 "to restore.", units=
"deg C", default=0.0)
1857 call get_param(param_file, mdl,
"SST_SOUTH", cs%T_south, &
1858 "With buoy_config linear, the sea surface temperature "//&
1859 "at the southern end of the domain toward which to "//&
1860 "to restore.", units=
"deg C", default=0.0)
1861 call get_param(param_file, mdl,
"SSS_NORTH", cs%S_north, &
1862 "With buoy_config linear, the sea surface salinity "//&
1863 "at the northern end of the domain toward which to "//&
1864 "to restore.", units=
"PSU", default=35.0)
1865 call get_param(param_file, mdl,
"SSS_SOUTH", cs%S_south, &
1866 "With buoy_config linear, the sea surface salinity "//&
1867 "at the southern end of the domain toward which to "//&
1868 "to restore.", units=
"PSU", default=35.0)
1871 call get_param(param_file, mdl,
"G_EARTH", cs%G_Earth, &
1872 "The gravitational acceleration of the Earth.", &
1873 units=
"m s-2", default = 9.80, scale=us%m_to_L**2*us%Z_to_m*us%T_to_s**2)
1875 call get_param(param_file, mdl,
"GUST_CONST", cs%gust_const, &
1876 "The background gustiness in the winds.", &
1877 units=
"Pa", default=0.0, scale=us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z)
1878 call get_param(param_file, mdl,
"FIX_USTAR_GUSTLESS_BUG", cs%fix_ustar_gustless_bug, &
1879 "If true correct a bug in the time-averaging of the gustless wind friction velocity", &
1881 call get_param(param_file, mdl,
"READ_GUST_2D", cs%read_gust_2d, &
1882 "If true, use a 2-dimensional gustiness supplied from "//&
1883 "an input file", default=.false.)
1884 if (cs%read_gust_2d)
then
1885 call get_param(param_file, mdl,
"GUST_2D_FILE", gust_file, &
1886 "The file in which the wind gustiness is found in "//&
1887 "variable gustiness.", fail_if_missing=.true.)
1888 call safe_alloc_ptr(cs%gust,g%isd,g%ied,g%jsd,g%jed)
1889 filename = trim(cs%inputdir) // trim(gust_file)
1890 call mom_read_data(filename,
'gustiness',cs%gust,g%domain, timelevel=1, &
1891 scale=us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z)
1896 if (trim(cs%wind_config) ==
"USER" .or. trim(cs%buoy_config) ==
"USER" )
then
1897 call user_surface_forcing_init(time, g, us, param_file, diag, cs%user_forcing_CSp)
1898 elseif (trim(cs%buoy_config) ==
"BFB" )
then
1899 call bfb_surface_forcing_init(time, g, us, param_file, diag, cs%BFB_forcing_CSp)
1900 elseif (trim(cs%buoy_config) ==
"dumbbell" )
then
1901 call dumbbell_surface_forcing_init(time, g, us, param_file, diag, cs%dumbbell_forcing_CSp)
1902 elseif (trim(cs%wind_config) ==
"MESO" .or. trim(cs%buoy_config) ==
"MESO" )
then
1903 call meso_surface_forcing_init(time, g, us, param_file, diag, cs%MESO_forcing_CSp)
1904 elseif (trim(cs%wind_config) ==
"ideal_hurr" .or.&
1905 trim(cs%wind_config) ==
"SCM_ideal_hurr")
then
1906 call idealized_hurricane_wind_init(time, g, us, param_file, cs%idealized_hurricane_CSp)
1907 elseif (trim(cs%wind_config) ==
"const")
then
1908 call get_param(param_file, mdl,
"CONST_WIND_TAUX", cs%tau_x0, &
1909 "With wind_config const, this is the constant zonal "//&
1910 "wind-stress", units=
"Pa", fail_if_missing=.true.)
1911 call get_param(param_file, mdl,
"CONST_WIND_TAUY", cs%tau_y0, &
1912 "With wind_config const, this is the constant meridional "//&
1913 "wind-stress", units=
"Pa", fail_if_missing=.true.)
1914 elseif (trim(cs%wind_config) ==
"SCM_CVmix_tests" .or. &
1915 trim(cs%buoy_config) ==
"SCM_CVmix_tests")
then
1916 call scm_cvmix_tests_surface_forcing_init(time, g, param_file, cs%SCM_CVmix_tests_CSp)
1919 call register_forcing_type_diags(time, diag, us, cs%use_temperature, cs%handles)
1922 call restart_init(param_file, cs%restart_CSp,
"MOM_forcing.res")
1925 call restart_init_end(cs%restart_CSp)
1927 if (
associated(cs%restart_CSp))
then
1928 call get_mom_input(dirs=dirs)
1931 if ((dirs%input_filename(1:1) ==
'n') .and. &
1932 (len_trim(dirs%input_filename) == 1)) new_sim = .true.
1933 if (.not.new_sim)
then
1934 call restore_state(dirs%input_filename, dirs%restart_input_dir, time_frc, &
1940 if (trim(cs%buoy_config) ==
"file")
then
1941 cs%SW_nlev = num_timelevels(cs%shortwave_file, cs%SW_var, min_dims=3)
1942 cs%LW_nlev = num_timelevels(cs%longwave_file, cs%LW_var, min_dims=3)
1943 cs%latent_nlev = num_timelevels(cs%latentheat_file, cs%latent_var, 3)
1944 cs%sens_nlev = num_timelevels(cs%sensibleheat_file, cs%sens_var, min_dims=3)
1946 cs%evap_nlev = num_timelevels(cs%evaporation_file, cs%evap_var, min_dims=3)
1947 cs%precip_nlev = num_timelevels(cs%rain_file, cs%rain_var, min_dims=3)
1948 cs%runoff_nlev = num_timelevels(cs%runoff_file, cs%lrunoff_var, 3)
1950 cs%SST_nlev = num_timelevels(cs%SSTrestore_file, cs%SST_restore_var, 3)
1951 cs%SSS_nlev = num_timelevels(cs%salinityrestore_file, cs%SSS_restore_var, 3)
1954 if (trim(cs%wind_config) ==
"file") &
1955 cs%wind_nlev = num_timelevels(cs%wind_file, cs%stress_x_var, min_dims=3)
1959 call user_revise_forcing_init(param_file, cs%urf_CS)
1961 call cpu_clock_end(id_clock_forcing)
1962 end subroutine surface_forcing_init
1966 subroutine surface_forcing_end(CS, fluxes)
1969 type(
forcing),
optional,
intent(inout) :: fluxes
1975 if (
present(fluxes))
call deallocate_forcing_type(fluxes)
1979 if (
associated(cs))
deallocate(cs)
1982 end subroutine surface_forcing_end