Can anyone explain in detail the method that gromacs uses to calculate
surface area? Which surface area is it calculating (e.g., connolley,
accessible)?

I am under the impression that the final result is given in nm^2. Is
this correct?


Carrying out the surface area calculation in gromacs using the OPLS all
atom force field I obtain a value of 4.43 nm^2  for methane (carbon =
opls_138, sigma = 0.35 and hydrogen = opls_140, sigma = 0.25) (CT-HC
ideal bond length = 0.109, e=284512). The command was 

g_sas -f conf.gro -s topol.tpr -oa -or -tv -probe 0.2 -ndots 200 -nopbc

and the result was taken from resarea.xvg


I have written a program to calculate the geometric accessible surface
area using the van der waals radius of the atoms and a probe of a given
size. The surface is defined from the centre of the probe. Using the
above values for sigma I obtain a value for methane of 1.15 nm^2 with a
probe diameter of 0.2 nm.

Obviously the two do not match and I am not sure why. 

The individual atom areas for my molecule are:
C 0.22
H 0.24
H 0.23
H 0.22
H 0.24

Attached are the files used in the gromacs calculation along with the
output and my own program. If any one has any suggestions please let me
know.
# This file was created Wed May 11 13:05:56 2011
# by the following command:
# g_sas -f conf.gro -s topol.tpr -oa -or -tv -probe 0.2 -ndots 200 -nopbc 
#
# g_sas is part of G R O M A C S:
#
# Gravel Rubs Often Many Awfully Cauterized Sores
#
@    title "Solvent Accessible Surface"
@    xaxis  label "Time (ps)"
@    yaxis  label "Area (nm\S2\N)"
@TYPE xy
@ view 0.15, 0.15, 0.75, 0.85
@ legend on
@ legend box on
@ legend loctype view
@ legend 0.78, 0.8
@ legend length 2
@ s0 legend "Hydrophobic"
@ s1 legend "Hydrophilic"
@ s2 legend "Total"
@ s3 legend "D Gsolv"
         0     4.43467           0     4.43467    -11.1328
# This file was created Wed May 11 13:05:56 2011
# by the following command:
# g_sas -f conf.gro -s topol.tpr -oa -or -tv -probe 0.2 -ndots 200 -nopbc 
#
# g_sas is part of G R O M A C S:
#
# Gravel Rubs Often Many Awfully Cauterized Sores
#
@    title "Area per atom"
@    xaxis  label "Atom #"
@    yaxis  label "Area (nm\S2\N)"
@TYPE xy
1 1.53938 0.000291734
2 0.723823 0
3 0.723823 0
4 0.723823 0
5 0.723823 0
methane
   5
    1DRG    CT    1   0.000   0.000   0.000
    1DRG    HC    2   0.000   0.000   1.089
    1DRG    HC    3   1.027   0.000  -0.363
    1DRG    HC    4  -0.513  -0.889  -0.363
    1DRG    HC    5  -0.513   0.889  -0.363
   7.90500   7.90500   7.90500
; VARIOUS PREPROCESSING OPTIONS
title                    = Yo
cpp                      = /usr/bin/cpp
include                  = 
define                   = 

; RUN CONTROL PARAMETERS
integrator               = steep
; Start time and timestep in ps
tinit                    = 0
dt                       = 0.002
nsteps                   = 1000
; For exact run continuation or redoing part of a run
init_step                = 0
; mode for center of mass motion removal
comm-mode                = Linear
; number of steps for center of mass motion removal
nstcomm                  = 1
; group(s) for center of mass motion removal
comm-grps                = 

; LANGEVIN DYNAMICS OPTIONS
; Temperature, friction coefficient (amu/ps) and random seed
;bd-temp                  = 300
bd-fric                  = 0
ld-seed                  = 1993

