aboutsummaryrefslogtreecommitdiff
path: root/wqflask/wqflask/my_pylmm/pyLMM/input.py
blob: f7b556a548395900ae8ca8dbe1da7a4007cbfeda (about) (plain)
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
# pylmm is a python-based linear mixed-model solver with applications to GWAS

# Copyright (C) 2013  Nicholas A. Furlotte (nick.furlotte@gmail.com)

#The program is free for academic use. Please contact Nick Furlotte
#<nick.furlotte@gmail.com> if you are interested in using the software for
#commercial purposes.

#The software must not be modified and distributed without prior
#permission of the author.

#THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
#"AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
#LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
#A PARTICULAR PURPOSE ARE DISCLAIMED.  IN NO EVENT SHALL THE COPYRIGHT OWNER OR
#CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
#EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
#PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
#PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
#LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
#NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
#SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.

import os
import sys
import numpy as np
import struct
import pdb

class plink:
    def __init__(self,fbase,kFile=None,phenoFile=None,type='b',normGenotype=True,readKFile=False):
        self.fbase = fbase
        self.type = type
        self.indivs = self.getIndivs(self.fbase,type)
        self.kFile = kFile
        self.phenos = None
        self.normGenotype = normGenotype
        self.phenoFile = phenoFile
        # Originally I was using the fastLMM style that has indiv IDs embedded.
        # NOW I want to use this module to just read SNPs so I'm allowing 
        # the programmer to turn off the kinship reading.
        self.readKFile = readKFile

        if self.kFile:
            self.K = self.readKinship(self.kFile)
        elif os.path.isfile("%s.kin" % fbase): 
            self.kFile = "%s.kin" %fbase
            if self.readKFile:
                self.K = self.readKinship(self.kFile)
        else: 
            self.kFile = None
            self.K = None

        self.getPhenos(self.phenoFile)

        self.fhandle = None
        self.snpFileHandle = None

    def __del__(self): 
        if self.fhandle: self.fhandle.close()
        if self.snpFileHandle: self.snpFileHandle.close()

    def getSNPIterator(self):
        if not self.type == 'b': 
            sys.stderr.write("Have only implemented this for binary plink files (bed)\n")
            return

        # get the number of snps
        file = self.fbase + '.bim'
        i = 0
        f = open(file,'r')
        for line in f: i += 1
        f.close()
        self.numSNPs = i
        self.have_read = 0
        self.snpFileHandle = open(file,'r')

        self.BytestoRead = self.N / 4 + (self.N % 4 and 1 or 0)
        self._formatStr = 'c'*self.BytestoRead

        file = self.fbase + '.bed'
        self.fhandle = open(file,'rb')
     
        magicNumber = self.fhandle.read(2)
        order = self.fhandle.read(1)
        if not order == '\x01': 
            sys.stderr.write("This is not in SNP major order - you did not handle this case\n")
            raise StopIteration

        return self

    def __iter__(self):
        return self.getSNPIterator()

    def next(self):
        if self.have_read == self.numSNPs:
            raise StopIteration
        X = self.fhandle.read(self.BytestoRead)
        XX = [bin(ord(x)) for x in struct.unpack(self._formatStr,X)]
        self.have_read += 1
        return self.formatBinaryGenotypes(XX,self.normGenotype),self.snpFileHandle.readline().strip().split()[1]
    
    def formatBinaryGenotypes(self,X,norm=True):
        D = { \
              '00': 0.0, \
              '10': 0.5, \
              '11': 1.0, \
              '01': np.nan \
           }
  
        D_tped = { \
              '00': '1 1', \
              '10': '1 2', \
              '11': '2 2', \
              '01': '0 0' \
           }
  
        #D = D_tped
              
        G = []
        for x in X:
            if not len(x) == 10:
                xx = x[2:]
                x = '0b' + '0'*(8 - len(xx)) + xx
            a,b,c,d = (x[8:],x[6:8],x[4:6],x[2:4]) 
            L = [D[y] for y in [a,b,c,d]]
            G += L
        # only take the leading values because whatever is left should be null
        G = G[:self.N]
        G = np.array(G)
        if norm:
            G = self.normalizeGenotype(G)
        return G
    
    def normalizeGenotype(self,G):
        # print "Before",G
        # print G.shape
        x = True - np.isnan(G)
        m = G[x].mean()
        s = np.sqrt(G[x].var())
        G[np.isnan(G)] = m
        if s == 0: G = G - m
        else: G = (G - m) / s
        # print "After",G
        return G
    
    def getPhenos(self,phenoFile=None):
        if not phenoFile:
            self.phenoFile = phenoFile = self.fbase+".phenos"
        if not os.path.isfile(phenoFile): 
            sys.stderr.write("Could not find phenotype file: %s\n" % (phenoFile))
            return
        f = open(phenoFile,'r')
        keys = []
        P = []
        for line in f:
            v = line.strip().split()
            keys.append((v[0],v[1]))
            P.append([(x == 'NA' or x == '-9') and np.nan or float(x) for x in v[2:]])
        f.close()
        P = np.array(P)
     
        # reorder to match self.indivs
        D = {}
        L = []
        for i in range(len(keys)):
            D[keys[i]] = i
        for i in range(len(self.indivs)):
            if not D.has_key(self.indivs[i]):
                continue 
            L.append(D[self.indivs[i]])
        P = P[L,:]
     
        self.phenos = P
        return P
    
    def getIndivs(self,base,type='b'):
        if type == 't':
            famFile = "%s.tfam" % base
        else:
            famFile = "%s.fam" % base
        keys = []
        i = 0
        f = open(famFile,'r')
        for line in f:
           v = line.strip().split()
           famId = v[0]
           indivId = v[1]
           k = (famId.strip(),indivId.strip())
           keys.append(k)
           i += 1
        f.close()
     
        self.N = len(keys)
        sys.stderr.write("Read %d individuals from %s\n" % (self.N, famFile))
     
        return keys
    
    def readKinship(self,kFile):
        # Assume the fastLMM style
        # This will read in the kinship matrix and then reorder it
        # according to self.indivs - additionally throwing out individuals 
        # that are not in both sets
        if self.indivs == None or len(self.indivs) == 0:
            sys.stderr.write("Did not read any individuals so can't load kinship\n")
            return 
     
        sys.stderr.write("Reading kinship matrix from %s\n" % (kFile) )
     
        f = open(kFile,'r')
        # read indivs 
        v = f.readline().strip().split("\t")[1:]
        keys = [tuple(y.split()) for y in v]
        D = {}
        for i in range(len(keys)): D[keys[i]] = i
     
        # read matrix
        K = []
        for line in f:
            K.append([float(x) for x in line.strip().split("\t")[1:]])
        f.close()
        K  = np.array(K)
     
        # reorder to match self.indivs
        L = []
        KK = []
        X = []
        for i in range(len(self.indivs)):
            if not D.has_key(self.indivs[i]):
                X.append(self.indivs[i])
            else: 
                KK.append(self.indivs[i])
                L.append(D[self.indivs[i]])
        K = K[L,:][:,L]
        self.indivs = KK
        self.indivs_removed = X
        if len(self.indivs_removed):
            sys.stderr.write("Removed %d individuals that did not appear in Kinship\n" % (len(self.indivs_removed)))
        return K 
     
    def getCovariates(self,covFile=None):
        if not os.path.isfile(covFile): 
            sys.stderr.write("Could not find covariate file: %s\n" % (phenoFile))
            return
        f = open(covFile,'r')
        keys = []
        P = []
        for line in f:
            v = line.strip().split()
            keys.append((v[0],v[1]))
            P.append([x == 'NA' and np.nan or float(x) for x in v[2:]])
        f.close()
        P = np.array(P)
     
        # reorder to match self.indivs
        D = {}
        L = []
        for i in range(len(keys)):
            D[keys[i]] = i
        for i in range(len(self.indivs)):
            if not D.has_key(self.indivs[i]): continue 
            L.append(D[self.indivs[i]])
        P = P[L,:]
     
        return P