Package 'TrenchR'

Title: Tools for Microclimate and Biophysical Ecology
Description: Tools for translating environmental change into organismal response. Microclimate models to vertically scale weather station data to organismal heights. The biophysical modeling tools include both general models for heat flows and specific models to predict body temperatures for a variety of ectothermic taxa. Additional functions model and temporally partition air and soil temperatures and solar radiation. Utility functions estimate the organismal and environmental parameters needed for biophysical ecology. 'TrenchR' focuses on relatively simple and modular functions so users can create transparent and flexible biophysical models. Many functions are derived from Gates (1980) <doi:10.1007/978-1-4612-6024-0> and Campbell and Norman (1988) <isbn:9780387949376>.
Authors: Lauren Buckley [aut, cre] , Bryan Briones Ortiz [aut], Isaac Caruso [aut], Aji John [aut] , Ofir Levy [aut] , Abigail Meyer [aut], Eric Riddell [aut] , Yutaro Sakairi [aut], Juniper Simonis [aut] , Brian Helmuth [ctb]
Maintainer: Lauren Buckley <[email protected]>
License: MIT + file LICENSE
Version: 1.1.1
Built: 2024-11-12 05:58:19 UTC
Source: https://github.com/trenchproject/trenchr

Help Index


Actual Vapor Pressure from Dewpoint Temperature

Description

The function calculates actual vapor pressure from dewpoint temperature based on Stull (2000); Riddell et al. (2018).

Usage

actual_vapor_pressure(T_dewpoint)

Arguments

T_dewpoint

numeric dewpoint temperature (K).

Value

numeric actual vapor pressure (kPa).

Author(s)

Eric Riddell

References

Riddell EA, Odom JP, Damm JD, Sears MW (2018). “Plasticity reveals hidden resistance to extinction under climate change in the global hotspot of salamander diversity.” Science Advances, 4(4). doi:10.1126/sciadv.aar5471.

Stull RB (2000). Meteorology for Scientists and Engineers. Brooks Cole. ISBN 978-0534372149.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

actual_vapor_pressure(T_dewpoint = 293)

Air Temperature Profile using MICRO Routine

Description

The function estimates air temperature (C) at a specified height (m). Estimates a single, unsegmented temperature profile using the MICRO routine from NicheMapR (Kearney and Porter 2017).

Usage

air_temp_profile(T_r, u_r, zr, z0, z, T_s)

Arguments

T_r

numeric air temperature (C) at reference height.

u_r

numeric windspeed (m s-1) at reference height.

zr

numeric initial reference height (m).

z0

numeric surface roughness (m).

z

numeric height to scale (m).

T_s

numeric surface temperatures (C).

Value

numeric air temperature (C).

References