; ENERGY MINIMIZATION OPTIONS
; Force tolerance and initial step-size
emtol                    = 100
emstep                   = 0.01
; Max number of iterations in relax_shells
niter                    = 20
; Step size (1/ps^2) for minimization of flexible constraints
fcstep                   = 0
; Frequency of steepest descents steps when doing CG
nstcgsteep               = 1000
nbfgscorr                = 10

; OUTPUT CONTROL OPTIONS
; Output frequency for coords (x), velocities (v) and forces (f)
nstxout                  = 0
nstvout                  = 0
nstfout                  = 0
; Checkpointing helps you continue after crashes
nstcheckpoint            = 1000
; Output frequency for energies to log file and energy file
nstlog                   = 50
nstenergy                = 50
; Output frequency and precision for xtc file
nstxtcout                = 50
xtc-precision            = 1000
; This selects the subset of atoms for the xtc file. You can
; select multiple groups. By default all atoms will be written.
xtc-grps                 = 
; Selection of energy groups
energygrps               = 

; NEIGHBORSEARCHING PARAMETERS
; nblist update frequency
nstlist                  = 10
; ns algorithm (simple or grid)
ns_type                  = grid
; Periodic boundary conditions: xyz (default), no (vacuum)
; or full (infinite systems only)
pbc                      = xyz
; nblist cut-off        
rlist                    = 0.9
domain-decomposition     = no

; OPTIONS FOR ELECTROSTATICS AND VDW
; Method for doing electrostatics
coulombtype              = Cut-off
rcoulomb-switch          = 0
rcoulomb                 = 0.9
; Dielectric constant (DC) for cut-off or DC of reaction field
epsilon-r                = 1
; Method for doing Van der Waals
vdw-type                 = Cut-off
; cut-off lengths       
rvdw-switch              = 0
rvdw                     = 0.9
; Apply long range dispersion corrections for Energy and Pressure
DispCorr                 = EnerPres
; Extension of the potential lookup tables beyond the cut-off
table-extension          = 1
; Spacing for the PME/PPPM FFT grid
fourierspacing           = 0.12
; FFT grid size, when a value is 0 fourierspacing will be used
fourier_nx               = 0
fourier_ny               = 0
fourier_nz               = 0
; EWALD/PME/PPPM parameters
pme_order                = 4
ewald_rtol               = 1e-05
ewald_geometry           = 3d
epsilon_surface          = 0
optimize_fft             = no

; GENERALIZED BORN ELECTROSTATICS
; Algorithm for calculating Born radii
gb_algorithm             = Still
; Frequency of calculating the Born radii inside rlist
nstgbradii               = 1
; Cutoff for Born radii calculation; the contribution from atoms
; between rlist and rgbradii is updated every nstlist steps
rgbradii                 = 2
; Salt concentration in M for Generalized Born models
gb_saltconc              = 0

; IMPLICIT SOLVENT (for use with Generalized Born electrostatics)
implicit_solvent         = No

; OPTIONS FOR WEAK COUPLING ALGORITHMS
; Temperature coupling  
Tcoupl                   = berendsen
; Groups to couple separately
tc-grps                  = System
; Time constant (ps) and reference temperature (K)
tau_t                    = 0.1
ref_t                    = 300
; Pressure coupling     
Pcoupl                   = berendsen
Pcoupltype               = isotropic
; Time constant (ps), compressibility (1/bar) and reference P (bar)
tau_p                    = 1.0
compressibility          = 4.5e-5
ref_p                    = 1.0
; Random seed for Andersen thermostat
andersen_seed            = 815131

; SIMULATED ANNEALING  
; Type of annealing for each temperature group (no/single/periodic)
annealing                = no
; Number of time points to use for specifying annealing in each group
annealing_npoints        = 
; List of times at the annealing points for each group
annealing_time           = 
; Temp. at each annealing point, for each group.
annealing_temp           = 

; GENERATE VELOCITIES FOR STARTUP RUN
gen_vel                  = yes
gen_temp                 = 300
gen_seed                 = 1993

