skesa-rs 0.2.1

Rust port of NCBI's SKESA genome assembler
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
use std::fs::File;
use std::io::{BufWriter, Write};

use clap::{Parser, Subcommand};

use crate::assembler::{self, AssemblerParams};
use crate::contig_output;
use crate::kmer_counter;
use crate::kmer_output;
use crate::reads_getter::ReadsGetter;

fn hardware_threads() -> usize {
    std::thread::available_parallelism()
        .map(|threads| threads.get())
        .unwrap_or(1)
}

fn resolve_cores(cores: i32) -> usize {
    let hardware = hardware_threads();
    if cores == 0 {
        hardware
    } else {
        (cores as usize).min(hardware)
    }
}

fn warn_if_cores_reduced(cores: i32) {
    let hardware = hardware_threads();
    if cores > hardware as i32 {
        eprintln!(
            "WARNING: number of cores was reduced to the hardware limit of {} cores",
            hardware
        );
    }
}

fn flush_file_or_report<W: Write>(writer: &mut W, path: &str) -> bool {
    match writer.flush() {
        Ok(()) => true,
        Err(e) => {
            eprintln!("Can't write to file {}: {}", path, e);
            false
        }
    }
}

fn flush_stdout_or_report<W: Write>(writer: &mut W) -> bool {
    match writer.flush() {
        Ok(()) => true,
        Err(e) => {
            eprintln!("Write failed: {}", e);
            false
        }
    }
}

#[derive(Parser)]
#[command(name = "skesa-rs", version = env!("CARGO_PKG_VERSION"), about = "Rust port of NCBI's SKESA genome assembler\n\nOriginal: SKESA 2.5.1 by Souvorov et al., Genome Biology 2018")]
struct Cli {
    #[command(subcommand)]
    command: Commands,
}

#[derive(Subcommand)]
enum Commands {
    /// De-novo sequence assembler for microbial genomes
    Skesa(SkesaArgs),
    /// Count k-mers in a read set
    Kmercounter(KmercounterArgs),
    /// Target-enriched de-novo assembler (unsupported; use C++ backend)
    Saute(SauteArgs),
    /// Protein-guided target-enriched assembler (unsupported; use C++ backend)
    SauteProt(SauteProtArgs),
    /// Connect contigs using reads and output GFA graph (unsupported; use C++ backend)
    GfaConnector(GfaConnectorArgs),
    /// Compute assembly statistics from a FASTA file
    Stats(StatsArgs),
}

#[derive(Parser)]
struct SkesaArgs {
    /// Input fasta/fastq file(s) (can be specified multiple times, can be gzipped)
    #[arg(long)]
    reads: Vec<String>,

    /// C++ SKESA compatibility spelling for FASTA input files
    #[arg(long)]
    fasta: Vec<String>,

    /// C++ SKESA compatibility spelling for FASTQ input files
    #[arg(long)]
    fastq: Vec<String>,

    /// C++ SKESA compatibility flag. Rust auto-detects gzip per input file.
    #[arg(long, default_value_t = false)]
    gz: bool,

    /// SRA run accession input is not supported by this Rust port
    #[arg(long, alias = "sra_run")]
    sra_run: Vec<String>,

    /// Indicates that single fasta/fastq files contain paired reads
    #[arg(long, alias = "use_paired_ends", default_value_t = false)]
    use_paired_ends: bool,

    /// Output file for contigs (stdout if not specified)
    #[arg(long, alias = "contigs_out")]
    contigs_out: Option<String>,

    /// Number of cores to use (0 = all)
    #[arg(long, default_value_t = 0)]
    cores: i32,

    /// Memory available in GB (only for sorted counter)
    #[arg(long, default_value_t = 32)]
    memory: i32,

    /// Use hash counter instead of sorted counter
    #[arg(long, alias = "hash_count", default_value_t = false)]
    hash_count: bool,

    /// C++ compatibility flag for hash-based read finding (default off)
    #[arg(long, alias = "hash_find", default_value_t = false)]
    hash_find: bool,