Kearney MR, Porter WP (2017). “NicheMapR - an R package for biophysical modelling: the microclimate model.” Ecography, 40, 664-674. doi:10.1111/ecog.02360.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), degree_days(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), monthly_solar_radiation(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

air_temp_profile(T_r = 20,
                   u_r = 0.1, 
                   zr  = 0.1, 
                   z0  = 0.2, 
                   z   = 0.15, 
                   T_s = 25)

Air Temperature at a Specified Height Under Neutral Conditions

Description

The function calculates air temperature (C) at a specified height (m) within a boundary layer near the surface. The velocity profile is the neutral profile described by Sellers (1965). Function is included as equations (2) and (3) of Porter et al. (1973).

Usage

air_temp_profile_neutral(T_r, zr, z0, z, T_s)

Arguments

T_r

numeric air temperature (C) at reference height.

zr

numeric initial reference height (m).

z0

numeric surface roughness (m).

z

numeric height to scale to (m).

T_s

numeric surface temperatures (C).

Value

numeric air temperature (C).

References

Porter WP, Mitchell JW, Bekman A, DeWitt CB (1973). “Behavioral implications of mechanistic ecology: thermal and behavioral modeling of desert ectotherms and their microenvironments.” Oecologia, 13, 1-54.

Sellers WD (1965). Physical climatology. University of Chicago Press, Chicago, IL, USA.

See Also

Other microclimate functions: air_temp_profile_segment(), air_temp_profile(), degree_days(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), monthly_solar_radiation(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

air_temp_profile_neutral(T_r = 20, 
                           zr  = 0.1, 
                           z0  = 0.2, 
                           z   = 0.15, 
                           T_s = 25)

Air Temperature at a Specified Height

Description

The function calculates air temperature (C) at a specified height (m). Estimates a three segment velocity and temperature profile based on user-specified, experimentally determined values for 3 roughness heights and reference heights. Multiple heights are appropriate in heterogenous areas with, for example, a meadow, bushes, and rocks. Implements the MICROSEGMT routine from NicheMapR as described in Kearney and Porter (2017).

Usage

air_temp_profile_segment(T_r, u_r, zr, z0, z, T_s)

Arguments

T_r

numeric a vector of air temperatures (C) at the 3 reference heights.

u_r

numeric a vector of wind speeds (m s-1) at the 3 reference heights.

zr

numeric a vector of 3 reference heights (meters).

z0

numeric a vector of 3 experimentally determined roughness heights (meters).

z

numeric height for air temperature estimation (meters).

T_s

numeric surface temperatures (C).

Value

numeric air temperature (C).

References

Kearney MR, Porter WP (2017). “NicheMapR - an R package for biophysical modelling: the microclimate model.” Ecography, 40, 664-674. doi:10.1111/ecog.02360.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile(), degree_days(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), monthly_solar_radiation(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

air_temp_profile_segment(T_r = c(25, 22, 20), 
                           u_r = c(0.01, 0.025, 0.05), 
                           zr  = c(0.05, 0.25, 0.5), 
                           z0  = c(0.01, 0.15, 0.2), 
                           z   = 0.3, 
                           T_s = 27)

Air Pressure from Elevation

Description

The function estimates air pressure (kPa) as a function of elevation (Engineering ToolBox 2003).

Usage

airpressure_from_elev(elev)

Arguments

elev

numeric elevation (meters).

Value

numeric air pressure (kPa).

References

Engineering ToolBox (2003). Atmospheric Pressure vs. Elevation above Sea Level. https://www.engineeringtoolbox.com/air-altitude-pressure-d_462.html.

See Also

Other utility functions: azimuth_angle(), day_of_year(), daylength(), dec_angle(), solar_noon(), temperature conversions, zenith_angle()

Examples

airpressure_from_elev(elev = 1500)

Convert Angles Between Radians and Degrees

Description

The function converts angles in radians to degrees or degrees to radians.

Usage

radians_to_degrees(rad)

degrees_to_radians(deg)

Arguments

rad

numeric angle (radians).

deg

numeric angle (degrees).

Value

numeric angle (degrees or radians).

Examples

radians_to_degrees(0.831)
  degrees_to_radians(47.608)

Azimuth Angle

Description

The function calculates the azimuth angle, the angle (degrees) from which the sunlight is coming measured from true north or south measured in the horizontal plane. The azimuth angle is measured with respect to due south, increasing in the counter clockwise direction so 90 degrees is east (Campbell and Norman 1998).

Usage

azimuth_angle(doy, lat, lon, hour, offset = NA)

Arguments

doy

numeric day of year (1-366). This can be obtained from standard date via day_of_year.

lat

numeric latitude (decimal degrees).

lon

numeric longitude (decimal degrees).

hour

numeric hour of the day.

offset

numeric number of hours to add to UTC (Coordinated Universal Time) to get local time (to improve accuracy but not always necessary). Optional. Defaults to NA.

Value

numeric azimuth angle (degrees).

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

See Also

Other utility functions: airpressure_from_elev(), day_of_year(), daylength(), dec_angle(), solar_noon(), temperature conversions, zenith_angle()

Examples

azimuth_angle(doy    = 112, 
                lat    = 47.61, 
                lon    = -122.33, 
                hour   = 12, 
                offset = -8)

Boundary Layer Resistance

Description

The function estimates boundary layer resistance under free convection based on the function in Riddell et al. (2018).

Usage

boundary_layer_resistance(T_a, e_s, e_a, elev, D, u = NA)

Arguments

T_a

numeric air temperature (K).

e_s

numeric saturation vapor pressure (kPa).

e_a

numeric actual vapor pressure (kPa).

elev

numeric elevation (m).

D

numeric characteristic dimension (e.g., body diameter) (m).

u

numeric wind speed (m s-1), if not provided assume free convection; if provided, use forced convection if appropriate.

Value

numeric boundary layer resistance (s cm-1).

Author(s)

Eric Riddell

References

Riddell EA, Odom JP, Damm JD, Sears MW (2018). “Plasticity reveals hidden resistance to extinction under climate change in the global hotspot of salamander diversity.” Science Advances, 4(4). doi:10.1126/sciadv.aar5471.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

boundary_layer_resistance(T_a  = 293,
                            e_s  = 2.5,
                            e_a  = 2.0,
                            elev = 500,
                            D    = 0.007,
                            u    = 2)

General Use Constants

Description

Basic functions for numerical constants for conversions.

Usage

specific_heat_h2o(units = "J_kg-1_K-1")

latent_heat_vaporization_h2o(units = "J_kg-1")

stefan_boltzmann_constant(units = "W_m-2_K-4")

von_karman_constant(units = "")

Arguments

units

character indicating units

  • specific_heat_h2o:

    • "J_kg-1_K-1": J kg-1 K-1

  • latent_heat_vaporization_h2o:

    • "J_kg-1": J kg-1

  • stefan_boltzmann_constant:

    • "W_m-2_K-4": W m-2 K-4

    • "mW_cm-2_K-4": mW cm-2 K-4

  • con_karman_constant:

    • "": dimensionless

Value

numeric values in units.

Examples

specific_heat_h2o()
  latent_heat_vaporization_h2o()
  stefan_boltzmann_constant()

Calendar Day from Date

Description

The function converts a date (day, month, year) to Calendar Day (day of year).

Usage

day_of_year(day, format = "%Y-%m-%d")

Arguments

day

character numerical date in standard format (e.g. "2017-01-02", "01-02", "01/02/2017" etc).

format

character date format following as.POSIXlt conventions. Default value = "%Y-%m-%d".

Value

numeric Julian day number, 1-366 (e.g. 1 for January 1st).

See Also

Other utility functions: airpressure_from_elev(), azimuth_angle(), daylength(), dec_angle(), solar_noon(), temperature conversions, zenith_angle()

Examples

day_of_year(day    = "2017-04-22", 
              format = "%Y-%m-%d")
  day_of_year(day    = "2017-04-22")
  day_of_year(day    = "04/22/2017", 
              format = "%m/%d/%Y")

Day Length

Description

The function calculates daylength in hours as a function of latitude and day of year. Uses the CMB model (Campbell and Norman 1998).

Usage

daylength(lat, doy)

Arguments

lat

numeric latitude (decimal degrees).

doy

numeric day of year (1-366). This can be obtained from standard date via day_of_year.

Value

numeric day length (hours).

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

See Also

Other utility functions: airpressure_from_elev(), azimuth_angle(), day_of_year(), dec_angle(), solar_noon(), temperature conversions, zenith_angle()

Examples

daylength(lat = 47.61, 
            doy = 112)

Solar Declination in Radians

Description

The function calculates solar declination, which is the angular distance of the sun north or south of the earth’s equator, based on the day of year (Campbell and Norman 1998).

Usage

dec_angle(doy)

Arguments

doy

numeric day of year (1-366). This can be obtained from standard date via day_of_year.

Value

numeric declination angle (radians).

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

See Also

Other utility functions: airpressure_from_elev(), azimuth_angle(), day_of_year(), daylength(), solar_noon(), temperature conversions, zenith_angle()

Examples

dec_angle(doy = 112)
  dec_angle(doy = 360)

Degree Days

Description

The function calculates degree days using the following approximations: single or double sine wave, single or double triangulation (University of California Integrated Pest Management Program 2016). Double approximation methods assume symmetry, such that a day's thermal minimum is equal to that of the previous day. Double sine wave approximation of degree days from Allen (1976).

Usage

degree_days(T_min, T_max, LDT = NA, UDT = NA, method = "single.sine")

Arguments

T_min

numeric Minimum temperature of the day (C).

T_max

numeric Maximum temperature of the day (C).

LDT

numeric lower developmental threshold (C).

UDT

numeric upper developmental threshold (C).

method

character type of method being used. Current choices: "single.sine", "double.sine", "single.triangulation", and "double.triangulation".

Value

numeric degree days (C).

References

Allen JC (1976). “A Modified Sine Wave Method for Calculating Degree Days.” Environmental Entomology, 5(3), 388-396. doi:10.1093/ee/5.3.388.

University of California Integrated Pest Management Program (2016). Degree Days: Methods. https://ipm.ucanr.edu/WEATHER/ddfigindex.html.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), air_temp_profile(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), monthly_solar_radiation(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

degree_days(T_min  = 7, 
              T_max  = 14, 
              LDT    = 12, 
              UDT    = 33, 
              method = "single.sine")
  degree_days(T_min  = 7, 
              T_max  = 14, 
              LDT    = 12, 
              UDT    = 33, 
              method = "single.triangulation")

Direct Solar Radiation

Description

The function estimates direct solar radiation (W/m2) based on latitude, day of year, elevation, and time. The function uses two methods (McCullough and Porter 1971; Campbell and Norman 1998) compiled in Tracy et al. (1983).

Usage

direct_solar_radiation(lat, doy, elev, t, t0, method = "Campbell 1977")

Arguments

lat

numeric latitude (degrees).

doy

numeric day of year (1-366).

elev

numeric elevation (m).

t

numeric local time (decimal hours).

t0

numeric time of local noon (decimal hours), can be estimated using solar_noon.

method

character method for estimating direct solar radiation, options: "Campbell 1977" (default), "Gates 1962".

Value

numeric direct solar radiation (W/m2).

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

McCullough EC, Porter WP (1971). “Computing Clear Day Solar Radiation Spectra for the Terrestrial Ecological Environment.” Ecology, 52(6), 1008-1015. doi:10.2307/1933806.

Tracy CR, Hammond KA, Lechleitner RA, II WJS, Thompson DB, Whicker AD, Williamson SC (1983). “Estimating clear-day solar radiation: an evaluation of three models.” Journal of Thermal Biology, 8(3), 247-251. doi:10.1016/0306-4565(83)90003-7.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), air_temp_profile(), degree_days(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), monthly_solar_radiation(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

direct_solar_radiation(lat    = 47.61, 
                         doy    = 112, 
                         elev   = 1500, 
                         t      = 9, 
                         t0     = 12, 
                         method = "Campbell 1977")

Hourly Solar Radiation

Description

The function estimates hourly solar radiation (W m-2 hr-1) as a function of daily global solar radiation (W m-2 d-1). Based on Tham et al. (2010) and Tham et al. (2011).

Usage

diurnal_radiation_variation(doy, S, hour, lon, lat)

Arguments

doy

numeric the day of year.

S

numeric solar radiation (W m-2 d-1).

hour

numeric hour (0-24).

lon

numeric longitude (degrees).

lat

numeric latitude (degrees).

Value

numeric hourly solar radiation (W m-2).

References

Tham Y, Muneer T, Davison B (2010). “Estimation of hourly averaged solar irradiation: evaluation of models.” Building Services Engingeering Research Technology, 31(1). doi:10.1177/0143624409350547.

Tham Y, Muneer T, Davison B (2011). “Prediction of hourly solar radiation on horizontal and inclined surfaces for Muscat/Oman.” The Journal of Engineering Research, 8(2), 19-31. doi:10.24200/tjer.vol8iss2pp19-31.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), air_temp_profile(), degree_days(), direct_solar_radiation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), monthly_solar_radiation(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

diurnal_radiation_variation(doy    = 112, 
                            S = 8000, 
                            hour   = 12, 
                            lon    = -122.33, 
                            lat    = 47.61)

Hourly Temperature Variation assuming a Sine Interpolation

Description

The function estimates temperature for a specified hour using the sine interpolation in Campbell and Norman (1998).

Usage

diurnal_temp_variation_sine(T_max, T_min, t)

Arguments

T_max, T_min

numeric maximum and minimum daily temperatures (C).

t

numeric time for temperature estimate (hour).

Value

numeric temperature (C) at a specified hour.

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), air_temp_profile(), degree_days(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), monthly_solar_radiation(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

diurnal_temp_variation_sine(T_max = 30, 
                              T_min = 10, 
                              t     = 11)

Hourly Temperature Variation assuming Sine and Exponential Components

Description

The function estimates temperature across hours using a diurnal temperature variation function incorporating sine and exponential components (Parton and Logan 1981).

Usage

diurnal_temp_variation_sineexp(
  T_max,
  T_min,
  t,
  t_r,
  t_s,
  alpha = 2.59,
  beta = 1.55,
  gamma = 2.2
)

Arguments

T_max, T_min

numeric maximum and minimum daily temperatures (C).

t

numeric time for temperature estimate (hour).

t_r, t_s

numeric times of sunrise and sunset (hour).

alpha

numeric time difference between t_x (time of maximum temperature) and noon (hour).

beta

numeric time difference between t_x and sunrise (hour).

gamma

numeric decay parameter for rate of t change from sunset to t_n (time of minimum temp).

Details

Default alpha, beta, and gamma values are the average of 5 North Carolina sites (Wann et al. 1985).

Other alpha, beta, and gamma parameterizations include values for Denver, Colorado from Parton and Logan (1981):

  • 150 cm air temperature: alpha = 1.86, beta = 2.20, gamma = -0.17

  • 10 cm air temperature: alpha = 1.52, beta = 2.00, gamma = -0.18

  • soil surface temperature: alpha = 0.50, beta = 1.81, gamma = 0.49

  • 10cm soil temperature: alpha = 0.45, beta = 2.28, gamma = 1.83

Value

numeric temperature (C) at a specified hour.

References

Parton WJ, Logan JA (1981). “A model for diurnal variation in soil and air temperature.” Agricultural Meteorology, 23, 205-216. doi:10.1016/0002-1571(81)90105-9, https://www.sciencedirect.com/science/article/abs/pii/0002157181901059.

Wann M, Yen D, Gold HJ (1985). “Evaluation and calibration of three models for daily cycle of air temperature.” Agricultural and Forest Meteorology, 34, 121-128. https://doi.org/10.1016/0168-1923(85)90013-9.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), air_temp_profile(), degree_days(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), monthly_solar_radiation(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

diurnal_temp_variation_sineexp(T_max = 30, 
                                 T_min = 10, 
                                 t     = 11, 
                                 t_r   = 6, 
                                 t_s   = 18, 
                                 alpha = 2.59, 
                                 beta  = 1.55, 
                                 gamma = 2.2)

Hourly Temperature Variation using Sine and Square Root Functions

Description

The function estimates temperature for a specified hour using sine and square root functions (Cesaraccio et al. 2001).

Usage

diurnal_temp_variation_sinesqrt(t, t_r, t_s, T_max, T_min, T_minp)

Arguments

t

numeric hour or hours for temperature estimate.

t_r, t_s

numeric sunrise and sunset hours (0-23).

T_max, T_min

numeric maximum and minimum temperatures of current day (C).

T_minp

numeric minimum temperature of following day (C).

Value

numeric temperature (C) at a specified hour.

References

Cesaraccio C, Spano D, Duce P, Snyder RL (2001). “An improved model for determining degree-day values from daily temperature data.” International Journal of Biometeorology, 45, 161-169. doi:10.1007/s004840100104.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), air_temp_profile(), degree_days(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sine(), monthly_solar_radiation(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

diurnal_temp_variation_sinesqrt(t      = 8, 
                                  t_r    = 6, 
                                  t_s    = 18, 
                                  T_max  = 30, 
                                  T_min  = 10, 
                                  T_minp = 12)

Properties of Dry Air

Description

The function calculates several properties of dry air and related characteristics shown as output variables below. The program is based on equations from List (1971) and code implementation from NicheMapR (Kearney and Porter 2017; Kearney and Porter 2020).

The user must supply values for the input variables db, bp, and alt. If alt is known (-1000 < alt < 20000) but not BP, then set bp = 0.

Usage

DRYAIR(db, bp = 0, alt = 0)

Arguments

db

numeric dry bulb temperature (C).

bp

numeric barometric pressure (Pa).

alt

numeric altitude (m).

Value

Named list with elements

  • patmos: numeric standard atmospheric pressure (Pa)

  • density: numeric density (kg m-3)

  • visdyn: numeric dynamic viscosity (kg m-1 s-1)

  • viskin: numeric kinematic viscosity (m2 s-1)

  • difvpr: numeric diffusivity of water vapor in air (m2 s-1)

  • thcond: numeric thermal conductivity (W K-1 m-1)

  • htovpr: numeric latent heat of vaporization of water (J kg-1)

  • tcoeff: numeric temperature coefficient of volume expansion (K-1)

  • ggroup: numeric group of variables in Grashof number (m-3 K-1)

  • bbemit: numeric black body emittance (W m-2)

  • emtmax: numeric wave length of maximum emittance (m)

References

Kearney MR, Porter WP (2017). “NicheMapR - an R package for biophysical modelling: the microclimate model.” Ecography, 40, 664-674. doi:10.1111/ecog.02360.

Kearney MR, Porter WP (2020). “NicheMapR - an R package for biophysical modelling: the ectotherm and Dynamic Energy Budget models.” Ecography, 43(1), 85-96. doi:10.1111/ecog.04680.

List RJ (1971). “Smithsonian Meteorological Tables.” Smithsonian Miscellaneous Collections, 114(1), 1-527. https://repository.si.edu/handle/10088/23746.

Examples

DRYAIR(db = 30, 
         bp = 100*1000, 
         alt = 0)

External Resistance to Water Vapor Transfer

Description

The function estimate external resistance to water vapor transfer using the Lewis rule relating heat and mass transport (Spotila et al. 1992)

Usage

external_resistance_to_water_vapor_transfer(H, ecp = 12000)

Arguments

H

numeric heat transfer (convection) coefficient (W m-2 K-1).

ecp

numeric aggregate parameter (J m-3 K-1) that is the product of the density of air (kg m-3) and the specific heat of air at constant pressure (J kg-1 K-1). Default of 12000 (J m-3 K-1) commonly assumed.

Value

numeric external resistance to water vapor transfer (s m-1).

References

Spotila JR, Feder ME, Burggren WW (1992). “Biophysics of Heat and Mass Transfer.” Environmental Physiology of the Amphibians. https://press.uchicago.edu/ucp/books/book/chicago/E/bo3636401.html.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

external_resistance_to_water_vapor_transfer(H = 20)

Determine if Convection is Free or Forced

Description

The function compares the Grashof and Reynolds numbers to determine whether convection is free or forced (Gates 1980).

Usage

free_or_forced_convection(Gr, Re)

Arguments

Gr

numeric Grashof Number (dimensionless).

Re

numeric Reynolds Number (dimensionless).

Value

character "free", "forced", or "intermediate".

References

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

free_or_forced_convection(Gr = 100,
                            Re = 5)

Grashof Number

Description

The function estimates the Grashof Number, which describes the ability of a parcel of fluid warmer or colder than the surrounding fluid to rise against or fall with the attractive force of gravity. The Grashof Number is estimated as the ratio of a buoyant force times an inertial force to the square of a viscous force (Campbell and Norman 1998).

Usage

Grashof_number(T_a, T_g, D, nu)

Arguments

T_a

numeric Air temperature (C).

T_g

numeric ground (surface) temperature (C).

D

numeric characteristic dimension (e.g., body diameter) (meters).

nu

numeric the kinematic viscosity, ratio of dynamic viscosity to density of the fluid (m2 s-1); can calculate from DRYAIR or WETAIR.

Value

numeric Grashof number.

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

See Also

Other biophysical models: Grashof_number_Gates(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Grashof_number(T_a = 30,
                 T_g = 35,
                 D  = 0.001,
                 nu = 1.2)

Grashof Number as in Gates (1980)

Description

The function estimates the Grashof Number, which describes the ability of a parcel of fluid warmer or colder than the surrounding fluid to rise against or fall with the attractive force of gravity (Gates 1980). The Grashof Number is estimated as the ratio of a buoyant force times an inertial force to the square of a viscous force.

Usage

Grashof_number_Gates(T_a, T_g, beta, D, nu)

Arguments

T_a

numeric Air temperature (C).

T_g

numeric Ground (surface) temperature (C).

beta

numeric coefficient of volumetric thermal expansion, beta = 3.67 x 10-3 C-1 in air and 41.9 x 10-4 C-1 in water.

D

numeric is characteristic dimension (e.g., body diameter) (m)

nu

numeric is the kinematic viscosity, the ratio of dynamic viscosity to density of the fluid (m2 s-1); can calculate from DRYAIR or WETAIR.

Value

numeric Grashof number.

References

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

See Also

Other biophysical models: Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Grashof_number_Gates(T_a   = 30,
                       T_g   = 35,
                       beta = 0.00367,
                       D    = 0.001,
                       nu   = 1.2)

Estimate the Heat Transfer Coefficient Empirically

Description

The function estimates the heat transfer coefficient for various taxa based on empirical measurements (Mitchell 1976).

Usage

heat_transfer_coefficient(u, D, K, nu, taxon = "cylinder")

Arguments

u

numeric air velocity (m s-1).

D

numeric characteristic dimension (e.g., diameter or snout-vent length) (meters).

K

numeric thermal conductivity of air (W K-1 m-1), can calculate using DRYAIR or WETAIR.

nu

numeric kinematic viscosity of air (m2 s-1), can calculate using DRYAIR or WETAIR.

taxon

character which class of organism, current choices: "sphere", "cylinder", "frog", "lizard_surface", "lizard_elevated", "flyinginsect", "spider".
"cylinder" assumes 40 < Re < 4000. "lizard_surface" and "lizard_elevated" assume the lizard is prostrate on and elevated above the surface, respectively. The values are the average for lizards parallel and perpendicular to the air flow.

Value

numeric heat transfer coefficient, H_L (W K-1 m-2).

References

Mitchell JW (1976). “Heat transfer from spheres and other animal forms.” Biophysical Journal, 16(6), 561-569. ISSN 0006-3495, doi:10.1016/S0006-3495(76)85711-6.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

heat_transfer_coefficient(u     = 0.5,
                            D     = 0.05,
                            K     = 25.7 * 10^(-3),
                            nu    = 15.3 * 10^(-6),
                            taxon = "cylinder")

Estimate the Heat Transfer Coefficient Using a Spherical Approximation

Description

The function estimates the heat transfer coefficient for various taxa. Approximates forced convective heat transfer for animal shapes using the convective relationship for a sphere (Mitchell 1976).

Usage

heat_transfer_coefficient_approximation(u, D, K, nu, taxon = "sphere")

Arguments

u

numeric air velocity (m s-1).

D

numeric characteristic dimension (e.g., diameter or snout-vent length) (meters).

K

numeric thermal conductivity of air (W m-2 K-1), can calculate using DRYAIR or WETAIR.

nu

numeric kinematic Viscosity of air (m2 s-1), can calculate using DRYAIR or WETAIR.

taxon

character of which class of organism, current choices: "sphere", "frog", "lizard", "flyinginsect", "spider".

Value

numeric heat transfer coefficient, H_L (W m-2 K-1).

References

Mitchell JW (1976). “Heat transfer from spheres and other animal forms.” Biophysical Journal, 16(6), 561-569. ISSN 0006-3495, doi:10.1016/S0006-3495(76)85711-6.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

heat_transfer_coefficient_approximation(u     = 3,
                                         D     = 0.05,
                                         K     = 25.7 * 10^(-3),
                                         nu    = 15.3 * 10^(-6),
                                         taxon = "sphere")

Estimate the Heat Transfer Coefficient using Simple Relationships

Description

The function estimates the heat transfer coefficient (Mitchell 1976) using either the relationship in Spotila et al. (1992) or that in Gates (1980).

Usage

heat_transfer_coefficient_simple(u, D, type)

Arguments

u

numeric air velocity (m s-1).

D

numeric characteristic dimension (e.g., diameter or snout-vent length) (m).

type

character choice between "Spotila" and "Gates" for equation to use.

Value

numeric heat transfer coefficient, H_L (W m-2 K-1).

References

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

Mitchell JW (1976). “Heat transfer from spheres and other animal forms.” Biophysical Journal, 16(6), 561-569. ISSN 0006-3495, doi:10.1016/S0006-3495(76)85711-6.

Spotila JR, Feder ME, Burggren WW (1992). “Biophysics of Heat and Mass Transfer.” Environmental Physiology of the Amphibians. https://press.uchicago.edu/ucp/books/book/chicago/E/bo3636401.html.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

heat_transfer_coefficient_simple(u    = 0.5,
                                   D    = 0.05,
                                   type = "Gates")

Organism Mass from Length

Description

The function estimates mass (g) from length (m) for a variety of taxa.

Usage

mass_from_length(l, taxon)

Arguments

l

numeric vector of length (m). Can be 1 or more values.
Snout-vent length is used for amphibians and reptiles, except turtles where length is carapace length.

taxon

character taxon of organism, current choices: "insect", "lizard", "salamander", "frog", "snake", "turtle".

Details

All models follow (m = a lb) with mass in grams and length in meters.

  • Lizards: Meiri (2010):

    a=16368.17a = 16368.17
    b=3.022b = 3.022

  • Salamanders: Pough (1980):

    a=13654.4a = 13654.4
    b=2.94b = 2.94

  • Frogs: Pough (1980):

    a=181197.1a = 181197.1
    b=3.24b = 3.24

  • Snakes: Pough (1980):

    a=723.6756a = 723.6756
    b=3.02b = 3.02

  • Turtles: Pough (1980):

    a=93554.48a = 93554.48
    b=2.69b = 2.69

  • Insects: Sample et al. (1993):

    a=806.0827a = 806.0827
    b=2.494b = 2.494

Value

numeric mass (g).

References

Meiri S (2010). “Length - weight allometries in lizards.” Journal of Zoology, 281(3), 218-226. doi:10.1111/j.1469-7998.2010.00696.x.

Pough FH (1980). “The Advantages of Ectothermy for Tetrapods.” The American Naturalist, 115(1), 92–112. ISSN 00030147, 15375323.

Sample BE, Cooper RJ, Greer RD, Whitmore RC (1993). “Estimation of Insect Biomass by Length and Width.” The American Midland Naturalist, 129(2), 234–240. ISSN 00030031, 19384238, doi:10.2307/2426503.

See Also

Other allometric functions: proportion_silhouette_area_shapes(), proportion_silhouette_area(), surface_area_from_length(), surface_area_from_mass(), surface_area_from_volume(), volume_from_length()

Examples

mass_from_length(l     = 0.04,
                   taxon = "insect")
  mass_from_length(l     = 0.04,
                   taxon = "lizard")
  mass_from_length(l     = 0.04,
                   taxon = "salamander")
  mass_from_length(l     = 0.04,
                   taxon = "frog")
  mass_from_length(l     = 0.04, 
                   taxon = "snake")
  mass_from_length(l     = 0.04, 
                   taxon = "turtle")

Average Monthly Solar Radiation

Description

The function estimates average monthly solar radiation (W m-2 d-1) using basic topographic and climatic information as input. Cloudiness is stochastically modeled, so output will vary between functional calls. Based on Nikolov and Zeller (1992).

Usage

monthly_solar_radiation(lat, lon, doy, elev, T_a, hp, P)

Arguments

lat

numeric latitude (degrees).

lon

numeric longitude (degrees).

doy

numeric day of year (1-366).

elev

numeric elevation (meters).

T_a

numeric mean monthly air temperature (C).

hp

numeric mean month relative humidity (percentage).

P

numeric total monthly precipitation (mm).

Value

numeric average monthly solar radiation (W m-2).

References

Nikolov NT, Zeller K (1992). “A solar radiation algorithm for ecosystem dynamic models.” Ecological Modelling, 61(3-4), 149-168. doi:10.24200/tjer.vol8iss2pp19-31.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), air_temp_profile(), degree_days(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

monthly_solar_radiation(lat  = 47.61, 
                          lon  = -122.33, 
                          doy  = 112, 
                          elev = 1500, 
                          T_a  = 15, 
                          hp   = 50, 
                          P   = 50)

Nusselt Number from the Grashof Number

Description

The function estimates the Nusselt number from the Grashof Number (Gates 1980).

Usage

Nusselt_from_Grashof(Gr)

Arguments

Gr

numeric Grashof Number (dimensionless).

Value

numeric Nusselt number (dimensionless).

References

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Nusselt_from_Grashof(Gr = 5)

Nusselt Number from the Reynolds Number

Description

The function estimates the Nusselt number from the Reynolds number for various taxa using Mitchell (1976) (Table 1: Convective Heat Transfer Relations for Animal Shapes).

Usage

Nusselt_from_Reynolds(Re, taxon = "cylinder")

Arguments

Re

numeric Reynolds Number (dimensionless).

taxon

character which class of organism. Current choices: "sphere", "cylinder", "frog", "lizard_traverse_to_air_flow", "lizard_parallel_to_air_flow", "lizard_surface", "lizard_elevated", "flyinginsect", "spider".

Value

numeric Nusselt number (dimensionless).

References

Mitchell JW (1976). “Heat transfer from spheres and other animal forms.” Biophysical Journal, 16(6), 561-569. ISSN 0006-3495, doi:10.1016/S0006-3495(76)85711-6.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Nusselt_from_Reynolds(Re    = 5,
                        taxon = "cylinder")

Nusselt Number

Description

The function estimates the Nusselt Number, which describes dimensionless conductance (Gates 1980).

Usage

Nusselt_number(H_L, D, K)

Arguments

H_L

numeric convective heat transfer coefficient (W m-2 K-1).

D

numeric characteristic dimension (e.g., body diameter) (m).

K

numeric thermal conductivity (W K-1 m-1).

Value

numeric Nusselt number.

References

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Nusselt_number(H_L = 20,
                 D   = 0.01,
                 K   = 0.5)

Diffuse Fraction for Partitioning Solar Radiation

Description

The function partitions solar radiation (W m-2) into direct and diffuse components by estimating the diffuse fraction (k_d). The function uses the models presented in Wong and Chow (2001).

Usage

partition_solar_radiation(method, kt, lat = NA, sol.elev = NA)

Arguments

method

character method to use for estimating the diffuse fraction, currently available: "Liu_Jordan", "Orgill_Hollands", "Erbs", "Olyphant", "Spencer", "Reindl-1", "Reindl-2", "Lam_Li".

kt

numeric the clearness index (dimensionless), which is the ratio of the global solar radiation measured at the surface to the total solar radiation at the top of the atmosphere. (0-1)

lat

numeric latitude (degrees). Needed only if method is "Spencer".

sol.elev

numeric the solar elevation angles (degrees). Needed only if method is "Reindl-2".

Value

numeric diffuse fraction.

References

Wong LT, Chow WK (2001). “Solar radiation model.” Applied Energy, 69(3), 191-224. ISSN 0306-2619, doi:10.1016/S0306-2619(01)00012-5, https://www.sciencedirect.com/science/article/pii/S0306261901000125.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), air_temp_profile(), degree_days(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), monthly_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

partition_solar_radiation(method   = "Erbs", 
                            kt       = 0.5, 
                            lat      = 40, 
                            sol.elev = 60)

Prandtl Number

Description

The function estimates the Prandtl Number, which describes the ratio of kinematic viscosity to thermal diffusivity (Gates 1980).

Usage

Prandtl_number(c_p, mu, K)

Arguments

c_p

numeric specific heat at constant pressure (J mol-1 K-1).

mu

numeric dynamic viscosity (mol s-1 m-1).

K

numeric thermal conductivity (W K-1 m-1).

Value

numeric Prandtl number.

References

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Prandtl_number(c_p = 29.3,
                 mu  = 0.00001,
                 K   = 0.5)

Ratio of Diffuse to Direct Solar Radiation

Description

The function estimates the ratio of diffuse to direct solar radiation based on the approximation of the SOLRAD model (McCullough and Porter 1971) described in Tracy et al. (1983).

Usage

proportion_diffuse_solar_radiation(psi, p_a, rho)

Arguments

psi

numeric Zenith angle of the sun (degrees).

p_a

numeric Atmospheric pressure (kPa).

rho

numeric albedo of the substrate (fraction of 1).

Value

numeric diffuse fraction.

References

McCullough EC, Porter WP (1971). “Computing Clear Day Solar Radiation Spectra for the Terrestrial Ecological Environment.” Ecology, 52(6), 1008-1015. doi:10.2307/1933806.

Tracy CR, Hammond KA, Lechleitner RA, II WJS, Thompson DB, Whicker AD, Williamson SC (1983). “Estimating clear-day solar radiation: an evaluation of three models.” Journal of Thermal Biology, 8(3), 247-251. doi:10.1016/0306-4565(83)90003-7.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), air_temp_profile(), degree_days(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), monthly_solar_radiation(), partition_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

proportion_diffuse_solar_radiation(psi = 60, 
                                     p_a = 86.1, 
                                     rho   = 0.25)

Organism Silhouette Area

Description

The function estimates the projected (silhouette) area as a portion of the surface area of the organism as a function of zenith angle. The function is useful for estimating absorbed solar radiation.

Usage

proportion_silhouette_area(psi, taxon, raz = 0, posture = "prostrate")

Arguments

psi

numeric zenith angle in degrees between 0 and 360.

taxon

character organism name. Current choices are "lizard", "frog", and "grasshopper".

raz

numeric relative solar azimuth angle (in degrees). Required if taxon = "lizard". This is the horizontal angle of the sun relative to the head and frontal plane of the lizard and options currently include 0 (in front), 90 (to side), and 180 (behind) degrees.

posture

character value describing posture. Required if taxon = "lizard". Options include "prostrate" (default) and "elevated".

Details

Relationships come from

  • Lizards: Muth (1977)

  • Frogs: Tracy (1976)

  • Grasshoppers: Anderson et al. (1979)

Value

numeric silhouette area as a proportion.

References

Anderson RV, Tracy CR, Abramsky Z (1979). “Habitat Selection in Two Species of Short-Horned Grasshoppers. The Role of Thermal and Hydric Stresses.” Oecologia, 38(3), 359–374. doi:10.1007/BF00345194.

Muth A (1977). “Thermoregulatory Postures and Orientation to the Sun: A Mechanistic Evaluation for the Zebra-Tailed Lizard, Callisaurus draconoides.” Copeia, 4, 710 - 720.

Tracy CR (1976). “A Model of the Dynamic Exchanges of Water and Energy between a Terrestrial Amphibian and Its Environment.” Ecological Monographs, 46(3), 293-326. doi:10.2307/1942256.

See Also

Other allometric functions: mass_from_length(), proportion_silhouette_area_shapes(), surface_area_from_length(), surface_area_from_mass(), surface_area_from_volume(), volume_from_length()

Examples

proportion_silhouette_area(psi     = 60,   
                             taxon = "frog")
  proportion_silhouette_area(psi     = 60, 
                             taxon = "grasshopper")
  proportion_silhouette_area(psi       = 60, 
                             taxon   = "lizard", 
                             posture = "prostrate", 
                             raz     = 90)
  proportion_silhouette_area(psi       = 60, 
                             taxon   = "lizard", 
                             posture = "elevated", 
                             raz     = 180)

Organism Silhouette Area using Shape Approximations

Description

The function estimates the projected (silhouette) area as a portion of the surface area of the organism. The function estimates the projected area as a function of the dimensions and the angle between the solar beam and the longitudinal axis of the solid, using Figure 11.6 in Campbell and Norman (1998). The function is useful for estimating absorbed solar radiation.

Usage

proportion_silhouette_area_shapes(shape, theta, h, d)

Arguments

shape

character Shape to use to approximate an organism. Shapes are assumed to be prolate or have the longest axis parallel with the ground. Current choices are "spheroid", "cylinder flat ends", and "cylinder hemisphere ends".

theta

numeric angle between the solar beam and the longitudinal axis (degrees).

h

numeric height (long axis in m). Cross section length for spheroid.

d

numeric diameter (short axis in m). Cross section length for spheroid.

Value

numeric silhouette area as a proportion.

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

See Also

Other allometric functions: mass_from_length(), proportion_silhouette_area(), surface_area_from_length(), surface_area_from_mass(), surface_area_from_volume(), volume_from_length()

Examples

proportion_silhouette_area_shapes(shape = "spheroid", 
                                    theta = 60,  
                                    h     = 0.01,  
                                    d     = 0.001)
  proportion_silhouette_area_shapes(shape = "cylinder flat ends",  
                                    theta = 60,  
                                    h     = 0.01, 
                                    d     = 0.001)
  proportion_silhouette_area_shapes(shape = "cylinder hemisphere ends",  
                                    theta = 60,  
                                    h     = 0.01, 
                                    d     = 0.001)

Conductance Assuming Animal Thermal Conductivity is Rate Limiting

Description

The function calculates conductance (W) of an ectothermic animal to its substrate. Method assumes the major resistance to conduction is within surface layers of the animal and that the interior of the animal is equal in temperature to its surface (thermally well mixed) (Spotila et al. 1992).

Usage

Qconduction_animal(T_g, T_b, d, K = 0.5, A, proportion)

Arguments

T_g

numeric ground surface temperature (K).

T_b

numeric body temperature (K).

d

numeric mean thickness of the animal skin (surface) in (meters). The function assumes a well mixed interior.

K

numeric thermal conductivity (W K-1 m-1). K = 0.5 for naked skin and K = 0.15 for insect cuticle (Galushko et al. 2005). The conductivity of the ground is generally greater than that of animal tissues, so animal thermal conductivity is generally the rate limiting step.

A

numeric surface area (m2).

proportion

numeric proportion of body in contact with the surface (0-1).

Value

numeric conductance (W).

References

Galushko D, Ermakov N, Karpovski M, Palevski A, Ishay JS, Bergman DJ (2005). “Electrical, thermoelectric and thermophysical properties of hornet cuticle.” Semiconductor Science and Technology, 20(3), 286–289. doi:10.1088/0268-1242/20/3/005.

Spotila JR, Feder ME, Burggren WW (1992). “Biophysics of Heat and Mass Transfer.” Environmental Physiology of the Amphibians. https://press.uchicago.edu/ucp/books/book/chicago/E/bo3636401.html.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Qconduction_animal(T_g        = 293,
                     T_b        = 303,
                     d          = 10^-6,
                     K          = 0.5,
                     A          = 10^-3,
                     proportion = 0.2)

Conductance Assuming Substrate Thermal Conductivity is Rate Limiting

Description

The function calculates conductance (W) of an ectothermic animal to its substrate. The method assumes the major resistance to conduction is the substrate and that the interior of the animal is equal in temperature to its surface (thermally well mixed) (Spotila et al. 1992).

Usage

Qconduction_substrate(T_g, T_b, D, K_g = 0.5, A, proportion)

Arguments

T_g

numeric surface temperature (K).

T_b

numeric body temperature (K).

D

numeric characteristic dimension of the animal (m).

K_g

numeric thermal conductivity of substrate (W K-1 m-1).

A

numeric surface area (m2).

proportion

numeric proportion in contact to the surface.

Value

numeric conductance (W).

References

Spotila JR, Feder ME, Burggren WW (1992). “Biophysics of Heat and Mass Transfer.” Environmental Physiology of the Amphibians. https://press.uchicago.edu/ucp/books/book/chicago/E/bo3636401.html.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Qconduction_substrate(T_g        = 293,
                        T_b        = 303,
                        D          = 0.01,
                        K_g        = 0.3,
                        A          = 10^-2,
                        proportion = 0.2)

Organismal Convection

Description

The function calculates convection from an organism to its environment as in Mitchell (1976). It includes an enhancement factor associated with outdoor environments.

Usage

Qconvection(T_a, T_b, A, proportion, H_L = 10.45, ef = 1.23)

Arguments

T_a

numeric air temperature (K).

T_b

numeric initial body temperature (K).

A

numeric surface area (m2).

proportion

numeric proportion of surface area exposed to air.

H_L

numeric convective heat transfer coefficient (W K-1 m-2).

ef

numeric is the enhancement factor, used to adjust H to field conditions. Approximated as mean value of 1.23 by default, but see Mitchell (1976) for further information.

Value

numeric convection (W).

References

Mitchell JW (1976). “Heat transfer from spheres and other animal forms.” Biophysical Journal, 16(6), 561-569. ISSN 0006-3495, doi:10.1016/S0006-3495(76)85711-6.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Qconvection(T_a        = 293,
              T_b        = 303,
              H_L        = 10.45,
              A          = 0.0025,
              proportion = 0.85)

Emitted Thermal Radiation

Description

The function estimates the net thermal radiation (W) emitted by the surface of an animal (Gates 1980; Spotila et al. 1992).

Usage

Qemitted_thermal_radiation(
  epsilon = 0.96,
  A,
  psa_dir,
  psa_ref,
  T_b,
  T_g,
  T_sky = NA,
  T_a,
  enclosed = FALSE
)

Arguments

epsilon

numeric longwave infrared emissivity of skin (proportion), 0.95 to 1 for most animals (Gates 1980).

A

numeric surface area ((m2).

psa_dir

numeric view factor indicating the proportion of surface area exposed to sky (or enclosure) (0-1)

psa_ref

numeric view factor indicating the proportion surface area exposed to ground (0-1).

T_b

numeric body surface temperature (K).

T_g

numeric ground surface temperature (K).

T_sky

numeric Estimate effective radiant temperature of sky (K)

T_a

numeric ambient air temperature (K), only required if the animal is in an enclosed environment.

enclosed

logical whether the animal is an enclosed environment or not.

Value

numeric emitted thermal radiation, Qemit (W).

References

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

Spotila JR, Feder ME, Burggren WW (1992). “Biophysics of Heat and Mass Transfer.” Environmental Physiology of the Amphibians. https://press.uchicago.edu/ucp/books/book/chicago/E/bo3636401.html.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Qemitted_thermal_radiation(epsilon  = 0.96,
                             A        = 1,
                             psa_dir  = 0.4,
                             psa_ref  = 0.5,
                             T_b      = 303,
                             T_g      = 293,
                             T_a      = 298,
                             enclosed = FALSE)

Heat Loss Associated with Evaporative Water Loss

Description

The function estimates heat loss associated with evaporative water loss for an amphibian (Spotila et al. 1992) or lizard. The lizard estimation is based on empirical measurements in Porter et al. (1973)).

