Dear Xiaohu, Thanks for bringing this up. The comment has been there for ever. Since I could not think of an application where one would not be using pbc at the time. However, your clusters prove me wrong.
In any case, as a work-around, you may want to use a bigger box, with long enough cut-offs, no PME. Hope this gets you started, Gerrit On 27 Jan 2011, at 21:02, Xiaohu Li wrote: > Dear GMX code writters, > Could anyone tell me why this comments in the code mdlib/qmmm.c > appear?(version 4.5.3) > at line ~ 714, in the beginning of subroutine update_QMMMrec > ====================================================== > /* updates the coordinates of both QM atoms and MM atoms and stores > * them in the QMMMrec. > * > * NOTE: is NOT yet working if there are no PBC. Also in ns.c, simple > * ns needs to be fixed! > */ > ====================================================== > First of all, I have a non-PBC job with electronic embedding and it > fails right in this routine. So I guess this is the reason, for non-PBC, it > crashes. But for oniom, it went through. Is there any simple fix to have the > non-PBC qmmm work? I'm also interested in doing non-PBC calculation(large > cluster). > > > Xiaohu > > On Thu, Jan 27, 2011 at 9:11 AM, Christoph Riplinger <c...@thch.uni-bonn.de> > wrote: > Hi Xiaohu, > > We are using gromacs/ORCA quite a lot. It works like any of the other > interfaced programs and whether you should use it depends on your needs to > the QM part of the QM/MM calculation. > > I am using gromacs 4.0.7, but I downloaded the 4.5.3 version to test it. > > I could not reproduce your problem with ONIOM. For the electrostatic > embedding bug, I agree. This didn't work in my test calculation, but my gdb > told me that the abort occured before the qm_orca.c was even called. I > suspect the bug is in the general qmmm.c code. I do not have access to the > other QM-programs and cannot test this, so if you have, could you please run > the same job with a different QM-program. > > Hope that helps > > Christoph > > > > > On 01/26/2011 06:48 AM, Xiaohu Li wrote: >> >> Hi, All, >> I'm trying to see if anybody has experience of using the interface of >> gromacs and ORCA(since it's free). I know that the following link gave >> information on how >> http://wwwuser.gwdg.de/~ggroenh/qmmm.html#code >> But.....But, the gromacs in the above link is quite old(3.2). I >> download the latest 4.5.3 and followed the instructions in the above link >> and I was trying to optimize an simple cluster(no pbc) where part of it are >> treated using QM. here is the example mdp file >> ========================================================================================================================= >> title = cpeptide >> integrator = steep ;integrator includes energy minimization >> algorithms >> dt = 0.002 ; ps ! >> nsteps = 10000 >> nstlist = 1 >> ns_type = simple >> rlist = 3.0 >> rcoulomb = 3.0 >> coulombtype = cut-off >> vdwtype = cut-off >> rvdw = 3.0 >> pbc = no >> periodic_molecules = no >> constraints = none >> energygrps = qm_part mm_part >> ; QM/MM calculation stuff >> QMMM = yes >> QMMM-grps = qm_part >> QMmethod = rhf >> QMbasis = 3-21G >> QMMMscheme = oniom >> QMcharge = 0 >> QMmult = 1 >> ; >> ; Energy minimizing stuff >> ; >> emtol = 60 ; minimization thresold (kj/mol.nm-1) 1 >> hartree/bohr= 49614.75241 kj/mol.nm-1 1 kj/mol.nm-1=2.01553e-5 hartree/bohr >> emstep = 0.01 ; minimization step in nm >> ========================================================================================================================= >> I set up the BASENAME and ORCA_PATH as told in the instruction. >> first of all, the normal electronic embedding just simply gave segmentation >> fault error right after the it prints information on number of steps of >> optimization. >> >> So I switch to ONIOM, this time, at least, orca is called and energy and >> gradient are both generated. However, when it comes to read the energy and >> gradient, it always crashes when tried to read gradient, this is at line 346 >> source code src/mdlib/qm_orca.c >> ============================================ >> sscanf(buf,"%lf\n", &QMgrad[k][XX]); >> ============================================ >> a segmentation fault error is printed. If I replace the &QMgrad[k][XX] by an >> temporary variable temp >> sscanf(buf,"%lf\n", &temp); >> temp gets the correct value and if I use, >> QMgrad[k][XX]=temp >> and tries to print QMgrad[k][XX], a bus error will be printed. >> I did some research online, seems that usually this implies an memory bug in >> the code which is the most difficult bug one can ever encounter. >> So has anyone successfully used gromacs and orca to do QMMM? >> Generally, would anyone recommend using gromacs to do QMMM? >> >> Cheers, >> Xiaohu > > > > > -- > ================================================= > Xiaohu Li > Post-doctoral Research Associate > Department of Chemistry > Northwestern University > 2145 Sheridan Road > Evanston IL 60208 > U.S.A > ================================================= >
-- 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