STORMM Source Documentation
Loading...
Searching...
No Matches
atomgraph_abstracts.h
1// -*-c++-*-
2#ifndef STORMM_ATOMGRAPH_ABSTRACTS_H
3#define STORMM_ATOMGRAPH_ABSTRACTS_H
4
5#include "copyright.h"
6#include "atomgraph_enumerators.h"
7
8namespace stormm {
9namespace topology {
10
12template <typename T> struct UreyBradleyTerm {
13 int i_atom;
14 int k_atom;
16 T keq;
17 T leq;
18};
19
21template <typename T> struct CharmmImprTerm {
22 int i_atom;
23 int j_atom;
24 int k_atom;
25 int l_atom;
27 T keq;
29};
30
32template <typename T> struct CmapTerm {
33 int i_atom;
34 int j_atom;
35 int k_atom;
36 int l_atom;
37 int m_atom;
41 T* surf;
42};
43
45template <typename T> struct BondTerm {
46 int i_atom;
47 int j_atom;
49 T keq;
50 T leq;
51};
52
54template <typename T> struct AngleTerm {
55 int i_atom;
56 int j_atom;
57 int k_atom;
59 T keq;
61};
62
67template <typename T> struct DihedralTerm {
68 int i_atom;
69 int j_atom;
70 int k_atom;
71 int l_atom;
80};
81
88template <typename T> struct ValenceKit {
89
96 explicit ValenceKit(int natom_in, int nbond_in, int nangl_in, int ndihe_in, int nbond_param_in,
97 int nangl_param_in, int ndihe_param_in, int ninfr14_in, int nattn14_param_in,
98 int nubrd_in, int ncimp_in, int ncmap_in, int nubrd_param_in,
99 int ncimp_param_in, int ncmap_surf_in, const T* bond_keq_in,
100 const T* bond_leq_in, const T* angl_keq_in,
101 const T* angl_theta_in, const T* dihe_amp_in, const T* dihe_freq_in,
102 const T* dihe_phi_in, const T* attn14_elec_in, const T* attn14_vdw_in,
103 const int* bond_i_atoms_in, const int* bond_j_atoms_in,
104 const int* bond_param_idx_in, const char4* bond_modifiers_in,
105 const uint* bond_relevance_in, const int* angl_i_atoms_in,
106 const int* angl_j_atoms_in, const int* angl_k_atoms_in,
107 const int* angl_param_idx_in, const char4* angl_modifiers_in,
108 const uint* angl_relevance_in, const int* dihe_i_atoms_in,
109 const int* dihe_j_atoms_in, const int* dihe_k_atoms_in,
110 const int* dihe_l_atoms_in, const int* dihe_param_idx_in,
111 const int* dihe14_param_idx_in, const char4* dihe_modifiers_in,
112 const int* infr14_i_atoms_in, const int* infr14_l_atoms_in,
113 const int* infr14_param_idx_in, const int* ubrd_i_atoms_in,
114 const int* ubrd_k_atoms_in, const int* ubrd_param_idx_in,
115 const int* cimp_i_atoms_in, const int* cimp_j_atoms_in,
116 const int* cimp_k_atoms_in, const int* cimp_l_atoms_in,
117 const int* cimp_param_idx_in, const int* cmap_i_atoms_in,
118 const int* cmap_j_atoms_in, const int* cmap_k_atoms_in,
119 const int* cmap_l_atoms_in, const int* cmap_m_atoms_in,
120 const int* cmap_dim_in, const int* cmap_surf_bounds_in,
121 const int* cmap_patch_bounds_in, const int* cmap_surf_idx_in,
122 const T* ubrd_keq_in, const T* ubrd_leq_in, const T* cimp_keq_in,
123 const T* cimp_phi_in, const T* cmap_surf_in, const T* cmap_dphi_in,
124 const T* cmap_dpsi_in, const T* cmap_dphi_dpsi_in, const T* cmap_patches_in,
125 const int* bond_asgn_atoms_in, const int* bond_asgn_index_in,
126 const int* bond_asgn_terms_in, const int* bond_asgn_bounds_in,
127 const int* angl_asgn_atoms_in, const int* angl_asgn_index_in,
128 const int* angl_asgn_terms_in, const int* angl_asgn_bounds_in,
129 const int* dihe_asgn_atoms_in, const int* dihe_asgn_index_in,
130 const int* dihe_asgn_terms_in, const int* dihe_asgn_bounds_in,
131 const int* ubrd_asgn_atoms_in, const int* ubrd_asgn_index_in,
132 const int* ubrd_asgn_terms_in, const int* ubrd_asgn_bounds_in,
133 const int* cimp_asgn_atoms_in, const int* cimp_asgn_index_in,
134 const int* cimp_asgn_terms_in, const int* cimp_asgn_bounds_in,
135 const int* cmap_asgn_atoms_in, const int* cmap_asgn_index_in,
136 const int* cmap_asgn_terms_in, const int* cmap_asgn_bounds_in);
137
141 ValenceKit(const ValenceKit &original) = default;
142 ValenceKit(ValenceKit &&other) = default;
144
145 // The purpose of this struct is to store a collection of pointers for HPC kernels. As such, it
146 // does not have any private member variables. As a provider of parameters, it also does not
147 // allow modification of the data that it points to.
148 const int natom;
151 const int nbond;
152 const int nangl;
153 const int ndihe;
155 const int nbond_param;
156 const int nangl_param;
157 const int ndihe_param;
158 const int ninfr14;
159 const int nattn14_param;
160 const int nubrd;
161 const int ncimp;
162 const int ncmap;
163 const int nubrd_param;
164 const int ncimp_param;
165 const int ncmap_surf;
166 const T* bond_keq;
167 const T* bond_leq;
168 const T* angl_keq;
169 const T* angl_theta;
170 const T* dihe_amp;
171 const T* dihe_freq;
172 const T* dihe_phi;
173 const T* attn14_elec;
176 const T* attn14_vdw;
179 const int* bond_i_atoms;
180 const int* bond_j_atoms;
181 const int* bond_param_idx;
183 const uint* bond_relevance;
186 const int* angl_i_atoms;
187 const int* angl_j_atoms;
188 const int* angl_k_atoms;
189 const int* angl_param_idx;
191 const uint* angl_relevance;
194 const int* dihe_i_atoms;
195 const int* dihe_j_atoms;
196 const int* dihe_k_atoms;
199 const int* dihe_l_atoms;
200 const int* dihe_param_idx;
201 const int* dihe14_param_idx;
206 const int* infr14_i_atoms;
209 const int* infr14_l_atoms;
210 const int* infr14_param_idx;
212 const int* ubrd_i_atoms;
213 const int* ubrd_k_atoms;
214 const int* ubrd_param_idx;
215 const int* cimp_i_atoms;
216 const int* cimp_j_atoms;
217 const int* cimp_k_atoms;
218 const int* cimp_l_atoms;
219 const int* cimp_param_idx;
220 const int* cmap_i_atoms;
221 const int* cmap_j_atoms;
223 const int* cmap_k_atoms;
225 const int* cmap_l_atoms;
227 const int* cmap_m_atoms;
228 const int* cmap_dim;
229 const int* cmap_surf_bounds;
230 const int* cmap_patch_bounds;
232 const int* cmap_surf_idx;
233 const T* ubrd_keq;
234 const T* ubrd_leq;
235 const T* cimp_keq;
236 const T* cimp_phi;
237 const T* cmap_surf;
238 const T* cmap_dphi;
239 const T* cmap_dpsi;
240 const T* cmap_dphi_dpsi;
241 const T* cmap_patches;
242
243 // The following arrays provide a different perspective on the valence terms. Each bond, angle,
244 // dihedral, Urey-Bradley, CMAP, and other valence term is assigned to one of its constituent
245 // atoms. Atoms are then responsible for executing any of their assigned terms. Atom a_idx
246 // is responsible for terms in each array indexed by [(term)_asgn_bounds[a_idx] ...
247 // (term)_asgn_bounds[a_idx] + 1), which provides the relevant elements of the associated
248 // (term)_asgn_index array (for the parameter index) or the (term)_assigned_atoms array
249 // (multiply by (number of atoms in the term - 1), but this lists all other atoms in the term).
250 // This allows instant lookup of all relevant bonds, angles, CMAPs, etc. based on any given atom
251 // number, although to get all bonds in a given molecule one still has to loop over all atoms in
252 // the molecule. To find all bonds associated with a given atom a_idx, look for all bonds that
253 // a_idx itself is responsible for and then cycle through its 1:2 exclusions neighbor list to see
254 // which of those atom are responsible for bonds back to a_idx. That will enumerate all bonds,
255 // with parameters, affecting atom a_idx.
256 const int* bond_asgn_atoms;
257 const int* bond_asgn_index;
258 const int* bond_asgn_terms;
260 const int* bond_asgn_bounds;
261 const int* angl_asgn_atoms;
262 const int* angl_asgn_index;
263 const int* angl_asgn_terms;
265 const int* angl_asgn_bounds;
266 const int* dihe_asgn_atoms;
267 const int* dihe_asgn_index;
268 const int* dihe_asgn_terms;
270 const int* dihe_asgn_bounds;
271 const int* ubrd_asgn_atoms;
272 const int* ubrd_asgn_index;
273 const int* ubrd_asgn_terms;
274 const int* ubrd_asgn_bounds;
275 const int* cimp_asgn_atoms;
276 const int* cimp_asgn_index;
277 const int* cimp_asgn_terms;
278 const int* cimp_asgn_bounds;
279 const int* cmap_asgn_atoms;
280 const int* cmap_asgn_index;
281 const int* cmap_asgn_terms;
282 const int* cmap_asgn_bounds;
283};
284
287template <typename T> struct NonbondedKit {
288
291 explicit NonbondedKit(int natom_in, int n_q_types_in, int n_lj_types_in,
292 const T coulomb_constant_in, const T* charge_in, const int* q_idx_in,
293 const int* lj_idx_in, const T* q_parameter_in, const T* lja_coeff_in,
294 const T* ljb_coeff_in, const T* ljc_coeff_in, const T* lja_14_coeff_in,
295 const T* ljb_14_coeff_in, const T* ljc_14_coeff_in, const T* lj_sigma,
296 const T* lj_14_sigma, const int* nb11x_in, const int* nb11_bounds_in,
297 const int* nb12x_in, const int* nb12_bounds_in, const int* nb13x_in,
298 const int* nb13_bounds_in, const int* nb14x_in, const int* nb14_bounds_in,
299 const T* lj_type_corr_in);
300
304 NonbondedKit(const NonbondedKit &original) = default;
305 NonbondedKit(NonbondedKit &&other) = default;
307
308 // Member variables again store a collection of atomic parameters
309 const int natom;
310 const int n_q_types;
311 const int n_lj_types;
313 const T* charge;
314 const int* q_idx;
315 const int* lj_idx;
316 const T* q_parameter;
326 const T* lja_coeff;
328 const T* ljb_coeff;
329 const T* ljc_coeff;
330 const T* lja_14_coeff;
332 const T* ljb_14_coeff;
333 const T* ljc_14_coeff;
334 const T* lj_sigma;
335 const T* lj_14_sigma;
336 const int* nb11x;
337 const int* nb11_bounds;
338 const int* nb12x;
339 const int* nb12_bounds;
340 const int* nb13x;
342 const int* nb13_bounds;
343 const int* nb14x;
344 const int* nb14_bounds;
345 const T* lj_type_corr;
347};
348
352template <typename T> struct ImplicitSolventKit {
353
354 explicit ImplicitSolventKit(int natom_in, ImplicitSolventModel igb_in, T dielectric_in,
355 T saltcon_in, const int* neck_gb_idx_in, const T* pb_radii_in,
356 const T* gb_screen_in, const T* gb_alpha_in,
357 const T* gb_beta_in, const T* gb_gamma_in);
358
361 ImplicitSolventKit(const ImplicitSolventKit &original) = default;
362 ImplicitSolventKit(ImplicitSolventKit &&other) = default;
364
365 // Member variables again store a collection of atomic parameters
366 const int natom;
367 const ImplicitSolventModel igb;
369 const T dielectric;
370 const T saltcon;
371 const int* neck_gb_idx;
373 const T* pb_radii;
374 const T* gb_screen;
375 const T* gb_alpha;
376 const T* gb_beta;
377 const T* gb_gamma;
378};
379
383
386 ChemicalDetailsKit(int natom_in, int nres_in, int nmol_in, int free_dof_in, int cnst_dof_in,
387 const char4* atom_names_in, const char4* res_names_in,
388 const char4* atom_types_in, const int* z_numbers_in, const int* res_limits_in,
389 const int* atom_numbers_in, const int* res_numbers_in, const int* mol_home_in,
390 const int* mol_contents_in, const int* mol_limits_in, const double* masses_in,
391 const float* sp_masses_in, const double* inv_masses_in,
392 const float* sp_inv_masses_in);
393
397 ChemicalDetailsKit(const ChemicalDetailsKit &original) = default;
398 ChemicalDetailsKit(ChemicalDetailsKit &&other) = default;
400
401 // Member variables store the atom count, residue count, and other ways to quantify the system
402 // in addition to many useful pointers
403 const int natom;
404 const int nres;
405 const int nmol;
406 const int free_dof;
408 const int cnst_dof;
413 const int* z_numbers;
414 const int* res_limits;
415 const int* atom_numbers;
416 const int* res_numbers;
418 const int* mol_home;
419 const int* mol_contents;
421 const int* mol_limits;
422 const double* masses;
423 const float* sp_masses;
424 const double* inv_masses;
425 const float* sp_inv_masses;
426};
427
430template <typename T> struct VirtualSiteKit {
431
434 explicit VirtualSiteKit(int nsite_in, int nframe_set_in, const int* vs_atoms_in,
435 const int* frame1_idx_in, const int* frame2_idx_in,
436 const int* frame3_idx_in, const int* frame4_idx_in,
437 const int* vs_param_idx_in, const int* vs_types_in, const T* dim1_in,
438 const T* dim2_in, const T* dim3_in);
439
443 VirtualSiteKit(const VirtualSiteKit &original) = default;
444 VirtualSiteKit(VirtualSiteKit &&other) = default;
446
447 const int nsite;
448 const int nframe_set;
449 const int* vs_atoms;
450 const int* frame1_idx;
451 const int* frame2_idx;
452 const int* frame3_idx;
453 const int* frame4_idx;
454 const int* vs_param_idx;
457 const int* vs_types;
458 const T* dim1;
459 const T* dim2;
460 const T* dim3;
461};
462
466template <typename T> struct ConstraintKit {
467
470 explicit ConstraintKit(const int nsettle_in, const int nsett_param_in, const int ngroup_in,
471 const int ncnst_param_in, const int* settle_ox_atoms_in,
472 const int* settle_h1_atoms_in, const int* settle_h2_atoms_in,
473 const int* settle_param_idx_in, const int* group_list_in,
474 const int* group_bounds_in, const int* group_param_idx_in,
475 const int* group_param_bounds_in, const T* settle_mo_in,
476 const T* settle_mh_in, const T* settle_moh_in, const T* settle_mormt_in,
477 const T* settle_mhrmt_in, const T* settle_ra_in, const T* settle_rb_in,
478 const T* settle_rc_in, const T* group_sq_lengths_in,
479 const T* group_inv_masses_in);
480
484 ConstraintKit(const ConstraintKit &original) = default;
485 ConstraintKit(ConstraintKit &&other) = default;
487
488 const int nsettle;
489 const int nsett_param;
490 const int ngroup;
491 const int ncnst_param;
492 const int* settle_ox_atoms;
493 const int* settle_h1_atoms;
494 const int* settle_h2_atoms;
495 const int* settle_param_idx;
496 const int* group_list;
501 const int* group_bounds;
504 const int* group_param_idx;
511 const T* settle_mo;
512 const T* settle_mh;
513 const T* settle_moh;
514 const T* settle_mormt;
518 const T* settle_mhrmt;
519 const T* settle_ra;
520 const T* settle_rb;
521 const T* settle_rc;
524};
525
526} // namespace topology
527} // namespace stormm
528
529#include "atomgraph_abstracts.tpp"
530
531#endif
ChemicalDetailsKit(int natom_in, int nres_in, int nmol_in, int free_dof_in, int cnst_dof_in, const char4 *atom_names_in, const char4 *res_names_in, const char4 *atom_types_in, const int *z_numbers_in, const int *res_limits_in, const int *atom_numbers_in, const int *res_numbers_in, const int *mol_home_in, const int *mol_contents_in, const int *mol_limits_in, const double *masses_in, const float *sp_masses_in, const double *inv_masses_in, const float *sp_inv_masses_in)
Simple constructor based on many detais. This, like other abstracts, will most likely be produced by ...
Definition atomgraph_abstracts.cpp:9
NonbondedKit(int natom_in, int n_q_types_in, int n_lj_types_in, const T coulomb_constant_in, const T *charge_in, const int *q_idx_in, const int *lj_idx_in, const T *q_parameter_in, const T *lja_coeff_in, const T *ljb_coeff_in, const T *ljc_coeff_in, const T *lja_14_coeff_in, const T *ljb_14_coeff_in, const T *ljc_14_coeff_in, const T *lj_sigma, const T *lj_14_sigma, const int *nb11x_in, const int *nb11_bounds_in, const int *nb12x_in, const int *nb12_bounds_in, const int *nb13x_in, const int *nb13_bounds_in, const int *nb14x_in, const int *nb14_bounds_in, const T *lj_type_corr_in)
As with most other astracts, the constructor is the only member function of significance....
ValenceKit(int natom_in, int nbond_in, int nangl_in, int ndihe_in, int nbond_param_in, int nangl_param_in, int ndihe_param_in, int ninfr14_in, int nattn14_param_in, int nubrd_in, int ncimp_in, int ncmap_in, int nubrd_param_in, int ncimp_param_in, int ncmap_surf_in, const T *bond_keq_in, const T *bond_leq_in, const T *angl_keq_in, const T *angl_theta_in, const T *dihe_amp_in, const T *dihe_freq_in, const T *dihe_phi_in, const T *attn14_elec_in, const T *attn14_vdw_in, const int *bond_i_atoms_in, const int *bond_j_atoms_in, const int *bond_param_idx_in, const char4 *bond_modifiers_in, const uint *bond_relevance_in, const int *angl_i_atoms_in, const int *angl_j_atoms_in, const int *angl_k_atoms_in, const int *angl_param_idx_in, const char4 *angl_modifiers_in, const uint *angl_relevance_in, const int *dihe_i_atoms_in, const int *dihe_j_atoms_in, const int *dihe_k_atoms_in, const int *dihe_l_atoms_in, const int *dihe_param_idx_in, const int *dihe14_param_idx_in, const char4 *dihe_modifiers_in, const int *infr14_i_atoms_in, const int *infr14_l_atoms_in, const int *infr14_param_idx_in, const int *ubrd_i_atoms_in, const int *ubrd_k_atoms_in, const int *ubrd_param_idx_in, const int *cimp_i_atoms_in, const int *cimp_j_atoms_in, const int *cimp_k_atoms_in, const int *cimp_l_atoms_in, const int *cimp_param_idx_in, const int *cmap_i_atoms_in, const int *cmap_j_atoms_in, const int *cmap_k_atoms_in, const int *cmap_l_atoms_in, const int *cmap_m_atoms_in, const int *cmap_dim_in, const int *cmap_surf_bounds_in, const int *cmap_patch_bounds_in, const int *cmap_surf_idx_in, const T *ubrd_keq_in, const T *ubrd_leq_in, const T *cimp_keq_in, const T *cimp_phi_in, const T *cmap_surf_in, const T *cmap_dphi_in, const T *cmap_dpsi_in, const T *cmap_dphi_dpsi_in, const T *cmap_patches_in, const int *bond_asgn_atoms_in, const int *bond_asgn_index_in, const int *bond_asgn_terms_in, const int *bond_asgn_bounds_in, const int *angl_asgn_atoms_in, const int *angl_asgn_index_in, const int *angl_asgn_terms_in, const int *angl_asgn_bounds_in, const int *dihe_asgn_atoms_in, const int *dihe_asgn_index_in, const int *dihe_asgn_terms_in, const int *dihe_asgn_bounds_in, const int *ubrd_asgn_atoms_in, const int *ubrd_asgn_index_in, const int *ubrd_asgn_terms_in, const int *ubrd_asgn_bounds_in, const int *cimp_asgn_atoms_in, const int *cimp_asgn_index_in, const int *cimp_asgn_terms_in, const int *cimp_asgn_bounds_in, const int *cmap_asgn_atoms_in, const int *cmap_asgn_index_in, const int *cmap_asgn_terms_in, const int *cmap_asgn_bounds_in)
The constructor takes a Fortran77-worthy list of arguments. It would not be any cleaner to make this ...
Definition stormm_vector_types.h:141
VirtualSiteKit(int nsite_in, int nframe_set_in, const int *vs_atoms_in, const int *frame1_idx_in, const int *frame2_idx_in, const int *frame3_idx_in, const int *frame4_idx_in, const int *vs_param_idx_in, const int *vs_types_in, const T *dim1_in, const T *dim2_in, const T *dim3_in)
The constructor follows other abstracts and is produced based on pointers from an AtomGraph object.
ConstraintKit(const int nsettle_in, const int nsett_param_in, const int ngroup_in, const int ncnst_param_in, const int *settle_ox_atoms_in, const int *settle_h1_atoms_in, const int *settle_h2_atoms_in, const int *settle_param_idx_in, const int *group_list_in, const int *group_bounds_in, const int *group_param_idx_in, const int *group_param_bounds_in, const T *settle_mo_in, const T *settle_mh_in, const T *settle_moh_in, const T *settle_mormt_in, const T *settle_mhrmt_in, const T *settle_ra_in, const T *settle_rb_in, const T *settle_rc_in, const T *group_sq_lengths_in, const T *group_inv_masses_in)
The constructor follows other abstracts and is produced based on pointers from an AtomGraph object.
Unguarded struct to store the ingredients of a single bond angle interaction.
Definition atomgraph_abstracts.h:54
int k_atom
Third atom in the term.
Definition atomgraph_abstracts.h:57
int j_atom
Second (central) atom in the term.
Definition atomgraph_abstracts.h:56
int param_idx
Parameter index (for tracing purposes–critical values are given next)
Definition atomgraph_abstracts.h:58
int i_atom
First atom in the term.
Definition atomgraph_abstracts.h:55
T theta_eq
Harmonic equilibrium angle (radians)
Definition atomgraph_abstracts.h:60
T keq
Harmonic angle stiffness.
Definition atomgraph_abstracts.h:59
Unguarded struct to store the ingredients of a single bond stretching interaction.
Definition atomgraph_abstracts.h:45
int j_atom
Atom j in the term.
Definition atomgraph_abstracts.h:47
T keq
Harmonic stretching stiffness.
Definition atomgraph_abstracts.h:49
int i_atom
Atom i in the term.
Definition atomgraph_abstracts.h:46
T leq
Harmonic equilibrium length (Angstroms)
Definition atomgraph_abstracts.h:50
int param_idx
Parameter index (for tracing purposes–critical values are given next)
Definition atomgraph_abstracts.h:48
Unguarded struct to store the ingredients of a single CHARMM improper dihedral.
Definition atomgraph_abstracts.h:21
int i_atom
First atom in the term.
Definition atomgraph_abstracts.h:22
int l_atom
Fourth atom in the term.
Definition atomgraph_abstracts.h:25
T keq
Harmonic planar angle stiffness.
Definition atomgraph_abstracts.h:27
T phi_eq
Harmonic equilibrium angle (radians)
Definition atomgraph_abstracts.h:28
int param_idx
Parameter index (for tracing purposes–critical values are given next)
Definition atomgraph_abstracts.h:26
int j_atom
Second (central) atom in the term.
Definition atomgraph_abstracts.h:23
int k_atom
Third atom in the term.
Definition atomgraph_abstracts.h:24
const int natom
The number of atoms in the system.
Definition atomgraph_abstracts.h:403
const int * mol_home
Molecule index to which each atom belongs.
Definition atomgraph_abstracts.h:418
const int * mol_limits
Molecule limits, the bounds by which to read mol_contents.
Definition atomgraph_abstracts.h:421
const int * res_numbers
Definition atomgraph_abstracts.h:416
ChemicalDetailsKit(int natom_in, int nres_in, int nmol_in, int free_dof_in, int cnst_dof_in, const char4 *atom_names_in, const char4 *res_names_in, const char4 *atom_types_in, const int *z_numbers_in, const int *res_limits_in, const int *atom_numbers_in, const int *res_numbers_in, const int *mol_home_in, const int *mol_contents_in, const int *mol_limits_in, const double *masses_in, const float *sp_masses_in, const double *inv_masses_in, const float *sp_inv_masses_in)
Simple constructor based on many detais. This, like other abstracts, will most likely be produced by ...
Definition atomgraph_abstracts.cpp:9
const char4 * atom_names
Names of all atoms in the system.
Definition atomgraph_abstracts.h:410
const int * z_numbers
Atomic numbers for all atoms in the system.
Definition atomgraph_abstracts.h:413
ChemicalDetailsKit(const ChemicalDetailsKit &original)=default
Take the default copy and move constructors. The assignment operators will get implicitly deleted as ...
const float * sp_inv_masses
Inverse masses of atoms in the system (single precision)
Definition atomgraph_abstracts.h:425
const double * inv_masses
Inverse masses of atoms in the system.
Definition atomgraph_abstracts.h:424
const double * masses
Masses of atoms in the system.
Definition atomgraph_abstracts.h:422
const int * atom_numbers
Structural atom numbers for every atom.
Definition atomgraph_abstracts.h:415
const char4 * res_names
Names of all residues in the system.
Definition atomgraph_abstracts.h:411
const int free_dof
Definition atomgraph_abstracts.h:406
const float * sp_masses
Masses of atoms in the system (single precision)
Definition atomgraph_abstracts.h:423
const int * mol_contents
Definition atomgraph_abstracts.h:419
const int nres
The number of residues in the system.
Definition atomgraph_abstracts.h:404
const int cnst_dof
Definition atomgraph_abstracts.h:408
const int * res_limits
Residue limits, a capped (exclusive) prefix sum.
Definition atomgraph_abstracts.h:414
const char4 * atom_types
Atom type names for all atoms in the system.
Definition atomgraph_abstracts.h:412
const int nmol
The number of molecules in the system.
Definition atomgraph_abstracts.h:405
Unguarded struct to store the ingredients of a single CMAP interaction.
Definition atomgraph_abstracts.h:32
int surf_idx
Definition atomgraph_abstracts.h:38
int k_atom
Third atom in the first dihedral term, or second atom in the second.
Definition atomgraph_abstracts.h:35
int l_atom
Fourth atom in the first dihedral term, or third atom in the second.
Definition atomgraph_abstracts.h:36
T * surf
Pointer to the surface potential values (stored in column-major format)
Definition atomgraph_abstracts.h:41
int surf_dim
Dimension of the (square, periodic) CMAP potential grid.
Definition atomgraph_abstracts.h:40
int j_atom
Second atom in the first dihedral term, or first atom in the second.
Definition atomgraph_abstracts.h:34
int m_atom
Fourth atom in the second dihedral term.
Definition atomgraph_abstracts.h:37
int i_atom
First atom in the first dihedral term.
Definition atomgraph_abstracts.h:33
const T * group_inv_masses
Inverse masses for the particles in each group.
Definition atomgraph_abstracts.h:523
const T * settle_ra
Internal distance measurement of SETTLE groups.
Definition atomgraph_abstracts.h:519
const T * settle_mhrmt
Proportional mass of "hydrogen" in SETTLE systems.
Definition atomgraph_abstracts.h:518
const int ngroup
Number of "hub and spoke" constrained groups of bonds.
Definition atomgraph_abstracts.h:490
const T * settle_mh
Mass of the hydrogen (or, lightest) atom.
Definition atomgraph_abstracts.h:512
const int ncnst_param
The number of distinct constraint group parameter sets.
Definition atomgraph_abstracts.h:491
const int nsettle
Number of SETTLE (analytic) constrained groups of bonds.
Definition atomgraph_abstracts.h:488
const int * group_bounds
Definition atomgraph_abstracts.h:501
const T * settle_mo
Mass of the oxygen (or, heaviest) atom.
Definition atomgraph_abstracts.h:511
const T * settle_rb
Internal distance measurement of SETTLE groups.
Definition atomgraph_abstracts.h:520
const int * settle_ox_atoms
Central 'oxygen' atoms for analytic SETTLE constraints.
Definition atomgraph_abstracts.h:492
const int * settle_param_idx
Parameter set indices for each SETLLE group.
Definition atomgraph_abstracts.h:495
ConstraintKit(const int nsettle_in, const int nsett_param_in, const int ngroup_in, const int ncnst_param_in, const int *settle_ox_atoms_in, const int *settle_h1_atoms_in, const int *settle_h2_atoms_in, const int *settle_param_idx_in, const int *group_list_in, const int *group_bounds_in, const int *group_param_idx_in, const int *group_param_bounds_in, const T *settle_mo_in, const T *settle_mh_in, const T *settle_moh_in, const T *settle_mormt_in, const T *settle_mhrmt_in, const T *settle_ra_in, const T *settle_rb_in, const T *settle_rc_in, const T *group_sq_lengths_in, const T *group_inv_masses_in)
The constructor follows other abstracts and is produced based on pointers from an AtomGraph object.
const int * group_param_idx
Definition atomgraph_abstracts.h:504
const T * settle_moh
Combined mass of the heaviest and one of the two light atoms.
Definition atomgraph_abstracts.h:513
const int * settle_h1_atoms
First 'hydrogen' atoms for analytic SETTLE constraints.
Definition atomgraph_abstracts.h:493
const int nsett_param
The number of distinct SETTLE parameter sets.
Definition atomgraph_abstracts.h:489
ConstraintKit(const ConstraintKit &original)=default
Take the default copy and move constructors. The assignment operators will get implicitly deleted as ...
const int * settle_h2_atoms
Second 'hydrogen' atoms for analytic SETTLE constraints.
Definition atomgraph_abstracts.h:494
const T * settle_rc
Internal distance measurement of SETTLE groups.
Definition atomgraph_abstracts.h:521
const T * group_sq_lengths
Bond length targets for the group atoms.
Definition atomgraph_abstracts.h:522
const int * group_list
Definition atomgraph_abstracts.h:496
const int * group_param_bounds
Bounds array for the constraint group parameter arrays below.
Definition atomgraph_abstracts.h:510
const T * settle_mormt
Definition atomgraph_abstracts.h:514
Unguarded struct to store the ingredients of a single dihedral (proper or improper) torsion interacti...
Definition atomgraph_abstracts.h:67
int k_atom
Third atom in the dihedral (central atom if this is an improper)
Definition atomgraph_abstracts.h:70
T periodicity
Periodicity of the dihedral (integral value, but given as a real number)
Definition atomgraph_abstracts.h:75
T phase
Phase angle for the dihedral.
Definition atomgraph_abstracts.h:74
int j_atom
Second atom in the dihedral.
Definition atomgraph_abstracts.h:69
int param_idx
Parameter index (for tracing purposes–critical values are given next)
Definition atomgraph_abstracts.h:72
int l_atom
Fourth atom in the dihedral.
Definition atomgraph_abstracts.h:71
T amplitude
Cosine series amplitude.
Definition atomgraph_abstracts.h:73
T elec_screen
Definition atomgraph_abstracts.h:76
int i_atom
First atom in the dihedral.
Definition atomgraph_abstracts.h:68
T vdw_screen
Definition atomgraph_abstracts.h:78
const T * gb_alpha
Generalized born alpha parameters (one parameter per atom)
Definition atomgraph_abstracts.h:375
const int * neck_gb_idx
Definition atomgraph_abstracts.h:371
const T saltcon
The salt concentration to use in GB calculations.
Definition atomgraph_abstracts.h:370
const T * gb_gamma
Generalized born gamma parameters (one parameter per atom)
Definition atomgraph_abstracts.h:377
const int natom
The number of atoms in the system.
Definition atomgraph_abstracts.h:366
const T dielectric
The dielectric constant to use (default 80.0)
Definition atomgraph_abstracts.h:369
const ImplicitSolventModel igb
Definition atomgraph_abstracts.h:367
const T * gb_beta
Generalized born beta parameters (one parameter per atom)
Definition atomgraph_abstracts.h:376
ImplicitSolventKit(const ImplicitSolventKit &original)=default
Take the default copy and move constructors as well as assignment operators.
const T * gb_screen
Generalized Born screening factors.
Definition atomgraph_abstracts.h:374
const T * pb_radii
Atomic PB radii.
Definition atomgraph_abstracts.h:373
const T * ljb_14_coeff
Lennard_jones B coefficients for 1:4 interactions.
Definition atomgraph_abstracts.h:332
const T * lja_coeff
Definition atomgraph_abstracts.h:326
const int * nb13_bounds
Non-bonded 1:3 exclusion array bounds for each atom.
Definition atomgraph_abstracts.h:342
const T * ljb_coeff
Lennard_jones B coefficients, again tabulated as a square matrix.
Definition atomgraph_abstracts.h:328
const int * nb12x
Non-bonded 1:2 exclusions, applicable to bonded atoms.
Definition atomgraph_abstracts.h:338
const T * lja_14_coeff
Definition atomgraph_abstracts.h:330
const T * lj_sigma
Lennard_jones sigma parameters.
Definition atomgraph_abstracts.h:334
const int * nb11_bounds
Non-bonded 1:1 exclusion array bounds for each atom.
Definition atomgraph_abstracts.h:337
const T coulomb_constant
Coulomb's constant in units of kcal-A/mol-e^2.
Definition atomgraph_abstracts.h:312
const T * ljc_coeff
Lennard_jones C coefficients, again tabulated as a square matrix.
Definition atomgraph_abstracts.h:329
const T * ljc_14_coeff
Lennard_jones C coefficients for 1:4 interactions.
Definition atomgraph_abstracts.h:333
const int * nb12_bounds
Non-bonded 1:2 exclusion array bounds for each atom.
Definition atomgraph_abstracts.h:339
const T * q_parameter
Definition atomgraph_abstracts.h:316
const int n_lj_types
The number of unique Lennard-Jones atom types in the system.
Definition atomgraph_abstracts.h:311
NonbondedKit(const NonbondedKit &original)=default
Take the default copy and move constructors. The assignment operators will get implicitly deleted as ...
const int * q_idx
Partial charge type indices for all atoms.
Definition atomgraph_abstracts.h:314
NonbondedKit(int natom_in, int n_q_types_in, int n_lj_types_in, const T coulomb_constant_in, const T *charge_in, const int *q_idx_in, const int *lj_idx_in, const T *q_parameter_in, const T *lja_coeff_in, const T *ljb_coeff_in, const T *ljc_coeff_in, const T *lja_14_coeff_in, const T *ljb_14_coeff_in, const T *ljc_14_coeff_in, const T *lj_sigma, const T *lj_14_sigma, const int *nb11x_in, const int *nb11_bounds_in, const int *nb12x_in, const int *nb12_bounds_in, const int *nb13x_in, const int *nb13_bounds_in, const int *nb14x_in, const int *nb14_bounds_in, const T *lj_type_corr_in)
As with most other astracts, the constructor is the only member function of significance....
const int n_q_types
The number of unique charge types in the system.
Definition atomgraph_abstracts.h:310
const T * charge
Partial atomic charges on all atoms.
Definition atomgraph_abstracts.h:313
const int * nb11x
Non-bonded 1:1 exclusions, applicable to virtual sites.
Definition atomgraph_abstracts.h:336
const int * nb14_bounds
Non-bonded 1:4 exclusion array bounds for each atom.
Definition atomgraph_abstracts.h:344
const int * nb14x
Non-bonded 1:4 exclusions, applicable to atoms in proper torsions.
Definition atomgraph_abstracts.h:343
const int * nb13x
Definition atomgraph_abstracts.h:340
const int natom
The number of atoms in the system.
Definition atomgraph_abstracts.h:309
const int * lj_idx
Lennard-Jones type indices of all atoms.
Definition atomgraph_abstracts.h:315
const T * lj_14_sigma
Lennard_jones sigma parameters for 1:4 interactions.
Definition atomgraph_abstracts.h:335
const T * lj_type_corr
Definition atomgraph_abstracts.h:345
Unguarded struct to store the ingredients of a single Urey-Bradley interaction.
Definition atomgraph_abstracts.h:12
int param_idx
Parameter index (for tracing purposes–critical values are given next)
Definition atomgraph_abstracts.h:15
int i_atom
Atom i in the term.
Definition atomgraph_abstracts.h:13
T leq
Harmonic equilibrium length (Angstroms)
Definition atomgraph_abstracts.h:17
int k_atom
Atom k in the term (there is no atom j in this harmonic two-body potential)
Definition atomgraph_abstracts.h:14
T keq
Harmonic stretching stiffness.
Definition atomgraph_abstracts.h:16
const int * dihe14_param_idx
Definition atomgraph_abstracts.h:201
const int ndihe
Definition atomgraph_abstracts.h:153
const int * cmap_patch_bounds
Definition atomgraph_abstracts.h:230
const T * dihe_freq
Periodicities for cosine-based dihedral / torsion terms.
Definition atomgraph_abstracts.h:171
const T * ubrd_leq
Array of unique Urey-Bradley 1:3 harmonic equilibrium lengths.
Definition atomgraph_abstracts.h:234
const int * dihe_asgn_index
Parameter indices for dihedral terms.
Definition atomgraph_abstracts.h:267
const int * bond_asgn_terms
Definition atomgraph_abstracts.h:258
const int * ubrd_asgn_terms
Urey-Bradley term indices assigned to each atom.
Definition atomgraph_abstracts.h:273
const T * cmap_surf
Array of CHARMM CMAP surface values.
Definition atomgraph_abstracts.h:237
const int nangl_param
Number of unqiue bond angle parameters.
Definition atomgraph_abstracts.h:156
const int * angl_j_atoms
Array of second atoms in each bond angle.
Definition atomgraph_abstracts.h:187
const int * ubrd_asgn_bounds
Bounds for each atom's assigned Urey-Bradley list.
Definition atomgraph_abstracts.h:274
const int * bond_asgn_atoms
Other atoms for bond terms (one atom per term)
Definition atomgraph_abstracts.h:256
const int * cimp_j_atoms
Second atoms in each CHARMM harmonic improper dihedral.
Definition atomgraph_abstracts.h:216
const int * dihe_param_idx
Parameter indices for each cosine-based dihedral in the system.
Definition atomgraph_abstracts.h:200
const int * cimp_l_atoms
Fourth atoms in each CHARMM harmonic improper dihedral.
Definition atomgraph_abstracts.h:218
const int * dihe_asgn_terms
Definition atomgraph_abstracts.h:268
const int * cmap_asgn_atoms
Other atoms for CMAP terms (four atoms per term)
Definition atomgraph_abstracts.h:279
const int * cmap_dim
Dimension of each CMAP surface (each surface is a square grid)
Definition atomgraph_abstracts.h:228
const int * bond_asgn_bounds
Bounds for each atom's assigned bond list.
Definition atomgraph_abstracts.h:260
const int * dihe_i_atoms
Array of first atoms in each cosine-based dihedral.
Definition atomgraph_abstracts.h:194
const int * angl_asgn_index
Parameter indices for bond angle terms.
Definition atomgraph_abstracts.h:262
const int * cmap_asgn_bounds
Bounds for each atom's assigned CMAP list.
Definition atomgraph_abstracts.h:282
const int * dihe_l_atoms
Array of fourth atoms in each cosine-based dihedral.
Definition atomgraph_abstracts.h:199
const int * ubrd_k_atoms
Second atoms in each Urey-Bradley harmonic angle interaction.
Definition atomgraph_abstracts.h:213
const int nbond
Number of bonds in the system.
Definition atomgraph_abstracts.h:151
const int * cimp_k_atoms
Third (center) atoms in each CHARMM harmonic improper dihedral.
Definition atomgraph_abstracts.h:217
const T * attn14_vdw
Definition atomgraph_abstracts.h:176
const int * ubrd_param_idx
Parameter indices for each Urey-Bradley interaction.
Definition atomgraph_abstracts.h:214
const int ncmap
Number of CHARMM CMAP surface spline-based interactions.
Definition atomgraph_abstracts.h:162
ValenceKit(int natom_in, int nbond_in, int nangl_in, int ndihe_in, int nbond_param_in, int nangl_param_in, int ndihe_param_in, int ninfr14_in, int nattn14_param_in, int nubrd_in, int ncimp_in, int ncmap_in, int nubrd_param_in, int ncimp_param_in, int ncmap_surf_in, const T *bond_keq_in, const T *bond_leq_in, const T *angl_keq_in, const T *angl_theta_in, const T *dihe_amp_in, const T *dihe_freq_in, const T *dihe_phi_in, const T *attn14_elec_in, const T *attn14_vdw_in, const int *bond_i_atoms_in, const int *bond_j_atoms_in, const int *bond_param_idx_in, const char4 *bond_modifiers_in, const uint *bond_relevance_in, const int *angl_i_atoms_in, const int *angl_j_atoms_in, const int *angl_k_atoms_in, const int *angl_param_idx_in, const char4 *angl_modifiers_in, const uint *angl_relevance_in, const int *dihe_i_atoms_in, const int *dihe_j_atoms_in, const int *dihe_k_atoms_in, const int *dihe_l_atoms_in, const int *dihe_param_idx_in, const int *dihe14_param_idx_in, const char4 *dihe_modifiers_in, const int *infr14_i_atoms_in, const int *infr14_l_atoms_in, const int *infr14_param_idx_in, const int *ubrd_i_atoms_in, const int *ubrd_k_atoms_in, const int *ubrd_param_idx_in, const int *cimp_i_atoms_in, const int *cimp_j_atoms_in, const int *cimp_k_atoms_in, const int *cimp_l_atoms_in, const int *cimp_param_idx_in, const int *cmap_i_atoms_in, const int *cmap_j_atoms_in, const int *cmap_k_atoms_in, const int *cmap_l_atoms_in, const int *cmap_m_atoms_in, const int *cmap_dim_in, const int *cmap_surf_bounds_in, const int *cmap_patch_bounds_in, const int *cmap_surf_idx_in, const T *ubrd_keq_in, const T *ubrd_leq_in, const T *cimp_keq_in, const T *cimp_phi_in, const T *cmap_surf_in, const T *cmap_dphi_in, const T *cmap_dpsi_in, const T *cmap_dphi_dpsi_in, const T *cmap_patches_in, const int *bond_asgn_atoms_in, const int *bond_asgn_index_in, const int *bond_asgn_terms_in, const int *bond_asgn_bounds_in, const int *angl_asgn_atoms_in, const int *angl_asgn_index_in, const int *angl_asgn_terms_in, const int *angl_asgn_bounds_in, const int *dihe_asgn_atoms_in, const int *dihe_asgn_index_in, const int *dihe_asgn_terms_in, const int *dihe_asgn_bounds_in, const int *ubrd_asgn_atoms_in, const int *ubrd_asgn_index_in, const int *ubrd_asgn_terms_in, const int *ubrd_asgn_bounds_in, const int *cimp_asgn_atoms_in, const int *cimp_asgn_index_in, const int *cimp_asgn_terms_in, const int *cimp_asgn_bounds_in, const int *cmap_asgn_atoms_in, const int *cmap_asgn_index_in, const int *cmap_asgn_terms_in, const int *cmap_asgn_bounds_in)
The constructor takes a Fortran77-worthy list of arguments. It would not be any cleaner to make this ...
const int nubrd
Number of Urey-Bradley 1:3 harmonic angle interactions.
Definition atomgraph_abstracts.h:160
const int * cmap_k_atoms
Definition atomgraph_abstracts.h:223
const int * cimp_asgn_bounds
Bounds for each atom's assigned CHARMM improper list.
Definition atomgraph_abstracts.h:278
const T * cmap_dpsi
Array of CHARMM CMAP surface "psi" derivative values.
Definition atomgraph_abstracts.h:239
ValenceKit(const ValenceKit &original)=default
Take the default copy and move constructors. The move assignment operator will get implicitly deleted...
const int ninfr14
The number of inferred 1:4 exclusions.
Definition atomgraph_abstracts.h:158
const int * bond_param_idx
Parameter indices for each bond in the system.
Definition atomgraph_abstracts.h:181
const char4 * bond_modifiers
Modifying details of each bond stretching term.
Definition atomgraph_abstracts.h:182
const int * cimp_i_atoms
First atoms in each CHARMM harmonic improper dihedral.
Definition atomgraph_abstracts.h:215
const T * attn14_elec
Definition atomgraph_abstracts.h:173
const T * cmap_dphi_dpsi
Array of CHARMM CMAP surface cross derivative values.
Definition atomgraph_abstracts.h:240
const int * cmap_asgn_terms
CMAP term indices assigned to each atom.
Definition atomgraph_abstracts.h:281
const char4 * angl_modifiers
Modifying details of each angle bending term.
Definition atomgraph_abstracts.h:190
const int * bond_asgn_index
Parameter indices for bond terms.
Definition atomgraph_abstracts.h:257
const T * bond_leq
Equilibrium lengths for all unique bonds.
Definition atomgraph_abstracts.h:167
const int * cmap_asgn_index
Parameter indices for CMAP terms.
Definition atomgraph_abstracts.h:280
const int * infr14_i_atoms
Definition atomgraph_abstracts.h:206
const T * dihe_phi
Phase angles for cosine-based dihedral / torsion terms.
Definition atomgraph_abstracts.h:172
const T * bond_keq
Equilibrium stiffness constants for all unique bonds.
Definition atomgraph_abstracts.h:166
const int * angl_asgn_bounds
Bounds for each atom's assigned bond angle list.
Definition atomgraph_abstracts.h:265
const int * angl_i_atoms
Array of first atoms in each bond angle.
Definition atomgraph_abstracts.h:186
const uint * bond_relevance
Definition atomgraph_abstracts.h:183
const T * cmap_patches
Array of interlaced CHARMM CMAP patch matrices.
Definition atomgraph_abstracts.h:241
const int * cmap_l_atoms
Definition atomgraph_abstracts.h:225
const int * ubrd_i_atoms
First atoms in each Urey-Bradley harmonic angle interaction.
Definition atomgraph_abstracts.h:212
const int * angl_k_atoms
Array of third atoms in each bond angle.
Definition atomgraph_abstracts.h:188
const int nubrd_param
Number of unique Urey-Bradley parameters.
Definition atomgraph_abstracts.h:163
const int * angl_asgn_atoms
Other atoms for bond angle terms (two atoms per term)
Definition atomgraph_abstracts.h:261
const int * infr14_l_atoms
Inferred 1:4 interactions' J atoms.
Definition atomgraph_abstracts.h:209
const int * ubrd_asgn_atoms
Other atoms for Urey-Bradley terms (one atom per term)
Definition atomgraph_abstracts.h:271
const T * ubrd_keq
Array of unique Urey-Bradley 1:3 harmonic stiffnesses.
Definition atomgraph_abstracts.h:233
const int * dihe_k_atoms
Definition atomgraph_abstracts.h:196
const T * angl_theta
Equilibrium lengths for all unique bond angles.
Definition atomgraph_abstracts.h:169
const int * cimp_asgn_atoms
Other atoms for CHARMM impropers (three atoms per term)
Definition atomgraph_abstracts.h:275
const int * cmap_i_atoms
First atoms of the first dihedral in each CHARMM CMAP term.
Definition atomgraph_abstracts.h:220
const int ncmap_surf
Number of unique CMAP surfaces.
Definition atomgraph_abstracts.h:165
const int ndihe_param
Number of unique cosine dihedral parameters.
Definition atomgraph_abstracts.h:157
const int * dihe_asgn_bounds
Bounds for each atom's assigned dihedral list.
Definition atomgraph_abstracts.h:270
const int nbond_param
Number of unique bond parameters.
Definition atomgraph_abstracts.h:155
const int * bond_i_atoms
Array of first atoms in each bond.
Definition atomgraph_abstracts.h:179
const int * ubrd_asgn_index
Parameter indices for Urey-Bradley terms.
Definition atomgraph_abstracts.h:272
const int * dihe_asgn_atoms
Other atoms for dihedral terms (three atoms per term)
Definition atomgraph_abstracts.h:266
const T * cimp_keq
Array of unique CHARMM harmonic improper dihedral stiffnesses.
Definition atomgraph_abstracts.h:235
const int * cmap_surf_bounds
Bounds array for individual CMAP surfaces.
Definition atomgraph_abstracts.h:229
const int nangl
Number of bond angles in the system.
Definition atomgraph_abstracts.h:152
const int * angl_asgn_terms
Definition atomgraph_abstracts.h:263
const T * angl_keq
Equilibrium stiffness constants for all unique bond angles.
Definition atomgraph_abstracts.h:168
const int nattn14_param
Number of 1:4 attenuation factor pairs (vdW, electrostatic)
Definition atomgraph_abstracts.h:159
const int * cimp_param_idx
Parameter indices for each CHARMM harmonic improper dihedral.
Definition atomgraph_abstracts.h:219
const int * cmap_m_atoms
Fourth atoms of the second dihedral in each CHARMM CMAP term.
Definition atomgraph_abstracts.h:227
const T * dihe_amp
Amplitudes for all unique cosine-based dihedrals.
Definition atomgraph_abstracts.h:170
const int * cmap_surf_idx
Parameter (surface) indices for CHARMM CMAP interactions.
Definition atomgraph_abstracts.h:232
const T * cimp_phi
Array of unique CHARMM harmonic improper dihedral phase angles.
Definition atomgraph_abstracts.h:236
const int natom
Definition atomgraph_abstracts.h:148
const int ncimp_param
Number of unique CHARMM harmonic improper parameters.
Definition atomgraph_abstracts.h:164
const int * dihe_j_atoms
Array of second atoms in each cosine-based dihedral.
Definition atomgraph_abstracts.h:195
const int * infr14_param_idx
Definition atomgraph_abstracts.h:210
const int * cimp_asgn_terms
CHARMM improper term indices assigned to each atom.
Definition atomgraph_abstracts.h:277
const int ncimp
Number of CHARMM harmonic improper interactions.
Definition atomgraph_abstracts.h:161
const int * bond_j_atoms
Array of second atoms in each bond.
Definition atomgraph_abstracts.h:180
const T * cmap_dphi
Array of CHARMM CMAP surface "phi" derivative values.
Definition atomgraph_abstracts.h:238
const int * cimp_asgn_index
Parameter indices for CHARMM improper dihedral terms.
Definition atomgraph_abstracts.h:276
const int * angl_param_idx
Parameter indices for each bond angle in the system.
Definition atomgraph_abstracts.h:189
const int * cmap_j_atoms
Definition atomgraph_abstracts.h:221
const char4 * dihe_modifiers
Definition atomgraph_abstracts.h:203
const uint * angl_relevance
Definition atomgraph_abstracts.h:191
const T * dim2
Frame dimension 2.
Definition atomgraph_abstracts.h:459
const int * frame1_idx
Topological indices of frame atom 1 (the parent atom)
Definition atomgraph_abstracts.h:450
const int * vs_types
Virtual site frame types.
Definition atomgraph_abstracts.h:457
const T * dim1
Frame dimension 1.
Definition atomgraph_abstracts.h:458
const int nsite
The number of virtual sites in the topology.
Definition atomgraph_abstracts.h:447
const int nframe_set
The number of unique frame parameter sets in the topology.
Definition atomgraph_abstracts.h:448
const T * dim3
Frame dimension 3.
Definition atomgraph_abstracts.h:460
const int * frame3_idx
Topological indices of frame atom 3.
Definition atomgraph_abstracts.h:452
VirtualSiteKit(const VirtualSiteKit &original)=default
Take the default copy and move constructors. The assignment operators will get implicitly deleted as ...
const int * frame2_idx
Topological indices of frame atom 2.
Definition atomgraph_abstracts.h:451
const int * vs_atoms
Topological indicies of frame atoms.
Definition atomgraph_abstracts.h:449
VirtualSiteKit(int nsite_in, int nframe_set_in, const int *vs_atoms_in, const int *frame1_idx_in, const int *frame2_idx_in, const int *frame3_idx_in, const int *frame4_idx_in, const int *vs_param_idx_in, const int *vs_types_in, const T *dim1_in, const T *dim2_in, const T *dim3_in)
The constructor follows other abstracts and is produced based on pointers from an AtomGraph object.
const int * vs_param_idx
Definition atomgraph_abstracts.h:454
const int * frame4_idx
Topological indices of frame atom 4.
Definition atomgraph_abstracts.h:453