    /// Estimated number of distinct k-mers for bloom filter (millions, only for hash counter)
    #[arg(long, alias = "estimated_kmers", default_value_t = 100)]
    estimated_kmers: i32,

    /// Skip bloom filter; use estimated_kmers as hash table size
    #[arg(long, alias = "skip_bloom_filter", default_value_t = false)]
    skip_bloom_filter: bool,

    /// Minimal k-mer length for assembly
    #[arg(long, default_value_t = 21)]
    kmer: i32,

    /// Maximal k-mer length for assembly (0 = auto)
    #[arg(long, alias = "max_kmer", default_value_t = 0)]
    max_kmer: i32,

    /// Number of assembly iterations from minimal to maximal k-mer length
    #[arg(long, default_value_t = 11)]
    steps: i32,

    /// Minimal count for k-mers retained for comparing alternate choices
    #[arg(long, alias = "min_count")]
    min_count: Option<i32>,

    /// Minimum acceptable average count for estimating maximal k-mer length
    #[arg(long, alias = "max_kmer_count")]
    max_kmer_count: Option<i32>,

    /// Percentage of reads containing 19-mer for adapter detection (1.0 disables)
    #[arg(long, alias = "vector_percent", default_value_t = 0.05)]
    vector_percent: f64,

    /// Expected insert size for paired reads (0 = auto)
    #[arg(long, alias = "insert_size", default_value_t = 0)]
    insert_size: i32,

    /// Maximum noise to signal ratio acceptable for extension
    #[arg(long, default_value_t = 0.1)]
    fraction: f64,

    /// Maximal SNP length
    #[arg(long, alias = "max_snp_len", default_value_t = 150)]
    max_snp_len: i32,

    /// Minimal contig length reported in output
    #[arg(long, alias = "min_contig", default_value_t = 200)]
    min_contig: i32,

    /// Allow additional step for SNP discovery
    #[arg(long, alias = "allow_snps", default_value_t = false)]
    allow_snps: bool,

    /// Don't use paired-end information
    #[arg(long, alias = "force_single_ends", default_value_t = false)]
    force_single_ends: bool,

    /// Input file with seeds
    #[arg(long)]
    seeds: Option<String>,

    /// Output FASTA for each iteration
    #[arg(long)]
    all: Option<String>,

    /// Output k-mer file (binary)
    #[arg(long, alias = "dbg_out")]
    dbg_out: Option<String>,

    /// File for histogram
    #[arg(long)]
    hist: Option<String>,

    /// File for connected paired reads
    #[arg(long, alias = "connected_reads")]
    connected_reads: Option<String>,

    /// Output assembly in GFA format
    #[arg(long, alias = "gfa_out")]
    gfa_out: Option<String>,

    /// Use the legacy single-pass k-mer counter (collects all raw k-mers in
    /// one Vec before dedup). Faster on small inputs but inflates peak RSS
    /// by the raw k-mer total — defaults to false to keep RSS comparable
    /// with the bundled C++ build, which uses the bucketed counter always.
    #[arg(long, alias = "single_pass_counter", default_value_t = false)]
    single_pass_counter: bool,
}

#[derive(Parser)]
struct KmercounterArgs {
    /// Input fasta/fastq file(s) (can be specified multiple times)
    #[arg(long)]
    reads: Vec<String>,

    /// SRA run accession input is not supported by this Rust port
    #[arg(long, alias = "sra_run")]
    sra_run: Vec<String>,

    /// K-mer length
    #[arg(long, default_value_t = 21)]
    kmer: i32,

    /// Minimal count for k-mers retained
    #[arg(long, alias = "min_count", default_value_t = 2)]
    min_count: i32,

    /// Percentage of reads containing 19-mer for adapter detection (1.0 disables)
    #[arg(long, alias = "vector_percent", default_value_t = 0.05)]
    vector_percent: f64,

    /// Disable directional filtering in graph
    #[arg(long, alias = "no_strand_info", default_value_t = false)]
    no_strand_info: bool,

