OSDN Git Service

Ruby: Molecule#make_front, set_name, get_view_rotation, get_view_scale, get_view_tran...
[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 #ifdef __cplusplus
25 extern "C" {
26 #endif
27         
28 #define ATOMS_MAX_CONNECTS 12
29 #define ATOMS_MAX_SYMMETRY 12
30 #define ATOMS_MAX_NUMBER 100000000  /*  Sufficiently large value  */
31
32 /*  Conversion between kcal/mol and internal energy unit (am*ang^2/fs^2, am = atomic mass) */
33 #define KCAL2INTERNAL (4.184e-4)
34 #define INTERNAL2KCAL (1.0/KCAL2INTERNAL)
35 #define J2INTERNAL (1e-4)
36 #define INTERNAL2J (1.0/J2INTERNAL)
37         
38 #define BOLTZMANN (8.31441e-3*J2INTERNAL)
39 #define PI 3.14159265358979
40                 
41 /*  Anisotropic thermal parameter  */
42 typedef struct Aniso {
43         Double  bij[6];    /*  b11, b22, b33, b12, b23, b31 (ORTEP type 0) */
44         char has_bsig;     /*  Has sigma values?  */
45         Double  bsig[6];   /*  sigma values  */
46         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.  */
47 } Aniso;
48
49 /*  Symmetry operation  */
50 /*  If periodic box is defined, dx/dy/dz denote multiples of the axes of the periodic box.
51     Otherwise, dx/dy/dz denote offset to the x/y/z coordinates of the atoms.  */
52 typedef struct Symop {
53         signed int dx : 4;
54         signed int dy : 4;
55         signed int dz : 4;
56         unsigned int sym : 8;
57         unsigned int alive: 1;
58 } Symop;
59
60 /*  Exflags  */
61 enum {
62         kAtomHiddenFlag = 1
63 };
64
65 /*  Atom record  */
66 typedef struct Atom {
67         Int    segSeq;
68         char   segName[4];
69         Int    resSeq;
70         char   resName[4];
71         char   aname[4];
72         UInt   type;
73         Double  charge;
74         Double  weight;
75         char   element[4];
76         Int    atomicNumber;
77         Int    nconnects;   /*  Number of connections (= bonds)  */
78         Int    connects[ATOMS_MAX_CONNECTS];
79         Vector r;  /*  position  */
80         Vector v;  /*  velocity  */
81         Vector f;  /*  force  */
82         Vector sigma;   /*  For crystallographic data only; sigma for each crystallographic coordinates  */
83                                         /*  (Unlike r, these are not converted to the cartesian system)  */
84         Double  occupancy;
85         Double  tempFactor;
86         Aniso  *aniso;
87         Int    intCharge;
88         Int    exflags;
89         Int    nframes;  /*  Multiple frames  */
90         Vector *frames;
91         Symop  symop;    /*  For symmetry-expanded atom  */
92         Int    symbase;  /*  The index of original atom for symmetry-expansion  */
93         Int    labelid;  /*  The label ID; 0 for no label  */
94         short  wrap_dx, wrap_dy, wrap_dz; /*  Calculated by md_wrap_coordinates; used only in wrapped output.  */
95         Double fix_force; /*  0: no fix, >0: fix at fix_pos with harmonic potential, <0: fix at fix_pos without force  */
96         Vector fix_pos;
97         Byte   mm_exclude;        /*  If nonzero, then this atom is excluded from MM/MD calculations  */
98         Byte   periodic_exclude;  /*  If nonzero, then this atom is excluded from periodic calculations  */
99 } Atom;
100
101 extern Int gSizeOfAtomRecord;
102
103 #define ATOM_AT_INDEX(p, i)  ((Atom *)((char *)(p) + (i) * gSizeOfAtomRecord))
104 #define ATOM_NEXT(p)         ((Atom *)((char *)(p) + gSizeOfAtomRecord))
105 #define ATOM_PREV(p)         ((Atom *)((char *)(p) - gSizeOfAtomRecord))
106 #define SYMOP_ALIVE(s) (s.dx || s.dy || s.dz || s.sym)
107 #define SYMOP_EQUAL(s1, s2) (s1.dx == s2.dx && s1.dy == s2.dy && s1.dz == s2.dz && s1.sym == s2.sym)
108
109 /*  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. */
110 extern Atom *AtomDuplicate(Atom *dst, const Atom *src);
111 /*  Clean the content of an atom record  */
112 extern void AtomClean(Atom *ap);
113
114 /*  MolEnumerable type code  */
115 enum {
116         kAtomKind = 0,
117         kBondKind = 1,
118         kAngleKind = 2,
119         kDihedralKind = 3,
120         kImproperKind = 4,
121         kResidueKind = 5,
122         kEndKind
123 };
124
125 /*  Enumerable class to access to atoms, bonds, etc.  */
126 typedef struct MolEnumerable {
127         struct Molecule *mol;
128         int    kind;
129 } MolEnumerable;
130
131 /*  Atom reference  */
132 typedef struct AtomRef {
133         struct Molecule *mol;
134         int idx;
135 } AtomRef;
136
137 /*  Crystallographic cell parameter (also used as periodic box in MD) */
138 typedef struct XtalCell {
139         Double  cell[6];     /*  a, b, c, alpha, beta, gamma (in degree)  */
140         Double  rcell[6];    /*  Reciprocal cell  */
141         Vector  axes[3];     /*  Cartesian unit vectors along the three axis  */
142         Vector  origin;      /*  Cartesian origin of the periodic box  */
143         char    flags[3];    /*  1 for periodic, 0 for non-periodic  */
144         char    has_sigma;   /*  Has sigma?  */
145         Transform tr;        /*  Crystal coord -> cartesian  */
146         Transform rtr;       /*  Cartesian -> crystal coord  */
147         Double  cellsigma[6];  /*  For crystallographic data; sigma for the cell parameters  */
148 } XtalCell;
149
150 /*  Periodic box parameter  */
151 #if 0
152 typedef struct PeriodicBox {
153         Vector  axes[3];     /*  Unit vectors along the three axis  */
154         Vector  origin;      /*  Origin of the periodic box  */
155         char    flags[3];    /*  1 for periodic, 0 for non-periodic  */
156         Transform tr;        /*  Internal coord -> cartesian  */
157         Transform rtr;       /*  Cartesian -> internal coord  */        
158 } PeriodicBox;
159 #endif
160
161 /*  Expanded atoms  */
162 typedef struct ExAtom {
163         Int    index;        /*  Base atom index  */
164         Vector dr;           /*  Translational offset  */
165         Int    symop;        /*  Symmetry operation  */
166         Int    labelid;      /*  Label ID; 0 for no label  */
167 } ExAtom;
168
169 /*  3-Dimensional distribution  */
170 typedef struct Cube {
171         Int idn;             /*  Integer identifier (such as MO number)  */
172         Vector origin;
173         Vector dx, dy, dz;
174         Int nx, ny, nz;
175         Double *dp;          /*  Value for point (ix, iy, iz) is in dp[(ix*ny+iy)*nz+iz]  */
176 } Cube;
177
178 /*  Gaussian orbital symmetry types  */
179 enum {
180         kGTOType_S,
181         kGTOType_SP,
182         kGTOType_P,
183         kGTOType_D,
184         kGTOType_D5,
185         kGTOtype_F,
186         kGTOType_F7,
187         kGTOType_UU
188 };
189
190 /*  Exponent/coefficient info for a single gaussian primitive  */
191 typedef struct PrimInfo {
192         Double A;            /*  Exponent  */
193         Double C;            /*  Contraction coefficient  */
194         Double Csp;          /*  P(S=P) contraction coefficient  */
195 } PrimInfo;
196
197 /*  Gaussian orbital shell information  */
198 typedef struct ShellInfo {
199         signed char sym;     /*  Symmetry of the basis; S, P, ... */
200         signed char ncomp;   /*  Number of components (S: 1, P: 3, SP: 4, etc.)  */
201         signed char nprim;   /*  Number of primitives for this shell  */
202         Int p_idx;           /*  Index to the PrimInfo (exponent/coefficient) table  */
203         Int cn_idx;          /*  Index to the normalized (cached) contraction coefficient table  */
204         Int a_idx;           /*  Index to the atom which this primitive belongs to */
205         Int m_idx;           /*  Index to the MO matrix  */
206 } ShellInfo;
207
208 /*  Basis set and MO information  */
209 typedef struct BasisSet {
210         Int nshells;         /*  Number of gaussian orbital shells  */
211         ShellInfo *shells;   /*  Gaussian orbital shells  */
212         Int npriminfos;      /*  Size of primitive information table  */
213         PrimInfo *priminfos; /*  Primitive information table  */
214         Int ncns;            /*  Number of normalized (cached) contraction coefficient values  */
215         Double *cns;         /*  Normalized (cached) contraction coefficients; (up to 10 values for each primitive)  */
216         Int natoms;          /*  Number of atoms; separately cached here because MO info should be invariant during editing */
217         Vector *pos;         /*  Positions of atoms; the unit is bohr, not angstrom  */
218         Double *nuccharges;  /*  Nuclear charges (for ECP atoms)  */
219         Int ne_alpha, ne_beta;  /*  Number of alpha/beta electrons  */
220         Int rflag;           /*  0: UHF, 1: RHF, 2:ROHF  */
221         Int ncomps;          /*  Number of AO components; equal to sum of shells[i].ncomp  */
222         Int nmos;            /*  Number of MOs; equal to ncomps if close shell, ncomps*2 if open shell */
223         Double *mo;          /*  MO matrix (mo[i][j] represents the j-th AO coefficient for the i-th MO)  */
224         Double *moenergies;  /*  MO energies  */
225         Double *scfdensities; /*  SCF densities; lower triangle of a symmetric matrix (size nmos*(nmos+1)/2)  */
226         Int ncubes;          /*  Number of calculated MOs  */
227         Cube **cubes;        /*  Calculated MOs (an array of pointers to Cubes)  */
228 } BasisSet;
229
230 /*  Electrostatic potential  */
231 typedef struct Elpot {
232         Vector pos;
233         Double esp;
234 } Elpot;
235
236 /*  Molecule record  */
237 typedef struct Molecule {
238         Object base;
239         Int    natoms;
240         Atom   *atoms;
241         Int    nbonds;
242         Int    *bonds;       /*  The size of array is 2*nbonds  */
243         Int    nangles;
244         Int    *angles;      /*  The size of array is 3*nangles  */
245         Int    ndihedrals;
246         Int    *dihedrals;   /*  The size of array is 4*ndihedrals  */
247         Int    nimpropers;
248         Int    *impropers;   /*  The size of array is 4*nimpropers  */
249         Int    nresidues;    /*  Number of residues; maximum residue number + 1 (because residue 0 is 'undefined residue')  */
250         char   (*residues)[4];
251         XtalCell   *cell;
252 /*      char   is_xtal_coord; *//*  True if the coordinates are measured in crystallographic units  */
253         Int    nsyms;        /*  Symmetry operations; syms are always described in crystallographic units (even when is_xtal_coord is false)  */
254         Transform *syms;
255 /*      PeriodicBox *box;    *//*  Periodic box  */
256         IntGroup *selection;
257         Int    nexatoms;
258         ExAtom *exatoms;
259         Int    nexbonds;
260         Int    *exbonds;     /*  The size of array is 2*nbonds; Atom index >= 0 : base atoms, < 0 : expanded atoms at index -exbonds[n]-1  */
261         Int    nframes;      /*  The number of frames (>= 1). This is a cached value, and should be
262                                                          recalculated from the atoms if it is -1  */
263         Int    cframe;       /*  The current frame number  */
264
265         Vector *frame_cells; /*  The cell vectors for frames; (nframes*3) array of Vectors  */
266
267         struct MainView *mview;  /*  Reference to the MainView object if present (no retain)  */
268         Int    modifyCount;  /*  Internal counter for modification. This value is not to be modified
269                                  manually; instead, call MoleculeIncrementModifyCount() whenever
270                                                      modification is done, which also takes care necessary notification
271                                                          to the other part of the application (system dependent)  */
272
273         struct MDArena *arena;  /*  Reference to the MDArena record during MM/MD run (no retain)  */
274
275         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) */
276
277         /*  Information from the dcd files  */
278         Int    startStep;     /*  the timestep for frame 0  */
279         Int    stepsPerFrame; /*  the number of timesteps between neighboring frames  */
280         Double psPerStep;     /*  picosecond per step  */
281
282         /*  Information for basis sets and MOs  */
283         BasisSet *bset;
284         
285         /*  Electrostatic potential  */
286         Int    nelpots;
287         Elpot  *elpots;
288
289         /*  Parameters specific for this molecule  */
290         struct Parameter *par;
291         
292         /*  Flag to request rebuilding MD internal information  */
293         Byte   needsMDRebuild;
294         
295         /*  Flag to clear selection of the parameter table  */
296         Byte   parameterTableSelectionNeedsClear;
297         
298         /*  Flag to request copying coordinates to MD arena  */
299         Byte   needsMDCopyCoordinates;
300
301         /*  Prohibit modification of the topology (to avoid interfering MD) */
302         Byte   noModifyTopology;
303         
304         /*  Flag to request aborting a subthread  */
305         Byte   requestAbortThread;
306
307         /*  Flag to signal that a subthread is terminated  */
308         Byte   threadTerminated;
309
310         /*  Mutex object. If non-NULL, it should be locked before modifying molecule  */
311         void *mutex;
312         
313         /*  Ruby pointer (see ruby_bind.c)  */
314         void *exmolobj;
315
316 } Molecule;
317
318 int strlen_limit(const char *s, int limit);
319
320 Molecule *MoleculeNew(void);
321 int MoleculeLoadFile(Molecule *mp, const char *fname, const char *ftype, char *errbuf, int errbufsize);
322 int MoleculeLoadPsfFile(Molecule *mp, const char *fname, char *errbuf, int errbufsize);
323 int MoleculeLoadTepFile(Molecule *mp, const char *fname, char *errbuf, int errbufsize);
324 int MoleculeLoadShelxFile(Molecule *mp, const char *fname, char *errbuf, int errbufsize);
325 int MoleculeLoadGaussianFchkFile(Molecule *mp, const char *fname, char *errbuf, int errbufsize);
326 int MoleculeLoadMbsfFile(Molecule *mp, const char *fname, char *errbuf, int errbufsize);
327 Molecule *MoleculeNewWithName(const char *name);
328 Molecule *MoleculeInitWithAtoms(Molecule *mp, const Atom *atoms, int natoms);
329 Molecule *MoleculeInitWithMolecule(Molecule *mp2, const Molecule *mp);
330 void MoleculeSetName(Molecule *par, const char *name);
331 const char *MoleculeGetName(Molecule *mp);
332 void MoleculeSetPath(Molecule *mol, const char *fname);
333 const char *MoleculeGetPath(Molecule *mol);
334 Molecule *MoleculeWithName(const char *name);
335 Molecule *MoleculeRetain(Molecule *mp);
336 void MoleculeRelease(Molecule *mp);
337 extern void MoleculeReleaseExternalHook(Molecule *mol);
338 void MoleculeExchange(Molecule *mp1, Molecule *mp2);
339
340 int MoleculeAddGaussianOrbitalShell(Molecule *mol, Int sym, Int nprims, Int a_idx);
341 int MoleculeAddGaussianPrimitiveCoefficients(Molecule *mol, Double exponent, Double contraction, Double contraction_sp);
342 int MoleculeSetMOCoefficients(Molecule *mol, Int idx, Double energy, Int ncomps, Double *coeffs);
343 int MoleculeAllocateBasisSetRecord(Molecule *mol, Int rflag, Int ne_alpha, Int ne_beta);
344
345 void MoleculeIncrementModifyCount(Molecule *mp);
346 void MoleculeClearModifyCount(Molecule *mp);
347
348 MolEnumerable *MolEnumerableNew(Molecule *mol, int kind);
349 void MolEnumerableRelease(MolEnumerable *mseq);
350 AtomRef *AtomRefNew(Molecule *mol, int idx);
351 void AtomRefRelease(AtomRef *aref);
352
353 void MoleculeSetCell(Molecule *mp, Double a, Double b, Double c, Double alpha, Double beta, Double gamma, int convertCoordinates);
354 void MoleculeSetAniso(Molecule *mp, int n1, int type, Double x11, Double x22, Double x33, Double x12, Double x23, Double x31, const Double *sigmap);
355 int MoleculeSetPeriodicBox(Molecule *mp, const Vector *ax, const Vector *ay, const Vector *az, const Vector *ao, const char *periodic);
356
357 int MoleculeReadCoordinatesFromFile(Molecule *mp, const char *fname, const char *ftype, char *errbuf, int errbufsize);
358 int MoleculeReadCoordinatesFromPdbFile(Molecule *mp, const char *fname, char *errbuf, int errbufsize);
359 int MoleculeReadCoordinatesFromDcdFile(Molecule *mp, const char *fname, char *errbuf, int errbufsize);
360
361 int MoleculeReadExtendedInfo(Molecule *mp, const char *fname, char *errbuf, int errbufsize);
362 int MoleculeWriteExtendedInfo(Molecule *mp, const char *fname, char *errbuf, int errbufsize);
363
364 int MoleculeWriteToFile(Molecule *mp, const char *fname, const char *ftype, char *errbuf, int errbufsize);
365 int MoleculeWriteToPsfFile(Molecule *mp, const char *fname, char *errbuf, int errbufsize);
366 int MoleculeWriteToPdbFile(Molecule *mp, const char *fname, char *errbuf, int errbufsize);
367 int MoleculeWriteToDcdFile(Molecule *mp, const char *fname, char *errbuf, int errbufsize);
368 int MoleculeWriteToTepFile(Molecule *mp, const char *fname, char *errbuf, int errbufsize);
369 void MoleculeDump(Molecule *mol);
370
371 int MoleculePrepareMDArena(Molecule *mol, int check_only, char **retmsg);
372
373 char *MoleculeSerialize(Molecule *mp, Int *outLength, Int *timep);
374 Molecule *MoleculeDeserialize(const char *data, Int length, Int *timep);
375
376 void MoleculeCleanUpResidueTable(Molecule *mp);
377 int MoleculeChangeNumberOfResidues(Molecule *mp, int nresidues);
378 int MoleculeChangeResidueNumberWithArray(Molecule *mp, IntGroup *group, Int *resSeqs);
379 int MoleculeChangeResidueNumber(Molecule *mp, IntGroup *group, int resSeq);
380 int MoleculeOffsetResidueNumbers(Molecule *mp, IntGroup *group, int offset, int nresidues);
381 int MoleculeChangeResidueNames(Molecule *mp, int argc, Int *resSeqs, char *names);
382 int MoleculeMaximumResidueNumber(Molecule *mp, IntGroup *group);
383 int MoleculeMinimumResidueNumber(Molecule *mp, IntGroup *group);
384
385 int MoleculeCreateAnAtom(Molecule *mp, const Atom *ap, int pos);
386 int MoleculeMerge(Molecule *dst, Molecule *src, IntGroup *where, int resSeqOffset);
387 int MoleculeUnmerge(Molecule *src, Molecule **dstp, IntGroup *where, int resSeqOffset);
388 int MoleculeExtract(Molecule *src, Molecule **dstp, IntGroup *where, int dummyFlag);
389 int MoleculeAddBonds(Molecule *mp, Int nbonds, const Int *bonds);
390 int MoleculeDeleteBonds(Molecule *mp, Int nbonds, const Int *bonds);
391 int MoleculeAddAngles(Molecule *mp, const Int *angles, IntGroup *where);
392 int MoleculeDeleteAngles(Molecule *mp, Int *angles, IntGroup *where);
393 int MoleculeAddDihedrals(Molecule *mp, const Int *dihedrals, IntGroup *where);
394 int MoleculeDeleteDihedrals(Molecule *mp, Int *dihedrals, IntGroup *where);
395 int MoleculeAddImpropers(Molecule *mp, const Int *impropers, IntGroup *where);
396 int MoleculeDeleteImpropers(Molecule *mp, Int *impropers, IntGroup *where);
397 int MoleculeLookupBond(Molecule *mp, Int n1, Int n2);
398 int MoleculeLookupAngle(Molecule *mp, Int n1, Int n2, Int n3);
399 int MoleculeLookupDihedral(Molecule *mp, Int n1, Int n2, Int n3, Int n4);
400 int MoleculeLookupImproper(Molecule *mp, Int n1, Int n2, Int n3, Int n4);
401
402 /*
403 Int     MoleculeReplaceAllAngles(Molecule *mol, Int nangles, const Int *angles, Int **outAngles);
404 Int MoleculeReplaceAllDihedrals(Molecule *mol, Int ndihedrals, const Int *dihedrals, Int **outDihedrals);
405 Int MoleculeReplaceAllImpropers(Molecule *mol, Int nimpropers, const Int *impropers, Int **outImpropers);
406
407 Int MoleculeFindAllAngles(Molecule *mol, Int **outAngles);
408 Int MoleculeFindAllDihedrals(Molecule *mol, Int **outDihedrals);
409 Int MoleculeFindAllImpropers(Molecule *mol, Int **outImpropers);
410 */
411
412 Int MoleculeFindMissingAngles(Molecule *mol, Int **outAngles);
413 Int MoleculeFindMissingDihedrals(Molecule *mol, Int **outDihedrals);
414 Int MoleculeFindMissingImpropers(Molecule *mol, Int **outImpropers);
415         
416 IntGroup *MoleculeSearchBondsIncludingAtoms(Molecule *mp, IntGroup *atomgroup);
417 IntGroup *MoleculeSearchAnglesIncludingAtoms(Molecule *mp, IntGroup *atomgroup);
418 IntGroup *MoleculeSearchDihedralsIncludingAtoms(Molecule *mp, IntGroup *atomgroup);
419 IntGroup *MoleculeSearchImpropersIncludingAtoms(Molecule *mp, IntGroup *atomgroup);
420
421 IntGroup *MoleculeSearchBondsAcrossAtomGroup(Molecule *mp, IntGroup *atomgroup);
422
423 IntGroup *MoleculeSearchAnglesIncludingBond(Molecule *mp, int n1, int n2);
424 IntGroup *MoleculeSearchDihedralsIncludingBond(Molecule *mp, int n1, int n2);
425 IntGroup *MoleculeSearchImpropersIncludingBond(Molecule *mp, int n1, int n2);
426
427 int MoleculeLookupAtomInResidue(Molecule *mp, int n1, int resno);
428 int MoleculeAnalyzeAtomName(const char *s, char *resName, int *resSeq, char *atomName);
429 int MoleculeAtomIndexFromString(Molecule *mp, const char *s);
430
431 int MoleculeFindCloseAtoms(Molecule *mp, Int index, Double limit, Int *outNbonds, Int **outBonds, Int triangle);
432 int MoleculeGuessBonds(Molecule *mp, Double limit, Int *outNbonds, Int **outBonds);
433 int MoleculeAreAtomsConnected(Molecule *mp, int n1, int n2);
434 int MoleculeRebuildTablesFromConnects(Molecule *mp);
435
436 void MoleculeGetAtomName(Molecule *mp, int index, char *buf, int bufsize);
437
438 void MoleculeSetSelection(Molecule *mp, IntGroup *select);
439 IntGroup *MoleculeGetSelection(Molecule *mp);
440 void MoleculeSelectAtom(Molecule *mp, int n1, int extending);
441 void MoleculeUnselectAtom(Molecule *mp, int n1);
442 void MoleculeToggleSelectionOfAtom(Molecule *mp, int n1);
443 int MoleculeIsAtomSelected(Molecule *mp, int n1);
444 int MoleculeIsBondSelected(Molecule *mp, int n1, int n2);
445 IntGroup *MoleculeModifySelectionByRemovingAtoms(Molecule *mp, IntGroup *selection, IntGroup *remove);
446
447 int MoleculeTransformBySymop(Molecule *mp, const Vector *vpin, Vector *vpout, Symop symop);
448 int MoleculeAddExpandedAtoms(Molecule *mp, Symop symop, IntGroup *group);
449 int MoleculeAmendBySymmetry(Molecule *mp, IntGroup *group, IntGroup **groupout, Vector **vpout);
450
451 int MoleculeShowAllAtoms(Molecule *mp);
452 int MoleculeShowReverse(Molecule *mp);
453 int MoleculeHideAtoms(Molecule *mp, IntGroup *ig);
454
455 int MoleculeRenumberAtoms(Molecule *mp, const Int *new2old, Int *old2new_out, Int isize);
456
457 void MoleculeTransform(Molecule *mp, Transform tr, IntGroup *group);
458 void MoleculeMove(Molecule *mp, Transform tr, IntGroup *group);
459 void MoleculeTranslate(Molecule *mp, const Vector *vp, IntGroup *group);
460 void MoleculeRotate(Molecule *mp, const Vector *axis, Double angle, const Vector *center, IntGroup *group);
461 void MoleculeReaxis(Molecule *mp, const Vector *xaxis, const Vector *yaxis, const Vector *zaxis, IntGroup *group);
462 int MoleculeCenterOfMass(Molecule *mp, Vector *center, IntGroup *group);
463 int MoleculeBounds(Molecule *mp, Vector *min, Vector *max, IntGroup *group);
464         
465 Int *MoleculeSearchEquivalentAtoms(Molecule *mol, IntGroup *ig);
466         
467 void MoleculeAddExpansion(Molecule *mp, Vector dr, Int symop, IntGroup *group, Double limit);
468 void MoleculeClearExpansion(Molecule *mp, IntGroup *group);
469 void MoleculeRemoveExpansion(Molecule *mp, Vector dr, Int symop, IntGroup *group);
470 void MoleculeAutoExpansion(Molecule *mp, const float *boxstart, const float *boxend, IntGroup *group, Double limit);
471
472 void MoleculeXtalToCartesian(Molecule *mp, Vector *dst, const Vector *src);
473 void MoleculeCartesianToXtal(Molecule *mp, Vector *dst, const Vector *src);
474 Double MoleculeMeasureBond(Molecule *mp, const Vector *vp1, const Vector *vp2);
475 Double MoleculeMeasureAngle(Molecule *mp, const Vector *vp1, const Vector *vp2, const Vector *vp3);
476 Double MoleculeMeasureDihedral(Molecule *mp, const Vector *vp1, const Vector *vp2, const Vector *vp3, const Vector *vp4);
477
478 IntGroup *MoleculeFragmentExcludingAtomGroup(Molecule *mp, int n1, IntGroup *exatoms);
479 IntGroup *MoleculeFragmentExcludingAtoms(Molecule *mp, int n1, int argc, int *argv);
480 IntGroup *MoleculeFragmentWithAtomGroups(Molecule *mp, IntGroup *inatoms, IntGroup *exatoms);
481 int MoleculeIsFragmentDetachable(Molecule *mp, IntGroup *group, int *n1, int *n2);
482 int MoleculeIsFragmentRotatable(Molecule *mp, IntGroup *group, int *n1, int *n2, IntGroup **rotGroup);
483
484 int MoleculeGetNumberOfFrames(Molecule *mp);
485 int MoleculeInsertFrames(Molecule *mp, IntGroup *group, const Vector *inFrame, const Vector *inFrameCell);
486 int MoleculeRemoveFrames(Molecule *mp, IntGroup *group, Vector *outFrame, Vector *outFrameCell);
487 /*int MoleculeInsertFrame(Molecule *mp, int index, const Vector *inFrame);
488 int MoleculeRemoveFrame(Molecule *mp, int frame, Vector *outFrame); */
489 int MoleculeSelectFrame(Molecule *mp, int frame, int copyback);
490
491 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);
492 int MoleculeGetDefaultMOGrid(Molecule *mp, Int npoints, Vector *op, Vector *xp, Vector *yp, Vector *zp, Int *nx, Int *ny, Int *nz);
493 const Cube *MoleculeGetCubeAtIndex(Molecule *mp, Int index);
494 int MoleculeLookUpCubeWithMONumber(Molecule *mp, Int mono);
495 int MoleculeClearCubeAtIndex(Molecule *mp, Int index);
496 int MoleculeOutputCube(Molecule *mp, Int index, const char *fname, const char *comment);
497
498 #define kMoleculePasteboardType "Molecule"
499 #define kParameterPasteboardType "Parameter"
500
501 STUB int MoleculeCallback_writeToPasteboard(const char *type, const void *data, int length);
502 STUB int MoleculeCallback_readFromPasteboard(const char *type, void **dptr, int *length);
503 STUB int MoleculeCallback_isDataInPasteboard(const char *type);
504
505 STUB Molecule *MoleculeCallback_openNewMolecule(const char *fname);
506 STUB void MoleculeCallback_notifyModification(Molecule *mp, int now_flag);
507 STUB Molecule *MoleculeCallback_currentMolecule(void);
508 STUB Molecule *MoleculeCallback_moleculeAtIndex(int idx);
509 STUB Molecule *MoleculeCallback_moleculeAtOrderedIndex(int idx);
510 STUB void MoleculeCallback_displayName(Molecule *mol, char *buf, int bufsize);
511 STUB void MoleculeCallback_pathName(Molecule *mol, char *buf, int bufsize);
512 STUB int MoleculeCallback_setDisplayName(Molecule *mol, const char *name);
513
514 STUB void MoleculeCallback_lockMutex(void *mutex);
515 STUB void MoleculeCallback_unlockMutex(void *mutex);
516 STUB void MoleculeCallback_cannotModifyMoleculeDuringMDError(Molecule *mol);
517
518 /*  This is also defined in Molby_extern.h, but it may be called from functions in Molecule.c  */
519 STUB int MyAppCallback_checkInterrupt(void);
520         
521 void MoleculeLock(Molecule *mol);
522 void MoleculeUnlock(Molecule *mol);
523
524 #if 0
525 #define __MoleculeLock(mol) MoleculeLock(mol)
526 #define __MoleculeUnlock(mol) MoleculeUnlock(mol)
527 #else
528 #define __MoleculeLock(mol)
529 #define __MoleculeUnlock(mol)
530 #endif
531         
532 #ifdef __cplusplus
533 }
534 #endif
535                 
536 #endif /* __Molecule_h__ */