; OPTIONS FOR BONDS    
constraints              = all-bonds
; Type of constraint algorithm
constraint-algorithm     = Lincs
; Do not constrain the start configuration
unconstrained-start      = no
; Use successive overrelaxation to reduce the number of shake iterations
Shake-SOR                = no
; Relative tolerance of shake
shake-tol                = 1e-04
; Highest order in the expansion of the constraint coupling matrix
lincs-order              = 4
; Number of iterations in the final step of LINCS. 1 is fine for
; normal simulations, but use 2 to conserve energy in NVE runs.
; For energy minimization with constraints it should be 4 to 8.
lincs-iter               = 1
; Lincs will write a warning to the stderr if in one step a bond
; rotates over more degrees than
lincs-warnangle          = 30
; Convert harmonic bonds to morse potentials
morse                    = no

; ENERGY GROUP EXCLUSIONS
; Pairs of energy groups for which all non-bonded interactions are excluded
energygrp_excl           = 


[moleculetype]
; Name nrexcl
DRG      0

[ atoms ]
;   nr      type  resnr resid  atom  cgnr   charge     mass
     1       opls_138       1  DRG     CT     1    0.000  12.0000  
     2       opls_140       1  DRG     HC     1    0.000  01.0000  
     3       opls_140       1  DRG     HC     1    0.000  01.0000   
     4       opls_140       1  DRG     HC     1    0.000  01.0000   
     5       opls_140       1  DRG     HC     1    0.000  01.0000   

[ bonds ]
; i    j  func       b0          kb
  1    2     1    0.10800   284512.0
  1    3     1    0.10800   284512.0
  1    4     1    0.10800   284512.0
  1    5     1    0.10800   284512.0

# This file was created Wed May 11 13:05:56 2011
# by the following command:
# g_sas -f conf.gro -s topol.tpr -oa -or -tv -probe 0.2 -ndots 200 -nopbc 
#
# g_sas is part of G R O M A C S:
#
# Gravel Rubs Often Many Awfully Cauterized Sores
#
@    title "Area per residue"
@    xaxis  label "Residue"
@    yaxis  label "Area (nm\S2\N)"
@TYPE xy
         1     4.43467          0
#include "ffoplsaa.itp"
#include "methane.itp"
;#include "methnb.itp"



[System]
DRG

[Molecules]
DRG 1
program surface

!-------------------------------------------------------
!Adapted from code by Tina Duren (18/05/10 - program accessible_surface_area)
!to calculate the accessible surface area using Monte Carlo Integeration
!http://www.see.ed.ac.uk/~tduren/research/surface_area/ortho/program/accessible_surface_area.f90
!----------------------------------------------------

implicit none

integer, parameter:: maxatoms = 6000	!max number of atoms to be tested
integer, parameter:: nmonomers = 50	! number of monomers in the chain 
integer, parameter:: npoints = 200	!number of points to try per atom
real, parameter:: probe	=2.0		!probe diameter in angstroms

!read in the xyz file
character(30) :: myfile(0:60)
integer:: natoms	
integer::n
character(20):: title
character(2):: id(maxatoms)
real:: xatom(maxatoms), yatom(maxatoms), zatom(maxatoms)

!calculate the surface
real:: ratom(maxatoms)			!the vdw sigma of the atom

Integer:: iseed				!random number
Real:: ran0				!for generating random number from the function
real:: phi, costheta, theta		!angles for generating random points
real, parameter:: pi = 3.14159

real:: xpoint, ypoint, zpoint		!co-ords of a random point on the surface
real:: d				!distance 
integer:: naccept 			!counts the number of points accepted
real :: totalpoints			!total number of points tried for the mol.

real:: atom_area, sum_area, mol_area	!surface areas

integer:: i, j, k, ifile=0, pim		!for do loops etc
integer :: flag             		!for counting overlaps
real:: wt(maxatoms)			!atomic mass
real:: mass				!mass of the chain