    /// Estimated number of distinct k-mers for bloom filter (millions)
    #[arg(long, alias = "estimated_kmers", default_value_t = 100)]
    estimated_kmers: i32,

    /// Skip bloom filter; use estimated_kmers as hash table size
    #[arg(long, alias = "skip_bloom_filter", default_value_t = false)]
    skip_bloom_filter: bool,

    /// De Bruijn graph binary output file
    #[arg(long, alias = "dbg_out")]
    dbg_out: Option<String>,

    /// Text k-mer output file
    #[arg(long, alias = "text_out")]
    text_out: Option<String>,

    /// Histogram output file
    #[arg(long)]
    hist: Option<String>,

    /// Number of cores to use (0 = all)
    #[arg(long, default_value_t = 0)]
    cores: i32,
}

#[derive(Parser)]
struct SauteArgs {
    /// Input fasta/fastq file(s) for reads
    #[arg(long, required = true)]
    reads: Vec<String>,

    /// Input file with reference/target sequences
    #[arg(long, required = true)]
    targets: String,

    /// Output GFA graph file
    #[arg(long)]
    gfa: Option<String>,

    /// K-mer length (default: automatic)
    #[arg(long)]
    kmer: Option<i32>,

    /// Minimal count for k-mers
    #[arg(long, alias = "min_count", default_value_t = 2)]
    min_count: i32,

    /// Number of cores (0 = all)
    #[arg(long, default_value_t = 0)]
    cores: i32,

    /// Percentage for adapter detection (1.0 disables)
    #[arg(long, alias = "vector_percent", default_value_t = 0.05)]
    vector_percent: f64,

    /// Estimated distinct k-mers (millions)
    #[arg(long, alias = "estimated_kmers", default_value_t = 1000)]
    estimated_kmers: i32,

    /// Indicates paired reads
    #[arg(long, alias = "use_paired_ends", default_value_t = false)]
    use_paired_ends: bool,

    /// Extend graph ends using de-novo assembly
    #[arg(long, alias = "extend_ends", default_value_t = false)]
    extend_ends: bool,
}

#[derive(Parser)]
struct SauteProtArgs {
    /// Input fasta/fastq file(s) for reads
    #[arg(long, required = true)]
    reads: Vec<String>,

    /// Input file with protein reference sequences
    #[arg(long, required = true)]
    targets: String,

    /// Output GFA graph file
    #[arg(long)]
    gfa: Option<String>,

    /// NCBI genetic code number (default: 1 = Standard)
    #[arg(long, alias = "genetic_code", default_value_t = 1)]
    genetic_code: u32,

    /// K-mer length (default: automatic)
    #[arg(long)]
    kmer: Option<i32>,

    /// Minimal count for k-mers
    #[arg(long, alias = "min_count", default_value_t = 2)]
    min_count: i32,

    /// Number of cores (0 = all)
    #[arg(long, default_value_t = 0)]
    cores: i32,

    /// Percentage for adapter detection (1.0 disables)
    #[arg(long, alias = "vector_percent", default_value_t = 0.05)]
    vector_percent: f64,

    /// Estimated distinct k-mers (millions)
    #[arg(long, alias = "estimated_kmers", default_value_t = 1000)]
    estimated_kmers: i32,
}

#[derive(Parser)]
struct GfaConnectorArgs {
    /// Input fasta/fastq file(s) for reads
    #[arg(long, required = true)]
    reads: Vec<String>,

    /// Input file with contigs
    #[arg(long, required = true)]
    contigs: String,

    /// Output GFA graph file (stdout if not specified)
    #[arg(long)]
    gfa: Option<String>,

    /// K-mer length
    #[arg(long)]
    kmer: Option<i32>,

    /// Minimal count for k-mers
    #[arg(long, alias = "min_count", default_value_t = 2)]
    min_count: i32,

    /// Number of cores (0 = all)
    #[arg(long, default_value_t = 0)]
    cores: i32,

