[gmx-developers] Calling do_steep from within do_md integrator
Mark Abraham
mark.j.abraham at gmail.com
Sun Feb 8 23:54:50 CET 2015
On Sun, Feb 8, 2015 at 4:57 PM, Andrew AbiMansour <anabiman at umail.iu.edu>
wrote:
> Hi David,
>
> Thanks for the reply. I am implementing an algorithm that modifies the
> coordinates during an MD run, but once the system is perturbed, it must be
> "fixed" before continuing MD (otherwise the system experiences high forces
> and GROMACS crashes). So calling energy minimization (SD) is one possible
> approach here. What does the function relax_shells_fc() do exactly?
>
> Another questions: does anyone know the difference between "nr" and
> "homenr" in the data struct "t_mdatoms" ?
>
One is probably the total number of atoms, and the other is the total
number "home" to this domain. We are slowly adding more more documentation,
but this area is still a wasteland until we finally get rid of the group
cut-off scheme and the maintenance burden drops from impossible to
"maniacal laughter..." Meanwhile, things like
git grep "md.*homenr"
are good allies for finding the code where such fields are set.
Mark
> Thanks!
>
>
> *Andrew Abi-Mansour*
>
> On Sun, Feb 8, 2015 at 5:12 AM, David van der Spoel <spoel at xray.bmc.uu.se>
> wrote:
>
>> On 2015-02-07 17:44, Andrew AbiMansour wrote:
>>
>>> Hi Mark,
>>>
>>> Thank you for the reply. I decided to just test my own simple SD
>>> integrator. After calling do_force() I update the positions. Example:
>>>
>>> int i = 0, dim = 3;
>>> real sd_step_size = 0.01;
>>>
>>> for(i = 0; i < mdatoms->nr; i++)
>>> for(d = 0; d < dim; d++) {
>>> state->x[i][d] += sd_step_size * f[i][d];
>>>
>>>
>>> The code *works* , but the output doesnt make much sense. When I print
>>> the force (f[i][d]), I get a lot of zeros. Am I doing this right? I'm
>>> afraid I might not be accessing the forces (stored locally on each
>>> processor) correctly. Any suggestions?
>>>
>>
>> It would be helpful to know what you are trying to achieve. Maybe the
>> functionality can be found in the routine relax_shells_fc already?
>>
>>
>>> Thanks!
>>>
>>>
>>> *Andrew Abi-Mansour*
>>>
>>> On Sat, Feb 7, 2015 at 6:54 AM, Mark Abraham <mark.j.abraham at gmail.com
>>> <mailto:mark.j.abraham at gmail.com>> wrote:
>>>
>>> Hi,
>>>
>>> Probably you are over- or under-managing periodic boundaries and/or
>>> domain decomposition. I would think you should refactor out of
>>> do_steep exactly the part you need, since it's not designed to be a
>>> callable utility routine. For example, you probably don't want it to
>>> ever call em_dd_partition_system()...
>>>
>>> Mark
>>>
>>> On Sat, Feb 7, 2015 at 7:15 AM, Andrew AbiMansour
>>> <anabiman at umail.iu.edu <mailto:anabiman at umail.iu.edu>> wrote:
>>>
>>> Hi,
>>>
>>> I would like to call the integrator function do_steep() from
>>> within do_md(). I modified some constraints parameters in
>>> t_inputrec before passing all the arguments supplied to do_md()
>>> by the user to do_steep(). This simple approach, however, does
>>> not work and GROMACS crashes complaining about a charge group
>>> moving too far. The coordinates seem to change by several nano
>>> meters over a single SD step. What is this happening? and Do you
>>> have any suggestions on how to implement this right?
>>>
>>> Thanks
>>>
>>> *Andrew Abi-Mansour*
>>>
>>> --
>>> Gromacs Developers mailing list
>>>
>>> * Please search the archive at
>>> http://www.gromacs.org/Support/Mailing_Lists/GMX-developers_List
>>> before posting!
>>>
>>> * Can't post? Read http://www.gromacs.org/Support/Mailing_Lists
>>>
>>> * For (un)subscribe requests visit
>>> https://maillist.sys.kth.se/mailman/listinfo/gromacs.org_
>>> gmx-developers
>>> or send a mail to gmx-developers-request at gromacs.org
>>> <mailto:gmx-developers-request at gromacs.org>.
>>>
>>>
>>>
>>> --
>>> Gromacs Developers mailing list
>>>
>>> * Please search the archive at
>>> http://www.gromacs.org/Support/Mailing_Lists/GMX-developers_List
>>> before posting!
>>>
>>> * Can't post? Read http://www.gromacs.org/Support/Mailing_Lists
>>>
>>> * For (un)subscribe requests visit
>>> https://maillist.sys.kth.se/mailman/listinfo/gromacs.org_
>>> gmx-developers
>>> or send a mail to gmx-developers-request at gromacs.org
>>> <mailto:gmx-developers-request at gromacs.org>.
>>>
>>>
>>>
>>>
>>>
>>
>> --
>> David van der Spoel, Ph.D., Professor of Biology
>> Dept. of Cell & Molec. Biol., Uppsala University.
>> Box 596, 75124 Uppsala, Sweden. Phone: +46184714205.
>> spoel at xray.bmc.uu.se http://folding.bmc.uu.se
>>
>> --
>> Gromacs Developers mailing list
>>
>> * Please search the archive at http://www.gromacs.org/
>> Support/Mailing_Lists/GMX-developers_List before posting!
>>
>> * Can't post? Read http://www.gromacs.org/Support/Mailing_Lists
>>
>> * For (un)subscribe requests visit
>> https://maillist.sys.kth.se/mailman/listinfo/gromacs.org_gmx-developers
>> or send a mail to gmx-developers-request at gromacs.org.
>>
>
>
> --
> Gromacs Developers mailing list
>
> * Please search the archive at
> http://www.gromacs.org/Support/Mailing_Lists/GMX-developers_List before
> posting!
>
> * Can't post? Read http://www.gromacs.org/Support/Mailing_Lists
>
> * For (un)subscribe requests visit
> https://maillist.sys.kth.se/mailman/listinfo/gromacs.org_gmx-developers
> or send a mail to gmx-developers-request at gromacs.org.
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://maillist.sys.kth.se/pipermail/gromacs.org_gmx-developers/attachments/20150208/324b622a/attachment-0001.html>
More information about the gromacs.org_gmx-developers
mailing list