!------------------------------------------------------------------------


! Seed for random number generator, change if you want to start from a 
! different random number.

iseed=-5213150

! to properly initialise the random number generator, a negative seed 
! is needeed

If (iseed > 0) iseed = -1 * iseed

!--------------------------------------------------------
!read in the xyz file

myfile(0)='PIM-PI-1/pim_1000.xyz'
myfile(1)='PIM-PI-1/pim_1001.xyz'
myfile(2)='PIM-PI-1/pim_1002.xyz'
myfile(3)='PIM-PI-1/pim_1003.xyz'
myfile(4)='PIM-PI-1/pim_1004.xyz'
myfile(5)='PIM-PI-1/pim_1005.xyz'
myfile(6)='PIM-PI-1/pim_1006.xyz'
myfile(7)='PIM-PI-1/pim_1007.xyz'
myfile(8)='PIM-PI-1/pim_1008.xyz'
myfile(9)='PIM-PI-1/pim_1009.xyz'
myfile(10)='PIM-PI-2/pim_1001.xyz'
myfile(11)='PIM-PI-2/pim_1002.xyz'
myfile(12)='PIM-PI-2/pim_1003.xyz'
myfile(13)='PIM-PI-2/pim_1004.xyz'
myfile(14)='PIM-PI-2/pim_1005.xyz'
myfile(15)='PIM-PI-2/pim_1006.xyz'
myfile(16)='PIM-PI-2/pim_1007.xyz'
myfile(17)='PIM-PI-2/pim_1008.xyz'
myfile(18)='PIM-PI-2/pim_1009.xyz'
myfile(19)='PIM-PI-2/pim_1010.xyz'
myfile(20)='PIM-PI-3/pim_1001.xyz'
myfile(21)='PIM-PI-3/pim_1002.xyz'
myfile(22)='PIM-PI-3/pim_1003.xyz'
myfile(23)='PIM-PI-3/pim_1004.xyz'
myfile(24)='PIM-PI-3/pim_1005.xyz'
myfile(25)='PIM-PI-3/pim_1006.xyz'
myfile(26)='PIM-PI-3/pim_1007.xyz'
myfile(27)='PIM-PI-3/pim_1008.xyz'
myfile(28)='PIM-PI-3/pim_1009.xyz'
myfile(29)='PIM-PI-3/pim_1010.xyz'
myfile(30)='PIM-PI-4/pim_1001.xyz'
myfile(31)='PIM-PI-4/pim_1002.xyz'
myfile(32)='PIM-PI-4/pim_1003.xyz'
myfile(33)='PIM-PI-4/pim_1004.xyz'
myfile(34)='PIM-PI-4/pim_1005.xyz'
myfile(35)='PIM-PI-4/pim_1006.xyz'
myfile(36)='PIM-PI-4/pim_1007.xyz'
myfile(37)='PIM-PI-4/pim_1008.xyz'
myfile(38)='PIM-PI-4/pim_1009.xyz'
myfile(39)='PIM-PI-4/pim_1010.xyz'
myfile(40)='PIM-PI-7/pim_1000.xyz'
myfile(41)='PIM-PI-7/pim_1001.xyz'
myfile(42)='PIM-PI-7/pim_1002.xyz'
myfile(43)='PIM-PI-7/pim_1003.xyz'
myfile(44)='PIM-PI-7/pim_1004.xyz'
myfile(45)='PIM-PI-7/pim_1005.xyz'
myfile(46)='PIM-PI-7/pim_1006.xyz'
myfile(47)='PIM-PI-7/pim_1007.xyz'
myfile(48)='PIM-PI-7/pim_1008.xyz'
myfile(49)='PIM-PI-7/pim_1009.xyz'
myfile(50)='PIM-PI-8/pim_1000.xyz'
myfile(51)='PIM-PI-8/pim_1001.xyz'
myfile(52)='PIM-PI-8/pim_1002.xyz'
myfile(53)='PIM-PI-8/pim_1003.xyz'
myfile(54)='PIM-PI-8/pim_1004.xyz'
myfile(55)='PIM-PI-8/pim_1005.xyz'
myfile(56)='PIM-PI-8/pim_1006.xyz'
myfile(57)='PIM-PI-8/pim_1007.xyz'
myfile(58)='PIM-PI-8/pim_1008.xyz'
myfile(59)='PIM-PI-8/pim_1009.xyz'


