From 10c8b56ecbdc1494d318dc0515ff4983e6f53758 Mon Sep 17 00:00:00 2001 From: Zilong-Li Date: Tue, 27 Aug 2024 09:42:14 +0200 Subject: [PATCH] copyHeader with specific samples --- tools/Makefile | 2 ++ tools/vcf_addDS.cpp | 4 ++-- tools/vcf_addINFO.cpp | 2 +- vcfpp.h | 5 +++-- 4 files changed, 8 insertions(+), 5 deletions(-) diff --git a/tools/Makefile b/tools/Makefile index d8d1ed1..980e204 100644 --- a/tools/Makefile +++ b/tools/Makefile @@ -21,4 +21,6 @@ clean: test: ./vcf_addDS -i test.vcf.gz -o out.vcf.gz + ./vcf_addDS -i test.vcf.gz -s Ind3,Ind1 -o out.vcf.gz ./vcf_addINFO -i test.vcf.gz -o out.vcf.gz + ./vcf_addINFO -i test.vcf.gz -s Ind0,Ind2 -o out.vcf.gz diff --git a/tools/vcf_addDS.cpp b/tools/vcf_addDS.cpp index 2e4774c..28bddf6 100644 --- a/tools/vcf_addDS.cpp +++ b/tools/vcf_addDS.cpp @@ -35,8 +35,8 @@ int main(int argc, char * argv[]) } // ========= core calculation part =========================================== BcfReader vcf(invcf, region, samples); - BcfWriter bw(outvcf, vcf.header); - bw.copyHeader(invcf); + BcfWriter bw(outvcf); + bw.copyHeader(invcf, samples); bw.header.addFORMAT("DS", "1", "Float", "Diploid Genotype Dosage"); // add DS tag into the header int nsamples = vcf.nsamples; vector gps, ds(nsamples); diff --git a/tools/vcf_addINFO.cpp b/tools/vcf_addINFO.cpp index ef5cbb3..aca7b26 100644 --- a/tools/vcf_addINFO.cpp +++ b/tools/vcf_addINFO.cpp @@ -37,7 +37,7 @@ int main(int argc, char * argv[]) // ========= core calculation part =========================================== BcfReader vcf(invcf, region, samples); BcfWriter bw(outvcf); - bw.copyHeader(invcf); + bw.copyHeader(invcf, samples); bw.header.addINFO("INFO", "1", "Float", "INFO score given genotype probability"); bw.header.addINFO("EAF", "1", "Float", "Estimated Allele Frequency"); int N = vcf.nsamples, i{0}; diff --git a/vcfpp.h b/vcfpp.h index a6a14fc..b70b997 100644 --- a/vcfpp.h +++ b/vcfpp.h @@ -1750,12 +1750,13 @@ class BcfWriter hp = &h; } - /// copy header of given VCF - void copyHeader(const std::string & vcffile) + /// copy header of given VCF and restrict on samples + void copyHeader(const std::string & vcffile, std::string samples = "-") { htsFile * fp2 = hts_open(vcffile.c_str(), "r"); if(!fp2) throw std::invalid_argument("I/O error: input file is invalid"); header.hdr = bcf_hdr_read(fp2); + header.setSamples(samples); hts_close(fp2); initalHeader(header); }