# file: generate_molmodel_source.py
# Run this script from the above src directory

import os
import sys
import re
from pyplusplus import module_builder
from pyplusplus.module_builder.call_policies import *
from simtk_wrap_utils import find_path, find_program
from simtk_wrap_utils import *
from pyplusplus import function_transformers as FT
from doxygen_doc_extractor import doxygen_doc_extractor

def generate_molmodel_source():
    """
    Use pyplusplus to generate boost.python source code for simtk molmodel module.
    See pyplusplus documentation at 
        http://www.language-binding.net/pyplusplus/pyplusplus.html
    """
    mb = create_module_builder(['molmodel/wrap_molmodel.h',])

    # Delegate wrapping details to other subroutines, to keep this subrouting from growing too large
    mb.namespace('SimTK').include()
    perform_usual_wrapping_tasks(mb)
    # disambiguate_atomindex(mb)
    disambiguate_typedef_aliases(mb)
    expose_amber99(mb)
    expose_clonable(mb)
    expose_compound(mb)
    expose_biopolymerresidue(mb)
    expose_residueinfo(mb)
    exclude_compoundmodeler(mb)
    expose_biotype(mb)
    exclude_grinpointers(mb)
    exclude_matrixhelpers(mb)
    expose_indexes(mb)
    expose_pdb(mb)
    expose_dumm(mb)
    expose_eventhandlers(mb)
    wrap_periodicpdbwriter(mb)
    wrap_periodicvmdreporter(mb)
    wrap_compoundsystem(mb)
    wrap_protein(mb)
    
    wrap_adaptbxostream(mb)
    exclude_rep_classes(mb)
    exclude_operators(mb)
    exclude_strange_classes(mb)
    exclude_isinstanceof(mb)
    wrap_properties(mb)

    # Don't rewrap anything already wrapped by simtk.common etc.
    # See http://www.language-binding.net/pyplusplus/documentation/multi_module_development.html
    mb.register_module_dependency('external/std/generated_code/')
    mb.register_module_dependency('simtkvecs/generated_code/')
    mb.register_module_dependency('simtkcommon/generated_code/')
    mb.register_module_dependency('simtkmath/generated_code/')
    mb.register_module_dependency('simbody/generated_code/')
    
    extractor = doxygen_doc_extractor()
    
    mb.build_code_creator(module_name='_molmodel', doc_extractor=extractor)
    mb.split_module(os.path.join(os.path.abspath('.'), 'molmodel', 'generated_code'))
    # If all succeeds, record this accomplishment by touching a particular file
    open(os.path.join(os.path.abspath('.'), 'molmodel', 'generated_code', 'generate_molmodel.stamp'), "w").close()

    # Create empty header file pyplusplus mistakenly thinks it needs, for value_traits specified in
    # simtk_indexing_helpers.hpp
    # 2>src\molmodel\generated_code\AtomTargetLocations.pypp.cpp(7) : fatal error C1083: Cannot open include file: '_Vec_less__3_comma__double_comma__1__greater___value_traits.pypp.hpp': No such file or directory
    # open(os.path.join(os.path.abspath('.'), 'molmodel', 'generated_code', '_Vec_less__3_comma__double_comma__1__greater___value_traits.pypp.hpp'), "a").close()    
    for fname in ['_Vec3__value_traits.pypp.hpp']:
        open(os.path.join(os.path.abspath('.'), 'molmodel', 'generated_code', fname), "a").close()    

def wrap_protein(mb):
    # Need constructor that takes a python file
    protein = mb.class_('Protein')
    protein.include_files.append('python_streambuf.h')
    protein.add_declaration_code("""//
        using boost_adaptbx::python::streambuf;
        static boost::shared_ptr<SimTK::Protein> makeProtein1(streambuf::istream& is) 
        {
            return boost::shared_ptr<SimTK::Protein>(new SimTK::Protein(is));
        }
        static boost::shared_ptr<SimTK::Protein> makeProtein2(bp::object& pyfile)
        {
            streambuf input(pyfile);
            streambuf::istream is(input);
            return boost::shared_ptr<SimTK::Protein>(new SimTK::Protein(is));
        }
        """)
    protein.add_registration_code('def("__init__", bp::make_constructor(makeProtein1))')
    protein.add_registration_code('def("__init__", bp::make_constructor(makeProtein2))', tail=False) # must be declared first

