[gmx-users] Software inconsistency error: update_coords called for velocity without VV integrator

Berk Hess gmx3 at hotmail.com
Tue Aug 24 11:33:05 CEST 2010


Hi,

This is something we forgot to update.
I fixed it for the next beta release.
If you can read/use diff's, you can fix it yourself with the diff below.

Berk

--- a/src/tools/gmx_membed.c
+++ b/src/tools/gmx_membed.c
@@ -3121,9 +3121,12 @@ double do_md_membed(FILE *fplog,t_commrec *cr,int nfile,const t_filenm fnm[],
                 update_pcouple(fplog,step,ir,state,pcoupl_mu,M,wcycle,
                                 upd,bInitStep);
 
-                /* velocity (for VV) */
-                update_coords(fplog,step,ir,mdatoms,state,f,fr->bTwinRange && bNStList,fr->f_twin,fcd,
-                              ekind,M,wcycle,upd,FALSE,etrtVELOCITY2,cr,nrnb,constr,&top->idef);
+               if (bVV)
+               {
+                   /* velocity half-step update */
+                   update_coords(fplog,step,ir,mdatoms,state,f,fr->bTwinRange && bNStList,fr->f_twin,fcd,
+                                 ekind,M,wcycle,upd,FALSE,etrtVELOCITY2,cr,nrnb,constr,&top->idef);
+               }
 
                 /* Above, initialize just copies ekinh into ekin,
                  * it doesn't copy position (for VV),


Date: Tue, 24 Aug 2010 16:41:35 +0800
From: zhongjin1000 at yahoo.com.cn
To: gmx-users at gromacs.org
Subject: [gmx-users] Software inconsistency error: update_coords called for	velocity without VV integrator






Hi ,
 When I use g_membed in GMX4.5.3 Beta 3 to insert a protein into a DMPC bilayer,  I came across such an error:
Program g_membed, VERSION 4.5-beta3
Source code file: update.c, line: 1595
Software inconsistency error:
update_coords called for velocity without VV integrator
For more information and tips for troubleshooting, please check the GROMACS

 
g_membed -f gmem.tpr -p 1b8l.top -xyinit 0.5 -xyend 1 -nxy 1000 -ndiff 11 -c OK.pdb
 
Reading file gmem.tpr, VERSION 4.5-beta3 (single precision)
Reading file gmem.tpr, VERSION 4.5-beta3 (single precision)
Select a group to embed in the membrane:
Group     0 (         System) has 16956 elements
Group     1 (        Protein) has  3708 elements
Group     2 (      Protein-H) has  2904 elements
Group     3 (        C-alpha) has   388 elements
Group     4 (       Backbone) has  1164 elements
Group     5 (      MainChain) has  1556 elements
Group     6 (   MainChain+Cb) has  1904 elements
Group     7 (    MainChain+H) has  1940 elements
Group     8 (      SideChain) has 
 1768 elements
Group     9 (    SideChain-H) has  1348 elements
Group    10 (    Prot-Masses) has  3708 elements
Group    11 (    non-Protein) has 13248 elements
Group    12 (          Other) has 13248 elements
Group    13 (           DMPC) has 13248 elements
Select a group: 1
Selected 1: 'Protein'
Select a group to embed Protein into (e.g. the membrane):
Group     0 (         System) has 16956 elements
Group     1 (        Protein) has  3708 elements
Group     2 (      Protein-H) has  2904 elements
Group     3 (        C-alpha) has   388 elements
Group     4 (       Backbone) has  1164 elements
Group     5 (      MainChain) has  1556 elements
Group     6 (   MainChain+Cb) has  1904 elements
Group     7 (    MainChain+H) has  1940 elements
Group     8 (     
 SideChain) has  1768 elements
Group     9 (    SideChain-H) has  1348 elements
Group    10 (    Prot-Masses) has  3708 elements
Group    11 (    non-Protein) has 13248 elements
Group    12 (          Other) has 13248 elements
Group    13 (           DMPC) has 13248 elements
Select a group: 13
Selected 13: 'DMPC'
The estimated area of the protein in the membrane is 15.930 nm^2
There are 130 lipids in the membrane part that overlaps the protein.
The area per lipid is 0.6809 nm^2.
Maximum number of lipids that will be removed is 46.
Will resize the protein by a factor of 0.500 in the xy plane and 1.000 in the z direction.
This resizing will be done with respect to the geometrical center of all protein atoms
that span the membrane region, i.e. z between 1.299 and 6.362
Embedding piece 0 with center of geometry: 4.933612 4.933665 3.830500
Will remove 0 Protein_chain_A molecules
Will remove 0 Protein_chain_B molecules
Will remove 0 Protein_chain_C molecules
Will remove 0 Protein_chain_D molecules
Will remove 37 DMPC molecules
Back Off! I just backed up 1b8l.top to ./#1b8l.top.2#
Back Off! I just backed up traj.trr to ./#traj.trr.1#
Back Off! I just backed up traj.xtc to ./#traj.xtc.1#
Back Off! I just backed up ener.edr to ./#ener.edr.1#
starting mdrun 'God Rules Over Mankind, Animals, Cosmos and Such'
1000 steps,      2.0 ps.
-------------------------------------------------------
Program g_membed, VERSION 4.5-beta3
Source code file: update.c, line: 1595
Software inconsistency error:
update_coords called for velocity without VV integrator
For more information and tips for troubleshooting, please check the GROMACS
website at http://www.gromacs.org/Documentation/Errors
-------------------------------------------------------
"Look at these, my work-strong arms" (P.J. Harvey)
 
How to solve such an problem? Thanks!
 
Zhongjin He







       
-- 
gmx-users mailing list    gmx-users at gromacs.org
http://lists.gromacs.org/mailman/listinfo/gmx-users
Please search the archive at http://www.gromacs.org/search before posting!
Please don't post (un)subscribe requests to the list. Use the 
www interface or send it to gmx-users-request at gromacs.org.
Can't post? Read http://www.gromacs.org/mailing_lists/users.php 		 	   		  
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://maillist.sys.kth.se/pipermail/gromacs.org_gmx-users/attachments/20100824/1d422e35/attachment.html>


More information about the gromacs.org_gmx-users mailing list