OSDN Git Service

Batch mode is implemented (still testing)
[molby/Molby.git] / MolLib / Molecule.h
1 /*
2  *  Molecule.h
3  *
4  *  Created by Toshi Nagata on 06/03/11.
5  *  Copyright 2006-2008 Toshi Nagata. All rights reserved.
6  *
7  This program is free software; you can redistribute it and/or modify
8  it under the terms of the GNU General Public License as published by
9  the Free Software Foundation version 2 of the License.
10  
11  This program is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
14  GNU General Public License for more details.
15 */
16
17 #ifndef __Molecule_h__
18 #define __Molecule_h__
19
20 #include "Types.h"
21 #include "Object.h"
22 #include "IntGroup.h"
23
24 #include <stdio.h>
25
26 #ifdef __cplusplus
27 extern "C" {
28 #endif
29         
30 #define ATOMS_MAX_SYMMETRY 12
31 #define ATOMS_MAX_NUMBER 100000000  /*  Sufficiently large value  */
32
33 /*  Conversion between kcal/mol and internal energy unit (am*ang^2/fs^2, am = atomic mass) */
34 #define KCAL2INTERNAL (4.184e-4)
35 #define INTERNAL2KCAL (1.0/KCAL2INTERNAL)
36 #define J2INTERNAL (1e-4)
37 #define INTERNAL2J (1.0/J2INTERNAL)
38         
39 #define BOLTZMANN (8.31441e-3*J2INTERNAL)
40 #define PI 3.14159265358979
41 #define PI2R 0.564189583547756    /*  1.0/sqrt(PI)  */
42         
43 /*  Anisotropic thermal parameter  */
44 typedef struct Aniso {
45         Double  bij[6];    /*  b11, b22, b33, b12, b13, b23 (ORTEP type 0) */
46         char has_bsig;     /*  Has sigma values?  */
47         Double  bsig[6];   /*  sigma values  */
48         Mat33  pmat;      /*  A 3x3 matrix whose three column vectors are the principal axes of the ellipsoid. Note: If the B matrix is not positive definite, the axis length corresponding to the negative eigenvalue is replaced with 0.001.  */
49 } Aniso;
50
51 /*  Symmetry operation  */
52 /*  If periodic box is defined, dx/dy/dz denote multiples of the axes of the periodic box.
53     Otherwise, dx/dy/dz denote offset to the x/y/z coordinates of the atoms.  */
54 typedef struct Symop {
55         signed int dx : 4;
56         signed int dy : 4;
57         signed int dz : 4;
58         unsigned int sym : 8;
59         unsigned int alive: 1;
60 } Symop;
61
62 /*  Exflags  */
63 enum {
64         kAtomHiddenFlag = 1
65 };
66
67 /*  Atom connection record  */
68 /*  If nconnects <= ATOM_CONNECT_LIMIT, data[] field is used. Otherwise,
69     memory is allocated by malloc().  */
70 #define ATOM_CONNECT_LIMIT 6
71 typedef struct AtomConnect {
72         Int    count;  /*  Number of connections  */
73         union {
74                 Int *ptr;
75                 Int data[ATOM_CONNECT_LIMIT];
76         } u;
77 } AtomConnect;
78
79 typedef struct PiAnchor {
80         AtomConnect connect;
81         Int ncoeffs;
82         Double *coeffs;
83 } PiAnchor;
84
85 /*  Atom record  */
86 typedef struct Atom {
87         Int    segSeq;
88         char   segName[4];
89         Int    resSeq;
90         char   resName[4];
91         char   aname[4];
92         UInt   type;
93         Double  charge;
94         Double  weight;
95         char   element[4];
96         Int    atomicNumber;
97         AtomConnect connect;
98         Vector r;  /*  position  */
99         Vector v;  /*  velocity  */
100         Vector f;  /*  force  */
101         Vector sigma;   /*  For crystallographic data only; sigma for each crystallographic coordinates  */
102                                         /*  (Unlike r, these are not converted to the cartesian system)  */
103         Double  occupancy;
104         Double  tempFactor;
105         Aniso  *aniso;
106         Int    intCharge;
107         Int    exflags;
108         Int    nframes;  /*  Multiple frames  */
109         Vector *frames;
110         Symop  symop;    /*  For symmetry-expanded atom  */
111         Int    symbase;  /*  The index of original atom for symmetry-expansion  */
112         PiAnchor *anchor;  /*  Non-NULL if this atom is a pi-anchor  */
113         Int    labelid;  /*  The label ID; 0 for no label  */
114         short  wrap_dx, wrap_dy, wrap_dz; /*  Calculated by md_wrap_coordinates; used only in wrapped output.  */
115         Double fix_force; /*  0: no fix, >0: fix at fix_pos with harmonic potential, <0: fix at fix_pos without force  */
116         Vector fix_pos;
117         Byte   mm_exclude;        /*  If nonzero, then this atom is excluded from MM/MD calculations  */
118         Byte   periodic_exclude;  /*  If nonzero, then this atom is excluded from periodic calculations  */
119         char   uff_type[6]; /*  UFF type string  */
120 } Atom;
121
122 extern Int gSizeOfAtomRecord;
123
124 #define ATOM_AT_INDEX(p, i)  ((Atom *)((char *)(p) + (i) * gSizeOfAtomRecord))
125 #define ATOM_NEXT(p)         ((Atom *)((char *)(p) + gSizeOfAtomRecord))
126 #define ATOM_PREV(p)         ((Atom *)((char *)(p) - gSizeOfAtomRecord))
127 #define SYMOP_ALIVE(s) ((s.dx || s.dy || s.dz || s.sym) != 0)
128 #define SYMOP_EQUAL(s1, s2) (s1.dx == s2.dx && s1.dy == s2.dy && s1.dz == s2.dz && s1.sym == s2.sym)
129 #define SYMMETRY_AT_INDEX(p, i) (*((i) == 0 ? &gIdentityTransform : &p[i]))
130
131 /*  atom.connects is a union entry, including direct data for nconnects <= ATOM_CONNECT_LIMIT
132     and malloc()'ed entry for nconnects > ATOM_CONNECT_LIMIT. The following functions
133         automatically take care of the memory allocation/deallocation.  */
134 Int *AtomConnectData(AtomConnect *ac);
135 void AtomConnectResize(AtomConnect *ac, Int nconnects);
136 void AtomConnectInsertEntry(AtomConnect *ac, Int idx, Int connect);
137 void AtomConnectDeleteEntry(AtomConnect *ac, Int idx);
138
139 #define ATOM_CONNECT_PTR(ac) ((ac)->count > ATOM_CONNECT_LIMIT ? (ac)->u.ptr : (ac)->u.data)
140
141 /*  Duplicate an atom. If dst is non-NULL, *src is copied to *dst and dst is returned. If dst is NULL, a new atom is allocated by malloc() and that atom is returned. It is the called's responsibility to release the returned memory. */
142 extern Atom *AtomDuplicate(Atom *dst, const Atom *src);
143         
144 /*  Duplicate an atom, except for the frame entry  */
145 extern Atom *AtomDuplicateNoFrame(Atom *dst, const Atom *src);
146         
147 /*  Clean the content of an atom record  */
148 extern void AtomClean(Atom *ap);
149
150 /*  MolEnumerable type code  */
151 enum {
152         kAtomKind = 0,
153         kBondKind = 1,
154         kAngleKind = 2,
155         kDihedralKind = 3,
156         kImproperKind = 4,
157         kResidueKind = 5,
158         kEndKind
159 };
160
161 /*  Enumerable class to access to atoms, bonds, etc.  */
162 typedef struct MolEnumerable {
163         struct Molecule *mol;
164         int    kind;
165 } MolEnumerable;
166
167 /*  Atom reference  */
168 typedef struct AtomRef {
169         struct Molecule *mol;
170         int idx;
171 } AtomRef;
172
173 /*  Crystallographic cell parameter (also used as periodic box in MD) */
174 typedef struct XtalCell {
175         Double  cell[6];     /*  a, b, c, alpha, beta, gamma (in degree)  */
176         Double  rcell[6];    /*  Reciprocal cell  */
177         Vector  axes[3];     /*  Cartesian unit vectors along the three axis  */
178         Vector  origin;      /*  Cartesian origin of the periodic box  */
179         char    flags[3];    /*  1 for periodic, 0 for non-periodic  */
180         char    has_sigma;   /*  Has sigma?  */
181         Transform tr;        /*  Crystal coord -> cartesian  */
182         Transform rtr;       /*  Cartesian -> crystal coord  */
183         Double  cellsigma[6];  /*  For crystallographic data; sigma for the cell parameters  */
184 } XtalCell;
185
186 /*  3-Dimensional distribution  */
187 typedef struct Cube {
188         Int idn;             /*  Integer identifier (such as MO number)  */
189         Vector origin;
190         Vector dx, dy, dz;
191         Int nx, ny, nz;
192         Double *dp;          /*  Value for point (ix, iy, iz) is in dp[(ix*ny+iy)*nz+iz]  */
193 } Cube;
194
195 /*  Gaussian orbital symmetry types  */
196 enum {
197         kGTOType_S,
198         kGTOType_SP,
199         kGTOType_P,
200         kGTOType_D,
201         kGTOType_D5,
202         kGTOType_F,
203         kGTOType_F7,
204         kGTOType_G,
205         kGTOType_G9,
206         kGTOType_UU
207 };
208
209 /*  Exponent/coefficient info for a single gaussian primitive  */
210 typedef struct PrimInfo {
211         Double A;            /*  Exponent  */
212         Double C;            /*  Contraction coefficient  */
213         Double Csp;          /*  P(S=P) contraction coefficient  */
214 } PrimInfo;
215
216 /*  Gaussian orbital shell information  */
217 typedef struct ShellInfo {
218         signed char sym;     /*  Symmetry of the basis; S, P, ... */
219         signed char ncomp;   /*  Number of components (S: 1, P: 3, SP: 4, etc.)  */
220         signed char nprim;   /*  Number of primitives for this shell  */
221         Int p_idx;           /*  Index to the PrimInfo (exponent/coefficient) table  */
222         Int cn_idx;          /*  Index to the normalized (cached) contraction coefficient table  */
223         Int a_idx;           /*  Index to the atom which this primitive belongs to */
224         Int m_idx;           /*  Index to the MO matrix  */
225 } ShellInfo;
226
227 /*  Basis set and MO information  */
228 typedef struct BasisSet {
229         Int nshells;         /*  Number of gaussian orbital shells  */
230         ShellInfo *shells;   /*  Gaussian orbital shells  */
231         Int npriminfos;      /*  Size of primitive information table  */
232         PrimInfo *priminfos; /*  Primitive information table  */
233         Int ncns;            /*  Number of normalized (cached) contraction coefficient values  */
234         Double *cns;         /*  Normalized (cached) contraction coefficients; (up to 10 values for each primitive)  */
235         Int natoms_bs;       /*  Number of atoms; separately cached here because MO info should be invariant during editing */
236         Double *nuccharges;  /*  Nuclear charges (for ECP atoms)  */
237         Int ne_alpha, ne_beta;  /*  Number of alpha/beta electrons  */
238         Int rflag;           /*  0: UHF, 1: RHF, 2:ROHF  */
239         Int ncomps;          /*  Number of AO components; equal to sum of shells[i].ncomp  */
240         Int nmos;            /*  Number of MOs; equal to ncomps if close shell, ncomps*2 if open shell */
241         Double *mo;          /*  MO matrix (mo[i][j] represents the j-th AO coefficient for the i-th MO)  */
242                                                  /*  Memory are allocated for (2*nmos+1) entries; the last entry is for displaying arbitrary vector  */
243         Double *moenergies;  /*  MO energies  */
244         Double *scfdensities; /*  SCF densities; lower triangle of a symmetric matrix (size nmos*(nmos+1)/2)  */
245         Int ncubes;          /*  Number of calculated MOs  */
246         Cube **cubes;        /*  Calculated MOs (an array of pointers to Cubes)  */
247 } BasisSet;
248
249 /*  Marching Cube (for drawing isosurface)  */
250 typedef struct MCubePoint {
251         Int key;       /*  key = ((ix*ny+iy)*nz+iz)*3+ii, ii=0/1/2 for x/y/z direction, respectively */
252         float d;       /*  offset toward the direction; 0 <= d < 1  */
253         float pos[3];  /*  cartesian coordinate of the point  */
254         float grad[3]; /*  gradient vector  */
255 } MCubePoint;
256         
257 typedef struct MCube {
258         char hidden;         /*  If non-zero, then this MCube is not drawn  */
259         Int idn;             /*  MO number  */
260         Vector origin;       /*  Cube origin */
261         Double dx, dy, dz;   /*  Cube steps */
262         Int nx, ny, nz;      /*  Cube dimension (must be multiples of 8)  */
263         Double thres;        /*  Threshold value  */
264         Double *dp;          /*  Value for point (ix, iy, iz) is in dp[(ix*ny+iy)*nz+iz]  */
265         Int nradii;
266         Double *radii;       /*  Estimated radius (with margin) for each atom  */
267         Double expand;       /*  Expand the estimated radius by this value (default: 1.0)  */
268         struct {
269                 /*  Flags for cube (ix, iy, iz)-(ix+1, iy+1, iz+1). It is an 8-bit */
270                 /*  integer representing whether the values at the 8 corners are */
271                 /*  larger than the threshold value or not. As special cases,  */
272                 /*  the values 0 and 255 (all corners are below or above the threshold) */
273                 /*  are represented as 255, and the value 0 is used to indicate "yet undefined". */
274                 unsigned char *fp;
275                 /*  Cube points and triangles: for positive and negative surfaces  */
276                 Int ncubepoints;
277                 MCubePoint *cubepoints;
278                 Int ntriangles;
279                 Int *triangles;  /*  Triangles; indices to cubepoints[]  */
280                 float rgba[4];   /*  Surface color  */
281         } c[2];
282 } MCube;
283
284 /*  Electrostatic potential  */
285 typedef struct Elpot {
286         Vector pos;
287         Double esp;
288 } Elpot;
289
290 /*  Properties (total energy etc.; specific for each frame)  */
291 typedef struct MolProp {
292         char *propname;
293         Double *propvals;
294 } MolProp;
295         
296 /*  Molecule record  */
297 typedef struct Molecule {
298         Object base;
299         Int    natoms;
300         Atom   *atoms;
301         Int    nbonds;
302         Int    *bonds;       /*  The size of array is 2*nbonds  */
303         Int    nangles;
304         Int    *angles;      /*  The size of array is 3*nangles  */
305         Int    ndihedrals;
306         Int    *dihedrals;   /*  The size of array is 4*ndihedrals  */
307         Int    nimpropers;
308         Int    *impropers;   /*  The size of array is 4*nimpropers  */
309         Int    nresidues;    /*  Number of residues; maximum residue number + 1 (because residue 0 is 'undefined residue')  */
310         char   (*residues)[4];
311         XtalCell   *cell;
312         Int    nsyms;        /*  Symmetry operations; syms are always described in crystallographic units (even when the unit cell is not defined)  */
313         Transform *syms;
314
315         IntGroup *selection;
316         Int    nframes;      /*  The number of frames (>= 1). This is a cached value, and should be
317                                                          recalculated from the atoms if it is -1  */
318         Int    cframe;       /*  The current frame number  */
319
320         Int    nframe_cells;
321         Vector *frame_cells; /*  The cell vectors for frames; (nframe_cells*4) array of Vectors  */
322
323         struct MainView *mview;  /*  Reference to the MainView object if present (no retain)  */
324         Int    modifyCount;  /*  Internal counter for modification. This value is not to be modified
325                                  manually; instead, call MoleculeIncrementModifyCount() whenever
326                                                      modification is done, which also takes care necessary notification
327                                                          to the other part of the application (system dependent)  */
328
329         struct MDArena *arena;  /*  Reference to the MDArena record during MM/MD run (no retain)  */
330
331         const char *path;     /*  The full path of the molecule, when this molecule is last saved to/loaded from file. Only used in the command-line version. (In GUI version, the path is obtained by the Document mechanism) */
332
333         /*  Information from the dcd files  */
334         Int    startStep;     /*  the timestep for frame 0  */
335         Int    stepsPerFrame; /*  the number of timesteps between neighboring frames  */
336         Double psPerStep;     /*  picosecond per step  */
337
338         /*  Information for basis sets and MOs  */
339         BasisSet *bset;
340         
341         /*  Marching cube  */
342         MCube *mcube;
343
344         /*  Electrostatic potential  */
345         Int    nelpots;
346         Elpot  *elpots;
347
348         /*  Properties  */
349         Int    nmolprops;
350         MolProp *molprops;
351
352         /*  Parameters specific for this molecule  */
353         struct Parameter *par;
354         
355         /*  Bond order (not to be used in MM/MD, but may be necessary to hold this info)  */
356         Int    nbondOrders;
357         Double *bondOrders;
358
359         /*  Flag to request rebuilding MD internal information  */
360         Byte   needsMDRebuild;
361         
362         /*  Flag to clear selection of the parameter table  */
363         Byte   parameterTableSelectionNeedsClear;
364         
365         /*  Flag to request copying coordinates to MD arena  */
366         Byte   needsMDCopyCoordinates;
367
368         /*  Prohibit modification of the topology (to avoid interfering MD) */
369         Byte   noModifyTopology;
370         
371         /*  Flag to request aborting a subthread  */
372         Byte   requestAbortThread;
373
374         /*  Flag to signal that a subthread is terminated  */
375         Byte   threadTerminated;
376
377         /*  Mutex object. If non-NULL, it should be locked before modifying molecule  */
378         void *mutex;
379         
380         /*  Flag to prohibit modification from user interface  */
381         Byte   dontModifyFromGUI;
382
383         /*  Ruby pointer (see ruby_bind.c)  */
384         void *exmolobj;
385         Byte exmolobjProtected;
386
387 } Molecule;
388
389 int strlen_limit(const char *s, int limit);
390
391 void BasisSetRelease(BasisSet *bset);
392
393 Molecule *MoleculeNew(void);
394 int MoleculeLoadFile(Molecule *mp, const char *fname, const char *ftype, char **errbuf);
395 int MoleculeLoadPsfFile(Molecule *mp, const char *fname, char **errbuf);
396 int MoleculeLoadTepFile(Molecule *mp, const char *fname, char **errbuf);
397 int MoleculeLoadShelxFile(Molecule *mp, const char *fname, char **errbuf);
398 int MoleculeLoadGaussianFchkFile(Molecule *mp, const char *fname, char **errbuf);
399 int MoleculeLoadMbsfFile(Molecule *mp, const char *fname, char **errbuf);
400 Molecule *MoleculeNewWithName(const char *name);
401 Molecule *MoleculeInitWithAtoms(Molecule *mp, const Atom *atoms, int natoms);
402 Molecule *MoleculeInitWithMolecule(Molecule *mp2, Molecule *mp);
403 void MoleculeSetName(Molecule *par, const char *name);
404 const char *MoleculeGetName(Molecule *mp);
405 void MoleculeSetPath(Molecule *mol, const char *fname);
406 const char *MoleculeGetPath(Molecule *mol);
407 Molecule *MoleculeWithName(const char *name);
408 Molecule *MoleculeRetain(Molecule *mp);
409 void MoleculeRelease(Molecule *mp);
410 void MoleculeExchange(Molecule *mp1, Molecule *mp2);
411
412 int MoleculeAddGaussianOrbitalShell(Molecule *mol, Int a_idx, Int sym, Int nprims);
413 int MoleculeAddGaussianPrimitiveCoefficients(Molecule *mol, Double exponent, Double contraction, Double contraction_sp);
414 int MoleculeGetGaussianComponentInfo(Molecule *mol, Int comp_idx, Int *outAtomIdx, char *outLabel, Int *outShellIdx);
415 int MoleculeSetMOCoefficients(Molecule *mol, Int idx, Double energy, Int ncomps, Double *coeffs);
416 int MoleculeGetMOCoefficients(Molecule *mol, Int idx, Double *energy, Int *ncoeffs, Double **coeffs);
417 int MoleculeSetMOInfo(Molecule *mol, Int rflag, Int ne_alpha, Int ne_beta);
418
419 void MoleculeIncrementModifyCount(Molecule *mp);
420 void MoleculeClearModifyCount(Molecule *mp);
421
422 MolEnumerable *MolEnumerableNew(Molecule *mol, int kind);
423 void MolEnumerableRelease(MolEnumerable *mseq);
424 AtomRef *AtomRefNew(Molecule *mol, int idx);
425 void AtomRefRelease(AtomRef *aref);
426
427 void MoleculeSetCell(Molecule *mp, Double a, Double b, Double c, Double alpha, Double beta, Double gamma, int convertCoordinates);
428 void MoleculeSetAniso(Molecule *mp, int n1, int type, Double x11, Double x22, Double x33, Double x12, Double x13, Double x23, const Double *sigmap);
429 void MoleculeSetAnisoBySymop(Molecule *mp, int idx);
430 int MoleculeSetPeriodicBox(Molecule *mp, const Vector *ax, const Vector *ay, const Vector *az, const Vector *ao, const char *periodic, int convertCoordinates);
431 int MoleculeCalculateCellFromAxes(XtalCell *cp, int calc_abc);
432
433 int MoleculeReadCoordinatesFromFile(Molecule *mp, const char *fname, const char *ftype, char **errbuf);
434 int MoleculeReadCoordinatesFromPdbFile(Molecule *mp, const char *fname, char **errbuf);
435 int MoleculeReadCoordinatesFromDcdFile(Molecule *mp, const char *fname, char **errbuf);
436
437 int MoleculeLoadGamessDatFile(Molecule *mol, const char *fname, char **errbuf);
438
439 int MoleculeReadExtendedInfo(Molecule *mp, const char *fname, char **errbuf);
440 int MoleculeWriteExtendedInfo(Molecule *mp, const char *fname, char **errbuf);
441
442 int MoleculeWriteToFile(Molecule *mp, const char *fname, const char *ftype, char **errbuf);
443 int MoleculeWriteToPsfFile(Molecule *mp, const char *fname, char **errbuf);
444 int MoleculeWriteToPdbFile(Molecule *mp, const char *fname, char **errbuf);
445 int MoleculeWriteToDcdFile(Molecule *mp, const char *fname, char **errbuf);
446 int MoleculeWriteToTepFile(Molecule *mp, const char *fname, char **errbuf);
447 int MoleculeWriteToMbsfFile(Molecule *mp, const char *fname, char **errbuf);
448 void MoleculeDump(Molecule *mol);
449
450 int MoleculePrepareMDArena(Molecule *mol, int check_only, char **retmsg);
451
452 char *MoleculeSerialize(Molecule *mp, Int *outLength, Int *timep);
453 Molecule *MoleculeDeserialize(const char *data, Int length, Int *timep);
454
455 void MoleculeCleanUpResidueTable(Molecule *mp);
456 int MoleculeChangeNumberOfResidues(Molecule *mp, int nresidues);
457 int MoleculeChangeResidueNumberWithArray(Molecule *mp, IntGroup *group, Int *resSeqs);
458 int MoleculeChangeResidueNumber(Molecule *mp, IntGroup *group, int resSeq);
459 int MoleculeOffsetResidueNumbers(Molecule *mp, IntGroup *group, int offset, int nresidues);
460 int MoleculeChangeResidueNames(Molecule *mp, int argc, Int *resSeqs, char *names);
461 int MoleculeMaximumResidueNumber(Molecule *mp, IntGroup *group);
462 int MoleculeMinimumResidueNumber(Molecule *mp, IntGroup *group);
463
464 struct MolAction;
465 #if defined(DEBUG)
466         int MoleculeCheckSanity(Molecule *mp);
467 #else
468 #define MoleculeCheckSanity(mp)
469 #endif
470
471 int MoleculeCreateAnAtom(Molecule *mp, const Atom *ap, int pos);
472 int MoleculeMerge(Molecule *dst, Molecule *src, IntGroup *where, int resSeqOffset, Int *nactions, struct MolAction ***actions, Int forUndo);
473 int MoleculeUnmerge(Molecule *src, Molecule **dstp, IntGroup *where, int resSeqOffset, Int *nactions, struct MolAction ***actions, Int forUndo);
474 int MoleculeExtract(Molecule *src, Molecule **dstp, IntGroup *where, int dummyFlag);
475 int MoleculeAddBonds(Molecule *mp, Int nbonds, const Int *bonds, IntGroup *where, Int autoGenerate);
476 int MoleculeDeleteBonds(Molecule *mp, Int *bonds, IntGroup *where, Int **outRemoved, IntGroup **outRemovedPos);
477 int MoleculeAssignBondOrders(Molecule *mp, const Double *orders, IntGroup *where);
478 int MoleculeGetBondOrders(Molecule *mp, Double *outOrders, IntGroup *where);
479 int MoleculeAddAngles(Molecule *mp, const Int *angles, IntGroup *where);
480 int MoleculeDeleteAngles(Molecule *mp, Int *angles, IntGroup *where);
481 int MoleculeAddDihedrals(Molecule *mp, const Int *dihedrals, IntGroup *where);
482 int MoleculeDeleteDihedrals(Molecule *mp, Int *dihedrals, IntGroup *where);
483 int MoleculeAddImpropers(Molecule *mp, const Int *impropers, IntGroup *where);
484 int MoleculeDeleteImpropers(Molecule *mp, Int *impropers, IntGroup *where);
485 int MoleculeLookupBond(Molecule *mp, Int n1, Int n2);
486 int MoleculeLookupAngle(Molecule *mp, Int n1, Int n2, Int n3);
487 int MoleculeLookupDihedral(Molecule *mp, Int n1, Int n2, Int n3, Int n4);
488 int MoleculeLookupImproper(Molecule *mp, Int n1, Int n2, Int n3, Int n4);
489
490 Int MoleculeFindMissingAngles(Molecule *mol, Int **outAngles);
491 Int MoleculeFindMissingDihedrals(Molecule *mol, Int **outDihedrals);
492 Int MoleculeFindMissingImpropers(Molecule *mol, Int **outImpropers);
493         
494 IntGroup *MoleculeSearchBondsIncludingAtoms(Molecule *mp, IntGroup *atomgroup);
495 IntGroup *MoleculeSearchAnglesIncludingAtoms(Molecule *mp, IntGroup *atomgroup);
496 IntGroup *MoleculeSearchDihedralsIncludingAtoms(Molecule *mp, IntGroup *atomgroup);
497 IntGroup *MoleculeSearchImpropersIncludingAtoms(Molecule *mp, IntGroup *atomgroup);
498
499 IntGroup *MoleculeSearchBondsAcrossAtomGroup(Molecule *mp, IntGroup *atomgroup);
500
501 IntGroup *MoleculeSearchAnglesIncludingBond(Molecule *mp, int n1, int n2);
502 IntGroup *MoleculeSearchDihedralsIncludingBond(Molecule *mp, int n1, int n2);
503 IntGroup *MoleculeSearchImpropersIncludingBond(Molecule *mp, int n1, int n2);
504
505 int MoleculeLookupAtomInResidue(Molecule *mp, int n1, int resno);
506 int MoleculeAnalyzeAtomName(const char *s, char *resName, int *resSeq, char *atomName);
507 int MoleculeAtomIndexFromString(Molecule *mp, const char *s);
508
509 int MoleculeFindCloseAtoms(Molecule *mp, const Vector *vp, Double radius, Double limit, Int *outNbonds, Int **outBonds, Int triangle);
510 int MoleculeGuessBonds(Molecule *mp, Double limit, Int *outNbonds, Int **outBonds);
511 int MoleculeRebuildTablesFromConnects(Molecule *mp);
512 int MoleculeAreAtomsConnected(Molecule *mol, int idx1, int idx2);
513         
514 void MoleculeGetAtomName(Molecule *mp, int index, char *buf, int bufsize);
515
516 void MoleculeSetSelection(Molecule *mp, IntGroup *select);
517 IntGroup *MoleculeGetSelection(Molecule *mp);
518 void MoleculeSelectAtom(Molecule *mp, int n1, int extending);
519 void MoleculeUnselectAtom(Molecule *mp, int n1);
520 void MoleculeToggleSelectionOfAtom(Molecule *mp, int n1);
521 int MoleculeIsAtomSelected(Molecule *mp, int n1);
522 int MoleculeIsBondSelected(Molecule *mp, int n1, int n2);
523 IntGroup *MoleculeModifySelectionByRemovingAtoms(Molecule *mp, IntGroup *selection, IntGroup *remove);
524
525 int MoleculeGetTransformForSymop(Molecule *mp, Symop symop, Transform *tf, int is_cartesian);
526 int MoleculeGetSymopForTransform(Molecule *mp, const Transform tf, Symop *symop, int is_cartesian);
527
528 int MoleculeTransformBySymop(Molecule *mp, const Vector *vpin, Vector *vpout, Symop symop);
529 int MoleculeAddExpandedAtoms(Molecule *mp, Symop symop, IntGroup *group, Int *indices, Int allowOverlap);
530 int MoleculeAmendBySymmetry(Molecule *mp, IntGroup *group, IntGroup **groupout, Vector **vpout);
531
532 int MoleculeShowAllAtoms(Molecule *mp);
533 int MoleculeShowReverse(Molecule *mp);
534 int MoleculeHideAtoms(Molecule *mp, IntGroup *ig);
535
536 int MoleculeRenumberAtoms(Molecule *mp, const Int *new2old, Int *old2new_out, Int isize);
537
538 void MoleculeTransform(Molecule *mp, Transform tr, IntGroup *group);
539 void MoleculeTranslate(Molecule *mp, const Vector *vp, IntGroup *group);
540 void MoleculeRotate(Molecule *mp, const Vector *axis, Double angle, const Vector *center, IntGroup *group);
541 int MoleculeCenterOfMass(Molecule *mp, Vector *center, IntGroup *group);
542 int MoleculeBounds(Molecule *mp, Vector *min, Vector *max, IntGroup *group);
543         
544 Int *MoleculeSearchEquivalentAtoms(Molecule *mol, IntGroup *ig);
545         
546 void MoleculeAddExpansion(Molecule *mp, Vector dr, Int symop, IntGroup *group, Double limit);
547 void MoleculeClearExpansion(Molecule *mp, IntGroup *group);
548 void MoleculeRemoveExpansion(Molecule *mp, Vector dr, Int symop, IntGroup *group);
549 void MoleculeAutoExpansion(Molecule *mp, const float *boxstart, const float *boxend, IntGroup *group, Double limit);
550
551 void MoleculeXtalToCartesian(Molecule *mp, Vector *dst, const Vector *src);
552 void MoleculeCartesianToXtal(Molecule *mp, Vector *dst, const Vector *src);
553 Double MoleculeMeasureBond(Molecule *mp, const Vector *vp1, const Vector *vp2);
554 Double MoleculeMeasureAngle(Molecule *mp, const Vector *vp1, const Vector *vp2, const Vector *vp3);
555 Double MoleculeMeasureDihedral(Molecule *mp, const Vector *vp1, const Vector *vp2, const Vector *vp3, const Vector *vp4);
556
557 IntGroup *MoleculeFragmentExcludingAtomGroup(Molecule *mp, int n1, IntGroup *exatoms);
558 IntGroup *MoleculeFragmentExcludingAtoms(Molecule *mp, int n1, int argc, int *argv);
559 IntGroup *MoleculeFragmentWithAtomGroups(Molecule *mp, IntGroup *inatoms, IntGroup *exatoms);
560 int MoleculeIsFragmentDetachable(Molecule *mp, IntGroup *group, int *n1, int *n2);
561 int MoleculeIsFragmentRotatable(Molecule *mp, IntGroup *group, int *n1, int *n2, IntGroup **rotGroup);
562
563 int MoleculeGetNumberOfFrames(Molecule *mp);
564 int MoleculeInsertFrames(Molecule *mp, IntGroup *group, const Vector *inFrame, const Vector *inFrameCell);
565 int MoleculeRemoveFrames(Molecule *mp, IntGroup *group, Vector *outFrame, Vector *outFrameCell);
566 int MoleculeSelectFrame(Molecule *mp, int frame, int copyback);
567 int MoleculeFlushFrames(Molecule *mp);
568 int MoleculeReorderFrames(Molecule *mp, const Int *old_idx);
569
570 int MoleculeCreateProperty(Molecule *mp, const char *name);
571 int MoleculeLookUpProperty(Molecule *mp, const char *name);
572 int MoleculeDeletePropertyAtIndex(Molecule *mp, int idx);
573 int MoleculeSetProperty(Molecule *mp, int idx, IntGroup *ig, const Double *values);
574 int MoleculeGetProperty(Molecule *mp, int idx, IntGroup *ig, Double *outValues);
575
576 void MoleculeUpdatePiAnchorPositions(Molecule *mol);
577 void MoleculeCalculatePiAnchorPosition(Molecule *mol, int idx);
578 int MoleculeSetPiAnchorList(Molecule *mol, Int idx, Int nentries, Int *entries, Double *weights, Int *nUndoActions, struct MolAction ***undoActions);
579         
580 int MoleculeCalcMO(Molecule *mp, Int mono, const Vector *op, const Vector *dxp, const Vector *dyp, const Vector *dzp, Int nx, Int ny, Int nz, int (*callback)(double progress, void *ref), void *ref);
581 int MoleculeGetDefaultMOGrid(Molecule *mp, Int npoints, Vector *op, Vector *xp, Vector *yp, Vector *zp, Int *nx, Int *ny, Int *nz);
582 const Cube *MoleculeGetCubeAtIndex(Molecule *mp, Int index);
583 int MoleculeLookUpCubeWithMONumber(Molecule *mp, Int mono);
584 int MoleculeClearCubeAtIndex(Molecule *mp, Int index);
585 int MoleculeOutputCube(Molecule *mp, Int index, const char *fname, const char *comment);
586
587 MCube *MoleculeClearMCube(Molecule *mol, Int nx, Int ny, Int nz, const Vector *origin, Double dx, Double dy, Double dz);
588 int MoleculeUpdateMCube(Molecule *mol, int idn);
589 void MoleculeDeallocateMCube(MCube *mcube);
590
591 extern char *gMoleculePasteboardType;
592 extern char *gParameterPasteboardType;
593 extern char *gLoadSaveErrorMessage;
594         
595 STUB void MoleculeRetainExternalObj(Molecule *mol);
596 STUB void MoleculeReleaseExternalObj(Molecule *mol);
597
598 STUB int MoleculeCallback_writeToPasteboard(const char *type, const void *data, int length);
599 STUB int MoleculeCallback_readFromPasteboard(const char *type, void **dptr, int *length);
600 STUB int MoleculeCallback_isDataInPasteboard(const char *type);
601
602 STUB Molecule *MoleculeCallback_openNewMolecule(const char *fname);
603 STUB void MoleculeCallback_notifyModification(Molecule *mp, int now_flag);
604 STUB Molecule *MoleculeCallback_currentMolecule(void);
605 STUB Molecule *MoleculeCallback_moleculeAtIndex(int idx);
606 STUB Molecule *MoleculeCallback_moleculeAtOrderedIndex(int idx);
607 STUB void MoleculeCallback_displayName(Molecule *mol, char *buf, int bufsize);
608 STUB void MoleculeCallback_pathName(Molecule *mol, char *buf, int bufsize);
609 STUB int MoleculeCallback_setDisplayName(Molecule *mol, const char *name);
610
611 STUB void MoleculeCallback_lockMutex(void *mutex);
612 STUB void MoleculeCallback_unlockMutex(void *mutex);
613 STUB void MoleculeCallback_disableModificationFromGUI(Molecule *mol);
614 STUB void MoleculeCallback_enableModificationFromGUI(Molecule *mol);
615         
616 STUB void MoleculeCallback_cannotModifyMoleculeDuringMDError(Molecule *mol);
617
618 STUB int MoleculeCallback_callSubProcessAsync(Molecule *mol, const char *cmd, int (*callback)(Molecule *, int), int (*timerCallback)(Molecule *, int), FILE *output, FILE *errout);
619
620 /*  This is also defined in Molby_extern.h, but it may be called from functions in Molecule.c  */
621 STUB int MyAppCallback_checkInterrupt(void);
622         
623 void MoleculeLock(Molecule *mol);
624 void MoleculeUnlock(Molecule *mol);
625
626 #if 0
627 #define __MoleculeLock(mol) MoleculeLock(mol)
628 #define __MoleculeUnlock(mol) MoleculeUnlock(mol)
629 #else
630 #define __MoleculeLock(mol)
631 #define __MoleculeUnlock(mol)
632 #endif
633         
634 #ifdef __cplusplus
635 }
636 #endif
637                 
638 #endif /* __Molecule_h__ */