advect1dRank1 Function

private pure function advect1dRank1(self, advection_method, halowidth, directional_type) result(adv)

Computes the advective tendency for the domain instance given the desired advection method as an input function and the number of halo cells. This function works only in cases where ndirs == 1.

This implementation accepts the methods that operate on spectrum arrays of rank 1 (omnidirectional) in 1-dimensional space:

  • advectUpwind1stOrder1dRank1
  • advectCentered2ndOrder1dRank1 function with the requested advection method

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

domain instance

public pure function advection_method(f, u, dx) result(tendency)
Arguments
Type IntentOptional AttributesName
real(kind=rk), intent(in), dimension(:,:):: f
real(kind=rk), intent(in), dimension(:,:):: u
real(kind=rk), intent(in), dimension(:):: dx
Return Value real(kind=rk), dimension(:,:), allocatable
integer(kind=ik), intent(in) :: halowidth

number of halo cells to use in the advection method

integer(kind=ik), intent(in), dimension(:):: directional_type

A global constant that helps resolve the interface of this specific prodedure

Return Value type(domain_type)

Calls

proc~~advect1drank1~~CallsGraph proc~advect1drank1 advect1dRank1 setspectrumarray setspectrumarray proc~advect1drank1->setspectrumarray proc~getgroupspeed getGroupSpeed proc~advect1drank1->proc~getgroupspeed proc~getspectrumarray getSpectrumArray proc~advect1drank1->proc~getspectrumarray ub ub proc~advect1drank1->ub lb lb proc~advect1drank1->lb proc~getgridspacingxwithhalo getGridSpacingXWithHalo proc~advect1drank1->proc~getgridspacingxwithhalo proc~getgroupspeed->proc~getgroupspeed proc~getgroupspeed->ub proc~getgroupspeed->lb spectrum spectrum proc~getgroupspeed->spectrum hw hw proc~getgroupspeed->hw proc~getspectrumarray->ub proc~getspectrumarray->lb proc~getspectrumarray->spectrum proc~getspectrumarray->hw proc~getspectrum getSpectrum proc~getspectrumarray->proc~getspectrum proc~getgridspacingxwithhalo->ub proc~getgridspacingxwithhalo->lb proc~getgridspacingxwithhalo->hw
Help

Source Code


Source Code

pure type(domain_type) function advect1dRank1(self,advection_method,halowidth,&
  directional_type) result(adv)
  !! Computes the advective tendency for the domain instance given the desired
  !! advection method as an input function and the number of halo cells. This 
  !! function works only in cases where `ndirs == 1`.
  !!
  !! This implementation accepts the methods that operate on spectrum arrays
  !! of rank 1 (omnidirectional) in 1-dimensional space:
  !!
  !!   * advectUpwind1stOrder1dRank1
  !!   * advectCentered2ndOrder1dRank1
  class(domain_type),intent(in) :: self
    !! `domain` instance
  interface
    pure function advection_method(f,u,dx) result(tendency)
      import :: rk
      real(kind=rk),dimension(:,:),intent(in) :: f
      real(kind=rk),dimension(:,:),intent(in) :: u
      real(kind=rk),dimension(:),intent(in) :: dx
      real(kind=rk),dimension(:,:),allocatable :: tendency
    endfunction advection_method
  endinterface
    !! function with the requested advection method
  integer(kind=ik),intent(in) :: halowidth
    !! number of halo cells to use in the advection method
  integer(kind=ik),dimension(:),intent(in) :: directional_type
    !! A global constant that helps resolve the interface of this specific
    !! prodedure
  integer(kind=ik) :: idm
  real(kind=rk),dimension(:,:),allocatable :: f
  real(kind=rk),dimension(:,:),allocatable :: cg
  real(kind=rk),dimension(:),allocatable :: dx
  associate(lb => self % lb,ub => self % ub,hw => halowidth)
  adv = self
  idm = ub(1)-lb(1)+1+2*hw
  f = reshape(self % getSpectrumArray([hw,0],.true.),[self % nfreqs,idm])
  cg = reshape(self % getGroupSpeed([hw,0],.true.),[self % nfreqs,idm])
  dx = reshape(self % getGridSpacingXWithHalo([hw,0],.true.),[idm])
  call adv % setSpectrumArray(advection_method(f,cg,dx))
  endassociate
endfunction advect1dRank1


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