PolyOrigin Documentation
Documentation for PolyOrigin.jl
's internal interface.
Contents
Index
PolyOrigin.polyMapRefine!
PolyOrigin.polyPhase
PolyOrigin.polyPhase!
PolyOrigin.polyReconstruct!
PolyOrigin.savegenodata
PolyOrigin.savegenoprob
PolyOrigin.setAbsPhase!
Internal Interface
PolyOrigin.polyPhase
— FunctionpolyPhase(genofile,pedfile, delimchar=',', missingstring="NA",
commentstring="#", keyargs...)
performs parental phasing and return phasedgeno::PolyGeno with phasedgeno.parentgeno being phased.
Positional arguments
genofile::AbstractString
: filename for genotypic data file.
pedfile::AbstractString
: filename for population pedigree.
see readPolyGeno
for the requirements of genofile
and pedfile
.
Keyword arguments
delimchar::AbstractChar=','
: text delimiter.
missingstring::AbstractString="NA"
: string code for missing value.
commentstring::AbstractString="#"
: rows that begins with commentstring will be ignored.
see keyargs in polyPhase!(polygeno::PolyGeno, keyargs...)
PolyOrigin.polyPhase!
— FunctionpolyPhase!(polygeno::PolyGeno, keyargs...)
performs parental phasing and return phasedgeno::PolyGeno with phasedgeno.parentgeno and polygeno.parentgeno being phased .
Positional arguments
polygeno::PolyGeno
: a struct that stores genotypic data and pedigree info.
Keyword arguments
doseerr::Real=0.01
: genotypic error probability.
seqerror::Real=0.001
: base sequencing error probability for GBS data.
chrpairing_phase::Integer=22
: chromosome pairing in parental phasing, with 22 being only bivalent formations and 44 being bi- and quadri-valent formations.
chrsubset::Union{Nothing,AbstractRange,AbstractVector}=nothing
: subset of chromosomes to be considered, with nothing denoting all chromosomes. Delete chromosome indices that are out of range.
snpsubset::Union{Nothing,AbstractRange,AbstractVector}=nothing
: subset of markers to be considered, with nothing denoting all markers. within a chromosome, marker index starts from 1, and marker indices that are larger than the number of markers within the chromosome are deleted.
isparallel::Bool=true
: if true, multicore computing over chromosomes.
delmarker::Bool=true
: if true, delete markers during parental phasing.
"delsiglevel::Real=0.05": significance level for deleting markers.
maxstuck::Integer=5
: the max number of consecutive iterations that are rejected in a phasing run.
maxiter::Integer=30
: the max number of iterations in a phasing run."
minrun::Integer=3
: if the min number of phasing runs that are at the same local maximimum or have the same parental phases reaches minrun, phasing algorithm will stop before reaching the maxrun.
maxrun::Integer=10
: the max number of phasing runs.
byparent::Union{Nothing,Bool}=nothing
: if true, update parental phases parent by parent; if false, update parental phases one subpopulation by subpopulation. The nothing denotes that it is true if a connected component is a simple F1 cross, and false otherwise.
byneighbor::Union{Nothing,Bool}=nothing
: if ture, udpate the combination of bivalent or multivalents in parents by their neighbors; if false, consider all the possible combinations. The nothing denotes thtat it is true if max ploidy>=6, and false otherwise.
refhapfile::Union{Nothing,AbstractString} = nothing
: reference haplotype file for setting absolute parental phases. It has the same format as the input genofile, except that parental genotypes are phased and offspring genotypes are ignored if they exist.
missingstring::AbstractString="NA"
: string code for missing value.
outstem::Union{Nothing,AbstractString}="outstem"
: stem of output filenames. If nothing, no output files.
logfile::Union{Nothing,AbstractString,IO}= (isnothing(outstem) ? nothing : string(outstem,".log"))
: log file or IO for writing log. If nothing, no log file.
workdir::AbstractString = pwd()
: directory for reading and writing files.
verbose::Bool=true
: if true, print messages on console.
PolyOrigin.polyMapRefine!
— FunctionpolyMapRefine!(polygeno::PolyGeno, keyargs...)
performs marker map refinning for polygeno with phased parent genotypes. Modifies polygeno.markermap into a refined genetic map.
Positional arguments
polygeno::PolyGeno
: a struct that stores genotypic data and pedigree info.
Keyword arguments
doseerr::Real=0.01
: genotypic error probability for offspring phasing.
chrpairing::Integer=44
: chromosome pairing in offspring decoding, with 22 being only bivalent formations and 44 being bivalent and quadrivalent formations.
chrsubset::Union{Nothing,AbstractRange,AbstractVector}=nothing
: subset of chromosomes to be considered, with nothing denoting all chromosomes. Delete chromosome indices that are out of range.
snpsubset::Union{Nothing,AbstractRange,AbstractVector}=nothing
: subset of markers to be considered, with nothing denoting all markers. within a chromosome, marker index starts from 1, and marker indices that are larger than the number of markers within the chromosome are deleted.
isparallel::Bool=true
: if true, multicore computing over chromosomes.
refineorder::Bool=false
: if true, refine marker mordering.
maxwinsize::Integer=50
: max size of sliding windown in map refinning.
inittemperature::Real=4
: initial temperature of simulated annealing in map refinning.
coolingrate::Real=0.5
: cooling rate of annealing temperature in map refinning.
stripdis::Real=20
: a chromosome end in map refinement is removed if it has a distance gap > stripdis (centiMorgan) and it contains less than 5% markers.
maxdoseerr::Real=0.5
: markers in map refinement are removed it they have error rates > maxdoseerr.
skeletonsize::Integer=50
: the number of markers in the skeleton map that is used to re-scale inter-map distances.
logfile::Union{AbstractString,IOStream}=string(outstem,".log")
: output filenames or stream for writing log.
workdir::AbstractString = pwd()
: directory for reading and writing files.
verbose::Bool=true
: if true, print messages on console.
PolyOrigin.polyReconstruct!
— FunctionpolyReconstruct!(polygeno::PolyGeno, keyargs...)
performs ancestral inference from polygeno and return polyancestry::PolyAncestry.
Positional arguments
polygeno::PolyGeno
: a struct that stores genotypic data and pedigree info.
Keyword arguments
doseerr::Real=0.01
: genotypic error probability for offspring phasing.
seqerror::Real=0.001
: base sequencing error probability for GBS data.
chrpairing::Integer=44
: chromosome pairing in offspring decoding, with 22 being only bivalent formations and 44 being bivalent and quadrivalent formations.
chrsubset::Union{Nothing,AbstractRange,AbstractVector}=nothing
: subset of chromosomes to be considered, with nothing denoting all chromosomes. Delete chromosome indices that are out of range.
snpsubset::Union{Nothing,AbstractRange,AbstractVector}=nothing
: subset of markers to be considered, with nothing denoting all markers. within a chromosome, marker index starts from 1, and marker indices that are larger than the number of markers within the chromosome are deleted.
isparallel::Bool=true
: if true, multicore computing over chromosomes.
correctthreshold::AbstractFloat=0.15
: a candidate marker is selected for parental error correction if the fraction of offspring genotypic error >= correctthreshold.
isinfererror::Bool=false
: if true, infer dosage error rate per marker.
isplot::Bool=false
: if true, plot haploprob for all offspring and save in the folder "outstem_plots".
outstem::Union{Nothing,AbstractString}="outstem"
: stem of output filenames. If nothing, no output files.
logfile::Union{Nothing,AbstractString,IO}= (isnothing(outstem) ? nothing : string(outstem,".log"))
: log file or IO for writing log. If nothing, no log file.
workdir::AbstractString = pwd()
: directory for reading and writing files.
verbose::Bool=true
: if true, print messages on console.
PolyOrigin.setAbsPhase!
— FunctionsetAbsPhase!(refhapfile,phasedgeno::PolyGeno,io=nothing,verbose=true)
set absolute parental phases, based on the reference haplotype file.
Positional arguments
refhapfile::AbstractString
: reference haplotype file has the same format as the input genofile, except that parental genotypes are phased and offspring genotypes are ignored if they exist.
phasedgeno::PolyGeno
: polygeno struct with phased parent genotypes .
Keyward arguments
workdir::AbstractString = pwd()
: directory for reading refhapfile
.
io::Union{Nothing,IOStream}
: stream for writing log.
verbose::Bool=true
: true if print messages on console.
setAbsPhase!(refhapfile, polyancestry,io=nothing,verbose=true)
set absolute parental phases and consistent polyancestry.genoprob and polyancestry.haploprob, based on the reference haplotype file.
Positional arguments
refhapfile::AbstractString
: reference haplotype file has the same format as the input genofile, except that parental genotypes are phased and offspring genotypes are ignored if they exist.
polyancestry::PolyAncestry
: polyancestry struct with phased parent genotypes and inference genoprob.
Keyward arguments
workdir::AbstractString = pwd()
: directory for reading trueancestryfile
.
io::Union{Nothing,IOStream}
: stream for writing log.
verbose::Bool=true
: true if print messages on console.
PolyOrigin.savegenodata
— Functionsavegenodata(outfile,polygeno,missingstring="NA",workdir=pwd())
saves genotypic data of the struct polygeno into outfile
Positional arguments
outfile::AbstractString
: filename for saving results.
polygeno::PolyGeno
: a struct returned by readPolyGeno
.
Keyward arguments
missingstring::AbstractString="NA"
: string code for missing value.
workdir::AbstractString = pwd()
: directory for writing outfile.
PolyOrigin.savegenoprob
— Functionsavegenoprob(outfile,polyancestry,missingstring="NA",workdir=pwd())
saves polyancestry.genoprob and polyancestry.parentgeno into outfile.
Positional arguments
outfile::AbstractString
: filename for saving results.
polyancestry::PolyAncestry
: a struct returned by polyOrigin
.
Keyward arguments
missingstring::AbstractString="NA"
: string code for missing value.
workdir::AbstractString = pwd()
: directory for writing outfile.