Public methods
Method new()
Starts the process of building a new simulation
by creating a new SimParam object and assigning a founder
population to the class. It is recommended that you save the
object with the name "SP", because subsequent functions will
check your global enviroment for an object of this name if
their simParam arguments are NULL. This allows you to call
these functions without explicitly supplying a simParam
argument with every call.
Usage
SimParam$new(founderPop)
Arguments
founderPopan object of MapPop-class
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
Method setTrackPed()
Sets pedigree tracking for the simulation.
By default pedigree tracking is turned off. When turned on,
the pedigree of all individuals created will be tracked,
except those created by hybridCross. Turning
off pedigree tracking will turn off recombination tracking
if it is turned on.
Usage
SimParam$setTrackPed(isTrackPed, force = FALSE)
Arguments
isTrackPedshould pedigree tracking be on.
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing.
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$setTrackPed(TRUE)
Method setTrackRec()
Sets recombination tracking for the simulation.
By default recombination tracking is turned off. When turned
on recombination tracking will also turn on pedigree tracking.
Recombination tracking keeps records of all individuals created,
except those created by hybridCross, because their
pedigree is not tracked.
Usage
SimParam$setTrackRec(isTrackRec, force = FALSE)
Arguments
isTrackRecshould recombination tracking be on.
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing.
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$setTrackRec(TRUE)
Method resetPed()
Resets the internal lastId, the pedigree
and recombination tracking (if in use) to the
supplied lastId. Be careful using this function because
it may introduce a bug if you use individuals from
the deleted portion of the pedigree.
Usage
SimParam$resetPed(lastId = 0L)
Arguments
lastIdlast ID to include in pedigree
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
#Create population
pop = newPop(founderPop, simParam=SP)
pop@id # 1:10
#Create another population after reseting pedigree
SP$resetPed()
pop2 = newPop(founderPop, simParam=SP)
pop2@id # 1:10
Method restrSegSites()
Sets restrictions on which segregating sites
can serve as SNP and/or QTL.
Usage
SimParam$restrSegSites(
minQtlPerChr = NULL,
minSnpPerChr = NULL,
overlap = FALSE,
minSnpFreq = NULL
)
Arguments
minQtlPerChrthe minimum number of segSites for QTLs.
Can be a single value or a vector values for each
chromosome.
minSnpPerChrthe minimum number of segSites for SNPs.
Can be a single value or a vector values for each
chromosome.
overlapshould SNP and QTL sites be allowed to overlap.
minSnpFreqminimum allowable frequency for SNP loci.
No minimum SNP frequency is used if value is NULL.
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$restrSegSites(minQtlPerChr=5, minSnpPerChr=5)
Method setGender()
Changes how gender is used in the simulation.
The default gender of a simulation is "no". To add gender
to the simulation, run this function with "yes_sys" or
"yes_rand". The value "yes_sys" will systematically assign
gender to newly created individuals as first male, then female.
Thus, odd numbers of individuals will have one more male than
female. The value "yes_rand" will randomly assign gender to
individuals.
Usage
SimParam$setGender(gender, force = FALSE)
Arguments
genderacceptable value are "no", "yes_sys", or
"yes_rand"
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing.
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$setGender("yes_sys")
Method addSnpChip()
Randomly assigns eligble SNPs to a SNP chip
Usage
SimParam$addSnpChip(nSnpPerChr, minSnpFreq = NULL, refPop = NULL)
Arguments
nSnpPerChrnumber of SNPs per chromosome.
Can be a single value or nChr values.
minSnpFreqminimum allowable frequency for SNP loci.
If NULL, no minimum frequency is used.
refPopreference population for calculating SNP
frequency. If NULL, the founder population is used.
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$addSnpChip(10)
Method addStructuredSnpChip()
Randomly selects the number of snps in structure and then
assigns them to chips based on structure
Usage
SimParam$addStructuredSnpChip(nSnpPerChr, structure, force = FALSE)
Arguments
nSnpPerChrnumber of SNPs per chromosome.
Can be a single value or nChr values.
structurea matrix. Rows are snp chips, columns are chips.
If value is true then that snp is on that chip.
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing.
Method addTraitA()
Randomly assigns eligble QTLs for one or more additive traits.
If simulating more than one trait, all traits will be pleiotrophic
with correlated additive effects.
Usage
SimParam$addTraitA(
nQtlPerChr,
mean = 0,
var = 1,
corA = NULL,
gamma = FALSE,
shape = 1,
force = FALSE
)
Arguments
nQtlPerChrnumber of QTLs per chromosome. Can be a single value or nChr values.
meana vector of desired mean genetic values for one or more traits
vara vector of desired genetic variances for one or more traits
corAa matrix of correlations between additive effects
gammashould a gamma distribution be used instead of normal
shapethe shape parameter for the gamma distribution
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing.
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$addTraitA(10)
Method addTraitAD()
Randomly assigns eligble QTLs for one or more traits with dominance.
If simulating more than one trait, all traits will be pleiotrophic
with correlated effects.
Usage
SimParam$addTraitAD(
nQtlPerChr,
mean = 0,
var = 1,
meanDD = 0,
varDD = 0,
corA = NULL,
corDD = NULL,
useVarA = TRUE,
gamma = FALSE,
shape = 1,
force = FALSE
)
Arguments
nQtlPerChrnumber of QTLs per chromosome. Can be a single value or nChr values.
meana vector of desired mean genetic values for one or more traits
vara vector of desired genetic variances for one or more traits
meanDDmean dominance degree
varDDvariance of dominance degree
corAa matrix of correlations between additive effects
corDDa matrix of correlations between dominance degrees
useVarAtune according to additive genetic variance if true. If
FALSE, tuning is performed according to total genetic variance.
gammashould a gamma distribution be used instead of normal
shapethe shape parameter for the gamma distribution
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing.
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$addTraitAD(10, meanDD=0.5)
Method addTraitAG()
Randomly assigns eligble QTLs for one ore more additive GxE traits.
If simulating more than one trait, all traits will be pleiotrophic
with correlated effects.
Usage
SimParam$addTraitAG(
nQtlPerChr,
mean = 0,
var = 1,
varGxE = 1e-06,
varEnv = 0,
corA = NULL,
corGxE = NULL,
gamma = FALSE,
shape = 1,
force = FALSE
)
Arguments
nQtlPerChrnumber of QTLs per chromosome. Can be a single value or nChr values.
meana vector of desired mean genetic values for one or more traits
vara vector of desired genetic variances for one or more traits
varGxEa vector of total genotype-by-environment variances for the traits
varEnva vector of environmental variances for one or more traits
corAa matrix of correlations between additive effects
corGxEa matrix of correlations between GxE effects
gammashould a gamma distribution be used instead of normal
shapethe shape parameter for the gamma distribution
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing.
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$addTraitAG(10, varGxE=2)
Method addTraitADG()
Randomly assigns eligble QTLs for a trait with dominance and GxE.
Usage
SimParam$addTraitADG(
nQtlPerChr,
mean = 0,
var = 1,
varEnv = 1e-06,
varGxE = 1e-06,
meanDD = 0,
varDD = 0,
corA = NULL,
corDD = NULL,
corGxE = NULL,
useVarA = TRUE,
gamma = FALSE,
shape = 1,
force = FALSE
)
Arguments
nQtlPerChrnumber of QTLs per chromosome. Can be a single
value or nChr values.
meana vector of desired mean genetic values for one or more traits
vara vector of desired genetic variances for one or more traits
varEnva vector of environmental variances for one or more traits
varGxEa vector of total genotype-by-environment variances for the traits
meanDDmean dominance degree
varDDvariance of dominance degree
corAa matrix of correlations between additive effects
corDDa matrix of correlations between dominance degrees
corGxEa matrix of correlations between GxE effects
useVarAtune according to additive genetic variance if true
gammashould a gamma distribution be used instead of normal
shapethe shape parameter for the gamma distribution
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing.
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$addTraitADG(10, meanDD=0.5, varGxE=2)
Method addTraitAE()
Randomly assigns eligble QTLs for one or more additive and epistasis
traits. If simulating more than one trait, all traits will be pleiotrophic
with correlated additive effects.
Usage
SimParam$addTraitAE(
nQtlPerChr,
mean = 0,
var = 1,
relAA = 0,
corA = NULL,
corAA = NULL,
useVarA = TRUE,
gamma = FALSE,
shape = 1,
force = FALSE
)
Arguments
nQtlPerChrnumber of QTLs per chromosome. Can be a single value or nChr values.
meana vector of desired mean genetic values for one or more traits
vara vector of desired genetic variances for one or more traits
relAAthe relative value of additive-by-additive variance compared
to additive variance in a diploid organism with allele frequency 0.5
corAa matrix of correlations between additive effects
corAAa matrix of correlations between additive-by-additive effects
useVarAtune according to additive genetic variance if true. If
FALSE, tuning is performed according to total genetic variance.
gammashould a gamma distribution be used instead of normal
shapethe shape parameter for the gamma distribution
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing.
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
Method addTraitADE()
Randomly assigns eligble QTLs for one or more traits with dominance and
epistasis. If simulating more than one trait, all traits will be pleiotrophic
with correlated effects.
Usage
SimParam$addTraitADE(
nQtlPerChr,
mean = 0,
var = 1,
meanDD = 0,
varDD = 0,
relAA = 0,
corA = NULL,
corDD = NULL,
corAA = NULL,
useVarA = TRUE,
gamma = FALSE,
shape = 1,
force = FALSE
)
Arguments
nQtlPerChrnumber of QTLs per chromosome. Can be a single value or nChr values.
meana vector of desired mean genetic values for one or more traits
vara vector of desired genetic variances for one or more traits
meanDDmean dominance degree
varDDvariance of dominance degree
relAAthe relative value of additive-by-additive variance compared
to additive variance in a diploid organism with allele frequency 0.5
corAa matrix of correlations between additive effects
corDDa matrix of correlations between dominance degrees
corAAa matrix of correlations between additive-by-additive effects
useVarAtune according to additive genetic variance if true. If
FALSE, tuning is performed according to total genetic variance.
gammashould a gamma distribution be used instead of normal
shapethe shape parameter for the gamma distribution
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing.
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$addTraitADE(10)
Method addTraitAEG()
Randomly assigns eligble QTLs for one or more additive and epistasis
GxE traits. If simulating more than one trait, all traits will be pleiotrophic
with correlated effects.
Usage
SimParam$addTraitAEG(
nQtlPerChr,
mean = 0,
var = 1,
relAA = 0,
varGxE = 1e-06,
varEnv = 0,
corA = NULL,
corAA = NULL,
corGxE = NULL,
useVarA = TRUE,
gamma = FALSE,
shape = 1,
force = FALSE
)
Arguments
nQtlPerChrnumber of QTLs per chromosome. Can be a single value or nChr values.
meana vector of desired mean genetic values for one or more traits
vara vector of desired genetic variances for one or more traits
relAAthe relative value of additive-by-additive variance compared
to additive variance in a diploid organism with allele frequency 0.5
varGxEa vector of total genotype-by-environment variances for the traits
varEnva vector of environmental variances for one or more traits
corAa matrix of correlations between additive effects
corAAa matrix of correlations between additive-by-additive effects
corGxEa matrix of correlations between GxE effects
useVarAtune according to additive genetic variance if true. If
FALSE, tuning is performed according to total genetic variance.
gammashould a gamma distribution be used instead of normal
shapethe shape parameter for the gamma distribution
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing.
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$addTraitAEG(10, varGxE=2)
Method addTraitADEG()
Randomly assigns eligble QTLs for a trait with dominance,
epistasis and GxE.
Usage
SimParam$addTraitADEG(
nQtlPerChr,
mean = 0,
var = 1,
varEnv = 1e-06,
varGxE = 1e-06,
meanDD = 0,
varDD = 0,
relAA = 0,
corA = NULL,
corDD = NULL,
corAA = NULL,
corGxE = NULL,
useVarA = TRUE,
gamma = FALSE,
shape = 1,
force = FALSE
)
Arguments
nQtlPerChrnumber of QTLs per chromosome. Can be a single
value or nChr values.
meana vector of desired mean genetic values for one or more traits
vara vector of desired genetic variances for one or more traits
varEnva vector of environmental variances for one or more traits
varGxEa vector of total genotype-by-environment variances for the traits
meanDDmean dominance degree
varDDvariance of dominance degree
relAAthe relative value of additive-by-additive variance compared
to additive variance in a diploid organism with allele frequency 0.5
corAa matrix of correlations between additive effects
corDDa matrix of correlations between dominance degrees
corAAa matrix of correlations between additive-by-additive effects
corGxEa matrix of correlations between GxE effects
useVarAtune according to additive genetic variance if true
gammashould a gamma distribution be used instead of normal
shapethe shape parameter for the gamma distribution
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing.
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$addTraitADEG(10, meanDD=0.5, varGxE=2)
Method manAddTrait()
Manually add a new trait to the simulation.
Usage
SimParam$manAddTrait(
lociMap,
varA = NA_real_,
varG = NA_real_,
varE = NA_real_,
force = FALSE
)
Arguments
lociMapa new object descended from
LociMap-class
varAthe value for varA in the base population, optional
varGthe value for varG in the base population, optional
varEdefault error variance for phenotype, optional
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing
Method switchTrait()
Switch a trait in the simulation.
Usage
SimParam$switchTrait(
traitPos,
lociMap,
varA = NA_real_,
varG = NA_real_,
varE = NA_real_,
force = FALSE
)
Arguments
traitPosan integer indicate which trait to switch
lociMapa new object descended from
LociMap-class
varAthe value for varA in the base population, optional
varGthe value for varG in the base population, optional
varEdefault error variance for phenotype, optional
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing
Method removeTrait()
Remove a trait from the simulation
Usage
SimParam$removeTrait(traits, force = FALSE)
Arguments
traitsan integer vector indicating which traits to remove
forceshould the check for a running simulation be
ignored. Only set to TRUE if you know what you are doing
Method setVarE()
Defines a default value for error
variances in the simulation.
Usage
SimParam$setVarE(h2 = NULL, H2 = NULL, varE = NULL)
Arguments
h2a vector of desired narrow-sense heritabilities
H2a vector of desired broad-sense heritabilities
varEa vector of error variances
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$addTraitA(10)
SP$setVarE(h2=0.5)
Method setCorE()
Defines a correlation structure for default
error variances. You must call setVarE first to define
the default error variances.
Usage
SimParam$setCorE(corE)
Arguments
corEa correlation matrix for the error variances
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$addTraitA(10, mean=c(0,0), var=c(1,1), corA=diag(2))
SP$setVarE(varE=c(1,1))
E = 0.5*diag(2)+0.5 #Positively correlated error
SP$setCorE(E)
Method rescaleTraits()
Linearly scales all traits to achieve desired
values of means and variances in the founder population.
Usage
SimParam$rescaleTraits(
mean = 0,
var = 1,
varEnv = 0,
varGxE = 1e-06,
useVarA = TRUE
)
Arguments
meana vector of new trait means
vara vector of new trait variances
varEnva vector of new environmental variances
varGxEa vector of new GxE variances
useVarAtune according to additive genetic variance if true
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$addTraitA(10)
#Create population
pop = newPop(founderPop, simParam=SP)
meanG(pop)
#Change mean to 1
SP$rescaleTraits(mean=1)
#Run resetPop for change to take effect
pop = resetPop(pop, simParam=SP)
meanG(pop)
Method setRecombRatio()
Set the relative recombination rates between males
and females. This allows for gender specific recombination rates,
under the assumption of equivalent recombination landscapes.
Usage
SimParam$setRecombRatio(femaleRatio)
Arguments
femaleRatiorelative ratio of recombination in females compared to
males. A value of 2 indicate twice as much recombination in females. The
value must be greater than 0. (default is 1)
Examples
#Create founder haplotypes
founderPop = quickHaplo(nInd=10, nChr=1, segSites=10)
#Set simulation parameters
SP = SimParam$new(founderPop)
SP$setRecombRatio(2) #Twice as much recombination in females
Method switchGenMap()
Replaces existing genetic map.
Usage
SimParam$switchGenMap(genMap, centromere = NULL)
Arguments
genMapa list of length nChr containing
numeric vectors for the position of each segregating
site on a chromosome.
centromerea numeric vector of centromere
positions. If NULL, the centromere are assumed to
be metacentric.
Method switchFemaleMap()
Replaces existing female genetic map.
Usage
SimParam$switchFemaleMap(genMap, centromere = NULL)
Arguments
genMapa list of length nChr containing
numeric vectors for the position of each segregating
site on a chromosome.
centromerea numeric vector of centromere
positions. If NULL, the centromere are assumed to
be metacentric.
Method switchMaleMap()
Replaces existing male genetic map.
Usage
SimParam$switchMaleMap(genMap, centromere = NULL)
Arguments
genMapa list of length nChr containing
numeric vectors for the position of each segregating
site on a chromosome.
centromerea numeric vector of centromere
positions. If NULL, the centromere are assumed to
be metacentric.
Method addToRec()
For internal use only.
Usage
SimParam$addToRec(hist)
Arguments
histnew recombination history
Method updateLastId()
For internal use only.
Usage
SimParam$updateLastId(lastId)
Arguments
lastIdlast ID assigned
Method addToPed()
For internal use only.
Usage
SimParam$addToPed(lastId, mother, father, isDH)
Arguments
lastIdID of last individual
mothervector of mother IDs
fathervector of father IDs
isDHvector of DH indicators
Method clone()
The objects of this class are cloneable with this method.
Usage
SimParam$clone(deep = FALSE)
Arguments
deepWhether to make a deep clone.