File: fingerprint.cpp

package info (click to toggle)
openbabel 2.2.3-1
  • links: PTS, VCS
  • area: main
  • in suites: squeeze
  • size: 36,644 kB
  • ctags: 33,717
  • sloc: cpp: 242,528; ansic: 87,037; sh: 10,280; perl: 5,518; python: 5,156; pascal: 793; makefile: 747; cs: 392; xml: 97; ruby: 54; java: 23
file content (690 lines) | stat: -rw-r--r-- 26,160 bytes parent folder | download
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
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
/**********************************************************************
fingerprint.cpp - Implementation of fingerpring base class and fastsearching

Copyright (C) 2005 by Chris Morley
 
This file is part of the Open Babel project.
For more information, see <http://openbabel.sourceforge.net/>

This program is free software; you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation version 2 of the License.
 
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
GNU General Public License for more details.
***********************************************************************/

#include <openbabel/babelconfig.h>

#include <vector>
#include <algorithm>
#include <iosfwd>
#include <cstring>
#include <fstream>

#include <openbabel/fingerprint.h>
#include <openbabel/oberror.h>

using namespace std;
namespace OpenBabel
{
#if defined(__CYGWIN32__) || defined(__MINGW32__)
  // macro to implement static OBPlugin::PluginMapType& Map()
  PLUGIN_CPP_FILE(OBFingerprint)
#endif

  const unsigned int OBFingerprint::bitsperint = 8 * sizeof(unsigned int);

  void OBFingerprint::SetBit(vector<unsigned int>& vec, const unsigned int n)
  {
    vec[n/Getbitsperint()] |= (1 << (n % Getbitsperint()));
  }

  bool OBFingerprint::GetBit(const vector<unsigned int>& vec, const unsigned int n)
  {
    unsigned int word =vec[n/Getbitsperint()];
    return (word &= (1 << (n % Getbitsperint())))!=0;
  }

  ////////////////////////////////////////	
  void OBFingerprint::Fold(vector<unsigned int>& vec, unsigned int nbits)
  {
    if(nbits<Getbitsperint())
    {
      stringstream ss;
      ss << "Can't fold to less than " << Getbitsperint() << "bits";
      obErrorLog.ThrowError(__FUNCTION__, ss.str(), obError);
      return;
    }
    while(vec.size()*Getbitsperint()/2 >= nbits) 
      vec.erase(transform(vec.begin(),vec.begin()+vec.size()/2,
                          vec.begin()+vec.size()/2, vec.begin(), bit_or()), vec.end());
  }

  ////////////////////////////////////////
/*  bool OBFingerprint::GetNextFPrt(std::string& id, OBFingerprint*& pFPrt)
  {
    Fptpos iter;
    if(id.empty())
      iter=FPtsMap().begin();
    else
      {
        iter=FPtsMap().find(id);
        if(iter!=FPtsMap().end())
          ++iter;
      }
    if(iter==FPtsMap().end())
      return false;
    id    = iter->first;
    pFPrt = iter->second;
    return true;
  }

  OBFingerprint* OBFingerprint::FindFingerprint(const string& ID)
  {
    if(ID.empty())
      return _pDefault;
    Fptpos iter = FPtsMap().find(ID);
    if(iter==FPtsMap().end())
      return NULL;
    else
      return iter->second;
  }
*/
  double OBFingerprint::Tanimoto(const vector<unsigned int>& vec1, const vector<unsigned int>& vec2) 
  {
    //Independent of sizeof(unsigned int)
    if(vec1.size()!=vec2.size())
      return -1; //different number of bits
    int andbits=0, orbits=0;
    for (unsigned i=0;i<vec1.size();++i)
      {
        int andfp = vec1[i] & vec2[i];
        int orfp = vec1[i] | vec2[i];
        //Count bits
      /* GCC 3.4 supports a "population count" builtin, which on many targets is
         implemented with a single instruction.  There is a fallback definition
         in libgcc in case a target does not have one, which should be just as
         good as the static function below.  */
#if __GNUC__ > 3 || (__GNUC__ == 3 && __GNUC_MINOR__ >= 4)
        andbits += __builtin_popcount(andfp);
        orbits += __builtin_popcount(orfp);
#else
        for(;andfp;andfp=andfp<<1)
          if(andfp<0) ++andbits;
        for(;orfp;orfp=orfp<<1)
          if(orfp<0) ++orbits;
#endif
      }
    if(orbits==0)
      return 0.0;
    return((double)andbits/(double)orbits);
  }

