One more detail: It's nothing to do with SMILES, because it does the exact
same thing if I modify the program to read SD files.

Craig

On Tue, Mar 19, 2013 at 3:12 PM, Noel O'Boyle <baoille...@gmail.com> wrote:

> I haven't looked into the details of your code (but I will if no-one
> else does), but regarding the relationship between different things,
> you may find the following notes I have made useful:
>
> In the Open Babel world, the "valence" methods of an atom refer to the
> number of bonds, rather than the sum of BOs (which is actual valence).
> In other words, they (a) are misnamed, and (b) our _impval should not
> be the underlying key data structure.
>
> GetValence() -> Num explicit bonds
> GetImplicitValence() -> Num explicit bonds plus num implicit bonds
> (i.e. total number of bonds) (_impval property on OBAtom)
> GetHvyAtomValence() -> Num explicit bonds to non-H atoms
>
> To find the number of implicit Hs, GetImplicitValence() - GetValence()
>
> BOSum() -> The sum of the BOs of explicit bonds (otherwise known as
> ExplicitValence)
>
> To find the actual valence: BOSum() + GetImplicitHydrogenCount()
>
> - Noel
>
> On 19 March 2013 19:08, Craig James <cja...@emolecules.com> wrote:
> > Below is a test program, highly stripped down, that illustrates a
> problem I
> > can't figure out.  The idea is to find ions that are neutral and show
> them
> > with a charge (don't worry about the chemistry behind this; this is
> highly
> > stripped down from the real code and just illustrates the inconsistency).
> > When I run it, I get this:
> >
> > [BH3+3]
> > [Li+]
> > [Na+]
> > [K+]
> > [Be+2]
> > [Mg+2]
> > [Ca+2]
> > [AlH3+3]
> > [SiH4+4]
> > [Fe+2]
> >
> > Note that B, Si, and Fe don't lose their hydrogens, whereas all the
> others
> > do ... in spite of the fact that the program explicitely removes all H
> atoms
> > and resets the "valence perceived" flag.
> >
> > Digging through smilesformat.cpp, atom.cpp, the typer and so forth, I'm
> > completely baffled as to what's going on.  (If anyone ever had time to
> > document the relationship and "synonym-ness" between bond order, number
> of
> > bonds, valence, implicit valence, charge, hydrogen count, and spin
> > multiplicity, that would sure be helpful!)
> >
> > I compiled this on Linux using:
> >
> > g++ -g -I/usr/local/openbabel/include/openbabel-2.0 -D_GNU_SOURCE \
> >    -L/usr/local/openbabel/lib -lopenbabel -ldl -o foo foo.cpp
> >
> > Thanks!
> > Craig
> >
> >
> > #include <sstream>
> >
> > #include <openbabel/babelconfig.h>
> > #include <openbabel/mol.h>
> > #include <openbabel/obconversion.h>
> >
> > using namespace std;
> > using namespace OpenBabel;
> >
> > static string test_smiles[] = {"[BH3]", "[LiH]", "[NaH]", "[KH]",
> "[BeH]",
> > "[MgH2]", "[CaH2]", "[AlH3]", "[SiH4]", "[FeH2]"};
> > #define NSMILES (sizeof(test_smiles) / sizeof(string))
> >
> > static void remove_explicit_hydrogen(OBMol *pmol)
> > {
> >   OBAtom *atom;
> >   vector<OBNodeBase*> delete_h;
> >   vector<OBNodeBase*>::iterator ai;
> >
> >   for (atom = pmol->BeginAtom(ai); atom; atom = pmol->NextAtom(ai)) {
> >     if (atom->GetAtomicNum() == 1)
> >       delete_h.push_back(atom);
> >   }
> >   for (ai = delete_h.begin(); ai != delete_h.end(); ai++)
> >     pmol->DeleteHydrogen((OBAtom*)*ai);
> > }
> >
> > int main(int argc,char **argv)
> > {
> >   OBMol *pmol = new OBMol;
> >   OBConversion conv(&cin, &cout);
> >   conv.SetInAndOutFormats("smi", "smi");
> >
> >   for (int i = 0; i < NSMILES; i++) {
> >     pmol->Clear();
> >     pmol->SetTitle("");
> >     conv.ReadString(pmol, test_smiles[i]);
> >
> >     if (pmol->Empty()) {
> >       cerr << "Couldn't read SMILES: '" << test_smiles[i] << "'\n";
> >       continue;
> >     }
> >     pmol->BeginModify();
> >     remove_explicit_hydrogen(pmol);
> >     OBAtom *a = pmol->GetAtom(1);
> >
> >     int anum = a->GetAtomicNum();
> >     int pos_charge = 0;
> >     switch (anum) {
> >     case 3: case 11: case 19: pos_charge = 1; break;    // +1 - Li, Na, K
> >     case 4: case 12: case 20: pos_charge = 2; break;    // +2 - Be, Mg,
> Ca
> >     case 26:                  pos_charge = 2; break;    // +2 - Fe
> >     case 5: case 13:          pos_charge = 3; break;    // +3 - B, Al
> >     case 14:                  pos_charge = 4; break;    // +4 - Si
> >     }
> >     a->SetFormalCharge(pos_charge);
> >     pmol->UnsetImplicitValencePerceived();
> >     atomtyper.AssignImplicitValence(*pmol);
> >     pmol->EndModify();
> >
> >     conv.Write(pmol);
> >   }
> > }
> >
> >
> >
> ------------------------------------------------------------------------------
> > Everyone hates slow websites. So do we.
> > Make your web apps faster with AppDynamics
> > Download AppDynamics Lite for free today:
> > http://p.sf.net/sfu/appdyn_d2d_mar
> > _______________________________________________
> > OpenBabel-Devel mailing list
> > OpenBabel-Devel@lists.sourceforge.net
> > https://lists.sourceforge.net/lists/listinfo/openbabel-devel
> >
>
------------------------------------------------------------------------------
Everyone hates slow websites. So do we.
Make your web apps faster with AppDynamics
Download AppDynamics Lite for free today:
http://p.sf.net/sfu/appdyn_d2d_mar
_______________________________________________
OpenBabel-Devel mailing list
OpenBabel-Devel@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/openbabel-devel

Reply via email to