Usage

Qevaporation(A, T_b, taxon, e_s = NA, e_a = NA, hp = NA, H = NA, r_i = NA)

Arguments

A

numeric surface area (m2).

T_b

numeric body temperature (K).

taxon

character organism type. Current choices: "lizard", "amphibian_wetskin" (fully wet skin), "amphibian" (not fully wet skin).

e_s

numeric saturation water vapor density at skin surface (kg m-3) (needed if amphibian).

e_a

numeric saturation water vapor density in ambient air (kg m-3) (needed if amphibian).

hp

numeric relative humidity (0-1) (needed if amphibian).

H

numeric convective heat transfer coefficient (W m-2 K-1) (needed if amphibian).

r_i

numeric internal (cutaneous) resistance to vapor transport (s m-1) (needed if amphibian).

Value

numeric evaporative heat loss (W).

References

Porter WP, Mitchell JW, Bekman A, DeWitt CB (1973). “Behavioral implications of mechanistic ecology: thermal and behavioral modeling of desert ectotherms and their microenvironments.” Oecologia, 13, 1-54.

Spotila JR, Feder ME, Burggren WW (1992). “Biophysics of Heat and Mass Transfer.” Environmental Physiology of the Amphibians. https://press.uchicago.edu/ucp/books/book/chicago/E/bo3636401.html.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Qevaporation(A      = 0.1,
               T_b    = 293,
               taxon = "amphibian",
               e_s = 0.003,
               e_a = 0.002,
               hp     = 0.5,
               H     = 20,
               r_i   = 50)
  Qevaporation(A     = 0.1,
               T_b   = 293,
               taxon = "lizard")

