RNAlib-2.1.9
data_structures.h
Go to the documentation of this file.
1 #ifndef __VIENNA_RNA_PACKAGE_DATA_STRUCTURES_H__
2 #define __VIENNA_RNA_PACKAGE_DATA_STRUCTURES_H__
3 
4 #include "energy_const.h"
10 /* to use floats instead of doubles in pf_fold() comment next line */
11 #define LARGE_PF
12 #ifdef LARGE_PF
13 #define FLT_OR_DBL double
14 #else
15 #define FLT_OR_DBL float
16 #endif
17 
18 #ifndef NBASES
19 #define NBASES 8
20 #endif
21 
22 #ifndef MAXALPHA
23 
26 #define MAXALPHA 20
27 #endif
28 
32 #define MAXDOS 1000
33 
34 #define VRNA_GQUAD_MAX_STACK_SIZE 7
35 #define VRNA_GQUAD_MIN_STACK_SIZE 2
36 #define VRNA_GQUAD_MAX_LINKER_LENGTH 15
37 #define VRNA_GQUAD_MIN_LINKER_LENGTH 1
38 #define VRNA_GQUAD_MIN_BOX_SIZE ((4*VRNA_GQUAD_MIN_STACK_SIZE)+(3*VRNA_GQUAD_MIN_LINKER_LENGTH))
39 #define VRNA_GQUAD_MAX_BOX_SIZE ((4*VRNA_GQUAD_MAX_STACK_SIZE)+(3*VRNA_GQUAD_MAX_LINKER_LENGTH))
40 
41 
42 /*
43 * ############################################################
44 * Here are the type definitions of various datastructures
45 * shared among the Vienna RNA Package
46 * ############################################################
47 */
48 
52 typedef struct plist {
53  int i;
54  int j;
55  float p;
56  int type;
57 } plist;
58 
62 typedef struct cpair {
63  int i,j,mfe;
64  float p, hue, sat;
65 } cpair;
66 
71 typedef struct {
72  float X; /* X coords */
73  float Y; /* Y coords */
74 } COORDINATE;
75 
79 typedef struct sect {
80  int i;
81  int j;
82  int ml;
83 } sect;
84 
88 typedef struct bondT {
89  unsigned int i;
90  unsigned int j;
91 } bondT;
92 
96 typedef struct bondTEn {
97  int i;
98  int j;
99  int energy;
100 } bondTEn;
101 
106 typedef struct{
107  int dangles;
114  int noLP;
115  int noGU;
117  int logML;
118  int circ;
119  int gquad;
122 
126 typedef struct{
127  int id;
128  int stack[NBPAIRS+1][NBPAIRS+1];
129  int hairpin[31];
130  int bulge[MAXLOOP+1];
131  int internal_loop[MAXLOOP+1];
132  int mismatchExt[NBPAIRS+1][5][5];
133  int mismatchI[NBPAIRS+1][5][5];
134  int mismatch1nI[NBPAIRS+1][5][5];
135  int mismatch23I[NBPAIRS+1][5][5];
136  int mismatchH[NBPAIRS+1][5][5];
137  int mismatchM[NBPAIRS+1][5][5];
138  int dangle5[NBPAIRS+1][5];
139  int dangle3[NBPAIRS+1][5];
140  int int11[NBPAIRS+1][NBPAIRS+1][5][5];
141  int int21[NBPAIRS+1][NBPAIRS+1][5][5][5];
142  int int22[NBPAIRS+1][NBPAIRS+1][5][5][5][5];
143  int ninio[5];
144  double lxc;
145  int MLbase;
146  int MLintern[NBPAIRS+1];
147  int MLclosing;
148  int TerminalAU;
149  int DuplexInit;
150  int Tetraloop_E[200];
151  char Tetraloops[1401];
152  int Triloop_E[40];
153  char Triloops[241];
154  int Hexaloop_E[40];
155  char Hexaloops[1801];
156  int TripleC;
157  int MultipleCA;
158  int MultipleCB;
159  int gquad [VRNA_GQUAD_MAX_STACK_SIZE + 1]
160  [3*VRNA_GQUAD_MAX_LINKER_LENGTH + 1];
161 
162  double temperature;
166 } paramT;
167 
171 typedef struct{
172  int id;
173  double expstack[NBPAIRS+1][NBPAIRS+1];
174  double exphairpin[31];
175  double expbulge[MAXLOOP+1];
176  double expinternal[MAXLOOP+1];
177  double expmismatchExt[NBPAIRS+1][5][5];
178  double expmismatchI[NBPAIRS+1][5][5];
179  double expmismatch23I[NBPAIRS+1][5][5];
180  double expmismatch1nI[NBPAIRS+1][5][5];
181  double expmismatchH[NBPAIRS+1][5][5];
182  double expmismatchM[NBPAIRS+1][5][5];
183  double expdangle5[NBPAIRS+1][5];
184  double expdangle3[NBPAIRS+1][5];
185  double expint11[NBPAIRS+1][NBPAIRS+1][5][5];
186  double expint21[NBPAIRS+1][NBPAIRS+1][5][5][5];
187  double expint22[NBPAIRS+1][NBPAIRS+1][5][5][5][5];
188  double expninio[5][MAXLOOP+1];
189  double lxc;
190  double expMLbase;
191  double expMLintern[NBPAIRS+1];
192  double expMLclosing;
193  double expTermAU;
194  double expDuplexInit;
195  double exptetra[40];
196  double exptri[40];
197  double exphex[40];
198  char Tetraloops[1401];
199  double expTriloop[40];
200  char Triloops[241];
201  char Hexaloops[1801];
202  double expTripleC;
203  double expMultipleCA;
204  double expMultipleCB;
205  double expgquad[VRNA_GQUAD_MAX_STACK_SIZE + 1]
206  [3*VRNA_GQUAD_MAX_LINKER_LENGTH + 1];
207 
208  double kT;
209  double pf_scale;
211  double temperature;
212  double alpha;
221 } pf_paramT;
222 
223 
224 
225 /*
226 * ############################################################
227 * SUBOPT data structures
228 * ############################################################
229 */
230 
231 
235 typedef struct {
236  int i;
237  int j;
238 } PAIR;
239 
243 typedef struct {
244  int i;
245  int j;
246  int array_flag;
247 } INTERVAL;
248 
252 typedef struct {
253  float energy;
254  char *structure;
255 } SOLUTION;
256 
257 /*
258 * ############################################################
259 * COFOLD data structures
260 * ############################################################
261 */
262 
266 typedef struct cofoldF {
267  /* free energies for: */
268  double F0AB;
269  double FAB;
270  double FcAB;
271  double FA;
272  double FB;
273 } cofoldF;
274 
278 typedef struct ConcEnt {
279  double A0;
280  double B0;
281  double ABc;
282  double AAc;
283  double BBc;
284  double Ac;
285  double Bc;
286 } ConcEnt;
287 
291 typedef struct pairpro{
292  struct plist *AB;
293  struct plist *AA;
294  struct plist *A;
295  struct plist *B;
296  struct plist *BB;
297 }pairpro;
298 
309 typedef struct {
310  unsigned i;
311  unsigned j;
312  float p;
313  float ent;
314  short bp[8];
315  char comp;
316 } pair_info;
317 
318 
319 /*
320 * ############################################################
321 * FINDPATH data structures
322 * ############################################################
323 */
324 
328 typedef struct move {
329  int i; /* i,j>0 insert; i,j<0 delete */
330  int j;
331  int when; /* 0 if still available, else resulting distance from start */
332  int E;
333 } move_t;
334 
338 typedef struct intermediate {
339  short *pt;
340  int Sen;
341  int curr_en;
344 
348 typedef struct path {
349  double en;
350  char *s;
351 } path_t;
352 
353 /*
354 * ############################################################
355 * RNAup data structures
356 * ############################################################
357 */
358 
362 typedef struct pu_contrib {
363  double **H;
364  double **I;
365  double **M;
366  double **E;
367  int length;
368  int w;
369 } pu_contrib;
370 
374 typedef struct interact {
375  double *Pi;
376  double *Gi;
377  double Gikjl;
379  double Gikjl_wo;
380  int i;
381  int k;
382  int j;
383  int l;
384  int length;
385 } interact;
386 
390 typedef struct pu_out {
391  int len;
392  int u_vals;
393  int contribs;
394  char **header;
395  double **u_values;
396 } pu_out;
397 
401 typedef struct constrain{
402  int *indx;
403  char *ptype;
404 } constrain;
405 
406 /*
407 * ############################################################
408 * RNAduplex data structures
409 * ############################################################
410 */
411 
415 typedef struct {
416  int i;
417  int j;
418  int end;
419  char *structure;
420  double energy;
421  double energy_backtrack;
422  double opening_backtrack_x;
423  double opening_backtrack_y;
424  int offset;
425  double dG1;
426  double dG2;
427  double ddG;
428  int tb;
429  int te;
430  int qb;
431  int qe;
432 } duplexT;
433 
434 /*
435 * ############################################################
436 * RNAsnoop data structures
437 * ############################################################
438 */
439 
443 typedef struct node {
444  int k;
445  int energy;
446  struct node *next;
447 } folden;
448 
452 typedef struct {
453  int i;
454  int j;
455  int u;
456  char *structure;
457  float energy;
458  float Duplex_El;
459  float Duplex_Er;
460  float Loop_E;
461  float Loop_D;
462  float pscd;
463  float psct;
464  float pscg;
465  float Duplex_Ol;
466  float Duplex_Or;
467  float Duplex_Ot;
468  float fullStemEnergy;
469 } snoopT;
470 
471 
472 
473 
474 
475 
476 
477 /*
478 * ############################################################
479 * PKplex data structures
480 * ############################################################
481 */
482 
486 typedef struct dupVar{
487  int i;
488  int j;
489  int end;
490  char *pk_helix;
491  char *structure;
492  double energy;
493  int offset;
494  double dG1;
495  double dG2;
496  double ddG;
497  int tb;
498  int te;
499  int qb;
500  int qe;
501  int inactive;
502  int processed;
503 } dupVar;
504 
505 
506 
507 /*
508 * ############################################################
509 * 2Dfold data structures
510 * ############################################################
511 */
512 
527 typedef struct{
528  int k;
529  int l;
530  float en;
531  char *s;
533 
539 typedef struct{
542  char *ptype;
543  char *sequence;
544  short *S, *S1;
545  unsigned int maxD1;
546  unsigned int maxD2;
549  unsigned int *mm1;
550  unsigned int *mm2;
552  int *my_iindx;
554  double temperature;
555 
556  unsigned int *referenceBPs1;
557  unsigned int *referenceBPs2;
558  unsigned int *bpdist;
560  short *reference_pt1;
561  short *reference_pt2;
562  int circ;
563  int dangles;
564  unsigned int seq_length;
565 
566  int ***E_F5;
567  int ***E_F3;
568  int ***E_C;
569  int ***E_M;
570  int ***E_M1;
571  int ***E_M2;
572 
573  int **E_Fc;
574  int **E_FcH;
575  int **E_FcI;
576  int **E_FcM;
577 
578  int **l_min_values;
579  int **l_max_values;
580  int *k_min_values;
581  int *k_max_values;
582 
583  int **l_min_values_m;
584  int **l_max_values_m;
585  int *k_min_values_m;
586  int *k_max_values_m;
587 
588  int **l_min_values_m1;
589  int **l_max_values_m1;
590  int *k_min_values_m1;
591  int *k_max_values_m1;
592 
593  int **l_min_values_f;
594  int **l_max_values_f;
595  int *k_min_values_f;
596  int *k_max_values_f;
597 
598  int **l_min_values_f3;
599  int **l_max_values_f3;
600  int *k_min_values_f3;
601  int *k_max_values_f3;
602 
603  int **l_min_values_m2;
604  int **l_max_values_m2;
605  int *k_min_values_m2;
606  int *k_max_values_m2;
607 
608  int *l_min_values_fc;
609  int *l_max_values_fc;
610  int k_min_values_fc;
611  int k_max_values_fc;
612 
613  int *l_min_values_fcH;
614  int *l_max_values_fcH;
615  int k_min_values_fcH;
616  int k_max_values_fcH;
617 
618  int *l_min_values_fcI;
619  int *l_max_values_fcI;
620  int k_min_values_fcI;
621  int k_max_values_fcI;
622 
623  int *l_min_values_fcM;
624  int *l_max_values_fcM;
625  int k_min_values_fcM;
626  int k_max_values_fcM;
627 
628  /* auxilary arrays for remaining set of coarse graining (k,l) > (k_max, l_max) */
629  int *E_F5_rem;
630  int *E_F3_rem;
631  int *E_C_rem;
632  int *E_M_rem;
633  int *E_M1_rem;
634  int *E_M2_rem;
635 
636  int E_Fc_rem;
637  int E_FcH_rem;
638  int E_FcI_rem;
639  int E_FcM_rem;
640 
641 #ifdef COUNT_STATES
642  unsigned long ***N_F5;
643  unsigned long ***N_C;
644  unsigned long ***N_M;
645  unsigned long ***N_M1;
646 #endif
647 } TwoDfold_vars;
648 
661 typedef struct{
662  int k;
663  int l;
664  FLT_OR_DBL q;
666 
673 typedef struct{
674 
675  unsigned int alloc;
676  char *ptype;
677  char *sequence;
678  short *S, *S1;
679  unsigned int maxD1;
680  unsigned int maxD2;
682  double temperature; /* temperature in last call to scale_pf_params */
683  double init_temp; /* temperature in last call to scale_pf_params */
684  FLT_OR_DBL *scale;
685  FLT_OR_DBL pf_scale;
686  pf_paramT *pf_params; /* holds all [unscaled] pf parameters */
687 
688  int *my_iindx;
689  int *jindx;
691  short *reference_pt1;
692  short *reference_pt2;
693 
694  unsigned int *referenceBPs1;
695  unsigned int *referenceBPs2;
696  unsigned int *bpdist;
698  unsigned int *mm1;
699  unsigned int *mm2;
701  int circ;
702  int dangles;
703  unsigned int seq_length;
704 
705  FLT_OR_DBL ***Q;
706  FLT_OR_DBL ***Q_B;
707  FLT_OR_DBL ***Q_M;
708  FLT_OR_DBL ***Q_M1;
709  FLT_OR_DBL ***Q_M2;
710 
711  FLT_OR_DBL **Q_c;
712  FLT_OR_DBL **Q_cH;
713  FLT_OR_DBL **Q_cI;
714  FLT_OR_DBL **Q_cM;
715 
716  int **l_min_values;
717  int **l_max_values;
718  int *k_min_values;
719  int *k_max_values;
720 
721  int **l_min_values_b;
722  int **l_max_values_b;
723  int *k_min_values_b;
724  int *k_max_values_b;
725 
726  int **l_min_values_m;
727  int **l_max_values_m;
728  int *k_min_values_m;
729  int *k_max_values_m;
730 
731  int **l_min_values_m1;
732  int **l_max_values_m1;
733  int *k_min_values_m1;
734  int *k_max_values_m1;
735 
736  int **l_min_values_m2;
737  int **l_max_values_m2;
738  int *k_min_values_m2;
739  int *k_max_values_m2;
740 
741  int *l_min_values_qc;
742  int *l_max_values_qc;
743  int k_min_values_qc;
744  int k_max_values_qc;
745 
746  int *l_min_values_qcH;
747  int *l_max_values_qcH;
748  int k_min_values_qcH;
749  int k_max_values_qcH;
750 
751  int *l_min_values_qcI;
752  int *l_max_values_qcI;
753  int k_min_values_qcI;
754  int k_max_values_qcI;
755 
756  int *l_min_values_qcM;
757  int *l_max_values_qcM;
758  int k_min_values_qcM;
759  int k_max_values_qcM;
760 
761  /* auxilary arrays for remaining set of coarse graining (k,l) > (k_max, l_max) */
762  FLT_OR_DBL *Q_rem;
763  FLT_OR_DBL *Q_B_rem;
764  FLT_OR_DBL *Q_M_rem;
765  FLT_OR_DBL *Q_M1_rem;
766  FLT_OR_DBL *Q_M2_rem;
767 
768  FLT_OR_DBL Q_c_rem;
769  FLT_OR_DBL Q_cH_rem;
770  FLT_OR_DBL Q_cI_rem;
771  FLT_OR_DBL Q_cM_rem;
772 
774 
775 #endif