    /// Maximum extension length
    #[arg(long, alias = "ext_len", default_value_t = 2000)]
    ext_len: i32,

    /// Noise-to-signal ratio threshold
    #[arg(long, default_value_t = 0.1)]
    fraction: f64,
}

#[derive(Parser)]
struct StatsArgs {
    /// Input FASTA file with contigs
    #[arg(required = true)]
    fasta: String,

    /// Minimum contig length to include
    #[arg(long, alias = "min_length", default_value_t = 0)]
    min_length: usize,
}

fn run_stats(args: &StatsArgs) -> i32 {
    match crate::assembly_stats::fasta_lengths(&args.fasta) {
        Ok(mut lengths) => {
            if args.min_length > 0 {
                lengths.retain(|&l| l >= args.min_length);
            }
            let stats = crate::assembly_stats::AssemblyStats::from_lengths(&lengths);
            stats.print();
            0
        }
        Err(e) => {
            eprintln!("{}", e);
            1
        }
    }
}

fn run_gfa_connector(_args: &GfaConnectorArgs) -> i32 {
    eprintln!("gfa_connector is not yet parity-supported in skesa-rs; use the bundled C++ gfa_connector for now");
    1
}

fn run_saute_prot(_args: &SauteProtArgs) -> i32 {
    eprintln!("saute_prot is not yet parity-supported in skesa-rs; use the bundled C++ saute_prot for now");
    1
}

fn run_saute(_args: &SauteArgs) -> i32 {
    eprintln!("saute is not yet parity-supported in skesa-rs; use the bundled C++ saute for now");
    1
}

fn skesa_input_files(args: &SkesaArgs) -> Vec<String> {
    let _gzip_compat_flag = args.gz;
    let mut files = Vec::new();
    files.extend(args.reads.iter().cloned());
    files.extend(args.fasta.iter().cloned());
    files.extend(args.fastq.iter().cloned());
    files.sort();
    let original_len = files.len();
    files.dedup();
    if files.len() != original_len {
        eprintln!("WARNING: duplicate input entries were removed from file list");
    }
    files
}