Metabolism as a Function of Mass

Description

The function estimates the field metabolic rate (W) of various taxa as a function of mass (g). The function does not account for temperature and is based on empirical relationships from Nagy (2005).

Usage

Qmetabolism_from_mass(m, taxon = "reptile")

Arguments

m

numeric mass (grams).

taxon

character taxon for calculation. Options: "reptile", "bird", "mammal".

Value

numeric metabolism (W).

References

Nagy KA (2005). “Field metabolic rate and body size.” Journal of Experimental Biology, 208, 1621-1625. doi:10.1242/jeb.01553, https://journals.biologists.com/jeb/article/208/9/1621/9364/Field-metabolic-rate-and-body-size.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Qmetabolism_from_mass(m     = 12,
                        taxon = "reptile")

Metabolism as a Function of Mass and Body Temperature

Description

The function estimates basal (or resting) metabolic rate (W) as a function of mass (g) and temperature (C). The function is based on empirical data and the metabolic theory of ecology (assumes a 3/4 scaling exponent) (Gillooly et al. 2001).

Usage

Qmetabolism_from_mass_temp(m, T_b, taxon)

Arguments

m

numeric mass (grams).

T_b

numeric body temperature (C).

taxon

character organism type. Options: "bird", "mammal", "reptile", "amphibian", "invertebrate".

Value

numeric basal metabolism (W).

References

Gillooly JF, Brown JH, West GB, Savage VM, Charnov EL (2001). “Effects of size and temperature on metabolic rate.” Science, 293, 2248-2251. doi:10.1126/science.1061967.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Qmetabolism_from_mass_temp(m     = 100,
                             T_b   = 30,
                             taxon = "reptile")

Net Energy Exchange Between an Animal and the Environment

Description

The function estimates the net energy exchange (W) between an animal and the environment. The function follows Gates (1980) and others.

Usage

Qnet_Gates(Qabs, Qemit, Qconv, Qcond, Qmet, Qevap)

Arguments

Qabs

numeric solar radiation absorbed (W).

Qemit

numeric thermal radiation emitted (W).

Qconv

numeric energy exchange due to convection; Energy exchange from an animal to its surrounding environment (air or water) (W).

Qcond

numeric energy exchange due to conduction; Energy exchange from animal to a surface if they are in contact (W).

Qmet

numeric energy emitted due to metabolism (W).

Qevap

numeric energy emitted due to evaporative water loss (W).

Value

numeric net energy exchange (W).

References

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Qnet_Gates(Qabs  = 500,
           Qemit = 10, 
           Qconv = 100, 
           Qcond = 100, 
           Qmet  = 10, 
           Qevap = 5)

Absorbed Solar and Thermal Radiation

Description

The function estimates solar and thermal radiation (W) absorbed by the surface of an animal following (Gates 1980) and (Spotila et al. 1992).

Usage

Qradiation_absorbed(
  a = 0.9,
  A,
  psa_dir,
  psa_dif,
  psa_ref,
  S_dir,
  S_dif,
  S_ref = NA,
  rho = NA
)

Arguments

a

numeric solar absorptivity of animal surface (proportion), default value is for reptiles (0-1).

A

numeric surface area (m2).

psa_dir

numeric proportion surface area exposed to direct solar radiation (0-1).

psa_dif

numeric proportion surface area exposed to diffuse solar radiation (0-1).

psa_ref

numeric proportion surface area exposed to reflected solar radiation (0-1).

S_dir

numeric direct solar radiation (W m-2).

S_dif

numeric diffuse solar radiation (W m-2).

S_ref

numeric reflected solar radiation (W m-2), either provided or estimated if surface albedo is provided instead

rho

numeric is surface albedo (proportion), optional (not used) if reflected radiation is provided. Values available in (Gates 1980) Table 8.2.

Value

numeric solar radiation absorbed (W)

References

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

Spotila JR, Feder ME, Burggren WW (1992). “Biophysics of Heat and Mass Transfer.” Environmental Physiology of the Amphibians. https://press.uchicago.edu/ucp/books/book/chicago/E/bo3636401.html.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Qradiation_absorbed(a       = 0.9,
                      A       = 1,
                      psa_dir = 0.5,
                      psa_dif = 0.5,
                      psa_ref = 0.5,
                      S_dir   = 1000,
                      S_dif   = 200,
                      rho     = 0.5)

Absorbed Thermal Radiation

Description

The function estimates longwave (thermal) radiation (W) absorbed from the sky and the ground (Campbell and Norman 1998; Riddell et al. 2018).

Usage

Qthermal_radiation_absorbed(
  T_a,
  T_g,
  epsilon_ground = 0.97,
  a_longwave = 0.965
)

Arguments

T_a

numeric air temperature (C).

T_g

numeric ground temperature (C).

epsilon_ground

numeric emissivity (proportion) for more soil types (Campbell and Norman 1998).

a_longwave

numeric absorptance (proportion) of organism to longwave radiation (Bartlett and Gates 1967; Buckley 2008).

Value

numeric thermal radiation absorbed (W).

Author(s)

Eric Riddell

References

Bartlett PN, Gates DM (1967). “The energy budget of a lizard on a tree trunk.” Ecology, 48, 316-322.

Buckley LB (2008). “Linking traits to energetics and population dynamics to predict lizard ranges in changing environments.” American Naturalist, 171(1), E1 - E19. doi:10.1086/523949, https://pubmed.ncbi.nlm.nih.gov/18171140/.

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

Riddell EA, Odom JP, Damm JD, Sears MW (2018). “Plasticity reveals hidden resistance to extinction under climate change in the global hotspot of salamander diversity.” Science Advances, 4(4). doi:10.1126/sciadv.aar5471.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Qthermal_radiation_absorbed(T_a            = 20,
                              T_g            = 25,
                              epsilon_ground = 0.97,
                              a_longwave     = 0.965)

Reynolds Number

Description

The function estimates the Reynolds Number, which describes the dynamic properties of the fluid surrounding the animal as the ratio of internal viscous forces (Gates 1980).

Usage

Reynolds_number(u, D, nu)

Arguments

u

numeric wind speed (m s-1).

D

numeric characteristic dimension (e.g., body diameter) (m)

nu

numeric the kinematic viscosity, ratio of dynamic viscosity to density of the fluid (m2 s-1); can calculate from DRYAIR or WETAIR.

Value

numeric Reynolds number.

References

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Reynolds_number(u  = 1,
                  D  = 0.001,
                  nu = 1.2)

Saturation Vapor Pressure

Description

The function calculates saturation vapor pressure (kPa) based on the Clausius-Clapeyron equation (Stull 2000; Riddell et al. 2018).

Usage

saturation_vapor_pressure(T_a)

Arguments

T_a

numeric air temperature (K).

Value

numeric saturation vapor pressure, e_s (kPa).

Author(s)

Eric Riddell

References

Riddell EA, Odom JP, Damm JD, Sears MW (2018). “Plasticity reveals hidden resistance to extinction under climate change in the global hotspot of salamander diversity.” Science Advances, 4(4). doi:10.1126/sciadv.aar5471.

Stull RB (2000). Meteorology for Scientists and Engineers. Brooks Cole. ISBN 978-0534372149.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_water_vapor_pressure()

Examples

saturation_vapor_pressure(T_a = 293)

Saturation Water Vapor Pressure

Description

The function approximates saturation water vapor pressure as a function of ambient temperature for temperatures from 0 to 40 C using Rosenberg (1974) in Spotila et al. (1992). See also NicheMapR WETAIR and DRYAIR (Kearney and Porter 2020).

Usage

saturation_water_vapor_pressure(T_a)

Arguments

T_a

numeric air temperature (C).

Value

numeric Saturation water vapor pressure, e_s (Pa).

References

Kearney MR, Porter WP (2020). “NicheMapR - an R package for biophysical modelling: the ectotherm and Dynamic Energy Budget models.” Ecography, 43(1), 85-96. doi:10.1111/ecog.04680.

Rosenberg NJ (1974). Microclimate: the biological environment. Wiley, New York.

Spotila JR, Feder ME, Burggren WW (1992). “Biophysics of Heat and Mass Transfer.” Environmental Physiology of the Amphibians. https://press.uchicago.edu/ucp/books/book/chicago/E/bo3636401.html.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure()

Examples

saturation_water_vapor_pressure(T_a = 20)

Soil Thermal Conductivity

Description

The function estimates soil thermal conductivity (W m-1 K-1) using the methods of deVries (1963).

Usage

soil_conductivity(x, lambda, g_a)

Arguments

x

numeric vector of volume fractions of soil constituents (e.g., clay, quartz, minerals other than quartz, organic matter, water, air). The volume fractions should sum to 1. Note that x and lambda values in the example correspond to these soil constituents.

lambda

numeric vector of the thermal conductivities (W m-1 K-1) of the soil constituents.

g_a

numeric shape factor on soil particles. The soil particles are assumed to be ellipsoids with axes g_a, g_b, and g_c, where g_a + g_b + g_c = 1 and g_a = g_b. deVries (1952) suggests g_a = g_b = 0.125.

Value

numeric soil thermal conductivity (W m-1 K-1).

Author(s)

Joseph Grigg

References

deVries DA (1952). “Thermal Conductivity of Soil.” Nature, 178, 1074. doi:10.1038/1781074a0.

deVries DA (1963). “Thermal Properties of Soils.” In Physics of Plant Environment. North Holland Publishing Company. doi:10.1002/qj.49709038628.

See Also

Other soil temperature functions: soil_specific_heat(), soil_temperature_equation(), soil_temperature_function(), soil_temperature_integrand(), soil_temperature()

Examples

soil_conductivity(x      = c(0.10, 0.40, 0.11, 0.01, 0.2, 0.18), 
                    lambda = c(0.10, 0.40, 0.11, 0.01, 0.2, 0.18), 
                    g_a     = 0.125)

Soil Specific Heat

Description

The function estimates soil specific heat (J kg-1 K-1) using the methods of deVries (1963). The function incorporates the volume fraction of organic material, minerals, and water in soil.

Usage

soil_specific_heat(x_o, x_m, x_w, rho_so)

Arguments

x_o

numeric volume fraction of organic material (0-1).

x_m

numeric volume fraction of minerals (0-1).

x_w

numeric volume fraction of water (0-1).

rho_so

numeric particle density of soil in (kg m-3) (bulk density).

Value

numeric soil specific heat (J kg-1 K-1).

Author(s)

Joseph Grigg

References

deVries DA (1963). “Thermal Properties of Soils.” In Physics of Plant Environment. North Holland Publishing Company. doi:10.1002/qj.49709038628.

See Also

Other soil temperature functions: soil_conductivity(), soil_temperature_equation(), soil_temperature_function(), soil_temperature_integrand(), soil_temperature()

Examples

soil_specific_heat(x_o    = 0.01, 
                     x_m    = 0.6, 
                     x_w    = 0.2, 
                     rho_so = 1620)

Calculate Soil Temperature using ODEs

Description

This function is called to calculate soil temperature (C) as in Beckman et al. (1973). This function calls soil_temperature_function, which uses ODEs to calculate a soil profile using equations from deVries (1963)

Usage

soil_temperature(
  z_r.intervals = 12,
  z_r,
  z,
  T_a,
  u,
  Tsoil0,
  z0,
  SSA,
  TimeIn,
  S,
  water_content = 0.2,
  air_pressure,
  rho_so = 1620,
  shade = FALSE
)

Arguments

z_r.intervals

numeric the number of intervals in the soil profile to calculate, defaults to 12.

z_r

numeric reference height (m).

z

numeric interval of the soil profile to return (1 to z_r.intervals).

T_a

numeric vector of air temperature (degrees C), Note: missing values will be linearly interpolated.

u

numeric vector of wind speeds (m s-1).

Tsoil0

numeric initial soil temperature (degrees C).

z0

numeric surface roughness (m).

SSA

numeric solar absorptivity of soil surface as a fraction.

TimeIn

numeric vector of time periods for the model.

S

numeric vector of solar radiation (W m-2).

water_content

numeric percent water content (percent).

air_pressure

numeric air pressure (kPa).

rho_so

numeric particle density of soil.

shade

logical whether or not soil temperature should be calculated in the shade.

Value

numeric soil temperature (C).

Author(s)

Joseph Grigg

References

Beckman WA, Mitchell JW, Porter WP (1973). “Thermal Model for Prediction of a Desert Iguana's Daily and Seasonal Behavior.” Journal of Heat Transfer, 95(2), 257-262. doi:10.1115/1.3450037.

deVries DA (1963). “Thermal Properties of Soils.” In Physics of Plant Environment. North Holland Publishing Company. doi:10.1002/qj.49709038628.

