trans_behav | R Documentation |
This model uses the transient heat budget models (i.e. accounting for heat storage and hence lag-effects of body mass) to simulate thermoregulatory behaviour of a diurnally active ectotherm. It uses a set of events to break out of the ordinary differential equation solver of the transient heat budget (onelump_var or twolump functions) to simulate thermoregulation around setpoints, specifically the transition from sitting in the shade overnight to basking in the sun perpendicular to the sun's rays (T_B_min), from basking to foraging in the open in an 'average' posture (T_F_min), from foraging in the open to moving into the shade (T_F_max), and then transitioning to basking in the afternoon and finally retreating to the shade for the evening. It also computes the operative temperature (Te) of a non-thermoregulating animal in the open, i.e. the steady state body temperature the animal would come to if it had zero heat capacity, and the body temperature of a non thermoregulating animal in the open accounting for the effect of thermal mass. The function computes a series of summary statistics about the length of basking and foraging bouts, activity time, time spent above high temperature thresholds (T_F_max and CT_max), and the extremes of body temperature.
trans_behav(
Tc_init = rep(20, 60),
Ts_init = Tc_init + 0.1,
To_init = Tc_init + 0.2,
Ww_g = 500,
T_F_min = 33,
T_F_max = 38,
T_B_min = 25,
CT_max = 43,
rho_body = 932,
x_shell = 0.001,
lump = 1,
q = 0,
c_body = 3073,
c_body_inner = c_body,
c_body_outer = c_body,
k_flesh = 0.5,
k_inner = k_flesh,
k_outer = k_flesh,
emis = 0.95,
alpha = 0.85,
geom = 2,
shape_b = 1/5,
shape_c = 1/5,
shape_coefs = c(10.4713, 0.688, 0.425, 0.85, 3.798, 0.683, 0.694, 0.743),
posture = "n",
orient = 1,
fatosk = 0.4,
fatosb = 0.4,
dyn_q = 0,
fluid = 0,
alpha_sub = 0.2,
pdif = 0.1,
shade = 90,
metout = metout,
shadmet = shadmet,
soil = soil,
shadsoil = shadsoil,
press = 101325
)
Tc_init |
= Tairf(1), initial temperature (°C) Organism shape, 0-5, Determines whether standard or custom shapes/surface area/volume relationships are used: 0=plate, 1=cyl, 2=ellips, 3=lizard (desert iguana), 4=frog (leopard frog), 5=custom (see details) |
Ts_init |
= Tc_init + 0.1, initial shell temperature (°C) |
To_init |
= Tc_init + 0.2, initial surface temperature (°C) |
Ww_g |
= 500, weight (g) |
T_F_min |
= 33, # minimum foraging body temperature threshold (°C) |
T_F_max |
= 38, # maximum foraging body temperature threshold (°C) |
T_B_min |
= 25, # basking body temperature threshold (°C) |
CT_max |
= 43, # critical thermal maximum (°C) |
rho_body |
= 932, animal density (kg/m3) |
x_shell |
= 0.001, shell thickness, m |
lump |
= 1, one lump (1) or two lump (2) model? |
q |
= 0, metabolic rate (W/m3) |
c_body |
= 3073, specific heat of flesh (J/kg-°C) |
c_body_inner |
= c_body, Specific heat of flesh J/(kg-°C) |
c_body_outer |
= c_body, Specific heat of outer shell J/(kg-°C) |
k_flesh |
= 0.5, conductivity of flesh (W/m-K) |
k_inner |
= k_flesh, Thermal conductivity of inner shell (W/m-K, range: 0.412-2.8) |
k_outer |
= k_flesh, Thermal conductivity of outer shell (W/m-K, range: 0.412-2.8) |
emis |
= 0.95, emissivity of skin (-) |
alpha |
= 0.85, animal solar absorptivity (-) |
geom |
= 2, Organism shape, 0-5, Determines whether standard or custom shapes/surface area/volume relationships are used: 0=plate, 1=cyl, 2=ellips, 3=lizard (desert iguana), 4=frog (leopard frog), 5=custom (see parameter 'shape_coeffs') |
shape_b |
= 1/5, Proportionality factor (-) for going from volume to area, represents ratio of width:height for a plate, length:diameter for cylinder, b axis:a axis for ellipsoid |
shape_c |
= 1/5, Proportionality factor (-) for going from volume to area, represents ratio of length:height for a plate, c axis:a axis for ellipsoid |
shape_coefs |
= c(10.4713,.688,0.425,0.85,3.798,.683,0.694,.743), Custom surface area coefficients. Operates if posture = 5, and consists of 4 pairs of values representing the parameters a and b of a relationship AREA=a*Ww_g^b, where AREA is in cm2 and Ww_g is in g. The first pair are a and b for total surface area, then a and b for ventral area, then for sillhouette area normal to the sun, then sillhouette area perpendicular to the sun |
posture |
= 'n', pointing normal 'n' or parallel 'p' to the sun's rays, or average 'a'? |
orient |
= 1, does the object orient toward the sun? (0,1) |
fatosk |
= 0.4, solar configuration factor to sky (-) |
fatosb |
= 0.4, solar configuration factor to substrate (-) |
dyn_q |
= 0, dynamic metabolic heat generation as a function of temperature, based on approxfun called qf (1) or constant based on parameter q (0) |
fluid |
= 0, fluid type, air (0) or water (1) |
alpha_sub |
= 0.2, substrate solar reflectivity, decimal percent |
pdif |
= 0.1, proportion of solar energy that is diffuse (rather than direct beam) |
shade |
= 90, maximum shade level (%) |
metout |
= metout, aboveground minimum shade microclimate output table from NicheMapR's microclimate model |
shadmet |
= shadmet, metout, aboveground maximum shademicroclimate output table from NicheMapR's microclimate model |
soil |
= soil, minimum shade soil temperature output table from NicheMapR's microclimate model |
shadsoil |
= shadsoil, maximum shade soil temperature output table from NicheMapR's microclimate model |
press |
= 101325, air pressure (Pa) Outputs: day_results variables:
act_window variables:
sum_stats variables:
|
library(NicheMapR)
# define animal biophysical functional traits
Ww_g <- 500 # wet weight (g)
Usrhyt <- 0.05 # height of animal (mid-point) above ground (m)
alpha <- 0.85 # solar absorptivity (-)
T_F_min <- 33 # minimum foraging Tb (deg C)
T_F_max <- 43 # maximum foraging Tb (deg C)
T_B_min <- 18 # basking Tb, moving from shade to sun (deg C)
CT_max <- 48 # critical thermal maximum (deg C)
shape_b <- 1/5 # shape coefficient a, -
shape_c <- 1/5 # shape coefficient b, -
rho_body <- 1000 # animal density, kg/m3
c_body <- 3073 # heat capacity (J/kg-C)
q <- 0 # metabolic rate, W/m3
k_flesh <- 0.5 # thermal conductivity of flesh, W/mK
geom <- 2 # shape, -
# get microclimate data
loc <- c(130, -25)
maxshade <- 90
micro <- micro_global(loc = loc, maxshade = maxshade, Usrhyt = Usrhyt) # run the model with default location and settings
metout <- as.data.frame(micro$metout) # above ground microclimatic conditions, min shade
soil <- as.data.frame(micro$soil) # soil temperatures, minimum shade
shadmet <- as.data.frame(micro$shadmet) # above ground microclimatic conditions, min shade
shadsoil <- as.data.frame(micro$shadsoil) # soil temperatures, minimum shade
# get air pressure
elevation <- micro$elev
press <- 101325 * ((1 - (0.0065 * elevation / 288)) ^ (1 / 0.190284))
mons <- c("January", "February", "March", "April", "May", "June", "July", "August", "September", "October", "November", "December")
DOYs <- unique(metout$DOY)
# loop through each month and run transient model with behaviour
for(i in 1:12){
# subset current month
metout_in <- subset(metout, DOY == DOYs[i])
shadmet_in <- subset(shadmet, DOY == DOYs[i])
soil_in <- subset(soil, DOY == DOYs[i])
shadsoil_in <- subset(shadmet, DOY == DOYs[i])
# run transient behavioural simulation
trans <- trans_behav(Ww_g = Ww_g, alpha = alpha, T_F_min = T_F_min, T_F_max = T_F_max,
CT_max = CT_max, T_B_min = T_B_min, geom = geom, shape_b = shape_b, shape_c = shape_c,
rho_body = rho_body, k_flesh = k_flesh, q = q, lump = 1,
metout = metout_in, shadmet = shadmet_in, soil = soil_in, shadsoil = shadsoil_in,
press = press, alpha_sub = 1 - micro$REFL, shade = micro$maxshade)
results <- as.data.frame(trans$day_results)
sum_stats <- as.data.frame(trans$sum_stats)
act_window <- as.data.frame(trans$act_window)
# collate
if(i == 1){
all_act_window <- act_window
}else{
all_act_window <- rbind(all_act_window, act_window)
}
results$hours <- results$time / 3600
# plot hourly results for the current day
plot(results$Tb_open ~ results$hours, type = 'l', ylim = c(-10, 80), col = 'grey', xaxs = 'i', ylab = "temperature, deg C", xlab = "time",
main = paste0(if(length(loc) == 2){paste("lon", loc[1], "lat", loc[2])}else{loc}, ", ", mons[i], ", ", Ww_g," g"), xlim = c(0, 23))
grid(nx = 23, ny = 0, col = "lightgray", lty = "dotted", lwd = par("lwd"), equilogs = TRUE)
abline(T_F_max, 0, col = 'red', lty = 2)
abline(T_F_min, 0, col = 'light blue', lty = 2)
abline(CT_max, 0, col = 'red')
points(results$T_air_shd ~ results$hours, type = 'l', col = 'blue')
points(results$Tb ~ results$hours, type = 'l', col = 'orange', lty = 1, lwd = 2)
text(3, 60, paste0("bouts ", round(sum_stats$bouts_sun, 0)), cex = 1)
text(3, 65, paste0("maximum bout ", round(sum_stats$max_foraging_bout_sun / 60, 1), " hrs"), cex = 1)
text(3, 70, paste0("total activity ", round(sum_stats$sum_activity_sun / 60, 1), " hrs"), cex = 1)
}
# make seasonal activity plot
all_act_window$ZEN <- metout$ZEN
all_act_window$DOY <- metout$DOY
foraging<-subset(all_act_window, forage_sun > 0)
night<-subset(all_act_window, ZEN==90)
with(night, plot(time ~ DOY, pch=15, cex = 2, xlim = c(1, 365), col = 'dark blue', xlab = 'day of year', ylab = 'hour of day', main = "Seasonal Activity Plot, Sun"))
with(foraging, points(time ~ DOY, pch = 15, cex = forage_sun / 30, col = 'orange'))
foraging<-subset(all_act_window, forage_shd > 0)
with(night, plot(time ~ DOY, pch=15, cex = 2, xlim = c(1, 365), col = 'dark blue', xlab = 'day of year', ylab = 'hour of day', main = "Seasonal Activity Plot, Shade"))
with(foraging, points(time ~ DOY, pch = 15, cex = forage_shd / 30, col = 'orange'))
mtext(text = paste0('Seasonal Activity Plot, ', if(length(loc) == 2){paste("lon", loc[1], "lat", loc[2])}else{loc}, " ", Ww_g," g"), outer = TRUE, side = 3, line = 0)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.