
		    AGBNP portable library (v1.1)
		    =============================

			  Emilio Gallicchio
		       emilio@hpcp.rutgers.edu

		   Last Modified: October 2005

 License Agreement
 -----------------

 Copyright (c) Emilio Gallicchio, Ronald M. Levy 2004

 This software (AGBNP v1.1) was developed by Emilio Gallicchio and
 Ronald M. Levy at Rutgers University, 610 Taylor Rd., Piscataway, NJ
 08854, U.S.A.
  
 Installation, use, reproduction, display, modification and
 redistribution of this software, with or without modification, in
 source and binary forms, is permitted only with the consent of the
 authors and copyright holders.

 Any academic report, publication, or other academic disclosure of
 results obtained with this software will acknowledge this software's
 use by the following citation:

    E. Gallicchio, R. M. Levy, AGBNP: An Analytic Implicit Solvent
    Model Suitable for Molecular Dynamics Simulations and
    High-Resolution Modeling, J. Comput. Chem., 25, 479-499 (2004).


 DISCLAIMER

 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
 "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
 LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY, OF
 SATISFACTORY QUALITY, AND FITNESS FOR A PARTICULAR PURPOSE OR USE ARE
 DISCLAIMED. THE COPYRIGHT HOLDERS AND CONTRIBUTORS MAKE NO
 REPRESENTATION THAT THE SOFTWARE, MODIFICATIONS, ENHANCEMENTS OR
 DERIVATIVE WORKS THEREOF, WILL NOT INFRINGE ANY PATENT, COPYRIGHT,
 TRADEMARK, TRADE SECRET OR OTHER PROPRIETARY RIGHT.

 LIMITATION OF LIABILITY

 THE COPYRIGHT HOLDERS AND CONTRIBUTORS AND ANY OTHER OFFICER, AGENT,
 OR EMPLOYEE OF THE UNIVERSITY SHALL HAVE NO LIABILITY TO LICENSEE OR
 OTHER PERSONS FOR DIRECT, INDIRECT, SPECIAL, INCIDENTAL,
 CONSEQUENTIAL, EXEMPLARY, OR PUNITIVE DAMAGES OF ANY CHARACTER
 INCLUDING, WITHOUT LIMITATION, PROCUREMENT OF SUBSTITUTE GOODS OR
 SERVICES, LOSS OF USE, DATA OR PROFITS, OR BUSINESS INTERRUPTION,
 HOWEVER CAUSED AND ON ANY THEORY OF CONTRACT, WARRANTY, TORT
 (INCLUDING NEGLIGENCE), PRODUCT LIABILITY OR OTHERWISE, ARISING IN
 ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
 POSSIBILITY OF SUCH DAMAGES.

 Introduction
 ------------

 AGBNP is an analytical implicit solvent model based on the pairwise
 descreening (PD) Generalized Born (GB) model and a non-polar
 solvation free energy (NP) estimator which takes into account
 independently the work of cavity formation and the solute-solvent van
 der Waals interaction energy. The model and its derivation are
 described in detail in the following paper:
 
 E. Gallicchio, R. M. Levy, AGBNP: An Analytic Implicit Solvent Model
 Suitable for Molecular Dynamics Simulations and High-Resolution
 Modeling, J. Comput. Chem., 25, 479-499 (2004),
 
 Only a brief overview is given here, just what is indispensable to
 better illustrate the library interface (API).
 
 In the Generalized Born model each atom is assigned a Born radius
 (Bi). The Born radius measures in some sense how buried an atom is in
 the solute. The more buried an atom the larger its Born radius. In
 the Pairwise Descreening (PD) scheme the Born radius (or to be precise
 its inverse) of atom i is computed as the sum of a kernel function Qji
 over the atom neighbors j:
 
 1/Bi = 1/Ri - Sum_j s_ji Qji                                   (1)
 
 Qji depends on the distance between atoms i and j and their van der
 Waals radii. In Eq.(1) Ri is the van der Waals radius of atom i and
 s_ji < 1 is a scaling coefficient that takes into account the volume
 of atom j that overlaps with other atoms and therefore should not be
 fully counted toward the amount of solvent volume displaced by atom
 j. AGBNP is unique among PD-GB models in that the overlap scaling
 coefficients depend on solute conformation and are computed from
 purely geometric considerations. As such they are not parameterized
 with respect to experimental and Poisson Boltzmann data. Hydrogen
 atoms do not contribute to descreening. That is Qji = 0 if j is an
 hydrogen atom.
 
 Once the Born radii are assigned the electrostatic solvation free
 energy (DGelec) is computed using the GB equation:
 
 DGelec = - (1/ein - 1/eout ) Sum_ij qi qj/fij                  (2)
 
 where eout and ein are the interior and exterior dielectric constants
 respectively, qi is the partial charge of atom i and fij is the GB
 distance function that depends on the distance, rij, between atoms i
 and j and the geometric average, Bij, of the Born radii of atoms i and
 j:
 
 fij = sqrt{ rij^2 + Bij^2 exp(-rij^2/Bij^2) }.                 (3)
 
 The non-polar hydration free energy estimator is composed of two
 terms. The first, related to the cavity hydration free energy
 (DGcav), is proportional to the solute surface area of each atom
 
 DGcav = Sum_i gamma_i A_i                                     (4)
 
 where gamma_i is a surface tension parameter that depend on the atom
 type of each atom, and A_i is the surface area of atom i. The surface
 area is defined as the van der Waals surface area obtained by
 increasing the van der Waals radius of each atom by 0.5
 Angstroms. The surface area of each atom is calculated using an
 analytical model based on the same method used to calculate overlap
 scaling factors. Hydrogen atoms do not contribute to the surface
 area, that is they can be thought as of atoms of zero radius in this
 respect.

 The second component of the non-polar hydration free energy estimator
 is a solute-solvent van der Waals interaction energy term
 (DGvdw). This is expressed as the sum of individual atomic terms
 (DGvdw_i). The solute-solvent van der Waals interaction energy of
 atom i, DGvdw_i, is expressed in terms of the Born radius Bi:

 DGvdw = Sum_i DGvdw_i = 
 Sum_i { alpha_i (-16 pi /3) rho epsi_i sigma_i^6/(Bi + Rw)^3 
                                                   + delta_i } (5)

 where alpha_i and delta_i are adjustable parameters for each atom
 type, rho=0.033428 1/Angstrom^3 is the number density of water at
 standard conditions, sigma_i and epsi_i are the Lennard-Jones
 parameters for the interaction energy of atom i with water, and
 Rw=1.4 Angstroms is the approximate radius of a water molecule.

 The AGBNP model is fully analytical with first derivatives. It
 requires the partial charge, van der Waals radius, and the non-polar
 parameters gamma_i and

 a_i = alpha_i (-16 pi /3) rho epsi_i sigma_i^6                (6)

 of each atom.

 Parameterization
 ----------------

 The AGBNP model is based on several parameters, some of which
 are derived from the underlying force field either without change
 (partial charges) or with small modifications (van der Waals radii,
 see original AGBNP reference). The non-polar parameters however need
 to be determined anew. The original AGBNP paper [J. Comp. Chem., 25,
 479-499 (2004)] describes the parameterization procedure of the AGBNP
 model with the OPLS-AA force field (1999 version). The list of the latest
 parameters of the OPLS-AA/AGBNP model is included in this
 distribution in the file "agbnp.param".

 The format of the "agbnp.param" file is as follows:

   Column   Content
  --------------------------------------------------------------------
     1       OPLS symbolic type
     2       Value of the van der Waals radius [Ang]
     3       Value of the pure gamma parameter [(kcal/mol)/Ang^2]
     4       Value of the pure alpha parameter [dimensionless]
     5       Value of the pure delta parameter [kcal/mol]
     6       Value of the correction gamma parameter [(kcal/mol)/Ang^2]
     7       Value of the correction alpha parameter [dimensionless]
     8       Value of the correction delta parameter [kcal/mol]
     9       Value of the screening parameter [dimensionless]

 lines that begin with '#' are comments.
 
 We do not support nor recommend the use of AGBNP with any force field
 other than OPLS-AA with the parameters listed in the included
 "agbnp.param" file. You are free to derive AGBNP parameters for other
 force fields. However, if you do, you are on your own.

 Pre-compiled binary distribution:
 --------------------------------

 Contents:

 The distribution includes two libraries and their corresponding
 header files:

 libagbnp.a agbnp.h
 libnblist.a libnblist.h

 plus a parameter file

 agbnp.param

 Platform supported:

 Linux-x86

 Installation:

 Unpack the distribution and copy the header files and libraries
 to a location of your choice.

 Usage
 -----

 Compilation:

 To link your code to the AGBNP add the locations to the search path
 of the compiler. For example if the header files are in
 /agbnp/include/ and the libraries are in /agbnp/lib/ this is an
 example of the compilation command:

 gcc -O -I/agbnp/include/ mycode.c -o mycode \
 -L/agbnp/lib/ -lagbnp -lnblist


 AGBNP C API
 -----------

 The header file agbnp.h must be included to access the AGBNP API functions.

 #include "agbnp.h"


 > int agbnp_initialize( void );
 
 Initializes libagbnp library. No library functions can be used until
 this call is made.
 
 Return values:
 AGBNP_OK - library initialized OK.
 AGBNP_ERR - error in initialization.
 
 > void agbnp_terminate( void );
 
 Terminates libagbnp library. Frees all associated storage.  Once this
 call is made then the library can no longer be used until it is
 initialized again by calling agbnp_initialize().
 
 > int agbnp_new(int *tag, int natoms, double *x, double *y, 
        double *z, double *r, double *charge, 
        double dielectric_in, double dielectric_out, 
        double *igamma, double *sgamma, 
        double *ia, double *sa, 
        double *idelta, double *sdelta, 
        double *ab, int nhydrogen, int *ihydrogen, 
        int ndummy, int *idummy, int *isfrozen, 
        NeighList *neigh_list, NeighList *excl_neigh_list);
 
 Creates a new public instance of an AGBNP structure. If successful
 this function returns in the 'tag' argument a non-negative integer id
 that identifies the instance. All subsequent operations on this
 instance should reference this tag. Any number of AGBNP instances can
 be generated calling the agbnp_new() function. Each instance is
 independent from the others, that is each may correspond to different
 systems, parameters, etc.
 
 Input parameters (units in []):
 
 natoms:   number of atoms of the system. Atoms are numbered
           from 0 to natoms-1.
 
 x, y, z:  Cartesian coordinates of each atom [Angstroms].
 
 charge:   partial charge of each atom [e].
 
 dielectric_in:  value of the interior (solute) relative dielectric
                 constant.
 
 dielectric_out: value of the exterior (solvent) relative dielectric
                 constant.
 
 igamma:   value of the gamma non-polar parameter for each atom 
           [(kcal/mol)/Ang^2]
 
 sgamma:   value of the gamma non-polar correction parameter for each
           atom [(kcal/mol)/Ang^2]
  
 ia:       value of the "a" van der Waals non-polar parameter for each
           atom [(kcal/mol) Ang^3] (Eq.6)
 
 sa:       value of the "a" van der Waals non-polar correction
           parameter for each atom [(kcal/mol) Ang^3]
 
 idelta:   value of the delta non-polar parameter for each atom [kcal/mol]
 
 sdelta:   value of the delta non-polar correction parameter for each
           atom [kcal/mol]
 
 Note:  The values of the non-polar parameters used internally are the
 sum of the ideal and correction values. However the non-polar energy
 derived from each is reported separately as a pure non-polar energy
 and a correction energy term. The correction energy term has the same
 expression as the non-polar estimator given above (this could change
 in the future) but it is calculated using the set of correction
 parameters rather than the ideal non-polar parameters. If you are
 confused about this just set the correction parameters to zero.
 
 ab:       value of the screening GB pair energy parameter for each
           atom. This feature is experimental (see [Felts, A.K.,
           Y. Harano, E. Gallicchio, R.M. Levy.  Free energy surfaces
           of beta-hairpin and alpha-helical peptides generated by
           replica exchange molecular dynamics with the AGBNP implicit
           solvent model. PROTEINS: Structure, Function, and
           Bioinformatics, 56, 310-321 (2004)]). The standard value is
           1 for each atom.
 
 nhydrogen: number of hydrogen atoms in the system.
 
 ihydrogen: atom index of each hydrogen atom. Again, atom indexes start
            from zero.
 
 ndummy:   number of dummy atoms in the system (currently ignored)
 
 idummy:   atom index of each dummy atom (currently ignored)
 
 isfrozen: a flag for each atom, 1 if the atom is frozen and 0
           otherwise (free atom). If no atoms are frozen can also set
           this pointer to NULL. CPU time-saving features are turned
           on if isfrozen is not NULL and at least one element is not
           zero. Also see agbnp_init_frozen().
 
 neigh_list: if not NULL instructs AGBNP to use this Verlet neighbor
             list to define the neighbors of each atom. neigh_list is
             a Verlet neighbor list as defined by the libnblist
             library - included in the AGBNP distribution. For
             documentation pertaining to libnblist see below. The
             neighbor list is used to define atom neighbors for the
             surface areas, Born radii, and GB pair energy
             calculation. If neigh_list is not used (set to NULL) a
             neighbor list is not used and every pair of atoms is
             considered a neighbor pair. Not using a neighbor list can
             increase significantly the CPU time and memory
             consumption especially for large systems. It is the
             responsibility of the caller to properly create and
             update the neighbor list.
 
 excl_neigh_list: if not NULL instructs AGBNP to use the Verlet
                  neighbor list pointed by excl_neigh_list as an
                  additional source of neighbors for each atom. In
                  molecular mechanics programs the non-bonded neighbor
                  list does not contain excluded atoms (such as pairs
                  of covalently linked atoms) and sometime does not
                  contain pairs of frozen atoms. Instead AGBNP
                  requires that the list of neighbors includes every
                  atom neighbor including non-bonded excluded atoms or
                  frozen atoms. The list of excluded atoms for each
                  atom is conveniently stored in a Verlet neighbor
                  list format and provided to AGBNP using this
                  parameter.

 Note: The arrays of parameters provided to AGBNP (such as 'x', 'r',
 'igamma', etc.) are copied by agbnp_new() to internal AGBNP data
 structures. The caller is allowed to free storage associated with
 these arrays. The exception to this rule are the neighbor
 lists. AGBNP saves the memory pointers to the neighbor lists
 structures and uses these pointers to access the neighbor
 list. Therefore the memory location of the neighbor lists should not
 be changed by the caller during the life of the AGBNP instance. Also
 note that the agbnp_terminate() and agbnp_delete() functions do not
 free memory storage of the Verlet neighbor lists. It is the
 responsibility of the caller to manage the memory associated with the
 neighbor lists.
 
 Return values:
 AGBNP_OK - AGBNP instance created, id tag is returned in 'tag'.
 AGBNP_ERR - error creating AGBNP instance. Consult error message
             on stderr.
 
 
 > int agbnp_delete(int tag);
 
 Deletes the AGBNP instance referenced by 'tag'. Associated memory is freed.
 Return values:
 AGBNP_OK - AGBNP instance deleted.
 AGBNP_ERR - error deleting AGBNP instance. Consult error message
             on stderr.generalized born and van der Waals energies
 
 
 > int agbnp_agb_energy(int tag, double *x, double *y, double *z,
        double *sp, double *br, double *egb, double (*dgbdr)[3],
        double *evdw, double (*dvwdr)[3], double *ecorr);
 
 Calculates the Generalized Born and solute-solvent van der Waals
 interaction energy components of the AGBNP instance referenced by
 'tag'. Current Cartesian coordinates (in Angstroms) for each atom are
 provided by the caller in arrays 'x', 'y' and 'z'. The parameters
 specified in the call to agbnp_new() (such as van der Waals radii and
 non-polar parameters) are used. The overlap scaling coefficient for
 each atom are returned in the array 'sp' and the Born radii in array
 'br' (that should be of size natoms or larger). The GB energy in
 kcal/mol is returned in 'egb', the solute-solvent van der Waals
 interaction energy in kcal/mol is returned in 'evdw', and the
 corresponding correction energy in 'ecorr'. The gradient of the GB
 energy for each atom is returned in 'dgbdr' and the gradient of the
 van der Waals energy in 'dvwdr' (these arrays should be of size natoms
 - 3*natoms*sizeof(double) bytes - or larger, dgbdr[i][0] is the
 x-component of the GB gradient with respect to atom i, the
 corresponding y and z components are dgbdr[i][1] and dgbdr[i][2],
 respectively), similarly for dvwdr. Gradients are given in units of
 (kcal/mol)/Angstrom. The van der Waals gradient includes the gradient
 of the correction energy. The force on atom i is the negative of the
 gradient.
 
 Return values:
 AGBNP_OK - AGBNP energy calculated.
 AGBNP_ERR - Error calculating energy. Consult error message
             on stderr.
 
 > int agbnp_cavity_energy(int tag, double *x, double *y, double *z,
        double *mol_volume, double *surf_area,
        double *ecav, double *ecorr, double (*decav)[3]);

 Calculates the free energy of cavity formation component of the AGBNP
 instance referenced by 'tag'. Current Cartesian coordinates (in
 Angstroms) for each atom are provided by the caller in arrays 'x', 'y'
 and 'z'. The parameters specified in the call to agbnp_new() (such as
 van der Waals radii and surface tension parameters are used). This
 routine uses the van der Waals radii augmented by 0.5 Angstroms. The
 molecular volume is returned in 'mol_volume' (in Angstrom^3). The
 surface area of each atom in square Angstrom is returned in the array
 surf_area (that should be of size natoms or larger). The cavity free
 energy in kcal/mol is returned in 'ecav', the corresponding correction
 energy in 'ecorr'. The total gradient in (kcal/mol)/Angstrom of
 ecav+ecorr is returned in 'decav' similarly as for the GB and van der
 Waals energy gradients described above.
 
 Return values:
 AGBNP_OK - AGBNP cavity energy calculated.
 AGBNP_ERR - Error calculating cavity energy. Consult error message
             on stderr.
 
 > int agbnp_bornr(int tag, double *x, double *y, double *z,
        double *sp, double *br);
 
 Same as agbnp_agb_energy() but returns Born radii and overlap scaling
 coefficients only, no energies, no derivatives.
 
 Return values:
 AGBNP_OK - AGBNP Born radii and overlap scaling coefficients calculated.
 AGBNP_ERR - Error calculating Born radii and overlap scaling
             coefficients. Consult error messages on stderr.
 
 > int agbnp_init_frozen(int tag, double *x, double *y, double *z, 
        int *isfrozen);
 
 Updates constant terms cache of the AGBNP instance referenced by
 'tag'. When using frozen atoms AGBNP attempts to save CPU time by
 saving the value of quantities that are deemed constant. The values of
 the constant quantities depends on the position of the free atoms
 relative to the frozen atoms and therefore, due to atomic motion, it
 is necessary to call this function periodically to update them. It is
 the responsibility of the caller to call this routine whenever
 needed. Calling this routine too infrequently would cause poor energy
 conservation in constant energy MD simulations at best, and in some cases
 it may lead to unphysical behavior (such as atom overlaps, etc.). On
 the other hand calling it too frequently wastes CPU time defying the
 purpose of the frozen atom scheme to save CPU time. A good
 general-purpose scheme is to call agbnp_init_frozen() whenever
 neighbor lists are updated. Current Cartesian coordinates (in
 Angstroms) for each atom are provided by the caller in arrays 'x', 'y'
 and 'z'. The list of frozen atoms, specified by the array 'isfrozen'
 as for agbnp_new(), can be different from the list provided to
 agbnp_new() or previous calls to agbnp_init_frozen(). Therefore this
 routine can be also used to update the list of frozen atoms.
 
 Return values:
 AGBNP_OK - Cache of constant quantities updated.
 AGBNP_ERR - Error updating cache of constant quantities.
             Consult error messages on stderr.
 
 
 
 Verlet Neighbor List Utility Functions (libnblist)
 =================================================
 
 Introduction
 ------------
 
 The libnblist library provides a data structure to hold a Verlet
 neighbor list and a set of utility functions to manage the memory
 associated with the Verlet neighbor list.
 
 To use this facility the libnblist header file must be included
 
 #include "libnblist.h"
 
 this file is automatically included by agbnp.h.
 
 Structure
 ---------
 
 libnblist defines the NeighList data type, a C structure composed of
 the following members:
 
 - Core members:
 
    int natoms; number of atoms.
 
    int neighl_size; Holds the allocated size of the Verlet neighbor
                     list, this is the maximum number of neighbors that
                     can be stored in the list.

    int *nne; an array of size natoms that holds the number of
              neighbors for each atom.
 
    int *neighl1; the pointer to the memory space where the neighbor
                  list is actually stored.
 
    int **neighl; an array of pointers into neighl1 of size natoms
                  pointing to the start of the list of neighbors for
                  each atom. So for example, neighl[iat][j], 0 <= j <
                  nne[iat], are the atom indexes of the neighbors of
                  atom iat.
 
 - Atom indexes remapping members:

    int idx_remap; whether (yes > 0, no = 0) to do atom indexes remapping.
 
    int *int2ext; an optional array of size natoms to translate from
                  internal to external atom indexes.
 
    int *ext2int; an optional array of size max(external index)+1 to
                  translate from external to internal atom indexes.
 
    Note: libnblist can be used to store Verlet neighbor lists for
    calling programs using different atom indexing schemes - starting
    from zero, starting from one, random, etc . To do so libnblist
    makes a distinction between internal and external atom
    indexes. This is best illustrated with an example. Let us say that
    in the caller program 5 atoms are defined with atomic indexes 7,
    9, 4, 6, and 2. The neighbor lists for these atoms are stored in
    the Verlet neighbor list in this order so that neighl[0][j]
    identifies the neighbors of atom 7, neighl[1][j] identifies the
    neighbors of atom 9, etc. 0 is the internal index for atom 7, 1 is
    the internal atom index for atom 9, etc. In order to reference and
    dereference the neighbor list mapping arrays are constructed. The
    content of the int2ext array in this case would be:
 
    int2ext = { 7, 9, 4, 6, 2 };
 
    and the content of the ext2int array would be:
 
    ext2int = { -, -, 4, -, 2, -, 3, 0, -, 1 };
 
    So that for example int2ext[4] = 2 and ext2int[9] = 1. The
    neighbors of atom 9 are neighl[ext2int[9]][j], 0 <= j <
    nne[ext2int[9]]. The Verlet neighbor list stores the external atom
    indexes of the neighbors so that for example in iat =
    neighl[ext2int[9]][1], iat may be 6.  AGBNP uses the libnblist
    indexes remapping arrays to convert from internal AGBNP atom
    indexes (which start from 0 and end at natoms-1) to external atom
    indexes (that could be anything) and viceversa. It is the
    responsibility of the caller to properly set up the remapping
    arrays.
 
 - Other data:
 
    The libnblist data structure has place holders for other types of
    data (Periodic Boundary Conditions translation vectors and general
    additional data). These data structures are not described here as
    they are not relevant to AGBNP usage.
 
 Libnblist C API
 ---------------
 
 libnblist provides simple memory management functions for Verlet
 neighbor lists.
 
 The caller usually creates a libnblist structure by allocating memory
 for it:
 
 NeighList *nl;
 nl = (NeighList *)malloc(sizeof(NeighList));
 
 nl is then used as an argument for the utility functions that follow.
 
 > void nblist_reset_neighbor_list(NeighList *neigh_list);
 
 resets structure to pristine state. All pointers NULL, variables set to zero.
 
 > int nblist_reallocate_neighbor_list(NeighList *nl, int natoms, int nlsize);
 
 Reallocates memory for the Verlet neighbor list pointed by
 'nl'. 'natoms' is the new number of atoms and 'nlsize' is the new list
 size. If natoms > nl->natoms the list is resized in terms of the
 number of atoms. If nlsize > nl->neighl_size the size of the Verlet
 neighbor list is reallocated as well. In either case the previous
 contents of the neighbor list structure are preserved, including index
 remapping arrays and additional data.
 Return values:
      NBLIST_OK - Neighbor list memory reallocated.
      NBLIST_ERR - Error reallocating memory for neighbor
                   list. Neighbor list may be unusable at this point. 
                   Consult error messages on stderr.
 
 > int nblist_delete_neighbor_list(NeighList *nl);
 
 Frees memory associated with the Verlet neighbor list referenced by
 'nl'. The memory associated with the neighbor list structure itself
 is not freed. The calling program is responsible for calling free(nl)
 if needed.
 Return values:
      NBLIST_OK - Neighbor list memory freed.
      NBLIST_ERR - Error freeing memory of neighbor list.
                   Consult error messages on stderr.
 
 The following is a code section that illustrates how to create and
 fill a Verlet neighbor list (without atom indexes remapping, external
 atom indexes go from 0 to natoms-1):
 
 
   int nnl; /* neighbor list counter */
   int iat,jat,nlsize;
   double dx, dy, dz, d2;
   double cutoff;
   double nlsize_increment = 1.2;
   int natoms = 100;
   NeighList *nl;
   double x[100],y[100],z[100];
 
   /* set coordinates */
   ...
 
   /* allocates and reset neighbor list structure */
   nl = (NeighList *)malloc(sizeof(NeighList));
   nblist_reset_neighbor_list(nl);
   
   /* intial allocation of neighbor list */
   nlsize = natoms;
   if(nblist_reallocate_neighbor_list(nl,natoms,nlsize) != NBLIST_OK){
      fprintf(stderr,"Error allocating neighbor list.\n");
      exit(0);
   }  
 
   /* constructs neighbor list */
   nnl = 0;
   for(iat=0;iat<natoms;iat++){
       while(nnl + natoms >= nl->neighl_size){
 	nlsize = nlsize_increment*nl->neighl_size;
 	if(nblist_reallocate_neighbor_list(nl,natoms,nlsize) != NBLIST_OK){
 	  fprintf(stderr,"Error allocating neighbor list.\n");
 	  exit(0);
 	}
       }
       /* resets number of neighbors for atom iat */
       nl->nne[iat] = 0; 
       /*set pntr to beg. of neigh. list for atom iat*/
       nl->neighl[iat] = &(nl->neighl1[nnl]); 
       for(jat=iat+1;jat<natoms;jat++){
 	  /* calculate distance between iat and jat */
 	  dx = x[jat] - x[iat];
 	  dy = y[jat] - y[iat];
 	  dz = z[jat] - z[iat];
 	  d2 = dx*dx + dy*dy + dz*dz;
 	  if(d2<cutoff*cutoff){
 	    /* jat is a neighbor */
 	    nl->neighl1[nnl] = jat; /* place jat in neighbor list */
 	    nnl += 1;               /* updates neighbor list counter */
 	    nl->nne[iat] += 1;      /* updates counter of iat neighbors */
 	  }
       }
   }
 
 
 
 
 