def wrap_compoundsystem(mb):
    system = mb.class_('CompoundSystem')
    system.member_function('adoptCompound').call_policies = with_custodian_and_ward(1, 2)

def wrap_periodicvmdreporter(mb):
    # PeriodicPdbWriter needs wrappers to pass python file objects as arguments in place of ostreams
    vmd = mb.class_('PeriodicVmdReporter')
    vmd.held_type = 'std::auto_ptr< %s >' % vmd.wrapper_alias
    vmd_ptr_type = 'std::auto_ptr< %s >' % vmd.decl_string
    vmd.add_registration_code(
        'bp::implicitly_convertible<%s, %s >();' % (vmd.held_type, vmd_ptr_type),
        works_on_instance=False)
    # Permit implicit conversion to immediate base class
    parent = vmd.bases[0].related_class
    vmd.add_registration_code(
        'bp::implicitly_convertible<%s, std::auto_ptr< %s > >();' % (
            vmd_ptr_type, parent.decl_string),
        works_on_instance=False)

def wrap_periodicpdbwriter(mb):
    # PeriodicPdbWriter needs wrappers to pass python file objects as arguments in place of ostreams
    ppw = mb.class_('PeriodicPdbWriter')
    ppw.include_files.append('python_streambuf.h')
    ppw.add_declaration_code("""//
        using boost_adaptbx::python::streambuf;
        // Need to return a pointer, not a boost::shared_ptr, because held_type
        // of PeriodicPdbWriter will be std::auto_ptr<>
        static SimTK::PeriodicPdbWriter* makePeriodicPdbWriter(
                const CompoundSystem& system, streambuf::ostream& os, Real interval) 
        {
            return new SimTK::PeriodicPdbWriter(system, os, interval);
        }""")
    ppw.add_registration_code("""
        def("__init__", 
            bp::make_constructor(makePeriodicPdbWriter)
                // , ( bp::arg("system"), bp::arg("os"), bp::arg("interval") )
                // ,  bp::with_custodian_and_ward_postcall< 0, 1, bp::with_custodian_and_ward_postcall< 0, 2 > >()
        )
    """, tail=False) # Must be declared first, for other constructors to be available
    ppw.held_type = 'std::auto_ptr< %s >' % ppw.wrapper_alias
    ppw_ptr_type = 'std::auto_ptr< %s >' % ppw.decl_string
    ppw.add_registration_code(
        'bp::implicitly_convertible<%s, %s >();' % (ppw.held_type, ppw_ptr_type),
        works_on_instance=False)
    # Permit implicit conversion to immediate base class
    parent = ppw.bases[0].related_class
    ppw.add_registration_code(
        'bp::implicitly_convertible<%s, std::auto_ptr< %s > >();' % (
            ppw_ptr_type, parent.decl_string),
        works_on_instance=False)

def wrap_adaptbxostream(mb):
    mb.add_declaration_code('#include "register_adaptbx_ostream.hpp"')
    mb.add_registration_code('register_adaptbx_ostream_class();', tail=False) # must register before Protein, PeriodicPdbWriter

