real do_optimize_poisson(FILE *log, bool bVerbose, t_inputrec *ir, int natoms, rvec x[], rvec f[], real charge[], rvec box, real phi[], t_commrec *cr, t_nrnb *nrnb, rvec f_ref[], real phi_ref[], rvec beta, bool bOld) { #define BMIN 1.6 #define DB 0.025 #define NB (1+(int)((2.1-BMIN)/DB)) static bool bFirst = TRUE; static bool bSecond = TRUE; static t_PSgrid *pot,*rho; static int maxnit; static real r1,rc; real rmsf[NB][NB][NB],rmsf_min,rrmsf; ivec minimum; const real tol = 1e-2; int i,m,bx,by,bz; char buf[128]; real ener; ener = 0.0; if (bFirst) { maxnit = ir->userint1; fprintf(log,"Will use Poisson Solver for long-range electrostatics\n"); fprintf(log,"Grid size is %d x %d x %d\n",ir->nkx,ir->nky,ir->nkz); if ((ir->nkx < 4) || (ir->nky < 4) || (ir->nkz < 4)) fatal_error(0,"Grid must be at least 4 points in all directions"); pot = mk_PSgrid(ir->nkx,ir->nky,ir->nkz); rho = mk_PSgrid(ir->nkx,ir->nky,ir->nkz); r1 = ir->rcoulomb_switch; rc = ir->rcoulomb; for(m=0; (m<DIM); m++) beta[m] = 4.0/3.0; bFirst = FALSE; } /* Make the grid empty */ clear_PSgrid(rho); spread_q_poisson(log,bVerbose,TRUE,natoms,x,charge,box,rc,rho,nrnb, bOld,r1); symmetrize_PSgrid(log,rho,0.0); if (bSecond) copy_PSgrid(pot,rho); /* Second step: solving the poisson equation in real space */ (void) solve_poisson(log,pot,rho,bVerbose,nrnb,maxnit,tol,box); symmetrize_PSgrid(log,pot,0.0); /* Third and last step: gather the forces, energies and potential * from the grid. */ #define BETA(n) (BMIN+n*DB) /* Optimization of beta in progress */ for(bx=0; (bx<NB); bx++) { beta[XX] = BETA(bx); for(by=0; (by<NB); by++) { beta[YY] = BETA(by); for(bz=0; (bz<NB); bz++) { beta[ZZ] = BETA(bz); for(i=0; (i<natoms); i++) { phi[i] = 0.0; clear_rvec(f[i]); } ener = ps_gather_f(log,bVerbose,natoms,x,f,charge,box, phi,pot,beta,nrnb); sprintf(buf,"Poisson, beta = %g\n",beta[XX]); rmsf[bx][by][bz] = analyse_diff(log,buf,natoms,f_ref,f,phi_ref,phi,NULL, /*"fcorr.xvg","pcorr.xvg"*/NULL,NULL,NULL,NULL); } } } rmsf_min = rmsf[0][0][0]; minimum[XX] = minimum[YY] = minimum[ZZ] = 0; for(bx=0; (bx<NB); bx++) { beta[XX] = BETA(bx); for(by=0; (by<NB); by++) { beta[YY] = BETA(by); for(bz=0; (bz<NB); bz++) { beta[ZZ] = BETA(bz); rrmsf = rmsf[bx][by][bz]; fprintf(log,"Beta: %6.3f %6.3f %6.3f RMSF: %8.3f\n", beta[XX],beta[YY],beta[ZZ],rrmsf); if (rrmsf < rmsf_min) { rmsf_min = rrmsf; minimum[XX] = bx; minimum[YY] = by; minimum[ZZ] = bz; } } } } beta[XX] = BETA(minimum[XX]); beta[YY] = BETA(minimum[YY]); beta[ZZ] = BETA(minimum[ZZ]); fprintf(log,"Minimum RMSF %8.3f at Beta = %6.3f %6.3f %6.3f\n", rmsf_min,beta[XX],beta[YY],beta[ZZ]); /* Computing optimum once more... */ for(i=0; (i<natoms); i++) { phi[i] = 0.0; clear_rvec(f[i]); } ener = ps_gather_f(log,bVerbose,natoms,x,f,charge,box,phi,pot,beta,nrnb); return ener; }
int main(int argc,char *argv[]) { static char *desc[] = { "testlr tests the PPPM and Ewald method for the", "long range electrostatics problem." }; static t_filenm fnm[] = { { efTPX, NULL, NULL, ffREAD }, { efHAT, "-g", "ghat", ffOPTRD }, { efOUT, "-o", "rho", ffOPTWR }, { efOUT, "-op", "lr-pb", ffOPTWR }, { efOUT, "-of", "lr-four", ffOPTWR }, { efOUT, "-opt", "tot-pb", ffOPTWR }, { efOUT, "-oft", "tot-four", ffOPTWR }, { efOUT, "-fin", "lr-four", ffOPTWR }, { efEPS, "-es", "sr", ffOPTWR }, { efEPS, "-elf", "lr-four", ffOPTWR }, { efEPS, "-etf", "tot-four", ffOPTWR }, { efEPS, "-qr", "qk-real", ffOPTWR }, { efEPS, "-qi", "qk-im", ffOPTWR }, { efEPS, "-elp", "lr-pb", ffOPTWR }, { efEPS, "-etp", "tot-pb", ffOPTWR }, { efEPS, "-rho", "rho", ffOPTWR }, { efEPS, "-qq", "charge", ffOPTWR }, { efXVG, "-gt", "gk-tab", ffOPTWR }, { efXVG, "-fcorr","fcorr", ffWRITE }, { efXVG, "-pcorr","pcorr", ffWRITE }, { efXVG, "-ftotcorr","ftotcorr", ffWRITE }, { efXVG, "-ptotcorr","ptotcorr", ffWRITE }, { efLOG, "-l", "fptest", ffWRITE }, { efXVG, "-gr", "spread", ffOPTWR }, { efPDB, "-pf", "pqr-four", ffOPTWR }, { efPDB, "-phitot", "pppm-phitot", ffOPTWR } }; #define NFILE asize(fnm) FILE *log; t_topology top; t_tpxheader stath; t_inputrec ir; t_block *excl; t_forcerec *fr; t_commrec *cr; t_mdatoms *mdatoms; t_graph *graph; int i,step,nre,natoms,nmol; rvec *x,*f_sr,*f_excl,*f_four,*f_pppm,*f_pois,box_size,hbox; matrix box; real t,lambda,vsr,*charge,*phi_f,*phi_pois,*phi_s,*phi_p3m,*rho; static bool bFour=FALSE,bVerbose=FALSE,bGGhat=FALSE,bPPPM=TRUE, bPoisson=FALSE,bOld=FALSE,bOldEwald=TRUE; static int nprocs = 1; static t_pargs pa[] = { { "-np", FALSE, etINT, &nprocs, "Do it in parallel" }, { "-ewald", FALSE, etBOOL, &bFour, "Do an Ewald solution"}, { "-pppm", FALSE, etBOOL, &bPPPM, "Do a PPPM solution" }, { "-poisson",FALSE, etBOOL, &bPoisson,"Do a Poisson solution" }, { "-v", FALSE, etBOOL, &bVerbose,"Verbose on"}, { "-ghat", FALSE, etBOOL, &bGGhat, "Generate Ghat function"}, { "-old", FALSE, etBOOL, &bOld, "Use old function types"}, { "-oldewald",FALSE,etBOOL, &bOldEwald,"Use old Ewald code"} }; CopyRight(stderr,argv[0]); parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_CAN_VIEW, NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL); if (nprocs > 1) { cr = init_par(&argc,argv); open_log(ftp2fn(efLOG,NFILE,fnm),cr); log = stdlog; } else { cr = init_par(&argc,argv); log = ftp2FILE(efLOG,NFILE,fnm,"w"); stdlog = log; } /* Read topology and coordinates */ read_tpxheader(ftp2fn(efTPX,NFILE,fnm),&stath,FALSE); snew(x,stath.natoms); snew(f_sr,stath.natoms); snew(f_excl,stath.natoms); snew(f_four,stath.natoms); snew(f_pppm,stath.natoms); snew(f_pois,stath.natoms); read_tpx(ftp2fn(efTPX,NFILE,fnm),&step,&t,&lambda,&ir, box,&natoms,x,NULL,NULL,&top); excl=&(top.atoms.excl); nmol=top.blocks[ebMOLS].nr; /* Allocate space for potential, charges and rho (charge density) */ snew(charge,stath.natoms); snew(phi_f,stath.natoms); snew(phi_p3m,stath.natoms); snew(phi_pois,stath.natoms); snew(phi_s,stath.natoms); snew(rho,stath.natoms); /* Set the charges */ for(i=0; (i<natoms); i++) charge[i]=top.atoms.atom[i].q; /* Make a simple box vector instead of tensor */ for(i=0; (i<DIM); i++) box_size[i]=box[i][i]; /* Set some constants */ fr = mk_forcerec(); mdatoms = atoms2md(&(top.atoms),FALSE,FALSE); set_LRconsts(log,ir.rcoulomb_switch,ir.rcoulomb,box_size,fr); init_forcerec(log,fr,&ir,&(top.blocks[ebMOLS]),cr, &(top.blocks[ebCGS]),&(top.idef),mdatoms,box,FALSE); calc_shifts(box,box_size,fr->shift_vec,FALSE); /* Periodicity stuff */ graph = mk_graph(&(top.idef),top.atoms.nr,FALSE,FALSE); shift_self(graph,fr->shift_vec,x); calc_LRcorrections(log,0,natoms,ir.rcoulomb_switch, ir.rcoulomb,charge,excl,x,f_excl,bOld); pr_f("f_excl.dat",natoms,f_excl); /* Compute the short range potential */ put_atoms_in_box(natoms,box,x); vsr=phi_sr(log,natoms,x,charge,ir.rcoulomb, ir.rcoulomb_switch,box_size,phi_s,excl,f_sr,bOld); pr_f("f_sr.dat",natoms,f_sr); /* Plot the short range potential in a matrix */ calc_ener(log,"Short Range",TRUE,nmol,natoms,phi_s,charge,excl); if (bFour) test_four(log,NFILE,fnm,&(top.atoms),&ir,x,f_four,box_size,charge,phi_f, phi_s,nmol,cr,bOld,bOldEwald); if (bPPPM) test_pppm(log,bVerbose,bGGhat,opt2fn("-g",NFILE,fnm), &(top.atoms),&ir,x,f_pppm,charge,box_size,phi_p3m,phi_s,nmol, cr,bOld,&(top.blocks[ebCGS])); if (bPoisson) test_poisson(log,bVerbose, &(top.atoms),&ir,x,f_pois,charge,box_size,phi_pois, phi_s,nmol,cr,bFour,f_four,phi_f,bOld); if (bPPPM && bFour) analyse_diff(log,"PPPM", top.atoms.nr,f_four,f_pppm,phi_f,phi_p3m,phi_s, opt2fn("-fcorr",NFILE,fnm), opt2fn("-pcorr",NFILE,fnm), opt2fn("-ftotcorr",NFILE,fnm), opt2fn("-ptotcorr",NFILE,fnm)); if (bPoisson && bFour) analyse_diff(log,"Poisson", top.atoms.nr,f_four,f_pois,phi_f,phi_pois,phi_s, opt2fn("-fcorr",NFILE,fnm), opt2fn("-pcorr",NFILE,fnm), opt2fn("-ftotcorr",NFILE,fnm), opt2fn("-ptotcorr",NFILE,fnm)); gmx_fio_fclose(log); thanx(stderr); return 0; }