source: tags/ms_r16q2/RNA3D/RNA3D_StructureData.cxx

Last change on this file was 12803, checked in by westram, 10 years ago
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 50.3 KB
Line 
1#include "RNA3D_GlobalHeader.hxx"
2#include "RNA3D_Global.hxx"
3#include "RNA3D_OpenGLGraphics.hxx"
4#include "RNA3D_Graphics.hxx"
5#include "RNA3D_StructureData.hxx"
6
7#include <cerrno>
8#include <string>
9#include <iostream>
10#include <fstream>
11
12#include <arbdbt.h>
13#include <aw_msg.hxx>
14#include <aw_awars.hxx>
15#include <aw_root.hxx>
16#include <ed4_extern.hxx>
17#include <ed4_plugins.hxx>
18#include <BI_basepos.hxx>
19
20
21#define COLORLINK (ED4_G_SBACK_0 - RNA3D_GC_SBACK_0)  // to link to the colors defined in primary editor ed4_defs.hxx
22#define SAICOLORS (ED4_G_CBACK_0 - RNA3D_GC_CBACK_0)
23
24using namespace std;
25
26static Struct3Dinfo   *start3D                 = NULL;
27static Struct2Dinfo   *start2D                 = NULL;
28static Struct2Dplus3D *start2D3D               = NULL;
29static HelixNrInfo    *start                   = NULL;
30static CurrSpecies    *startSp                 = NULL;
31static bool            bOldSpeciesDataExists   = false;
32static Insertions     *startIns                = NULL;
33static bool            bOldInsertionDataExists = false;
34static bool            bStartPosStored         = false;
35static bool            bEndPosStored           = false;
36
37static char *find_data_file(const char *name) {
38    char *fname = GB_lib_file(false, "rna3d/", name);
39    if (!fname) throw string("file not found: ")+name;
40    return fname;
41}
42
43static void throw_IO_error(const char *filename) __ATTR__NORETURN;
44static void throw_IO_error(const char *filename) {
45    string error = string("IO-Error: ")+strerror(errno)+" ('"+filename+"')";
46    throw error;
47}
48
49GBDATA *Structure3D::gb_main = 0;
50
51Structure3D::Structure3D(ED4_plugin_host& host_)
52    : strCen(new Vector3(0.0, 0.0, 0.0)),
53      iInterval(25),
54      iMapSAI(0),
55      iMapSearch(0),
56      iMapEnable(0),
57      iStartPos(0),
58      iEndPos(0),
59      iEColiStartPos(0),
60      iEColiEndPos(0),
61      iTotalSubs(0),
62      iTotalDels(0),
63      iTotalIns(0),
64      LSU_molID(3), // default molecule to load in case of 23S rRNA : 1C2W
65      HelixBase(500),
66      EColiRef(NULL),
67      Host(host_),
68      GRAPHICS(new OpenGLGraphics)
69{
70}
71
72Structure3D::~Structure3D() {
73    delete GRAPHICS;
74    delete strCen;
75}
76
77
78void Structure3D::StoreCoordinates(float x, float y, float z, char base, unsigned int pos) {
79    Struct3Dinfo *data, *temp;
80    data = new Struct3Dinfo;
81    data->x    = x;
82    data->y    = y;
83    data->z    = z;
84    data->base = base;
85    data->pos  = pos;
86    data->next = NULL;
87
88    if (start3D == NULL)
89        start3D = data;
90    else {
91        temp = start3D;
92        // We know this is not NULL - list not empty!
93        while (temp->next != NULL) {
94            temp = temp->next;  // Move to next link in chain
95        }
96        temp->next = data;
97    }
98}
99
100// ------------------Selecting rRNA type--------------------------//
101
102int Structure3D::FindTypeOfRNA() {
103    int rnaType = 0;
104    GB_push_transaction(gb_main);  // opening a transaction
105
106    GBDATA *gbTemplate = GBT_find_SAI(gb_main, "ECOLI");
107    if (!gbTemplate) {
108        aw_message("SAI:ECOLI not found");
109    }
110    else {
111        char   *ali_name          = GBT_get_default_alignment(gb_main);
112        GBDATA *gbAlignment       = GB_entry(gbTemplate, ali_name);
113        GBDATA *gbTemplateSeqData = gbAlignment ? GB_entry(gbAlignment, "data") : 0;
114
115        if (!gbTemplateSeqData) {
116          aw_message("Could not find species in the database!");
117        }
118        else {
119            const char *pTemplateSeqData  = GB_read_char_pntr(gbTemplateSeqData);
120            int iSeqLen = strlen(pTemplateSeqData);
121            int iBaseCount = 0;
122            for (int i = 0; i<iSeqLen; i++) {
123                if ((pTemplateSeqData[i] != '.') && (pTemplateSeqData[i] != '-')) {
124                    iBaseCount++;
125                }
126            }
127
128            if (iBaseCount > 2000)      rnaType = LSU_23S;
129            else if (iBaseCount > 300)  rnaType = SSU_16S;
130            else                        rnaType = LSU_5S;
131        }
132    }
133    GB_pop_transaction(gb_main);
134
135    return rnaType;
136}
137
138// ----------Delete old molecule data-------------------------------//
139// Struct2Dinfo, Struct2Dplus3D, Struct3Dinfo, HelixNrInfo, Insertions
140
141void Structure3D::DeleteOldMoleculeData() {
142    // Struct2Dinfo -> start2D
143    {
144        Struct2Dinfo *tmp, *data;
145        for (data = start2D; data != NULL; data = tmp) {
146            tmp = data->next;
147            delete data;
148        }
149        start2D = NULL;
150    }
151
152    // Struct2Dplus3D ->start2D3D
153    {
154        Struct2Dplus3D *tmp, *data;
155        for (data = start2D3D; data != NULL; data = tmp) {
156            tmp = data->next;
157            delete data;
158        }
159        start2D3D = NULL;
160    }
161
162    // Struct3Dinfo ->start3D
163    {
164        Struct3Dinfo *tmp, *data;
165        for (data = start3D; data != NULL; data = tmp) {
166            tmp = data->next;
167            delete data;
168        }
169        start3D = NULL;
170    }
171
172    // HelixNrInfo -> start
173    {
174        HelixNrInfo *tmp, *data;
175        for (data = start; data != NULL; data = tmp) {
176            tmp = data->next;
177            delete data;
178        }
179        start = NULL;
180    }
181
182    // Delete insertions data
183    DeleteOldInsertionData();
184
185    // Delete mapped species data
186    DeleteOldSpeciesData ();
187}
188
189// =========== Reading 3D Coordinates from PDB file ====================//
190
191void Structure3D::ReadCoOrdinateFile() {
192    static char *DataFile = 0;
193    int          rnaType  = FindTypeOfRNA();
194
195    RNA3D->bDisplayComments = true; // displaying comments in main window
196
197    switch (rnaType) {
198    case LSU_23S:
199        switch (LSU_molID) {
200        case _1PNU:
201            DataFile = find_data_file("Ecoli_1PNU_23S_rRNA.pdb");
202            sprintf(globalComment, "The 3D molecule rendered from PDB entry : 1PNU at 8.7 Angstrom.");
203            break;
204        case _1VOR:
205            DataFile = find_data_file("Ecoli_1VOR_23S_rRNA.pdb");
206            sprintf(globalComment, "The 3D molecule is rendered from PDB entry [1VOR] with 11.5 Angstrom resolution.");
207            break;
208        case _1C2W:
209            DataFile = find_data_file("Ecoli_1C2W_23S_rRNA.pdb");
210            sprintf(globalComment, "The 3D molecule is rendered from PDB entry [1C2W] with 7.5 Angstrom resolution.");
211            break;
212        }
213        break;
214    case LSU_5S:
215        DataFile = find_data_file("Ecoli_1C2X_5S_rRNA.pdb");
216        sprintf(globalComment, "The 3D molecule is rendered from PDB entry [1C2X] with 7.5 Angstrom resolution.");
217        break;
218    case SSU_16S:
219        DataFile = find_data_file("Ecoli_1M5G_16S_rRNA.pdb");
220        sprintf(globalComment, "The 3D molecule is rendered from PDB entry [1M5G] with 11.5 Angstrom resolution.");
221        break;
222    }
223
224    ifstream readData;
225    readData.open(DataFile, ios::in);
226    if (!readData.is_open()) {
227        throw_IO_error(DataFile);
228    }
229
230    int          cntr                 = 0;
231    unsigned int last3Dpos            = 0;
232    static bool  bEColiStartPosStored = false;
233
234    while (!readData.eof()) {
235        char buf[256];
236        readData.getline(buf, 100);
237
238        string line(buf);
239        if ((line.find("ATOM") != string::npos) || (line.find("HETATM") != string::npos)) {
240            string atom = line.substr(77, 2);
241            if (atom.find("P") != string::npos) {
242                char     base = line[19];
243                unsigned pos  = atoi(line.substr(22, 4).c_str());
244
245                // special filter for 23S rRNA structure (IVOR/IPNU)
246                // IVOR/IPNU contains artifacts and are not mentioned in any of the
247                // remarks of PDB file
248                if (last3Dpos != pos && !(pos >= 3093)) {
249                    float X = atof(line.substr(31, 8).c_str());
250                    float Y = atof(line.substr(39, 8).c_str());
251                    float Z = atof(line.substr(47, 8).c_str());
252
253                    StoreCoordinates(X, Y, Z, base, pos);
254                    last3Dpos  = pos;
255                    strCen->x += X; strCen->y += Y; strCen->z += Z;
256                    cntr++;
257                }
258                if (!bEColiStartPosStored) {
259                    iEColiStartPos = pos;
260                    bEColiStartPosStored = true;
261                }
262                iEColiEndPos = pos;
263            }
264        }
265    }
266
267    cout<<"--------------------------------------------------"<<endl
268        <<globalComment<<endl
269        <<"Structure contains "<<iEColiEndPos-iEColiStartPos<<"( "<<iEColiStartPos<<" - "
270        <<iEColiEndPos<<" ) positions."<<endl
271        <<"---------------------------------------------------"<<endl;
272
273    strCen->x = strCen->x/cntr;
274    strCen->y = strCen->y/cntr;
275    strCen->z = strCen->z/cntr;
276
277    readData.close();
278    free(DataFile);
279}
280
281
282void Structure3D::Store2Dinfo(char *info, int pos, int helixNr) {
283    Struct2Dinfo *data = new Struct2Dinfo;
284
285    data->base    = info[0];
286    data->mask    = info[1];
287    data->code    = info[2];
288    data->pos     = pos;
289    data->helixNr = helixNr;
290
291    data->next = NULL;
292
293    if (start2D == NULL) {
294        start2D = data;
295    }
296    else {
297        Struct2Dinfo *temp = start2D;
298        // We know this is not NULL - list not empty!
299        while (temp->next != NULL) {
300            temp = temp->next;  // Move to next link in chain
301        }
302        temp->next = data;
303    }
304}
305
306// =========== Reading Secondary Structure Data from Ecoli Secondary Structure Mask file ====================//
307
308void Structure3D::GetSecondaryStructureInfo() {
309    static char *DataFile = 0;
310    int rnaType    = FindTypeOfRNA();
311
312    switch (rnaType) {
313        case LSU_23S:
314            DataFile = find_data_file("SecondaryStructureModel_23SrRNA.data");
315            break;
316        case LSU_5S:
317            DataFile = find_data_file("SecondaryStructureModel_5SrRNA.data");
318            break;
319        case SSU_16S:
320            DataFile = find_data_file("SecondaryStructureModel_16SrRNA.data");
321            break;
322    }
323
324    char  buf[256];
325
326    int  pos         = 0;
327    int  helixNr     = 0;
328    int  lastHelixNr = 0;
329    char info[4];
330    info[3]          = '\0';
331    bool insideHelix = false;
332    bool skip        = false;
333
334    ifstream readData;
335    readData.open(DataFile, ios::in);
336    if (!readData.is_open()) {
337        throw_IO_error(DataFile);
338    }
339
340    while (!readData.eof()) {
341        readData.getline(buf, 100);
342        char *tmp;
343        tmp = strtok(buf, " ");
344        for (int i = 0; tmp != NULL; tmp = strtok(NULL, " "), i++) {
345            switch (i) {
346                case 0: pos = atoi(tmp);      break;
347                case 1: info[0] = tmp[0];     break;
348                case 2: info[1] = tmp[0];     break;
349                case 3: info[2] = tmp[0];     break;
350                case 4: helixNr = atoi(tmp); break;
351            }
352        }
353
354        bool is_SE = (info[2] == 'S') || (info[2] == 'E'); // 'S' = helix start 'E' = helix end
355        if (is_SE) {
356            if (helixNr>0) lastHelixNr = helixNr;
357
358            if (!insideHelix) insideHelix = true;
359            else {
360                Store2Dinfo(info, pos, lastHelixNr);
361                skip        = true;
362                insideHelix = false;
363            }
364        }
365        if (insideHelix) Store2Dinfo(info, pos, lastHelixNr);
366        else {
367            if (skip) skip = false;
368            else Store2Dinfo(info, pos, 0);
369        }
370    }
371    readData.close();
372    free(DataFile);
373}
374
375void Structure3D::Store2D3Dinfo(Struct2Dinfo *s2D, Struct3Dinfo *s3D) {
376    Struct2Dplus3D *data, *temp;
377    data = new Struct2Dplus3D;
378    data->base    = s2D->base;
379    data->mask    = s2D->mask;
380    data->code    = s2D->code;
381    data->pos     = s2D->pos;
382    data->helixNr = s2D->helixNr;
383    data->x       = s3D->x;
384    data->y       = s3D->y;
385    data->z       = s3D->z;
386    data->next    = NULL;
387
388    if (start2D3D == NULL)
389        start2D3D = data;
390    else {
391        temp = start2D3D;
392        while (temp->next != NULL) {
393            temp = temp->next;
394        }
395        temp->next = data;
396    }
397}
398
399// =========== Combining Secondary Structure Data with 3D Coordinates =======================//
400
401void Structure3D::Combine2Dand3DstructureInfo() {
402    Struct3Dinfo *temp3D;
403    Struct2Dinfo *temp2D;
404    int cntr = 0;
405
406    cout<<"Missing Base Positions : "<<endl;
407
408    temp3D = start3D;
409    temp2D = start2D;
410    while ((temp3D != NULL) &&  (temp2D != NULL)) {
411        if (temp2D->pos == temp3D->pos) {
412            Store2D3Dinfo(temp2D, temp3D);
413        }
414        else {
415            while (temp2D->pos != temp3D->pos) {
416                cout<<temp2D->pos<<", "; // missing base positions
417                cntr++;
418                temp2D = temp2D->next;
419            }
420            Store2D3Dinfo(temp2D, temp3D);
421        }
422
423        temp3D = temp3D->next;
424        temp2D = temp2D->next;
425    }
426    cout<<endl<<"Total No. of bases missing = "<<cntr<<endl;
427
428    // printing comments in the main window
429    {
430        RNA3D->bDisplayComments = true;
431        char buf[256];
432        sprintf(buf, "Total No. of bases missing = %d. See the console messages for the actual missing base positions.", cntr);
433        strcat(globalComment, buf);
434    }
435}
436
437void Structure3D::PointsToQuads(float x, float y, float z) {
438
439    if (RNA3D->bPointSpritesSupported) {
440        glVertex3f(x, y, z);
441    }
442    else {
443        glBegin(GL_QUADS);
444        // Front Face
445        glTexCoord2f(0, 0); glVertex3f(x - 1, y + 1, z + 1);
446        glTexCoord2f(1, 0); glVertex3f(x + 1, y + 1, z + 1);
447        glTexCoord2f(1, 1); glVertex3f(x + 1, y - 1, z + 1);
448        glTexCoord2f(0, 1); glVertex3f(x - 1, y - 1, z + 1);
449
450        // Back Face
451        glTexCoord2f(0, 0); glVertex3f(x + 1, y + 1, z - 1);
452        glTexCoord2f(1, 0); glVertex3f(x - 1, y + 1, z - 1);
453        glTexCoord2f(1, 1); glVertex3f(x - 1, y - 1, z - 1);
454        glTexCoord2f(0, 1); glVertex3f(x + 1, y - 1, z - 1);
455
456        // Top Face
457        glTexCoord2f(0, 0); glVertex3f(x + 1, y + 1, z + 1);
458        glTexCoord2f(1, 0); glVertex3f(x - 1, y + 1, z + 1);
459        glTexCoord2f(1, 1); glVertex3f(x - 1, y + 1, z - 1);
460        glTexCoord2f(0, 1); glVertex3f(x + 1, y + 1, z - 1);
461
462        // Bottom Face
463        glTexCoord2f(0, 0); glVertex3f(x + 1, y - 1, z - 1);
464        glTexCoord2f(1, 0); glVertex3f(x - 1, y - 1, z - 1);
465        glTexCoord2f(1, 1); glVertex3f(x - 1, y - 1, z + 1);
466        glTexCoord2f(0, 1); glVertex3f(x + 1, y - 1, z + 1);
467
468        // Left Face
469        glTexCoord2f(0, 0); glVertex3f(x + 1, y + 1, z + 1);
470        glTexCoord2f(1, 0); glVertex3f(x + 1, y + 1, z - 1);
471        glTexCoord2f(1, 1); glVertex3f(x + 1, y - 1, z - 1);
472        glTexCoord2f(0, 1); glVertex3f(x + 1, y - 1, z + 1);
473
474        // Right Face
475        glTexCoord2f(0, 0); glVertex3f(x - 1, y + 1, z - 1);
476        glTexCoord2f(1, 0); glVertex3f(x - 1, y + 1, z + 1);
477        glTexCoord2f(1, 1); glVertex3f(x - 1, y - 1, z + 1);
478        glTexCoord2f(0, 1); glVertex3f(x - 1, y - 1, z - 1);
479
480        glEnd();
481    }
482}
483
484void Structure3D::PositionsToCoordinatesDispList(int listID, int *pos, int len) {
485    Struct2Dplus3D *t;
486    int tmpPos = 0;
487
488    glNewList(listID, GL_COMPILE);
489    {
490        if (RNA3D->bPointSpritesSupported) {
491            glBegin(GL_POINTS);
492        }
493        for (int i = 0; i < len; i++)
494            {
495                tmpPos = pos[i];
496                t = start2D3D;
497                while (t != NULL) {
498                    if (t->pos == tmpPos) {
499                        PointsToQuads(t->x, t->y, t->z);
500                        break;
501                    }
502                    t = t->next;
503                }
504            }
505        if (RNA3D->bPointSpritesSupported) {
506            glEnd();
507        }
508    }
509    glEndList();
510}
511
512void Structure3D::GenerateSecStructureNonHelixRegions() {
513    Struct2Dplus3D *t;
514    const int MAX_BASE = 1000;
515    int baseA[MAX_BASE], baseG[MAX_BASE], baseC[MAX_BASE], baseU[MAX_BASE];
516    int a, g, c, u; a=g=c=u=0;
517
518    {
519        t = start2D3D;
520        while (t != NULL) {
521            if (t->helixNr == 0) {
522                switch (t->base)
523                    {
524                    case 'A': baseA[a++] = t->pos; break;
525                    case 'G': baseG[g++] = t->pos; break;
526                    case 'C': baseC[c++] = t->pos; break;
527                    case 'U': baseU[u++] = t->pos; break;
528                    }
529            }
530            t = t->next;
531        }
532    }
533
534    PositionsToCoordinatesDispList(NON_HELIX_A, baseA, a);
535    PositionsToCoordinatesDispList(NON_HELIX_G, baseG, g);
536    PositionsToCoordinatesDispList(NON_HELIX_C, baseC, c);
537    PositionsToCoordinatesDispList(NON_HELIX_U, baseU, u);
538}
539
540void Structure3D::GenerateSecStructureHelixRegions() {
541    Struct2Dplus3D *t;
542    const int MAX_BASE = 1000;
543    int baseA[MAX_BASE], baseG[MAX_BASE], baseC[MAX_BASE], baseU[MAX_BASE];
544    int a, g, c, u; a=g=c=u=0;
545
546    {
547        t = start2D3D;
548        while (t != NULL) {
549            if (t->helixNr > 0) {
550                if ((t->mask == '[') || (t->mask == ']') ||
551                    (t->mask == '<') || (t->mask == '>'))
552                    {
553                        switch (t->base)
554                            {
555                            case 'A': baseA[a++] = t->pos; break;
556                            case 'G': baseG[g++] = t->pos; break;
557                            case 'C': baseC[c++] = t->pos; break;
558                            case 'U': baseU[u++] = t->pos; break;
559                            }
560                    }
561            }
562            t = t->next;
563        }
564    }
565
566    PositionsToCoordinatesDispList(HELIX_A, baseA, a);
567    PositionsToCoordinatesDispList(HELIX_G, baseG, g);
568    PositionsToCoordinatesDispList(HELIX_C, baseC, c);
569    PositionsToCoordinatesDispList(HELIX_U, baseU, u);
570}
571
572void Structure3D::GenerateSecStructureUnpairedHelixRegions() {
573    Struct2Dplus3D *t;
574    const int MAX_BASE = 500;
575    int baseA[MAX_BASE], baseG[MAX_BASE], baseC[MAX_BASE], baseU[MAX_BASE];
576    int a, g, c, u; a=g=c=u=0;
577
578    {
579        t = start2D3D;
580        while (t != NULL) {
581            if (t->helixNr > 0) {
582                if (t->mask == '.')
583                    {
584                        switch (t->base)
585                            {
586                            case 'A': baseA[a++] = t->pos; break;
587                            case 'G': baseG[g++] = t->pos; break;
588                            case 'C': baseC[c++] = t->pos; break;
589                            case 'U': baseU[u++] = t->pos; break;
590                            }
591                    }
592            }
593            t = t->next;
594        }
595    }
596
597    PositionsToCoordinatesDispList(UNPAIRED_HELIX_A, baseA, a);
598    PositionsToCoordinatesDispList(UNPAIRED_HELIX_G, baseG, g);
599    PositionsToCoordinatesDispList(UNPAIRED_HELIX_C, baseC, c);
600    PositionsToCoordinatesDispList(UNPAIRED_HELIX_U, baseU, u);
601}
602
603// ==============================================================================
604// Tertiary Interactions of 16S ribosomal RNA model of E.coli.
605// Reference : http://www.rna.icmb.utexas.edu/
606// Year of Last Update : 2001.
607// Pseudoknots and Triple Base pairs are extracted and displayed in
608// the 3D model.
609// ==============================================================================
610
611void Structure3D::GenerateTertiaryInteractionsDispLists() {
612    Struct2Dplus3D *t;
613    static char *DataFile = 0;
614    int rnaType    = FindTypeOfRNA();
615
616    switch (rnaType) {
617    case LSU_23S:
618        DataFile = find_data_file("TertiaryInteractions_23SrRNA.data");
619        break;
620    case LSU_5S:
621        cout<<"There are no tertiary interactions observed in 5S rRNA model!"<<endl;
622        return;
623    case SSU_16S:
624        DataFile = find_data_file("TertiaryInteractions_16SrRNA.data");
625        break;
626    }
627
628    cout<<"Tertiary Interactions are fetched from comparative RNA website [http://www.rna.icmb.utexas.edu]."<<endl
629        <<"The same are located in the file \""<<DataFile<<"\"."<<endl;
630
631    char  buf[256];
632
633    ifstream readData;
634    readData.open(DataFile, ios::in);
635    if (!readData.is_open()) {
636        throw_IO_error(DataFile);
637    }
638
639    int K[100];
640    int R[100];
641    int k, r; k = r = 0;
642
643    while (!readData.eof()) {
644        readData.getline(buf, 100);
645        char *tmp;
646        tmp = strtok(buf, " ");
647        if (tmp != NULL)
648            {
649                if (strcmp(tmp, "KNOT") == 0) {
650                    tmp = strtok(NULL, ":");
651                    while (tmp != NULL) {
652                        K[k++] = atoi(tmp);
653                        tmp = strtok(NULL, ":");
654                    }
655                }
656                else if (strcmp(tmp, "TRIPLE") == 0) {
657                    tmp = strtok(NULL, ":");
658                    while (tmp != NULL) {
659                        R[r++] = atoi(tmp);
660                        tmp = strtok(NULL, ":");
661                    }
662                }
663            }
664    }
665    readData.close();
666
667    glNewList(ECOLI_TERTIARY_INTRACTION_PSEUDOKNOTS, GL_COMPILE);
668    {
669        for (int i = 0; i < k;) {
670            glBegin(GL_LINES);
671            for (int j = 0; j < 2; j++)
672                {
673                    t = start2D3D;
674                    while (t != NULL) {
675                        if (t->pos == K[i]) {
676                            glVertex3f(t->x, t->y, t->z);
677                            i++;
678                            break;
679                        }
680                        t = t->next;
681                    }
682                }
683            glEnd();
684        }
685    }
686    glEndList();
687
688    glNewList(ECOLI_TERTIARY_INTRACTION_TRIPLE_BASES, GL_COMPILE);
689    {
690        for (int i = 0; i < r;) {
691            glBegin(GL_LINE_STRIP);
692            for (int j = 0; j < 3; j++)
693                {
694                    t = start2D3D;
695                    while (t != NULL) {
696                        if (t->pos == R[i]) {
697                            glVertex3f(t->x, t->y, t->z);
698                            i++;
699                            break;
700                        }
701                        t = t->next;
702                    }
703                }
704            glEnd();
705        }
706    }
707    glEndList();
708    free(DataFile);
709}
710
711// ==============================================================================
712
713void Structure3D::StoreHelixNrInfo(float x, float y, float z, int helixNr) {
714    HelixNrInfo *data, *temp;
715    data = new HelixNrInfo;
716    data->helixNr = helixNr;
717    data->x       = x;
718    data->y       = y;
719    data->z       = z;
720    data->next    = NULL;
721
722    if (start == NULL)
723        start = data;
724    else {
725        temp = start;
726        while (temp->next != NULL) {
727            temp = temp->next;
728        }
729        temp->next = data;
730    }
731}
732
733void Structure3D::GenerateHelixDispLists(int HELIX_NR_ID, int HELIX_NR) {
734    Struct2Dinfo *temp2D;
735    Struct2Dplus3D *temp2D3D;
736
737    const int MAX = 500;
738
739    int thisStrandPos[MAX], otherStrandPos[MAX];
740    int i, j; i = j = 0;
741    {
742        temp2D = start2D;
743        while (temp2D != NULL) {
744            if (temp2D->helixNr == HELIX_NR) {
745                if ((temp2D->mask == '[') || (temp2D->mask == '<'))
746                    thisStrandPos[i++]  = temp2D->pos;
747                if ((temp2D->mask == ']') || (temp2D->mask == '>'))
748                    otherStrandPos[j++] = temp2D->pos;
749            }
750            temp2D = temp2D->next;
751        }
752    }
753
754    int tempPos = 0;
755    float x1, x2, y1, y2, z1, z2; x1=x2=y1=y2=z1=z2=0.0;
756
757    bool bThisStrand, bOtherStrand;
758    bThisStrand = bOtherStrand = false;
759
760    bool bMissingBasePair = false;
761
762    glNewList(HELIX_NR_ID, GL_COMPILE);
763    {
764        for (int k = 0, l = j-1; k < i && l >= 0; k++, l--)
765            {
766                tempPos = thisStrandPos[k];
767                {
768                    temp2D3D = start2D3D;
769                    while (temp2D3D != NULL) {
770                        if (temp2D3D->pos == tempPos) {
771                            bThisStrand = true;
772                            x1=temp2D3D->x; y1=temp2D3D->y; z1=temp2D3D->z;
773                            break;
774                        }
775                        temp2D3D = temp2D3D->next;
776                    }
777                }
778                tempPos = otherStrandPos[l];
779                {
780                    temp2D3D = start2D3D;
781                    while (temp2D3D != NULL) {
782                        if (temp2D3D->pos == tempPos) {
783                            bOtherStrand = true;
784                            x2=temp2D3D->x; y2=temp2D3D->y; z2=temp2D3D->z;
785                            break;
786                        }
787                        temp2D3D = temp2D3D->next;
788                    }
789                }
790                // if bases present in both the strands then draw a bond
791                // and store the helix number information
792                if (bThisStrand && bOtherStrand) {
793                    glVertex3f(x1, y1, z1);
794                    glVertex3f(x2, y2, z2);
795
796                    x1 = (x1+x2)/2; y1 = (y1+y2)/2; z1 = (z1+z2)/2;
797                    StoreHelixNrInfo(x1, y1, z1, HELIX_NR);
798
799                    bThisStrand = bOtherStrand = false;
800                }
801                else {
802                    bMissingBasePair = true;
803                    cout<<thisStrandPos[k]<<"-"<<tempPos<<"("<<HELIX_NR_ID<<"), ";
804                }
805            }
806        if (bMissingBasePair) cout<<endl;
807    }
808    glEndList();
809}
810
811void Structure3D::GenerateHelixNrDispList(int startHx, int endHx) {
812    HelixNrInfo *t;
813
814    glNewList(HELIX_NUMBERS, GL_COMPILE);
815    {
816        char POS[50];
817        for (int i = startHx; i <= endHx; i++) {
818            t = start;
819            while (t != NULL) {
820                if (t->helixNr == i) {
821                    sprintf(POS, "%d", t->helixNr);
822                    GRAPHICS->PrintString(t->x, t->y, t->z, POS, GLUT_BITMAP_8_BY_13);
823                }
824                t = t->next;
825            }
826        }
827    }
828    glEndList();
829
830    glNewList(HELIX_NUMBERS_POINTS, GL_COMPILE);
831    {
832        if (RNA3D->bPointSpritesSupported) {
833            glBegin(GL_POINTS);
834        }
835        for (int i = startHx; i <= endHx; i++) {
836            t = start;
837            while (t != NULL) {
838                if (t->helixNr == i) {
839                    PointsToQuads(t->x, t->y, t->z);
840                }
841                t = t->next;
842            }
843        }
844        if (RNA3D->bPointSpritesSupported) {
845            glEnd();
846        }
847    }
848    glEndList();
849}
850
851void Structure3D::GenerateDisplayLists() {
852
853    GenerateMoleculeSkeleton();
854    ComputeBasePositions();
855
856    int rnaType = FindTypeOfRNA();
857    switch (rnaType) {
858    case SSU_16S:
859        if (glIsList(HelixBase) != GL_TRUE)
860            {
861                for (int i = 1; i <= 50; i++) {
862                    GenerateHelixDispLists(HelixBase+i, i);
863                }
864            }
865        break;
866    case LSU_23S:
867        cout<<"=========================================================="<<endl
868            <<"Missing base pairs (bracket number indicates helix number) "<<endl;
869        if (glIsList(HelixBase) != GL_TRUE)
870            {
871                for (int i = 1; i <= 101; i++) {
872                    GenerateHelixDispLists(HelixBase+i, i);
873                }
874            }
875        cout<<"=========================================================="<<endl;
876        break;
877    case LSU_5S:
878        if (glIsList(HelixBase) != GL_TRUE)
879            {
880                for (int i = 1; i <= 5; i++) {
881                    GenerateHelixDispLists(HelixBase+i, i);
882                }
883            }
884        break;
885    }
886
887    GenerateSecStructureHelixRegions();
888    GenerateSecStructureNonHelixRegions();
889    GenerateSecStructureUnpairedHelixRegions();
890
891    GenerateTertiaryInteractionsDispLists();
892}
893
894void Structure3D::GenerateMoleculeSkeleton() {
895    Struct2Dplus3D *t;
896
897    glNewList(STRUCTURE_BACKBONE, GL_COMPILE);
898    {
899        glBegin(GL_LINE_STRIP);
900        t = start2D3D;
901        while (t != NULL) {
902            glVertex3f(t->x, t->y, t->z);
903            t = t->next;
904        }
905        glEnd();
906    }
907    glEndList();
908
909    glNewList(STRUCTURE_BACKBONE_CLR, GL_COMPILE);
910    {
911        glBegin(GL_LINE_STRIP);
912        t = start2D3D;
913        while (t != NULL)
914            {
915                if (t->helixNr > 0) {
916                    if ((t->mask == '[') || (t->mask == ']') ||
917                        (t->mask == '<') || (t->mask == '>'))
918                        {
919                            GRAPHICS->SetColor(RNA3D_GC_BASES_HELIX);
920                            glVertex3f(t->x, t->y, t->z);
921                        }
922                    if (t->mask == '.') {
923                        GRAPHICS->SetColor(RNA3D_GC_BASES_UNPAIRED_HELIX);
924                        glVertex3f(t->x, t->y, t->z);
925                    }
926                }
927                 if (t->helixNr == 0) {
928                     GRAPHICS->SetColor(RNA3D_GC_BASES_NON_HELIX);
929                     glVertex3f(t->x, t->y, t->z);
930                 }
931                 t = t->next;
932            }
933        glEnd();
934    }
935    glEndList();
936}
937
938void Structure3D::GenerateCursorPositionDispList(long pos) {
939    Struct3Dinfo *temp;
940
941    glNewList(ECOLI_CURSOR_POSITION, GL_COMPILE);
942    {
943        glBegin(GL_POINTS);
944        temp = start3D;
945        while (temp != NULL) {
946            if (temp->pos == pos) {
947#ifdef DEBUG
948                cout<<"Cursor Position : "<<pos<<endl;
949#endif
950                glVertex3f(temp->x, temp->y, temp->z);
951                break;
952            }
953            temp = temp->next;
954        }
955        glEnd();
956    }
957    glEndList();
958}
959
960void Structure3D::ComputeBasePositions() {
961    Struct3Dinfo *temp;
962
963    char POS[50];
964    float spacer = 1.5;
965    int posSkip = iInterval;
966    if (posSkip <= 0) posSkip = 25; // default interval value
967
968    glNewList(STRUCTURE_POS, GL_COMPILE);
969    {
970        temp = start3D;
971        while (temp != NULL) {
972            if (temp->pos%posSkip == 0) {
973                sprintf(POS, "%d", temp->pos);
974                GRAPHICS->PrintString(temp->x-spacer, temp->y, temp->z-spacer, POS, GLUT_BITMAP_HELVETICA_10);
975            }
976            temp = temp->next;
977        }
978    }
979    glEndList();
980
981    glNewList(STRUCTURE_POS_ANCHOR, GL_COMPILE);
982    {
983        glLineWidth(0.5);
984        glBegin(GL_LINES);
985        temp = start3D;
986        while (temp != NULL) {
987            if (temp->pos%posSkip == 0) {
988                glVertex3f(temp->x, temp->y, temp->z);
989                glVertex3f(temp->x-spacer, temp->y, temp->z-spacer);
990            }
991            temp = temp->next;
992        }
993        glEnd();
994    }
995    glEndList();
996}
997
998void Structure3D::PrepareSecondaryStructureData() {
999//     const char
1000//         outFile[]      = "/nfshome/yadhu/DataBase/Struct3D/model16SrRNA/test.data",
1001//         EcoliFile[]    = "/nfshome/yadhu/DataBase/Struct3D/model16SrRNA/ECOLI_GAPS",
1002//         HelixNrFile[]  = "/nfshome/yadhu/DataBase/Struct3D/model16SrRNA/HELIX_NR",
1003//         HelixGapFile[] = "/nfshome/yadhu/DataBase/Struct3D/model16SrRNA/HELIX_GAP",
1004//         ErrorMsg[] = "\n *** Error Opening File : ";
1005
1006    const char
1007        outFile[]      = "/home/Yadhu/DataBase/Struct3D/model23SrRNA/test.data",
1008        EcoliFile[]    = "/home/Yadhu/DataBase/Struct3D/model23SrRNA/Ecoli_23SrRNA_with_gaps.seq",
1009        HelixNrFile[]  = "/home/Yadhu/DataBase/Struct3D/model23SrRNA/HelixNumbers.msk",
1010        HelixGapFile[] = "/home/Yadhu/DataBase/Struct3D/model23SrRNA/SecStructureMask_with_gaps.msk",
1011        ErrorMsg[] = "\n *** Error Opening File : ";
1012
1013//     const char
1014//         outFile[]      = "/home/Yadhu/DataBase/Struct3D/model5SrRNA/test.data",
1015//         EcoliFile[]    = "/home/Yadhu/DataBase/Struct3D/model5SrRNA/Ecoli_5SrRNA_with_gaps.seq",
1016//         HelixNrFile[]  = "/home/Yadhu/DataBase/Struct3D/model5SrRNA/HelixNumbers.msk",
1017//         HelixGapFile[] = "/home/Yadhu/DataBase/Struct3D/model5SrRNA/SecStructureMask_with_gaps.msk",
1018//         ErrorMsg[] = "\n *** Error Opening File : ";
1019
1020    int fileLen = 0;
1021    char *helixNrBuf, *helixGapBuf, *ecoliBuf;
1022
1023    ofstream out;
1024    out.open(outFile, ios::out);
1025    if (!out.is_open())  cerr<<ErrorMsg<<outFile<<endl;
1026
1027    ifstream inFile;
1028    {
1029        {
1030            inFile.open(EcoliFile, ios::binary);
1031            if (!inFile.is_open())  cerr<<ErrorMsg<<HelixGapFile<<endl;
1032
1033            inFile.seekg (0, ios::end);  // get length of file
1034            fileLen = inFile.tellg();
1035            inFile.seekg (0, ios::beg);
1036
1037            ecoliBuf = new char[fileLen];    // allocate memory:
1038
1039            inFile.read (ecoliBuf, fileLen);    // read data as a block:
1040            inFile.close();
1041        }
1042
1043        {
1044            inFile.open(HelixNrFile, ios::binary);
1045            if (!inFile.is_open())  cerr<<ErrorMsg<<HelixGapFile<<endl;
1046
1047            inFile.seekg (0, ios::end);  // get length of file
1048            fileLen = inFile.tellg();
1049            inFile.seekg (0, ios::beg);
1050
1051            helixNrBuf = new char[fileLen];    // allocate memory:
1052
1053            inFile.read (helixNrBuf, fileLen);    // read data as a block:
1054            inFile.close();
1055        }
1056
1057        {
1058            inFile.open(HelixGapFile, ios::binary);
1059            if (!inFile.is_open())  cerr<<ErrorMsg<<HelixNrFile<<endl;
1060
1061            inFile.seekg (0, ios::end);  // get length of file
1062            fileLen = inFile.tellg();
1063            inFile.seekg (0, ios::beg);
1064
1065            helixGapBuf = new char[fileLen];    // allocate memory:
1066
1067            inFile.read (helixGapBuf, fileLen);    // read data as a block:
1068            inFile.close();
1069        }
1070
1071        char helixNr[4]; helixNr[3] = '\0';
1072        int skip, gaps, k;
1073        skip = gaps = k = 0;
1074
1075        for (unsigned int i = 0; i < strlen(ecoliBuf); i++)
1076            {
1077                if (ecoliBuf[i] == '\n')    skip++;
1078                else {
1079                    if ((ecoliBuf[i] == '.') || (ecoliBuf[i] == '-'))  gaps++;
1080                    else
1081                        {
1082                            int pos = (i - (skip + gaps)) + 1;
1083                            out<<pos<<"  "<<ecoliBuf[i]<<"  "<<helixGapBuf[i]<<"  ";
1084                            switch (helixGapBuf[i])
1085                                {
1086                                case '.': out<<"N  "; break;
1087                                case '[': out<<"S  "; break;
1088                                case ']': out<<"E  "; break;
1089                                case '<': out<<"H  "; break;
1090                                case '>': out<<"H  "; break;
1091                                }
1092                            for (unsigned int j = i; helixNrBuf[j] != '.'; j++) {
1093                                helixNr[k++] = helixNrBuf[j];
1094                            }
1095                            for (int l = k; l < 4; l++)  helixNr[l] = '\0';
1096                            k = 0;
1097                            if ((helixGapBuf[i] == '[') || (helixGapBuf[i] == ']')) {
1098                                out<<helixNr;
1099                            }
1100                            out<<endl;
1101                        }
1102                }
1103            }
1104        delete [] helixNrBuf;
1105        delete [] ecoliBuf;
1106        delete [] helixGapBuf;
1107        out.close();
1108    }
1109}
1110
1111void Structure3D::StoreCurrSpeciesDifference(char base, int pos) {
1112    Struct2Dplus3D *st;
1113    st = start2D3D;
1114    while (st != NULL) {
1115        if (st->pos == pos) {
1116            CurrSpecies *data, *temp;
1117            data = new CurrSpecies;
1118            data->base = base;
1119            data->pos  = pos;
1120            data->x = st->x;
1121            data->y = st->y;
1122            data->z = st->z;
1123            data->next = NULL;
1124
1125            if (startSp == NULL) {
1126                startSp = data;
1127                bOldSpeciesDataExists = true;
1128            }
1129            else {
1130                temp = startSp;
1131                while (temp->next != NULL) {
1132                    temp = temp->next;
1133                }
1134                temp->next = data;
1135            }
1136            break;
1137        }
1138        st = st->next;
1139    }
1140}
1141
1142void Structure3D::DeleteOldSpeciesData() {
1143    CurrSpecies *tmp, *data;
1144
1145    for (data = startSp; data != NULL; data = tmp) {
1146        tmp = data->next;
1147        delete data;
1148    }
1149    startSp = NULL;
1150
1151    bOldSpeciesDataExists = false;
1152}
1153
1154void Structure3D::GenerateBaseDifferencePositionDisplayList() {
1155    CurrSpecies *t;
1156    Struct3Dinfo *temp;
1157
1158    char POS[50];
1159    float spacer = 1.5;
1160
1161    glNewList(MAP_SPECIES_BASE_DIFFERENCE_POS, GL_COMPILE);
1162    {
1163        for (t = startSp; t != NULL; t = t->next)
1164            {
1165                for (temp = start3D; temp != NULL; temp = temp->next)
1166                    {
1167                        if (temp->pos == t->pos) {
1168                            sprintf(POS, "%d", temp->pos);
1169                            GRAPHICS->PrintString(temp->x-spacer, temp->y, temp->z-spacer, POS, GLUT_BITMAP_8_BY_13);
1170                        }
1171                    }
1172            }
1173    }
1174    glEndList();
1175
1176    glNewList(MAP_SPECIES_BASE_DIFFERENCE_POS_ANCHOR, GL_COMPILE);
1177    {
1178        glLineWidth(1.0);
1179        glBegin(GL_LINES);
1180
1181        for (t = startSp; t != NULL; t = t->next)
1182            {
1183                for (temp = start3D; temp != NULL; temp = temp->next)
1184                    {
1185                        if (temp->pos == t->pos) {
1186                            glVertex3f(temp->x, temp->y, temp->z);
1187                            glVertex3f(temp->x-spacer, temp->y, temp->z-spacer);
1188                        }
1189                    }
1190            }
1191        glEnd();
1192    }
1193    glEndList();
1194}
1195
1196void Structure3D::BuildDisplayList(int listID, int *pos, int len) {
1197    CurrSpecies *t;
1198    int tmpPos = 0;
1199
1200    glNewList(listID, GL_COMPILE);
1201    {
1202        if (RNA3D->bPointSpritesSupported) {
1203            glBegin(GL_POINTS);
1204        }
1205        for (int i = 0; i < len; i++)
1206            {
1207                tmpPos = pos[i];
1208                t = startSp;
1209                while (t != NULL) {
1210                    if (t->pos == tmpPos) {
1211                        PointsToQuads(t->x, t->y, t->z);
1212                        break;
1213                    }
1214                    t = t->next;
1215                }
1216            }
1217        if (RNA3D->bPointSpritesSupported) {
1218            glEnd();
1219        }
1220    }
1221    glEndList();
1222}
1223
1224void Structure3D::GenerateBaseDifferenceDisplayList() {
1225    CurrSpecies *t;
1226
1227    glNewList(MAP_SPECIES_BASE_DIFFERENCE, GL_COMPILE);
1228    {
1229        if (RNA3D->bPointSpritesSupported) {
1230            glBegin(GL_POINTS);
1231        }
1232        t = startSp;
1233        while (t != NULL) {
1234            PointsToQuads(t->x, t->y, t->z);
1235            t = t->next;
1236        }
1237        if (RNA3D->bPointSpritesSupported) {
1238            glEnd();
1239        }
1240    }
1241    glEndList();
1242
1243    const int MAX_BASE = 400;
1244    int baseA[MAX_BASE], baseG[MAX_BASE], baseC[MAX_BASE],
1245        baseU[MAX_BASE], deletion[MAX_BASE], miss[MAX_BASE];
1246    int a, g, c, u, d, m; a=g=c=u=d=m=0;
1247
1248    {
1249        t = startSp;
1250        while (t != NULL) {
1251            switch (t->base)
1252                {
1253                case 'A': baseA[a++]     = t->pos; break;
1254                case 'G': baseG[g++]     = t->pos; break;
1255                case 'C': baseC[c++]     = t->pos; break;
1256                case 'U': baseU[u++]     = t->pos; break;
1257                case '-': deletion[d++] = t->pos; break;
1258                case '.': miss[m++] = t->pos; break;
1259                }
1260            t = t->next;
1261        }
1262        BuildDisplayList(MAP_SPECIES_BASE_A, baseA, a);
1263        BuildDisplayList(MAP_SPECIES_BASE_U, baseU, u);
1264        BuildDisplayList(MAP_SPECIES_BASE_G, baseG, g);
1265        BuildDisplayList(MAP_SPECIES_BASE_C, baseC, c);
1266        BuildDisplayList(MAP_SPECIES_DELETION, deletion, d);
1267        BuildDisplayList(MAP_SPECIES_MISSING, miss, m);
1268        GenerateBaseDifferencePositionDisplayList();
1269
1270        iTotalSubs = a+g+c+u;  // Summing up the substitutions/mutations
1271        iTotalDels = d;        // Storing total number of deletions
1272#ifdef DEBUG
1273        cout<<"Substitutions = "<<iTotalSubs<<endl;
1274        cout<<"Deletions     = "<<iTotalDels<<endl;
1275#endif
1276    }
1277}
1278
1279void Structure3D::StoreInsertions(char base, int pos) {
1280    Insertions *data, *temp;
1281    data = new Insertions;
1282    data->base = base;
1283    data->pos  = pos;
1284    data->next = NULL;
1285
1286    if (startIns == NULL) {
1287        startIns = data;
1288        bOldInsertionDataExists = true;
1289    }
1290    else {
1291        temp = startIns;
1292        while (temp->next != NULL) {
1293            temp = temp->next;
1294        }
1295        temp->next = data;
1296    }
1297}
1298
1299void Structure3D::DeleteOldInsertionData() {
1300    Insertions *tmp, *data;
1301
1302    for (data = startIns; data != NULL; data = tmp) {
1303        tmp = data->next;
1304        delete data;
1305    }
1306    startIns = NULL;
1307
1308    bOldInsertionDataExists = false;
1309}
1310
1311void Structure3D::GenerateInsertionDisplayList() {
1312    Insertions   *ins;
1313    Struct3Dinfo *str;
1314    char inserts[500];
1315    int i, cntr;
1316    float spacer = 2.0;
1317    iTotalIns = 0;
1318
1319    glNewList(MAP_SPECIES_INSERTION_BASES, GL_COMPILE);
1320    {
1321        for (str = start3D; str != NULL; str = str->next)
1322        {
1323            i = cntr = 0;
1324            for (ins = startIns; ins != NULL; ins = ins->next)
1325                {
1326                    if (str->pos == ins->pos) {
1327                        inserts[i++] = ins->base;
1328                        cntr++;
1329                    }
1330                }
1331            if (cntr>0) {
1332                inserts[i] = '\0';
1333                char buffer[strlen(inserts) + 10];
1334                sprintf(buffer, "%d:%s", cntr, inserts);
1335                GRAPHICS->PrintString(str->x, str->y+spacer, str->z,
1336                                      buffer, GLUT_BITMAP_8_BY_13);
1337
1338                iTotalIns += cntr; // Summing up the insertions
1339            }
1340        }
1341    }
1342    glEndList();
1343
1344#ifdef DEBUG
1345    cout<<"Insertions = "<<iTotalIns<<endl;
1346#endif
1347
1348    glNewList(MAP_SPECIES_INSERTION_BASES_ANCHOR, GL_COMPILE);
1349    {
1350        glLineWidth(1.0);
1351        glBegin(GL_LINES);
1352
1353        for (str = start3D; str != NULL; str = str->next)
1354        {
1355            for (ins = startIns; ins != NULL; ins = ins->next)
1356                {
1357                    if (str->pos == ins->pos) {
1358                        glVertex3f(str->x, str->y, str->z);
1359                        glVertex3f(str->x, str->y+spacer, str->z);
1360                    }
1361                }
1362        }
1363        glEnd();
1364    }
1365    glEndList();
1366
1367    glNewList(MAP_SPECIES_INSERTION_POINTS, GL_COMPILE);
1368    {
1369        if (RNA3D->bPointSpritesSupported) {
1370            glBegin(GL_POINTS);
1371        }
1372        for (str = start3D; str != NULL; str = str->next)
1373        {
1374            for (ins = startIns; ins != NULL; ins = ins->next)
1375                {
1376                    if (str->pos == ins->pos) {
1377                        PointsToQuads(str->x, str->y, str->z);
1378                        break;
1379                    }
1380                }
1381        }
1382        if (RNA3D->bPointSpritesSupported) {
1383            glEnd();
1384        }
1385    }
1386    glEndList();
1387}
1388
1389void Structure3D::MapCurrentSpeciesToEcoliTemplate(AW_root *awr) {
1390
1391    GB_push_transaction(gb_main);
1392
1393    GBDATA *gbTemplate = GBT_find_SAI(gb_main, "ECOLI");
1394
1395    if (!gbTemplate) {
1396        aw_message("SAI:ECOLI not found");
1397    }
1398    else {
1399        char *pSpeciesName = awr->awar(AWAR_SPECIES_NAME)->read_string();
1400        if (pSpeciesName) {
1401            Host.announce_current_species(pSpeciesName);
1402            GBDATA *gbSpecies = GBT_find_species(gb_main, pSpeciesName);
1403            if (gbSpecies) {
1404                char   *ali_name          = GBT_get_default_alignment(gb_main);
1405                GBDATA *gbAlignment       = GB_entry(gbTemplate, ali_name);
1406                GBDATA *gbTemplateSeqData = gbAlignment ? GB_entry(gbAlignment, "data") : 0;
1407
1408                if (!gbTemplateSeqData) {
1409                    aw_message(GBS_global_string("Mapping impossible, since SAI:ECOLI has no data in alignment '%s'", ali_name));
1410                }
1411                else {
1412                    const char *pTemplateSeqData  = GB_read_char_pntr(gbTemplateSeqData);
1413
1414                    if (!RNA3D->bEColiRefInitialized) {
1415                        EColiRef = new BI_ecoli_ref;
1416                        EColiRef->init(gb_main);
1417                        RNA3D->bEColiRefInitialized = true;
1418                    }
1419
1420                    char buf[100];
1421                    char *pSpFullName = GB_read_string(GB_entry(gbSpecies, "full_name"));
1422                    sprintf(buf, "%s : %s", pSpeciesName, pSpFullName);
1423                    awr->awar(AWAR_3D_SELECTED_SPECIES)->write_string(buf);
1424
1425                    GBDATA *gbSeqData    = GBT_find_sequence(gbSpecies, ali_name);
1426                    const char *pSeqData = GB_read_char_pntr(gbSeqData);
1427
1428                    if (pSeqData && pTemplateSeqData) {
1429                        int iSeqLen = strlen(pTemplateSeqData);
1430
1431                        if (bOldSpeciesDataExists) {
1432                            DeleteOldSpeciesData();
1433                        }
1434
1435                        for (int i = 0, iSeqCount = 0; i<iSeqLen; i++) {
1436                            if ((pTemplateSeqData[i] != '.') && (pTemplateSeqData[i] != '-'))
1437                            {
1438                                if (!bStartPosStored) {
1439                                    iStartPos = i;
1440                                    bStartPosStored = true;
1441                                }
1442                                if (pTemplateSeqData[i] != pSeqData[i]) {
1443                                    StoreCurrSpeciesDifference(pSeqData[i], iSeqCount);
1444                                }
1445                                iSeqCount++;
1446                            }
1447                        }
1448
1449                        for (int i = iSeqLen; i>0; i--) {
1450                            if ((pTemplateSeqData[i] != '.') && (pTemplateSeqData[i] != '-')) {
1451                                if (!bEndPosStored) {
1452                                    iEndPos = i;
1453                                    bEndPosStored = true;
1454                                    break;
1455                                }
1456                            }
1457                        }
1458
1459                        if (bOldInsertionDataExists) {
1460                            DeleteOldInsertionData();
1461                        }
1462
1463                        for (int i = iStartPos, iSeqCount = 0; i < iEndPos; i++) {
1464                            if ((pTemplateSeqData[i] != '.') &&
1465                               (pTemplateSeqData[i] != '-'))
1466                                { // Store EColi base positions : Insertion point !
1467                                    iSeqCount++;
1468                                }
1469
1470                            if ((pTemplateSeqData[i] == '-') &&
1471                               (pSeqData[i]         != '-') &&
1472                               (pSeqData[i]         != '.'))
1473                                {
1474                                    StoreInsertions(pSeqData[i], iSeqCount);
1475                                }
1476                        }
1477                    }
1478                    free(pSpFullName);
1479                }
1480                free(ali_name);
1481            }
1482            GenerateBaseDifferenceDisplayList();
1483            GenerateInsertionDisplayList();
1484        }
1485        free(pSpeciesName);
1486    }
1487    GB_pop_transaction(gb_main);
1488}
1489
1490static bool ValidSearchColor(int iColor, int mode) {
1491    bool isValid = false;
1492
1493    switch (mode) {
1494        case SAI:    isValid = (iColor >= RNA3D_GC_CBACK_0) && (iColor < RNA3D_GC_MAX); break;
1495        case SEARCH: isValid = (iColor >= RNA3D_GC_SBACK_0) && (iColor < RNA3D_GC_MAX); break;
1496    }
1497    return isValid;
1498}
1499
1500void Structure3D::MapSearchStringsToEcoliTemplate(AW_root * /* awr */) {
1501    const char *pSearchColResults = 0;
1502
1503    if (iMapSearch) pSearchColResults = Host.get_search_background(iStartPos, iEndPos);
1504
1505    if (pSearchColResults) {
1506        int iColor = 0;
1507
1508        glNewList(MAP_SEARCH_STRINGS_TO_STRUCTURE, GL_COMPILE);
1509        {
1510            if (RNA3D->bPointSpritesSupported) glBegin(GL_POINTS);
1511            for (int i = iStartPos; i < iEndPos; i++) {
1512                if (RNA3D->bEColiRefInitialized) {
1513                    long absPos   = (long) i;
1514                    long EColiPos = EColiRef->abs_2_rel(absPos);
1515
1516                    for (Struct3Dinfo *t = start3D; t != NULL; t = t->next) {
1517                        if ((t->pos == EColiPos) && (pSearchColResults[i] >= 0)) {
1518                            iColor = pSearchColResults[i] - COLORLINK;
1519
1520                            if (ValidSearchColor(iColor, SEARCH)) {
1521                                RNA3D->cGraphics->SetColor(iColor);
1522                                PointsToQuads(t->x, t->y, t->z);
1523                            }
1524                            break;
1525                        }
1526                    }
1527                }
1528            }
1529            if (RNA3D->bPointSpritesSupported) glEnd();
1530        }
1531        glEndList();
1532
1533        glNewList(MAP_SEARCH_STRINGS_BACKBONE, GL_COMPILE);
1534        {
1535            int iLastClr = 0; int iLastPos = 0; Vector3 vLastPt;
1536            glBegin(GL_LINES);
1537            for (int i = iStartPos; i < iEndPos; i++) {
1538                if (RNA3D->bEColiRefInitialized) {
1539                    long absPos   = (long) i;
1540                    long EColiPos = EColiRef->abs_2_rel(absPos);
1541
1542                    for (Struct3Dinfo *t = start3D; t != NULL; t = t->next) {
1543                        if ((t->pos == EColiPos) && (pSearchColResults[i] >= 0)) {
1544                            iColor = pSearchColResults[i] - COLORLINK;
1545
1546                            if (ValidSearchColor(iColor, SEARCH)) {
1547                                if ((iLastClr == iColor) && (iLastPos == EColiPos-1)) {
1548                                    RNA3D->cGraphics->SetColor(iColor);
1549                                    glVertex3f(vLastPt.x, vLastPt.y, vLastPt.z);
1550                                    glVertex3f(t->x, t->y, t->z);
1551                                }
1552                                iLastPos = EColiPos;
1553                                iLastClr = iColor;
1554                                vLastPt.x = t->x; vLastPt.y = t->y; vLastPt.z = t->z;
1555                            }
1556                            break;
1557                        }
1558                    }
1559                }
1560            }
1561            glEnd();
1562        }
1563        glEndList();
1564
1565        RNA3D->bMapSearchStringsDispListCreated = true;
1566    }
1567    else cout<<"BuildColorString did not get the colors : SAI cannot be Visualized!"<<endl;
1568}
1569
1570void Structure3D::MapSaiToEcoliTemplate() {
1571    const char *pSearchColResults = 0;
1572
1573    if (iMapSAI && Host.SAIs_visualized()) {
1574        pSearchColResults = Host.get_SAI_background(iStartPos, iEndPos); // returns 0 if sth went wrong
1575    }
1576
1577    if (pSearchColResults) {
1578        int iColor = 0;
1579
1580        glNewList(MAP_SAI_TO_STRUCTURE, GL_COMPILE);
1581        {
1582            if (RNA3D->bPointSpritesSupported) glBegin(GL_POINTS);
1583
1584            for (int i = iStartPos; i < iEndPos; i++) {
1585                if (RNA3D->bEColiRefInitialized) {
1586                    long absPos   = (long) i;
1587                    long EColiPos = EColiRef->abs_2_rel(absPos);
1588
1589                    for (Struct3Dinfo *t = start3D; t != NULL; t = t->next) {
1590                        if ((t->pos == EColiPos) && (pSearchColResults[i] >= 0)) {
1591                            iColor = pSearchColResults[i-1] - SAICOLORS;
1592
1593                            if (ValidSearchColor(iColor, SAI)) {
1594                                RNA3D->cGraphics->SetColor(iColor);
1595                                PointsToQuads(t->x, t->y, t->z);
1596                            }
1597                            break;
1598                        }
1599                    }
1600                }
1601            }
1602
1603            if (RNA3D->bPointSpritesSupported) glEnd();
1604        }
1605        glEndList();
1606
1607        RNA3D->bMapSaiDispListCreated = true;
1608    }
1609    else cout<<"ED4_getSaiColorString did not get the colors : SAI cannot be Visualized!"<<endl;
1610}
Note: See TracBrowser for help on using the repository browser.