Gromacs
5.1
|
#include "gromacs/legacyheaders/typedefs.h"
#include "gromacs/legacyheaders/types/simple.h"
Implementation of the AdResS method.
Functions | |
real | adress_weight (rvec x, int adresstype, real adressr, real adressw, rvec *ref, struct t_pbc *pbc, t_forcerec *fr) |
calculates the AdResS weight of a particle More... | |
void | update_adress_weights_com (FILE *fplog, int cg0, int cg1, t_block *cgs, rvec x[], t_forcerec *fr, t_mdatoms *mdatoms, struct t_pbc *pbc) |
update the weight of all coarse-grained particles in several charge groups for com vsites More... | |
void | update_adress_weights_cog (t_iparams ip[], t_ilist ilist[], rvec x[], t_forcerec *fr, t_mdatoms *mdatoms, struct t_pbc *pbc) |
update the weight of all coarse-grained particles for cog vsites More... | |
void | update_adress_weights_atom (int cg0, int cg1, t_block *cgs, rvec x[], t_forcerec *fr, t_mdatoms *mdatoms, struct t_pbc *pbc) |
update the weight of all coarse-grained particles in several charge groups for atom vsites More... | |
void | update_adress_weights_atom_per_atom (int cg0, int cg1, t_block *cgs, rvec x[], t_forcerec *fr, t_mdatoms *mdatoms, struct t_pbc *pbc) |
update the weight on per atom basis of all coarse-grained particles in several charge groups for atom vsites More... | |
void | adress_thermo_force (int cg0, int cg1, t_block *cgs, rvec x[], rvec f[], t_forcerec *fr, t_mdatoms *mdatoms, struct t_pbc *pbc) |
add AdResS IC thermodynamic force to f_novirsum More... | |
void | adress_set_kernel_flags (int n_ex, int n_hyb, int n_cg, t_mdatoms *mdatoms) |
checks weather a cpu calculates only coarse-grained or explicit interactions More... | |
gmx_bool | egp_explicit (t_forcerec *fr, int egp_nr) |
looks up if a energy group is explicit More... | |
gmx_bool | egp_coarsegrained (t_forcerec *fr, int egp_nr) |
looks up if a energy group is coarse-grained More... | |
void adress_set_kernel_flags | ( | int | n_ex, |
int | n_hyb, | ||
int | n_cg, | ||
t_mdatoms * | mdatoms | ||
) |
checks weather a cpu calculates only coarse-grained or explicit interactions
[in] | n_ex | number of explicit particles |
[in] | n_hyb | number of hybrid particles |
[in] | n_cg | number of coarse-grained particles |
[in,out] | mdatoms | the struct containing all the atoms properties |
void adress_thermo_force | ( | int | cg0, |
int | cg1, | ||
t_block * | cgs, | ||
rvec | x[], | ||
rvec | f[], | ||
t_forcerec * | fr, | ||
t_mdatoms * | mdatoms, | ||
struct t_pbc * | pbc | ||
) |
add AdResS IC thermodynamic force to f_novirsum
[in] | cg0 | first charge group to update |
[in] | cg1 | last+1 charge group to update |
[in] | cgs | block containing the cg index |
[in] | x | array with all the particle positions |
[in,out] | f | the force array pointing at f_novirsum from sim_util.c |
[in] | fr | the forcerec containing all the parameters |
[in] | mdatoms | the struct containing all the atoms properties |
[in] | pbc | for shortest distance to fr->adress_refs |
real adress_weight | ( | rvec | x, |
int | adresstype, | ||
real | adressr, | ||
real | adressw, | ||
rvec * | ref, | ||
struct t_pbc * | pbc, | ||
t_forcerec * | fr | ||
) |
calculates the AdResS weight of a particle
[in] | x | position of the particle |
[in] | adresstype | type of address weight function eAdressOff - explicit simulation eAdressConst - constant weight all over the box eAdressXSplit - split in x direction with ref as center eAdressSphere - spherical splitting with ref as center else - weight = 1 - explicit simulation |
[in] | adressr | radius/size of the explicit zone |
[in] | adressw | size of the hybrid zone |
[in] | ref | center of the explicit zone for adresstype 1 - unused for adresstype 2 - only ref[0] is used |
[in] | pbc | pbc struct for calculating shortest distance |
[in] | fr | the forcerec containing all the parameters |
gmx_bool egp_coarsegrained | ( | t_forcerec * | fr, |
int | egp_nr | ||
) |
looks up if a energy group is coarse-grained
[in] | fr | the forcerec containing all the parameters |
[in] | egp_nr | energy group number |
gmx_bool egp_explicit | ( | t_forcerec * | fr, |
int | egp_nr | ||
) |
looks up if a energy group is explicit
[in] | fr | the forcerec containing all the parameters |
[in] | egp_nr | energy group number |
void update_adress_weights_atom | ( | int | cg0, |
int | cg1, | ||
t_block * | cgs, | ||
rvec | x[], | ||
t_forcerec * | fr, | ||
t_mdatoms * | mdatoms, | ||
struct t_pbc * | pbc | ||
) |
update the weight of all coarse-grained particles in several charge groups for atom vsites
[in] | cg0 | first charge group to update |
[in] | cg1 | last+1 charge group to update |
[in] | cgs | block containing the cg index |
[in] | x | array with all the particle positions |
[in] | fr | the forcerec containing all the parameters |
[in,out] | mdatoms | the struct containing all the atoms properties |
[in] | pbc | for shortest distance in adress_weight |
void update_adress_weights_atom_per_atom | ( | int | cg0, |
int | cg1, | ||
t_block * | cgs, | ||
rvec | x[], | ||
t_forcerec * | fr, | ||
t_mdatoms * | mdatoms, | ||
struct t_pbc * | pbc | ||
) |
update the weight on per atom basis of all coarse-grained particles in several charge groups for atom vsites
[in] | cg0 | first charge group to update |
[in] | cg1 | last+1 charge group to update |
[in] | cgs | block containing the cg index |
[in] | x | array with all the particle positions |
[in] | fr | the forcerec containing all the parameters |
[in,out] | mdatoms | the struct containing all the atoms properties |
[in] | pbc | for shortest distance in adress_weight |
void update_adress_weights_cog | ( | t_iparams | ip[], |
t_ilist | ilist[], | ||
rvec | x[], | ||
t_forcerec * | fr, | ||
t_mdatoms * | mdatoms, | ||
struct t_pbc * | pbc | ||
) |
update the weight of all coarse-grained particles for cog vsites
[in] | ip | contains interaction parameters, in this case the number of constructing atoms n for vsitesn |
[in] | ilist | list of interaction types, in this case the virtual site types are what's important |
[in] | x | array with all the particle positions |
[in] | fr | the forcerec containing all the parameters |
[in,out] | mdatoms | the struct containing all the atoms properties |
[in] | pbc | for shortest distance in adress_weight |
void update_adress_weights_com | ( | FILE * | fplog, |
int | cg0, | ||
int | cg1, | ||
t_block * | cgs, | ||
rvec | x[], | ||
t_forcerec * | fr, | ||
t_mdatoms * | mdatoms, | ||
struct t_pbc * | pbc | ||
) |
update the weight of all coarse-grained particles in several charge groups for com vsites
[in,out] | fplog | log file in case of debug |
[in] | cg0 | first charge group to update |
[in] | cg1 | last+1 charge group to update |
[in] | cgs | block containing the cg index |
[in] | x | array with all the particle positions |
[in] | fr | the forcerec containing all the parameters |
[in,out] | mdatoms | the struct containing all the atoms properties |
[in] | pbc | for shortest distance in adress_weight |