See Also

Other soil temperature functions: soil_conductivity(), soil_specific_heat(), soil_temperature_equation(), soil_temperature_function(), soil_temperature_integrand()

Examples

set.seed(123)
  temp_vector       <- runif(48, min = -10, max = 10)
  wind_speed_vector <- runif(48, min = 0, max = 0.4)
  time_vector       <- rep(1:24, 2)
  solrad_vector     <- rep(c(rep(0, 6), 
                             seq(10, 700, length.out = 6), 
                             seq(700, 10, length.out = 6), 
                             rep(0, 6)),
                           2)

  soil_temperature(z_r.intervals = 12, 
                   z_r           = 1.5, 
                   z             = 2, 
                   T_a           = temp_vector, 
                   u             = wind_speed_vector, 
                   Tsoil0        = 20, 
                   z0            = 0.02, 
                   SSA           = 0.7, 
                   TimeIn        = time_vector, 
                   S             = solrad_vector, 
                   water_content = 0.2, 
                   air_pressure  = 85, 
                   rho_so        = 1620, 
                   shade         = FALSE)

Core Function Called to Solve Equation for Soil Temperature

Description

The function called by soil_temperature_function to solve equation for soil temperature from Beckman et al. (1973).

Usage

soil_temperature_equation(L, rho_a, c_a, u_inst, z_r, z0, T_inst, T_s)

Arguments

L

numeric Monin-Obukhov length, a measure of the instability of heat flow (see Beckman et al. (1973)).

rho_a

numeric density of air (kg m-3).

c_a

numeric specific heat of air (J kg-1 C-1).

u_inst

numeric instantaneous wind speed (m s-1).

z_r

numeric reference height (m).

z0

numeric surface roughness (m).

T_inst

numeric instantaneous air temperature (C).

T_s

numeric initial soil surface temperature (C).

Value

numeric soil temperature (C).

Author(s)

Joseph Grigg

References

Beckman WA, Mitchell JW, Porter WP (1973). “Thermal Model for Prediction of a Desert Iguana's Daily and Seasonal Behavior.” Journal of Heat Transfer, 95(2), 257-262. doi:10.1115/1.3450037.

See Also

Other soil temperature functions: soil_conductivity(), soil_specific_heat(), soil_temperature_function(), soil_temperature_integrand(), soil_temperature()

Examples

soil_temperature_equation(L      = -10, 
                            rho_a  = 1.177, 
                            c_a    = 1006, 
                            u_inst = 0.3, 
                            z_r    = 1.5, 
                            z0     = 0.02, 
                            T_inst = 8, 
                            T_s    = 20)

Core Function for Calculating Soil Temperature

Description

This function is called to calculate soil temperature as in Beckman et al. (1973). Parameters are passed as a list to facilitating solving the equations. This function is not intended to be called directly. The energy balance equations are from Porter et al. (1973) and Kingsolver (1979)

Usage

soil_temperature_function(j, T_so, params)

Arguments

j

numeric the number of the iteration of running the model.

T_so

numeric the initial soil temperature profile in C.

params

list containing the following param, which are described or calculated in soil_temperature: SSA, epsilon_so, k_so, c_so, dz, z_r, z0, S, T_a, u, rho_a, rho_so, c_a, TimeIn, dt, shade.

Value

Soil temperature profile as a list.

Author(s)

Joseph Grigg

References

Beckman WA, Mitchell JW, Porter WP (1973). “Thermal Model for Prediction of a Desert Iguana's Daily and Seasonal Behavior.” Journal of Heat Transfer, 95(2), 257-262. doi:10.1115/1.3450037.

Kingsolver JG (1979). “Thermal and hydric aspects of environmental hetergeneity in the pitcher plant mosquito.” Ecological Monographs, 49, 357-376.

Porter WP, Mitchell JW, Bekman A, DeWitt CB (1973). “Behavioral implications of mechanistic ecology: thermal and behavioral modeling of desert ectotherms and their microenvironments.” Oecologia, 13, 1-54.

See Also

Other soil temperature functions: soil_conductivity(), soil_specific_heat(), soil_temperature_equation(), soil_temperature_integrand(), soil_temperature()

Examples

set.seed(123)
  temp_vector       <- runif(96, min = -10, max = 10)
  wind_speed_vector <- runif(96, min = 0, max = 0.4)
  time_vector       <- rep(1:24, 4)
  solrad_vector     <- rep(c(rep(0, 6), 
                             seq(10, 700, length.out = 6), 
                             seq(700, 10, length.out = 6), 
                             rep(0, 6)),
                           4)
  params            <- list(SSA        = 0.7, 
                            epsilon_so = 0.98, 
                            k_so       = 0.293, 
                            c_so       = 800, 
                            dz         = 0.05, 
                            z_r        = 1.5, 
                            z0         = 0.02, 
                            S          = solrad_vector, 
                            T_a        = temp_vector, 
                            u          = wind_speed_vector, 
                            rho_a      = 1.177, 
                            rho_so     = 1620,
                            c_a        = 1006, 
                            TimeIn     = time_vector, 
                            dt         = 60 * 60, 
                            shade      = FALSE)

soil_temperature_function(j      = 1, 
                          T_so   = rep(20,13), 
                          params = params)

Solve Equation for Soil Temperature

Description

This function is called by soil_temperature_equation to solve the equation for soil temperature from Beckman et al. (1973). The function represents the integrand in the equation. It is not intended to be called directly.

Usage

soil_temperature_integrand(x, L, z0)

Arguments

x

numeric vector of volume fractions of soil constituents (e.g., clay, quartz, minerals other than quartz, organic matter, water, air). The volume fractions should sum to 1. Note that x and lambda values in the example correspond to these soil constituents.

L

numeric Monin-Obukhov length, a measure of the instability of heat flow (Beckman et al. 1973).

z0

numeric surface roughness (m).

Value

numeric integrand for soil temperature function.

Author(s)

Joseph Grigg

References

Beckman WA, Mitchell JW, Porter WP (1973). “Thermal Model for Prediction of a Desert Iguana's Daily and Seasonal Behavior.” Journal of Heat Transfer, 95(2), 257-262. doi:10.1115/1.3450037.

See Also

Other soil temperature functions: soil_conductivity(), soil_specific_heat(), soil_temperature_equation(), soil_temperature_function(), soil_temperature()

Examples

soil_temperature_integrand(x  = c(0.10, 0.40, 0.11, 0.01, 0.2, 0.18), 
                             L  = -10, 
                             z0 = 0.2)

Time of Solar Noon

Description

The function calculates the time of solar noon in hours as a function of the day of year and longitude (Campbell and Norman 1998).

Usage

solar_noon(lon, doy, offset = NA)

Arguments

lon

numeric longitude (decimal degrees).

doy

numeric day of year (1-366). This can be obtained from a standard date via day_of_year.

offset

numeric number of hours to add to UTC (Coordinated Universal Time) to get local time (improves accuracy but not always necessary). Defaults to NA.

Value

numeric time of solar noon (hours).

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

See Also

Other utility functions: airpressure_from_elev(), azimuth_angle(), day_of_year(), daylength(), dec_angle(), temperature conversions, zenith_angle()

Examples

solar_noon(lon = -122.335,
             doy = 112)

Estimate the Three Components of Solar Radiation (Direct, Diffuse and Reflected)

Description

The function estimate direct, diffuse, and reflected components of solar radiation (W m-2) as a function of day of year using the model in Campbell and Norman (1998).

Usage

solar_radiation(doy, psi, tau, elev, rho = 0.7)

Arguments

doy

numeric the day of year; day_of_year.

psi

numeric zenith angle (radians).

tau

numeric atmospheric transmissivity (proportion), which is ratio of global solar radiation at ground level to extra-terrestrial solar radiation.

elev

numeric elevation (meters).

rho

numeric albedo as a proportion (0-1).

Value