fn run_skesa(args: &SkesaArgs) -> i32 {
    let output = crate::output::StdioOutput;
    let input_files = skesa_input_files(args);

    if !args.sra_run.is_empty() {
        eprintln!("SRA input is not supported; use --reads with local FASTA/FASTQ files");
        return 1;
    }
    if args.hash_count {
        eprintln!(
            "--hash_count assembly mode is not yet parity-supported in skesa-rs; use the bundled C++ skesa for this mode"
        );
        return 1;
    }
    let _hash_find = args.hash_find;
    if input_files.is_empty() {
        eprintln!("Provide some input reads");
        return 1;
    }
    if args.cores < 0 {
        eprintln!("Value of --cores must be >= 0");
        return 1;
    }
    if args.steps <= 0 {
        eprintln!("Value of --steps must be > 0");
        return 1;
    }
    if args.fraction >= 1.0 {
        eprintln!("Value of --fraction must be < 1 (more than 0.25 is not recommended)");
        return 1;
    }
    if args.fraction < 0.0 {
        eprintln!("Value of --fraction must be >= 0");
        return 1;
    }
    if args.max_snp_len < 0 {
        eprintln!("Value of --max_snp_len must be >= 0");
        return 1;
    }
    if args.kmer <= 0 {
        eprintln!("Value of --kmer must be > 0");
        return 1;
    }
    if args.kmer < 21 || args.kmer % 2 == 0 {
        eprintln!("Kmer must be an odd number >= 21");
        return 1;
    }
    if args.max_kmer < 0 {
        eprintln!("Value of --max_kmer must be > 0");
        return 1;
    }
    if let Some(min_count) = args.min_count {
        if min_count <= 0 {
            eprintln!("Value of --min_count must be > 0");
            return 1;
        }
    }
    if let Some(max_kmer_count) = args.max_kmer_count {
        if max_kmer_count <= 0 {
            eprintln!("Value of --max_kmer_count must be > 0");
            return 1;
        }
    }
    if args.kmer as usize > crate::kmer::MAX_KMER {
        eprintln!("Not supported kmer length");
        return 1;
    }
    if args.max_kmer > 0 && args.max_kmer as usize > crate::kmer::MAX_KMER {
        eprintln!("Not supported kmer length");
        return 1;
    }
    if args.insert_size < 0 {
        eprintln!("Value of --insert_size must be >= 0");
        return 1;
    }
    if args.vector_percent > 1.0 {
        eprintln!("Value of --vector_percent  must be <= 1");
        return 1;
    }
    if args.vector_percent <= 0.0 {
        eprintln!("Value of --vector_percent  must be > 0");
        return 1;
    }
    if args.estimated_kmers <= 0 {
        eprintln!("Value of --estimated_kmers must be > 0");
        return 1;
    }
    if args.memory <= 0 {
        eprintln!("Value of --memory must be > 0");
        return 1;
    }
    if args.memory <= 2 {
        eprintln!("Memory provided is insufficient to do runs in 10 cycles for the read coverage. We find that 16 Gb for 20x coverage of a 5 Mb genome is usually sufficient");
        return 1;
    }
    if args.min_contig <= 0 {
        eprintln!("Value of --min_contig must be > 0");
        return 1;
    }

    // Load reads
    let ncores = resolve_cores(args.cores);
    let rg = match ReadsGetter::new_with_ncores_and_output(
        &input_files,
        args.use_paired_ends,
        ncores,
        &output,
    ) {
        Ok(rg) => rg,
        Err(e) => {
            eprintln!("{}", e);
            return 1;
        }
    };

    // Clip adapters (only clone if clipping is needed)
    let clipped;
    let reads_ref = if args.vector_percent < 1.0 {
        clipped = {
            let mut v = rg.reads().to_vec();
            crate::reads_getter::clip_adapters_with_output(
                &mut v,
                args.vector_percent,
                100,
                &output,
            );
            v
        };
        &clipped[..]
    } else {
        eprintln!("Adapters clip is disabled");
        rg.reads()
    };

    let raw_kmer_num: usize = reads_ref
        .iter()
        .map(|read_pair| {
            read_pair[0].kmer_num(args.kmer as usize) + read_pair[1].kmer_num(args.kmer as usize)
        })
        .sum();
    if let Err(e) = crate::sorted_counter::sorted_counter_plan(
        raw_kmer_num,
        reads_ref.len(),
        args.kmer as usize,
        args.memory as usize,
    ) {
        eprintln!("{}", e);
        return 1;
    }

    let min_count = args.min_count.unwrap_or(2) as usize;
    let max_kmer_count = args.max_kmer_count.unwrap_or(10) as usize;
    let estimate_min_count = args.min_count.is_none() && args.max_kmer_count.is_none();

    crate::sorted_counter::set_single_pass_counter(args.single_pass_counter);

    let params = AssemblerParams {
        min_kmer: args.kmer as usize,
        max_kmer: args.max_kmer as usize,
        steps: args.steps as usize,
        fraction: args.fraction,
        max_snp_len: args.max_snp_len as usize,
        min_count,
        estimate_min_count,
        max_kmer_count,
        force_single_reads: args.force_single_ends,
        insert_size: args.insert_size as usize,
        allow_snps: args.allow_snps,
        ncores,
        memory_gb: args.memory as usize,
        retain_all_iterations: args.all.is_some(),
        retain_all_graphs: args.hist.is_some()
            || args.dbg_out.is_some()
            || args.all.is_some()
            || args.allow_snps,
    };

    let seeds = match load_seed_fasta(args.seeds.as_deref()) {
        Ok(seeds) => seeds,
        Err(e) => {
            eprintln!("{}", e);
            return 1;
        }
    };

    let result = assembler::run_assembly_with_output(reads_ref, &params, &seeds, &output);

    // Output contigs
    let min_contig = args.min_contig as usize;
    let first_graph = result.graphs.first();

    if let Some(ref path) = args.contigs_out {
        let file = match File::create(path) {
            Ok(f) => f,
            Err(e) => {
                eprintln!("Can't open file {}: {}", path, e);
                return 1;
            }
        };
        let mut writer = BufWriter::new(file);
        if let Some((kmer_len, ref kmers)) = first_graph {
            if let Err(e) = contig_output::write_contigs_with_abundance(
                &mut writer,
                &result.contigs,
                kmers,
                *kmer_len,
                min_contig,
            ) {
                eprintln!("Can't write to file {}: {}", path, e);
                return 1;
            }
        } else if let Err(e) =
            contig_output::write_contigs(&mut writer, &result.contigs, min_contig)
        {
            eprintln!("Can't write to file {}: {}", path, e);
            return 1;
        }
        if !flush_file_or_report(&mut writer, path) {
            return 1;
        }
    } else {
        let stdout = std::io::stdout();
        let mut writer = BufWriter::new(stdout.lock());
        if let Some((kmer_len, ref kmers)) = first_graph {
            if let Err(e) = contig_output::write_contigs_with_abundance(
                &mut writer,
                &result.contigs,
                kmers,
                *kmer_len,
                min_contig,
            ) {
                eprintln!("Write failed: {}", e);
                return 1;
            }
        } else if let Err(e) =
            contig_output::write_contigs(&mut writer, &result.contigs, min_contig)
        {
            eprintln!("Write failed: {}", e);
            return 1;
        }
        if !flush_stdout_or_report(&mut writer) {
            return 1;
        }
    }

    // GFA output
    if let Some(ref path) = args.gfa_out {
        let file = match File::create(path) {
            Ok(f) => f,
            Err(e) => {
                eprintln!("Can't open file {}: {}", path, e);
                return 1;
            }
        };
        let mut writer = BufWriter::new(file);
        if let Err(e) = crate::gfa::write_gfa(&mut writer, &result.contigs, min_contig) {
            eprintln!("Can't write GFA to {}: {}", path, e);
            return 1;
        }
        if !flush_file_or_report(&mut writer, path) {
            return 1;
        }
    }

    if let Some(ref path) = args.all {
        let file = match File::create(path) {
            Ok(f) => f,
            Err(e) => {
                eprintln!("Can't open file {}: {}", path, e);
                return 1;
            }
        };
        let mut writer = BufWriter::new(file);
        if let Err(e) = contig_output::write_all_iterations(
            &mut writer,
            &result,
            args.allow_snps,
            !seeds.is_empty(),
        ) {
            eprintln!("Can't write to file {}: {}", path, e);
            return 1;
        }
        if !flush_file_or_report(&mut writer, path) {
            return 1;
        }
    }

    if let Some(ref path) = args.hist {
        let file = match File::create(path) {
            Ok(f) => f,
            Err(e) => {
                eprintln!("Can't open file {}: {}", path, e);
                return 1;
            }
        };
        let mut writer = BufWriter::new(file);
        for (kmer_len, kmers) in &result.graphs {
            let bins = crate::sorted_counter::get_bins(kmers);
            for (count, freq) in bins {
                if let Err(e) = writeln!(writer, "{}\t{}\t{}", kmer_len, count, freq) {
                    eprintln!("Can't write to file {}: {}", path, e);
                    return 1;
                }
            }
        }
        if !flush_file_or_report(&mut writer, path) {
            return 1;
        }
    }

    if let Some(ref path) = args.connected_reads {
        let file = match File::create(path) {
            Ok(f) => f,
            Err(e) => {
                eprintln!("Can't open file {}: {}", path, e);
                return 1;
            }
        };
        let mut writer = BufWriter::new(file);
        let mut num = 0usize;
        let mut read_iter = result.connected_reads.string_iter();
        while !read_iter.at_end() {
            num += 1;
            if let Err(e) = writeln!(writer, ">ConnectedRead_{}\n{}", num, read_iter.get()) {
                eprintln!("Can't write to file {}: {}", path, e);
                return 1;
            }
            read_iter.advance();
        }
        if !flush_file_or_report(&mut writer, path) {
            return 1;
        }
    }

    // DBG output: SKESA's assembler writes sorted graph records.
    if let Some(ref path) = args.dbg_out {
        let file = match File::create(path) {
            Ok(f) => f,
            Err(e) => {
                eprintln!("Can't open file {}: {}", path, e);
                return 1;
            }
        };
        let mut writer = BufWriter::new(file);
        for (_, kmers) in &result.graphs {
            if let Err(e) = crate::graph_io::write_sorted_graph(&mut writer, kmers, true) {
                eprintln!("Can't write to file {}: {}", path, e);
                return 1;
            }
        }
        if !flush_file_or_report(&mut writer, path) {
            return 1;
        }
    }

    eprintln!("DONE");
    0
}

