snl_DCCM2012 Function

public pure elemental function snl_DCCM2012(spectrum, sds_tendency, snl_coefficient) result(tendency)

proc~~snl_dccm2012~~UsesGraph proc~snl_dccm2012 snl_DCCM2012 module~mod_spectrum mod_spectrum module~mod_spectrum->proc~snl_dccm2012 module~mod_const mod_const module~mod_const->proc~snl_dccm2012 module~mod_const->module~mod_spectrum json_module json_module json_module->module~mod_spectrum module~mod_utility mod_utility module~mod_utility->module~mod_spectrum module~mod_precision mod_precision module~mod_precision->module~mod_spectrum module~mod_precision->module~mod_const module~mod_precision->module~mod_utility datetime_module datetime_module datetime_module->module~mod_spectrum iso_c_binding iso_c_binding iso_c_binding->module~mod_spectrum iso_fortran_env iso_fortran_env iso_fortran_env->module~mod_precision
Help

Returns a spectrum instance with the non-linear wave-wave energy transfer ($S_{nl}$) tendency formulated by Donelan et al. (2012).

References:

Donelan, M. A., M. Curcic, S. S. Chen, and A. K. Magnusson, 2012: Modeling waves and wind stress, J. Geophys. Res. Oceans, 117, C00J23, doi:10.1029/2011JC007787.

Arguments

Type IntentOptional AttributesName
type(spectrum_type), intent(in) :: spectrum

Spectrum instance

type(spectrum_type), intent(in) :: sds_tendency

Spectral dissipation tendency instance

real(kind=rk), intent(in) :: snl_coefficient

Linear coefficient of the dissipation function

Return Value type(spectrum_type)

Result tendency instance

Calls

proc~~snl_dccm2012~~CallsGraph proc~snl_dccm2012 snl_DCCM2012 getfrequency getfrequency proc~snl_dccm2012->getfrequency getspectrum getspectrum proc~snl_dccm2012->getspectrum getwavenumber getwavenumber proc~snl_dccm2012->getwavenumber getwavenumberspacing getwavenumberspacing proc~snl_dccm2012->getwavenumberspacing getdirections getdirections proc~snl_dccm2012->getdirections float float proc~snl_dccm2012->float
Help

Source Code


Source Code

pure elemental function snl_DCCM2012(spectrum,sds_tendency,snl_coefficient)&
  result(tendency)

  !! Returns a spectrum instance with the non-linear wave-wave energy transfer
  !! ($S_{nl}$) tendency formulated by Donelan et al. (2012).
  !!
  !! References:
  !!
  !! Donelan, M. A., M. Curcic, S. S. Chen, and A. K. Magnusson, 2012: Modeling
  !! waves and wind stress, *J. Geophys. Res. Oceans*, **117**, C00J23,
  !! doi:10.1029/2011JC007787.

  use mod_spectrum,only:spectrum_type
  use mod_const,only:twopi

  type(spectrum_type),intent(in) :: spectrum
    !! Spectrum instance
  type(spectrum_type),intent(in) :: sds_tendency
    !! Spectral dissipation tendency instance
  real(kind=rk),intent(in) :: snl_coefficient
    !! Linear coefficient of the dissipation function
  type(spectrum_type) :: tendency
    !! Result tendency instance

  real(kind=rk),dimension(:,:),allocatable :: sds_spectrum
  real(kind=rk),dimension(:,:),allocatable :: s_nl

  real(kind=rk),dimension(:),allocatable :: f
  real(kind=rk),dimension(:),allocatable :: th
  real(kind=rk),dimension(:),allocatable :: k
  real(kind=rk),dimension(:),allocatable :: dk

  real(kind=rk),dimension(:),allocatable :: w1,w2
  real(kind=rk) :: bf1,bf1a,bf2,dlnf

  integer :: nfreq,nfreqs
  integer :: ndir,ndirs

  f = spectrum % getFrequency()
  k = spectrum % getWavenumber()
  dk = spectrum % getWavenumberSpacing()
  th = spectrum % getDirections()

  nfreqs = size(f)
  ndirs = size(th)

  tendency = spectrum

  sds_spectrum = sds_tendency % getSpectrum()*spectrum % getSpectrum()

  allocate(s_nl(nfreqs,ndirs))
  s_nl = 0

  allocate(w1(nfreqs),w2(nfreqs))
  w1 = 0
  w2 = 0

  dlnf = (log(f(size(f)))-log(f(1)))/float(nfreqs-1)
  bf1 = exp(-16*dlnf**2)
  bf2 = exp(-64*dlnf**2)
  bf1a = bf1/(bf1+bf2)
  bf2  = bf2/(bf1+bf2)
  bf1  = bf1a

  do nfreq = 1,nfreqs-2
    w1(nfreq) = snl_coefficient*bf1*k(nfreq+1)*dk(nfreq+1)/(k(nfreq)*dk(nfreq))
    w2(nfreq) = snl_coefficient*bf2*k(nfreq+2)*dk(nfreq+2)/(k(nfreq)*dk(nfreq))
  enddo

  do concurrent(ndir=1:ndirs)
    do nfreq = 1,nfreqs-2
      s_nl(nfreq,ndir) = w1(nfreq)*sds_spectrum(nfreq+1,ndir)&
                       + w2(nfreq)*sds_spectrum(nfreq+2,ndir)&
                       - snl_coefficient*sds_spectrum(nfreq,ndir)
    enddo
  enddo

  tendency = s_nl
  deallocate(w1,w2,s_nl,sds_spectrum)