  //*****************************************************************
  bool FastSearch::Find(OBBase* pOb, vector<unsigned int>& SeekPositions,
                        unsigned int MaxCandidates)
  {
    ///Finds chemical objects in datafilename (which must previously have been indexed)
    ///that have all the same bits set in their fingerprints as in the fingerprint of
    ///a pattern object. (Usually a substructure search in molecules.)
    ///The type of fingerprint and its degree of folding does not have to be specified
    ///here because the values in the index file are used.
    ///The positions of the candidate matching molecules in the original datafile are returned.
	
    vector<unsigned int> vecwords;
    _pFP->GetFingerprint(pOb,vecwords, _index.header.words * OBFingerprint::Getbitsperint());
	
    vector<unsigned int>candidates; //indices of matches from fingerprint screen
    candidates.reserve(MaxCandidates);

    unsigned int dataSize = _index.header.nEntries;
    //	GetFingerprint(mol, vecwords, _index.header.words, _index.header.fptype); 
	
    unsigned int words = _index.header.words;
    unsigned int* nextp = &_index.fptdata[0];
    unsigned int* ppat0 = &vecwords[0];
    register unsigned int* p;
    register unsigned int* ppat;
    register unsigned int a;
    unsigned int i; // need address of this, can't be register
    for(i=0;i<dataSize; ++i) //speed critical section
      {
        p=nextp;
        nextp += words;
        ppat=ppat0;
        a=0;
        while(p<nextp)
          {
            if ( (a=((*ppat) & (*p++)) ^ (*ppat++)) ) break;
          }
        if(!a)
          {
            candidates.push_back(i);
            if(candidates.size()>=MaxCandidates)
              break;
          }
      }

    if(i<_index.header.nEntries) //premature end to search
      {
        stringstream errorMsg;
        errorMsg << "Stopped looking after " << i << " molecules." << endl;
        obErrorLog.ThrowError(__FUNCTION__, errorMsg.str(), obWarning);
      }

    vector<unsigned int>::iterator itr;
    for(itr=candidates.begin();itr!=candidates.end();++itr)
      {
        SeekPositions.push_back(_index.seekdata[*itr]);
      }
    return true;
  }

////////////////////////////////////////////////////////////
 bool FastSearch::FindMatch(OBBase* pOb, vector<unsigned int>& SeekPositions,
                            unsigned int MaxCandidates)
{
//Similar to FastSearch::Find() except that successful candidates have all bits the same as the target
  vector<unsigned int> vecwords;
  _pFP->GetFingerprint(pOb,vecwords, _index.header.words * OBFingerprint::Getbitsperint());

  vector<unsigned int>candidates; //indices of matches from fingerprint screen

  unsigned int dataSize = _index.header.nEntries;
  unsigned int words = _index.header.words;
  unsigned int* nextp = &_index.fptdata[0]; // start of next FP in index
  unsigned int* ppat0 = &vecwords[0];       // start of target FP
  register unsigned int* p;                 // current position in index
  register unsigned int* ppat;              // current position in target FP
  unsigned int i; // need address of this, can't be register
  for(i=0;i<dataSize; ++i) //speed critical section
  {
    p=nextp;
    nextp += words;
    ppat=ppat0;

    while((*p++ == *ppat++ ) && (p<nextp));

    if(p==nextp)
    {
      candidates.push_back(i);
      if(candidates.size()>=MaxCandidates)
        break;
    }
  }
  
  vector<unsigned int>::iterator itr;
  for(itr=candidates.begin();itr!=candidates.end();++itr)
    {
      SeekPositions.push_back(_index.seekdata[*itr]);
    }
  return true;
}

  /////////////////////////////////////////////////////////
  bool FastSearch::FindSimilar(OBBase* pOb, multimap<double, unsigned int>& SeekposMap,
                               double MinTani)
  {
    vector<unsigned int> targetfp;
    _pFP->GetFingerprint(pOb,targetfp, _index.header.words * OBFingerprint::Getbitsperint());

    unsigned int words = _index.header.words;
    unsigned int dataSize = _index.header.nEntries;
    unsigned int* nextp = &_index.fptdata[0];
    register unsigned int* p;
    register unsigned int i;
    for(i=0;i<dataSize; ++i) //speed critical section
      {
        p=nextp;
        nextp += words;
        double tani = OBFingerprint::Tanimoto(targetfp,p);
        if(tani>MinTani)
          SeekposMap.insert(pair<const double, unsigned int>(tani,_index.seekdata[i]));
      }	
    return true;
  }

