1 module dhtslib.memory; 2 3 import std.traits : isPointer, isSomeFunction, ReturnType, isSafe; 4 import core.stdc.stdlib : calloc, malloc, free; 5 import core.lifetime : move; 6 import std.typecons : RefCounted, RefCountedAutoInitialize; 7 import core.atomic : atomicOp; 8 import htslib; 9 10 /// can we use @live for scope checking? 11 enum dip1000Enabled = isSafe!((int x) => *&x); 12 13 static if(dip1000Enabled) 14 pragma(msg, "Using -dip1000 for scope checking and safety"); 15 16 pragma(inline, true): 17 /**! Logs an event with severity HTS_LOG_ERROR and compile-time context. */ 18 private void hts_log_errorNoGC(const(char)[] ctx)( string msg) @trusted @nogc nothrow 19 { 20 string open_error_color = "\x1b[0;31m"; 21 string close_color = "\x1b[0m"; 22 static newCtx = ctx ~ '\0'; 23 hts_log(htsLogLevel.HTS_LOG_ERROR, newCtx.ptr, "%.*s%.*s%.*s", 24 cast(int)open_error_color.length, open_error_color.ptr, 25 cast(int)msg.length, msg.ptr, 26 cast(int)close_color.length, close_color.ptr); 27 } 28 29 /// Template struct that wraps an htslib 30 /// pointer and reference counts it and then 31 /// destroys with destroyFun when it goes 32 /// truly out of scope 33 struct SafeHtslibPtr(T, alias destroyFun) 34 if(!isPointer!T && isSomeFunction!destroyFun) 35 { 36 @safe @nogc nothrow: 37 38 /// data pointer 39 T * ptr; 40 /// reference counting 41 shared int* refct; 42 /// initialized? 43 bool initialized; 44 45 /// ctor that respects scope 46 this(T * rawPtr) @trusted return scope 47 { 48 this.ptr = rawPtr; 49 this.refct = cast(shared int *) calloc(int.sizeof,1); 50 (*this.refct) = 1; 51 this.initialized = true; 52 } 53 54 /// postblit that respects scope 55 this(this) @trusted return scope 56 { 57 if(initialized)atomicOp!"+="(*this.refct, 1); 58 } 59 60 /// allow SafeHtslibPtr to be used as 61 /// underlying ptr type 62 alias getRef this; 63 64 /// get underlying data pointer 65 @property nothrow pure @nogc 66 ref inout(T*) getRef() inout return 67 { 68 return ptr; 69 } 70 71 /// take ownership of underlying data pointer 72 @property nothrow pure @nogc 73 T* moveRef() 74 { 75 T * ptr; 76 move(this.getRef, ptr); 77 return ptr; 78 } 79 80 /// dtor that respects scope 81 ~this() @trusted return scope 82 { 83 84 if(!this.initialized) return; 85 if(atomicOp!"-="(*this.refct, 1)) return; 86 if(this.ptr){ 87 free(cast(int*)this.refct); 88 /// if destroy function return is void 89 /// just destroy 90 /// else if return is int 91 /// destroy then check return value 92 /// else don't compile 93 static if(is(ReturnType!destroyFun == void)) 94 destroyFun(this.ptr); 95 else static if(is(ReturnType!destroyFun == int)) 96 { 97 auto err = destroyFun(this.ptr); 98 if(err != 0) 99 hts_log_errorNoGC!__FUNCTION__("Couldn't destroy/close "~T.stringof~" * data using function "~__traits(identifier, destroyFun)); 100 }else{ 101 static assert(0, "HtslibPtr doesn't recognize destroy function return type"); 102 } 103 } 104 } 105 } 106 107 /// reference counted bam1_t wrapper 108 /// can be used directly as a bam1_t * 109 alias Bam1 = SafeHtslibPtr!(bam1_t, bam_destroy1); 110 111 /// reference counted bam_hdr_t wrapper 112 /// can be used directly as a bam_hdr_t * 113 alias BamHdr = SafeHtslibPtr!(bam_hdr_t, bam_hdr_destroy); 114 115 /// reference counted bcf1_t wrapper 116 /// can be used directly as a bcf1_t * 117 alias Bcf1 = SafeHtslibPtr!(bcf1_t, bcf_destroy); 118 119 /// reference counted bcf_hdr_t wrapper 120 /// can be used directly as a bcf_hdr_t * 121 alias BcfHdr = SafeHtslibPtr!(bcf_hdr_t, bcf_hdr_destroy); 122 123 /// reference counted htsFile wrapper 124 /// can be used directly as a htsFile * 125 alias HtsFile = SafeHtslibPtr!(htsFile, hts_close); 126 127 /// reference counted htsFile wrapper 128 /// can be used directly as a htsFile * 129 alias HtsIdx = SafeHtslibPtr!(hts_idx_t, hts_idx_destroy); 130 131 /// reference counted htsFile wrapper 132 /// can be used directly as a htsFile * 133 alias HtsItr = SafeHtslibPtr!(hts_itr_t, hts_itr_destroy); 134 135 /// reference counted htsFile wrapper 136 /// can be used directly as a htsFile * 137 alias VcfFile = HtsFile; 138 139 /// reference counted tbx_t wrapper 140 /// can be used directly as a tbx_t * 141 alias Tbx = SafeHtslibPtr!(tbx_t, tbx_destroy); 142 143 /// reference counted BGZF wrapper 144 /// can be used directly as a BGZF * 145 alias Bgzf = SafeHtslibPtr!(BGZF, bgzf_close); 146 147 /// reference counted faidx_t wrapper 148 /// can be used directly as a faidx_t * 149 alias Faidx = SafeHtslibPtr!(faidx_t, fai_destroy); 150 151 /// reference counted Kstring wrapper 152 /// can be used directly as a kstring_t * 153 alias Kstring = SafeHtslibPtr!(kstring_t, ks_free); 154 155 alias HtsItrMulti = HtsItr; 156 157 debug(dhtslib_unittest) unittest 158 { 159 auto rc1 = Bam1(bam_init1); 160 assert(rc1.core.pos == 0); 161 // No more allocation, add just one extra reference count 162 auto rc2 = rc1; 163 // Reference semantics 164 rc2.core.pos = 42; 165 assert(rc1.core.pos == 42); 166 } 167 168 static if(dip1000Enabled){ 169 debug(dhtslib_unittest) unittest 170 { 171 auto rc1 = Bam1(bam_init1); 172 assert(rc1.core.pos == 0); 173 // No more allocation, add just one extra reference count 174 auto rc2 = rc1; 175 // Reference semantics 176 rc2.core.pos = 42; 177 assert(rc1.core.pos == 42); 178 } 179 180 debug(dhtslib_unittest) unittest 181 { 182 auto testfun(bool noScope = false)() 183 { 184 // auto rc = Bam1(bam_init); 185 auto rc1 = Bam1(bam_init1); 186 assert(rc1.core.pos == 0); 187 // No more allocation, add just one extra reference count 188 static if(noScope) 189 auto rc2 = rc1; 190 else 191 scope rc2 = rc1; 192 // Reference semantics 193 rc2.core.pos = 42; 194 assert(rc1.core.pos == 42); 195 return rc2.getRef; 196 } 197 static assert(!__traits(compiles,testfun!true)); 198 static assert(!__traits(compiles,testfun!false)); 199 } 200 201 debug(dhtslib_unittest) unittest 202 { 203 auto testfun(bool noScope = false)() 204 { 205 struct T 206 { 207 Bam1 rc1; 208 } 209 T test; 210 test.rc1 = Bam1(bam_init1); 211 212 assert(test.rc1.core.pos == 0); 213 // No more allocation, add just one extra reference count 214 static if(noScope) 215 auto rc2 = test.rc1; 216 else 217 scope rc2 = rc1; 218 // Reference semantics 219 rc2.core.pos = 42; 220 assert(test.rc1.core.pos == 42); 221 return rc2.getRef; 222 } 223 static assert(!__traits(compiles,testfun!true)); 224 static assert(!__traits(compiles,testfun!false)); 225 } 226 }