open(40,file='output.txt')
write(40,*) 'nmonomers:', nmonomers, '    npoints:', npoints, '    probe diameter:', probe

do pim=0,1
!pim = 5
!  if (pim==0) then
!	 write (40,*)'PIM-PI-1'
!	natoms =93*nmonomers
!  else if (pim==1) then
!	write (40,*)'PIM-PI-2'
!	natoms =106*nmonomers
!  else if (pim==2) then
!	write (40,*)'PIM-PI-3'
!	natoms =100*nmonomers
!  else if (pim==3) then
!	 write (40,*)'PIM-PI-4'
!	natoms =97*nmonomers
!  else if (pim==4) then
!	write (40,*)'PIM-PI-7'
!	natoms =87*nmonomers
!  else ! (pim==5) then
!	write (40,*)'PIM-PI-8'
!	natoms =109*nmonomers
!  end if

if (pim==0) then
  write(40,*)'Argon'
  natoms = 1
  open(10,file='Argon.xyz')
  print*,'opened argon'
else if(pim==1)then
  write(40,*)'Methane'
  natoms = 5
  open(10,file='methane.xyz')
    print*,'opened methane'
end if

!do ifile=0,9
	!ifile = 0
!	open(10,file=myfile(pim*10+ifile))
read(10,*) n			!1st line in the file is the no. of atoms
	!if(n > maxatoms) then
	!  Print*, 'Too many atoms found in xyz file'
	!  stop
	!end if
    
read (10,*) title			!2nd line of the file is the molecule name
print*,title
if(pim==0)then
   read(10,*) id(1), xatom(1), yatom(1), zatom(1)
else
   do i=1,natoms					!reads in the co-ordinates
	read(10,*) id(i), xatom(i), yatom(i), zatom(i)
	!print*, i, id(i), xatom(i), yatom(i), zatom(i)
   end do
end if

close (10)
  
!----------------calculate the surface-----------------------------------

!assign Van der Waals radii to the atoms
mass = 0.0

do i=1,natoms
  if (id(i)== 'C') then
    ratom(i)= 3.5
    wt(i)=12.0
  else if(id(i) == 'Ar') then
    ratom(i)= 3.5	!OPLS_135 (CT)
    wt(i) = 40.0
  else !id == H
    ratom(i) = 2.5
    wt(i) = 1.0
  end if
  print*, id(i), ratom(i)
  mass = mass + wt(i)
end do
print*,mass


!place a random point on an atom's surface and 
!find the distance to all other atoms
!count the no. of points that do not lie wthin the radius of another atom

sum_area = 0.0

do i=1,natoms		!loop over all atoms in chain
  naccept = 0

  do j=1,npoints	!generate random points on surface of atom