numeric radiation components - direct, diffused and reflected (W m-2).

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), air_temp_profile(), degree_days(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), monthly_solar_radiation(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), surface_roughness(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

solar_radiation(doy  = 112, 
                  psi  = 1, 
                  tau  = 0.6, 
                  elev = 1500, 
                  rho  = 0.7)

Organism Surface Area from Length

Description

This function estimates surface area (m2) from length (m) by approximating the animal's body as a rotational ellipsoid with half the body length as the semi-major axis.

Usage

surface_area_from_length(l)

Arguments

l

numeric length (m).

Details

Following Samietz et al. (2005) and Lactin and Johnson (1998).

Value

numeric surface area (m2).

References

Lactin DJ, Johnson DL (1998). “Convective heat loss and change in body temperature of grasshopper and locust nymphs: Relative importance of wind speed, insect size and insect orientation.” Journal of Thermal Biology, 23(1), 5-13. ISSN 0306-4565, doi:10.1016/S0306-4565(97)00037-5, https://www.sciencedirect.com/science/article/pii/S0306456597000375.

Samietz J, Salser MA, Dingle H (2005). “Altitudinal variation in behavioural thermoregulation: local adaptation vs. plasticity in California grasshoppers.” Journal of Evolutionary Biology, 18(4), 1087-1096. doi:10.1111/j.1420-9101.2005.00893.x.

See Also

Other allometric functions: mass_from_length(), proportion_silhouette_area_shapes(), proportion_silhouette_area(), surface_area_from_mass(), surface_area_from_volume(), volume_from_length()

Examples

surface_area_from_length(l = 0.04)

Organism Surface Area from Mass

Description

The function estimates surface area (m2) from mass (g) for one of a variety of taxa.

Usage

surface_area_from_mass(m, taxon)

Arguments

m

numeric vector of mass (g).

taxon

character taxonomic classification of organism, current choices: "lizard", "salamander", "frog", "insect".

Details

All models follow (SA = a Mb) with mass in grams and surface area in meters2.

  • Lizards (Norris 1965; Porter and James 1979; Roughgarden 1981; O'Connor 1999; Fei et al. 2012):

    a=0.000314πa = 0.000314 \pi
    b=2/3b = 2/3

  • Salamanders (Whitford and Hutchison 1967; Riddell et al. 2017):

    a=0.000842a = 0.000842
    b=0.694b = 0.694

  • Frogs (McClanahan and Baldwin 1969):

    a=0.00099a = 0.00099
    b=0.56b = 0.56

  • Insects (Lactin and Johnson 1997):

    a=0.0013a = 0.0013
    b=0.8b = 0.8

Value

numeric surface area (m2).

References

Fei T, Skidmore AK, Venus V, Wang T, Schlerf M, Toxopeus B, van Overjijk S, Bian M, Liu Y (2012). “A body temperature model for lizards as estimated from the thermal environment.” Journal of Thermal Biology, 37(1), 56-64. ISSN 0306-4565, doi:10.1016/j.jtherbio.2011.10.013, https://www.sciencedirect.com/science/article/pii/S0306456511001513.

Lactin DJ, Johnson DL (1997). “Response of body temperature to solar radiation in restrained nymphal migratory grasshoppers (Orthoptera: Acrididae): influences of orientation and body size.” Physiological Entomology, 22(2), 131-139. doi:10.1111/j.1365-3032.1997.tb01150.x.

McClanahan L, Baldwin R (1969). “Rate of water uptake through the integument of the desert toad, Bufo punctatus.” Comparative Biochemistry and Physiology, 28(1), 381-389. ISSN 0010-406X, doi:10.1016/0010-406X(69)91351-6.

Norris KS (1965). “Color adaptation in desert reptiles and its thermal relationships.” In Symposium on Lizard Ecology, 162- 229. U. Missouri Press.

O'Connor M (1999). “Physiological and ecological implications of a simple model of heating and cooling in reptiles.” Journal of Thermal Biology, 24, 113-136.

Porter WP, James FC (1979). “Behavioral Implications of Mechanistic Ecology II: The African Rainbow Lizard, Agama agama.” Copeia, 1979(4), 594–619. ISSN 00458511, 19385110, doi:10.2307/1443867.

Riddell EA, Apanovitch EK, Odom JP, Sears MW (2017). “Physical calculations of resistance to water loss improve predictions of species range models.” Ecological Monographs, 87(1), 21-33. doi:10.1002/ecm.1240.

Roughgarden J (1981). “Resource partitioning of space and its relationship to body temperature in Anolis lizard populations.” Oecologia, 50, 256 – 264. https://link.springer.com/article/10.1007/BF00348048.

Whitford WG, Hutchison VH (1967). “Body Size and Metabolic Rate in Salamanders.” Physiological Zoology, 40(2), 127-133. doi:10.1086/physzool.40.2.30152447.

See Also

Other allometric functions: mass_from_length(), proportion_silhouette_area_shapes(), proportion_silhouette_area(), surface_area_from_length(), surface_area_from_volume(), volume_from_length()

Examples

surface_area_from_mass(m     = 1:50, 
                         taxon = "lizard")
  surface_area_from_mass(m     = 1:50,  
                         taxon = "salamander")
  surface_area_from_mass(m     = 1:50,  
                         taxon = "frog")
  surface_area_from_mass(m     = seq(0.1, 5, 0.1),  
                         taxon = "insect")

Organism Surface Area from Volume

Description

The function estimates surface area (m2) from volume (m3) for a variety of taxa following Mitchell (1976).

Usage

surface_area_from_volume(V, taxon)

Arguments

V

numeric vector of volume (m3). Can be one or more values.

taxon

character taxon of organism, current choices: "lizard", "frog", "sphere".

Details

All models follow (SA = Ka V2/3) with surface area and volume in meters.

  • Lizards: Norris (1965):

    Ka=11.0Ka = 11.0

  • Frogs: Tracy (1972):

    Ka=11.0Ka = 11.0

  • Sphere: Mitchell (1976):

    Ka=4.83Ka = 4.83

Value

numeric surface area (m2).

References

Mitchell JW (1976). “Heat transfer from spheres and other animal forms.” Biophysical Journal, 16(6), 561-569. ISSN 0006-3495, doi:10.1016/S0006-3495(76)85711-6.

Norris KS (1965). “Color adaptation in desert reptiles and its thermal relationships.” In Symposium on Lizard Ecology, 162- 229. U. Missouri Press.

Tracy CR (1972). “Newton's Law: Its Application for Expressing Heat Losses from Homeotherms.” BioScience, 22(11), 656-659. ISSN 0006-3568, doi:10.2307/1296267.

See Also

Other allometric functions: mass_from_length(), proportion_silhouette_area_shapes(), proportion_silhouette_area(), surface_area_from_length(), surface_area_from_mass(), volume_from_length()

Examples

surface_area_from_volume(V     = 0.001, 
                           taxon = "lizard")
  surface_area_from_volume(V     = 0.001,  
                           taxon = "frog")
  surface_area_from_volume(V     = 0.001,  
                           taxon = "sphere")

Surface Roughness from Empirical Measurements

Description

The function estimates surface roughness (m) from empirical wind speed (m s-1) data collected at a vector of heights (m) (Kingsolver and Buckley 2015; Campbell and Norman 1998; Porter and James 1979).

Usage

surface_roughness(u_r, zr)

Arguments

u_r

numeric wind velocity (m s-1) at a vector of reference heights.

zr

numeric vector of reference heights (m).

Value

numeric surface roughness (m).

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

Kingsolver JG, Buckley LB (2015). “Climate variability slows evolutionary responses of Colias butterflies to recent climate change.” Proceedings of the Royal Society B, 282(1802). doi:10.1098/rspb.2014.2470.

Porter WP, James FC (1979). “Behavioral Implications of Mechanistic Ecology II: The African Rainbow Lizard, Agama agama.” Copeia, 1979(4), 594–619. ISSN 00458511, 19385110, doi:10.2307/1443867.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), air_temp_profile(), degree_days(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), monthly_solar_radiation(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), wind_speed_profile_neutral(), wind_speed_profile_segment()

Examples

surface_roughness(u_r = c(0.01, 0.025, 0.05, 0.1, 0.2), 
                    zr  = c(0.05, 0.25, 0.5, 0.75, 1))

Effective radiant temperature of sky (K)

Description

This function estimates the effective radiant temperature of the sky (K) using either the Brunt or Swinbank formulations (Gates 1980).

Usage

T_sky(T_a, formula)

Arguments

T_a

numeric ambient air temperature (K)

formula

character Formula to use, current choices are: "Brunt", "Swinbank".

Value

numeric T_sky (K).

References

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

T_sky(T_a=293, formula="Swinbank")

Operative Environmental Temperature of a Butterfly

Description

The function estimates body temperatures (C, operative environmental temperatures) of a butterfly based on Kingsolver (1983) and Buckley and Kingsolver (2012). The function is designed for butterflies that bask with closed wings such as Colias.

Usage

Tb_butterfly(
  T_a,
  T_g,
  T_sh,
  u,
  S_sdir,
  S_sdif,
  z,
  D,
  delta,
  alpha,
  r_g = 0.3,
  shade = FALSE
)

Arguments

T_a

numeric air temperature (C).

T_g

numeric surface temperature (C) in the sunlight.

T_sh

numeric surface temperature (C) in the shade.

u

numeric wind speed (m s-1).

S_sdir

numeric direct solar radiation flux (W m-2).

S_sdif

numeric diffuse solar radiation flux (W m-2).

z

numeric solar zenith angle (degrees).

D

numeric thoracic diameter (cm).

delta

numeric thoracic fur thickness (mm).

alpha

numeric wing solar absorptivity (proportion). The range for Colias butterflies is 0.4 to 0.7.

r_g

numeric substrate solar reflectivity (proportion). See Kingsolver (1983).

shade

logical whether body temperature should be calculated in sun (FALSE) or shade (TRUE).

Details

Thermal radiative flux is calculated following Gates (1980) based on Swinbank (1960). Kingsolver (1983) estimates using the Brunt equation with black body sky temperature from Swinbank (1963).

Value

numeric predicted body (operative environmental) temperature (C).

References

Buckley LB, Kingsolver JG (2012). “The demographic impacts of shifts in climate means and extremes on alpine butterflies.” Functional Ecology, 26(4), 969-977. doi:10.1111/j.1365-2435.2012.01969.x.

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

Kingsolver JG (1983). “Thermoregulation and Flight in Colias Butterflies: Elevational Patterns and Mechanistic Limitations.” Ecology, 64(3), 534-545. doi:10.2307/1939973.

Swinbank WC (1960). “Wind profile in thermally stratified flow.” Nature, 186, 463-464.

Swinbank WC (1963). “Long-wave radiation from clear skies.” Quarterly Journal of the Royal Meteorological Society, 89, 339-348.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tb_butterfly(T_a    = 25, 
               T_g    = 25, 
               T_sh   = 20, 
               u      = 0.4, 
               S_sdir = 300, 
               S_sdif = 100, 
               z      = 30, 
               D      = 0.36, 
               delta  = 1.46, 
               alpha  = 0.6, 
               r_g    = 0.3)

Operative Environmental Temperature of an Ectotherm based on Campbell and Norman (1988)

Description

The function estimates body temperatures (C, operative environmental temperature) of an ectotherm using an approximation based on Campbell and Norman (1998) and Mitchell (1976).

Usage

Tb_CampbellNorman(
  T_a,
  T_g,
  S,
  a_s = 0.7,
  a_l = 0.96,
  epsilon = 0.96,
  c_p = 29.3,
  D,
  u
)

Arguments

T_a

numeric air temperature (C).

T_g

numeric ground temperature (C).

S

numeric flux density of solar radiation (W m-2), combining direct, diffuse, and reflected radiation accounting for view factors.

a_s

numeric organismal solar absorptivity (proportion).

a_l

numeric organismal thermal absorptivity (proportion); 0.965 for lizards (Bartlett and Gates 1967).

epsilon

numeric longwave infrared emissivity of skin (proportion), 0.95 to 1 for most animals (Gates 1980).

c_p

numeric specific heat of air (J mol-1 C-1).

D

numeric characteristic dimension of the animal (m).

u

numeric wind speed (m s-1).

Details

Boundary conductance uses a factor of 1.4 to account for increased convection (Mitchell 1976). The function assumes forced conduction.

Value

numeric operative environmental temperature, T_e (C).

References

Bartlett PN, Gates DM (1967). “The energy budget of a lizard on a tree trunk.” Ecology, 48, 316-322.

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

Mitchell JW (1976). “Heat transfer from spheres and other animal forms.” Biophysical Journal, 16(6), 561-569. ISSN 0006-3495, doi:10.1016/S0006-3495(76)85711-6.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tb_CampbellNorman (T_a     = 30, 
                   T_g     = 30, 
                   S       = 823, 
                   a_s = 0.7, 
                   a_l = 0.96, 
                   epsilon = 0.96, 
                   c_p     = 29.3, 
                   D       = 0.17, 
                   u       = 1)

Operative Environmental Temperature of an Ectotherm Based on Gates (1980)

Description

The function predicts body temperatures (C, operative environmental temperature) of an ectotherm using the approximation in Gates (1980). The functions omits evaporative and metabolic heat loss (Mitchell 1976; Kingsolver 1983).

Usage

Tb_Gates(
  A,
  D,
  psa_dir,
  psa_ref,
  psa_air,
  psa_g,
  T_g,
  T_a,
  Qabs,
  epsilon,
  H_L,
  ef = 1.3,
  K
)

Arguments

A

numeric surface area (m2).

D

numeric characteristic dimension for conduction (meters).

psa_dir

numeric view factor for the proportion surface area exposed to direct radiation from the sky (or enclosure) (0-1).

psa_ref

numeric view factor for proportion surface area exposed to reflected radiation from the ground (0-1).

psa_air

numeric proportion surface area exposed to air (0-1).

psa_g

numeric proportion surface area in contact with substrate (0-1).

T_g

numeric ground surface temperature (C).

T_a

numeric ambient air temperature (C).

Qabs

numeric Solar radiation absorbed (W).

epsilon

numeric longwave infrared emissivity of skin (proportion), 0.95 to 1 for most animals (Gates 1980).

H_L

numeric Convective heat transfer coefficient (W m-2 K-1).

ef

numeric enhancement factor used to adjust H_L to field conditions using h_L approximation from Mitchell (1976). Approximated as 1.23 by default, but see Mitchell (1976) for relationship.

K

numeric Thermal conductivity (W K-1 m-1), K = 0.5 W K-1 m-1 for naked skin, K = 0.15 W K-1 m-1for insect cuticle Galushko et al. (2005); conductivity of the ground is generally greater than that of animal tissues, so animal thermal conductivity is generally the rate limiting step.

Value

numeric operative environmental temperature, T_e (C).

References

Galushko D, Ermakov N, Karpovski M, Palevski A, Ishay JS, Bergman DJ (2005). “Electrical, thermoelectric and thermophysical properties of hornet cuticle.” Semiconductor Science and Technology, 20(3), 286–289. doi:10.1088/0268-1242/20/3/005.

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

Kingsolver JG (1983). “Thermoregulation and Flight in Colias Butterflies: Elevational Patterns and Mechanistic Limitations.” Ecology, 64(3), 534-545. doi:10.2307/1939973.

Mitchell JW (1976). “Heat transfer from spheres and other animal forms.” Biophysical Journal, 16(6), 561-569. ISSN 0006-3495, doi:10.1016/S0006-3495(76)85711-6.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tb_Gates (A       = 0.1, 
           D       = 0.025, 
           psa_dir = 0.6, 
           psa_ref = 0.4, 
           psa_air = 0.5, 
           psa_g   = 0.1, 
           T_g     = 30, 
           T_a     = 37, 
           Qabs    = 2, 
           epsilon = 0.95, 
           H_L     = 10, 
           ef      = 1.23, 
           K       = 0.5)

Operative Environmental Temperature of an Ectotherm Based on a Variant of Gates (1980)

Description

The function predicts body temperatures (C, operative environmental temperature) of an ectotherm using the approximation in Gates (1980). The function omits evaporative and metabolic heat loss.

Usage

Tb_Gates2(A, D, T_g, T_a, Qabs, u, epsilon)

Arguments

A

numeric surface area (m2).

D

numeric characteristic dimension for conduction (meters).

T_g

numeric ground surface temperature (C).

T_a

numeric ambient air temperature (C).

Qabs

numeric Solar radiation absorbed (W).

u

numeric Wind speed (m s-1).

epsilon

numeric longwave infrared emissivity of skin (proportion), 0.95 to 1 for most animals (Gates 1980).

Value

numeric operative environmental temperature (C).

References

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tb_Gates2(A       = 1, 
            D       = 0.001, 
            T_g     = 27, 
            T_a     = 37, 
            Qabs    = 2, 
            u       = 0.1, 
            epsilon = 1)

Operative Environmental Temperature of a Grasshopper

Description

The function estimates body temperatures (C, operative environmental temperatures) of a grasshopper based on Lactin and Johnson (1998). Part of the model is based on Swinbank (1963), following Gates (1962) in Kingsolver (1983).

Usage

Tb_grasshopper(
  T_a,
  T_g,
  u,
  S,
  K_t,
  psi,
  l,
  Acondfact = 0.25,
  z = 0.001,
  abs = 0.7,
  r_g = 0.3
)

Arguments

T_a

numeric air temperature (C).

T_g

numeric surface temperature (C). Kingsolver (1983) assumes T_g - T_a = 8.4.

u

numeric wind speed (m s-1).

S

numeric total (direct + diffuse) solar radiation flux (W m-2).

K_t

numeric clearness index (dimensionless), which is the ratio of the global solar radiation measured at the surface to the total solar radiation at the top of the atmosphere.

psi

numeric solar zenith angle (degrees).

l

numeric grasshopper length (m).

Acondfact

numeric the proportion of the grasshopper surface area that is in contact with the ground.

z

numeric distance from the ground to the grasshopper (m).

abs

numeric absorptivity of the grasshopper to solar radiation (proportion). See Anderson et al. (1979).

r_g

numeric substrate solar reflectivity (proportion). See Kingsolver (1983).

Details

Total radiative flux is calculated as thermal radiative heat flux plus convective heat flux, following Kingsolver (1983), with the Erbs et al. (1982) model from Wong and Chow (2001).

Energy balance is based on Kingsolver (1983).

Radiation is calculated without area dependence (Anderson et al. 1979).

The body of a grasshopper female is approximated by a rotational ellipsoid with half the body length as the semi-major axis (Samietz et al. 2005).

The diffuse fraction is corrected following Olyphant (1984).

Value

numeric predicted body (operative environmental) temperature (C).

References

Anderson RV, Tracy CR, Abramsky Z (1979). “Habitat Selection in Two Species of Short-Horned Grasshoppers. The Role of Thermal and Hydric Stresses.” Oecologia, 38(3), 359–374. doi:10.1007/BF00345194.

Erbs D, Klein S, Duffie J (1982). “Estimation of the diffuse radiation fraction for hourly, daily and monthly-average global radiation.” Solar Energy, 28, 293-302.

Gates DM (1962). “Leaf temperature and energy exchange.” Archiv fur Meteorologie, Geophysik und Bioklimatologie, Serie B volume, 12, 321-336.

Kingsolver JG (1983). “Thermoregulation and Flight in Colias Butterflies: Elevational Patterns and Mechanistic Limitations.” Ecology, 64(3), 534-545. doi:10.2307/1939973.

Lactin DJ, Johnson DL (1998). “Convective heat loss and change in body temperature of grasshopper and locust nymphs: Relative importance of wind speed, insect size and insect orientation.” Journal of Thermal Biology, 23(1), 5-13. ISSN 0306-4565, doi:10.1016/S0306-4565(97)00037-5, https://www.sciencedirect.com/science/article/pii/S0306456597000375.

Olyphant G (1984). “Insolation Topoclimates and Potential Ablation in Alpine Snow Accumulation Basins: Front Range, Colorado.” Water Resources Research, 20(4), 491-498.

Samietz J, Salser MA, Dingle H (2005). “Altitudinal variation in behavioural thermoregulation: local adaptation vs. plasticity in California grasshoppers.” Journal of Evolutionary Biology, 18(4), 1087-1096. doi:10.1111/j.1420-9101.2005.00893.x.

Swinbank WC (1963). “Long-wave radiation from clear skies.” Quarterly Journal of the Royal Meteorological Society, 89, 339-348.

Wong LT, Chow WK (2001). “Solar radiation model.” Applied Energy, 69(3), 191-224. ISSN 0306-2619, doi:10.1016/S0306-2619(01)00012-5, https://www.sciencedirect.com/science/article/pii/S0306261901000125.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tb_grasshopper(T_a       = 25, 
                 T_g       = 25,      
                 u         = 0.4, 
                 S         = 400, 
                 K_t       = 0.7, 
                 psi       = 30, 
                 l         = 0.02, 
                 Acondfact = 0.25, 
                 z         = 0.001, 
                 abs       = 0.7, 
                 r_g       = 0.3)

Operative Environmental Temperature of a Limpet

Description

The function estimates body temperatures (C, operative environmental temperatures) of a limpet based on Denny and Harley (2006).

Usage

Tb_limpet(T_a, T_r, l, h, S, u, psi, c, position = "anterior")

Arguments

T_a

numeric air temperature (C).

T_r

numeric rock surface temperature (C) in the sunlight.

l

numeric limpet length (anterior/posterior axis, m).

h

numeric limpet height (dorsal/ventral axis, m).

S

numeric solar irradiance (W m-2).

u

numeric wind speed (m s-1).

psi

numeric solar zenith angle (degrees). Can be calculated from zenith_angle function.

c

numeric fraction of the sky covered by cloud (proportion).

position

character direction of the limpet that is facing upwind. Options are "anterior", "posterior", and "broadside".

Details

The original equation uses a finite-difference approach where they divide the rock into series of chunks, and calculate the temperature at each node to derive the conductive heat. For simplification, here it takes the rock temperature as a parameter, and conductive heat is calculated as a product of the area, thermal conductivity of rock and the temperature difference between the rock and the body.

Limpets are simulated as cones following and using solar emissivity values from Campbell and Norman (1998).

The area of the limpet's shell (m2) is projected according to the direction at which sunlight strikes the organism (Pennell and Deignan 1989).

Air conductivity values (W m-1 K-1) are calculated following Denny and Harley (2006).

Value

numeric predicted body (operative environmental) temperature (C).

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

Denny MW, Harley CDG (2006). “Hot limpets: predicting body temperature in a conductance-mediated thermal system.” Journal of Experimental Biology, 209(13), 2409-2419. ISSN 0022-0949, doi:10.1242/jeb.02356.

