
4. Pedigree Analysis and Population Genetics
Source:vignettes/pedigree-analysis.Rmd
pedigree-analysis.RmdThis vignette summarizes a practical workflow for pedigree analysis
with visPedigree, with emphasis on the interpretation of
the main indicators used in breeding and conservation genetics.
The discussion is organized around five questions:
- How complete and deep is the pedigree?
- How long is the generation cycle?
- Has genetic diversity been eroded by unequal founder use, bottlenecks, or drift?
- How large is the effective population size under different definitions?
- Are relationship, inbreeding, subpopulation structure, and gene flow under control?
1. Setup and Data Preparation
Different package datasets are used for different analytical tasks.
-
deep_pedis useful for pedigree depth and diversity summaries. -
big_family_size_pedcontains aYearcolumn and is suitable for generation intervals. -
small_pedis convenient for relationship examples. -
inbred_pedis convenient for inbreeding and classification examples.
library(visPedigree)
library(data.table)
data(deep_ped, package = "visPedigree")
data(big_family_size_ped, package = "visPedigree")
data(small_ped, package = "visPedigree")
data(inbred_ped, package = "visPedigree")
tp_deep <- tidyped(deep_ped)
tp_small <- tidyped(small_ped)
tp_inbred <- tidyped(inbred_ped)2. Pedigree Overview with pedstats()
pedstats() provides a compact structural summary with
three components:
-
$summary: pedigree size and parental structure. -
$ecg: pedigree completeness and ancestral depth. -
$gen_intervals: generation intervals, only when a usable time column is available.
Here deep_ped has no explicit birth-date column.
Accordingly, pedstats() returns $summary and
$ecg, whereas $gen_intervals remains
NULL.
stats_deep <- pedstats(tp_deep)
stats_deep$summary
#> N NSire NDam NFounder MaxGen
#> <int> <int> <int> <int> <int>
#> 1: 4399 483 554 138 13
tail(stats_deep$ecg)
#> Ind ECG FullGen MaxGen
#> <char> <num> <num> <num>
#> 1: K110997Q 7.616211 5 12
#> 2: K110997Z 6.722656 5 12
#> 3: K110998Q 4.606445 2 12
#> 4: K110998Z 6.417969 5 12
#> 5: K110999Q 7.345215 5 12
#> 6: K110999Z 6.875977 5 12
stats_deep$gen_intervals
#> NULLPedigree depth and pedigree time are related but distinct quantities.
The former is addressed by pedecg(), whereas the latter is
addressed by pedgenint().
3. Pedigree Completeness with pedecg()
Equivalent complete generations (ECG) summarize the amount of ancestral information available for each individual:
where is the number of generations between individual and ancestor . ECG increases with both pedigree depth and pedigree completeness.
In practice:
-
ECGcombines depth and completeness. -
FullGencounts how many fully known ancestral generations exist. -
MaxGenrecords the deepest known ancestral path.
ECG is especially useful because it provides the depth
adjustment required by realized effective population size estimators
based on inbreeding or coancestry.
ecg_deep <- pedecg(tp_deep)
ecg_deep[order(-ECG)][1:10]
#> Ind ECG FullGen MaxGen
#> <char> <num> <num> <num>
#> 1: K110034Q 9.078125 6 12
#> 2: K110052L 9.078125 6 12
#> 3: K110060H 9.078125 6 12
#> 4: K110069Z 9.078125 6 12
#> 5: K110097Q 9.078125 6 12
#> 6: K110118M 9.078125 6 12
#> 7: K110131M 9.078125 6 12
#> 8: K110138M 9.078125 6 12
#> 9: K110155M 9.078125 6 12
#> 10: K110165M 9.078125 6 124. Generation Intervals with pedgenint()
Generation interval is the age of a parent at the birth of its offspring:
pedgenint() estimates this quantity for seven pathway
summaries:
-
SS,SD,DS,DD: sex-specific gametic pathways. -
SO,DO: sex-independent sire-offspring and dam-offspring pathways. -
Average: all parent-offspring pairs combined.
The function accepts Date/POSIXct columns,
date strings, or numeric years. When only an integer year is available,
it is converted internally to YYYY-07-01 as a mid-year
approximation.
tp_time <- tidyped(big_family_size_ped)
genint_year <- pedgenint(tp_time, timevar = "Year", unit = "year")
#> Numeric time column detected. Converting to Date (YYYY-07-01). For finer precision, convert to Date beforehand.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
genint_year
#> Pathway N Mean SD
#> <char> <int> <num> <num>
#> 1: Average 280512 1.001093 0.03770707
#> 2: DD 607 1.164398 0.37080261
#> 3: DO 140256 1.001093 0.03770714
#> 4: DS 507 1.196959 0.39776787
#> 5: SD 607 1.164398 0.37080261
#> 6: SO 140256 1.001093 0.03770714
#> 7: SS 507 1.196959 0.39776787The optional cycle parameter adds GenEquiv,
which compares the observed mean interval with a target breeding
cycle:
Values larger than 1 indicate that the realized generation interval exceeds the target cycle.
genint_cycle <- pedgenint(tp_time, timevar = "Year", unit = "year", cycle = 1.2)
#> Numeric time column detected. Converting to Date (YYYY-07-01). For finer precision, convert to Date beforehand.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
genint_cycle[Pathway %in% c("SS", "SD", "DS", "DD", "Average")]
#> Pathway N Mean SD GenEquiv
#> <char> <int> <num> <num> <num>
#> 1: Average 280512 1.001093 0.03770707 0.8342445
#> 2: DD 607 1.164398 0.37080261 0.9703321
#> 3: DS 507 1.196959 0.39776787 0.9974660
#> 4: SD 607 1.164398 0.37080261 0.9703321
#> 5: SS 507 1.196959 0.39776787 0.99746605. Subpopulation Structure with pedsubpop()
Before interpreting diversity or relationship metrics, it is useful to check whether the pedigree forms a single connected population or a mixture of separate components.
pedsubpop() has two modes:
-
by = NULL: summarize disconnected pedigree components. -
by = "...": summarize the pedigree by an existing grouping variable.
ped_demo <- data.table(
Ind = c("A", "B", "C", "D", "E", "F", "G", "H"),
Sire = c(NA, NA, "A", NA, NA, "E", NA, NA),
Dam = c(NA, NA, "B", NA, NA, NA, NA, NA),
Sex = c("male", "female", "male", "female", "male", "female", "male", "female"),
Batch = c("L1", "L1", "L1", "L1", "L2", "L2", "L3", "L3")
)
tp_demo <- tidyped(ped_demo)
pedsubpop(tp_demo)
#> Group N N_Sire N_Dam N_Founder
#> <char> <int> <num> <num> <int>
#> 1: GP1 3 1 1 2
#> 2: GP2 2 1 0 1
#> 3: Isolated 3 0 0 3
pedsubpop(tp_demo, by = "Batch")
#> Group N N_Sire N_Dam N_Founder
#> <char> <int> <int> <int> <int>
#> 1: L1 4 1 1 3
#> 2: L3 2 0 0 2
#> 3: L2 2 1 0 1This is a compact summary tool. When the actual split pedigree
objects are needed for downstream analysis, use
splitped().
6. Diversity Indicators with pediv()
pediv() is the integrated diversity summary. It
combines:
- founder and ancestor contributions (
fe,fa), - founder genome equivalents (
fg), - three effective population size estimates (
Ne).
All of these quantities depend on the definition of the
reference population. In the present example, the
reference population is defined as the most recent two generations in
deep_ped.
div_res <- pediv(tp_deep, reference = ref_pop, top = 10, seed = 42L)
#> Calculating founder and ancestor contributions...
#> Calculating founder contributions...
#> Calculating ancestor contributions (Boichard's iterative algorithm)...
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Calculating Ne (coancestry) and fg...
#> Calculating Ne (inbreeding)...
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Calculating Ne (demographic)...
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
div_res$summary
#> NRef NFounder feH fe NAncestor faH fa fafe fg
#> <int> <int> <num> <num> <int> <num> <num> <num> <num>
#> 1: 3471 157 92.0943 64.73344 94 60.1444 44.12033 0.681569 19.17965
#> MeanCoan NSampledCoan NeCoancestry NeInbreeding NeDemographic
#> <num> <int> <num> <num> <num>
#> 1: 0.0260693 1000 124.5124 98.00425 374.2021
div_res$ancestors
#> Ind Contrib CumContrib Rank
#> <char> <num> <num> <int>
#> 1: K500I804 0.05248848 0.05248848 1
#> 2: K60GXQ91 0.04525893 0.09774741 2
#> 3: K60NXQ91 0.04525893 0.14300634 3
#> 4: K600532I 0.04445765 0.18746399 4
#> 5: K700S069 0.03927182 0.22673581 5
#> 6: K900D251 0.03622875 0.26296456 6
#> 7: K700S011 0.02906223 0.29202679 7
#> 8: K700U416 0.02890017 0.32092697 8
#> 9: K700U650 0.02890017 0.34982714 9
#> 10: K500I869 0.02831497 0.37814211 106.1 fe, fa, and fg: what do
they measure?
These three indicators describe diversity loss from complementary angles.
Effective number of founders (fe)
where
is the expected contribution of founder
to the reference population. fe answers the question: if
all founders had contributed equally, how many founders would generate
the same diversity?
Use fe when the main concern is unequal founder
representation.
Effective number of ancestors (fa)
where
is the marginal contribution of ancestor
after removing the contributions already explained by more influential
ancestors. fa is usually smaller than fe when
bottlenecks occurred.
Use fa when the main concern is genetic
bottlenecks caused by a limited set of influential
ancestors.
Founder genome equivalents (fg)
In its simplest interpretation,
where
is the mean coancestry of the reference population. In
visPedigree, fg is estimated from a
diagonal-corrected mean coancestry:
where is the size of the full reference cohort, is the mean off-diagonal additive relationship among sampled individuals, and is their mean inbreeding coefficient.
Use fg when the main concern is the amount of
founder genome still surviving after unequal use, bottlenecks, and
drift. In practice, fg is often the most
conservative diversity indicator.
6.2 Shannon-entropy effective numbers: feH and
faH
The classical fe and fa are based on
(Hill number order
),
which is disproportionately influenced by a few large contributions.
visPedigree also reports the information-theoretic
counterpart at order
,
derived from Shannon entropy:
$$ f_e^{(H)} = \exp\!\Bigl(-\sum_{i=1}^{f} p_i \ln p_i\Bigr), \qquad f_a^{(H)} = \exp\!\Bigl(-\sum_{j=1}^{a} q_j \ln q_j\Bigr) $$
These are the exponentials of Shannon entropy, which give a diversity
measure that counts all contributors — including rare ones — more
equitably than the classical
metric. They appear in the pediv() summary as
feH and faH, and in the
pedcontrib() summary as f_e_H and
f_a_H.
The three diversity orders satisfy a monotone ordering:
where is the total number of founders and is the number of ancestors considered. When is much larger than 1, many minor founders still carry genetic material into the reference population but are invisible to the classical metric.
# Shannon metrics are included in pediv() output
div_res$summary[, .(NFounder, feH, fe, NAncestor, faH, fa)]
#> NFounder feH fe NAncestor faH fa
#> <int> <num> <num> <int> <num> <num>
#> 1: 157 92.0943 64.73344 94 60.1444 44.12033
# The ratio feH/fe reveals long-tail founder diversity
div_res$summary[, .(rho_founder = feH / fe, rho_ancestor = faH / fa)]
#> rho_founder rho_ancestor
#> <num> <num>
#> 1: 1.42267 1.36319
# pedcontrib() provides the same metrics via its summary
contrib_res <- pedcontrib(tp_deep, reference = ref_pop, mode = "both")
#> Calculating founder contributions...
#> Calculating ancestor contributions (Boichard's iterative algorithm)...
contrib_res$summary[c("n_founder", "f_e_H", "f_e", "n_ancestor", "f_a_H", "f_a")]
#> $n_founder
#> [1] 157
#>
#> $f_e_H
#> [1] 92.0943
#>
#> $f_e
#> [1] 64.73344
#>
#> $n_ancestor
#> [1] 94
#>
#> $f_a_H
#> [1] 60.1444
#>
#> $f_a
#> [1] 44.120337. Effective Population Size with pedne() and
pediv()
pediv() reports three complementary effective population
size definitions. Each addresses a distinct biological question.
7.1 Demographic Ne
where and are the numbers of contributing males and females. This is the easiest estimate to understand, but it ignores realized pedigree structure, inbreeding, and drift. It is therefore often optimistic.
7.2 Inbreeding Ne
This estimator uses the realized rate of inbreeding. ECG
standardizes individuals with different pedigree depths. Use this
estimate when the primary concern is the rate of inbreeding
accumulation.
7.3 Coancestry Ne
This estimator is based on the rate of coancestry among members of the reference population. Because it captures the accumulation of relatedness before it is fully expressed as realized inbreeding, it is often the strictest and most sensitive warning signal in managed breeding populations.
8. Average Relationship Trends with pedrel()
pedrel() summarizes the mean off-diagonal additive
relationship within groups:
where
.
A higher MeanRel means that individuals within the group
are, on average, more related by descent.
This summary is useful for tracking relatedness by generation or year.
tp_small$BirthYear <- 2010 + tp_small$Gen
rel_by_gen <- pedrel(tp_small, by = "Gen")
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
rel_by_gen
#> Gen NTotal NUsed MeanRel
#> <int> <int> <int> <num>
#> 1: 1 9 9 0.0000000
#> 2: 2 5 5 0.2000000
#> 3: 3 7 7 0.1547619
#> 4: 4 3 3 0.1145833
#> 5: 5 2 2 0.0468750
#> 6: 6 2 2 0.5507812The reference argument lets you focus on a subset of
interest inside each group, such as candidate breeders.
ref_ids <- c("Z1", "Z2", "X", "Y")
pedrel(tp_small, by = "Gen", reference = ref_ids)
#> Warning in FUN(X[[i]], ...): Group '1' has less than 2 individuals after
#> applying 'reference', returning NA_real_.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning in FUN(X[[i]], ...): Group '2' has less than 2 individuals after
#> applying 'reference', returning NA_real_.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning in FUN(X[[i]], ...): Group '3' has less than 2 individuals after
#> applying 'reference', returning NA_real_.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning in FUN(X[[i]], ...): Group '4' has less than 2 individuals after
#> applying 'reference', returning NA_real_.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Warning: Subsetting removed parent records. Result is a plain data.table, not a tidyped.
#> Use tidyped(tp, cand = ids, trace = "up") to extract a valid sub-pedigree.
#> Gen NTotal NUsed MeanRel
#> <int> <int> <int> <num>
#> 1: 1 9 0 NA
#> 2: 2 5 0 NA
#> 3: 3 7 0 NA
#> 4: 4 3 0 NA
#> 5: 5 2 2 0.0468750
#> 6: 6 2 2 0.5507812MeanRel and coancestry-based Ne are
conceptually linked, but they are not identical summaries.
MeanRel is an average additive relationship within a group,
whereas coancestry-based Ne is derived from the rate of
increase in coancestry across pedigree depth.
9. Inbreeding Trends with inbreed() and
pedfclass()
Wright’s inbreeding coefficient is the probability that the two alleles of an individual are identical by descent. A practical starting point is to inspect the mean trend by generation.
tp_inbred_f <- inbreed(tp_inbred)
f_trend <- tp_inbred_f[, .(MeanF = mean(f, na.rm = TRUE)), by = Gen]
f_trend
#> Gen MeanF
#> <int> <num>
#> 1: 1 0.0000
#> 2: 2 0.0000
#> 3: 3 0.2500
#> 4: 4 0.2500
#> 5: 5 0.4375For classification by inbreeding severity, pedfclass()
can be applied directly to the pedigree object. If the inbreeding
coefficient is not yet present, the function computes it internally.
pedfclass(tp_inbred)
#> Calculating inbreeding coefficients...
#> Key: <FClass>
#> FClass Count Percentage
#> <ord> <int> <num>
#> 1: F = 0 4 57.14286
#> 2: 0 < F <= 0.0625 0 0.00000
#> 3: 0.0625 < F <= 0.125 0 0.00000
#> 4: 0.125 < F <= 0.25 2 28.57143
#> 5: F > 0.25 1 14.28571Custom reporting thresholds can be supplied through
breaks.
pedfclass(tp_inbred, breaks = c(0.03125, 0.0625, 0.125, 0.25))
#> Calculating inbreeding coefficients...
#> Key: <FClass>
#> FClass Count Percentage
#> <ord> <int> <num>
#> 1: F = 0 4 57.14286
#> 2: 0 < F <= 0.03125 0 0.00000
#> 3: 0.03125 < F <= 0.0625 0 0.00000
#> 4: 0.0625 < F <= 0.125 0 0.00000
#> 5: 0.125 < F <= 0.25 2 28.57143
#> 6: F > 0.25 1 14.28571The default thresholds correspond approximately to familiar pedigree scenarios:
- : half-sib mating.
- : avuncular or grandparent-grandchild mating.
- : full-sib or parent-offspring mating.
10. Gene Flow and Partial Inbreeding
10.1 pedancestry(): founder-line proportions
pedancestry() tracks how founder groups contribute to
later descendants. This is useful when founders are labeled by line,
strain, or geographic origin.
ped_line <- data.table(
Ind = c("A", "B", "C", "D", "E", "F", "G"),
Sire = c(NA, NA, NA, NA, "A", "C", "E"),
Dam = c(NA, NA, NA, NA, "B", "D", "F"),
Sex = c("male", "female", "male", "female", "male", "female", "male"),
Line = c("Line1", "Line1", "Line2", "Line2", NA, NA, NA)
)
tp_line <- tidyped(ped_line)
anc <- pedancestry(tp_line, foundervar = "Line")
anc
#> Ind Line1 Line2
#> <char> <num> <num>
#> 1: A 1.0 0.0
#> 2: B 1.0 0.0
#> 3: C 0.0 1.0
#> 4: D 0.0 1.0
#> 5: E 1.0 0.0
#> 6: F 0.0 1.0
#> 7: G 0.5 0.510.2 pedpartial(): which ancestors explain
inbreeding?
pedpartial() decomposes total inbreeding into
contributions from selected ancestors. It is useful for identifying
which ancestors are most responsible for the observed inbreeding.
partial_f <- pedpartial(tp_inbred, ancestors = c("A", "B"))
#> Calculating partial inbreeding for 2 ancestors...
partial_f
#> Ind A B
#> <char> <num> <num>
#> 1: A 0.00000 0.00000
#> 2: B 0.00000 0.00000
#> 3: C 0.00000 0.00000
#> 4: D 0.00000 0.00000
#> 5: E 0.12500 0.12500
#> 6: F 0.00000 0.25000
#> 7: G 0.15625 0.2812511. Practical Interpretation
One useful interpretation order is:
- Check pedigree depth and completeness with
pedecg(). - Check generation timing with
pedgenint(). - Quantify diversity loss with
fe,fa, andfg; compare withfeHandfaHto assess long-tail founder value. - Compare demographic, inbreeding, and coancestry-based
Ne. - Monitor
MeanRelandMeanFover time. - Use
pedsubpop(),pedancestry(), andpedpartial()to diagnose structure, introgression, and bottlenecks.
References
- Boichard, D., Maignel, L., & Verrier, É. (1997). The value of using probabilities of gene origin to measure genetic variability in a population. Genetics Selection Evolution, 29(1), 5-23.
- Caballero, A., & Toro, M. A. (2000). Interrelations between effective population size and other pedigree tools for the management of conserved populations. Genetical Research, 75(3), 331-343.
- Cervantes, I., Goyache, F., Molina, A., Valera, M., & Gutiérrez, J. P. (2011). Estimation of effective population size from the rate of coancestry in pedigreed populations. Journal of Animal Breeding and Genetics, 128(1), 56-63.
- Gutiérrez, J. P., Cervantes, I., Molina, A., Valera, M., & Goyache, F. (2008). Individual increase in inbreeding allows estimating effective sizes from pedigrees. Genetics Selection Evolution, 40(4), 359-378.
- Gutiérrez, J. P., Cervantes, I., & Goyache, F. (2009). Improving the estimation of realized effective population sizes in farm animals. Journal of Animal Breeding and Genetics, 126(4), 327-332.
- Hill, M. O. (1973). Diversity and evenness: a unifying notation and its consequences. Ecology, 54(2), 427-432.
- Jost, L. (2006). Entropy and diversity. Oikos, 113(2), 363-375.
- Lacy, R. C. (1989). Analysis of founder representation in pedigrees: founder equivalents and founder genome equivalents. Zoo Biology, 8(2), 111-123.
- Wright, S. (1922). Coefficients of inbreeding and relationship. The American Naturalist, 56(645), 330-338.
- Wright, S. (1931). Evolution in Mendelian populations. Genetics, 16(2), 97-159.