OSDN Git Service

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