  /////////////////////////////////////////////////////////
  bool FastSearch::FindSimilar(OBBase* pOb, multimap<double, unsigned int>& SeekposMap,
                               int nCandidates)
  {
    ///If nCandidates is zero or omitted the original size of the multimap is used
    if(nCandidates)
      {
        //initialise the multimap with nCandidate zero entries
        SeekposMap.clear();
        int i;
        for(i=0;i<nCandidates;++i)
          SeekposMap.insert(pair<const double, unsigned int>(0,0));
      }
    else if(SeekposMap.size()==0)
      return false;

    vector<unsigned int> targetfp;
    _pFP->GetFingerprint(pOb,targetfp, _index.header.words * OBFingerprint::Getbitsperint());

    unsigned int words = _index.header.words;
    unsigned int dataSize = _index.header.nEntries;
    unsigned int* nextp = &_index.fptdata[0];
    register unsigned int* p;
    register unsigned int i;
    for(i=0;i<dataSize; ++i) //speed critical section
      {
        p=nextp;
        nextp += words;
        double tani = OBFingerprint::Tanimoto(targetfp,p);
        if(tani>SeekposMap.begin()->first)
          {	
            SeekposMap.insert(pair<const double, unsigned int>(tani,_index.seekdata[i]));
            SeekposMap.erase(SeekposMap.begin());
          }
      }	
    return true;
  }

  /////////////////////////////////////////////////////////
  string FastSearch::ReadIndex(istream* pIndexstream)
  {
    //Reads fs index from istream into member variables
    _index.Read(pIndexstream);

    _pFP = _index.CheckFP();	
    if(!_pFP)
      *(_index.header.datafilename) = '\0';

    return _index.header.datafilename; //will be empty on error
  }

  //////////////////////////////////////////////////////////
  string FastSearch::ReadIndexFile(string IndexFilename)
  {
    ifstream ifs(IndexFilename.c_str(),ios::binary);
    if(ifs)
      return ReadIndex(&ifs);
    else
    {
      string dum;
      return dum;
    }
  }

  //////////////////////////////////////////////////////////
  bool FptIndex::Read(istream* pIndexstream)
  {
    pIndexstream->read((char*)&(header), sizeof(FptIndexHeader));
    pIndexstream->seekg(header.headerlength);//allows header length to be changed
    if(pIndexstream->fail() || header.headerlength != sizeof(FptIndexHeader))
      {
        *(header.datafilename) = '\0';
        return false;
      }

    unsigned int nwords = header.nEntries * header.words;
    fptdata.resize(nwords);
    seekdata.resize(header.nEntries);

    pIndexstream->read((char*)&(fptdata[0]), sizeof(unsigned int) * nwords);
    pIndexstream->read((char*)&(seekdata[0]), sizeof(unsigned int) * header.nEntries);
	
    if(pIndexstream->fail())
      {
        *(header.datafilename) = '\0';
        return false;
      }
    return true;
  }

  //////////////////////////////////////////////////////////
  OBFingerprint* FptIndex::CheckFP()
  {
    //check that fingerprint type is available
    OBFingerprint* pFP = OBFingerprint::FindFingerprint(header.fpid);
    if(!pFP)
      {
        stringstream errorMsg;
        errorMsg << "Index has Fingerprints of type '" << header.fpid 
                 << " which is not currently loaded." << endl;
        obErrorLog.ThrowError(__FUNCTION__, errorMsg.str(), obError);
      }
    return pFP; //NULL if not available
  }

  //*******************************************************
  FastSearchIndexer::FastSearchIndexer(string& datafilename, ostream* os, 
                                       std::string& fpid, int FptBits)
  {
    ///Starts indexing process
    _indexstream = os;
    _nbits=FptBits;
    _pindex= new FptIndex;
    _pindex->header.headerlength = sizeof(FptIndexHeader);
    strncpy(_pindex->header.fpid,fpid.c_str(),15);
    _pindex->header.fpid[15]='\0'; //ensure fpid is terminated at 15 characters.
    strncpy(_pindex->header.datafilename, datafilename.c_str(), 255);

    //check that fingerprint type is available
    _pFP = _pindex->CheckFP();	
  }

