1 |
/********************************************************************** |
2 |
Copyright (C) 2005 by Chris Morley |
3 |
|
4 |
This program is free software; you can redistribute it and/or modify |
5 |
it under the terms of the GNU General Public License as published by |
6 |
the Free Software Foundation version 2 of the License. |
7 |
|
8 |
This program is distributed in the hope that it will be useful, |
9 |
but WITHOUT ANY WARRANTY; without even the implied warranty of |
10 |
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
11 |
GNU General Public License for more details. |
12 |
***********************************************************************/ |
13 |
#include "mol.hpp" |
14 |
#include "obconversion.hpp" |
15 |
#include "obmolecformat.hpp" |
16 |
#include "fingerprint.hpp" |
17 |
#include <vector> |
18 |
#include <string> |
19 |
#include <iomanip> |
20 |
|
21 |
using namespace std; |
22 |
namespace OpenBabel |
23 |
{ |
24 |
|
25 |
/// \brief Constructs and displays fingerprints. For details see OBFingerprint class |
26 |
class FingerprintFormat : public OBMoleculeFormat |
27 |
{ |
28 |
public: |
29 |
//Register this format type ID |
30 |
FingerprintFormat() {OBConversion::RegisterFormat("fpt",this);} |
31 |
|
32 |
virtual const char* Description() //required |
33 |
{ return |
34 |
"Fingerprint format\n \ |
35 |
Constructs and displays fingerprints and (for multiple input objects)\n \ |
36 |
the Tanimoto coefficient and whether a superstructure of the first object\n \ |
37 |
Options e.g. -xfFP3 -xn128\n \ |
38 |
f<id> fingerprint type\n \ |
39 |
N# fold to specified number of bits, 32, 64, 128, etc.\n \ |
40 |
h hex output when multiple molecules\n \ |
41 |
F displays the available fingerprint types\n \ |
42 |
"; |
43 |
}; |
44 |
|
45 |
virtual unsigned int Flags(){return NOTREADABLE;}; |
46 |
virtual bool WriteMolecule(OBBase* pOb, OBConversion* pConv); |
47 |
|
48 |
private: |
49 |
vector<unsigned int> firstfp; |
50 |
string firstname; |
51 |
bool IsPossibleSubstructure(vector<unsigned int>Mol, vector<unsigned int>Frag); |
52 |
}; |
53 |
|
54 |
//////////////////////////////////////////////////// |
55 |
//Make an instance of the format class |
56 |
FingerprintFormat theFingerprintFormat; |
57 |
|
58 |
//******************************************************************* |
59 |
bool FingerprintFormat::WriteMolecule(OBBase* pOb, OBConversion* pConv) |
60 |
{ |
61 |
ostream &ofs = *pConv->GetOutStream(); |
62 |
|
63 |
OBFingerprint* pFP; |
64 |
string id; |
65 |
if(pConv->IsOption("F")) |
66 |
{ |
67 |
while(OBFingerprint::GetNextFPrt(id, pFP)) |
68 |
ofs << id << " -- " << pFP->Description() << endl; |
69 |
return true; |
70 |
} |
71 |
|
72 |
|
73 |
bool hexoutput=false; |
74 |
if(pConv->IsOption("h") || (pConv->GetOutputIndex()==1 && pConv->IsLast())) |
75 |
hexoutput=true; |
76 |
|
77 |
string fpid; |
78 |
int nbits=0; |
79 |
const char* p=pConv->IsOption("f"); |
80 |
if(p) |
81 |
{ |
82 |
fpid=p; |
83 |
fpid = fpid.substr(0,fpid.find('"')); |
84 |
} |
85 |
|
86 |
pFP = OBFingerprint::FindFingerprint(fpid); |
87 |
if(!pFP) |
88 |
{ |
89 |
#ifdef HAVE_SSTREAM |
90 |
stringstream errorMsg; |
91 |
#else |
92 |
strstream errorMsg; |
93 |
#endif |
94 |
errorMsg << "Fingerprint type '" << fpid << "' not available" << endl; |
95 |
obErrorLog.ThrowError(__func__, errorMsg.str(), obError); |
96 |
return false; |
97 |
} |
98 |
|
99 |
p=pConv->IsOption("N"); |
100 |
if(p) |
101 |
nbits = atoi(p); |
102 |
|
103 |
vector<unsigned int> fptvec; |
104 |
if(!pFP->GetFingerprint(pOb, fptvec, nbits)) |
105 |
return false; |
106 |
|
107 |
OBMol* pmol = dynamic_cast<OBMol*>(pOb); |
108 |
if(pmol) |
109 |
ofs << ">" << pmol->GetTitle(); |
110 |
|
111 |
if(hexoutput) |
112 |
{ |
113 |
unsigned int i, bitsset=0; |
114 |
for (i=0;i<fptvec.size();++i) |
115 |
{ |
116 |
int wd = fptvec[i]; |
117 |
for(;wd;wd=wd<<1)//count bits set by shifting into sign bit until word==0 |
118 |
if(wd<0) ++bitsset; |
119 |
} |
120 |
ofs << " " << bitsset << " bits set. "; |
121 |
} |
122 |
|
123 |
if(pConv->GetOutputIndex()==1) |
124 |
{ |
125 |
//store the fingerprint and name of first molecule |
126 |
firstfp=fptvec; |
127 |
if(pmol) |
128 |
firstname=pmol->GetTitle(); |
129 |
if(firstname.empty()) |
130 |
firstname = "first mol"; |
131 |
} |
132 |
else |
133 |
{ |
134 |
ofs << " Tanimoto from " << firstname << " = " << OBFingerprint::Tanimoto(firstfp, fptvec); |
135 |
if(IsPossibleSubstructure(fptvec,firstfp)) |
136 |
ofs << "\nPossible superstructure of " << firstname; |
137 |
} |
138 |
ofs << endl; |
139 |
|
140 |
int i; |
141 |
|
142 |
if(hexoutput) |
143 |
{ |
144 |
for(i=fptvec.size()-1;i>=0;i--) |
145 |
{ |
146 |
ofs << hex << setfill('0') << setw(8) << fptvec[i] << " " ; |
147 |
if((fptvec.size()-i)%6==0) |
148 |
ofs <<endl; |
149 |
} |
150 |
ofs << dec << endl; |
151 |
} |
152 |
return true; |
153 |
} |
154 |
|
155 |
bool FingerprintFormat::IsPossibleSubstructure(vector<unsigned int>Mol, vector<unsigned int>Frag) |
156 |
{ |
157 |
//Returns false if Frag is definitely NOT a substructure of Mol |
158 |
unsigned int i; |
159 |
for (i=0;i<Mol.size();++i) |
160 |
if((Mol[i] & Frag[i]) ^ Frag[i]) return false; |
161 |
return true; |
162 |
} |
163 |
|
164 |
} //namespace OpenBabel |