Pennell S, Deignan J (1989). “Computing the Projected Area of a Cone.” SIAM Review, 31, 299-302. doi:10.1137/1031052.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tb_limpet(T_a      = 25, 
            T_r      = 30, 
            l        = 0.0176, 
            h        = 0.0122, 
            S        = 1300, 
            u        = 1, 
            psi      = 30, 
            c        = 1, 
            position = "anterior")

Operative Environmental Temperature of a Limpet Based on a Model by Helmuth

Description

The function predicts body temperatures (C, operative environmental temperatures) of a limpet. The function was provided by Brian Helmuth – although radiation and convection are altered from his original model – and based on Denny and Harley (2006).

Usage

Tb_limpetBH(T_a, T_r, l, h, S, u, s_aspect, s_slope, c)

Arguments

T_a

numeric air temperature (C).

T_r

numeric rock surface temperature (C) in the sunlight.

l

numeric limpet length (anterior/posterior axis) (m).

h

numeric limpet height (dorsal/ventral axis) (m).

S

numeric solar irradiance (W m-2).

u

numeric wind speed (m s-1).

s_aspect

numeric solar aspect angle (degree), the angle between the limpet's length dimension and the vector to the Sun. Generally between 70 and 110 degrees.

s_slope

numeric solar elevation angle (degree), the altitude of the sun, which is the angle between the horizon and the sun.

c

numeric fraction of the sky covered by clouds.

Details

The original equation uses a finite-difference approach where they divide the rock into series of chunks, and calculate the temperature at each node to derive the conductive heat. For simplification, here it takes the rock temperature as a parameter, and conductive heat is calculated by the product of the area, thermal conductivity of rock and the difference in temperatures of the rock and the body.

Limpets are simulated as cones following and using solar emissivity values from Campbell and Norman (1998).

The area of the limpet's shell (m2) is projected in the direction at which sunlight strikes the organism Pennell and Deignan (1989).

Air conductivity values (W m-1 K-1) are calculated following Denny and Harley (2006).

Value

numeric predicted body (operative environmental) temperature (C).

Author(s)

Brian Helmuth et al.

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

Denny MW, Harley CDG (2006). “Hot limpets: predicting body temperature in a conductance-mediated thermal system.” Journal of Experimental Biology, 209(13), 2409-2419. ISSN 0022-0949, doi:10.1242/jeb.02356.

Pennell S, Deignan J (1989). “Computing the Projected Area of a Cone.” SIAM Review, 31, 299-302. doi:10.1137/1031052.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tb_limpetBH(T_a = 25,
              T_r = 30,
              l = 0.0176,
              h = 0.0122,
              S = 1300,
              u = 1,
              s_aspect = 90,
              s_slope = 60,
              c = 1)

Operative Environmental Temperature of a Lizard

Description

The function estimates body temperature (C, operative environmental temperature) of a lizard based on Campbell and Norman (1998). The function was designed for Sceloporus lizards and described in Buckley (2008).

Usage

Tb_lizard(
  T_a,
  T_g,
  u,
  svl,
  m,
  psi,
  rho_s,
  elev,
  doy,
  sun = TRUE,
  surface = TRUE,
  a_s = 0.9,
  a_l = 0.965,
  epsilon_s = 0.965,
  F_d = 0.8,
  F_r = 0.5,
  F_a = 0.5,
  F_g = 0.5
)

Arguments

T_a

numeric air temperature (C).

T_g

numeric surface temperature (C).

u

numeric wind speed (m s-1).

svl

numeric lizard snout vent length (mm).

m

numeric lizard mass (g); note that it can be estimated as in mass_from_length: 3.55 x 10-5 x length3

psi

numeric solar zenith angle (degrees).

rho_s

numeric surface albedo (proportion). ~ 0.25 for grass, ~ 0.1 for dark soil, > 0.75 for fresh snow (Campbell and Norman 1998).

elev

numeric elevation (m).

doy

numeric day of year (1-366).

sun

logical indicates whether lizard is in sun (TRUE) or shade (FALSE).

surface

logical indicates whether lizard is on ground surface (TRUE) or above the surface (FALSE, e.g. in a tree).

a_s

numeric lizard solar absorptivity (proportion), a_s = 0.9 (Gates 1980) (Table 11.4).

a_l

numeric lizard thermal absorptivity (proportion), a_l = 0.965 (Bartlett and Gates 1967).

epsilon_s

numeric surface emissivity of lizards (proportion), epsilon_s = 0.965 (Bartlett and Gates 1967).

F_d

numeric the view factor between the surface of the lizard and diffuse solar radiation (proportion). i.e., the portion of the lizard surface that is exposed to diffuse solar radiation (Bartlett and Gates 1967).

F_r

numeric the view factor between the surface of the lizard and reflected solar radiation (proportion).

F_a

numeric the view factor between the surface of the lizard and atmospheric radiation (proportion).

F_g

numeric the view factor between the surface of the lizard and ground thermal radiation (proportion).

Details

The proportion of radiation that is direct is determined following Sears et al. (2011).

Boundary conductance uses a factor of 1.4 to account for increased convection (Mitchell 1976).

Value

T_e numeric predicted body (operative environmental) temperature (C).

References

Bartlett PN, Gates DM (1967). “The energy budget of a lizard on a tree trunk.” Ecology, 48, 316-322.

Buckley LB (2008). “Linking traits to energetics and population dynamics to predict lizard ranges in changing environments.” American Naturalist, 171(1), E1 - E19. doi:10.1086/523949, https://pubmed.ncbi.nlm.nih.gov/18171140/.

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

Gates DM (1980). Biophysical Ecology. Springer-Verlag, New York, NY, USA.

Mitchell JW (1976). “Heat transfer from spheres and other animal forms.” Biophysical Journal, 16(6), 561-569. ISSN 0006-3495, doi:10.1016/S0006-3495(76)85711-6.

Sears MW, Raskin E, Angilletta Jr. MJ (2011). “The World Is not Flat: Defining Relevant Thermal Landscapes in the Context of Climate Change.” Integrative and Comparative Biology, 51(5), 666-675. ISSN 1540-7063, doi:10.1093/icb/icr111, https://academic.oup.com/icb/article-pdf/51/5/666/1757893/icr111.pdf.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tb_lizard(T_a       = 25, 
            T_g       = 30, 
            u         = 0.1, 
            svl       = 60, 
            m         = 10, 
            psi       = 34, 
            rho_s     = 0.24, 
            elev      = 500, 
            doy       = 200, 
            sun       = TRUE, 
            surface   = TRUE, 
            a_s   = 0.9, 
            a_l   = 0.965, 
            epsilon_s = 0.965, 
            F_d       = 0.8, 
            F_r       = 0.5, 
            F_a       = 0.5, 
            F_g       = 0.5)

Operative Temperature of a Lizard Using Fei et al. (2012)

Description

The function predicts body temperature (C, operative environmental temperature) of a lizard based on Fei et al. (2012).

Usage

Tb_lizard_Fei(T_a, T_g, S, lw, shade, m, Acondfact, Agradfact)

Arguments

T_a

numeric air temperature at lizard height (C).

T_g

numeric surface temperature (C).

S

numeric total (direct + diffuse) solar radiation flux (W m-2).

lw

numeric downward flux of near-infrared radiation (W m-2).

shade

numeric proportion of shade.

m

numeric lizard mass (g).

Acondfact

numeric proportion of the lizard projected area that is in contact with the ground. Acondfact = 0.1 for standing and Acondfact = 0.4 for lying on ground.

Agradfact

numeric proportion of the lizard projected area exposed to radiation from the ground. Agradfact = 0.3 for standing and Agradfact = 0.0 for lying on ground.

Details

Thermal radiative flux is calculated following Fei et al. (2012) based on Bartlett and Gates (1967) and Porter et al. (1973).

Value

numeric predicted body (operative environmental) temperature (K).

Author(s)

Ofir Levy

References

Bartlett PN, Gates DM (1967). “The energy budget of a lizard on a tree trunk.” Ecology, 48, 316-322.

Fei T, Skidmore AK, Venus V, Wang T, Schlerf M, Toxopeus B, van Overjijk S, Bian M, Liu Y (2012). “A body temperature model for lizards as estimated from the thermal environment.” Journal of Thermal Biology, 37(1), 56-64. ISSN 0306-4565, doi:10.1016/j.jtherbio.2011.10.013, https://www.sciencedirect.com/science/article/pii/S0306456511001513.

Porter WP, Mitchell JW, Bekman A, DeWitt CB (1973). “Behavioral implications of mechanistic ecology: thermal and behavioral modeling of desert ectotherms and their microenvironments.” Oecologia, 13, 1-54.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tb_lizard_Fei(T_a       = 20,
                T_g       = 27,
                S         = 1300, 
                lw        = 60, 
                shade     = 0.5, 
                m         = 10.5, 
                Acondfact = 0.1, 
                Agradfact = 0.3)

Operative Environmental Temperature of a Mussel

Description

The function estimates body temperature (C, operative environmental temperature) of a mussel. The function implements a steady-state model, which assumes unchanging environmental conditions.

Usage

Tb_mussel(l, h, T_a, T_g, S, k_d, u, psi, cl, evap = FALSE, group = "solitary")

Arguments

l

numeric mussel length (anterior/posterior axis, m).

h

numeric mussel height (dorsal/ventral axis, m). It is reasonable to assume h = 0.5 * l.

T_a

numeric air temperature (C).

T_g

numeric ground temperature (C).

S

numeric direct solar flux density (W m-2).

k_d

numeric diffuse fraction, proportion of solar radiation that is diffuse.

u

numeric wind speed (m s-1).

psi

numeric solar zenith angle (degrees): can be calculated from zenith_angle.

cl

numeric fraction of the sky covered by cloud.

evap

logical Whether mussel is gaping to evaporatively cool. If TRUE, the function assumes a constant mass loss rate of 5 percent of the initial body mass per hour.

group

character; options are "aggregated": mussels living in beds; "solitary": solitary individual, anterior or posterior end facing upwind; and "solitary_valve": solitary individual, valve facing upwind.

Details

Thermal radiative flux is calculated following Helmuth (1998), Helmuth (1999), and Idso and Jackson (1969).

Value

numeric predicted body (operative environmental) temperature (C).

References

Helmuth B (1999). “Thermal biology of rocky intertidal mussels: quantifying body temperatures using climatological data.” Ecology, 80(1), 15-34. doi:10.2307/176977.

Helmuth BST (1998). “Intertidal Mussel Microclimates: Predicting the Body Temperature of a Sessile Invertebrate.” Ecological Monographs, 68(1), 51–74. ISSN 00129615, doi:10.2307/2657143.

Idso SB, Jackson RD (1969). “Thermal radiation from the atmosphere.” Journal of Geophysical Research (1896-1977), 74(23), 5397-5403. doi:10.1029/JC074i023p05397.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tb_mussel(l     = 0.1, 
          h     = 0.05, 
          T_a   = 25, 
          T_g   = 30, 
          S     = 500, 
          k_d   = 0.2, 
          u     = 2, 
          psi   = 30, 
          evap  = FALSE, 
          cl    = 0.5, 
          group = "solitary")

Humid Operative Environmental Temperature of a Salamander

Description

The function estimates the humid body temperature (C, operative environmental temperature) using an adaptation of Campbell and Norman (1998) described in Riddell et al. (2018).

Usage

Tb_salamander_humid(r_i, r_b, D, T_a, elev, e_a, e_s, Qabs, epsilon = 0.96)

Arguments

r_i

numeric internal (skin) resistance (s cm-1).

r_b

numeric boundary layer resistance (s cm-1).

D

numeric body diameter (m); can estimate as diameter = 0.0016*log(mass) + 0.0061 for mass(g).

T_a

numeric ambient temperature (C).

elev

numeric elevation (m).

e_a

numeric actual vapor pressure (kPa).

e_s

numeric saturation vapor pressure (kPa).

Qabs

numeric Solar and thermal radiation absorbed (W).

epsilon

numeric emissivity of salamander skin (proportion).

Value

numeric humid operative temperature (C).

Author(s)

Eric Riddell

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

