Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

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

# Copyright 2006-2013 by Peter Cock.  All rights reserved. 

# 

# This code is part of the Biopython distribution and governed by its 

# license.  Please see the LICENSE file that should have been included 

# as part of this package. 

"""Bio.AlignIO support for "stockholm" format (used in the PFAM database). 

 

You are expected to use this module via the Bio.AlignIO functions (or the 

Bio.SeqIO functions if you want to work directly with the gapped sequences). 

 

For example, consider a Stockholm alignment file containing the following:: 

 

    # STOCKHOLM 1.0 

    #=GC SS_cons       .................<<<<<<<<...<<<<<<<........>>>>>>>.. 

    AP001509.1         UUAAUCGAGCUCAACACUCUUCGUAUAUCCUC-UCAAUAUGG-GAUGAGGGU 

    #=GR AP001509.1 SS -----------------<<<<<<<<---..<<-<<-------->>->>..-- 

    AE007476.1         AAAAUUGAAUAUCGUUUUACUUGUUUAU-GUCGUGAAU-UGG-CACGA-CGU 

    #=GR AE007476.1 SS -----------------<<<<<<<<-----<<.<<-------->>.>>---- 

 

    #=GC SS_cons       ......<<<<<<<.......>>>>>>>..>>>>>>>>............... 

    AP001509.1         CUCUAC-AGGUA-CCGUAAA-UACCUAGCUACGAAAAGAAUGCAGUUAAUGU 

    #=GR AP001509.1 SS -------<<<<<--------->>>>>--->>>>>>>>--------------- 

    AE007476.1         UUCUACAAGGUG-CCGG-AA-CACCUAACAAUAAGUAAGUCAGCAGUGAGAU 

    #=GR AE007476.1 SS ------.<<<<<--------->>>>>.-->>>>>>>>--------------- 

    // 

 

This is a single multiple sequence alignment, so you would probably load this 

using the Bio.AlignIO.read() function: 

 

    >>> from Bio import AlignIO 

    >>> align = AlignIO.read("Stockholm/simple.sth", "stockholm") 

    >>> print(align) 

    SingleLetterAlphabet() alignment with 2 rows and 104 columns 

    UUAAUCGAGCUCAACACUCUUCGUAUAUCCUC-UCAAUAUGG-G...UGU AP001509.1 

    AAAAUUGAAUAUCGUUUUACUUGUUUAU-GUCGUGAAU-UGG-C...GAU AE007476.1 

    >>> for record in align: 

    ...     print("%s %i" % (record.id, len(record))) 

    AP001509.1 104 

    AE007476.1 104 

 

This example file is clearly using RNA, so you might want the alignment object 

(and the SeqRecord objects it holds) to reflect this, rather than simple using 

the default single letter alphabet as shown above.  You can do this with an 

optional argument to the Bio.AlignIO.read() function: 

 

    >>> from Bio import AlignIO 

    >>> from Bio.Alphabet import generic_rna 

    >>> align = AlignIO.read("Stockholm/simple.sth", "stockholm", 

    ...                      alphabet=generic_rna) 

    >>> print(align) 

    RNAAlphabet() alignment with 2 rows and 104 columns 

    UUAAUCGAGCUCAACACUCUUCGUAUAUCCUC-UCAAUAUGG-G...UGU AP001509.1 

    AAAAUUGAAUAUCGUUUUACUUGUUUAU-GUCGUGAAU-UGG-C...GAU AE007476.1 

 

In addition to the sequences themselves, this example alignment also includes 

some GR lines for the secondary structure of the sequences.  These are 

strings, with one character for each letter in the associated sequence: 

 

    >>> for record in align: 

    ...     print(record.id) 

    ...     print(record.seq) 

    ...     print(record.letter_annotations['secondary_structure']) 

    AP001509.1 

    UUAAUCGAGCUCAACACUCUUCGUAUAUCCUC-UCAAUAUGG-GAUGAGGGUCUCUAC-AGGUA-CCGUAAA-UACCUAGCUACGAAAAGAAUGCAGUUAAUGU 

    -----------------<<<<<<<<---..<<-<<-------->>->>..---------<<<<<--------->>>>>--->>>>>>>>--------------- 

    AE007476.1 

    AAAAUUGAAUAUCGUUUUACUUGUUUAU-GUCGUGAAU-UGG-CACGA-CGUUUCUACAAGGUG-CCGG-AA-CACCUAACAAUAAGUAAGUCAGCAGUGAGAU 

    -----------------<<<<<<<<-----<<.<<-------->>.>>----------.<<<<<--------->>>>>.-->>>>>>>>--------------- 

 

Any general annotation for each row is recorded in the SeqRecord's annotations 

dictionary.  You can output this alignment in many different file formats 

using Bio.AlignIO.write(), or the MultipleSeqAlignment object's format method: 

 

    >>> print(align.format("fasta")) 

    >AP001509.1 

    UUAAUCGAGCUCAACACUCUUCGUAUAUCCUC-UCAAUAUGG-GAUGAGGGUCUCUAC-A 

    GGUA-CCGUAAA-UACCUAGCUACGAAAAGAAUGCAGUUAAUGU 

    >AE007476.1 

    AAAAUUGAAUAUCGUUUUACUUGUUUAU-GUCGUGAAU-UGG-CACGA-CGUUUCUACAA 

    GGUG-CCGG-AA-CACCUAACAAUAAGUAAGUCAGCAGUGAGAU 

    <BLANKLINE> 

 

Most output formats won't be able to hold the annotation possible in a 

Stockholm file: 

 

    >>> print(align.format("stockholm")) 

    # STOCKHOLM 1.0 

    #=GF SQ 2 

    AP001509.1 UUAAUCGAGCUCAACACUCUUCGUAUAUCCUC-UCAAUAUGG-GAUGAGGGUCUCUAC-AGGUA-CCGUAAA-UACCUAGCUACGAAAAGAAUGCAGUUAAUGU 

    #=GS AP001509.1 AC AP001509.1 

    #=GS AP001509.1 DE AP001509.1 

    #=GR AP001509.1 SS -----------------<<<<<<<<---..<<-<<-------->>->>..---------<<<<<--------->>>>>--->>>>>>>>--------------- 

    AE007476.1 AAAAUUGAAUAUCGUUUUACUUGUUUAU-GUCGUGAAU-UGG-CACGA-CGUUUCUACAAGGUG-CCGG-AA-CACCUAACAAUAAGUAAGUCAGCAGUGAGAU 

    #=GS AE007476.1 AC AE007476.1 

    #=GS AE007476.1 DE AE007476.1 

    #=GR AE007476.1 SS -----------------<<<<<<<<-----<<.<<-------->>.>>----------.<<<<<--------->>>>>.-->>>>>>>>--------------- 

    // 

    <BLANKLINE> 

 

Note that when writing Stockholm files, AlignIO does not break long sequences 

up and interleave them (as in the input file shown above).  The standard 

allows this simpler layout, and it is more likely to be understood by other 

tools. 

 

Finally, as an aside, it can sometimes be useful to use Bio.SeqIO.parse() to 

iterate over the alignment rows as SeqRecord objects - rather than working 

with Alignnment objects. Again, if you want to you can specify this is RNA: 

 

    >>> from Bio import SeqIO 

    >>> from Bio.Alphabet import generic_rna 

    >>> for record in SeqIO.parse("Stockholm/simple.sth", "stockholm", 

    ...                           alphabet=generic_rna): 

    ...     print(record.id) 

    ...     print(record.seq) 

    ...     print(record.letter_annotations['secondary_structure']) 

    AP001509.1 

    UUAAUCGAGCUCAACACUCUUCGUAUAUCCUC-UCAAUAUGG-GAUGAGGGUCUCUAC-AGGUA-CCGUAAA-UACCUAGCUACGAAAAGAAUGCAGUUAAUGU 

    -----------------<<<<<<<<---..<<-<<-------->>->>..---------<<<<<--------->>>>>--->>>>>>>>--------------- 

    AE007476.1 

    AAAAUUGAAUAUCGUUUUACUUGUUUAU-GUCGUGAAU-UGG-CACGA-CGUUUCUACAAGGUG-CCGG-AA-CACCUAACAAUAAGUAAGUCAGCAGUGAGAU 

    -----------------<<<<<<<<-----<<.<<-------->>.>>----------.<<<<<--------->>>>>.-->>>>>>>>--------------- 

 

Remember that if you slice a SeqRecord, the per-letter-annotions like the 

secondary structure string here, are also sliced: 

 

    >>> sub_record = record[10:20] 

    >>> print(sub_record.seq) 

    AUCGUUUUAC 

    >>> print(sub_record.letter_annotations['secondary_structure']) 

    -------<<< 

""" 