endfunction snl_DCCM2012


advect1dRank1 advect1dRank2 advect2dRank2 advectCentered2ndOrder advectCentered2ndOrder1dRank0 advectCentered2ndOrder1dRank1 advectCentered2ndOrder1dRank2 advectCentered2ndOrder2dRank0 advectCentered2ndOrder2dRank1 advectCentered2ndOrder2dRank2 advectUpwind1stOrder advectUpwind1stOrder1dRank0 advectUpwind1stOrder1dRank1 advectUpwind1stOrder1dRank2 advectUpwind1stOrder2dRank0 advectUpwind1stOrder2dRank1 advectUpwind1stOrder2dRank2 assign_array_1d assign_array_2d assign_spectrum_array_1d assign_spectrum_array_2d backward_euler backward_euler_domain backward_euler_spectrum constructor constructor constructor_1d constructor_2d diff diff_1d diff_2d diff_periodic diff_periodic_1d diff_periodic_2d domain_add_domain domain_add_real domain_div_domain domain_div_real domain_mult_domain domain_mult_real domain_sub_domain domain_sub_real domain_type domain_unary_minus donelanHamiltonHui donelanHamiltonHuiDirectionalSpectrum donelanHamiltonHuiDirectionalSpreading elevation eq eq exact_exponential exact_exponential_domain exact_exponential_spectrum forward_euler forward_euler_domain forward_euler_spectrum frequencyMoment frequencyMoment ge getAirDensity getAirDensity getAmplitude getAxisX getAxisY getCurrent_u getCurrent_u getCurrent_v getCurrent_v getDepth getDepth getDepthLevels getDirections getDirections getDirections2d getElevation getElevation getFrequency getFrequency getFrequency2d getGravity getGravity getGrid getGridRotation getGridSpacingX getGridSpacingXWithHalo getGridSpacingY getGridSpacingYWithHalo getGroupSpeed getGroupSpeed getGroupSpeed2d getLatitude getLongitude getLowerBounds getLowerBounds getPhaseSpeed getPhaseSpeed getPhaseSpeed2d getSpectrum getSpectrum getSpectrumArray getSurfaceTension getSurfaceTension getUpperBounds getUpperBounds getWaterDensity getWaterDensity getWaveAction getWavelength getWavenumber getWavenumber2d getWavenumberSpacing gravityClairaut grid_type gt horizontalAcceleration horizontalVelocity integrate integrate_domain integrate_spectrum isAllocated isAllocated isMonochromatic isOmnidirectional jonswap jonswapPeakFrequency le lt meanPeriod meanPeriod meanPeriodZeroCrossing meanPeriodZeroCrossing meanSquareSlope meanSquareSlopeDirectional momentum_x momentum_y momentumFlux_xx momentumFlux_xy momentumFlux_yy neq neq nondimensionalDepth nondimensionalEnergy nondimensionalFetch nondimensionalFrequency nondimensionalRoughness_H1986 nondimensionalRoughness_S1974 nondimensionalTime omnidirectionalSpectrum ones ones_int ones_real peakedness peakFrequency peakFrequencyDiscrete phillips piersonMoskowitz piersonMoskowitzPeakFrequency pressure range range_int range_real readJSON real2d_mult_spectrum real_add_domain real_add_spectrum real_div_domain real_div_spectrum real_mult_domain real_mult_spectrum real_sub_domain real_sub_spectrum saturationSpectrum sbf_DCCM2012 sbf_JONSWAP sds_DCCM2012 sdt_DCCM2012 setAirDensity setAirDensity setCurrent1d setCurrent2d setDepth setDepth setElevation setElevation setGravity setGravity setSpectrum1d setSpectrum1d setSpectrum2d setSpectrum2d setSpectrumArray1d1d setSpectrumArray1d2d setSpectrumArray2d2d setSurfaceTension setSurfaceTension setWaterDensity setWaterDensity significantSurfaceOrbitalVelocity significantWaveHeight significantWaveHeight sin_DCCM2012 snl_DCCM2012 spectrum_add_real spectrum_add_spectrum spectrum_div_real spectrum_div_spectrum spectrum_mult_real spectrum_mult_real2d spectrum_mult_spectrum spectrum_sub_real spectrum_sub_spectrum spectrum_type spectrum_unary_minus stokesDrift stokesDrift2d tile tile_1d_int tile_1d_real tile_2d_int tile_2d_real tile_3d_int tile_3d_real ursellNumber verticalAcceleration verticalVelocity waveAge wavenumber wavenumberMoment wavenumberMoment wavenumberSpectrum writeJSON writeJSON zeros zeros_int zeros_real