mapply_calc_post_prob_presence: Mapply version of calc_post_prob_presence()

Description Usage Arguments Details Value Note Author(s) References See Also Examples

View source: R/BioGeoBEARS_detection_v1.R

Description

This function applies calc_post_prob_presence to all cells of the input matrices obs_target_species and obs_all_species. These matrices obviously must have the same dimensions.

Usage

1
2
3
4
  mapply_calc_post_prob_presence(prior_prob_presence = 0.01,
    obs_target_species, obs_all_species,
    mean_frequency = 0.1, dp = 1, fdp = 0,
    print_progress = "")

Arguments

prior_prob_presence

The prior probability of presence, i.e. when no detection or taphonomic control data whatsoever is available. Default is set to 0.01 which expresses my totally uninformed bias that in whatever your data is, your species of interest probably doesn't live in the typical area you are looking at.

obs_target_species

A scalar or column/vector/matrix of detection counts, e.g. as produced from the output from read_detections.

obs_all_species

A scalar or column/vector/matrix of detection counts, e.g. as produced from the output from read_controls.

mean_frequency

This is the proportion of samples from the taphonomic control group that will truly be from this OTU, GIVEN that the OTU is present. This could be estimated, but a decent first guess is (total # samples of OTU of interest / total # of samples in the taphonomic control group where the OTU is known to be present). All that is really needed is some reasonable value, such that more sampling without detection lowers the likelihood of the data on the hypothesis of true presence, and vice versa. This value can only be 1 when the number of detections = the number of taphonomic control detections, for every OTU and area. This is the implicit assumption in e.g. standard historical biogeography analyses in LAGRANGE or BioGeoBEARS.

dp

The detection probability. This is the per-sample probability that you will correctly detect the OTU in question, when you are looking at it. Default is 1, which is the implicit assumption in standard analyses.

fdp

The false detection probability. This is probability of falsely concluding a detection of the OTU of interest occurred, when in fact the specimen was of something else. The default is 0, which assumes zero error rate, i.e. the assumption being made in all historical biogeography analyses that do not take into account detection probability. These options are being included for completeness, but it may not be wise to try to infer mean_frequency, dp and fdp all at once due to identifiability issues (and estimation of fdp may take a very large amount of data). However, fixing some of these parameters to reasonable values can allow the user to effectively include beliefs about the uncertainty of the input data into the analysis, if desired.

print_progress

If not the default (""), print whatever is in print_progress, followed by a space (for error checking/surveying results).

Details

The inputs are the same as for calc_post_prob_presence, except that obs_target_species and obs_all_species can be matrices.

Value

pp_df A matrix of the posterior probability of presence, given the prior probability, the model parameters, and the data.

Note

Go BEARS!

Author(s)

Nicholas J. Matzke matzke@berkeley.edu

References

http://phylo.wikidot.com/matzke-2013-international-biogeography-society-poster http://en.wikipedia.org/wiki/Bayes'_theorem

Matzke_2012_IBS

Bottjer_Jablonski_1988

Bayes_1763

See Also

calc_obs_like, calc_post_prob_presence, mapply_calc_obs_like Pdata_given_rangerow, mapply, tiplikes_wDetectionModel

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
# Calculate posterior probability of presence in an area,
# given a dp (detection probability) and detection model.

# soft-coded input files
extdata_dir = np(system.file("extdata", package="BioGeoBEARS"))
detects_fn = np(paste(extdata_dir, "/Psychotria_detections_v1.txt", sep=""))
controls_fn = np(paste(extdata_dir, "/Psychotria_controls_v1.txt", sep=""))

OTUnames=NULL
areanames=NULL
tmpskip=0

detects_df = read_detections(detects_fn, OTUnames=NULL, areanames=NULL, tmpskip=0)
controls_df = read_controls(controls_fn, OTUnames=NULL, areanames=NULL, tmpskip=0)

detects_df
controls_df
detects_df / controls_df


# Calculate data likelihoods, and posterior probability of presence=TRUE
mean_frequency=0.1
dp=1
fdp=0

mapply_calc_obs_like(truly_present=TRUE, obs_target_species=detects_df,
obs_all_species=controls_df, mean_frequency, dp, fdp)

mapply_calc_obs_like(truly_present=FALSE, obs_target_species=detects_df,
obs_all_species=controls_df, mean_frequency, dp, fdp)

mapply_calc_post_prob_presence(prior_prob_presence=0.01,
obs_target_species=detects_df,
obs_all_species=controls_df, mean_frequency, dp, fdp)

BioGeoBEARS documentation built on May 29, 2017, 8:36 p.m.