Pipelines#
- class pyhmmer.plan7.Pipeline#
An HMMER3 accelerated sequence/profile comparison pipeline.
The Plan7 pipeline handles the platform-accelerated comparison of sequences to profile HMMs. It performs either a search (comparing a single query profile to a target sequence database) or a scan (comparing a single query sequence to a target profile database). The two methods are yielding equivalent results: if you have a collection of \(M\) sequences and \(N\) HMMs to compare, doing a search or a scan should give the same raw scores. The E-values will however be different if
Z
anddomZ
where not set manually: \(Z\) will be set to \(M\) for a search, and to \(N\) for a scan.The main reason for which you should choose search or scan is the relative size of the sequences and HMMs databases. In the original HMMER3 code, the memory was managed in a way that you never had to load the entirety of the target sequences in memory. In PyHMMER, the methods accept both reading the target database from a file, or loading it entirely into memory. A scan is always slower than a search because of the overhead introduced when reconfiguring a profile for a new sequence.
- background#
The null background model to use to compute scores.
- Type:
- randomness#
The random number generator being used by the pipeline.
- Type:
Changed in version 0.4.2: Added the
randomness
attribute.- __init__(alphabet, background=None, *, bias_filter=True, null2=True, seed=42, Z=None, domZ=None, F1=0.02, F2=0.001, F3=1e-05, E=10.0, T=None, domE=10.0, domT=None, incE=0.01, incT=None, incdomE=0.01, incdomT=None, bit_cutoffs=None)#
Instantiate and configure a new accelerated comparison pipeline.
- Parameters:
alphabet (
Alphabet
) – The biological alphabet the of the HMMs and sequences that are going to be compared. Used to build the background model.background (
Background
, optional) – The background model to use with the pipeline, orNone
to create and use a default one. The pipeline needs ownership of the background model, so any background model passed there will be copied.
- Keyword Arguments:
bias_filter (
bool
) – Whether or not to enable composition bias filter. Defaults toTrue
.null2 (
bool
) – Whether or not to compute biased composition score corrections. Defaults toTrue
.seed (
int
, optional) – The seed to use with the random number generator. Pass 0 to use a one-time arbitrary seed, orNone
to keep the default seed from HMMER.Z (
int
, optional) – The effective number of comparisons done, for E-value calculation. Leave asNone
to auto-detect by counting the number of sequences queried.domZ (
int
, optional) – The number of significant sequences found, for domain E-value calculation. Leave asNone
to auto-detect by counting the number of sequences reported.F1 (
float
) – The MSV filter threshold.F2 (
float
) – The Viterbi filter threshold.F3 (
float
) – The uncorrected Forward filter threshold.E (
float
) – The per-target E-value threshold for reporting a hit.T (
float
, optional) – The per-target bit score threshold for reporting a hit. If given, takes precedence overE
.domE (
float
) – The per-domain E-value threshold for reporting a domain hit.domT (
float
, optional) – The per-domain bit score threshold for reporting a domain hit. If given, takes precedence overdomE
.incE (
float
) – The per-target E-value threshold for including a hit in the resultingTopHits
.incT (
float
, optional) – The per-target bit score threshold for including a hit in the resultingTopHits
. If given, takes precedence overincE
.incdomE (
float
) – The per-domain E-value threshold for including a domain in the resultingTopHits
.incdomT (
float
, optional) – The per-domain bit score thresholds for including a domain in the resultingTopHits
. If given, takes precedence overincdomE
.bit_cutoffs (
str
, optional) – The model-specific thresholding option to use for reporting hits. WithNone
(the default), use global pipeline options; otherwise pass one of"noise"
,"gathering"
or"trusted"
to use the appropriate cutoffs.
Hint
In order to run the pipeline in slow/max mode, disable the bias filter, and set the three filtering parameters to \(1.0\):
>>> dna = easel.Alphabet.dna() >>> max_pli = Pipeline(dna, bias_filter=False, F1=1.0, F2=1.0, F3=1.0)
Changed in version 0.4.6: Added keywords arguments to set the E-value thresholds.
- arguments()#
Generate an
argv
-like list with pipeline options as CLI flags.Example
>>> alphabet = easel.Alphabet.amino() >>> plan7.Pipeline(alphabet).arguments() [] >>> plan7.Pipeline(alphabet, F1=0.01).arguments() ['--F1', '0.01']
Added in version 0.6.0.
- clear()#
Reset the pipeline configuration to its default state.
- iterate_hmm(query, sequences, builder=None, select_hits=None)#
Run the pipeline to find homologous sequences to a query HMM.
- Parameters:
query (
HMM
) – The sequence object to use to query the sequence database.sequences (
DigitalSequenceBlock
) – The target sequences to query with the HMM.builder (
Builder
, optional) – A HMM builder to use to convert the query and subsequent alignments to aHMM
. IfNone
is given, this method will create one with the default parameters.select_hits (callable, optional) – A function or callable object for manually selecting hits during each iteration. It should take a single
TopHits
argument and change the inclusion of every individual hits by setting theincluded
anddropped
flags of eachHit
manually.
- Returns:
IterativeSearch
– An iterator object yielding the hits, sequence alignment, and HMM for each iteration.- Raises:
AlphabetMismatch – When the alphabet of the current pipeline does not match the alphabet of the given query or database sequences.
Hint
This method corresponds to running
jackhmmer
with thequery
sequence against thesequences
database.Caution
Default values used for
jackhmmer
do not correspond to the default parameters used for creating a pipeline in the other cases. To have truly identical results to thejackhmmer
results in default mode, create thePipeline
object withincE=0.001
andincdomE=0.001
.See also
The
iterate_seq
, which does the same operation with a query sequence instead of a query HMM, and contains more details and examples.Changed in version 0.7.0: Targets must now be inside a
DigitalSequenceBlock
.
- iterate_seq(query, sequences, builder=None, select_hits=None)#
Run the pipeline to find homologous sequences to a query sequence.
This method implements an iterative search over a database to find all sequences homologous to a query sequence. It is very sensitive to the pipeline inclusion thresholds (
incE
andincdomE
).Since this method returns an iterator, the local results of each iteration will be available for inspection before starting the next one. The
select_hits
callback in particular can be used for manually including/excluding hits in each iteration, which is not supported in the originaljackhmmer
, but available on the HMMER web client.- Parameters:
query (
DigitalSequence
) – The sequence object to use to query the sequence database.sequences (
DigitalSequenceBlock
) – The target sequences to query with the query sequence.builder (
Builder
, optional) – A HMM builder to use to convert the query and subsequent alignments to aHMM
. IfNone
is given, this method will create one with the default parameters.select_hits (callable, optional) – A function or callable object for manually selecting hits during each iteration. It should take a single
TopHits
argument and change the inclusion of individual hits with theHit.include
andHit.drop
methods.
- Returns:
IterativeSearch
– An iterator object yielding the hits, sequence alignment, and HMM for each iteration.- Raises:
AlphabetMismatch – When the alphabet of the current pipeline does not match the alphabet of the given query or database sequences.
Hint
This method corresponds to running
jackhmmer
with thequery
sequence against thesequences
database.Caution
Default values used for
jackhmmer
do not correspond to the default parameters used for creating a pipeline in the other cases. To have truly identical results to thejackhmmer
results in default mode, create thePipeline
object withincE=0.001
andincdomE=0.001
.Example
Starting from a pipeline and a query sequence, let’s first obtain the iterator over the successive results:
>>> abc = easel.Alphabet.amino() >>> pli = plan7.Pipeline(abc, incE=1e-3, incdomE=1e-3) >>> iterator = pli.iterate_seq(reductase, proteins)
Once this is ready, we can keep iterating until we converge:
>>> converged = False >>> while not converged: ... _, hits, _, converged, _ = next(iterator) ... print(f"Hits: {len(hits)} Converged: {converged}") Hits: 1 Converged: False Hits: 2 Converged: False Hits: 2 Converged: True
To prevent diverging searches from running infinitely, you could wrap the search in a
for
loop instead, using a number of maximum iterations as the upper boundary:>>> iterator = pli.iterate_seq(reductase, proteins) >>> max_iterations = 10 >>> for n in range(max_iterations): ... iteration = next(iterator) ... if iteration.converged: ... break
Added in version 0.6.0.
Changed in version 0.7.0: Targets must now be inside a
DigitalSequenceBlock
.
- search_msa(query, sequences, builder=None)#
Run the pipeline using a query alignment against a sequence database.
- Parameters:
query (
DigitalMSA
) – The multiple sequence alignment to use to query the sequence database.sequences (
DigitalSequenceBlock
orSequenceFile
) – The target sequences to query with the alignment, either pre-loaded in memory inside apyhmmer.easel.DigitalSequenceBlock
, or to be read iteratively from aSequenceFile
opened in digital mode.builder (
Builder
, optional) – A HMM builder to use to convert the query to aHMM
. IfNone
is given, it will use a default one.
- Returns:
TopHits
– the hits found in the sequence database.- Raises:
AlphabetMismatch – When the alphabet of the current pipeline does not match the alphabet of the given query.
ValueError – When the
Builder
fails to create an HMM from the givenDigitalMSA
query.
Hint
This method corresponds to running
phmmer
with thequery
multiple sequence alignment against thesequences
database.Caution
Internally, this method will create a new HMM from the query MSA using the
Builder.build_msa
method. HMMER requires that every HMM has a name, so theBuilder
will attempt to use the name of the query MSA to name the HMM. Passing an MSA without a name will result in an error.Added in version 0.3.0.
Changed in version 0.7.0: Targets can be inside a
DigitalSequenceBlock
or aSequenceFile
.
- search_seq(query, sequences, builder=None)#
Run the pipeline using a query sequence against a sequence database.
- Parameters:
query (
DigitalSequence
) – The sequence object to use to query the sequence database.sequences (
DigitalSequenceBlock
orSequenceFile
) – The target sequences to query with the query sequence, either pre-loaded in memory inside apyhmmer.easel.DigitalSequenceBlock
, or to be read iteratively from aSequenceFile
opened in digital mode.builder (
Builder
, optional) – A HMM builder to use to convert the query to aHMM
. IfNone
is given, it will use a default one.
- Returns:
TopHits
– the hits found in the sequence database.- Raises:
AlphabetMismatch – When the alphabet of the current pipeline does not match the alphabet of the given query.
ValueError – When the
Builder
fails to create an HMM from the givenDigitalSequence
query.
Hint
This method corresponds to running
phmmer
with thequery
sequence against thesequences
database.Added in version 0.2.0.
Changed in version 0.7.0: Targets can be inside a
DigitalSequenceBlock
or aSequenceFile
.
- T#
The per-target score threshold for reporting a hit.
If set to a non-
None
value, this threshold takes precedence over the per-target E-value threshold (Pipeline.E
).Added in version 0.4.8.
- Z#
The number of effective targets searched.
It is used to compute the independent e-value for each domain, and for an entire hit. If
None
, the parameter number will be set automatically after all the comparisons have been done. Otherwise, it can be set to an arbitrary number.
- bias_filter#
Whether or not to enable the biased comp HMM filter.
Added in version 0.4.1.
- Type:
- bit_cutoffs#
The model-specific thresholding option, if any.
Added in version 0.4.6.
- domT#
The per-domain score threshold for reporting a hit.
If set to a non-
None
value, this threshold takes precedence over the per-domain E-value threshold (Pipeline.domE
).Added in version 0.4.8.
- domZ#
The number of significant targets.
It is used to compute the conditional e-value for each domain. If
None
, the parameter number will be set automatically after all the comparisons have been done, and all hits have been thresholded. Otherwise, it can be set to an arbitrary number.
- incT#
The per-target score threshold for including a hit.
If set to a non-
None
value, this threshold takes precedence over the per-target E-value inclusion threshold (Pipeline.incE
).Added in version 0.4.8.
- incdomT#
The per-domain score threshold for including a hit.
If set to a non-
None
value, this threshold takes precedence over the per-domain E-value inclusion threshold (Pipeline.incdomE
).Added in version 0.4.8.
- scan_seq#
Run the pipeline using a query sequence against a profile database.
- Parameters:
query (
DigitalSequence
) – The sequence object to use to query the profile database.targets (
OptimizedProfileBlock
orHMMPressedFile
) – The optimized profiles to query, either pre-loaded in memory as anOptimizedProfileBlock
, or to be read iteratively from a pressed HMM file.
- Returns:
TopHits
– the hits found in the profile database.- Raises:
AlphabetMismatch – When the alphabet of the current pipeline does not match the alphabet of the given query or profiles.
Caution
In the current version, this method is not optimized to use the pressed database, even if it exists. This will cause the MSV and SSV filters to be rebuilt at each iteration, which could be slow. Consider at least pre-fetching the HMM database if calling this method several times in a row.
Hint
This method corresponds to running
hmmscan
with thequery
sequence against thehmms
database.Added in version 0.4.0.
Changed in version 0.7.0: Require optimized profiles to be inside an
OptimizedProfileBlock
.
- search_hmm#
Run the pipeline using a query HMM against a sequence database.
- Parameters:
query (
HMM
,Profile
orOptimizedProfile
) – The object to use to query the sequence database.sequences (
DigitalSequenceBlock
orSequenceFile
) – The target sequences to query with the HMM, either pre-loaded in memory inside aDigitalSequenceBlock
, or to be read iteratively from aSequenceFile
opened in digital mode.
- Returns:
TopHits
– the hits found in the sequence database.- Raises:
ValueError – When the pipeline is configured to use model-specific reporting thresholds but the
HMM
query doesn’t have the right cutoffs available.AlphabetMismatch – When the alphabet of the current pipeline does not match the alphabet of the given HMM.
Hint
This method corresponds to running
hmmsearch
with thequery
HMM against thesequences
database.Added in version 0.2.0.
Changed in version 0.4.9: Query can now be a
Profile
or anOptimizedProfile
.Changed in version 0.7.0: Targets can be inside a
DigitalSequenceBlock
or aSequenceFile
.
- class pyhmmer.plan7.LongTargetsPipeline(Pipeline)#
An HMMER3 pipeline tuned for long targets.
The default HMMER3 pipeline is configured not to accept target sequences longer than 100,000 residues. Although there is no strong limitation for this threshold, comparing a sequence of \(L\) residues to a profile with \(M\) nodes requires the allocation of a \(L imes M\) dynamic programming matrix.
For sequences too long, it’s actually more efficient memory-wise to use a sliding window to match the profile to the sequence. The usual comparison pipeline is then used to perform the comparison on each window, and results are merged once the entire sequence is done being processed. The context size \(C\) is large enough to accommodate for the entire profile, so that there is no risk of missing a hit in the overlaps between windows. The window size \(W\) can be changed with the
block_length
argument when instantiating a newLongTargetsPipeline
object.Added in version 0.4.9.
Added in version 0.10.8: The
window_length
andwindow_beta
keyword arguments.- __init__(alphabet, background=None, *, F1=0.02, F2=0.003, F3=3e-05, strand=None, B1=100, B2=240, B3=1000, block_length=262144, window_length=None, window_beta=None, **kwargs)#
Instantiate and configure a new long targets pipeline.
- Parameters:
alphabet (
Alphabet
) – The biological alphabet the of the HMMs and sequences that are going to be compared. Used to build the background model. A nucleotide alphabet is expected.background (
Background
, optional) – The background model to use with the pipeline, orNone
to create and use a default one. The pipeline needs ownership of the background model, so any background model passed there will be copied.
- Keyword Arguments:
strand (
str
, optional) – The strand to use when processing nucleotide sequences. Use"watson"
to use only the coding strand,"crick"
to use only the reverse strand, or leave asNone
to process both strands.B1 (
int
) – The window length to use for the biased-composition modifier of the MSV filter.B2 (
int
) – The window length to use for the biased-composition modifier of the Viterbi filter.B3 (
int
) – The window length to use for the biased-composition modifier of the Forward filter.block_length (
int
) – The number of residues to use as the window size \(W\) when reading blocks from the long target sequences.window_length (
int
) – The window length to use to compute E-values.window_beta (
float
) – The tail mass at which window length is determined.**kwargs – Any additional parameter will be passed to the
Pipeline
constructor.
- arguments()#
Generate an
argv
-like list with pipeline options as CLI flags.Example
>>> alphabet = easel.Alphabet.dna() >>> plan7.LongTargetsPipeline(alphabet).arguments() [] >>> plan7.LongTargetsPipeline(alphabet, B1=200).arguments() ['--B1', '200']
Added in version 0.6.0.
- search_msa(query, sequences, builder=None)#
Run the pipeline using a query alignment against a sequence database.
- Parameters:
query (
DigitalMSA
) – The multiple sequence alignment to use to query the sequence database.sequences (
DigitalSequenceBlock
) – The target sequences to query with the query alignment.builder (
Builder
, optional) – A HMM builder to use to convert the query to aHMM
. IfNone
is given, it will use a default one.
- Returns:
TopHits
– the hits found in the sequence database.- Raises:
AlphabetMismatch – When the alphabet of the current pipeline does not match the alphabet of the given query.
Hint
This method corresponds to running
nhmmer
with thequery
multiple sequence alignment against thesequences
database.Changed in version 0.7.0: Targets must now be inside a
DigitalSequenceBlock
.
- search_seq(query, sequences, builder=None)#
Run the pipeline using a query sequence against a sequence database.
- Parameters:
query (
DigitalSequence
) – The sequence object to use to query the sequence database.sequences (
DigitalSequenceBlock
) – The target sequences to query with the query sequence.builder (
Builder
, optional) – A HMM builder to use to convert the query to aHMM
. IfNone
is given, it will use a default one.
- Returns:
TopHits
– the hits found in the sequence database.- Raises:
AlphabetMismatch – When the alphabet of the current pipeline does not match the alphabet of the given query.
Hint
This method corresponds to running
nhmmer
with thequery
sequence against thesequences
database.Changed in version 0.7.0: Targets must now be inside a
DigitalSequenceBlock
.
- scan_seq#
Run the pipeline using a query sequence against a profile database.
This is currently unsupported for
LongTargetsPipeline
.
- search_hmm#
Run the pipeline using a query HMM against a sequence database.
- Parameters:
query (
HMM
,Profile
orOptimizedProfile
) – The object to use to query the sequence database.sequences (
DigitalSequenceBlock
) – The target sequences to query with the HMM.
- Returns:
TopHits
– the hits found in the sequence database.- Raises:
ValueError – When the pipeline is configured to use model-specific reporting thresholds but the
HMM
query doesn’t have the right cutoffs available.AlphabetMismatch – When the alphabet of the current pipeline does not match the alphabet of the given HMM.
Hint
This method corresponds to running
nhmmer
with thequery
HMM against thesequences
database.Changed in version 0.7.0: Targets must now be inside a
DigitalSequenceBlock
.
- class pyhmmer.plan7.Builder#
A factory for constructing new HMMs from raw sequences.
- alphabet#
The alphabet the builder is configured to use to convert sequences to HMMs.
- Type:
- randomness#
The random number generator being used by the builder.
- Type:
- score_matrix#
The name of the substitution matrix used to build HMMs for single sequences.
- Type:
Added in version 0.2.0.
Changed in version 0.4.2: Added the
randomness
attribute.- __init__(alphabet, *, architecture='fast', weighting='pb', effective_number='entropy', prior_scheme='alphabet', symfrac=0.5, fragthres=0.5, wid=0.62, esigma=45.0, eid=0.62, EmL=200, EmN=200, EvL=200, EvN=200, EfL=100, EfN=200, Eft=0.04, seed=42, ere=None, popen=None, pextend=None, score_matrix=None, window_length=None, window_beta=None)#
Create a new sequence builder with the given configuration.
- Parameters:
alphabet (
Alphabet
) – The alphabet the builder expects the sequences to be in.- Keyword Arguments:
architecture (
str
) – The algorithm to use to determine the model architecture, either"fast"
(the default), or"hand"
.weighting (
str
) – The algorithm to use for relative sequence weighting, either"pb"
(the default),"gsc"
,"blosum"
,"none"
, or"given"
.effective_number (
str
,int
, orfloat
) – The algorithm to use to determine the effective sequence number, either"entropy"
(the default),"exp"
,"clust"
,"none"
. A number can also be given directly to set the effective sequence number manually.prior_scheme (
str
) – The choice of mixture Dirichlet prior when parameterizing from counts, either"laplace"
or"alphabet"
(the default).symfrac (
float
) – The residue occurrence threshold for fast architecture determination.fragthresh (
float
) – A threshold such that a sequence is called a fragment when \(L \le fragthresh imes alen\).wid (
double
) – The percent identity threshold for BLOSUM relative weighting.esigma (
float
) – The minimum total relative entropy parameter for effective number entropy weights.eid (
float
) – The percent identity threshold for effective number clustering.EmL (
int
) – The length of sequences generated for MSV fitting.EmN (
int
) – The number of sequences generated for MSV fitting.EvL (
int
) – The lenght of sequences generated for Viterbi fitting.EvN (
int
) – The number of sequences generated for Viterbi fitting.EfL (
int
) – The lenght of sequences generated for Forward fitting.EfN (
int
) – The number of sequences generated for Forward fitting.Eft (
float
) – The tail mass used for Forward fitting.seed (
int
) – The seed to use to initialize the internal random number generator. If0
is given, an arbitrary seed will be chosen based on the current time.ere (
double
, optional) – The relative entropy target for effective number weighting, orNone
.popen (
float
, optional) – The gap open probability to use when building HMMs from single sequences. The default value depends on the alphabet: 0.02 for proteins, 0.03125 for nucleotides.pextend (
float
, optional) – The gap extend probability to use when building HMMs from single sequences. Default depends on the alphabet: 0.4 for proteins, 0.75 for nucleotides.score_matrix (
str
, optional) – The name of the score matrix to use when building HMMs from single sequences. The only allowed value for nucleotide alphabets is DNA1. For proteins, PAM30, PAM70, PAM120, PAM240, BLOSUM45, BLOSUM50, BLOSUM62 (the default), BLOSUM80 or BLOSUM90 can be used.window_length (
float
, optional) – The window length for nucleotide sequences, essentially the max expected hit length. If given, takes precedence overwindow_beta
.window_beta (
float
, optional) – The tail mass at which window length is determined for nucleotide sequences.
- build(sequence, background)#
Build a new HMM from
sequence
using the builder configuration.- Parameters:
sequence (
DigitalSequence
) – A single biological sequence in digital mode to build a HMM with.background (
Background
) – The background model to use to create the HMM.
- Returns:
(
HMM
,Profile
,OptimizedProfile
) – A tuple containing the new HMM as well as profiles to be used directly in aPipeline
.- Raises:
AlphabetMismatch – When either
sequence
orbackground
have the wrong alphabet for this builder.ValueError – When the HMM cannot be created successfully from the input; the error message contains more details.
Hint
The score matrix and the gap probabilities used here can be set when initializing the builder, or changed by setting a new value to the right attribute:
>>> alphabet = easel.Alphabet.amino() >>> background = plan7.Background(alphabet) >>> builder = plan7.Builder(alphabet, popen=0.05) >>> builder.score_matrix = "BLOSUM90" >>> hmm, _, _ = builder.build(proteins[0], background)
Changed in version 0.4.6: Sets the
HMM.creation_time
attribute with the current time.Changed in version 0.4.10: The score system is now cached between calls to
Builder.build
.
- build_msa(msa, background)#
Build a new HMM from
msa
using the builder configuration.- Parameters:
msa (
DigitalMSA
) – A multiple sequence alignment in digital mode to build a HMM with.background (
Background
) – The background model to use to create the HMM.
- Returns:
(
HMM
,Profile
,OptimizedProfile
) – A tuple containing the new HMM as well as profiles to be used directly in aPipeline
.- Raises:
AlphabetMismatch – When either
msa
orbackground
have the wrong alphabet for this builder.ValueError – When the HMM cannot be created successfully from the input; the error message contains more details.
Caution
HMMER requires that every HMM has a name, so the
Builder
will attempt to use the name of the sequence to name the HMM. Passing an MSA without a name will result in an error:>>> alphabet = easel.Alphabet.amino() >>> msa = easel.TextMSA(sequences=[ ... easel.TextSequence(name=b"ustiA", sequence="YAIG"), ... easel.TextSequence(name=b"ustiB", sequence="YVIG") ... ]) >>> builder = plan7.Builder(alphabet) >>> background = plan7.Background(alphabet) >>> builder.build_msa(msa.digitize(alphabet), background) Traceback (most recent call last): ... ValueError: Could not build HMM: Unable to name the HMM.
Added in version 0.3.0.
Changed in version 0.4.6: Sets the
HMM.creation_time
attribute with the current time.
- seed#
The seed given at builder initialization.
Setting this attribute to a different value will cause the internal random number generator to be reseeded immediately.
Changed in version 0.4.2: Avoid shadowing initial null seeds given on builder initialization.
- Type:
- class pyhmmer.plan7.Background#
The null background model of HMMER.
Changed in version 0.4.0: Added the
residue_frequencies
attribute.- __init__(alphabet, uniform=False)#
Create a new background model for the given
alphabet
.
- copy()#
Create a copy of the null model with the same parameters.