Dependencies | Installation and execution | Usage | Notes | References | Citations
fq2dna (FASTQ files to de novo assembly) is a command line tool written in Bash to ease the de novo assembly of archaea, bacteria or virus genomes from raw high-throughput sequencing (HTS) paired-end (PE) reads.
Every data pre- and post-processing step is managed by fq2dna (e.g. HTS read filtering and enhancing, well-tuned de novo assemblies, scaffold sequence accuracy assessment). The main purpose of fq2dna is to efficiently use different methods, programs and tools to quickly infer accurate genome assemblies (e.g. from 5 to 20 minutes to deal with a bacteria HTS sample using 12 threads). This mini-workflow can therefore be very useful to deal with large batches of whole-genome shotgun sequencing data.
fq2dna runs on UNIX, Linux and most OS X operating systems.
program | package | version | sources |
---|---|---|---|
gawk | - | > 4.0.0 | ftp.gnu.org/gnu/gawk |
bwa-mem2 | - | ≥ 2.2.1 | gitlab.pasteur.fr/GIPhy/contig_info |
contig_info | - | > 2.0 | gitlab.pasteur.fr/GIPhy/contig_info |
FASTA2AGP | - | ≥ 2.0 | gitlab.pasteur.fr/GIPhy/FASTA2AGP |
fqCleanER | - | ≥ 23.12 | gitlab.pasteur.fr/GIPhy/fqCleanER |
fqstats | fqtools | ≥ 1.2 | ftp.pasteur.fr/pub/gensoft/projects/fqtools |
ntCard | - | > 1.2 | github.com/bcgsc/ntCard |
Platon | - | > 1.5 | github.com/oschwengers/platon |
Prokka | - | ≥ 1.14.5 | github.com/tseemann/prokka |
samtools | - | ≥ 1.18 | github.com/samtools/samtools sourceforge.net/projects/samtools |
SAM2MAP | SAM2MSA | ≥ 0.4.3.1 | gitlab.pasteur.fr/GIPhy/SAM2MSA |
SPAdes | - | ≥ 3.15.5 | github.com/ablab/spades |
A. Clone this repository with the following command line:
git clone https://gitlab.pasteur.fr/GIPhy/fq2dna.git
B. Go to the created directory and give the execute
permission to the file fq2dna.sh
:
cd fq2dna/
chmod +x fq2dna.sh
C. Check the dependencies (and their version) using the following command line:
./fq2dna.sh -d
$PATH
variable (or if one compiled binary has a different
default name), it should be manually specified. To specify the location
of a specific binary, edit the file fq2dna.sh
and indicate
the local path to the corresponding binary(ies) within the code block
REQUIREMENTS
(approximately lines 80-200). For each
required program, the table below reports the corresponding variable
assignment instruction to edit (if needed) within the code block
REQUIREMENTS
program | variable assignment | program | variable assignment | |
---|---|---|---|---|
bwa-mem2 | BWAMEM2_BIN=bwa-mem2; |
ntcard | NTCARD_BIN=ntcard; |
|
contig_info | CONTIG_INFO_BIN=contig_info; |
Platon | PLATON_BIN=platon; |
|
FASTA2AGP | FASTA2AGP_BIN=FASTA2AGP; |
Prokka | PROKKA_BIN=prokka; |
|
fqCleanER | FQCLEANER_BIN=fqCleanER; |
samtools | SAMTOOLS_BIN=samtools; |
|
fqstats | FQSTATS_BIN=fqstats; |
SAM2MAP | SAM2MAP_BIN=SAM2MAP; |
|
gawk | GAWK_BIN=gawk; |
SPAdes | SPADES_BIN=spades.py; |
Note that depending on the installation of some required programs,
the corresponding variable can be assigned with complex commands. For
example, as FASTA2AGP is a Java tool that can be run using a
Java virtual machine, the executable jar file FASTA2AGP.jar
can be used by fq2dna by editing the corresponding variable
assignment instruction as follows:
FASTA2AGP_BIN="java -jar FASTA2AGP.jar"
.
E. Execute fq2dna with the following command line model:
./fq2dna.sh [options]
Run fq2dna without option to read the following documentation:
USAGE: fq2dna.sh [options]
Processing and assembling high-throughput sequencing (HTS) paired-end (PE) reads:
+ processing HTS reads using different steps: deduplicating [D], trimming/clipping [T], error
correction [E], contaminant removal [C], merging [M], and/or digital normalization [N]
+ de novo assembly [dna] of the whole genome from processed HTS reads, followed by polishing
and annotation steps (depending on the specified strategy)
OPTIONS:
-1 <infile> fwd (R1) FASTQ input file name from PE library 1 || input files can be compressed
-2 <infile> rev (R2) FASTQ input file name from PE library 1 || using either gzip (file
-3 <infile> fwd (R1) FASTQ input file name from PE library 2 || extension .gz), bzip2 (.bz or
-4 <infile> rev (R2) FASTQ input file name from PE library 2 || .bz2), or DSRC (.dsrc or
-5 <infile> fwd (R1) FASTQ input file name from PE library 3 || .dsrc2), or uncompressed (.fq
-6 <infile> rev (R2) FASTQ input file name from PE library 3 || or .fastq)
-o <outdir> path and name of the output directory (mandatory option)
-b <string> base name for output files (mandatory option)
-s <char> to set a predefined strategy for processing HTS reads among the following ones:
A Archaea: DT(C)E + [N|MN] + dna + polishing + annotation
B Bacteria: DT(C)E + [N|MN] + dna + polishing + annotation
E Eukaryote: DT(C) + [N|MN] + dna
P Prokaryote: DT(C)E + [N|MN] + dna + polishing
S Standard: DT(C) + [N|MN] + dna + polishing
V Virus: DT(C)E + [N|MN] + dna + polishing + annotation
(default: S)
-L <int> minimum required length for a contig (default: 300)
-T <"G S I"> Genus (G), Species (S) and Isolate (I) names to be used during annotation step;
should be set between quotation marks and separated by a blank space; only with
options -s A, -s B or -s V (default: "Genus sp. STRAIN")
-q <int> quality score threshold; all bases with Phred score below this threshold are
considered as non-confident during step [T] (default: 20)
-l <int> minimum required length for a read (default: half the average read length)
-p <int> maximum allowed percentage of non-confident bases (as ruled by option -q) per
read (default: 50)
-c <int> minimum allowed coverage depth during step [N] (default: 3)
-C <int> maximum allowed coverage depth during step [N] (default: 60)
-a <infile> to set a file containing every alien oligonucleotide sequence (one per line) to
be clipped during step [T] (see below)
-a <string> one or several key words (separated with commas), each corresponding to a set of
alien oligonucleotide sequences to be clipped during step [T]:
POLY nucleotide homopolymers
NEXTERA Illumina Nextera index Kits
IUDI Illumina Unique Dual index Kits
AMPLISEQ AmpliSeq for Illumina Panels
TRUSIGHT_PANCANCER Illumina TruSight RNA Pan-Cancer Kits
TRUSEQ_UD Illumina TruSeq Unique Dual index Kits
TRUSEQ_CD Illumina TruSeq Combinatorial Dual index Kits
TRUSEQ_SINGLE Illumina TruSeq Single index Kits
TRUSEQ_SMALLRNA Illumina TruSeq Small RNA Kits
Note that these sets of alien sequences are not exhaustive and will never replace
the exact oligos used for library preparation (default: "POLY")
-a AUTO to (try to) infer 3' alien oligonucleotide sequence(s) for step [T]; inferred
oligo(s) are completed with those from "POLYS" (see above)
-A <infile> to set sequence or k-mer model file(s) to carry out contaminant read removal
during step [C]; several comma-separated file names can be specified; allowed
file extensions: .fa, .fasta, .fna, .kmr or .kmz
-t <int> number of threads (default: 12)
-w <dir> path to tmp directory (default: $TMPDIR, otherwise /tmp)
-x to not remove (after completing) the tmp directory inside the one set with option
-w (default: not set)
-d checks dependencies and exit
-h prints this help and exit
EXAMPLES:
fq2dna.sh -1 r1.fastq -2 r2.fastq -o out -b ek12 -s P -T "Escherichia coli K12" -a AUTO
fq2dna.sh -1 rA.1.fq -2 rA.2.fq -3 rB.1.fq.gz -4 rB.2.fq.gz -o out -b name -a NEXTERA
fq2dna is able to consider up to three paired-ends
libraries (PE; options -1
to -6
). Input files
should be in FASTQ format and can be compressed using gzip, bzip2 or DSRC (Roguski
and Deorowicz 2014).
In brief, fq2dna first performs initial basic HTS read
preprocessing (i.e. step I) using fqCleanER,
i.e. deduplication (D), trimming/clipping (T; as ruled by options
-a
, -q
, -l
and -p
)
and error correction (E); when specified (option -A
), a
contaminating HTS read removal step (C ) can also be performed.
Next, fq2dna creates two distinct datasets from the
preprocessed HTS reads: a first one (i.e. step N)
obtained using a digital normalization procedure (N; as
ruled by options -C
and -c
), and a second one
(i.e. step M) by merging the PE HTS reads with short
insert size (M) followed by a digital normalization
procedure. Each of these two HTS datasets is used to infer a de
novo genome assembly (dnaN and dnaM, respectively) using SPAdes
(Bankevich et al. 2012).
Among the two genome assemblies, the less
fragmented one is retained, and next polished (i.e. correcting putative
local assembly errors such as mismatches and short indels) using samtools on the aligned step
I HTS reads. Polished scaffolds are then used together
with their generating (step N or M)
HTS reads to infer a genome coverage profile (e.g. Lindner et al. 2013);
based on this coverage profile, sufficiently long and well-covered
scaffold sequences are finally selected.
Depending on the specified
strategy (option -s
), selected scaffold sequences are
classified as chromosome/plasmid/undetermined (i.e. CPU) using Platon
(Schwengers et al. 2020) and/or annotated using Prokka (Seemann
2014).
Output files are all defined by the same specified prefix
(mandatory option -b
) and written into a specified output
directory (mandatory option -o
). Each output file content
is determined by its file extension:
output file | file content |
---|---|
<prefix>.stepI.log | fqCleanER log file of the step I |
<prefix>.stepM.log | fqCleanER log file of the step M |
<prefix>.stepN.log | fqCleanER log file of the step N |
<prefix>.all.fasta | the less fragmented SPAdes assembly (FASTA format) |
<prefix>.cov.info.txt | coverage profile summary generated using SAM2MAP |
<prefix>.scf.fasta | selected and polished scaffold sequences (FASTA format) |
<prefix>.scf.info.txt | residue content of the selected scaffold sequences (tab-delimited) |
<prefix>.scf.amb.txt | ambiguously assembled bases (tab-delimited) |
<prefix>.agp.fasta | contigs derived from the selected and polished scaffold sequences (FASTA format) |
<prefix>.agp | scaffolding information associated to the contigs (AGP format) |
<prefix>.dna.info.txt | descriptive statistics of each FASTA file content (tab-delimited) |
<prefix>.isd.txt | descriptive statistics of the insert size distribution |
<prefix>.gbk | assembled genome annotation (GenBank flat file format) |
<prefix>.gbk.info.txt | annotation statistics generated by Prokka |
Temporary files are written into a dedicated directory created
into the $TMPDIR
directory (when defined, otherwise
/tmp/
). When possible, it is highly recommended to set a
temp directory with large capacity (option -w
).
Default options lead to a de novo assembly inferred from
a subset of well-suited HTS read (digital normalization) corresponding
to 60x coverage depth (option -C
), which is a good tradeoff
to observe accurate results with fast running times. It is therefore
expected that (i) the genome coverage profile (output file
<prefix>._cov.info.txt_) is a unimodal distribution with
both mean and mode close to 60x, and (ii) most of the selected scaffold
sequences (<prefix>._scf.fasta_) are supported by 60x
coverage depth on average (covr values in FASTA headers).
Strong deviation from the specified coverage depth (option
-C
) can be indicative of a sequencing problem, e.g. low
coverage depth (small mode value), locally insufficient coverage depth
(very negative skewness), sample contamination (bimodal distribution).
To ease the identification of contaminated samples, an original
mono-modality (momo) index is available in the output file
<prefix>._cov.info.txt_, e.g. non-monomodal distributions
generally lead to momo < 0.99.
Of important note is that the genome coverage profile is also used to fit a theoretical distribution using SAM2MAP (i.e. Negative Binomial or Generalized Poisson, depending on the variance). Such a theoretical distribution enables to determine a coverage depth cutoff, under which any observed coverage depth is considered as significantly low (p-value ≤ 0.000001). Every assembled base with coverage depth lower than the estimated coverage cutoff (specified in FASTA headers) is written in lowercase (into both files <prefix>._all.fasta_ and <prefix>._scf.fasta_) and should be considered with caution. If the selected scaffold sequences (i.e. with average coverage depth greater than the coverage cutoff; output file <prefix>._scf.fasta_) do not seems to represent the whole genome, try to consider the entire set of assembled scaffolds (output file <prefix>._all.fasta_).
An ambiguously assembled base is defined as a (uppercase) character state N that is sufficiently covered by its generating HTS reads (e.g. non-monoallelic position). Every ambiguously assembled base is listed in the tab-delimited output file <prefix>._scf.amb.txt_ (one per line) together with an estimate of the sequenced content corresponding to this assembled position (i.e. proportions of A, C, G, T and gaps, respectively, estimated from the alignment of the generating HTS reads).
For more details about the different assembly and scaffold sequence statistics (output files <prefix>._dna.info.txt_ and <prefix>._scf.info.txt_, respectively), see the documentation of contig_info.
For more details about the chromosome/plasmid/undetermined
classification of each scaffold sequence (i.e. column CPU in
<prefix>._scf.info.txt_ when using option
-s B
), see the documentation of Platon.
For more details about the annotation procedure when using
strategies A, B or V (option -s
; output files
<prefix>._gbk_ and
<prefix>._gbk.info.txt_), see the documentation of Prokka.
In order to illustrate the usefulness of fq2dna and to better describe its output files, the following use case example describes its usage for (re)assembling the draft genome of Listeria monocytogenes 2HF33 (Duru et al. 2020).
All output files are available in the directory example/ (the four large sequence files are compressed using gzip), as well as the version of every used tool and program (see program.versions.txt).
Downloading input files
Paired-end sequencing of this genome was performed using Illumina Miseq, and the resulting pair of (compressed) FASTQ files (112 Mb and 128 Mb, respectively) can be downloaded using the following command lines:
wget ftp://ftp.sra.ebi.ac.uk/vol1/fastq/ERR403/006/ERR4032786/ERR4032786_1.fastq.gz
wget ftp://ftp.sra.ebi.ac.uk/vol1/fastq/ERR403/006/ERR4032786/ERR4032786_2.fastq.gz
As phiX genome was used as spike-in during Illumina sequencing (Duru et al. 2020), this putative contaminating sequence (5.4 kb) can be downloded using the following command line:
wget -O phiX.fasta "https://eutils.ncbi.nlm.nih.gov/entrez/eutils/efetch.fcgi?db=nuccore&rettype=fasta&id=NC_001422.1"
Running fq2dna
Below is a typical command line to run fq2dna on PE FASTQ files to deal with bacteria HTS data:
./fq2dna.sh -1 ERR4032786_1.fastq.gz -2 ERR4032786_2.fastq.gz -a NEXTERA -A phiX.fasta \\
-s B -T "Listeria monocytogenes 2HF33" -o 2HF33 -b Lm.2HF33 -t 12 -w /tmp
Of note, option -a NEXTERA
was set to clip the Nextera
XT technical oligonucleotides used during library preparation (see Duru
et al. 2020), the putative contaminating sequence file
phiX.fasta was set using option -A
, and the
strategy ‘Bacteria’ was set (option -s B
).
Different log outputs can be observed depending on the OS and/or installed program versions (see Dependencies). Below is the one observed when using fq2dna together with program and tool versions listed into example/program.versions.txt:
# fq2dna v25.03
# Copyright (C) 2016-2025 Institut Pasteur
+ https://gitlab.pasteur.fr/GIPhy/fq2dna
> System: x86_64-redhat-linux-gnu
> Bash: 4.4.20(1)-release
# input file(s):
> PE lib 1
+ FQ11: [111.9 Mb] ERR4032786_1.fastq.gz
+ FQ12: [127.8 Mb] ERR4032786_2.fastq.gz
# output directory
+ OUTDIR=2HF33
# tmp directory
+ TMP_DIR=/tmp/fq2dna.Lm.2HF33.Iqb7XNY1Q
# STRATEGY B: Bacteria
[00:00] Deduplicating, Trimming, Clipping, Decontaminating and Correcting reads ... [ok]
[03:00] Merging and/or Normalizing reads ... [ok]
[04:07] Approximating genome size ... [ok]
> 3112890 bps (N=3112026 M=3113754)
[04:12] Assembling genome with/without PE read merging (dnaM/dnaN) ... [ok]
> [dnaN] covr=67 arl=278 k=21,33,55,77,99,121
> [dnaM] covr=60 arl=279 k=21,33,55,77,99,121
[06:31] Comparing genome assemblies ... [ok]
> [dnaN] Nseq=32 Nres=3099402 NG50=361897 auGN=462403
> [dnaM] Nseq=33 Nres=3099184 NG50=452690 auGN=346675
> selecting dnaN
> [dna] Nseq=32 Nres=3099402 N50=361897 auN=464415
[06:33] Aligning PE reads against scaffolds .... [ok]
[06:51] Polishing scaffolds ...
> [dna] Nseq=32 Nres=3099547 NG50=361909 auGN=462419
> [dna] Nseq=31 Nres=3099603 NG50=361915 auGN=462429
> [dna] Nseq=31 Nres=3099603 N50=361915 auN=464406
[07:38] Processing scaffolds ... [ok]
> expect. cov. 60
> avg. cov. 65.3006
> cov. mode 65
> momo index 0.999494
> accuracy 0.994893
> min. cov. cutoff 40
> ambiguous pos. 19
> [dna] Nseq=26 Nres=3097515 N50=361915 auN=464718
[08:10] Annotating scaffold sequences ... [ok]
> taxon: Listeria monocytogenes 2HF33
> locus tag: LISMON2HF33
# output files:
+ de novo assembly (all scaffolds): 2HF33/Lm.2HF33.all.fasta
+ coverage profile summary: 2HF33/Lm.2HF33.cov.info.txt
+ de novo assembly (selected scaffolds): 2HF33/Lm.2HF33.scf.fasta
+ sequence stats (selected scaffolds): 2HF33/Lm.2HF33.scf.info.txt
+ ambiguous positions (selected scaffolds): 2HF33/Lm.2HF33.scf.amb.txt
+ de novo assembly (selected contigs): 2HF33/Lm.2HF33.agp.fasta
+ scaffolding info: 2HF33/Lm.2HF33.agp
+ descriptive statistics: 2HF33/Lm.2HF33.dna.info.txt
+ insert size statistics: 2HF33/Lm.2HF33.isd.txt
+ annotation (selected scaffolds): 2HF33/Lm.2HF33.scf.gbk
+ annotation info: 2HF33/Lm.2HF33.scf.gbk.info.txt
[12:17] exit
This log output shows that the complete analysis was performed on 12 threads in less than 15 minutes (HTS read processing in < 5 minutes, de novo assemblies in ~2 minutes, and scaffold sequence post-processing in ~4 minutes). One can also observe that the selected de novo assembly was dnaN (i.e. greater NG50 and auGN metrics), therefore showing that the PE HTS read merging step (M) does not always enable to obtain better genome assemblies.
Output files
The genome coverage profile summary file
(Lm.2HF33.cov.info.txt; see excerpt below) shows an average
coverage depth of ~65x, corresponding to a symmetrical coverage depth
distribution with mode 65×, close to the expected value (default option
-C 60
). The momo index is 0.999245 (above the
mono-modality cutoff, e.g. 0.99), therefore suggesting that no
contamination occurred in the HTS reads selected for assembly. The
overall accuracy is assessed by the following formula:
momo × min( ocov/ecov , mode/ocov ), where ocov and ecov are the
observed and expected coverage depths, respectively. As ocov = 65.3006,
ecov = 60 (option -C
) and mode = 65, this leads to
accuracy = 0.994893, therefore suggesting that the sequencing data are
of good quality (i.e. accuracy ≥ 0.99).
= observed coverage distribution: no.pos=3100533
avg=65.3006 mode=65 momo=0.999494 accuracy=0.994893
# Poisson(l) coverage tail distribution: l=0.17783505 w=0.00012514
* GP(l',r) coverage distribution: l'=96.67031669 r=-0.47633017 1-w=0.99987486
min.cov.cutoff=40 (p-value<=0.000001)
These different statistics therefore assess that the assembled sequences are accurate and well-supported by sufficient data. Finally, a generalized Poisson (GP) distribution (that fits the observed genome coverage profile) enables to determine a minimum coverage cutoff (i.e. 40×) under which any assembled base is considered as putatively not trustworthy (lowercase characters in the scaffold sequence file Lm.2HF33.scf.fasta).
The tab-delimited file Lm.2HF33.dna.info.txt shows that the de novo assembly procedure led to quite few scaffold and contig sequences:
#File Nseq Nres A C G T N %A %C %G %T %N %AT %GC Min Q25 Med Q75 Max Avg auN N50 N75 N90 L50 L75 L90
Lm.2HF33.all.fasta 34 3100965 943191 604057 567949 985468 300 30.41% 19.47% 18.31% 31.77% 0.00% 62.21% 37.79% 237 440 1664 99066 1091349 91204.85 559175 532745 125907 97509 2 6 10
Lm.2HF33.scf.fasta 24 3097730 942310 603413 567209 984585 213 30.41% 19.47% 18.31% 31.78% 0.00% 62.21% 37.79% 407 1251 46823 108604 1091380 129072.08 559792 532745 125913 97524 2 6 10
Lm.2HF33.agp.fasta 26 3097550 942310 603412 567209 984584 35 30.42% 19.48% 18.31% 31.78% 0.00% 62.21% 37.79% 407 1664 46823 125913 935180 119136.53 464718 361915 125913 90129 3 7 11
Of note, 10 small (min. length = 300 bps; option -L
)
and/or insufficiently covered sequences (cutoff = 40×; see
Lm.2HF33.cov.info.txt) were not selected for the final scaffold
sequence set (i.e. Lm.2HF33.scf.fasta), therefore leading to a
final assembled genome of total size 3,097,550 bps (Nres),
represented in 26 contigs (Nseq), with 37.79% GC-content.
The tab-delimited file Lm.2HF33.scf.info.txt (displayed below) shows different residue statistics for each final scaffold sequence. Note that each FASTA header (info files Lm.2HF33.all.fasta and Lm.2HF33.scf.fasta) contains the k-mer coverage depth returned by SPAdes (covk), the base coverage depth estimated by fq2dna (covr), and the (low) coverage cutoff derived from the coverage profile analysis (cutoff).
#Seq Nres A C G T N %A %C %G %T %N %AT %GC Pval CPU
NODE_0001 1091380 330565 212930 195683 352102 100 30.28% 19.51% 17.92% 32.26% 0.00% 62.56% 37.44% 0.0245 C
NODE_0002 532745 154855 110898 91560 175432 0 29.06% 20.81% 17.18% 32.92% 0.00% 62.00% 38.00% 0.5029 C
NODE_0003 361915 113234 65662 74541 108478 0 31.28% 18.14% 20.59% 29.97% 0.00% 61.27% 38.73% 0.0013 C
NODE_0004 206465 68180 35149 42251 60880 5 33.02% 17.02% 20.46% 29.48% 0.00% 62.52% 37.48% 0.0087 C
NODE_0005 126584 37676 26293 21299 41316 0 29.76% 20.77% 16.82% 32.63% 0.00% 62.41% 37.59% 0.2969 C
NODE_0006 125913 41257 22119 25733 36804 0 32.76% 17.56% 20.43% 29.22% 0.00% 62.00% 38.00% 0.9984 C
NODE_0007 108604 32364 22660 17840 35740 0 29.80% 20.86% 16.42% 32.90% 0.00% 62.71% 37.29% 0.0347 C
NODE_0008 103527 29748 21599 16744 35344 92 28.73% 20.86% 16.17% 34.13% 0.08% 62.94% 37.06% 0.0106 C
NODE_0009 99087 32746 17156 20513 28672 0 33.04% 17.31% 20.70% 28.93% 0.00% 61.99% 38.01% 0.8111 C
NODE_0010 97524 28264 20285 17249 31726 0 28.98% 20.80% 17.68% 32.53% 0.00% 61.52% 38.48% 0.0114 C
NODE_0011 64151 21653 11031 13320 18147 0 33.75% 17.19% 20.76% 28.28% 0.00% 62.05% 37.95% 0.9605 C
NODE_0012 61185 17687 12639 9717 21142 0 28.90% 20.65% 15.88% 34.55% 0.00% 63.47% 36.53% 0.0000 P
NODE_0013 46823 14085 9566 7968 15204 0 30.08% 20.43% 17.01% 32.47% 0.00% 62.56% 37.44% 0.1593 C
NODE_0014 35062 10004 7708 6005 11345 0 28.53% 21.98% 17.12% 32.35% 0.00% 60.89% 39.11% 0.0134 C
NODE_0015 22536 6000 4813 3643 8080 0 26.62% 21.35% 16.16% 35.85% 0.00% 62.48% 37.52% 0.1791 U
NODE_0016 4992 1394 1052 1469 1071 6 27.92% 21.07% 29.42% 21.45% 0.12% 49.44% 50.56% 0.0000 U
NODE_0017 2806 724 661 386 1035 0 25.80% 23.55% 13.75% 36.88% 0.00% 62.69% 37.31% 0.5037 U
NODE_0018 1664 524 292 266 582 0 31.49% 17.54% 15.98% 34.97% 0.00% 66.47% 33.53% 0.0080 U
NODE_0019 1251 317 215 219 498 2 25.33% 17.18% 17.50% 39.80% 0.15% 65.26% 34.74% 0.0523 C
NODE_0020 1126 398 199 265 263 1 35.34% 17.67% 23.53% 23.35% 0.08% 58.76% 41.24% 0.1237 U
NODE_0021 846 181 200 234 227 4 21.39% 23.64% 27.65% 26.83% 0.47% 48.46% 51.54% 0.0150 U
NODE_0022 628 239 86 114 186 3 38.05% 13.69% 18.15% 29.61% 0.47% 68.00% 32.00% 0.0648 U
NODE_0023 509 135 101 120 153 0 26.52% 19.84% 23.57% 30.05% 0.00% 56.59% 43.41% 0.0251 U
NODE_0027 407 80 99 70 158 0 19.65% 24.32% 17.19% 38.82% 0.00% 58.48% 41.52% 0.1986 C
The last column CPU also shows the classification of each scaffold
sequence into the category ‘Chromosome’, ‘Plasmid’ or ‘Undetermined’
(i.e. C, P, U, respectively). Such a classification shows that the 12th
scaffold sequence (i.e. NODE_0012
) is likely a plasmid
(subsequent BLAST searches indeed show that it is almost identical to pLM6179).
The tab-delimited file Lm.2HF33.scf.amb.txt (displayed below) lists the different assembled positions that are likely non mono-allelic.
#seq lgt pos base %A %C %G %T %gap
NODE_0001 1091380 481973 N 48.9 0 51.1 0 0
NODE_0001 1091380 481988 N 51.1 0 0 48.9 0
NODE_0004 206465 203562 N 26.1 0 0 73.9 0
NODE_0004 206465 203837 N 0 40.0 0 60.0 0
NODE_0004 206465 203838 N 0 40.8 0 59.2 0
NODE_0004 206465 203840 N 62.5 0 37.5 0 0
NODE_0004 206465 204068 N 38.3 0 61.7 0 0
NODE_0008 103527 13568 N 73.5 0 26.5 0 0
NODE_0008 103527 13622 N 28.8 0 71.2 0 0
NODE_0008 103527 13673 N 0 69.8 0 30.2 0
NODE_0016 4992 374 N 77.8 0 22.2 0 0
NODE_0016 4992 376 N 20.3 0 79.7 0 0
NODE_0019 1251 392 N 60.6 0 0 39.4 0
NODE_0021 846 314 N 0 0 0 59.7 40.3
NODE_0021 846 495 N 75.5 0 0 0 24.5
NODE_0021 846 497 N 0 0 76.0 24.0 0
NODE_0021 846 498 N 23.1 0 75.0 1.9 0
NODE_0022 628 468 N 0 71.7 0 28.3 0
NODE_0022 628 475 N 0 71.1 0 28.9 0
As the overall sequencing seems to not suffer from contamination (see above), (some of) these different ambiguous positions can be indicative of repeat regions of the genome.
Bankevich A, Nurk S, Antipov D, Gurevich A, Dvorkin M, Kulikov AS, Lesin V, Nikolenko S, Pham S, Prjibelski A, Pyshkin A, Sirotkin A, Vyahhi N, Tesler G, Alekseyev MA, Pevzner PA (2012) SPAdes: A New Genome Assembly Algorithm and Its Applications to Single-Cell Sequencing. Journal of Computational Biology, 19(5):455-477. doi:10.1089/cmb.2012.0021.
Duru IC, Andreevskaya M, Laine P, Rode TN, Ylinen A, Løvdal T, Bar N, Crauwels P, Riedel CU, Bucur FI, Nicolau AI, Auvinen P (2020) Genomic characterization of the most barotolerant Listeria monocytogenes RO15 strain compared to reference strains used to evaluate food high pressure processing. BMC Genomics, 21:455. doi:10.1186/s12864-020-06819-0.
Lindner MS, Kollock M, Zickmann F, Renard BY (2013) Analyzing genome coverage profiles with applications to quality control in metagenomics. Bioinformatics, 29(10):1260-1267. doi:10.1093/bioinformatics/btt147.
Roguski L, Deorowicz S (2014) DSRC 2: Industry-oriented compression of FASTQ files. Bioinformatics, 30(15):2213-2215. doi:10.1093/bioinformatics/btu208.
Schwengers O, Barth P, Falgenhauer L, Hain T, Chakraborty T, Goesmann A (2020) Platon: identification and characterization of bacterial plasmid contigs in short-read draft assemblies exploiting protein sequence-based replicon distribution scores. Microbial Genomics, 6(10):mgen000398. doi:10.1099/mgen.0.000398.
Seemann T (2014) Prokka: rapid prokaryotic genome annotation. Bioinformatics, 30(14):2068-2069. doi:10.1093/bioinformatics/btu153.
Abou Fayad A, Rafei R, Njamkepo E, Ezzeddine J, Hussein H, Sinno S, Gerges J-R, Barada S, Sleiman A, Assi M, Baakliny M, Hamedeh L, Mahfouz R, Dabboussi F, Feghali R, Mohsen Z, Rady A, Ghosn N, Abiad F, Abubakar A, Barakat A, Wauquier N, Quilici M-L, Hamze M, Weill F-X, Matar GM (2024) An unusual two-strain cholera outbreak in Lebanon, 2022-2023: a genomic epidemiology study. Nature Communications, 15:6963. doi:10.1038/s41467-024-51428-0
Crippa C, Pasquali F, Rodrigues C, De Cesare A, Lucchi A, Gambi L, Manfreda G, Brisse S, Palma F (2023) Genomic features of Klebsiella isolates from artisanal ready-to-eat food production facilities. Scientific Reports, 13:10957. doi:10.1038/s41598-023-37821-7
Hawkey J, Frézal L, Tran Dien A, Zhukova A, Brown D, Chattaway MA, Simon S, Izumiya H, Fields PI, De Lappe N, Kaftyreva L, Xu X, Isobe J, Clermont D, Njamkepo E, Akeda Y, Issenhuth-Jeanjean S, Makarova M, Wang Y, Hunt M, Jenkins BM, Ravel M, Guibert V, Serre E, Matveeva Z, Fabre L, Cormican M, Yue M, Zhu B, Morita M, Iqbal Z, Nodari CS, Pardos de la Gandara M, Weill F-X (2024) Genomic perspective on the bacillus causing paratyphoid B fever. Nature Communications, 15:10143. doi:10.1038/s41467-024-54418-4
Kämpfer P, Glaeser SP, Busse H-J, McInroy JA, Clermont D, Criscuolo A (2022) Pseudoneobacillus rhizosphaerae gen. nov., sp. nov., isolated from maize root rhizosphere. International Journal of Systematic and Evolutionary Biology, 72:5. doi:10.1099/ijsem.0.005367
Kämpfer P, Lipski A, McInroy JA, Clermont D, Criscuolo A, Glaeser SP (2022) Bacillus rhizoplanae sp. nov. from maize roots. International Journal of Systematic and Evolutionary Biology, 72:7. doi:10.1099/ijsem.0.005450
Kämpfer P, Lipski A, Lamothe L, Clermont D, Criscuolo A, McInroy JA, Glaeser SP (2023) Paenibacillus plantiphilus sp. nov. from the plant environment of Zea mays. Antonie van Leeuwenhoek, 116:883-892. doi:10.1007/s10482-023-01852-x
Kämpfer P, Glaeser SP, McInroy JA, Busse H-J, Clermont D, Criscuolo A (2024) Description of Cohnella rhizoplanae sp. nov., isolated from the root surface of soybean (Glycine max). Antonie van Leeuwenhoek, 118:41. doi:10.1007/s10482-024-02051-y
Palma F, Mangone I, Janowicz A, Moura A, Chiaverini A, Torresi M, Garofolo G, Criscuolo A, Brisse S, Di Pasquale A, Cammà C, Radomski N (2022) In vitro and in silico parameters for precise cgMLST typing of Listeria monocytogenes. BMC Genomics, 23:235. doi:10.1186/s12864-022-08437-4
Rahi P, Mühle E, Scandola C, Touak G, Clermont D (2024) Genome sequence-based identification of Enterobacter strains and description of Enterobacter pasteurii sp. nov. Microbiology Spectrum, 12:e03150-23. doi:10.1128/spectrum.03150-23
Shelomi M, Han C-J, Chen W-M, Chen H-K, Liaw S-J, Mühle E, Clermont5 D (2023) Chryseobacterium oryctis sp. nov., isolated from the gut of the beetle Oryctes rhinoceros, and Chryseobacterium kimseyorum sp. nov., isolated from a stick insect rearing cage. International Journal of Systematic and Evolutionary Biology, 73:4. doi:10.1099/ijsem.0.005813
Vautrin N, Alexandre K, Pestel-Caron M, Bernard E, Fabre R, Leoz M, Dahyot S, Caron F (2023) Contribution of Antibiotic Susceptibility Testing and CH Typing Compared to Next-Generation Sequencing for the Diagnosis of Recurrent Urinary Tract Infections Due to Genetically Identical Escherichia coli Isolates: a Prospective Cohort Study of Cystitis in Women. Microbiology Spectrum, 11(4):e02785-22. doi:10.1128/spectrum.02785-22