OSDN Git Service

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