Merge branch 'master' into pygromacs
[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) 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 %PostInitialisationCode
37     import_array();
38 %End
39
40 struct t_trxframe {
41
42 %TypeHeaderCode
43 #include <gromacs/fileio/trx.h>
44 #include "numpy_conv.h"
45 %End
46
47     int flags;
48     int not_ok;
49     bool bDouble;
50     int natoms;
51     double t0;
52     double tf;
53     double tpf;
54     double tppf;
55     bool bTitle;
56     const char *title;
57     bool bStep;
58     int step;
59     bool bTime;
60     double time;
61     bool bLambda;
62     bool bFepState;
63     double lambda;
64     int fep_state;
65     bool bAtoms;
66     t_atoms *atoms;
67     bool bPrec;
68     double prec;
69     bool bX;
70     SIP_PYOBJECT x {
71     %GetCode
72         sipPy = array2dToNumpy(sipCpp->natoms, 3, sipCpp->x);
73     %End
74     %SetCode
75
76     %End
77     };
78     bool bV;
79     SIP_PYOBJECT v {
80     %GetCode
81         sipPy = array2dToNumpy(sipCpp->natoms, 3, sipCpp->v);
82     %End
83     %SetCode
84
85     %End
86     };
87     bool bF;
88     SIP_PYOBJECT f {
89     %GetCode
90         sipPy = array2dToNumpy(sipCpp->natoms, 3, sipCpp->f);
91     %End
92     %SetCode
93
94     %End
95     };
96     bool bBox;
97     SIP_PYOBJECT box {
98     %GetCode
99         sipPy = array2dToNumpy(3, 3, sipCpp->box);
100     %End
101     %SetCode
102
103     %End
104     };
105     bool bPBC;
106     int ePBC;
107     // TODO
108     /*t_gmxvmdplugin* *vmdplugin; */
109 };
110
111 struct t_pbc {
112
113 %TypeHeaderCode
114 #include <gromacs/pbcutil/pbc.h>
115 #include "numpy_conv.h"
116 %End
117     int ePBC;
118     int ndim_ePBC;
119     int ePBCDX;
120     int dim;
121     SIP_PYOBJECT box {
122     %GetCode
123         sipPy = array2dToNumpy(3, 3, sipCpp->box);
124     %End
125     %SetCode
126
127     %End
128     };
129     SIP_PYOBJECT fbox_diag {
130     %GetCode
131         sipPy = array1dToNumpy(3, sipCpp->fbox_diag);
132     %End
133     %SetCode
134
135     %End
136     };
137     SIP_PYOBJECT hbox_diag {
138     %GetCode
139         sipPy = array1dToNumpy(3, sipCpp->hbox_diag);
140     %End
141     %SetCode
142
143     %End
144     };
145     SIP_PYOBJECT mhbox_diag {
146     %GetCode
147         sipPy = array1dToNumpy(3, sipCpp->mhbox_diag);
148     %End
149     %SetCode
150
151     %End
152     };
153     double max_cutoff2;
154     bool bLimitDistance;
155     double limit_distance2;
156     int  ntric_vec;
157     /*ivec tric_shift[MAX_NTRICVEC];*/
158     /*rvec tric_vec[MAX_NTRICVEC];*/
159 };
160
161 class TrajectoryAnalysisModuleData /NoDefaultCtors/ {
162
163 %TypeHeaderCode
164 #include <gromacs/trajectoryanalysis/analysismodule.h>
165 using namespace gmx;
166 %End
167
168 public:
169     virtual void finish() = 0;
170 };
171
172 %ModuleHeaderCode
173 #include "gromacs/utility/exceptions.h"
174 %End
175
176 %VirtualErrorHandler vehandler
177     SIP_RELEASE_GIL(sipGILState);
178
179     GMX_THROW(gmx::InternalError("Python virtual overload raised an exception, see traceback"));
180 %End
181
182 class TrajectoryAnalysisModule /VirtualErrorHandler=vehandler/ {
183
184 %TypeHeaderCode
185 #include <gromacs/trajectoryanalysis/analysismodule.h>
186 using namespace gmx;
187 %End
188 public:
189     virtual void initOptions(Options *options, TrajectoryAnalysisSettings *settings) = 0;
190     virtual void optionsFinished(Options *options, TrajectoryAnalysisSettings *settings);
191     virtual void initAnalysis(const TrajectoryAnalysisSettings &settings, const TopologyInformation &top) = 0;
192     virtual void initAfterFirstFrame(const TrajectoryAnalysisSettings &settings, const t_trxframe &fr);
193 //    virtual TrajectoryAnalysisModuleDataPointer       startFrames (const AnalysisDataParallelOptions &opt, const SelectionCollection &selections);
194     virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc, TrajectoryAnalysisModuleData *pdata ) = 0;
195     virtual void finishFrames(TrajectoryAnalysisModuleData *pdata);
196     virtual void finishAnalysis(int nframes) = 0;
197     virtual void writeOutput() = 0;
198     const char* name() const;
199     const char* description() const;
200     int datasetCount() const;
201 protected:
202     TrajectoryAnalysisModule(const char *name, const char *description);
203 private:
204     TrajectoryAnalysisModule(const TrajectoryAnalysisModule &other);
205 };
206
207 %Exception gmx::InconsistentInputError(SIP_Exception) {
208
209 %TypeHeaderCode
210 #include <gromacs/utility/exceptions.h>
211 using namespace std;
212 %End
213
214 %RaiseCode
215     const char *detail = sipExceptionRef.what();
216
217     SIP_BLOCK_THREADS
218     PyErr_SetString(sipException_gmx_InconsistentInputError, detail);
219     SIP_UNBLOCK_THREADS
220 %End
221 };