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 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208
|
//
// Copyright (C) 2003-2016 Greg Landrum and Rational Discovery LLC
//
// @@ All Rights Reserved @@
// This file is part of the RDKit.
// The contents are covered by the terms of the BSD license
// which is included in the file license.txt, found at the root
// of the RDKit source tree.
//
#include "rdmolops.h"
#include <RDGeneral/BoostStartInclude.h>
#include <RDBoost/python.h>
#include <boost/dynamic_bitset.hpp>
#include <RDGeneral/BoostStartInclude.h>
#include <RDGeneral/types.h>
#include <GraphMol/RDKitBase.h>
#include <GraphMol/MolOps.h>
#include <GraphMol/new_canon.h>
#include <GraphMol/SmilesParse/SmilesParse.h>
#include <GraphMol/SmilesParse/SmilesWrite.h>
#include <GraphMol/SmilesParse/SmartsWrite.h>
#include <GraphMol/FileParsers/FileParsers.h>
#include <GraphMol/FileParsers/SequenceParsers.h>
#include <GraphMol/FileParsers/SequenceWriters.h>
#include <RDGeneral/BadFileException.h>
#include <RDGeneral/FileParseException.h>
#include <RDBoost/Wrap.h>
#include <RDGeneral/Exceptions.h>
#include <RDGeneral/BadFileException.h>
#include <GraphMol/SanitException.h>
namespace python = boost::python;
using namespace RDKit;
void rdSanitExceptionTranslator(RDKit::MolSanitizeException const &x) {
std::ostringstream ss;
ss << "Sanitization error: " << x.message();
PyErr_SetString(PyExc_ValueError, ss.str().c_str());
}
void rdBadFileExceptionTranslator(RDKit::BadFileException const &x) {
std::ostringstream ss;
ss << "File error: " << x.message();
PyErr_SetString(PyExc_IOError, ss.str().c_str());
}
namespace RDKit {
std::string pyObjectToString(python::object input) {
python::extract<std::string> ex(input);
if (ex.check()) return ex();
std::wstring ws = python::extract<std::wstring>(input);
return std::string(ws.begin(), ws.end());
}
ROMol *MolFromSmiles(python::object ismiles, bool sanitize,
python::dict replDict) {
std::map<std::string, std::string> replacements;
for (unsigned int i = 0;
i < python::extract<unsigned int>(replDict.keys().attr("__len__")());
++i) {
replacements[python::extract<std::string>(replDict.keys()[i])] =
python::extract<std::string>(replDict.values()[i]);
}
RWMol *newM;
std::string smiles = pyObjectToString(ismiles);
try {
newM = SmilesToMol(smiles, 0, sanitize, &replacements);
} catch (...) {
newM = nullptr;
}
return static_cast<ROMol *>(newM);
}
ROMol *MolFromSmarts(python::object ismarts, bool mergeHs,
python::dict replDict) {
std::map<std::string, std::string> replacements;
for (unsigned int i = 0;
i < python::extract<unsigned int>(replDict.keys().attr("__len__")());
++i) {
replacements[python::extract<std::string>(replDict.keys()[i])] =
python::extract<std::string>(replDict.values()[i]);
}
std::string smarts = pyObjectToString(ismarts);
RWMol *newM;
try {
newM = SmartsToMol(smarts, 0, mergeHs, &replacements);
} catch (...) {
newM = nullptr;
}
return static_cast<ROMol *>(newM);
}
ROMol *MolFromTPLFile(const char *filename, bool sanitize = true,
bool skipFirstConf = false) {
RWMol *newM;
try {
newM = TPLFileToMol(filename, sanitize, skipFirstConf);
} catch (RDKit::BadFileException &e) {
PyErr_SetString(PyExc_IOError, e.message());
throw python::error_already_set();
} catch (...) {
newM = nullptr;
}
return static_cast<ROMol *>(newM);
}
ROMol *MolFromTPLBlock(python::object itplBlock, bool sanitize = true,
bool skipFirstConf = false) {
std::istringstream inStream(pyObjectToString(itplBlock));
unsigned int line = 0;
RWMol *newM;
try {
newM = TPLDataStreamToMol(&inStream, line, sanitize, skipFirstConf);
} catch (...) {
newM = nullptr;
}
return static_cast<ROMol *>(newM);
}
ROMol *MolFromMolFile(const char *molFilename, bool sanitize, bool removeHs,
bool strictParsing) {
RWMol *newM = nullptr;
try {
newM = MolFileToMol(molFilename, sanitize, removeHs, strictParsing);
} catch (RDKit::BadFileException &e) {
PyErr_SetString(PyExc_IOError, e.message());
throw python::error_already_set();
} catch (RDKit::FileParseException &e) {
BOOST_LOG(rdWarningLog) << e.message() << std::endl;
} catch (...) {
}
return static_cast<ROMol *>(newM);
}
ROMol *MolFromMolBlock(python::object imolBlock, bool sanitize, bool removeHs,
bool strictParsing) {
std::istringstream inStream(pyObjectToString(imolBlock));
unsigned int line = 0;
RWMol *newM = nullptr;
try {
newM =
MolDataStreamToMol(inStream, line, sanitize, removeHs, strictParsing);
} catch (RDKit::FileParseException &e) {
BOOST_LOG(rdWarningLog) << e.message() << std::endl;
} catch (...) {
}
return static_cast<ROMol *>(newM);
}
ROMol *MolFromSVG(python::object imolBlock, bool sanitize, bool removeHs) {
RWMol *res = nullptr;
res = RDKitSVGToMol(pyObjectToString(imolBlock), sanitize, removeHs);
return static_cast<ROMol *>(res);
}
ROMol *MolFromMol2File(const char *molFilename, bool sanitize = true,
bool removeHs = true, bool cleanupSubstructures = true) {
RWMol *newM;
try {
newM = Mol2FileToMol(molFilename, sanitize, removeHs, CORINA,
cleanupSubstructures);
} catch (RDKit::BadFileException &e) {
PyErr_SetString(PyExc_IOError, e.message());
throw python::error_already_set();
} catch (...) {
newM = nullptr;
}
return static_cast<ROMol *>(newM);
}
ROMol *MolFromMol2Block(std::string mol2Block, bool sanitize = true,
bool removeHs = true,
bool cleanupSubstructures = true) {
std::istringstream inStream(mol2Block);
RWMol *newM;
try {
newM = Mol2DataStreamToMol(inStream, sanitize, removeHs, CORINA,
cleanupSubstructures);
} catch (...) {
newM = nullptr;
}
return static_cast<ROMol *>(newM);
}
ROMol *MolFromPDBFile(const char *filename, bool sanitize, bool removeHs,
unsigned int flavor, bool proximityBonding) {
RWMol *newM = nullptr;
try {
newM = PDBFileToMol(filename, sanitize, removeHs, flavor, proximityBonding);
} catch (RDKit::BadFileException &e) {
PyErr_SetString(PyExc_IOError, e.message());
throw python::error_already_set();
} catch (RDKit::FileParseException &e) {
BOOST_LOG(rdWarningLog) << e.message() << std::endl;
} catch (...) {
}
return static_cast<ROMol *>(newM);
}
ROMol *MolFromPDBBlock(python::object molBlock, bool sanitize, bool removeHs,
unsigned int flavor, bool proximityBonding) {
std::istringstream inStream(pyObjectToString(molBlock));
RWMol *newM = nullptr;
try {
newM = PDBDataStreamToMol(inStream, sanitize, removeHs, flavor,
proximityBonding);
} catch (RDKit::FileParseException &e) {
BOOST_LOG(rdWarningLog) << e.message() << std::endl;
} catch (...) {
}
return static_cast<ROMol *>(newM);
}
ROMol *MolFromSequence(python::object seq, bool sanitize, int flavor) {
RWMol *newM = nullptr;
try {
newM = SequenceToMol(pyObjectToString(seq), sanitize, flavor);
} catch (RDKit::FileParseException &e) {
BOOST_LOG(rdWarningLog) << e.message() << std::endl;
} catch (...) {
}
return static_cast<ROMol *>(newM);
}
ROMol *MolFromFASTA(python::object seq, bool sanitize, int flavor) {
RWMol *newM = nullptr;
try {
newM = FASTAToMol(pyObjectToString(seq), sanitize, flavor);
} catch (RDKit::FileParseException &e) {
BOOST_LOG(rdWarningLog) << e.message() << std::endl;
} catch (...) {
}
return static_cast<ROMol *>(newM);
}
ROMol *MolFromHELM(python::object seq, bool sanitize) {
RWMol *newM = nullptr;
try {
newM = HELMToMol(pyObjectToString(seq), sanitize);
} catch (RDKit::FileParseException &e) {
BOOST_LOG(rdWarningLog) << e.message() << std::endl;
} catch (...) {
}
return static_cast<ROMol *>(newM);
}
std::string MolFragmentToSmilesHelper(
const ROMol &mol, python::object atomsToUse, python::object bondsToUse,
python::object atomSymbols, python::object bondSymbols,
bool doIsomericSmiles, bool doKekule, int rootedAtAtom, bool canonical,
bool allBondsExplicit, bool allHsExplicit) {
std::unique_ptr<std::vector<int>> avect =
pythonObjectToVect(atomsToUse, static_cast<int>(mol.getNumAtoms()));
if (!avect.get() || !(avect->size())) {
throw_value_error("atomsToUse must not be empty");
}
std::unique_ptr<std::vector<int>> bvect =
pythonObjectToVect(bondsToUse, static_cast<int>(mol.getNumBonds()));
std::unique_ptr<std::vector<std::string>> asymbols =
pythonObjectToVect<std::string>(atomSymbols);
std::unique_ptr<std::vector<std::string>> bsymbols =
pythonObjectToVect<std::string>(bondSymbols);
if (asymbols.get() && asymbols->size() != mol.getNumAtoms()) {
throw_value_error("length of atom symbol list != number of atoms");
}
if (bsymbols.get() && bsymbols->size() != mol.getNumBonds()) {
throw_value_error("length of bond symbol list != number of bonds");
}
std::string res = MolFragmentToSmiles(
mol, *avect.get(), bvect.get(), asymbols.get(), bsymbols.get(),
doIsomericSmiles, doKekule, rootedAtAtom, canonical, allBondsExplicit,
allHsExplicit);
return res;
}
std::vector<unsigned int> CanonicalRankAtoms(const ROMol &mol,
bool breakTies = true,
bool includeChirality = true,
bool includeIsotopes = true) {
std::vector<unsigned int> ranks(mol.getNumAtoms());
Canon::rankMolAtoms(mol, ranks, breakTies, includeChirality, includeIsotopes);
return ranks;
}
std::vector<int> CanonicalRankAtomsInFragment(const ROMol &mol,
python::object atomsToUse,
python::object bondsToUse,
python::object atomSymbols,
python::object bondSymbols,
bool breakTies = true)
{
std::unique_ptr<std::vector<int>> avect =
pythonObjectToVect(atomsToUse, static_cast<int>(mol.getNumAtoms()));
if (!avect.get() || !(avect->size())) {
throw_value_error("atomsToUse must not be empty");
}
std::unique_ptr<std::vector<int>> bvect =
pythonObjectToVect(bondsToUse, static_cast<int>(mol.getNumBonds()));
std::unique_ptr<std::vector<std::string>> asymbols =
pythonObjectToVect<std::string>(atomSymbols);
std::unique_ptr<std::vector<std::string>> bsymbols =
pythonObjectToVect<std::string>(bondSymbols);
if (asymbols.get() && asymbols->size() != mol.getNumAtoms()) {
throw_value_error("length of atom symbol list != number of atoms");
}
if (bsymbols.get() && bsymbols->size() != mol.getNumBonds()) {
throw_value_error("length of bond symbol list != number of bonds");
}
boost::dynamic_bitset<> atoms(mol.getNumAtoms());
for (size_t i = 0; i < avect->size(); ++i) atoms[(*avect)[i]] = true;
boost::dynamic_bitset<> bonds(mol.getNumBonds());
for (size_t i = 0; bvect.get() && i < bvect->size(); ++i)
bonds[(*bvect)[i]] = true;
std::vector<unsigned int> ranks(mol.getNumAtoms());
Canon::rankFragmentAtoms(mol, ranks, atoms, bonds, asymbols.get(),
bsymbols.get(), breakTies);
std::vector<int> resRanks(mol.getNumAtoms());
// set unused ranks to -1 for the Python interface
for (size_t i = 0; i < atoms.size(); ++i) {
if (!atoms[i]) {
resRanks[i] = -1;
} else {
resRanks[i] = static_cast<int>(ranks[i]);
}
}
return resRanks;
}
ROMol *MolFromSmilesHelper(python::object ismiles,
const SmilesParserParams ¶ms) {
std::string smiles = pyObjectToString(ismiles);
return SmilesToMol(smiles, params);
}
} // namespace RDKit
// MolSupplier stuff
#ifdef SUPPORT_COMPRESSED_SUPPLIERS
void wrap_compressedsdsupplier();
#endif
void wrap_sdsupplier();
void wrap_forwardsdsupplier();
void wrap_tdtsupplier();
void wrap_smisupplier();
#ifdef RDK_BUILD_COORDGEN_SUPPORT
void wrap_maesupplier();
#endif
// mol writer stuff
void wrap_smiwriter();
void wrap_sdwriter();
void wrap_tdtwriter();
void wrap_pdbwriter();
BOOST_PYTHON_MODULE(rdmolfiles) {
std::string docString;
python::scope().attr("__doc__") =
"Module containing RDKit functionality for working with molecular file "
"formats.";
python::register_exception_translator<RDKit::BadFileException>(
&rdBadFileExceptionTranslator);
docString =
"Construct a molecule from a TPL file.\n\n\
ARGUMENTS:\n\
\n\
- fileName: name of the file to read\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to True.\n\
\n\
- skipFirstConf: (optional) skips reading the first conformer.\n\
Defaults to False.\n\
This should be set to True when reading TPLs written by \n\
the CombiCode.\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def("MolFromTPLFile", RDKit::MolFromTPLFile,
(python::arg("fileName"), python::arg("sanitize") = true,
python::arg("skipFirstConf") = false),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Construct a molecule from a TPL block.\n\n\
ARGUMENTS:\n\
\n\
- fileName: name of the file to read\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to True.\n\
\n\
- skipFirstConf: (optional) skips reading the first conformer.\n\
Defaults to False.\n\
This should be set to True when reading TPLs written by \n\
the CombiCode.\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def("MolFromTPLBlock", RDKit::MolFromTPLBlock,
(python::arg("tplBlock"), python::arg("sanitize") = true,
python::arg("skipFirstConf") = false),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Construct a molecule from a Mol file.\n\n\
ARGUMENTS:\n\
\n\
- fileName: name of the file to read\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to true.\n\
\n\
- removeHs: (optional) toggles removing hydrogens from the molecule.\n\
This only make sense when sanitization is done.\n\
Defaults to true.\n\
\n\
- strictParsing: (optional) if this is false, the parser is more lax about.\n\
correctness of the content.\n\
Defaults to true.\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def(
"MolFromMolFile", RDKit::MolFromMolFile,
(python::arg("molFileName"), python::arg("sanitize") = true,
python::arg("removeHs") = true, python::arg("strictParsing") = true),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Construct a molecule from a Mol block.\n\n\
ARGUMENTS:\n\
\n\
- molBlock: string containing the Mol block\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to True.\n\
\n\
- removeHs: (optional) toggles removing hydrogens from the molecule.\n\
This only make sense when sanitization is done.\n\
Defaults to true.\n\
\n\
- strictParsing: (optional) if this is false, the parser is more lax about.\n\
correctness of the content.\n\
Defaults to true.\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def(
"MolFromMolBlock", RDKit::MolFromMolBlock,
(python::arg("molBlock"), python::arg("sanitize") = true,
python::arg("removeHs") = true, python::arg("strictParsing") = true),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Construct a molecule from an RDKit-generate SVG string.\n\n\
ARGUMENTS:\n\
\n\
- svg: string containing the SVG data (must include molecule metadata)\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to True.\n\
\n\
- removeHs: (optional) toggles removing hydrogens from the molecule.\n\
This only make sense when sanitization is done.\n\
Defaults to true.\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n\
NOTE: this functionality should be considered beta.\n\
\n";
python::def("MolFromRDKitSVG", RDKit::MolFromSVG,
(python::arg("molBlock"), python::arg("sanitize") = true,
python::arg("removeHs") = true),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Construct a molecule from a Tripos Mol2 file.\n\n\
NOTE:\n \
The parser expects the atom-typing scheme used by Corina.\n\
Atom types from Tripos' dbtranslate are less supported.\n\
Other atom typing schemes are unlikely to work.\n\
\n\
ARGUMENTS:\n \
\n\
- fileName: name of the file to read\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to true.\n\
\n\
- removeHs: (optional) toggles removing hydrogens from the molecule.\n\
This only make sense when sanitization is done.\n\
Defaults to true.\n\
\n\
- cleanupSubstructures: (optional) toggles standardizing some \n\
substructures found in mol2 files.\n\
Defaults to true.\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def("MolFromMol2File", RDKit::MolFromMol2File,
(python::arg("molFileName"), python::arg("sanitize") = true,
python::arg("removeHs") = true,
python::arg("cleanupSubstructures") = true),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Construct a molecule from a Tripos Mol2 block.\n\n\
NOTE:\n \
The parser expects the atom-typing scheme used by Corina.\n\
Atom types from Tripos' dbtranslate are less supported.\n\
Other atom typing schemes are unlikely to work.\n\
\n\
ARGUMENTS:\n\
\n\
- mol2Block: string containing the Mol2 block\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to True.\n\
\n\
- removeHs: (optional) toggles removing hydrogens from the molecule.\n\
This only make sense when sanitization is done.\n\
Defaults to true.\n\
\n\
- cleanupSubstructures: (optional) toggles standardizing some \n\
substructures found in mol2 files.\n\
Defaults to true.\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def("MolFromMol2Block", RDKit::MolFromMol2Block,
(python::arg("molBlock"), python::arg("sanitize") = true,
python::arg("removeHs") = true,
python::arg("cleanupSubstructures") = true),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Construct a molecule from a Mol file.\n\n\
ARGUMENTS:\n\
\n\
- fileName: name of the file to read\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to true.\n\
\n\
- removeHs: (optional) toggles removing hydrogens from the molecule.\n\
This only make sense when sanitization is done.\n\
Defaults to true.\n\
\n\
- strictParsing: (optional) if this is false, the parser is more lax about.\n\
correctness of the content.\n\
Defaults to true.\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def(
"MolFromMolFile", RDKit::MolFromMolFile,
(python::arg("molFileName"), python::arg("sanitize") = true,
python::arg("removeHs") = true, python::arg("strictParsing") = true),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Construct a molecule from a Mol block.\n\n\
ARGUMENTS:\n\
\n\
- molBlock: string containing the Mol block\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to True.\n\
\n\
- removeHs: (optional) toggles removing hydrogens from the molecule.\n\
This only make sense when sanitization is done.\n\
Defaults to true.\n\
\n\
- strictParsing: (optional) if this is false, the parser is more lax about.\n\
correctness of the content.\n\
Defaults to true.\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def(
"MolFromMolBlock", RDKit::MolFromMolBlock,
(python::arg("molBlock"), python::arg("sanitize") = true,
python::arg("removeHs") = true, python::arg("strictParsing") = true),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Returns a Mol block for a molecule\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
- includeStereo: (optional) toggles inclusion of stereochemical\n\
information in the output\n\
- confId: (optional) selects which conformation to output (-1 = default)\n\
- kekulize: (optional) triggers kekulization of the molecule before it's written,\n\
as suggested by the MDL spec.\n\
- forceV3000 (optional) force generation a V3000 mol block (happens automatically with \n\
more than 999 atoms or bonds)\n\
\n\
RETURNS:\n\
\n\
a string\n\
\n";
python::def("MolToMolBlock", RDKit::MolToMolBlock,
(python::arg("mol"), python::arg("includeStereo") = true,
python::arg("confId") = -1, python::arg("kekulize") = true,
python::arg("forceV3000") = false),
docString.c_str());
docString =
"Writes a Mol file for a molecule\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
- filename: the file to write to\n\
- includeStereo: (optional) toggles inclusion of stereochemical\n\
information in the output\n\
- confId: (optional) selects which conformation to output (-1 = default)\n\
- kekulize: (optional) triggers kekulization of the molecule before it's written,\n\
as suggested by the MDL spec.\n\
- forceV3000 (optional) force generation a V3000 mol block (happens automatically with \n\
more than 999 atoms or bonds)\n\
\n\
RETURNS:\n\
\n\
a string\n\
\n";
python::def(
"MolToMolFile", RDKit::MolToMolFile,
(python::arg("mol"), python::arg("filename"),
python::arg("includeStereo") = true, python::arg("confId") = -1,
python::arg("kekulize") = true, python::arg("forceV3000") = false),
docString.c_str());
python::class_<RDKit::SmilesParserParams, boost::noncopyable>(
"SmilesParserParams", "Parameters controlling SMILES Parsing")
.def_readwrite("maxIterations", &RDKit::SmilesParserParams::debugParse,
"controls the amount of debugging information produced")
.def_readwrite("parseName", &RDKit::SmilesParserParams::parseName,
"controls whether or not the molecule name is also parsed")
.def_readwrite(
"allowCXSMILES", &RDKit::SmilesParserParams::allowCXSMILES,
"controls whether or not the CXSMILES extensions are parsed")
.def_readwrite("sanitize", &RDKit::SmilesParserParams::sanitize,
"controls whether or not the molecule is sanitized before "
"being returned")
.def_readwrite("removeHs", &RDKit::SmilesParserParams::removeHs,
"controls whether or not Hs are removed before the "
"molecule is returned");
docString =
"Construct a molecule from a SMILES string.\n\n\
ARGUMENTS:\n\
\n\
- SMILES: the smiles string\n\
\n\
- params: used to provide optional parameters for the SMILES parsing\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def("MolFromSmiles", MolFromSmilesHelper,
(python::arg("SMILES"), python::arg("params")), docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Construct a molecule from a SMILES string.\n\n\
ARGUMENTS:\n\
\n\
- SMILES: the smiles string\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to True.\n\
\n\
- replacements: (optional) a dictionary of replacement strings (see below)\n\
Defaults to {}.\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n\
The optional replacements dict can be used to do string substitution of abbreviations \n\
in the input SMILES. The set of substitutions is repeatedly looped through until \n\
the string no longer changes. It is the responsiblity of the caller to make sure \n\
that substitutions results in legal and sensible SMILES. \n\
\n\
Examples of replacements: \n\
\n\
CC{Q}C with {'{Q}':'OCCO'} -> CCOCCOC \n\
C{A}C{Q}C with {'{Q}':'OCCO', '{A}':'C1(CC1)'} -> CC1(CC1)COCCOC \n\
C{A}C{Q}C with {'{Q}':'{X}CC{X}', '{A}':'C1CC1', '{X}':'N'} -> CC1CC1CCNCCNC \n\
\n";
python::def("MolFromSmiles", RDKit::MolFromSmiles,
(python::arg("SMILES"), python::arg("sanitize") = true,
python::arg("replacements") = python::dict()),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString = "Construct an atom from a SMILES string";
python::def("AtomFromSmiles", SmilesToAtom, python::arg("SMILES"),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString = "Construct a bond from a SMILES string";
python::def("BondFromSmiles", SmilesToBond, python::arg("SMILES"),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Construct a molecule from a SMARTS string.\n\n\
ARGUMENTS:\n\
\n\
- SMARTS: the smarts string\n\
\n\
- mergeHs: (optional) toggles the merging of explicit Hs in the query into the attached\n\
atoms. So, for example, 'C[H]' becomes '[C;!H0]'.\n\
Defaults to 0.\n\
\n\
- replacements: (optional) a dictionary of replacement strings (see below)\n\
Defaults to {}. See the documentation for MolFromSmiles for an explanation.\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def("MolFromSmarts", RDKit::MolFromSmarts,
(python::arg("SMARTS"), python::arg("mergeHs") = false,
python::arg("replacements") = python::dict()),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString = "Construct an atom from a SMARTS string";
python::def("AtomFromSmarts", SmartsToAtom, python::arg("SMARTS"),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString = "Construct a bond from a SMARTS string";
python::def("BondFromSmarts", SmartsToBond, python::arg("SMILES"),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Returns the canonical SMILES string for a molecule\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
- isomericSmiles: (optional) include information about stereochemistry in\n\
the SMILES. Defaults to true.\n\
- kekuleSmiles: (optional) use the Kekule form (no aromatic bonds) in\n\
the SMILES. Defaults to false.\n\
- rootedAtAtom: (optional) if non-negative, this forces the SMILES \n\
to start at a particular atom. Defaults to -1.\n\
- canonical: (optional) if false no attempt will be made to canonicalize\n\
the molecule. Defaults to true.\n\
- allBondsExplicit: (optional) if true, all bond orders will be explicitly indicated\n\
in the output SMILES. Defaults to false.\n\
- allHsExplicit: (optional) if true, all H counts will be explicitly indicated\n\
in the output SMILES. Defaults to false.\n\
\n\
RETURNS:\n\
\n\
a string\n\
\n";
python::def(
"MolToSmiles", RDKit::MolToSmiles,
(python::arg("mol"), python::arg("isomericSmiles") = true,
python::arg("kekuleSmiles") = false, python::arg("rootedAtAtom") = -1,
python::arg("canonical") = true, python::arg("allBondsExplicit") = false,
python::arg("allHsExplicit") = false, python::arg("doRandom") = false),
docString.c_str());
docString =
"Returns the canonical SMILES string for a fragment of a molecule\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
- atomsToUse : a list of atoms to include in the fragment\n\
- bondsToUse : (optional) a list of bonds to include in the fragment\n\
if not provided, all bonds between the atoms provided\n\
will be included.\n\
- atomSymbols : (optional) a list with the symbols to use for the atoms\n\
in the SMILES. This should have be mol.GetNumAtoms() long.\n\
- bondSymbols : (optional) a list with the symbols to use for the bonds\n\
in the SMILES. This should have be mol.GetNumBonds() long.\n\
- isomericSmiles: (optional) include information about stereochemistry in\n\
the SMILES. Defaults to true.\n\
- kekuleSmiles: (optional) use the Kekule form (no aromatic bonds) in\n\
the SMILES. Defaults to false.\n\
- rootedAtAtom: (optional) if non-negative, this forces the SMILES \n\
to start at a particular atom. Defaults to -1.\n\
- canonical: (optional) if false no attempt will be made to canonicalize\n\
the molecule. Defaults to true.\n\
- allBondsExplicit: (optional) if true, all bond orders will be explicitly indicated\n\
in the output SMILES. Defaults to false.\n\
- allHsExplicit: (optional) if true, all H counts will be explicitly indicated\n\
in the output SMILES. Defaults to false.\n\
- doRandom: (optional) if true, randomized the DFS transversal graph,\n\
so we can generate random smiles. Defaults to false.\n\
\n\
RETURNS:\n\
\n\
a string\n\
\n";
python::def(
"MolFragmentToSmiles", MolFragmentToSmilesHelper,
(python::arg("mol"), python::arg("atomsToUse"),
python::arg("bondsToUse") = 0, python::arg("atomSymbols") = 0,
python::arg("bondSymbols") = 0, python::arg("isomericSmiles") = true,
python::arg("kekuleSmiles") = false, python::arg("rootedAtAtom") = -1,
python::arg("canonical") = true, python::arg("allBondsExplicit") = false,
python::arg("allHsExplicit") = false),
docString.c_str());
docString =
"Returns a SMARTS string for a molecule\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
- isomericSmiles: (optional) include information about stereochemistry in\n\
the SMARTS. Defaults to true.\n\
\n\
RETURNS:\n\
\n\
a string\n\
\n";
python::def("MolToSmarts", RDKit::MolToSmarts,
(python::arg("mol"), python::arg("isomericSmiles") = true),
docString.c_str());
docString =
"Writes a molecule to a TPL file.\n\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
- fileName: name of the file to write\n\
- partialChargeProp: name of the property to use for partial charges\n\
Defaults to '_GasteigerCharge'.\n\
- writeFirstConfTwice: Defaults to False.\n\
This should be set to True when writing TPLs to be read by \n\
the CombiCode.\n\
\n";
python::def("MolToTPLFile", RDKit::MolToTPLFile,
(python::arg("mol"), python::arg("fileName"),
python::arg("partialChargeProp") = "_GasteigerCharge",
python::arg("writeFirstConfTwice") = false),
docString.c_str());
docString =
"Returns the Tpl block for a molecule.\n\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
- partialChargeProp: name of the property to use for partial charges\n\
Defaults to '_GasteigerCharge'.\n\
- writeFirstConfTwice: Defaults to False.\n\
This should be set to True when writing TPLs to be read by \n\
the CombiCode.\n\
\n\
RETURNS:\n\
\n\
a string\n\
\n";
python::def("MolToTPLBlock", RDKit::MolToTPLText,
(python::arg("mol"),
python::arg("partialChargeProp") = "_GasteigerCharge",
python::arg("writeFirstConfTwice") = false),
docString.c_str());
docString =
"Construct a molecule from a PDB file.\n\n\
ARGUMENTS:\n\
\n\
- fileName: name of the file to read\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to true.\n\
\n\
- removeHs: (optional) toggles removing hydrogens from the molecule.\n\
This only make sense when sanitization is done.\n\
Defaults to true.\n\
\n\
- flavor: (optional) \n\
\n\
- proximityBonding: (optional) toggles automatic proximity bonding\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def("MolFromPDBFile", RDKit::MolFromPDBFile,
(python::arg("molFileName"), python::arg("sanitize") = true,
python::arg("removeHs") = true, python::arg("flavor") = 0,
python::arg("proximityBonding") = true),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Construct a molecule from a PDB block.\n\n\
ARGUMENTS:\n\
\n\
- molBlock: string containing the PDB block\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to True.\n\
\n\
- removeHs: (optional) toggles removing hydrogens from the molecule.\n\
This only make sense when sanitization is done.\n\
Defaults to true.\n\
\n\
- flavor: (optional) \n\
\n\
- proximityBonding: (optional) toggles automatic proximity bonding\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def("MolFromPDBBlock", RDKit::MolFromPDBBlock,
(python::arg("molBlock"), python::arg("sanitize") = true,
python::arg("removeHs") = true, python::arg("flavor") = 0,
python::arg("proximityBonding") = true),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Returns a PDB block for a molecule\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
- confId: (optional) selects which conformation to output (-1 = default)\n\
- flavor: (optional) \n\
flavor & 1 : Write MODEL/ENDMDL lines around each record \n\
flavor & 2 : Don't write any CONECT records \n\
flavor & 4 : Write CONECT records in both directions \n\
flavor & 8 : Don't use multiple CONECTs to encode bond order \n\
flavor & 16 : Write MASTER record \n\
flavor & 32 : Write TER record \n\
\n\
RETURNS:\n\
\n\
a string\n\
\n";
python::def("MolToPDBBlock", RDKit::MolToPDBBlock,
(python::arg("mol"), python::arg("confId") = -1,
python::arg("flavor") = 0),
docString.c_str());
docString =
"Writes a PDB file for a molecule\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
- filename: name of the file to write\n\
- confId: (optional) selects which conformation to output (-1 = default)\n\
- flavor: (optional) \n\
flavor & 1 : Write MODEL/ENDMDL lines around each record \n\
flavor & 2 : Don't write any CONECT records \n\
flavor & 4 : Write CONECT records in both directions \n\
flavor & 8 : Don't use multiple CONECTs to encode bond order \n\
flavor & 16 : Write MASTER record \n\
flavor & 32 : Write TER record \n\
\n\
RETURNS:\n\
\n\
a string\n\
\n";
python::def("MolToPDBFile", RDKit::MolToPDBFile,
(python::arg("mol"), python::arg("filename"),
python::arg("confId") = -1, python::arg("flavor") = 0),
docString.c_str());
docString =
"Construct a molecule from a sequence string (currently only supports peptides).\n\n\
ARGUMENTS:\n\
\n\
- text: string containing the sequence\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to True.\n\
\n\
- flavor: (optional)\n\
0 Protein, L amino acids (default)\n\
1 Protein, D amino acids\n\
2 RNA, no cap\n\
3 RNA, 5' cap\n\
4 RNA, 3' cap\n\
5 RNA, both caps\n\
6 DNA, no cap\n\
7 DNA, 5' cap\n\
8 DNA, 3' cap\n\
9 DNA, both caps\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def("MolFromSequence", RDKit::MolFromSequence,
(python::arg("text"), python::arg("sanitize") = true,
python::arg("flavor") = 0),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Returns the sequence string for a molecule\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
\n\
NOTE: the molecule should contain monomer information in AtomMonomerInfo structures \n\
\n\
RETURNS:\n\
\n\
a string\n\
\n";
python::def("MolToSequence", RDKit::MolToSequence, (python::arg("mol")),
docString.c_str());
docString =
"Construct a molecule from a FASTA string (currently only supports peptides).\n\n\
ARGUMENTS:\n\
\n\
- text: string containing the FASTA\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to True.\n\
\n\
- flavor: (optional)\n\
0 Protein, L amino acids (default)\n\
1 Protein, D amino acids\n\
2 RNA, no cap\n\
3 RNA, 5' cap\n\
4 RNA, 3' cap\n\
5 RNA, both caps\n\
6 DNA, no cap\n\
7 DNA, 5' cap\n\
8 DNA, 3' cap\n\
9 DNA, both caps\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def("MolFromFASTA", RDKit::MolFromFASTA,
(python::arg("text"), python::arg("sanitize") = true,
python::arg("flavor") = 0),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Returns the FASTA string for a molecule\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
\n\
NOTE: the molecule should contain monomer information in AtomMonomerInfo structures \n\
\n\
RETURNS:\n\
\n\
a string\n\
\n";
python::def("MolToFASTA", RDKit::MolToFASTA, (python::arg("mol")),
docString.c_str());
docString =
"Construct a molecule from a HELM string (currently only supports peptides).\n\n\
ARGUMENTS:\n\
\n\
- text: string containing the HELM\n\
\n\
- sanitize: (optional) toggles sanitization of the molecule.\n\
Defaults to true.\n\
\n\
RETURNS:\n\
\n\
a Mol object, None on failure.\n\
\n";
python::def("MolFromHELM", RDKit::MolFromHELM,
(python::arg("text"), python::arg("sanitize") = true),
docString.c_str(),
python::return_value_policy<python::manage_new_object>());
docString =
"Returns the HELM string for a molecule\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
\n\
NOTE: the molecule should contain monomer information in AtomMonomerInfo structures \n\
\n\
RETURNS:\n\
\n\
a string\n\
\n";
python::def("MolToHELM", RDKit::MolToHELM, (python::arg("mol")),
docString.c_str());
docString =
"Returns the canonical atom ranking for each atom of a molecule fragment.\n\
If breakTies is False, this returns the symmetry class for each atom. The symmetry\n\
class is used by the canonicalization routines to type each atom based on the whole\n\
chemistry of the molecular graph. Any atom with the same rank (symmetry class) is\n\
indistinguishable. For example:\n\
\n\
> mol = MolFromSmiles('C1NCN1')\n\
> list(CanonicalRankAtoms(mol, breakTies=False))\n\
[0,1,0,1]\n\
\n\
In this case the carbons have the same symmetry class and the nitrogens have the same\n\
symmetry class. From the perspective of the Molecular Graph, they are indentical.\n\
\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
- breakTies: (optional) force breaking of ranked ties [default=True]\n\
- includeChirality: (optional) use chiral information when computing rank [default=True]\n\
- includeIsotopes: (optional) use isotope information when computing rank [default=True]\n\
\n\
RETURNS:\n\
\n\
a string\n\
\n";
python::def("CanonicalRankAtoms", CanonicalRankAtoms,
(python::arg("mol"), python::arg("breakTies") = true,
python::arg("includeChirality") = true,
python::arg("includeIsotopes") = true),
docString.c_str());
docString =
"Returns the canonical atom ranking for each atom of a molecule fragment\n\
See help(CanonicalRankAtoms) for more information.\n\
\n\
> mol = MolFromSmiles('C1NCN1.C1NCN1')\n\
> list(CanonicalRankAtomsInFragment(mol, atomsToUse=range(0,4), breakTies=False))\n\
[0,1,0,1,-1,-1,-1,-1]\n\
> list(CanonicalRankAtomsInFragment(mol, atomsToUse=range(4,8), breakTies=False))\n\
[-1,-1,-1,-1,0,1,0,1]\n\
\n\
ARGUMENTS:\n\
\n\
- mol: the molecule\n\
- atomsToUse : a list of atoms to include in the fragment\n\
- bondsToUse : (optional) a list of bonds to include in the fragment\n\
if not provided, all bonds between the atoms provided\n\
will be included.\n\
- atomSymbols : (optional) a list with the symbols to use for the atoms\n\
in the SMILES. This should have be mol.GetNumAtoms() long.\n\
- bondSymbols : (optional) a list with the symbols to use for the bonds\n\
in the SMILES. This should have be mol.GetNumBonds() long.\n\
- breakTies: (optional) force breaking of ranked ties\n\
\n\
RETURNS:\n\
\n\
a string\n\
\n";
python::def("CanonicalRankAtomsInFragment", CanonicalRankAtomsInFragment,
(python::arg("mol"), python::arg("atomsToUse"),
python::arg("bondsToUse") = 0, python::arg("atomSymbols") = 0,
python::arg("bondSymbols") = 0, python::arg("breakTies") = true),
docString.c_str());
/********************************************************
* MolSupplier stuff
*******************************************************/
#ifdef SUPPORT_COMPRESSED_SUPPLIERS
wrap_compressedsdsupplier();
#endif
wrap_sdsupplier();
wrap_forwardsdsupplier();
wrap_tdtsupplier();
wrap_smisupplier();
#ifdef RDK_BUILD_COORDGEN_SUPPORT
wrap_maesupplier();
#endif
// wrap_pdbsupplier();
/********************************************************
* MolWriter stuff
*******************************************************/
wrap_smiwriter();
wrap_sdwriter();
wrap_tdtwriter();
wrap_pdbwriter();
}
|