Add python api to gromacs build
[alexxy/gromacs.git] / src / python / sip / trajectoryanalysis / analysismodule.sip
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2011,2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
36 struct t_trxframe {
37
38 %TypeHeaderCode
39 #include <gromacs/fileio/trx.h>
40 #include "vec.h"
41 %End
42
43     int flags;
44     int not_ok;
45     bool bDouble;
46     int natoms;
47     double t0;
48     double tf;
49     double tpf;
50     double tppf;
51     bool bTitle;
52     const char *title;
53     bool bStep;
54     int step;
55     bool bTime;
56     double time;
57     bool bLambda;
58     bool bFepState;
59     double lambda;
60     int fep_state;
61     bool bAtoms;
62     // TODO
63     /*void *atoms;*/
64     bool bPrec;
65     double prec;
66     bool bX;
67     rvec *x {
68     %GetCode
69         rvecs *arr = new rvecs(sipCpp->x, sipCpp->natoms);
70         sipPy = sipConvertFromNewType(arr, sipType_rvecs, NULL);
71     %End
72     %SetCode
73
74     %End
75     };
76     bool bV;
77     rvec *v {
78     %GetCode
79         rvecs *arr = new rvecs(sipCpp->v, sipCpp->natoms);
80         sipPy = sipConvertFromNewType(arr, sipType_rvecs, NULL);
81     %End
82     %SetCode
83
84     %End
85     };
86     bool bF;
87     rvec *f {
88     %GetCode
89         rvecs *arr = new rvecs(sipCpp->f, sipCpp->natoms);
90         sipPy = sipConvertFromNewType(arr, sipType_rvecs, NULL);
91     %End
92     %SetCode
93
94     %End
95     };
96     bool bBox;
97     matrix box {
98     %GetCode
99         Matrix *mat = new Matrix(sipCpp->box);
100         sipPy = sipConvertFromNewType(mat, sipType_Matrix, NULL);
101     %End
102     %SetCode
103
104     %End
105     };
106     bool bPBC;
107     int ePBC;
108     // TODO
109     /*t_gmxvmdplugin* *vmdplugin; */
110 };
111
112 struct t_pbc {
113
114 %TypeHeaderCode
115 #include <gromacs/pbcutil/pbc.h>
116 #include "vec.h"
117 %End
118     int ePBC;
119     int ndim_ePBC;
120     int ePBCDX;
121     int dim;
122     matrix box {
123     %GetCode
124         Matrix *mat = new Matrix(sipCpp->box);
125         sipPy = sipConvertFromNewType(mat, sipType_Matrix, NULL);
126     %End
127     %SetCode
128
129     %End
130     };
131     rvec fbox_diag {
132     %GetCode
133         py_rvec *vec = new py_rvec(sipCpp->fbox_diag);
134         sipPy = sipConvertFromNewType(vec, sipType_py_rvec, NULL);
135     %End
136     %SetCode
137
138     %End
139     };
140     rvec hbox_diag {
141     %GetCode
142         py_rvec *vec = new py_rvec(sipCpp->fbox_diag);
143         sipPy = sipConvertFromNewType(vec, sipType_py_rvec, NULL);
144     %End
145     %SetCode
146
147     %End
148     };
149     rvec mhbox_diag {
150     %GetCode
151         py_rvec *vec = new py_rvec(sipCpp->mhbox_diag);
152         sipPy = sipConvertFromNewType(vec, sipType_py_rvec, NULL);
153     %End
154     %SetCode
155
156     %End
157     };
158     double max_cutoff2;
159     bool bLimitDistance;
160     double limit_distance2;
161     int  ntric_vec;
162     /*ivec tric_shift[MAX_NTRICVEC];*/
163     /*rvec tric_vec[MAX_NTRICVEC];*/
164 };
165
166 class TrajectoryAnalysisModuleData /NoDefaultCtors/ {
167
168 %TypeHeaderCode
169 #include <gromacs/trajectoryanalysis/analysismodule.h>
170 using namespace gmx;
171 %End
172
173 public:
174     virtual void finish() = 0;
175 };
176
177 class TrajectoryAnalysisModule {
178
179 %TypeHeaderCode
180 #include <gromacs/trajectoryanalysis/analysismodule.h>
181 using namespace gmx;
182 %End
183 public:
184     virtual void initOptions(Options *options, TrajectoryAnalysisSettings *settings) = 0;
185     virtual void optionsFinished(Options *options, TrajectoryAnalysisSettings *settings);
186     virtual void initAnalysis(const TrajectoryAnalysisSettings &settings, const TopologyInformation &top) = 0;
187     virtual void initAfterFirstFrame(const TrajectoryAnalysisSettings &settings, const t_trxframe &fr);
188 //    virtual TrajectoryAnalysisModuleDataPointer       startFrames (const AnalysisDataParallelOptions &opt, const SelectionCollection &selections);
189     virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc, TrajectoryAnalysisModuleData *pdata ) = 0;
190     virtual void finishFrames(TrajectoryAnalysisModuleData *pdata);
191     virtual void finishAnalysis(int nframes) = 0;
192     virtual void writeOutput() = 0;
193     const char* name() const;
194     const char* description() const;
195     int datasetCount() const;
196 protected:
197     TrajectoryAnalysisModule(const char *name, const char *description);
198 private:
199     TrajectoryAnalysisModule(const TrajectoryAnalysisModule &other);
200 };
201
202 %Exception gmx::InconsistentInputError(SIP_Exception) {
203
204 %TypeHeaderCode
205 #include <gromacs/utility/exceptions.h>
206 using namespace std;
207 %End
208
209 %RaiseCode
210     printf("lALLA\n");
211     const char *detail = sipExceptionRef.what();
212
213     SIP_BLOCK_THREADS
214     PyErr_SetString(sipException_gmx_InconsistentInputError, detail);
215     SIP_UNBLOCK_THREADS
216 %End
217 };
218
219 class TrajectoryAnalysisCommandLineRunner {
220 %TypeHeaderCode
221 #include <gromacs/trajectoryanalysis/cmdlinerunner.h>
222 %End
223
224 public:
225     TrajectoryAnalysisCommandLineRunner(TrajectoryAnalysisModule *module);
226     int run(SIP_PYLIST) throw (gmx::InconsistentInputError);
227     %MethodCode
228     int argc = PyList_GET_SIZE(a0);
229
230     char **argv = new char *[argc + 1];
231
232     // Convert the list.
233     for (int a = 0; a < argc; ++a)
234     {
235         PyObject *arg_obj = PyList_GET_ITEM(a0, a);
236         const char *arg = sipString_AsLatin1String(&arg_obj);
237
238         if (arg)
239         {
240             arg = strdup(arg);
241             Py_DECREF(arg_obj);
242         }
243         else
244         {
245             arg = "unknown";
246         }
247
248         argv[a] = const_cast<char *>(arg);
249     }
250
251     argv[argc] = NULL;
252
253     try {
254         sipCpp->run(argc, argv);
255     } catch (gmx::InconsistentInputError &e) {
256         sipIsErr = 1;
257         PyErr_SetString(sipException_gmx_InconsistentInputError, e.what());
258     }
259     %End
260 private:
261     TrajectoryAnalysisCommandLineRunner(const TrajectoryAnalysisCommandLineRunner &other);
262 };
263