hidive/
main.rs

1//! Hidive is a targeted genome co-assembler for biobank-scale long-read and short-read data.
2//! It is designed to assemble high-diversity loci from long-read data and to co-assemble these
3//! loci with short-read data. Hidive is built on top of the Skydive library, which provides
4//! a series-parallel graph representation of long-read data and a linked de Bruijn graph
5//! representation of short-read data.
6//!
7//! # Quick Start
8//!
9//! The following commands will compile the Rust and Python codebase.
10//!
11//! ## Install Rust
12//!
13//! ```sh
14//! curl --proto '=https' --tlsv1.2 -sSf https://sh.rustup.rs | sh
15//! ```
16//!
17//! ## Download and Build Hidive
18//!
19//! ```sh
20//! git clone https://github.com/broadinstitute/hidive.git
21//! cd hidive
22//! cargo build --release
23//! ```
24//!
25//! ## Configure Python Environment
26//!
27//! ```sh
28//! python -m venv venv
29//! . venv/bin/activate
30//! pip install -r dev-requirements.txt
31//! ```
32//!
33//! ## Build Hidive's Python Codebase, Pydive
34//!
35//! ```sh
36//! cd src/pydive/
37//! maturin develop --release
38//! ```
39//! # Prerequisites
40//! Hidive is designed to access local files or data in Google Cloud Storage (GCS).
41//! Within certain cloud-computing environments (i.e. Terra, All of Us Researcher Workbench),
42//! access to GCS is already configured. For accessing files in GCS on your local machine,
43//! you will also need to install the [Google Cloud CLI](https://cloud.google.com/sdk/docs/install-sdk).
44//! Then, configure your [Application Default Credentials (ADC)](https://cloud.google.com/docs/authentication/provide-credentials-adc#local-dev).
45//! If accessing [requester pays buckets](https://cloud.google.com/storage/docs/requester-pays)
46//! set the following environment variable before running hidive commands:
47//! ```sh
48//!  export GCS_REQUESTER_PAYS_PROJECT=<Google Project ID>
49//! ```
50
51use crate::rescue::SearchOption;
52use std::path::PathBuf;
53
54use clap::{Parser, Subcommand};
55
56mod assemble;
57mod build;
58mod cluster;
59mod coassemble;
60mod consensus;
61mod correct;
62mod eval_model;
63mod fetch;
64mod filter;
65mod genotype;
66mod impute;
67mod phase;
68mod recruit;
69mod rescue;
70mod train;
71mod trim;
72
73#[derive(Debug, Parser)]
74#[clap(name = "hidive")]
75#[clap(about = "Analysis of high-diversity loci through genome co-assembly of long/short reads.", long_about = None)]
76// #[clap(author = "Kiran V Garimella (kiran@broadinstitute.org)")]
77struct Cli {
78    #[clap(subcommand)]
79    command: Commands,
80}
81
82const DEFAULT_KMER_SIZE: usize = 17;
83const DEFAULT_WINDOW_SIZE: usize = 1000;
84
85#[derive(Debug, Subcommand)]
86enum Commands {
87    /// Train a error correction model with reads and ground truth assemblies.
88    #[clap(arg_required_else_help = true)]
89    Train {
90        /// Output path for trained model.
91        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
92        output: PathBuf,
93
94        /// Kmer-size.
95        #[clap(short, long, value_parser, default_value_t = DEFAULT_KMER_SIZE)]
96        kmer_size: usize,
97
98        // /// Test split.
99        // #[clap(short, long, value_parser, default_value_t = 0.2)]
100        // test_split: f32,
101        /// Number of training iterations.
102        #[clap(short, long, value_parser, default_value_t = 200)]
103        iterations: usize,
104
105        /// Indexed WGS BAM, CRAM, or FASTA files from which to extract relevant sequences.
106        #[clap(short, long, value_parser, required = true)]
107        long_read_seq_paths: Vec<PathBuf>,
108
109        /// Indexed WGS BAM, CRAM, or FASTA files from which to extract relevant sequences.
110        #[clap(short, long, value_parser, required = false)]
111        short_read_seq_paths: Vec<PathBuf>,
112
113        /// Indexed BAM files to use as ground truth (usually from ultra-high-quality assemblies).
114        #[clap(long, value_parser, required = true, num_args = 2)]
115        truth_seq_paths: Vec<PathBuf>,
116
117        /// Indexed BAM files to use as ground truth (usually from ultra-high-quality assemblies).
118        /// for testing
119        #[clap(long, value_parser, required = false)]
120        test_long_read_seq_paths: Vec<PathBuf>,
121
122        /// Indexed BAM files to use as ground truth (usually from ultra-high-quality assemblies).
123        /// for testing
124        #[clap(long, value_parser, required = false)]
125        test_short_read_seq_paths: Vec<PathBuf>,
126
127        /// Indexed BAM files to use as ground truth (usually from ultra-high-quality assemblies).
128        /// for testing
129        #[clap(long, value_parser, required = false, num_args = 2)]
130        test_truth_seq_paths: Vec<PathBuf>,
131
132        /// Turn on debug mode.
133        #[clap(short, long, value_parser)]
134        debug: bool,
135    },
136
137    /// Evaluate a trained model using long- and (optionally) short-read data with ground truth assemblies.
138    #[clap(arg_required_else_help = true)]
139    EvalModel {
140        /// Output path for trained model.
141        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
142        output: PathBuf,
143
144        /// Kmer-size.
145        #[clap(short, long, value_parser, default_value_t = DEFAULT_KMER_SIZE)]
146        kmer_size: usize,
147
148        /// Indexed WGS BAM, CRAM, or FASTA files from which to extract relevant sequences.
149        #[clap(short, long, value_parser, required = true)]
150        long_read_seq_paths: Vec<PathBuf>,
151
152        /// Indexed WGS BAM, CRAM, or FASTA files from which to extract relevant sequences.
153        #[clap(short, long, value_parser, required = false)]
154        short_read_seq_paths: Vec<PathBuf>,
155
156        /// Indexed BAM files to use as ground truth (usually from ultra-high-quality assemblies).
157        #[clap(long, value_parser, required = true, num_args = 2)]
158        truth_seq_paths: Vec<PathBuf>,
159
160        /// Path to trained model.
161        #[clap(short, long, value_parser, required = true)]
162        model_path: PathBuf,
163    },
164
165    /// Stream loci from CRAM/BAM/FASTA files stored locally or in Google Cloud Storage.
166    #[clap(arg_required_else_help = true)]
167    Fetch {
168        /// Output path for FASTQ file with reads spanning locus of interest.
169        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
170        output: PathBuf,
171
172        /// One or more genomic loci ("contig:start-stop[|name]", or BED format) to extract from WGS BAM files.
173        #[clap(short, long, value_parser, required = true)]
174        loci: Vec<String>,
175
176        /// Include unmapped reads.
177        #[clap(short, long, value_parser, default_value_t = 0)]
178        padding: u64,
179
180        /// Include reads from only one haplotagged haplotype.
181        #[clap(long, value_parser)]
182        haplotype: Option<u8>,
183
184        /// Indexed WGS BAM, CRAM, or FASTA files from which to extract relevant sequences.
185        #[clap(required = true, value_parser)]
186        seq_paths: Vec<PathBuf>,
187    },
188
189    /// Find more sequences (both aligned and unaligned) overlapping previously fetched reads.
190    #[clap(arg_required_else_help = true)]
191    Rescue {
192        /// Output path for FASTQ file with reads spanning locus of interest.
193        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
194        output: PathBuf,
195
196        /// Kmer-size.
197        #[clap(short, long, value_parser, default_value_t = DEFAULT_KMER_SIZE)]
198        kmer_size: usize,
199
200        /// Minimum percentage of k-mers to require before examining a read more carefully.
201        #[clap(short, long, value_parser, default_value_t = 70)]
202        min_kmers_pct: usize,
203
204        /// Option to search for reads based on alignment status or regions of interest.
205        #[clap(short, long, default_value = "contig-and-interval")]
206        search_option: SearchOption,
207
208        /// Reference FASTA (for guessing where reads mapped based on input FASTA filter files).
209        #[clap(short, long, value_parser, required = true)]
210        ref_path: Option<PathBuf>,
211
212        /// One or more genomic loci ("contig:start-stop[|name]", or BED format) to extract from WGS BAM files.
213        #[clap(short, long, value_parser, required = false)]
214        loci: Option<Vec<String>>,
215
216        /// FASTA/FASTQ files with reads to use as a filter for finding more reads.
217        #[clap(short, long, value_parser, required = true)]
218        fastx_paths: Vec<PathBuf>,
219
220        /// Indexed WGS BAM, CRAM, or FASTA files from which to find relevant sequences.
221        #[clap(required = true, value_parser)]
222        seq_paths: Vec<PathBuf>,
223    },
224
225    /// From rescued reads, recruit subset reads that are similar to input long reads.
226    #[clap(arg_required_else_help = true)]
227    Recruit {
228        /// Output path for FASTA file with reads spanning locus of interest.
229        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
230        output: PathBuf,
231
232        /// Kmer-size.
233        #[clap(short, long, value_parser, default_value_t = DEFAULT_KMER_SIZE)]
234        kmer_size: usize,
235
236        /// Minimum percentage of k-mers to require before examining a read more carefully.
237        #[clap(short, long, value_parser, default_value_t = 70)]
238        min_kmers_pct: usize,
239
240        /// FASTA files with reads to use as a filter for finding more reads.
241        #[clap(short, long, value_parser, required = true)]
242        fasta_paths: Vec<PathBuf>,
243
244        /// FASTA files from which to extract relevant sequences.
245        #[clap(required = true, value_parser)]
246        seq_paths: Vec<PathBuf>,
247    },
248
249    /// Further filter rescued reads to those most closely matching a long-read draft assembly.
250    #[clap(arg_required_else_help = true)]
251    Filter {
252        /// Output path for filtered short-read sequences.
253        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
254        output: PathBuf,
255
256        /// FASTA files with short-read sequences (may contain one or more samples).
257        #[clap(short, long, value_parser)]
258        gfa_path: PathBuf,
259
260        /// FASTA files with short-read sequences (may contain one or more samples).
261        #[clap(required = true, value_parser)]
262        short_read_fasta_paths: Vec<PathBuf>,
263    },
264
265    /// Cluster sequences based on k-mer presence/absence.
266    #[clap(arg_required_else_help = true)]
267    Cluster {
268        /// Output path for clustered sequences.
269        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
270        output: PathBuf,
271
272        /// Sample name.
273        #[clap(short, long, required = true, value_parser)]
274        sample_name: String,
275
276        /// One or more genomic loci ("contig:start-stop[|name]", or BED format) to extract from WGS BAM files.
277        #[clap(short, long, value_parser, required = true)]
278        from_loci: Vec<String>,
279
280        /// One or more genomic loci ("contig:start-stop[|name]", or BED format) to which extracted reads should be aligned.
281        #[clap(short, long, value_parser, required = true)]
282        to_loci: Vec<String>,
283
284        /// Reference FASTA (for guessing where reads mapped based on input FASTA filter files).
285        #[clap(short, long, value_parser, required = true)]
286        ref_path: PathBuf,
287
288        /// BAM file with integrated long- and short-read data.
289        #[clap(required = true, value_parser)]
290        bam_path: PathBuf,
291    },
292
293    /// Compute consensus sequences from aligned reads.
294    #[clap(arg_required_else_help = true)]
295    Consensus {
296        /// Output path for clustered sequences.
297        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
298        output: PathBuf,
299
300        /// One or more genomic loci ("contig:start-stop[|name]", or BED format) to extract from WGS BAM files.
301        #[clap(short, long, value_parser, required = true)]
302        loci: Vec<String>,
303
304        // /// Reference FASTA (for guessing where reads mapped based on input FASTA filter files).
305        // #[clap(short, long, value_parser, required = true)]
306        // ref_path: PathBuf,
307
308        // /// FASTA files with short-read sequences (may contain one or more samples).
309        // #[clap(required = true, value_parser)]
310        // short_read_fasta_path: PathBuf,
311        /// BAM file with integrated long- and short-read data.
312        #[clap(required = true, value_parser)]
313        bam_path: PathBuf,
314    },
315
316    /// Trim reads to a specific window around locus.
317    #[clap(arg_required_else_help = true)]
318    Trim {
319        /// Output path for trimmed BAM.
320        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
321        output: PathBuf,
322
323        /// One or more genomic loci ("contig:start-stop") to extract from WGS BAM files.
324        #[clap(short, long, value_parser)]
325        loci: Vec<String>,
326
327        /// Multi-sample BAM file with reads spanning locus of interest.
328        #[clap(required = true, value_parser)]
329        bam_path: PathBuf,
330    },
331
332    /// Build series-parallel graph from long-read data in multi-sample FASTA file.
333    #[clap(arg_required_else_help = true)]
334    Build {
335        /// Output path for series-parallel graph.
336        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
337        output: PathBuf,
338
339        /// Kmer-size
340        #[clap(short, long, value_parser, default_value_t = DEFAULT_KMER_SIZE)]
341        kmer_size: usize,
342
343        /// Name of sequence to use as reference.
344        #[clap(short, long, value_parser, required = true)]
345        reference_name: String,
346
347        /// Multi-sample FASTA file with reads spanning locus of interest.
348        #[clap(required = true, value_parser)]
349        fasta_path: PathBuf,
350    },
351
352    /// Cluster edge matrix and impute missing edges.
353    #[clap(arg_required_else_help = true)]
354    Impute {
355        /// Output path for series-parallel graph with imputed missing edges.
356        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
357        output: PathBuf,
358
359        /// Series-parallel graph.
360        #[clap(required = true, value_parser)]
361        graph: PathBuf,
362    },
363
364    /// Error-correct long reads using a linked de Bruijn graph.
365    #[clap(arg_required_else_help = true)]
366    Assemble {
367        /// Output path for corrected reads.
368        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
369        output: PathBuf,
370
371        /// Output path for corrected reads.
372        #[clap(short, long, value_parser, required = false)]
373        gfa_output: Option<PathBuf>,
374
375        /// Kmer-size
376        #[clap(short, long, value_parser, default_value_t = DEFAULT_KMER_SIZE)]
377        kmer_size: usize,
378
379        /// Trained error-cleaning model.
380        #[clap(short, long, required = true, value_parser)]
381        model_path: PathBuf,
382
383        /// FASTA files with long-read sequences (may contain one or more samples).
384        #[clap(required = true, value_parser)]
385        long_read_fasta_path: PathBuf,
386
387        /// FASTA files with short-read sequences (may contain one or more samples).
388        #[clap(required = true, value_parser)]
389        short_read_fasta_path: PathBuf,
390    },
391
392    /// Error-clean long reads using a linked de Bruijn graph.
393    #[clap(arg_required_else_help = true)]
394    Correct {
395        /// Output path for corrected reads.
396        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
397        output: PathBuf,
398
399        /// One or more genomic loci ("contig:start-stop[|name]", or BED format) to extract from WGS BAM files.
400        #[clap(short, long, value_parser, required = true)]
401        loci: Vec<String>,
402
403        /// Kmer-size
404        #[clap(short, long, value_parser, default_value_t = DEFAULT_KMER_SIZE)]
405        kmer_size: usize,
406
407        /// Window size
408        #[clap(short, long, value_parser, default_value_t = DEFAULT_WINDOW_SIZE)]
409        window_size: usize,
410
411        /// Trained error-cleaning model.
412        #[clap(short, long, required = true, value_parser)]
413        model_path: PathBuf,
414
415        /// FASTA files with long-read sequences (may contain one or more samples).
416        #[clap(required = true, value_parser)]
417        long_read_fasta_path: PathBuf,
418
419        /// FASTA files with short-read sequences (may contain one or more samples).
420        #[clap(required = true, value_parser)]
421        short_read_fasta_path: PathBuf,
422    },
423
424    /// Co-assemble target locus from long- and short-read data using a linked de Bruijn graph.
425    #[clap(arg_required_else_help = true)]
426    Coassemble {
427        /// Output path for assembled short-read sequences.
428        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
429        output: PathBuf,
430
431        /// Kmer-size
432        #[clap(short, long, value_parser, default_value_t = DEFAULT_KMER_SIZE)]
433        kmer_size: usize,
434
435        /// Trained error-cleaning model.
436        #[clap(short, long, required = true, value_parser)]
437        model_path: PathBuf,
438
439        /// FASTA files with reference subsequences.
440        #[clap(short, long, required = true, value_parser)]
441        reference_fasta_paths: Vec<PathBuf>,
442
443        /// FASTA files with long-read sequences (may contain one or more samples).
444        #[clap(required = true, value_parser)]
445        long_read_fasta_path: PathBuf,
446
447        /// FASTA files with short-read sequences (may contain one or more samples).
448        #[clap(required = false, value_parser)]
449        short_read_fasta_path: PathBuf,
450    },
451
452    /// Phase reads.
453    #[clap(arg_required_else_help = true)]
454    Phase {
455        /// Output path for assembled short-read sequences.
456        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
457        output: PathBuf,
458
459        /// Sample name.
460        #[clap(short, long, required = true, value_parser)]
461        sample_name: String,
462
463        /// One or more genomic loci ("contig:start-stop[|name]", or BED format) to extract from WGS BAM files.
464        #[clap(short, long, value_parser, required = true)]
465        loci: Vec<String>,
466
467        /// FASTA reference sequence.
468        #[clap(short, long, required = true, value_parser)]
469        reference_fasta_path: PathBuf,
470
471        /// BAM file with integrated long- and short-read data.
472        #[clap(required = true, value_parser)]
473        bam_path: PathBuf,
474    },
475
476    /// Re-genotype variants.
477    #[clap(arg_required_else_help = true)]
478    Genotype {
479        /// Output path for assembled short-read sequences.
480        #[clap(short, long, value_parser, default_value = "/dev/stdout")]
481        output: PathBuf,
482
483        /// Sample name.
484        #[clap(short, long, required = true, value_parser)]
485        sample_name: String,
486
487        /// One or more genomic loci ("contig:start-stop[|name]", or BED format) to extract from WGS BAM files.
488        #[clap(short, long, value_parser, required = true)]
489        loci: Vec<String>,
490
491        /// FASTA reference sequence.
492        #[clap(short, long, required = true, value_parser)]
493        reference_fasta_path: PathBuf,
494
495        /// VCF file
496        #[clap(required = true, value_parser)]
497        vcf_path: PathBuf,
498
499        /// BAM file with integrated long- and short-read data.
500        #[clap(required = true, value_parser)]
501        bam_path: PathBuf,
502    },
503}
504
505fn main() {
506    let args = Cli::parse();
507
508    skydive::elog!("Hidive version {}", env!("CARGO_PKG_VERSION"));
509    skydive::elog!("{:?}", args);
510
511    let start_time = std::time::Instant::now();
512
513    match args.command {
514        Commands::Train {
515            output,
516            kmer_size,
517            // test_split,
518            iterations,
519            long_read_seq_paths,
520            short_read_seq_paths,
521            truth_seq_paths,
522            test_long_read_seq_paths,
523            test_short_read_seq_paths,
524            test_truth_seq_paths,
525            debug,
526        } => {
527            train::start(
528                &output,
529                kmer_size,
530                iterations,
531                // test_split,
532                &long_read_seq_paths,
533                &short_read_seq_paths,
534                &truth_seq_paths,
535                &test_long_read_seq_paths,
536                &test_short_read_seq_paths,
537                &test_truth_seq_paths,
538                debug,
539            );
540        }
541        Commands::EvalModel {
542            output,
543            kmer_size,
544            long_read_seq_paths,
545            short_read_seq_paths,
546            truth_seq_paths,
547            model_path,
548        } => {
549            eval_model::start(
550                &output,
551                kmer_size,
552                &long_read_seq_paths,
553                &short_read_seq_paths,
554                &truth_seq_paths,
555                &model_path,
556            );
557        }
558        Commands::Fetch {
559            output,
560            loci,
561            padding,
562            haplotype,
563            seq_paths,
564        } => {
565            fetch::start(&output, &loci, padding, haplotype, &seq_paths);
566        }
567        Commands::Rescue {
568            output,
569            kmer_size,
570            min_kmers_pct,
571            search_option,
572            ref_path,
573            loci,
574            fastx_paths,
575            seq_paths,
576        } => {
577            rescue::start(
578                &output,
579                kmer_size,
580                min_kmers_pct,
581                search_option,
582                ref_path,
583                loci,
584                &fastx_paths,
585                &seq_paths,
586            );
587        }
588        Commands::Recruit {
589            output,
590            kmer_size,
591            min_kmers_pct,
592            fasta_paths,
593            seq_paths,
594        } => {
595            recruit::start(&output, kmer_size, min_kmers_pct, &fasta_paths, &seq_paths);
596        }
597        Commands::Filter {
598            output,
599            gfa_path,
600            short_read_fasta_paths,
601        } => {
602            filter::start(&output, &gfa_path, &short_read_fasta_paths);
603        }
604        Commands::Cluster {
605            output,
606            sample_name,
607            from_loci,
608            to_loci,
609            ref_path,
610            bam_path,
611        } => {
612            cluster::start(
613                &output,
614                &sample_name,
615                &from_loci,
616                &to_loci,
617                &ref_path,
618                &bam_path,
619            );
620        }
621        Commands::Consensus {
622            output,
623            loci,
624            // ref_path,
625            // short_read_fasta_path,
626            bam_path,
627        } => {
628            consensus::start(&output, &loci, &bam_path);
629        }
630        Commands::Trim {
631            output,
632            loci,
633            bam_path,
634        } => {
635            trim::start(&output, &loci, &bam_path);
636        }
637        Commands::Build {
638            output,
639            kmer_size,
640            fasta_path,
641            reference_name,
642        } => {
643            build::start(&output, kmer_size, &fasta_path, reference_name);
644        }
645        Commands::Impute { output, graph } => {
646            impute::start(&output, &graph);
647        }
648        Commands::Assemble {
649            output,
650            gfa_output,
651            kmer_size,
652            model_path,
653            long_read_fasta_path,
654            short_read_fasta_path,
655        } => {
656            assemble::start(
657                &output,
658                gfa_output,
659                kmer_size,
660                &model_path,
661                &long_read_fasta_path,
662                &short_read_fasta_path,
663            );
664        }
665        Commands::Correct {
666            output,
667            loci,
668            kmer_size,
669            window_size,
670            model_path,
671            long_read_fasta_path,
672            short_read_fasta_path,
673        } => {
674            correct::start(
675                &output,
676                &loci,
677                kmer_size,
678                window_size,
679                &model_path,
680                &long_read_fasta_path,
681                &short_read_fasta_path,
682            );
683        }
684        Commands::Coassemble {
685            output,
686            kmer_size,
687            model_path,
688            reference_fasta_paths,
689            long_read_fasta_path,
690            short_read_fasta_path,
691        } => {
692            coassemble::start(
693                &output,
694                kmer_size,
695                &model_path,
696                &reference_fasta_paths,
697                long_read_fasta_path,
698                short_read_fasta_path,
699            );
700        }
701        Commands::Phase {
702            output,
703            sample_name,
704            loci,
705            reference_fasta_path,
706            bam_path,
707        } => {
708            phase::start(
709                &output,
710                &sample_name,
711                &loci,
712                &reference_fasta_path,
713                &bam_path,
714            );
715        }
716        Commands::Genotype {
717            output,
718            sample_name,
719            loci,
720            reference_fasta_path,
721            vcf_path,
722            bam_path,
723        } => {
724            genotype::start(
725                &output,
726                &sample_name,
727                &loci,
728                &reference_fasta_path,
729                &vcf_path,
730                &bam_path,
731            );
732        }
733    }
734
735    skydive::elog!("Complete. Elapsed time: {}.", elapsed_time(start_time));
736}
737
738fn elapsed_time(start_time: std::time::Instant) -> String {
739    let end_time = std::time::Instant::now();
740    let elapsed_time = end_time.duration_since(start_time);
741
742    let elapsed_secs = elapsed_time.as_secs_f64();
743
744    if elapsed_secs < 60.0 {
745        format!("{:.2} seconds", elapsed_secs)
746    } else if elapsed_secs < 3600.0 {
747        format!("{:.2} minutes", elapsed_secs / 60.0)
748    } else if elapsed_secs < 86400.0 {
749        format!("{:.2} hours", elapsed_secs / 3600.0)
750    } else {
751        format!("{:.2} days", elapsed_secs / 86400.0)
752    }
753}