  /////////////////////////////////////////////////////////////
  FastSearchIndexer::FastSearchIndexer(FptIndex* pindex, std::ostream* os)
  {
    //Uses existing index
    _indexstream = os;
    _pindex = pindex;
    _nbits  = _pindex->header.words * OBFingerprint::Getbitsperint();

    //check that fingerprint type is available
    _pFP = _pindex->CheckFP();	
  }

  /////////////////////////////////////////////////////////////
  FastSearchIndexer::~FastSearchIndexer()
  {
    ///Saves index file
    _pindex->header.nEntries = _pindex->seekdata.size();
    _indexstream->write((const char*)&_pindex->header, sizeof(FptIndexHeader));
    _indexstream->write((const char*)&_pindex->fptdata[0], _pindex->fptdata.size()*sizeof(unsigned int));
    _indexstream->write((const char*)&_pindex->seekdata[0], _pindex->seekdata.size()*sizeof(unsigned int));
    if(!_indexstream)
      obErrorLog.ThrowError(__FUNCTION__,
                            "Difficulty writing index", obWarning);
    delete _pindex;
  }

  ///////////////////////////////////////////////////////////////
  bool FastSearchIndexer::Add(OBBase* pOb, std::streampos seekpos)
  {
    ///Adds a fingerprint
	
    vector<unsigned int> vecwords;
    if(!_pFP)
      return false;
    if(_pFP->GetFingerprint(pOb, vecwords, _nbits))
      {
        _pindex->header.words = vecwords.size(); //Use size as returned from fingerprint
        for(unsigned int i=0;i<_pindex->header.words;++i)
          _pindex->fptdata.push_back(vecwords[i]);
        _pindex->seekdata.push_back(seekpos);
        return true;	
      }
    obErrorLog.ThrowError(__FUNCTION__, "Failed to make a fingerprint", obWarning);
    return false;
  }

  /*!
    \class OBFingerprint fingerprint.h <openbabel/fingerprint.h>
    These fingerprints are condensed representation of molecules (or other objects)
    as a list of boolean values (actually bits in a vector<unsigned>) with length 
    of a power of 2. The main motivation is for fast searching of data sources
    containing large numbers of molecules (up to several million). Open Babel
    provides some routines which can search text files containing lists of molecules
    in any format. See the documentation on the class FastSearch.

    There are descriptions of molecular fingerprints at <br>
    http://www.daylight.com/dayhtml/doc/theory/theory.finger.html) and <br>
    http://www.mesaac.com/Fingerprint.htm <br>
    Many methods of preparing fingerprints have been described, but the type supported
    currently in OpenBabel has each bit representing a substructure (or other
    molecular property). If a substructure is present in the molecule, then a 
    particular bit is set to 1. But because the hashing method may also map other 
    substructures to the same bit, a match does not guarantee that a particular 
    substructure is present; there may be false positives.  However, with proper design,
    a large fraction of irrelevant molecules in a data set can be eliminated in a 
    fast search with boolean methods on the fingerprints. 
    It then becomes feasible to make a definitive substructure search by
    conventional methods on this reduced list even if it is slow. 

    OpenBabel provides a framework for applying new types of fingerprints without
    changing any existing code. They are derived from OBFingerprint and the
    source file is just compiled with the rest of OpenBabel. Alternatively,
    they can be separately compiled as a DLL or shared library and discovered
    when OpenBabel runs.

    For more on these specific implementations of fingerprints in Open
    Babel, please take a look at the developer's wiki:
    http://openbabel.org/wiki/Fingerprints
 
    Fingerprints derived from this abstract base class OBFingerprint can be for any 
    object derived from OBBase (not just for OBMol). 
    Each derived class provides an ID as a string and OBFingerprint keeps a map of 
    these to provides a pointer to the class when requested in FindFingerprint.
   
    <h4>-- To define a fingerprint type --</h4>
    The classes derived form OBFingerprint are required to provide
    a GetFingerprint() routine and a Description() routine
    \code
    class MyFpType : OBFingerprint 
    {
       MyFpType(const char* id) : OBFingerprint(id){};

       virtual bool GetFingerprint(OBBase* pOb, vector<unsigned int>& fp, int nbits) 
       {
          //Convert pOb to the required type, usually OBMol
          OBMol* pmol = dynamic_cast<OBMol*>(pOb);
          fp.resize(required_number_of_words);
          ... 
          use SetBit(fp,n); to set the nth bit

          if(nbits)
             Fold(fp, nbits);
       }
       
       virtual const char* Description(){ return "Some descriptive text";}
       ...
    };
    \endcode
    
    Declare a global instance with the ID you will use in -f options to specify
    its use.
    \code
    MyFpType theMyFpType("myfpID");
    \endcode

    <h4>-- To obtain a fingerprint --</h4>
    \code
    OBMol mol;
    ...
    vector<unsigned int> fp;
    OBFingerprint::GetDefault()->GetFingerprint(&mol, fp); //gets default size of fingerprint
    \endcode	
    or
    \code
    vector<unsigned int> fp;
    OBFingerPrint* pFP = OBFingerprint::FindFingerprint("myfpID");
    ...and maybe...
    pFP->GetFingerprint(&mol,fp, 128); //fold down to 128bits if was originally larger
    \endcode

    <h4>-- To print a list of available fingerprint types --</h4>
    \code
    std::string id;
    OBFingerPrint* pFPrt=NULL;
    while(OBFingerprint::GetNextFPrt(id, pFPrt))
    {
       cout << id << " -- " << pFPrt->Description() << endl;
    }
    \endcode

    Fingerprints are handled as vector<unsigned int> so that the number of bits
    in this vector and their order will be platform and compiler
    dependent, because of size of int types and endian differences.
    Use fingerprints (and fastsearch indexes containing them) only 
    for comparing with other fingerprints prepared on the same machine.

    The FingerprintFormat class is an output format which displays fingerprints
    as hexadecimal. When multiple molecules are supplied it will calculate the
    Tanimoto coefficient from the first molecule to each of the others. It also
    shows whether the first molecule is a possible substructure to all the others,
    i.e. whether all the bits set in the fingerprint for the first molecule are
    set in the fingerprint of the others. To display hexadecimal information when
    multiple molecules are provided it is necessay to use the -xh option.

    To see a list of available format types, type babel -F on the command line.
    The -xF option of the FingerprintFormat class also provides this output, but due
    to a quirk in the way the program works, it is necessary to have a valid input
    molecule for this option to work.
  */
	