def expose_eventhandlers(mb):
    # runtime error with unmodified VelocityRescalingThermostat
    for cls_name in ['VelocityRescalingThermostat', 'MassCenterMotionRemover']:
        cls = mb.class_(cls_name)
        cls.held_type = 'std::auto_ptr< %s >' % cls.wrapper_alias
        cls_ptr_type = 'std::auto_ptr< %s >' % cls.decl_string
        cls.add_registration_code(
            'bp::implicitly_convertible<%s, %s >();' % (cls.held_type, cls_ptr_type),
            works_on_instance=False)
        # Permit implicit conversion to all base classes
        parent = cls.bases[0].related_class
        cls.add_registration_code(
            'bp::implicitly_convertible<%s, std::auto_ptr< %s > >();' % (
                cls_ptr_type, parent.decl_string),
            works_on_instance=False)
        # Transform "handleEvent" method
        cls.member_function("handleEvent").add_transformation(
            FT.inout("lowestModified"))
        for ctor in cls.constructors():
            # don't want copy constructor or default constructor
            if len(ctor.argument_types) < 2:
                continue
            # First argument should be a system
            if ctor.argument_types[0].decl_string.find("System") >= 0:
                ctor.call_policies = with_custodian_and_ward(1, 2)        

def expose_dumm(mb):
    dumm = mb.class_('DuMMForceFieldSubsystem')
    #1>WARNING: SimTK::DuMMForceFieldSubsystem [class]
    #1>> warning W1031: `Py++` will generate class wrapper - user asked to
    #1>> expose non - public member function "defineIncompleteAtomClass"
    for fn in dumm.member_functions(lambda f: f.access_type != declarations.ACCESS_TYPES.PUBLIC):
        fn.exclude()
    system_t = declarations.declarated_t(mb.class_('MolecularMechanicsSystem'))
    system_ref_t = declarations.reference_t(system_t)
    dumm.constructor(arg_types=[system_ref_t]).call_policies = with_custodian_and_ward(1, 2)
    # There is a mysterious runtime problem with defineAtomClass_KA
    # Use defineAtomClass() instead
    dumm.member_functions('defineAtomClass_KA').exclude()

def expose_pdb(mb):
    #1>WARNING: std::ostream & SimTK::PdbChain::write(std::ostream & os, int & nextAtomSerialNumber, SimTK::Transform const & transform=SimTK::Transform_<double>()) const [member function]
    #1>> execution error W1009: The function takes as argument
    #1>> (name=nextAtomSerialNumber, pos=1) non-const reference to Python
    #1>> immutable type - function could not be called from Python. Take a look
    #1>> on "Function Transformation" functionality and define the
    #1>> transformation.
    mb.class_('PdbChain').member_function('write').add_transformation(
            FT.inout('nextAtomSerialNumber'))
    #1>WARNING: SimTK::PdbChain [class]
    #1>> warning W1031: `Py++` will generate class wrapper - user asked to
    #1>> expose non - public member function "parsePdbLine"
    mb.class_('PdbChain').member_function('parsePdbLine').exclude()
    #
    expose_pdbatomlocation(mb)
    
def expose_pdbatomlocation(mb):
    loc = mb.class_('PdbAtomLocation')
    loc.include()
    # loc.alias = 'UnitlessPdbAtomLocation'
    # Wrap PDB writing code
    # std::ostream & writePdb(std::ostream &os, const Transform  &transform) const 
    loc.add_declaration_code("""
        // Create operator<< for boost.python to do its __str__ magic with, if it has any.
        namespace SimTK {
            std::ostream& operator<<(std::ostream& os, const SimTK::PdbAtomLocation& loc) {
                loc.writePdb(os, SimTK::Transform());
                return os;
            }
        }
        std::string loc_write_str(const SimTK::PdbAtomLocation& loc, const Transform& transform = Transform()) {
            std::ostringstream str;
            loc.writePdb(str, transform);
            return str.str();
        }
    """)
    loc.add_registration_code('def( bp::self_ns::str(bp::self) )')
    loc.add_registration_code('def( "writePdb", &loc_write_str )')
    
def exclude_isinstanceof(mb):
    # run time error during "from simbios.simtk._molmodel import *"
    # KeyError: 'isInstanceOf'
    for fn in mb.member_functions('isInstanceOf'):
        # Staticness somehow leaves remains in the generated wrapper code.
        fn.has_static = False
        fn.exclude()
    