Riddell EA, Odom JP, Damm JD, Sears MW (2018). “Plasticity reveals hidden resistance to extinction under climate change in the global hotspot of salamander diversity.” Science Advances, 4(4). doi:10.1126/sciadv.aar5471.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_snail(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tb_salamander_humid(r_i     = 4,
                      r_b     = 1,
                      D       = 0.01,
                      T_a     = 20,
                      elev    = 500,
                      e_a     = 2.0,
                      e_s     = 2.5,
                      Qabs    = 400,
                      epsilon = 0.96)

Operative Environmental Temperature of a Marine Snail

Description

The function estimates body temperature (C, operative environmental temperature) of a marine snail. The function implements a steady-state model, which assumes unchanging environmental conditions and is based on (Iacarella and Helmuth 2012). Body temperature and desiccation constrain the activity of Littoraria irrorata within the Spartina alterniflora canopy. The function was provided by Brian Helmuth and is a simplified version of the published model.

Usage

Tb_snail(temp, l, S, u, CC, WL, WSH)

Arguments

temp

numeric air temperature (C).

l

numeric snail length (m).

S

numeric direct solar flux density (W m-2).

u

numeric wind speed (m s-1).

CC

numeric fraction of the sky covered by cloud (0-1).

WL

numeric water loss rate (kg s-1), 5 percent loss of body mass over one hour is a reasonable maximum level (Helmuth 1999).

WSH

numeric wind sensor height (m).

Details

Thermal radiative flux is calculated following Helmuth (1998), Helmuth (1999), and Idso and Jackson (1969).

Value

numeric predicted body (operative environmental) temperature (C).

Author(s)

Brian Helmuth et al.

References

Helmuth B (1999). “Thermal biology of rocky intertidal mussels: quantifying body temperatures using climatological data.” Ecology, 80(1), 15-34. doi:10.2307/176977.

Helmuth BST (1998). “Intertidal Mussel Microclimates: Predicting the Body Temperature of a Sessile Invertebrate.” Ecological Monographs, 68(1), 51–74. ISSN 00129615, doi:10.2307/2657143.

Iacarella J, Helmuth B (2012). “Body temperature and desiccation constrain the activity of Littoraria irrorata within the Spartina alterniflora canopy.” Journal of Thermal Biology, 37(1). doi:10.1016/j.jtherbio.2011.10.003.

Idso SB, Jackson RD (1969). “Thermal radiation from the atmosphere.” Journal of Geophysical Research (1896-1977), 74(23), 5397-5403. doi:10.1029/JC074i023p05397.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tbed_mussel(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tb_snail(temp  = 25, 
           l     = 0.012, 
           S = 800, 
           u    = 1, 
           CC    = 0.5, 
           WL    = 0, 
           WSH   = 10)

Operative Environmental Temperature of a Mussel Bed

Description

The function estimates body temperature of a mussel (C). The function implements a steady-state model, which assumes unchanging environmental conditions. Based on Helmuth (1999).

Usage

Tbed_mussel(l, T_a, S, k_d, u, evap = FALSE, cl = NA)

Arguments

l

numeric mussel length (anterior/posterior axis) (m).

T_a

numeric air temperature at 4 m above ground (C).

S

numeric direct solar flux density (W m-2).

k_d

numeric diffuse fraction, proportion of solar radiation that is diffuse.

u

numeric wind speed at 4 m above ground (m s-1).

evap

logical Are mussels gaping to evaporatively cool? If TRUE, assumes constant mass loss rate of 5 percent of initial body mass per hour.

cl

numeric fraction of the sky covered by cloud, optional.

Details

Conduction is considered negligible due to a small area of contact.

Thermal radiative flux is calculated following Helmuth (1998), Helmuth (1999), and Idso and Jackson (1969).

Value

numeric predicted body (operative environmental) temperature (C).

References

Helmuth B (1999). “Thermal biology of rocky intertidal mussels: quantifying body temperatures using climatological data.” Ecology, 80(1), 15-34. doi:10.2307/176977.

Helmuth BST (1998). “Intertidal Mussel Microclimates: Predicting the Body Temperature of a Sessile Invertebrate.” Ecological Monographs, 68(1), 51–74. ISSN 00129615, doi:10.2307/2657143.

Idso SB, Jackson RD (1969). “Thermal radiation from the atmosphere.” Journal of Geophysical Research (1896-1977), 74(23), 5397-5403. doi:10.1029/JC074i023p05397.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tsoil(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tbed_mussel(l    = 0.1, 
             T_a  = 25, 
             S    = 500, 
             k_d  = 0.2, 
             u    = 1, 
             evap = FALSE)

Convert Among Temperature Scales

Description

The function converts temperatures among Celsius, Fahrenheit, and Kelvin (J. Blischak et al. 2016).

Usage

fahrenheit_to_kelvin(temperature)

kelvin_to_celsius(temperature)

celsius_to_kelvin(temperature)

fahrenheit_to_celsius(temperature)

Arguments

temperature

numeric temperature (Celsius, Fahrenheit, or Kelvin).

Value

numeric temperature (Celsius, Fahrenheit, or Kelvin).

References

J. Blischak, D. Chen, H. Dashnow, Haine D (2016). Software Carpentry: Programming with R. doi:10.5281/zenodo.57541, Version 2016.06, June 2016.

See Also

Other utility functions: airpressure_from_elev(), azimuth_angle(), day_of_year(), daylength(), dec_angle(), solar_noon(), zenith_angle()

Examples

kelvin_to_celsius(temperature = 270)
  fahrenheit_to_kelvin(temperature = 85)
  fahrenheit_to_celsius(temperature = 85)
  celsius_to_kelvin(temperature = -10)

Gaussian-Quadratic Function Thermal Performance Curve

Description

The function constructs a thermal performance curve by combining as a Gaussian function to describe the rise in performance up to the optimal temperature and a quadratic decline to zero performance at critical thermal maxima and higher temperatures (Deutsch et al. 2008).

Usage

TPC(T_b, T_opt, CT_min, CT_max)

Arguments

T_b

numeric vector of temperature range (C).

T_opt

numeric thermal optima (C), the temperature at which peak performance occurs.

CT_min, CT_max

numeric critical thermal minimum and maximum (C), the lower and upper temperature limits for performance.

Value

performance

References

Deutsch CA, Tewksbury JJ, Huey RB, Sheldon KS, Ghalambor CK, Haak DC, Martin PR (2008). “Impacts of climate warming on terrestrial ectotherms across latitude.” Proceedings of the National Academy of Science of the United States of America, 105, 6668-6672. doi:10.1073/pnas.0709472105.

Examples

TPC(T_b    = 0:60, 
      T_opt  = 30, 
      CT_min = 10, 
      CT_max = 40)

Beta Function Thermal Performance Curve

Description

The function constructs a thermal performance curve based on a beta function (Asbury and Angilletta 2010).

Usage

TPC_beta(T_b, shift = -1, breadth = 0.1, aran = 0, tolerance = 43, skew = 0.7)

Arguments

T_b

numeric temperature (C).

shift

numeric mode of the thermal performance curve.

breadth

numeric breadth of the thermal performance curve.

aran

numeric scale performance value. If 0, no scaling; if 1, include a thermodynamic effect on mean performance.

tolerance

numeric maximal breath (C) of the thermal performance curve.

skew

numeric skewness of the thermal performance curve (0-1).

Value

numeric performance.

References

Asbury DA, Angilletta MJ (2010). “Thermodynamic effects on the evolution of performance curves.” American Naturalist, 176(2), E40-E49. doi:10.1086/653659.

Examples

TPC_beta(T_b       = 0:60, 
           shift     = -1, 
           breadth   = 0.1, 
           aran      = 0, 
           tolerance = 43, 
           skew      = 0.7)

Approximate Soil Temperature

Description

The function estimates soil temperature (C) at a given depth and hour by approximating diurnal variation as sinusoidal. The function is adapted from Campbell and Norman (1998) and described in Riddell et al. (2018).

Usage

Tsoil(T_g_max, T_g_min, hour, depth)

Arguments

T_g_max

numeric daily maximum soil surface temperature (C).

T_g_min

numeric daily minimum soil surface temperature (C).

hour

numeric hour of the day.

depth

numeric depth (cm).

Value

numeric soil temperature (C).

Author(s)

Eric Riddell

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

Riddell EA, Odom JP, Damm JD, Sears MW (2018). “Plasticity reveals hidden resistance to extinction under climate change in the global hotspot of salamander diversity.” Science Advances, 4(4). doi:10.1126/sciadv.aar5471.

See Also

Other biophysical models: Grashof_number_Gates(), Grashof_number(), Nusselt_from_Grashof(), Nusselt_from_Reynolds(), Nusselt_number(), Prandtl_number(), Qconduction_animal(), Qconduction_substrate(), Qconvection(), Qemitted_thermal_radiation(), Qevaporation(), Qmetabolism_from_mass_temp(), Qmetabolism_from_mass(), Qnet_Gates(), Qradiation_absorbed(), Qthermal_radiation_absorbed(), Reynolds_number(), T_sky(), Tb_CampbellNorman(), Tb_Gates2(), Tb_Gates(), Tb_butterfly(), Tb_grasshopper(), Tb_limpetBH(), Tb_limpet(), Tb_lizard_Fei(), Tb_lizard(), Tb_mussel(), Tb_salamander_humid(), Tb_snail(), Tbed_mussel(), actual_vapor_pressure(), boundary_layer_resistance(), external_resistance_to_water_vapor_transfer(), free_or_forced_convection(), heat_transfer_coefficient_approximation(), heat_transfer_coefficient_simple(), heat_transfer_coefficient(), saturation_vapor_pressure(), saturation_water_vapor_pressure()

Examples

Tsoil(T_g_max = 30,
        T_g_min = 15,
        hour    = 12,
        depth   = 5)

Saturation Vapor Pressure

Description

The function calculates saturation vapor pressure for a given air temperature. The program is based on equations from List (1971) and code implementation from NicheMapR (Kearney and Porter 2017; Kearney and Porter 2020).

Usage

VAPPRS(db)

Arguments

db

numeric Dry bulb temperature (C)

Value

esat numeric Saturation vapor pressure (Pa)

References

Kearney MR, Porter WP (2017). “NicheMapR - an R package for biophysical modelling: the microclimate model.” Ecography, 40, 664-674. doi:10.1111/ecog.02360.

Kearney MR, Porter WP (2020). “NicheMapR - an R package for biophysical modelling: the ectotherm and Dynamic Energy Budget models.” Ecography, 43(1), 85-96. doi:10.1111/ecog.04680.

List RJ (1971). “Smithsonian Meteorological Tables.” Smithsonian Miscellaneous Collections, 114(1), 1-527. https://repository.si.edu/handle/10088/23746.

Examples

VAPPRS(db = 30)

Organism Volume from Length

Description

The function estimates volume (m3) from length (m) for a variety of taxa following Mitchell (1976).

Usage

volume_from_length(l, taxon)

Arguments

l

numeric length (m).
Use snout-vent length for lizards and frogs.

taxon

character taxon of organism, current choices: "lizard", "frog", "sphere".

Details

Relationships come from

  • Lizards: Norris (1965)

  • Frogs: Tracy (1972)

  • Sphere: Mitchell (1976)

Value

numeric volume (m3).

References

Mitchell JW (1976). “Heat transfer from spheres and other animal forms.” Biophysical Journal, 16(6), 561-569. ISSN 0006-3495, doi:10.1016/S0006-3495(76)85711-6.

Norris KS (1965). “Color adaptation in desert reptiles and its thermal relationships.” In Symposium on Lizard Ecology, 162- 229. U. Missouri Press.

Tracy CR (1972). “Newton's Law: Its Application for Expressing Heat Losses from Homeotherms.” BioScience, 22(11), 656-659. ISSN 0006-3568, doi:10.2307/1296267.

See Also

Other allometric functions: mass_from_length(), proportion_silhouette_area_shapes(), proportion_silhouette_area(), surface_area_from_length(), surface_area_from_mass(), surface_area_from_volume()

Examples

volume_from_length(l     = 0.05,  
                     taxon = "lizard")
  volume_from_length(l     = 0.05,   
                     taxon = "frog")
  volume_from_length(l     = 0.05,   
                     taxon = "sphere")

Properties of Wet Air

Description

The function calculates several properties of humid air described as output variables below. The program is based on equations from List (1971) and code implementation from NicheMapR (Kearney and Porter 2017; Kearney and Porter 2020).
WETAIR must be used in conjunction with VAPPRS. Input variables are shown below. See Details.

Usage

WETAIR(db, wb = db, rh = 0, dp = 999, bp = 101325)

Arguments

db

numeric dry bulb temperature (C).

wb

numeric wet bulb temperature (C).

rh

numeric relative humidity (%).

dp

numeric dew point temperature (C).

bp

numeric barometric pressure (Pa).

Details

The user must supply known values for DB and BP (BP at one standard atmosphere is 101,325 pascals). Values for the remaining variables are determined by whether the user has either (1) psychrometric data (WB or RH), or (2) hygrometric data (DP):

  • Psychrometric data: If WB is known but not RH, then set RH = -1 and DP = 999. If RH is known but not WB then set WB = 0 and DP = 999

  • Hygrometric data: If DP is known, set WB = 0 and RH = 0

Value

Named list with elements

  • e: numeric saturation vapor pressure (Pa)

  • vd: numeric vapor density (kg m-3)

  • rw: numeric mixing ratio (kg kg-1)

  • tvir: numeric virtual temperature (K)

  • tvinc: numeric virtual temperature increment (K)

  • denair: numeric density of the air (kg m-3)

  • cp: numeric specific heat of air at constant pressure (J kg-1 K-1)

  • wtrpot: numeric water potential (Pa)

  • rh: numeric relative humidity (%)

References

Kearney MR, Porter WP (2017). “NicheMapR - an R package for biophysical modelling: the microclimate model.” Ecography, 40, 664-674. doi:10.1111/ecog.02360.

Kearney MR, Porter WP (2020). “NicheMapR - an R package for biophysical modelling: the ectotherm and Dynamic Energy Budget models.” Ecography, 43(1), 85-96. doi:10.1111/ecog.04680.

List RJ (1971). “Smithsonian Meteorological Tables.” Smithsonian Miscellaneous Collections, 114(1), 1-527. https://repository.si.edu/handle/10088/23746.

Examples

WETAIR(db = 30, 
         wb = 28, 
         rh = 60, 
         bp = 100 * 1000)

Wind Speed at a Specific Height Under Neutral Conditions

Description

The function calculates wind speed (m s-1) at a specified height (m) within a boundary layer near the surface. The profile assumes neutral conditions. The velocity profile is the neutral profile described by Sellers (1965). Function is equations (2) and (3) of Porter et al. (1973).

Usage

wind_speed_profile_neutral(u_r, zr, z0, z)

Arguments

u_r

numeric wind velocity (m s-1) at reference height.

zr

numeric initial reference height (m).

z0

numeric surface roughness (m).

z

numeric height to scale (m).

Value

numeric windspeed (m s-1).

References

Porter WP, Mitchell JW, Bekman A, DeWitt CB (1973). “Behavioral implications of mechanistic ecology: thermal and behavioral modeling of desert ectotherms and their microenvironments.” Oecologia, 13, 1-54.

Sellers WD (1965). Physical climatology. University of Chicago Press, Chicago, IL, USA.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), air_temp_profile(), degree_days(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), monthly_solar_radiation(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_segment()

Examples

wind_speed_profile_neutral(u_r = 0.1, 
                             zr  = 0.1, 
                             z0  = 0.2, 
                             z   = 0.15)

Wind Speed at a Specified Height

Description

The function calculates wind speed (m s-1) at a specified height (m). The function estimates a three segment velocity and temperature profile based on user-specified, experimentally determined values for 3 roughness heights and reference heights. Multiple heights are appropriate in heterogenous areas with, for example, a meadow, bushes, and rocks. Implements the MICROSEGMT routine from NicheMapR as described in Kearney and Porter (2017).

Usage

wind_speed_profile_segment(u_r, zr, z0, z)

Arguments

u_r

numeric a vector of wind speeds (m s-1) at the 3 reference heights.

zr

numeric a vector of 3 reference heights (m).

z0

numeric a vector of 3 experimentally determined roughness heights (m).

z

numeric height to scale (m).

Value

numeric wind speed (m s-1).

References

Kearney MR, Porter WP (2017). “NicheMapR - an R package for biophysical modelling: the microclimate model.” Ecography, 40, 664-674. doi:10.1111/ecog.02360.

See Also

Other microclimate functions: air_temp_profile_neutral(), air_temp_profile_segment(), air_temp_profile(), degree_days(), direct_solar_radiation(), diurnal_radiation_variation(), diurnal_temp_variation_sineexp(), diurnal_temp_variation_sinesqrt(), diurnal_temp_variation_sine(), monthly_solar_radiation(), partition_solar_radiation(), proportion_diffuse_solar_radiation(), solar_radiation(), surface_roughness(), wind_speed_profile_neutral()

Examples

wind_speed_profile_segment(u_r = c(0.01, 0.025, 0.05), 
                             zr  = c(0.05, 0.25, 0.5), 
                             z0  = c(0.01, 0.15, 0.2), 
                             z   = 0.3)

Zenith Angle

Description

The function calculates the zenith angle, the location of the sun as an angle (in degrees) measured from vertical (Campbell and Norman 1998).

Usage

zenith_angle(doy, lat, lon, hour, offset = NA)

Arguments

doy

numeric day of year (1-366). This can be obtained from standard date via day_of_year.

lat

numeric latitude (decimal degrees).

lon

numeric longitude (decimal degrees).

hour

numeric hour of the day (0-24).

offset

numeric the number of hours to add to UTC (Coordinated Universal Time) to get local time (improves accuracy but not always necessary). Optional. Defaults to NA.

Value

numeric zenith angle (degrees)

References

Campbell GS, Norman JM (1998). Introduction to environmental biophysics, 2nd ed. edition. Springer, New York. ISBN 0387949372.

See Also

Other utility functions: airpressure_from_elev(), azimuth_angle(), day_of_year(), daylength(), dec_angle(), solar_noon(), temperature conversions

Examples

zenith_angle(doy  = 112, 
               lat  = 47.61, 
               lon  = -122.33, 
               hour = 12)