18 use data_override_mod,
only : data_override_init, data_override
20 implicit none ;
private
22 #include <MOM_memory.h>
24 public mom_wave_interface_init
25 public mom_wave_interface_init_lite
26 public update_surface_waves
28 public update_stokes_drift
30 public get_langmuir_number
50 logical,
public :: usewaves
51 logical,
public :: lagrangianmixing
56 logical,
public :: stokesmixing
62 logical,
public :: coriolisstokes
64 integer,
public :: stklevelmode=1
71 real,
allocatable,
dimension(:),
public :: &
73 real,
allocatable,
dimension(:),
public :: &
75 real,
allocatable,
dimension(:),
public :: &
77 real,
allocatable,
dimension(:),
public :: &
79 real,
allocatable,
dimension(:,:,:),
public :: &
83 real,
allocatable,
dimension(:,:,:),
public :: &
87 real,
allocatable,
dimension(:,:),
public :: &
88 la_sl,& !< SL Langmuir number (directionality factored later)
92 real,
allocatable,
dimension(:,:),
public :: &
95 real,
allocatable,
dimension(:,:),
public :: &
98 real,
allocatable,
dimension(:,:,:),
public :: &
102 real,
allocatable,
dimension(:,:,:),
public :: &
106 real,
allocatable,
dimension(:,:,:),
public :: &
110 type(time_type),
pointer,
public :: time
119 real :: la_min = 0.05
122 integer,
public :: id_surfacestokes_x = -1 , id_surfacestokes_y = -1
123 integer,
public :: id_3dstokes_x = -1 , id_3dstokes_y = -1
124 integer,
public :: id_la_turb = -1
131 integer :: wavemethod=-99
141 integer,
public :: numbands =0
145 integer,
public :: partitionmode
149 integer :: datasource
157 character(len=40) :: surfbandfilename
159 logical :: dataoverrideisinitialized
165 logical :: la_misalignment = .false.
169 #include "version_variable.h"
171 character(len=40) :: mdl =
"MOM_wave_interface"
177 integer,
parameter :: testprof = 0, surfbands = 1, &
178 dhh85 = 2, lf17 = 3, null_wavemethod=-99, &
179 dataovr = 1, coupler = 2, input = 3
182 Real :: tp_stkx0, tp_stky0, tp_wvl
183 logical :: waveagepeakfreq
184 logical :: staticwaves, dhh85_is_set
185 real :: waveage, wavewind
192 subroutine mom_wave_interface_init(time, G, GV, US, param_file, CS, diag )
193 type(time_type),
target,
intent(in) :: time
199 type(
diag_ctrl),
target,
intent(inout) :: diag
202 character*(13) :: tmpstring1,tmpstring2
203 character*(5),
parameter :: null_string =
"EMPTY"
204 character*(12),
parameter :: testprof_string =
"TEST_PROFILE"
205 character*(13),
parameter :: surfbands_string =
"SURFACE_BANDS"
206 character*(5),
parameter :: dhh85_string =
"DHH85"
207 character*(4),
parameter :: lf17_string =
"LF17"
208 character*(12),
parameter :: dataovr_string =
"DATAOVERRIDE"
209 character*(7),
parameter :: coupler_string =
"COUPLER"
210 character*(5),
parameter :: input_string =
"INPUT"
213 if (
associated(cs))
then
214 call mom_error(fatal,
"wave_interface_init called with an associated"//&
215 "control structure.")
228 dataoverrideisinitialized = .false.
237 call get_param(param_file, mdl,
"LAGRANGIAN_MIXING", cs%LagrangianMixing, &
238 "Flag to use Lagrangian Mixing of momentum", units=
"", &
240 if (cs%LagrangianMixing)
then
242 call mom_error(fatal,
"Should you be enabling Lagrangian Mixing? Code not ready.")
244 call get_param(param_file, mdl,
"STOKES_MIXING", cs%StokesMixing, &
245 "Flag to use Stokes Mixing of momentum", units=
"", &
247 if (cs%StokesMixing)
then
249 call mom_error(fatal,
"Should you be enabling Stokes Mixing? Code not ready.")
251 call get_param(param_file, mdl,
"CORIOLIS_STOKES", cs%CoriolisStokes, &
252 "Flag to use Coriolis Stokes acceleration", units=
"", &
254 if (cs%CoriolisStokes)
then
256 call mom_error(fatal,
"Should you be enabling Coriolis-Stokes? Code not ready.")
260 call get_param(param_file,mdl,
"WAVE_METHOD",tmpstring1, &
261 "Choice of wave method, valid options include: \n"// &
262 " TEST_PROFILE - Prescribed from surface Stokes drift \n"// &
263 " and a decay wavelength.\n"// &
264 " SURFACE_BANDS - Computed from multiple surface values \n"// &
265 " and decay wavelengths.\n"// &
266 " DHH85 - Uses Donelan et al. 1985 empirical \n"// &
267 " wave spectrum with prescribed values. \n"// &
268 " LF17 - Infers Stokes drift profile from wind \n"// &
269 " speed following Li and Fox-Kemper 2017.\n", &
270 units=
'', default=null_string)
271 select case (trim(tmpstring1))
273 call mom_error(fatal,
"wave_interface_init called with no specified "//&
275 case (testprof_string)
276 wavemethod = testprof
277 call get_param(param_file,mdl,
"TP_STKX_SURF",tp_stkx0,&
278 'Surface Stokes (x) for test profile',&
279 units=
'm/s',default=0.1)
280 call get_param(param_file,mdl,
"TP_STKY_SURF",tp_stky0,&
281 'Surface Stokes (y) for test profile',&
282 units=
'm/s',default=0.0)
283 call get_param(param_file,mdl,
"TP_WVL",tp_wvl,&
284 units=
'm', default=50.0, scale=us%m_to_Z)
285 case (surfbands_string)
286 wavemethod = surfbands
287 call get_param(param_file, mdl,
"SURFBAND_SOURCE",tmpstring2, &
288 "Choice of SURFACE_BANDS data mode, valid options include: \n"// &
289 " DATAOVERRIDE - Read from NetCDF using FMS DataOverride. \n"// &
290 " COUPLER - Look for variables from coupler pass \n"// &
291 " INPUT - Testing with fixed values.", &
292 units=
'', default=null_string)
293 select case (trim(tmpstring2))
295 call mom_error(fatal,
"wave_interface_init called with SURFACE_BANDS"//&
296 " but no SURFBAND_SOURCE.")
297 case (dataovr_string)
299 call get_param(param_file, mdl,
"SURFBAND_FILENAME", surfbandfilename, &
300 "Filename of surface Stokes drift input band data.", default=
"StkSpec.nc")
301 case (coupler_string)
305 call get_param(param_file,mdl,
"SURFBAND_NB",numbands, &
306 "Prescribe number of wavenumber bands for Stokes drift. "// &
307 "Make sure this is consistnet w/ WAVENUMBERS, STOKES_X, and "// &
308 "STOKES_Y, there are no safety checks in the code.", &
310 allocate( cs%WaveNum_Cen(1:numbands) )
311 cs%WaveNum_Cen(:) = 0.0
312 allocate( cs%PrescribedSurfStkX(1:numbands))
313 cs%PrescribedSurfStkX(:) = 0.0
314 allocate( cs%PrescribedSurfStkY(1:numbands))
315 cs%PrescribedSurfStkY(:) = 0.0
316 allocate( cs%STKx0(g%isdB:g%iedB,g%jsd:g%jed,1:numbands))
317 cs%STKx0(:,:,:) = 0.0
318 allocate( cs%STKy0(g%isd:g%ied,g%jsdB:g%jedB,1:numbands))
319 cs%STKy0(:,:,:) = 0.0
321 call get_param(param_file,mdl,
"SURFBAND_WAVENUMBERS",cs%WaveNum_Cen, &
322 "Central wavenumbers for surface Stokes drift bands.",units=
'rad/m', &
324 call get_param(param_file,mdl,
"SURFBAND_STOKES_X",cs%PrescribedSurfStkX, &
325 "X-direction surface Stokes drift for bands.",units=
'm/s', &
327 call get_param(param_file,mdl,
"SURFBAND_STOKES_Y",cs%PrescribedSurfStkY, &
328 "Y-direction surface Stokes drift for bands.",units=
'm/s', &
331 call mom_error(fatal,
'Check WAVE_METHOD.')
336 call mom_error(warning,
"DHH85 only ever set-up for uniform cases w/"//&
337 " Stokes drift in x-direction.")
338 call get_param(param_file,mdl,
"DHH85_AGE_FP",waveagepeakfreq, &
339 "Choose true to use waveage in peak frequency.", &
340 units=
'', default=.false.)
341 call get_param(param_file,mdl,
"DHH85_AGE",waveage, &
342 "Wave Age for DHH85 spectrum.", &
343 units=
'', default=1.2)
344 call get_param(param_file,mdl,
"DHH85_WIND",wavewind, &
345 "Wind speed for DHH85 spectrum.", &
346 units=
'', default=10.0)
347 call get_param(param_file,mdl,
"STATIC_DHH85",staticwaves, &
348 "Flag to disable updating DHH85 Stokes drift.", &
353 call mom_error(fatal,
'Check WAVE_METHOD.')
357 call get_param(param_file, mdl,
"LA_DEPTH_RATIO", la_frachbl, &
358 "The depth (normalized by BLD) to average Stokes drift over in "//&
359 "Langmuir number calculation, where La = sqrt(ust/Stokes).", &
360 units=
"nondim",default=0.04)
361 call get_param(param_file, mdl,
"LA_MISALIGNMENT", la_misalignment, &
362 "Flag (logical) if using misalignment bt shear and waves in LA",&
364 call get_param(param_file, mdl,
"MIN_LANGMUIR", cs%La_min, &
365 "A minimum value for all Langmuir numbers that is not physical, "//&
366 "but is likely only encountered when the wind is very small and "//&
367 "therefore its effects should be mostly benign.",units=
"nondim",&
372 allocate(cs%Us_x(g%isdB:g%IedB,g%jsd:g%jed,g%ke))
374 allocate(cs%Us_y(g%isd:g%Ied,g%jsdB:g%jedB,g%ke))
377 allocate(cs%US0_x(g%isdB:g%iedB,g%jsd:g%jed))
379 allocate(cs%US0_y(g%isd:g%ied,g%jsdB:g%jedB))
382 allocate(cs%La_SL(g%isc:g%iec,g%jsc:g%jec))
383 allocate(cs%La_turb(g%isc:g%iec,g%jsc:g%jec))
385 cs%La_turb (:,:) = 0.0
387 if (cs%StokesMixing)
then
388 allocate(cs%KvS(g%isd:g%Ied,g%jsd:g%jed,g%ke))
393 cs%id_surfacestokes_y = register_diag_field(
'ocean_model',
'surface_stokes_y', &
394 cs%diag%axesCu1,time,
'Surface Stokes drift (y)',
'm s-1')
395 cs%id_surfacestokes_x = register_diag_field(
'ocean_model',
'surface_stokes_x', &
396 cs%diag%axesCv1,time,
'Surface Stokes drift (x)',
'm s-1')
397 cs%id_3dstokes_y = register_diag_field(
'ocean_model',
'3d_stokes_y', &
398 cs%diag%axesCvL,time,
'3d Stokes drift (y)',
'm s-1')
399 cs%id_3dstokes_x = register_diag_field(
'ocean_model',
'3d_stokes_x', &
400 cs%diag%axesCuL,time,
'3d Stokes drift (y)',
'm s-1')
401 cs%id_La_turb = register_diag_field(
'ocean_model',
'La_turbulent',&
402 cs%diag%axesT1,time,
'Surface (turbulent) Langmuir number',
'nondim')
405 end subroutine mom_wave_interface_init
409 subroutine mom_wave_interface_init_lite(param_file)
411 character*(5),
parameter :: null_string =
"EMPTY"
412 character*(4),
parameter :: lf17_string =
"LF17"
413 character*(13) :: tmpstring1
414 logical :: statisticalwaves
417 call get_param(param_file, mdl,
"LA_DEPTH_RATIO", la_frachbl, &
418 "The depth (normalized by BLD) to average Stokes drift over in "//&
419 "Langmuir number calculation, where La = sqrt(ust/Stokes).", &
420 units=
"nondim",default=0.04)
423 call get_param(param_file,mdl,
"USE_LA_LI2016",statisticalwaves, &
424 do_not_log=.true.,default=.false.)
425 if (statisticalwaves)
then
429 wavemethod = null_wavemethod
433 end subroutine mom_wave_interface_init_lite
436 subroutine update_surface_waves(G, GV, US, Day, dt, CS)
441 type(time_type),
intent(in) :: day
442 type(time_type),
intent(in) :: dt
444 integer :: ii, jj, kk, b
445 type(time_type) :: day_center
448 day_center = day + dt/2
450 if (wavemethod == testprof)
then
452 elseif (wavemethod==surfbands)
then
453 if (datasource==dataovr)
then
454 call surface_bands_by_data_override(day_center, g, gv, us, cs)
455 elseif (datasource==coupler)
then
457 elseif (datasource==input)
then
461 cs%STKx0(ii,jj,b) = cs%PrescribedSurfStkX(b)
466 cs%STKY0(ii,jj,b) = cs%PrescribedSurfStkY(b)
474 end subroutine update_surface_waves
478 subroutine update_stokes_drift(G, GV, US, CS, h, ustar)
483 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
485 real,
dimension(SZI_(G),SZJ_(G)), &
488 real :: top, midpoint, bottom, one_cm
490 real :: cmn_fac, wn, ustokes
492 integer :: ii, jj, kk, b, iim1, jjm1
494 one_cm = 0.01*us%m_to_Z
498 if (wavemethod==testprof)
then
499 decayscale = 4.*pi / tp_wvl
500 do ii = g%isdB,g%iedB
507 midpoint = bottom - gv%H_to_Z*0.25*(h(ii,jj,kk)+h(iim1,jj,kk))
508 bottom = bottom - gv%H_to_Z*0.5*(h(ii,jj,kk)+h(iim1,jj,kk))
509 cs%Us_x(ii,jj,kk) = tp_stkx0*exp(midpoint*decayscale)
514 do jj = g%jsdB,g%jedB
520 midpoint = bottom - gv%H_to_Z*0.25*(h(ii,jj,kk)+h(ii,jjm1,kk))
521 bottom = bottom - gv%H_to_Z*0.5*(h(ii,jj,kk)+h(ii,jjm1,kk))
522 cs%Us_y(ii,jj,kk) = tp_stky0*exp(midpoint*decayscale)
529 elseif (wavemethod==surfbands)
then
535 do ii = g%isdB,g%iedB
540 if (partitionmode==0)
then
542 cmn_fac = (1.0-exp(-one_cm*2.*cs%WaveNum_Cen(b))) / &
543 (one_cm*2.*cs%WaveNum_Cen(b))
544 elseif (partitionmode==1)
then
548 cs%US0_x(ii,jj) = cs%US0_x(ii,jj) + cs%STKx0(ii,jj,b)*cmn_fac
555 midpoint = bottom - gv%H_to_Z*0.25*(h(ii,jj,kk)+h(iim1,jj,kk))
556 bottom = bottom - gv%H_to_Z*0.5*(h(ii,jj,kk)+h(iim1,jj,kk))
558 if (partitionmode==0)
then
560 cmn_fac = (exp(top*2.*cs%WaveNum_Cen(b))-exp(bottom*2.*cs%WaveNum_Cen(b)))&
561 / ((top-bottom)*(2.*cs%WaveNum_Cen(b)))
562 elseif (partitionmode==1)
then
563 if (cs%StkLevelMode==0)
then
565 cmn_fac = exp(midpoint*2.*(2.*pi*cs%Freq_Cen(b)*us%T_to_s)**2/(us%L_to_Z**2*gv%g_Earth))
566 elseif (cs%StkLevelMode==1)
then
569 wn = (2.*pi*cs%Freq_Cen(b)*us%T_to_s)**2 / (us%L_to_Z**2*gv%g_Earth)
570 cmn_fac = (exp(2.*wn*top)-exp(2.*wn*bottom)) / (2.*wn*(top-bottom))
573 cs%US_x(ii,jj,kk) = cs%US_x(ii,jj,kk) + cs%STKx0(ii,jj,b)*cmn_fac
580 do jj = g%jsdB,g%jedB
583 if (partitionmode==0)
then
585 cmn_fac = (1.0-exp(-one_cm*2.*cs%WaveNum_Cen(b))) / &
586 (one_cm*2.*cs%WaveNum_Cen(b))
587 elseif (partitionmode==1)
then
591 cs%US0_y(ii,jj) = cs%US0_y(ii,jj) + cs%STKy0(ii,jj,b)*cmn_fac
598 midpoint = bottom - gv%H_to_Z*0.25*(h(ii,jj,kk)+h(ii,jjm1,kk))
599 bottom = bottom - gv%H_to_Z*0.5*(h(ii,jj,kk)+h(ii,jjm1,kk))
601 if (partitionmode==0)
then
603 cmn_fac = (exp(top*2.*cs%WaveNum_Cen(b)) - &
604 exp(bottom*2.*cs%WaveNum_Cen(b))) / &
605 ((top-bottom)*(2.*cs%WaveNum_Cen(b)))
606 elseif (partitionmode==1)
then
607 if (cs%StkLevelMode==0)
then
609 cmn_fac = exp(midpoint*2.*(2.*pi*cs%Freq_Cen(b)*us%T_to_s)**2/(us%L_to_Z**2*gv%g_Earth))
610 elseif (cs%StkLevelMode==1)
then
613 wn = (2.*pi*cs%Freq_Cen(b)*us%T_to_s)**2 / (us%L_to_Z**2*gv%g_Earth)
614 cmn_fac = (exp(2.*wn*top)-exp(2.*wn*bottom)) / (2.*wn*(top-bottom))
617 cs%US_y(ii,jj,kk) = cs%US_y(ii,jj,kk) + cs%STKy0(ii,jj,b)*cmn_fac
622 elseif (wavemethod==dhh85)
then
623 if (.not.(staticwaves .and. dhh85_is_set))
then
624 do ii = g%isdB,g%iedB
630 midpoint = bottom - gv%H_to_Z*0.25*(h(ii,jj,kk)+h(iim1,jj,kk))
631 bottom = bottom - gv%H_to_Z*0.5*(h(ii,jj,kk)+h(iim1,jj,kk))
636 call dhh85_mid(gv, us, midpoint, ustokes)
638 cs%US_x(ii,jj,kk) = ustokes
643 do jj = g%jsdB,g%jedB
648 midpoint = bottom - gv%H_to_Z*0.25*(h(ii,jj,kk)+h(ii,jjm1,kk))
649 bottom = bottom - gv%H_to_Z*0.5*(h(ii,jj,kk)+h(ii,jjm1,kk))
656 cs%US_y(ii,jj,kk) = 0.0
664 dhh85_is_set = .true.
668 do ii = g%isdB,g%iedB
670 cs%Us_x(ii,jj,kk) = 0.
674 do jj = g%jsdB,g%jedB
675 cs%Us_y(ii,jj,kk) = 0.
686 top = h(ii,jj,1)*gv%H_to_Z
687 call get_langmuir_number( la, g, gv, us, top, ustar(ii,jj), ii, jj, &
688 h(ii,jj,:),override_ma=.false.,waves=cs)
689 cs%La_turb(ii,jj) = la
694 if (cs%id_surfacestokes_y>0) &
695 call post_data(cs%id_surfacestokes_y, cs%us0_y, cs%diag)
696 if (cs%id_surfacestokes_x>0) &
697 call post_data(cs%id_surfacestokes_x, cs%us0_x, cs%diag)
698 if (cs%id_3dstokes_y>0) &
699 call post_data(cs%id_3dstokes_y, cs%us_y, cs%diag)
700 if (cs%id_3dstokes_x>0) &
701 call post_data(cs%id_3dstokes_x, cs%us_x, cs%diag)
702 if (cs%id_La_turb>0) &
703 call post_data(cs%id_La_turb, cs%La_turb, cs%diag)
705 end subroutine update_stokes_drift
709 subroutine surface_bands_by_data_override(day_center, G, GV, US, CS)
711 type(time_type),
intent(in) :: day_center
717 real :: temp_x(SZI_(G),SZJ_(G))
718 real :: temp_y(SZI_(G),SZJ_(G))
719 real :: Top, MidPoint
722 integer,
dimension(4) :: start, counter, dims, dim_id
723 character(len=12) :: dim_name(4)
724 character(20) :: varname, varread1, varread2
725 integer :: rcode_fr, rcode_wn, ncid, varid_fr, varid_wn, id, ndims
727 if (.not.dataoverrideisinitialized)
then
728 call data_override_init(ocean_domain_in=g%Domain%mpp_domain)
729 dataoverrideisinitialized = .true.
733 varread1 =
'wavenumber'
734 varread2 =
'frequency'
735 rcode_wn = nf90_open(trim(surfbandfilename), nf90_nowrite, ncid)
736 if (rcode_wn /= 0)
then
737 call mom_error(fatal,
"error opening file "//trim(surfbandfilename)//&
738 " in MOM_wave_interface.")
742 rcode_wn = nf90_inq_varid(ncid, varread1, varid_wn)
743 rcode_fr = nf90_inq_varid(ncid, varread2, varid_fr)
745 if (rcode_wn /= 0 .and. rcode_fr /= 0)
then
746 call mom_error(fatal,
"error finding variable "//trim(varread1)//&
747 " or "//trim(varread2)//
" in file "//trim(surfbandfilename)//
" in MOM_wave_interface.")
749 elseif (rcode_wn == 0)
then
752 rcode_wn = nf90_inquire_variable(ncid, varid_wn, ndims=ndims, &
754 if (rcode_wn /= 0)
then
755 call mom_error(fatal, &
756 'error inquiring dimensions MOM_wave_interface.')
758 rcode_wn = nf90_inquire_dimension(ncid, dims(1), dim_name(1), len=id)
759 if (rcode_wn /= 0)
then
760 call mom_error(fatal,
"error reading dimension 1 data for "// &
761 trim(varread1)//
" in file "// trim(surfbandfilename)// &
762 " in MOM_wave_interface.")
764 rcode_wn = nf90_inq_varid(ncid, dim_name(1), dim_id(1))
765 if (rcode_wn /= 0)
then
766 call mom_error(fatal,
"error finding variable "//trim(dim_name(1))//&
767 " in file "//trim(surfbandfilename)//
" in MOM_wave_interace.")
770 allocate( cs%WaveNum_Cen(1:id) )
771 cs%WaveNum_Cen(:) = 0.0
772 elseif (rcode_fr == 0)
then
775 rcode_fr = nf90_inquire_variable(ncid, varid_fr, ndims=ndims, &
777 if (rcode_fr /= 0)
then
778 call mom_error(fatal,&
779 'error inquiring dimensions MOM_wave_interface.')
781 rcode_fr = nf90_inquire_dimension(ncid, dims(1), dim_name(1), len=id)
782 if (rcode_fr /= 0)
then
783 call mom_error(fatal,
"error reading dimension 1 data for "// &
784 trim(varread2)//
" in file "// trim(surfbandfilename)// &
785 " in MOM_wave_interface.")
787 rcode_fr = nf90_inq_varid(ncid, dim_name(1), dim_id(1))
788 if (rcode_fr /= 0)
then
789 call mom_error(fatal,
"error finding variable "//trim(dim_name(1))//&
790 " in file "//trim(surfbandfilename)//
" in MOM_wave_interace.")
793 allocate( cs%Freq_Cen(1:id) )
796 allocate( cs%WaveNum_Cen(1:id) )
797 cs%WaveNum_Cen(:) = 0.0
798 allocate( cs%STKx0(g%isdB:g%iedB,g%jsd:g%jed,1:id))
799 cs%STKx0(:,:,:) = 0.0
800 allocate( cs%STKy0(g%isd:g%ied,g%jsdB:g%jedB,1:id))
801 cs%STKy0(:,:,:) = 0.0
808 if (partitionmode==0)
then
809 rcode_wn = nf90_get_var(ncid, dim_id(1), cs%WaveNum_Cen, start, counter)
810 if (rcode_wn /= 0)
then
811 call mom_error(fatal,&
812 "error reading dimension 1 values for var_name "// &
813 trim(varread1)//
",dim_name "//trim(dim_name(1))// &
814 " in file "// trim(surfbandfilename)//
" in MOM_wave_interface")
817 do b = 1,numbands ; cs%WaveNum_Cen(b) = us%Z_to_m*cs%WaveNum_Cen(b) ;
enddo
818 elseif (partitionmode==1)
then
819 rcode_fr = nf90_get_var(ncid, dim_id(1), cs%Freq_Cen, start, counter)
820 if (rcode_fr /= 0)
then
821 call mom_error(fatal,&
822 "error reading dimension 1 values for var_name "// &
823 trim(varread2)//
",dim_name "//trim(dim_name(1))// &
824 " in file "// trim(surfbandfilename)//
" in MOM_wave_interface")
828 cs%WaveNum_Cen(b) = (2.*pi*cs%Freq_Cen(b)*us%T_to_s)**2 / (us%L_to_Z**2*gv%g_Earth)
838 write(varname,
"(A3,I0)")
'Usx',b
839 call data_override(
'OCN',trim(varname), temp_x, day_center)
841 write(varname,
'(A3,I0)')
'Usy',b
842 call data_override(
'OCN',trim(varname), temp_y, day_center)
844 call pass_vector(temp_x, temp_y, g%Domain, to_all, agrid)
848 if (abs(temp_x(i,j)) > 10. .or. abs(temp_y(i,j)) > 10.)
then
859 cs%STKx0(i,j,b) = 0.5 * (temp_x(i,j) + temp_x(i+1,j))
864 cs%STKy0(i,j,b) = 0.5 * (temp_y(i,j) + temp_y(i,j+1))
868 call pass_vector(cs%STKx0(:,:,b),cs%STKy0(:,:,b), g%Domain, to_all)
871 end subroutine surface_bands_by_data_override
878 subroutine get_langmuir_number( LA, G, GV, US, HBL, ustar, i, j, &
879 H, U_H, V_H, Override_MA, Waves )
883 integer,
intent(in) :: i
884 integer,
intent(in) :: j
885 real,
intent(in) :: ustar
886 real,
intent(in) :: hbl
887 logical,
optional,
intent(in) :: override_ma
891 real,
dimension(SZK_(GV)),
optional, &
893 real,
dimension(SZK_(GV)),
optional, &
895 real,
dimension(SZK_(GV)),
optional, &
900 real,
intent(out) :: la
903 real :: top, bottom, midpoint
904 real :: dpt_lasl, sheardirection, wavedirection
905 real :: la_stkx, la_stky, la_stk
906 logical :: continueloop, use_ma
907 real,
dimension(SZK_(G)) :: us_h, vs_h
908 real,
dimension(NumBands) :: stkband_x, stkband_y
912 dpt_lasl = min(-0.1*us%m_to_Z, -la_frachbl*hbl)
914 use_ma = la_misalignment
915 if (
present(override_ma)) use_ma = override_ma
919 if (.not.(
present(h).and.
present(u_h).and.
present(v_h)))
then
920 call mom_error(fatal,
'Get_LA_waves requested to consider misalignment.')
922 continueloop = .true.
926 midpoint = bottom + gv%H_to_Z*0.5*h(kk)
927 bottom = bottom + gv%H_to_Z*h(kk)
928 if (midpoint > dpt_lasl .and. kk > 1 .and. continueloop)
then
929 sheardirection = atan2(v_h(1)-v_h(kk),u_h(1)-u_h(kk))
930 continueloop = .false.
935 if (wavemethod==testprof)
then
937 us_h(kk) = 0.5*(waves%US_X(i,j,kk)+waves%US_X(i-1,j,kk))
938 vs_h(kk) = 0.5*(waves%US_Y(i,j,kk)+waves%US_Y(i,j-1,kk))
940 call get_sl_average_prof( gv, dpt_lasl, h, us_h, la_stkx)
941 call get_sl_average_prof( gv, dpt_lasl, h, vs_h, la_stky)
942 la_stk = sqrt(la_stkx*la_stkx+la_stky*la_stky)
943 elseif (wavemethod==surfbands)
then
945 stkband_x(bb) = 0.5*(waves%STKx0(i,j,bb)+waves%STKx0(i-1,j,bb))
946 stkband_y(bb) = 0.5*(waves%STKy0(i,j,bb)+waves%STKy0(i,j-1,bb))
948 call get_sl_average_band(gv, dpt_lasl, numbands, waves%WaveNum_Cen, stkband_x, la_stkx )
949 call get_sl_average_band(gv, dpt_lasl, numbands, waves%WaveNum_Cen, stkband_y, la_stky )
950 la_stk = sqrt(la_stkx**2 + la_stky**2)
951 elseif (wavemethod==dhh85)
then
954 us_h(kk) = 0.5*(waves%US_X(i,j,kk)+waves%US_X(i-1,j,kk))
955 vs_h(kk) = 0.5*(waves%US_Y(i,j,kk)+waves%US_Y(i,j-1,kk))
957 call get_sl_average_prof( gv, dpt_lasl, h, us_h, la_stkx)
958 call get_sl_average_prof( gv, dpt_lasl, h, vs_h, la_stky)
959 la_stk = sqrt(la_stkx**2 + la_stky**2)
960 elseif (wavemethod==lf17)
then
961 call get_stokessl_lifoxkemper(ustar, hbl*la_frachbl, gv, us, la_stk, la)
962 elseif (wavemethod==null_wavemethod)
then
963 call mom_error(fatal,
"Get_Langmuir_number called without defining a WaveMethod. "//&
964 "Suggest to make sure USE_LT is set/overridden to False or "//&
965 "choose a wave method (or set USE_LA_LI2016 to use statistical "//&
969 if (.not.(wavemethod==lf17))
then
976 la = max(waves%La_min, sqrt(us%Z_to_m*us%s_to_T*ustar / (la_stk+1.e-10)))
980 wavedirection = atan2(la_stky, la_stkx)
981 la = la / sqrt(max(1.e-8, cos( wavedirection - sheardirection)))
985 end subroutine get_langmuir_number
1003 subroutine get_stokessl_lifoxkemper(ustar, hbl, GV, US, UStokes_SL, LA)
1004 real,
intent(in) :: ustar
1005 real,
intent(in) :: hbl
1008 real,
intent(out) :: UStokes_SL
1009 real,
intent(out) :: LA
1012 real,
parameter :: &
1014 u19p5_to_u10 = 1.075, &
1017 fm_into_fp = 1.296, &
1019 us_to_u10 = 0.0162, &
1022 real :: UStokes, hm0, fm, fp, vstokes, kphil, kstar
1023 real :: z0, z0i, r1, r2, r3, r4, tmp, lasl_sqr_i
1026 if (ustar > 0.0)
then
1028 call ust_2_u10_coare3p5(us%Z_to_m*us%s_to_T*ustar*sqrt(us%R_to_kg_m3*gv%Rho0/1.225), u10, gv, us)
1030 ustokes = us_to_u10*u10
1034 hm0 = 0.0246 *u10**2
1037 tmp = 2.0 * pi * u19p5_to_u10 * u10
1038 fp = 0.877 * us%L_T_to_m_s**2*us%m_to_Z * gv%g_Earth / tmp
1041 fm = fm_into_fp * fp
1047 vstokes = 0.125 * pi * r_loss * fm * hm0**2
1051 kphil = 0.176 * ustokes / vstokes
1057 kstar = kphil * 2.56
1059 z0 = abs(us%Z_to_m*hbl)
1062 r1 = ( 0.151 / kphil * z0i -0.84 ) * &
1063 ( 1.0 - exp(-2.0 * kphil * z0) )
1064 r2 = -( 0.84 + 0.0591 / kphil * z0i ) * &
1065 sqrt( 2.0 * pi * kphil * z0 ) * &
1066 erfc( sqrt( 2.0 * kphil * z0 ) )
1067 r3 = ( 0.0632 / kstar * z0i + 0.125 ) * &
1068 (1.0 - exp(-2.0 * kstar * z0) )
1069 r4 = ( 0.125 + 0.0946 / kstar * z0i ) * &
1070 sqrt( 2.0 * pi *kstar * z0) * &
1071 erfc( sqrt( 2.0 * kstar * z0 ) )
1072 ustokes_sl = ustokes * (0.715 + r1 + r2 + r3 + r4)
1073 la = sqrt(us%Z_to_m*us%s_to_T*ustar / ustokes_sl)
1079 end subroutine get_stokessl_lifoxkemper
1082 subroutine get_sl_average_prof( GV, AvgDepth, H, Profile, Average )
1085 real,
intent(in) :: AvgDepth
1086 real,
dimension(SZK_(GV)), &
1088 real,
dimension(SZK_(GV)), &
1089 intent(in) :: Profile
1091 real,
intent(out) :: Average
1094 real :: top, midpoint, bottom
1105 midpoint = bottom - gv%H_to_Z * 0.5*h(kk)
1106 bottom = bottom - gv%H_to_Z * h(kk)
1107 if (avgdepth < bottom)
then
1108 sum = sum + profile(kk) * (gv%H_to_Z * h(kk))
1109 elseif (avgdepth < top)
then
1110 sum = sum + profile(kk) * (top-avgdepth)
1118 if (abs(avgdepth) <= abs(bottom))
then
1119 average = sum / abs(avgdepth)
1120 elseif (abs(bottom) > 0.0)
then
1121 average = sum / abs(bottom)
1126 end subroutine get_sl_average_prof
1129 subroutine get_sl_average_band( GV, AvgDepth, NB, WaveNumbers, SurfStokes, Average )
1132 real,
intent(in) :: AvgDepth
1133 integer,
intent(in) :: NB
1134 real,
dimension(NB), &
1135 intent(in) :: WaveNumbers
1136 real,
dimension(NB), &
1137 intent(in) :: SurfStokes
1138 real,
intent(out) :: Average
1148 average = average + surfstokes(bb) * &
1149 (1.-exp(-abs(avgdepth * 2.0 * wavenumbers(bb)))) / &
1150 abs(avgdepth * 2.0 * wavenumbers(bb))
1154 end subroutine get_sl_average_band
1162 subroutine dhh85_mid(GV, US, zpt, UStokes)
1165 real,
intent(in) :: zpt
1166 real,
intent(out) :: UStokes
1168 real :: ann, Bnn, Snn, Cnn, Dnn
1169 real :: omega_peak, omega, u10, WA, domega
1170 real :: omega_min, omega_max, wavespec, Stokes
1172 integer :: Nomega, OI
1176 g_earth = us%L_T_to_m_s**2*us%m_to_Z * gv%g_Earth
1183 domega = (omega_max-omega_min)/real(nomega)
1186 if (waveagepeakfreq)
then
1187 omega_peak = g_earth / (wa * u10)
1189 omega_peak = 2. * pi * 0.13 * g_earth / u10
1192 ann = 0.006 * waveage**(-0.55)
1194 snn = 0.08 * (1.0 + 4.0 * waveage**3)
1197 cnn = cnn - 6.0*log10(wa)
1201 omega = omega_min + 0.5*domega
1203 dnn = exp( -0.5 * (omega-omega_peak)**2 / (snn**2 * omega_peak**2) )
1205 wavespec = (ann * g_earth**2 / (omega_peak*omega**4 ) ) * &
1206 exp(-bnn*(omega_peak/omega)**4)*cnn**dnn
1208 stokes = 2.0 * wavespec * omega**3 * &
1209 exp( 2.0 * omega**2 * us%Z_to_m*zpt / g_earth) / g_earth
1210 ustokes = ustokes + stokes*domega
1211 omega = omega + domega
1215 end subroutine dhh85_mid
1219 subroutine stokesmixing(G, GV, dt, h, u, v, Waves )
1224 real,
intent(in) :: dt
1225 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),&
1227 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1229 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1234 real :: dtauup, dtaudn
1243 do i = g%iscB, g%iecB
1244 h_lay = gv%H_to_Z*0.5*(h(i,j,k)+h(i+1,j,k))
1247 dtauup = 0.5*(waves%Kvs(i,j,k)+waves%Kvs(i+1,j,k)) * &
1248 (waves%us_x(i,j,k-1)-waves%us_x(i,j,k)) / &
1249 (0.5*(h_lay + gv%H_to_Z*0.5*(h(i,j,k-1)+h(i+1,j,k-1)) ))
1252 dtaudn = 0.5*(waves%Kvs(i,j,k+1)+waves%Kvs(i+1,j,k+1)) * &
1253 (waves%us_x(i,j,k)-waves%us_x(i,j,k+1)) / &
1254 (0.5*(h_lay + gv%H_to_Z*0.5*(h(i,j,k+1)+h(i+1,j,k+1)) ))
1255 u(i,j,k) = u(i,j,k) + dt * (dtauup-dtaudn) / h_lay
1261 do j = g%jscB, g%jecB
1263 h_lay = gv%H_to_Z*0.5*(h(i,j,k)+h(i,j+1,k))
1266 dtauup = 0.5*(waves%Kvs(i,j,k)+waves%Kvs(i,j+1,k)) * &
1267 (waves%us_y(i,j,k-1)-waves%us_y(i,j,k)) / &
1268 (0.5*(h_lay + gv%H_to_Z*0.5*(h(i,j,k-1)+h(i,j+1,k-1)) ))
1271 dtaudn =0.5*(waves%Kvs(i,j,k+1)+waves%Kvs(i,j+1,k+1)) * &
1272 (waves%us_y(i,j,k)-waves%us_y(i,j,k+1)) / &
1273 (0.5*(h_lay + gv%H_to_Z*0.5*(h(i,j,k+1)+h(i,j+1,k+1)) ))
1274 v(i,j,k) = v(i,j,k) + dt * (dtauup-dtaudn) / h_lay
1279 end subroutine stokesmixing
1287 subroutine coriolisstokes(G, GV, DT, h, u, v, WAVES, US)
1292 real,
intent(in) :: dt
1293 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1295 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1297 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1308 do i = g%iscB, g%iecB
1309 dvel = 0.25*(waves%us_y(i,j+1,k)+waves%us_y(i-1,j+1,k))*g%CoriolisBu(i,j+1) + &
1310 0.25*(waves%us_y(i,j,k)+waves%us_y(i-1,j,k))*g%CoriolisBu(i,j)
1311 u(i,j,k) = u(i,j,k) + dvel*us%s_to_T*dt
1317 do j = g%jscB, g%jecB
1319 dvel = 0.25*(waves%us_x(i+1,j,k)+waves%us_x(i+1,j-1,k))*g%CoriolisBu(i+1,j) + &
1320 0.25*(waves%us_x(i,j,k)+waves%us_x(i,j-1,k))*g%CoriolisBu(i,j)
1321 v(i,j,k) = v(i,j,k) - dvel*us%s_to_T*dt
1325 end subroutine coriolisstokes
1331 subroutine ust_2_u10_coare3p5(USTair, U10, GV, US)
1332 real,
intent(in) :: USTair
1333 real,
intent(out) :: U10
1338 real,
parameter :: vonkar = 0.4
1339 real,
parameter :: nu=1e-6
1340 real :: z0sm, z0, z0rough, u10a, alpha, CD
1349 z0sm = 0.11 * nu * us%m_to_Z / ustair
1350 u10 = ustair/sqrt(0.001)
1354 do while (abs(u10a/u10-1.) > 0.001)
1357 alpha = min(0.028, 0.0017 * u10 - 0.005)
1358 z0rough = alpha * (us%m_s_to_L_T*ustair)**2 / gv%g_Earth
1360 cd = ( vonkar / log(10.*us%m_to_Z / z0) )**2
1361 u10 = ustair/sqrt(cd)
1368 u10 = ustair/sqrt(0.0015)
1374 end subroutine ust_2_u10_coare3p5
1377 subroutine waves_end(CS)
1380 if (
allocated(cs%WaveNum_Cen))
deallocate( cs%WaveNum_Cen )
1381 if (
allocated(cs%Freq_Cen))
deallocate( cs%Freq_Cen )
1382 if (
allocated(cs%Us_x))
deallocate( cs%Us_x )
1383 if (
allocated(cs%Us_y))
deallocate( cs%Us_y )
1384 if (
allocated(cs%La_SL))
deallocate( cs%La_SL )
1385 if (
allocated(cs%La_turb))
deallocate( cs%La_turb )
1386 if (
allocated(cs%STKx0))
deallocate( cs%STKx0 )
1387 if (
allocated(cs%STKy0))
deallocate( cs%STKy0 )
1388 if (
allocated(cs%KvS))
deallocate( cs%KvS )
1389 if (
allocated(cs%Us0_y))
deallocate( cs%Us0_y )
1390 if (
allocated(cs%Us0_x))
deallocate( cs%Us0_x )
1395 end subroutine waves_end