Dear Christoph,

   Did you also change original bond and angle tables to gromacs format?

If so, how do I determine the min and max values for bond and angle tables?



Thank you.

Best regards,
Changwoon Jang

On Tue, Aug 30, 2016 at 2:04 PM, Christoph Junghans <[email protected]>
wrote:

> Hi,
>
> thanks for sending the potentials.
> Here is what I did:
> $ ls table_d*.txt
> table_d1.txt  table_d2.txt  table_d3.txt  table_d4.txt
> # These are the potentials out of csg_boltzmann
> $ cat table_d.xml
> <cg>
>   <bonded>
>     <!-- name of the interaction -->
>     <name>dihedral</name>
>     <min>-3.141</min>
>     <max>3.141</max>
>     <step>0.001</step>
>   </bonded>
>   <inverse>
>     <gromacs>
>       <pot_max>1e8</pot_max>
>       <table_bins>0.001</table_bins>
>     </gromacs>
>   </inverse>
> </cg>
> $ for i in table_d*.txt; do csg_call --options table_d.xml --ia-type
> dihedral --ia-name dihedral --sloppy-tables convert_potential gromacs
> $i ${i%.txt}.xvg; done
> This basically does:
> csg_call --options table_d.xml --ia-type dihedral --ia-name dihedral
> --sloppy-tables convert_potential gromacs table_d1.txt table_d1.xvg
> csg_call --options table_d.xml --ia-type dihedral --ia-name dihedral
> --sloppy-tables convert_potential gromacs table_d2.txt table_d2.xvg
> csg_call --options table_d.xml --ia-type dihedral --ia-name dihedral
> --sloppy-tables convert_potential gromacs table_d3.txt table_d3.xvg
> csg_call --options table_d.xml --ia-type dihedral --ia-name dihedral
> --sloppy-tables convert_potential gromacs table_d4.txt table_d4.xvg
> Please note you need the "--sloppy-tables" option as you want to
> ignore the third column of  table_d1.txt
> Alternatively you could replace the 3rd column by "i":
> $ awk '/^[^#]/{print $1,$2,"i"}' table_d1.txt > fixed_table_d1.txt
> Then a I ran gromacs
> $ gmx mdrun -tableb table_[abd]*.xvg -v
> Please mind that this is gromacs 2016. which is why you need the
> "-tableb" option. (Gromacs 5.0 was made obsolete a while back by its
> developers.)
>
> Seems to work for at least 10k steps.
>
> Cheers,
>
> Christoph
>
>
>
> 2016-08-29 20:52 GMT-06:00 Christoph Junghans <[email protected]>:
> > 2016-08-29 15:19 GMT-06:00 Chang Woon Jang <[email protected]>:
> >> Dear Christoph,
> >>
> >>    Can I run ibi without dihedrals? Actually, I do not need dihedral
> >> potentials even though CG molecules have dihedrals. Or, any suggestions
> for
> >> debugging the dihedrals. Actually, the warning signs disappeared but I
> have
> >> still the same Fatal error.
> > Please send me the table and the xml off-list and I will have a look.
> >
> > And of course you can run IBI on a subset of the interactions, see
> > <https://github.com/votca/csg-tutorials/tree/master/hexane/ibi>
> > In this example only the non-bonded interactions are iterated.
> >
> > Christoph
> >>
> >> Thank you.
> >>
> >> Best regards,
> >> Changwoon Jang
> >>
> >> On Mon, Aug 29, 2016 at 4:30 PM, Christoph Junghans <[email protected]
> >
> >> wrote:
> >>>
> >>> 2016-08-29 13:32 GMT-06:00 Chang Woon Jang <[email protected]>:
> >>> > Dear Votca Users,
> >>> >
> >>> >    I have four dihedral distributions directly from csg_boltzmann.
> >>> > However,
> >>> > only one dihedral distribution gives waning sign during IBI run as
> >>> > follows.
> >>> >
> >>> > WARNING: For the 359 non-zero entries for table 0 in table_d1.xvg the
> >>> > forces
> >>> > deviate on average 196% from minus the numerical derivative of the
> >>> > potential
> >>> That means something went wrong in the table conversion!
> >>> >
> >>> >
> >>> >
> >>> > Then, the IBI fails with Fatal error: 1 of the 1720 bonded
> interactions
> >>> > could not be calculated because some atoms involved moved further
> apart
> >>> > than
> >>> > the multi-body cut-off distance (2.02391 nm) or the two-body cut-off
> >>> > distance (2.02391 nm), see option -rdd, for pairs and tabulated bonds
> >>> > also
> >>> > see option -ddcheck
> >>> That is most likely due to the above warning.
> >>> >
> >>> >
> >>> > I am not sure why only one dihedral gives this warning and prevent
> IBI
> >>> > run.
> >>> > Do you have any comments on that?
> >>> IBI for dihedral interactions is a bit tricky and not very well
> >>> testing in VOTCA, so you will need to do a bit a debugging on your
> >>> side.
> >>>
> >>> Christoph
> >>>
> >>> >
> >>> > I have converted dihedral distribution radian to angle (only first
> >>> > column)
> >>> > due to the ibi requirement just using unit convert     x*180/3.14 (x
> is
> >>> > radian).
> >>> >
> >>> > Thank you.
> >>> >
> >>> > Best regards,
> >>> > Changwoon Jang,
> >>> >
> >>> > --
> >>> > You received this message because you are subscribed to the Google
> >>> > Groups
> >>> > "votca" group.
> >>> > To unsubscribe from this group and stop receiving emails from it,
> send
> >>> > an
> >>> > email to [email protected].
> >>> > To post to this group, send email to [email protected].
> >>> > Visit this group at https://groups.google.com/group/votca.
> >>> > For more options, visit https://groups.google.com/d/optout.
> >>>
> >>>
> >>>
> >>> --
> >>> Christoph Junghans
> >>> Web: http://www.compphys.de
> >>>
> >>> --
> >>> You received this message because you are subscribed to the Google
> Groups
> >>> "votca" group.
> >>> To unsubscribe from this group and stop receiving emails from it, send
> an
> >>> email to [email protected].
> >>> To post to this group, send email to [email protected].
> >>> Visit this group at https://groups.google.com/group/votca.
> >>> For more options, visit https://groups.google.com/d/optout.
> >>
> >>
> >>
> >>
> >> --
> >> Best regards,
> >> Changwoon Jang,
> >>
> >> Postdoctoral Research Fellow
> >> Department of Chemical & Biological Engineering, Drexel University
> >> 3141 Chestnut Street, Philadelphia, PA 19104
> >>
> >> Voice: (662) 617-2267
> >> E-mail: [email protected]
> >>
> >> --
> >> You received this message because you are subscribed to the Google
> Groups
> >> "votca" group.
> >> To unsubscribe from this group and stop receiving emails from it, send
> an
> >> email to [email protected].
> >> To post to this group, send email to [email protected].
> >> Visit this group at https://groups.google.com/group/votca.
> >> For more options, visit https://groups.google.com/d/optout.
> >
> >
> >
> > --
> > Christoph Junghans
> > Web: http://www.compphys.de
>
>
>
> --
> Christoph Junghans
> Web: http://www.compphys.de
>
> --
> You received this message because you are subscribed to the Google Groups
> "votca" group.
> To unsubscribe from this group and stop receiving emails from it, send an
> email to [email protected].
> To post to this group, send email to [email protected].
> Visit this group at https://groups.google.com/group/votca.
> For more options, visit https://groups.google.com/d/optout.
>



-- 
Best regards,
Changwoon Jang,

Postdoctoral Research Fellow
Department of Chemical & Biological Engineering, Drexel University
3141 Chestnut Street, Philadelphia, PA 19104

Voice: (662) 617-2267
E-mail: [email protected]

-- 
You received this message because you are subscribed to the Google Groups 
"votca" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to [email protected].
To post to this group, send email to [email protected].
Visit this group at https://groups.google.com/group/votca.
For more options, visit https://groups.google.com/d/optout.

Reply via email to