fn load_seed_fasta(path: Option<&str>) -> Result<Vec<String>, String> {
    let Some(path) = path else {
        return Ok(Vec::new());
    };

    let text =
        std::fs::read_to_string(path).map_err(|e| format!("Can't open file {}: {}", path, e))?;
    if text.is_empty() {
        eprintln!("Empty fasta file for seeds");
        return Ok(Vec::new());
    }
    if !text.starts_with('>') {
        return Err(format!("Invalid fasta file format in {}", path));
    }

    let mut seeds = Vec::new();
    for record in text[1..].split('>') {
        let Some(first_ret) = record.find('\n') else {
            return Err(format!("Invalid fasta file format in {}", path));
        };
        let mut sequence = record[first_ret + 1..].replace('\n', "");
        if sequence.ends_with('\r') {
            sequence.pop();
        }
        if sequence.chars().any(|c| !"ACGTYRWSKMDVHBN".contains(c)) {
            return Err(format!("Invalid fasta file format in {}", path));
        }
        seeds.push(sequence);
    }

    Ok(seeds)
}

fn run_kmercounter(args: &KmercounterArgs) -> i32 {
    let output = crate::output::StdioOutput;
    if !args.sra_run.is_empty() {
        eprintln!("SRA input is not supported; use --reads with local FASTA/FASTQ files");
        return 1;
    }
    if args.reads.is_empty() {
        eprintln!("Provide some input reads");
        return 1;
    }
    if args.cores < 0 {
        eprintln!("Value of --cores must be >= 0");
        return 1;
    }
    if args.kmer <= 0 {
        eprintln!("Value of --kmer must be > 0");
        return 1;
    }
    if args.kmer as usize > crate::kmer::MAX_KMER {
        eprintln!("Not supported kmer length");
        return 1;
    }
    if args.min_count <= 0 {
        eprintln!("Value of --min_count must be > 0");
        return 1;
    }
    if args.vector_percent > 1.0 {
        eprintln!("Value of --vector_percent  must be <= 1");
        return 1;
    }
    if args.vector_percent <= 0.0 {
        eprintln!("Value of --vector_percent  must be > 0");
        return 1;
    }
    if args.estimated_kmers <= 0 {
        eprintln!("Value of --estimated_kmers must be > 0");
        return 1;
    }

    // Load reads
    let rg = match ReadsGetter::new_with_ncores_and_output(
        &args.reads,
        false,
        resolve_cores(args.cores),
        &output,
    ) {
        Ok(rg) => rg,
        Err(e) => {
            eprintln!("{}", e);
            return 1;
        }
    };

    // Clip adapters (only clone reads if clipping is needed)
    let clipped;
    let reads_ref = if args.vector_percent < 1.0 {
        clipped = {
            let mut v = rg.reads().to_vec();
            crate::reads_getter::clip_adapters_with_output(
                &mut v,
                args.vector_percent,
                15,
                &output,
            );
            v
        };
        &clipped[..]
    } else {
        eprintln!("Adapters clip is disabled");
        rg.reads()
    };

    // Count k-mers
    let mb: usize = 1_000_000;
    let hash_table = kmer_counter::count_kmers(
        reads_ref,
        args.kmer as usize,
        args.min_count as usize,
        args.estimated_kmers as usize * mb,
        true, // is_stranded
        args.skip_bloom_filter,
    );

    // Text output
    if let Some(ref path) = args.text_out {
        let file = match File::create(path) {
            Ok(f) => f,
            Err(e) => {
                eprintln!("Can't open file {}: {}", path, e);
                return 1;
            }
        };
        let mut writer = BufWriter::new(file);
        if let Err(e) = kmer_output::write_text_output(&mut writer, &hash_table, args.kmer as usize)
        {
            eprintln!("Can't write to file {}: {}", path, e);
            return 1;
        }
        if !flush_file_or_report(&mut writer, path) {
            return 1;
        }
    }

    // Histogram output
    if let Some(ref path) = args.hist {
        let file = match File::create(path) {
            Ok(f) => f,
            Err(e) => {
                eprintln!("Can't open file {}: {}", path, e);
                return 1;
            }
        };
        let bins = hash_table.get_bins();
        let mut writer = BufWriter::new(file);
        if let Err(e) = kmer_output::write_histogram(&mut writer, &bins) {
            eprintln!("Can't write to file {}: {}", path, e);
            return 1;
        }
        if !flush_file_or_report(&mut writer, path) {
            return 1;
        }
    }

    // DBG output: write SKESA's hash graph file format.
    if let Some(ref path) = args.dbg_out {
        let mut sorted_kmers = crate::sorted_counter::count_kmers_sorted_with_output(
            reads_ref,
            args.kmer as usize,
            args.min_count as usize,
            32,
            &output,
        );
        crate::sorted_counter::get_branches_with_output(
            &mut sorted_kmers,
            args.kmer as usize,
            &output,
        );

        let file = match File::create(path) {
            Ok(f) => f,
            Err(e) => {
                eprintln!("Can't open file {}: {}", path, e);
                return 1;
            }
        };
        let mut writer = BufWriter::new(file);

        if let Err(e) = crate::hash_graph_output::write_hash_graph(
            &mut writer,
            reads_ref,
            &sorted_kmers,
            args.kmer as usize,
            !args.no_strand_info,
        ) {
            eprintln!("Can't write to file {}: {}", path, e);
            return 1;
        }
        if !flush_file_or_report(&mut writer, path) {
            return 1;
        }
    }

    eprintln!("DONE");
    0
}