!	print*, 'calculating point ', j, ' of atom ', i, NACCEPT
	! Generate random vector of length 1
	! First generate phi 0:+2pi

     	 phi = ran0(iseed)*2.0*pi
	
	! Generate cosTheta -1:1 to allow for even distribution of random vectors
	! on the unit sphere.See http://mathworld.wolfram.com/SpherePointPicking.html
	! for further explanations
	
     	costheta = 1 - ran0(iseed) * 2.0
      	theta = Acos(costheta)
     	xpoint = sin(theta)*cos(phi)
       	ypoint = sin(theta)*sin(phi)
       	zpoint = costheta

	! Make this vector of length (sigma+probe_diameter)/2.0 

     	 xpoint = xpoint*(ratom(i)+probe)/2.0
     	 ypoint = ypoint*(ratom(i)+probe)/2.0
     	 zpoint = zpoint*(ratom(i)+probe)/2.0

	! Translate the center of coordinate to the particle i center and apply PBC

     	 xpoint = xpoint + xatom(i)
     	 ypoint = ypoint + yatom(i)
     	 zpoint = zpoint + zatom(i)
      flag=0
      do k=1, natoms		!find distance to surface of all other atoms
        if(flag==0) then
        if(k == i) cycle
	       d = sqrt((xpoint-xatom(k))**2+(ypoint-yatom(k))**2+(zpoint-zatom(k))**2)
	       if ( d < (ratom(k)+probe)/2. ) flag=1
        endif
      end do
	if(flag==0) naccept=naccept+1
  end do	!j loop -npoints
!PRINT*, 'FINISH LOOP'

!use the fraction of points accepted to find the surface area of the chain
  totalpoints = REAL(naccept)/REAL(npoints)			!fraction of points accepted for atom i

 ! PRINT*, 'TOTALPOINTS',TOTALPOINTS,NACCEPT,NPOINTS
 ! print*, 'i',i, ratom(i)

  atom_area = 4.0*pi*((ratom(i)+probe)/2.)**2 * totalpoints	!surface area of atom
print*,naccept, atom_area
  sum_area = sum_area + atom_area		!surface area of all the atoms
end do   	!i loop -natoms

write(40,*) ifile+1,'::  ', sum_area, ' Asq.    ', sum_area*(6.023e23/1.e20)/mass , ' msq/gram'

print*, '---------- Chain', ifile+1, 'of PIM-PI-',pim+1,'----------'
print*, 'The surface area of the chain is:', sum_area, ' Angstroms squared'
print*, 'This is:', sum_area*(6.023e23/1.e20)/mass , ' metres squared per gram'

!end do !myfile
end do !pim

close(40)

end program surface


!----------------FUNCTIONS-------------------------------------

FUNCTION ran0(idum)
!
! Random number generator from W.H. Press et al, Numerical Recipes in
! FORTRAN, Cambridge University Press, 1992
!
 INTEGER idum,IA,IM,IQ,IR,NTAB,NDIV
 REAL  ran0,AM,EPS,RNMX
 PARAMETER (IA=16807,IM=2147483647,AM=1./IM,IQ=127773,IR=2836, &
   NTAB=32,NDIV=1+(IM-1)/NTAB,EPS=1.2e-7,RNMX=1.-EPS)
 INTEGER j,k,iv(NTAB),iy
 
 SAVE iv,iy
 DATA iv /NTAB*0/, iy /0/
 if (idum.le.0.or.iy.eq.0) then
     idum=max(-idum,1)
     do 11 j=NTAB+8,1,-1
        k=idum/IQ
        idum=IA*(idum-k*IQ)-IR*k
          if (idum.lt.0) idum=idum+IM
          if (j.le.NTAB) iv(j)=idum
11   continue
     iy=iv(1)
endif
k=idum/IQ
idum=IA*(idum-k*IQ)-IR*k
if (idum.lt.0) idum=idum+IM
j=1+iy/NDIV
iy=iv(j)
iv(j)=idum
ran0=min(AM*iy,RNMX)
return

END FUNCTION ran0

-- 
gmx-users mailing list    gmx-users@gromacs.org
http://lists.gromacs.org/mailman/listinfo/gmx-users
Please search the archive at 
http://www.gromacs.org/Support/Mailing_Lists/Search before posting!
Please don't post (un)subscribe requests to the list. Use the 
www interface or send it to gmx-users-requ...@gromacs.org.
Can't post? Read http://www.gromacs.org/Support/Mailing_Lists

Reply via email to