Re: [AMBER] parmed question regarding modifying Lennard-Jones parameters

From: InSuk Joung <i.joung.gmail.com>
Date: Fri, 8 Jan 2016 23:00:01 +0900

Unfortunately, changeLJSingleType() does not work.

In omm_nonbonded_force() in structure.py, it seems that sigma and epsilon
are actually used rathter than values in low-level parm arrays.

So, I tried

mineps = 0.01
minsig = 0.1
for a in prmtop.atoms:
    if a.epsilon_14 == 0.0:
        a.epsilon_14 = mineps
        print 'EPS14', a.epsilon_14
    if a.sigma_14 == 0.0:
        a.sigma_14 = minsig
        print 'SIG14', a.sigma_14
    if a.epsilon == 0.0:
        a.epsilon = mineps
        print 'EPS', a.epsilon
    if a.sigma == 0.0:
        a.sigma = minsig
        print 'SIG', a.sigma

For some reason, a new value cannot be assigned to epsilon_14. The print
line tells that a.epsilon_14 is not changed. (a.epsilon and a.sigma are
replaced) I thought it is the source of the problem.

Then I found this:
In line 678, parmed/topologyobjects.py
self._rmin14 = value
should be replaced with
self._epsilon14 = value

Now, it works fine!

On Fri, Jan 8, 2016 at 8:50 PM, Jason Swails <jason.swails.gmail.com> wrote:

> On Fri, Jan 8, 2016 at 3:45 AM, InSuk Joung <i.joung.gmail.com> wrote:
>
> > I am using prmed to load topology files and running simulations with
> > openMM. However, I want to change Lennard-Jones parameters for some atoms
> > (especially zero-radius hydrogens). In case of wild conformational
> search,
> > these zero-radius hydrogens cause nasty problems due to electrostatic
> > interactions.
> >
> > Before create openMM's system, I guess I need to modify them by giving
> tiny
> > values to both epsilons and sigma.
> > This is the part of my python script.
> >
> > prmtop = pmd.load_file('model.prmtop')
> > # something for changing lennard-jones parameters
> > system = prmtop.createSystem()
> >
> > So far, I tried to modify atom.epsilon and atom.sigma in prmtop.atoms,
> but
> > it does not seem to work.
> >
>
> ​Possibly because the information is not being flushed back to the raw data
> arrays... That's a complication with the AmberParm class -- there are two
> places where the data is held and it's difficult to keep them in-sync.
>
> I know that I can use interactive parmed and modify the topology file. I
> > just want to know if there is a possible way in script-level without
> > modifying the topology file.
> >
>
> Anything in the interactive interpreter can be used at the script level.
> So what you want to use is "changeLJSingleType", and you can use that
> inside your script as follows:
>
> import parmed as pmd
>
> prmtop = pmd.load_file('model.prmtop')
> pmd.tools.changeLJSingleType(prmtop, '.%HO', rmin=0.8,
> epsilon=0.05).execute()
> system = prmtop.createSystem() ...
>
> If that doesn't work, let me know.
>
> I hope everything's going well!
> Jason
>
> --
> Jason M. Swails
> BioMaPS,
> Rutgers University
> Postdoctoral Researcher
> _______________________________________________
> AMBER mailing list
> AMBER.ambermd.org
> http://lists.ambermd.org/mailman/listinfo/amber
>



-- 
Best,
InSuk Joung
_______________________________________________
AMBER mailing list
AMBER.ambermd.org
http://lists.ambermd.org/mailman/listinfo/amber
Received on Fri Jan 08 2016 - 06:30:03 PST
Custom Search