  /*! \class FastSearch fingerprint.h <openbabel/fingerprint.h>
    The FastSearch class searches an index to a datafile containing a list of molecules
    (or other objects) prepared by FastSearchIndexer.

    OpenBabel can also search files for molecules containing a substructure specified
    by a SMARTS string, using OBSmartsPattern or from the command line:
    \code
    babel datafile.xxx -outfile.yyy -sSMARTS
    \endcode
    But when the data file contains more than about 10,000 molecules this becomes
    rather too slow to be used interactively. To do it more quickly, an index
    of the molecules containing their fingerprints (see OBFingerprint) is prepared using 
    FastSearchIndexer. The indexing process may take a long time but only has to 
    be done once. The index can be searched very quickly with FastSearch. Because 
    of the nature of fingerprints a match to a bit does not guarantee 
    the presence of a particular substructure or other molecular property, so that
    a definitive answer may require a subsequent search of the (much reduced) list
    of candidates by another method (like OBSmartsPattern).

    Note that the index files are not portable. They need to be prepared on the
    computer that will access them.

    <h4>Using FastSearch and FastSearchIndexer in a program</h4>

    The index has two tables:
    - an array of fingerprints of the molecules
    - an array of the seek positions in the datasource file of all the molecules

    <h4>To prepare an fastsearch index file:</h4>
    - Open an ostream to the index file.
    - Make a FastSearchIndexer object on the heap or the stack, passing in as parameters:
    - the datafile name, the indexfile stream,
    - the id of the fingerprint type to be used,
    -  the number of bits the fingerprint is to be folded down to, If it is to be left
    unfolded, set fpid to 0 or do not specify it.
    .
    - For each molecule, call Add() with its pointer and its position in the datafile.<br>
    Currently the std::streampos value is implicitly cast to unsigned int so that
    for 32bit machines the datafile has to be no longer than about 2Gbyte.
    - The index file is written when the FastSearchIndexer object is deleted or goes out
    of scope.

    <h4>To search in a fastsearch index file</h4>
 
    - Open a std::istream to the indexfile (in binary mode on some systems)
    - Make a FastSearch object, read the index and open the datafile whose 
    name it provides
    \code	
    ifstream ifs(indexname,ios::binary);
    FastSearch fs;
    string datafilename = fs.ReadIndex(&ifs);
    if(datafilename.empty()
       return false;

    ifstream datastream(datafilename);
    if(!datastream)
       return false;
    \endcode

    <strong>To do a search for molecules which have all the substructure bits the
    OBMol object, patternMol</strong>
    \code
    vector<unsigned int>& SeekPositions;
    if(!fs.Find(patternMol, SeekPositions, MaxCandidates))
	    for(itr=SeekPositions.begin();itr!=SeekPositions.end();++itr)
      {
         datastream.seekg(*itr);
         ... read the candidate molecule
         and subject to more rigorous test if necessary
      }
    \endcode

    <strong>To do a similarity search based on the Tanimoto coefficient</strong>
    This is defined as: <br>
    <em>Number of bits set in (patternFP & targetFP)/Number of bits in (patternFP | targetFP)</em><br>
    where the boolean operations between the fingerprints are bitwise<br>
    The Tanimoto coefficient has no absolute meaning and depends on
    the design of the fingerprint.
    \code
    multimap<double, unsigned int> SeekposMap;
    // to find n best molecules
    fs.FindSimilar(&patternMol, SeekposMap, n);
    ...or
    // to find molecules with Tanimoto coefficient > MinTani
    fs.FindSimilar(&patternMol, SeekposMap, MinTani);

    multimap<double, unsigned int>::reverse_iterator itr;
    for(itr=SeekposMap.rbegin();itr!=SeekposMap.rend();++itr)
    {
       datastream.seekg(itr->second);
       // ... read the candidate molecule
       double tani = itr->first;
    }
    \endcode

    The FastSearchFormat class facilitates the use of these routine from the
    command line or other front end program. For instance:

    <strong>Prepare an index:</strong>
    \code
    babel datafile.xxx index.fs
    \endcode
    or
    \code
    babel datafile.xxx -ofs namedindex
    \endcode
    With options you can specify:
    - which type of fingerprint to be used, e.g. -xfFP2, 
    -	whether it is folded to a specified number of bits, e.g. -xn128 
    (which should be a power of2)
    - whether to pre-select the molecules which are indexed:
    - by structure e.g only ethers and esters, -sCOC
    - by excluding molecules with bezene rings, -vc1ccccc1
    - by position in the datafile e.g. only mols 10 to 90, -f10 -l90
    .
    <strong>Substructure search</strong> in a previously prepared index file
    \code
    babel index.fs outfile.yyy -sSMILES
    \endcode
    The datafile can also be used as the input file, provided the input format is specified as fs 
    \code
    babel datafile.xxx outfile.yyy -sSMILES -ifs
    \endcode
    A "slow" search not using fingerprints would be done on the same data by omitting -ifs.
    A "slow" search can use SMARTS strings, but the fast search is restricted
    to the subset which is valid SMILES.

    With the -S option, the target can be specified as a molecule in a file of any format 
    \code
    babel datafile.xxx outfile.yyy -Spattern_mol.zzz -ifs
    \endcode
    These searches have two stages: a fingerprint search which produces
    a number of candidate molecules and a definitive search which selects
    from these using SMARTS pattern matching. The maximum number of candidates
    is 4000 by default but you can change this with an option
    e.g. -al 8000  (Note that you need the space between l and the number.)
    If the fingerprint search reaches the maximum number it will not
    look further and will tell you at what stage it stopped.

    <strong>Similarity search</strong> in a previously prepared index file<br>
    This rather done (rather than a substructure search) if the -at option is used,
    \code
    babel datafile.xxx outfile.yyy -sSMILES -ifs -at0.7
    \endcode
    for instance
    - -at0.7 will recover all molecules with Tanimoto greater than 0.7 
    - -at15 (no decimal point) will recover the 15 molecules with largest coefficients.
    - -aa will add the Tanimoto coefficient to the titles of the output molecules.

    All stages, the indexing, the interpretation of the SMILES string in the -s option,
    the file in the -S option and the final SMARTS filter convert to OBMol and apply
    ConvertDativeBonds(). This ensures thatforms such as[N+]([O-])=O  and N(=O)=O 
    are recognized as chemically identical.     
  */

}//Openbabel

//! \file fingerprint.cpp
//! \brief Definitions for OBFingerprint base class and fastsearch classes