pub fn run_cli_from<I, T>(args: I) -> i32
where
    I: IntoIterator<Item = T>,
    T: Into<std::ffi::OsString> + Clone,
{
    // Honor SKESA_RS_RLIMIT_GB before any allocation-heavy work, so a tight
    // cap fires early during dev/testing instead of after we've already
    // grabbed memory.
    crate::rlimit::apply_from_env();

    let cli = Cli::parse_from(args);

    // Set rayon thread pool based on --cores if available
    let cores = match &cli.command {
        Commands::Skesa(args) => args.cores,
        Commands::Kmercounter(args) => args.cores,
        Commands::Saute(args) => args.cores,
        Commands::SauteProt(args) => args.cores,
        Commands::GfaConnector(args) => args.cores,
        Commands::Stats(_) => 0,
    };
    if cores >= 0 {
        warn_if_cores_reduced(cores);
        rayon::ThreadPoolBuilder::new()
            .num_threads(resolve_cores(cores))
            .build_global()
            .ok(); // ignore if already initialized
    }

    let exit_code = match &cli.command {
        Commands::Skesa(args) => run_skesa(args),
        Commands::Kmercounter(args) => run_kmercounter(args),
        Commands::Saute(args) => run_saute(args),
        Commands::SauteProt(args) => run_saute_prot(args),
        Commands::GfaConnector(args) => run_gfa_connector(args),
        Commands::Stats(args) => run_stats(args),
    };

    exit_code
}