from __future__ import print_function 

 

__docformat__ = "epytext en"  # not just plaintext 

from Bio.Seq import Seq 

from Bio.SeqRecord import SeqRecord 

from Bio.Align import MultipleSeqAlignment 

from .Interfaces import AlignmentIterator, SequentialAlignmentWriter 

 

 

class StockholmWriter(SequentialAlignmentWriter): 

    """Stockholm/PFAM alignment writer.""" 

 

    #These dictionaries should be kept in sync with those 

    #defined in the StockholmIterator class. 

    pfam_gr_mapping = {"secondary_structure": "SS", 

                       "surface_accessibility": "SA", 

                       "transmembrane": "TM", 

                       "posterior_probability": "PP", 

                       "ligand_binding": "LI", 

                       "active_site": "AS", 

                       "intron": "IN"} 

    #Following dictionary deliberately does not cover AC, DE or DR 

    pfam_gs_mapping = {"organism": "OS", 

                       "organism_classification": "OC", 

                       "look": "LO"} 

 

    def write_alignment(self, alignment): 

        """Use this to write (another) single alignment to an open file. 

 

        Note that sequences and their annotation are recorded 

        together (rather than having a block of annotation followed 

        by a block of aligned sequences). 

        """ 

        count = len(alignment) 

 

        self._length_of_sequences = alignment.get_alignment_length() 

        self._ids_written = [] 

 

        #NOTE - For now, the alignment object does not hold any per column 

        #or per alignment annotation - only per sequence. 

 

        if count == 0: 

            raise ValueError("Must have at least one sequence") 

        if self._length_of_sequences == 0: 

            raise ValueError("Non-empty sequences are required") 

 

        self.handle.write("# STOCKHOLM 1.0\n") 

        self.handle.write("#=GF SQ %i\n" % count) 

        for record in alignment: 

            self._write_record(record) 

        self.handle.write("//\n") 

 

    def _write_record(self, record): 

        """Write a single SeqRecord to the file""" 

        if self._length_of_sequences != len(record.seq): 

            raise ValueError("Sequences must all be the same length") 

 

        #For the case for stockholm to stockholm, try and use record.name 

        seq_name = record.id 

        if record.name is not None: 

            if "accession" in record.annotations: 

                if record.id == record.annotations["accession"]: 

                    seq_name = record.name 

 

        #In the Stockholm file format, spaces are not allowed in the id 

        seq_name = seq_name.replace(" ", "_") 

 

        if "start" in record.annotations \ 

        and "end" in record.annotations: 

            suffix = "/%s-%s" % (str(record.annotations["start"]), 

                                 str(record.annotations["end"])) 

            if seq_name[-len(suffix):] != suffix: 

                seq_name = "%s/%s-%s" % (seq_name, 

                                        str(record.annotations["start"]), 

                                        str(record.annotations["end"])) 

 

        if seq_name in self._ids_written: 

            raise ValueError("Duplicate record identifier: %s" % seq_name) 

        self._ids_written.append(seq_name) 

        self.handle.write("%s %s\n" % (seq_name, str(record.seq))) 

 

        #The recommended placement for GS lines (per sequence annotation) 

        #is above the alignment (as a header block) or just below the 

        #corresponding sequence. 

        # 

        #The recommended placement for GR lines (per sequence per column 

        #annotation such as secondary structure) is just below the 

        #corresponding sequence. 

        # 

        #We put both just below the corresponding sequence as this allows 

        #us to write the file using a single pass through the records. 

 

        #AC = Accession 

        if "accession" in record.annotations: 

            self.handle.write("#=GS %s AC %s\n" 

                % (seq_name, self.clean(record.annotations["accession"]))) 

        elif record.id: 

            self.handle.write("#=GS %s AC %s\n" 

                % (seq_name, self.clean(record.id))) 

 

        #DE = description 

        if record.description: 

            self.handle.write("#=GS %s DE %s\n" 

                % (seq_name, self.clean(record.description))) 

 

        #DE = database links 

        for xref in record.dbxrefs: 

            self.handle.write("#=GS %s DR %s\n" 

                % (seq_name, self.clean(xref))) 

 

        #GS = other per sequence annotation 

        for key, value in record.annotations.items(): 

            if key in self.pfam_gs_mapping: 

                data = self.clean(str(value)) 

                if data: 

                    self.handle.write("#=GS %s %s %s\n" 

                                      % (seq_name, 

                                         self.clean(self.pfam_gs_mapping[key]), 

                                         data)) 

            else: 

                #It doesn't follow the PFAM standards, but should we record 

                #this data anyway? 

                pass 

 

        #GR = per row per column sequence annotation 

        for key, value in record.letter_annotations.items(): 

            if key in self.pfam_gr_mapping and len(str(value)) == len(record.seq): 

                data = self.clean(str(value)) 

                if data: 

                    self.handle.write("#=GR %s %s %s\n" 

                                      % (seq_name, 

                                         self.clean(self.pfam_gr_mapping[key]), 

                                         data)) 

            else: 

                #It doesn't follow the PFAM standards, but should we record 

                #this data anyway? 

                pass 

 

 