def expose_residueinfo(mb):
    # link error
    mb.class_('ResidueInfo').constructor(arg_types=[None,None,None,None,None]).exclude()
    
def expose_biotype(mb):
    biotype = mb.class_('Biotype')
    biotype.member_function('SerineHB').exclude() # link error
    biotype.member_function('SerineCB').exclude() # link error

def expose_biopolymerresidue(mb):
    res = mb.class_('BiopolymerResidue')
    # link errors
    res.member_function('getParentCompoundAtomIndex').exclude()
    
def expose_compound(mb):
    # link errors
    cmpd = mb.class_('Compound')
    cmpd.member_function('updSubcompound').exclude()
    cmpd.member_function('getSubcompound').exclude()
    cmpd.member_function('getSubcompoundFrameInParentFrame').exclude()
    cmpd.member_function('calcAtomLocationInCompoundFrame').exclude()
    # Use transformer with FilestarOstream to convert python file object
    # to C++ ostream reference
    #for fn in cmpd.member_functions('writeDefaultPdb', arg_types=[None,None]):
        ## need to resolve ambiguity ostream vs. char* overloads
        #if declarations.is_reference(fn.argument_types[0]):
            #fn.add_transformation(FT.input_ostream("os"))
    #for fn in cmpd.member_functions('writePdb', arg_types=[None,None,None]):
        ## need to resolve ambiguity ostream vs. char* overloads
        #if declarations.is_reference(fn.argument_types[1]):
            #fn.add_transformation(FT.input_ostream("os"))
    ## adaptbx_ostream version
    cmpd.include_files.append("python_streambuf.h")
    cmpd.add_declaration_code("""
        using boost_adaptbx::python::streambuf;
        // returning ostream& causes crash...
        void cmpd_writeDefaultPdb_wrapper(const Compound& compound, bp::object& pyfile, const Transform& transform=Transform())
        {
            streambuf output(pyfile);
            streambuf::ostream os(output);
            compound.writeDefaultPdb(os, transform);
        }
        void cmpd_writePdb_wrapper(const Compound& compound, const State& state, bp::object& pyfile, const Transform& transform=Transform())
        {
            streambuf output(pyfile);
            streambuf::ostream os(output);
            compound.writePdb(state, os, transform);
        }
    """)
    cmpd.add_registration_code("""def("writeDefaultPdb", 
            &cmpd_writeDefaultPdb_wrapper,
            ( bp::arg("inst"), bp::arg("pyfile"), bp::arg("transform")=SimTK::Transform_<double>() ))
            """, tail=False)
    cmpd.add_registration_code("""def("writePdb", 
            &cmpd_writePdb_wrapper,
            ( bp::arg("inst"), bp::arg("state"), bp::arg("pyfile"), bp::arg("transform")=SimTK::Transform_<double>() ))
            """, tail=False)

def exclude_compoundmodeler(mb):
    # link errors
    mb.class_('CompoundModeler').exclude()

def exclude_grinpointers(mb):
    mb.classes(lambda c: c.name.startswith('GrinPointer<')).exclude()

def expose_clonable(mb):
    mb.class_('Clonable').member_function('clone').call_policies = return_value_policy(manage_new_object)

def expose_atomtargetlocations(mb):
    cls = mb.class_('AtomTargetLocations')
    cls.include_files.append("simtk_indexing_helpers.hpp")
    for inc in cls.include_files():
        print inc
        print dir(inc)

def disambiguate_atomindex(mb):
    count = 1
    for cls in mb.classes(lambda c: c.alias=='AtomIndex'):
        cls.alias = 'AtomIndex%d' % count
        count += 1

def expose_amber99(mb):
    amber = mb.class_('Amber99ForceSubsystem')
    mb.member_functions('updDowncast').exclude()
    mb.member_functions('downcast').exclude() # link error

if __name__ == "__main__":
    generate_molmodel_source()