class StockholmIterator(AlignmentIterator): 

    """Loads a Stockholm file from PFAM into MultipleSeqAlignment objects. 

 

    The file may contain multiple concatenated alignments, which are loaded 

    and returned incrementally. 

 

    This parser will detect if the Stockholm file follows the PFAM 

    conventions for sequence specific meta-data (lines starting #=GS 

    and #=GR) and populates the SeqRecord fields accordingly. 

 

    Any annotation which does not follow the PFAM conventions is currently 

    ignored. 

 

    If an accession is provided for an entry in the meta data, IT WILL NOT 

    be used as the record.id (it will be recorded in the record's 

    annotations).  This is because some files have (sub) sequences from 

    different parts of the same accession (differentiated by different 

    start-end positions). 

 

    Wrap-around alignments are not supported - each sequences must be on 

    a single line.  However, interlaced sequences should work. 

 

    For more information on the file format, please see: 

    http://www.bioperl.org/wiki/Stockholm_multiple_alignment_format 

    http://www.cgb.ki.se/cgb/groups/sonnhammer/Stockholm.html 

 

    For consistency with BioPerl and EMBOSS we call this the "stockholm" 

    format. 

    """ 

 

    #These dictionaries should be kept in sync with those 

    #defined in the PfamStockholmWriter class. 

    pfam_gr_mapping = {"SS": "secondary_structure", 

                       "SA": "surface_accessibility", 

                       "TM": "transmembrane", 

                       "PP": "posterior_probability", 

                       "LI": "ligand_binding", 

                       "AS": "active_site", 

                       "IN": "intron"} 

    #Following dictionary deliberately does not cover AC, DE or DR 

    pfam_gs_mapping = {"OS": "organism", 

                       "OC": "organism_classification", 

                       "LO": "look"} 

 

    def __next__(self): 

        try: 

            line = self._header 

            del self._header 

        except AttributeError: 

            line = self.handle.readline() 

        if not line: 

            #Empty file - just give up. 

            raise StopIteration 

        if not line.strip() == '# STOCKHOLM 1.0': 

            raise ValueError("Did not find STOCKHOLM header") 

 

        # Note: If this file follows the PFAM conventions, there should be 

        # a line containing the number of sequences, e.g. "#=GF SQ 67" 

        # We do not check for this - perhaps we should, and verify that 

        # if present it agrees with our parsing. 

 

        seqs = {} 

        ids = [] 

        gs = {} 

        gr = {} 

        gf = {} 

        passed_end_alignment = False 

        while True: 

            line = self.handle.readline() 

            if not line: 

                break  # end of file 

            line = line.strip()  # remove trailing \n 

            if line == '# STOCKHOLM 1.0': 

                self._header = line 

                break 

            elif line == "//": 

                #The "//" line indicates the end of the alignment. 

                #There may still be more meta-data 

                passed_end_alignment = True 

            elif line == "": 

                #blank line, ignore 

                pass 

            elif line[0] != "#": 

                #Sequence 

                #Format: "<seqname> <sequence>" 

                assert not passed_end_alignment 

                parts = [x.strip() for x in line.split(" ", 1)] 

                if len(parts) != 2: 

                    #This might be someone attempting to store a zero length sequence? 

                    raise ValueError("Could not split line into identifier " 

                                      + "and sequence:\n" + line) 

                id, seq = parts 

                if id not in ids: 

                    ids.append(id) 

                seqs.setdefault(id, '') 

                seqs[id] += seq.replace(".", "-") 

            elif len(line) >= 5: 

                #Comment line or meta-data 

                if line[:5] == "#=GF ": 

                    #Generic per-File annotation, free text 

                    #Format: #=GF <feature> <free text> 

                    feature, text = line[5:].strip().split(None, 1) 

                    #Each feature key could be used more than once, 

                    #so store the entries as a list of strings. 

                    if feature not in gf: 

                        gf[feature] = [text] 

                    else: 

                        gf[feature].append(text) 

                elif line[:5] == '#=GC ': 

                    #Generic per-Column annotation, exactly 1 char per column 

                    #Format: "#=GC <feature> <exactly 1 char per column>" 

                    pass 

                elif line[:5] == '#=GS ': 

                    #Generic per-Sequence annotation, free text 

                    #Format: "#=GS <seqname> <feature> <free text>" 

                    id, feature, text = line[5:].strip().split(None, 2) 

                    #if id not in ids: 

                    #    ids.append(id) 

                    if id not in gs: 

                        gs[id] = {} 

                    if feature not in gs[id]: 

                        gs[id][feature] = [text] 

                    else: 

                        gs[id][feature].append(text) 

                elif line[:5] == "#=GR ": 

                    #Generic per-Sequence AND per-Column markup 

                    #Format: "#=GR <seqname> <feature> <exactly 1 char per column>" 

                    id, feature, text = line[5:].strip().split(None, 2) 

                    #if id not in ids: 

                    #    ids.append(id) 

                    if id not in gr: 

                        gr[id] = {} 

                    if feature not in gr[id]: 

                        gr[id][feature] = "" 

                    gr[id][feature] += text.strip()  # append to any previous entry 

                    #TODO - Should we check the length matches the alignment length? 

                    #       For iterlaced sequences the GR data can be split over 

                    #       multiple lines 

            #Next line... 

 

        assert len(seqs) <= len(ids) 

        #assert len(gs)   <= len(ids) 

        #assert len(gr)   <= len(ids) 

 

        self.ids = ids 

        self.sequences = seqs 

        self.seq_annotation = gs 

        self.seq_col_annotation = gr 

 

        if ids and seqs: 

 

            if self.records_per_alignment is not None \ 

            and self.records_per_alignment != len(ids): 

                raise ValueError("Found %i records in this alignment, told to expect %i" 

                                 % (len(ids), self.records_per_alignment)) 

 

            alignment_length = len(list(seqs.values())[0]) 

            records = []  # Alignment obj will put them all in a list anyway 

            for id in ids: 

                seq = seqs[id] 

                if alignment_length != len(seq): 

                    raise ValueError("Sequences have different lengths, or repeated identifier") 

                name, start, end = self._identifier_split(id) 

                record = SeqRecord(Seq(seq, self.alphabet), 

                                   id=id, name=name, description=id, 

                                   annotations={"accession": name}) 

                #Accession will be overridden by _populate_meta_data if an explicit 

                #accession is provided: 

                record.annotations["accession"] = name 

 

                if start is not None: 

                    record.annotations["start"] = start 

                if end is not None: 

                    record.annotations["end"] = end 

 

                self._populate_meta_data(id, record) 

                records.append(record) 

            alignment = MultipleSeqAlignment(records, self.alphabet) 

 

            #TODO - Introduce an annotated alignment class? 

            #For now, store the annotation a new private property: 

            alignment._annotations = gr 

 

            return alignment 

        else: 

            raise StopIteration 

 

    def _identifier_split(self, identifier): 

        """Returns (name, start, end) string tuple from an identier.""" 

        if '/' in identifier: 

            name, start_end = identifier.rsplit("/", 1) 

            if start_end.count("-") == 1: 

                try: 

                    start, end = start_end.split("-") 

                    return name, int(start), int(end) 

                except ValueError: 

                    # Non-integers after final '/' - fall through 

                    pass 

        return identifier, None, None 

 

    def _get_meta_data(self, identifier, meta_dict): 

        """Takes an itentifier and returns dict of all meta-data matching it. 

 

        For example, given "Q9PN73_CAMJE/149-220" will return all matches to 

        this or "Q9PN73_CAMJE" which the identifier without its /start-end 

        suffix. 

 

        In the example below, the suffix is required to match the AC, but must 

        be removed to match the OS and OC meta-data:: 

 

            # STOCKHOLM 1.0 

            #=GS Q9PN73_CAMJE/149-220  AC Q9PN73 

            ... 

            Q9PN73_CAMJE/149-220               NKA... 

            ... 

            #=GS Q9PN73_CAMJE OS Campylobacter jejuni 

            #=GS Q9PN73_CAMJE OC Bacteria 

 

        This function will return an empty dictionary if no data is found.""" 

        name, start, end = self._identifier_split(identifier) 

        if name == identifier: 

            identifier_keys = [identifier] 

        else: 

            identifier_keys = [identifier, name] 

        answer = {} 

        for identifier_key in identifier_keys: 

            try: 

                for feature_key in meta_dict[identifier_key]: 

                    answer[feature_key] = meta_dict[identifier_key][feature_key] 

            except KeyError: 

                pass 

        return answer 

 

    def _populate_meta_data(self, identifier, record): 

        """Adds meta-date to a SecRecord's annotations dictionary. 

 

        This function applies the PFAM conventions.""" 

 

        seq_data = self._get_meta_data(identifier, self.seq_annotation) 

        for feature in seq_data: 

            #Note this dictionary contains lists! 

            if feature == "AC":  # ACcession number 

                assert len(seq_data[feature]) == 1 

                record.annotations["accession"] = seq_data[feature][0] 

            elif feature == "DE":  # DEscription 

                record.description = "\n".join(seq_data[feature]) 

            elif feature == "DR":  # Database Reference 

                #Should we try and parse the strings? 

                record.dbxrefs = seq_data[feature] 

            elif feature in self.pfam_gs_mapping: 

                record.annotations[self.pfam_gs_mapping[feature]] = ", ".join(seq_data[feature]) 

            else: 

                #Ignore it? 

                record.annotations["GS:" + feature] = ", ".join(seq_data[feature]) 

 

        #Now record the per-letter-annotations 

        seq_col_data = self._get_meta_data(identifier, self.seq_col_annotation) 

        for feature in seq_col_data: 

            #Note this dictionary contains strings! 

            if feature in self.pfam_gr_mapping: 

                record.letter_annotations[self.pfam_gr_mapping[feature]] = seq_col_data[feature] 

            else: 

                #Ignore it? 

                record.letter_annotations["GR:" + feature] = seq_col_data[feature] 

 

 

if __name__ == "__main__": 

    from Bio._utils import run_doctest 

    run_doctest()