STORMM Source Documentation
Loading...
Searching...
No Matches
atomgraph_synthesis.h
1// -*-c++-*-
2#ifndef STORMM_ATOMGRAPH_SYNTHESIS_H
3#define STORMM_ATOMGRAPH_SYNTHESIS_H
4
5#include "copyright.h"
6#include "Accelerator/hybrid.h"
7#include "Accelerator/gpu_details.h"
8#include "Constants/behavior.h"
9#include "Constants/generalized_born.h"
10#include "DataTypes/stormm_vector_types.h"
11#include "Math/reduction_enumerators.h"
12#include "Math/reduction_workunit.h"
13#include "Potential/energy_enumerators.h"
14#include "Restraints/restraint_apparatus.h"
15#include "Structure/structure_enumerators.h"
16#include "Topology/atomgraph.h"
17#include "Topology/atomgraph_enumerators.h"
18#include "Topology/topology_util.h"
19#include "UnitTesting/stopwatch.h"
20#include "static_mask_synthesis.h"
21#include "synthesis_abstracts.h"
22#include "valence_workunit.h"
23
24namespace stormm {
25namespace synthesis {
26
27using card::GpuDetails;
28using card::Hybrid;
29using card::HybridTargetLevel;
30using constants::ExceptionResponse;
31using energy::ValenceKernelSize;
32using stmath::RdwuPerSystem;
33using stmath::ReductionWorkUnit;
34using restraints::RestraintApparatus;
35using restraints::RestraintKit;
36using structure::ApplyConstraints;
37using topology::AtomGraph;
38using topology::AtomicRadiusSet;
39using topology::ChemicalDetailsKit;
40using topology::getRealParameters;
41using topology::ImplicitSolventModel;
42using topology::MassForm;
43using topology::UnitCellType;
44using testing::StopWatch;
45using namespace generalized_born_defaults;
46
56public:
57
84 AtomGraphSynthesis(const std::vector<AtomGraph*> &topologies_in,
85 const std::vector<RestraintApparatus*> &restraints_in,
86 const std::vector<int> &topology_indices_in,
87 const std::vector<int> &restraint_indices_in,
88 ExceptionResponse policy_in = ExceptionResponse::WARN,
89 const GpuDetails &gpu = null_gpu, StopWatch *timer_in = nullptr);
90
91 AtomGraphSynthesis(const std::vector<AtomGraph*> &topologies_in,
92 const std::vector<RestraintApparatus*> &restraints_in,
93 ExceptionResponse policy_in = ExceptionResponse::WARN,
94 const GpuDetails &gpu = null_gpu,
95 StopWatch *timer_in = nullptr);
96
97 AtomGraphSynthesis(const std::vector<AtomGraph*> &topologies_in,
98 const std::vector<int> &topology_indices_in,
99 ExceptionResponse policy_in = ExceptionResponse::WARN,
100 const GpuDetails &gpu = null_gpu, StopWatch *timer_in = nullptr);
101
102 AtomGraphSynthesis(const std::vector<AtomGraph*> &topologies_in,
103 ExceptionResponse policy_in = ExceptionResponse::WARN,
104 const GpuDetails &gpu = null_gpu,
105 StopWatch *timer_in = nullptr);
107
114 AtomGraphSynthesis(const AtomGraphSynthesis &original);
116 AtomGraphSynthesis& operator=(const AtomGraphSynthesis &original);
117 AtomGraphSynthesis& operator=(AtomGraphSynthesis &&original);
119
121 int getUniqueTopologyCount() const;
122
131 const AtomGraph* getSystemTopologyPointer(int system_index) const;
132 const std::vector<AtomGraph*> getSystemTopologyPointer() const;
134
137 const std::vector<AtomGraph*>& getUniqueTopologies() const;
138
141 std::vector<int> getTopologyIndices() const;
142
147 RestraintApparatus* getSystemRestraintPointer(int system_index) const;
148
150 int getSystemCount() const;
151
161 int getAtomCount() const;
162 int getAtomCount(int system_index) const;
164
166 int getPaddedAtomCount() const;
167
171 int getAtomOffset(int system_index) const;
172
174 int getVirtualSiteCount() const;
175
177 int getBondTermCount() const;
178
180 int getAngleTermCount() const;
181
183 int getDihedralTermCount() const;
184
186 int getUreyBradleyTermCount() const;
187
189 int getCharmmImproperTermCount() const;
190
192 int getCmapTermCount() const;
193
196 int getLJTypeCount() const;
197
199 int getChargeTypeCount() const;
200
202 int getBondParameterCount() const;
203
205 int getAngleParameterCount() const;
206
208 int getDihedralParameterCount() const;
209
212
215
217 int getCmapSurfaceCount() const;
218
221 const Hybrid<int>& getSystemAtomCounts() const;
222
225 const Hybrid<int>& getSystemAtomOffsets() const;
226
228 // synthesis, as a const reference to the internal Hybrid array member variable.
229 const Hybrid<int>& getDegreesOfFreedom() const;
230
232 // synthesis, as a const reference to the internal Hybrid array member variable.
234
237 UnitCellType getUnitCellType() const;
238
240 ImplicitSolventModel getImplicitSolventModel() const;
241
243 double getDielectricConstant() const;
244
246 double getSaltConcentration() const;
247
249 double getCoulombConstant() const;
250
258 std::vector<AtomicRadiusSet> getPBRadiiSet() const;
259 std::vector<AtomicRadiusSet> getPBRadiiSet(int low_limit, int high_limit) const;
260 AtomicRadiusSet getPBRadiiSet(int index) const;
262
275 template <typename T>
276 std::vector<T> getPartialCharges(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
277
278 template <typename T>
279 std::vector<T> getPartialCharges(int system_index,
280 HybridTargetLevel tier = HybridTargetLevel::HOST) const;
281
282 template <typename T>
283 std::vector<T> getPartialCharges(int system_index, int low_index, int high_index,
284 HybridTargetLevel tier = HybridTargetLevel::HOST) const;
286
299 template <typename T> std::vector<T> getAtomicMasses(HybridTargetLevel tier) const;
300 template <typename T> std::vector<T> getAtomicMasses(HybridTargetLevel tier,
301 int system_index) const;
302 template <typename T> std::vector<T> getAtomicMasses(HybridTargetLevel tier, int system_index,
303 int low_index, int high_index) const;
305
308 int getValenceWorkUnitCount() const;
309
314
316 ValenceKernelSize getValenceThreadBlockSize() const;
317
321
323 NbwuKind getNonbondedWorkType() const;
324
326 int getNonbondedWorkUnitCount() const;
327
331 int getRandomCacheDepth() const;
332
334 int getReductionWorkUnitCount() const;
335
338 RdwuPerSystem getRdwuPerSystem() const;
339
346
352 getDoublePrecisionValenceKit(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
353
359 getSinglePrecisionValenceKit(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
360
366 getDoublePrecisionRestraintKit(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
367
373 getSinglePrecisionRestraintKit(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
374
380 getDoublePrecisionNonbondedKit(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
381
387 getSinglePrecisionNonbondedKit(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
388
394 getDoublePrecisionAtomUpdateKit(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
395
401 getSinglePrecisionAtomUpdateKit(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
402
405 const AtomGraphSynthesis* getSelfPointer() const;
406
407#ifdef STORMM_USE_HPC
409 void upload();
410
412 void download();
413#endif
414
430 InitializationTask init_request = InitializationTask::NONE,
431 int random_cache_depth = 0, const GpuDetails &gpu = null_gpu);
432
454
455 void setImplicitSolventModel(ImplicitSolventModel igb_in,
457 const std::vector<AtomicRadiusSet> &radii_sets_in,
458 double dielectric_in = 80.0, double saltcon_in = 0.0,
459 ExceptionResponse policy = ExceptionResponse::WARN);
460
461 void setImplicitSolventModel(ImplicitSolventModel igb_in,
463 AtomicRadiusSet radii_set = AtomicRadiusSet::NONE,
464 double dielectric_in = 80.0, double saltcon_in = 0.0,
465 ExceptionResponse policy = ExceptionResponse::WARN);
466
467 void setImplicitSolventModel(ImplicitSolventModel igb_in,
468 const NeckGeneralizedBornTable &ngb_tab,
469 const std::vector<AtomicRadiusSet> &radii_sets_in,
470 double dielectric_in = 80.0, double saltcon_in = 0.0,
471 ExceptionResponse policy = ExceptionResponse::WARN);
472
473 void setImplicitSolventModel(ImplicitSolventModel igb_in,
474 const NeckGeneralizedBornTable &ngb_tab,
475 AtomicRadiusSet radii_set = AtomicRadiusSet::NONE,
476 double dielectric_in = 80.0, double saltcon_in = 0.0,
477 ExceptionResponse policy = ExceptionResponse::WARN);
479
480private:
481
482 // The first member variables pertain to totals across all systems: atoms, potential function
483 // terms, and sizes of composite parameter arrays that span all topologies
484 ExceptionResponse policy;
485 int topology_count;
486 int restraint_network_count;
487 int system_count;
489 int total_atoms;
490 int total_virtual_sites;
491 int total_bond_terms;
492 int total_angl_terms;
493 int total_dihe_terms;
494 int total_ubrd_terms;
495 int total_cimp_terms;
496 int total_cmap_terms;
497 int total_lj_types;
499 int total_charge_types;
501 int total_bond_params;
502 int total_angl_params;
503 int total_dihe_params;
504 int total_ubrd_params;
505 int total_cimp_params;
506 int total_cmap_surfaces;
507 int total_attn14_params;
509 int total_vste_params;
511 int total_position_restraints;
512 int total_distance_restraints;
513 int total_angle_restraints;
514 int total_dihedral_restraints;
515
516 // The presence of periodic boundaries or an implicit solvent must be common to all systems, as
517 // must the cutoff and other run parameters, although these are contained in other data
518 // structres. If an implicit solvent model is in effect, its parameters are taken as common to
519 // all systems. Other facets of the the system behavior, such as bond constraints, must also
520 // be common to all systems.
521 UnitCellType periodic_box_class;
525 ImplicitSolventModel gb_style;
530 double dielectric_constant;
532 double is_kappa;
534 double salt_concentration;
535 double gb_offset;
536 double gb_neckscale;
537 double gb_neckcut;
538 double coulomb_constant;
541 ApplyConstraints use_shake;
542 ApplyConstraints use_settle;
543 int largest_constraint_group;
545 char4 water_residue_name;
546
551 std::vector<AtomicRadiusSet> pb_radii_sets;
552
556 std::vector<AtomGraph*> topologies;
557
561 std::vector<RestraintApparatus*> restraint_networks;
562
566 std::vector<RestraintApparatus> restraint_dummies;
567
570 Hybrid<int> topology_indices;
571
574 Hybrid<int> restraint_indices;
575
576 // The following arrays are POINTER-kind objects, each directed at a segment of a data array
577 // for storing critical topology sizes. There is one of each of these numbers for every
578 // system in this work plan: each of the Hybrid arrays that follows is system_count in length.
579 // All of the descriptions below follow from the eponymous member variables in the AtomGraph
580 // object and should be taken to mean "... for each topology in the work plan." The pointers
581 // target the ARRAY-kind object int_system_data, at intervals of the number of systems in the
582 // work plan, with no padding for the warp size. If these POINTER-kind objects list bounds, the
583 // bounds themselves will reflect padding for the HPC warp size.
584 Hybrid<int> atom_counts;
585 Hybrid<int> residue_counts;
587 Hybrid<int> molecule_counts;
588 Hybrid<int> largest_residue_sizes;
589 Hybrid<int> last_solute_residues;
592 Hybrid<int> last_solute_atoms;
595 Hybrid<int> first_solvent_molecules;
598
599 // Valence term and off-center particle quantities
600 Hybrid<int> ubrd_term_counts;
601 Hybrid<int> cimp_term_counts;
602 Hybrid<int> cmap_term_counts;
603 Hybrid<int> bond_term_counts;
604 Hybrid<int> angl_term_counts;
605 Hybrid<int> dihe_term_counts;
606 Hybrid<int> virtual_site_counts;
608
609 // Restraint tallies for each system
610 Hybrid<int> posn_restraint_counts;
611 Hybrid<int> bond_restraint_counts;
612 Hybrid<int> angl_restraint_counts;
613 Hybrid<int> dihe_restraint_counts;
614
615 // Information relevant to non-bonded calculations
616 Hybrid<int> lj_type_counts;
617 Hybrid<int> total_exclusion_counts;
618
619 // Information relevant to the MD propagation algorithm
620 Hybrid<int> rigid_water_counts;
621 Hybrid<int> bond_constraint_counts;
622 Hybrid<int> degrees_of_freedom;
624 Hybrid<int> cnst_degrees_of_freedom;
626 Hybrid<int> nonrigid_particle_counts;
629
630 // The following are bounds arrays for the larger data arrays below: more segments of memory,
631 // each the number of systems in length. The contents of these arrays are, with some
632 // exceptions, prefix sums over the associated counts arrays above, taking into account zero
633 // padding for the stuff they are providing indices into. A series of systems with atom counts
634 // 31, 59, and 78 with an HPC warp size of 32 would generate an atom bounds array stating 0, 32,
635 // and 96. There is no capping value, and the (k+1)th bound minus the kth bounds does not
636 // provide information on the exact number associated with the kth system. Programs should read
637 // the exact number of atoms or force field terms per system and the starting bound from one of
638 // the arrays below. The arrays for valence terms spanning all systems provide the global
639 // indices of atoms in unified (but not re-ordered) arrays of atom descriptors and also
640 // coordinates. In this sense, these arrays are intermediates between the valence work units
641 // that follow and the original AtomGraph objects used to construct this AtomGraphSynthesis.
642 Hybrid<int> atom_offsets;
645 Hybrid<int> atom_bit_offsets;
646 Hybrid<int> residue_offsets;
648 Hybrid<int> molecule_offsets;
650 Hybrid<int> ubrd_term_offsets;
652 Hybrid<int> cimp_term_offsets;
653 Hybrid<int> cmap_term_offsets;
654 Hybrid<int> bond_term_offsets;
655 Hybrid<int> angl_term_offsets;
656 Hybrid<int> dihe_term_offsets;
657 Hybrid<int> virtual_site_offsets;
659 Hybrid<int> posn_restraint_offsets;
661 Hybrid<int> bond_restraint_offsets;
663 Hybrid<int> angl_restraint_offsets;
665 Hybrid<int> dihe_restraint_offsets;
667 Hybrid<int> sett_group_offsets;
669 Hybrid<int> cnst_group_offsets;
671 Hybrid<int> nb_exclusion_offsets;
673 Hybrid<int> lennard_jones_abc_offsets;
675
677 Hybrid<int> int_system_data;
678
679 // Maps of each system and its distinguishable parts. From here onwards, the AtomGraphSynthesis
680 // keeps a combination of POINTER-kind and ARRAY-kind Hybrid Object member variables, whatever is
681 // most convenient for grouping and storing data as it is produced. As above, many of the
682 // names match what is in the AtomGraph object. The arrays are laid out such that each system's
683 // details occur in one contiguous stretch, padded by blank elements up to a multiple of the
684 // HPC warp size. Some of these arrays are indices themselves, but will still need a bounds
685 // array, i.e. the point at which to start reading residue limits for the kth system. Bounds
686 // arrays for individual systems within each of these objects are written into the collected
687 // integer system data array above, accessed through one of the preceding Hybrid POINTER-kind
688 // objects.
689
690 // Atom and residue details
691 Hybrid<int2> residue_limits;
694 Hybrid<int> atom_struc_numbers;
696 Hybrid<int> residue_numbers;
697 Hybrid<int2> molecule_limits;
700 Hybrid<int> atomic_numbers;
701 Hybrid<int> mobile_atoms;
702 Hybrid<int> molecule_membership;
704 Hybrid<int> molecule_contents;
707 Hybrid<double> atomic_charges;
708 Hybrid<double> atomic_masses;
709 Hybrid<double> inverse_atomic_masses;
710 Hybrid<float> sp_atomic_charges;
711 Hybrid<float> sp_atomic_masses;
712 Hybrid<float> sp_inverse_atomic_masses;
713 Hybrid<char4> atom_names;
714 Hybrid<char4> atom_types;
715 Hybrid<char4> residue_names;
716 Hybrid<int> chem_int_data;
717 Hybrid<int2> chem_int2_data;
718 Hybrid<double> chem_double_data;
719 Hybrid<float> chem_float_data;
720 Hybrid<char4> chem_char4_data;
721
722 // Valence parameter maps and bounds, showing how the parameters of each unique topology map to
723 // those of the consensus tables in the synthesis.
724 Hybrid<int> ubrd_param_map;
725 Hybrid<int2> ubrd_param_map_bounds;
730 Hybrid<int> cimp_param_map;
731 Hybrid<int2> cimp_param_map_bounds;
732 Hybrid<int> cmap_param_map;
733 Hybrid<int2> cmap_param_map_bounds;
734 Hybrid<int> bond_param_map;
735 Hybrid<int2> bond_param_map_bounds;
736 Hybrid<int> angl_param_map;
737 Hybrid<int2> angl_param_map_bounds;
738 Hybrid<int> dihe_param_map;
739 Hybrid<int2> dihe_param_map_bounds;
740 Hybrid<int> attn14_param_map;
741 Hybrid<int2> attn14_param_map_bounds;
742 Hybrid<int> vste_param_map;
743 Hybrid<int2> vste_param_map_bounds;
744 Hybrid<int> sett_param_map;
745 Hybrid<int2> sett_param_map_bounds;
746 Hybrid<int> cnst_param_map;
747 Hybrid<int2> cnst_param_map_bounds;
749
750 // CHARMM and basic force field valence term details. Each of these objects points into one of
751 // the data arrays at the bottom of the section.
752 Hybrid<double> ubrd_stiffnesses;
753 Hybrid<double> ubrd_equilibria;
754 Hybrid<double> cimp_stiffnesses;
755 Hybrid<double> cimp_phase_angles;
756 Hybrid<int> cmap_surface_dimensions;
757 Hybrid<int> cmap_surface_bounds;
758 Hybrid<int> cmap_patch_bounds;
759 Hybrid<double> cmap_surfaces;
761 Hybrid<double> cmap_patches;
764 Hybrid<float> sp_ubrd_stiffnesses;
766 Hybrid<float> sp_ubrd_equilibria;
768 Hybrid<float> sp_cimp_stiffnesses;
769 Hybrid<float> sp_cimp_phase_angles;
770 Hybrid<float> sp_cmap_surfaces;
772 Hybrid<float> sp_cmap_patches;
775
776 // Basic force field valence term details, consensus tables form all topologies
777 Hybrid<double> bond_stiffnesses;
778 Hybrid<double> bond_equilibria;
779 Hybrid<double> angl_stiffnesses;
780 Hybrid<double> angl_equilibria;
781 Hybrid<double> dihe_amplitudes;
782 Hybrid<double> dihe_periodicities;
783 Hybrid<double> dihe_phase_angles;
784 Hybrid<double> attn14_elec_factors;
786 Hybrid<double> attn14_vdw_factors;
788 Hybrid<float> sp_bond_stiffnesses;
789 Hybrid<float> sp_bond_equilibria;
790 Hybrid<float> sp_angl_stiffnesses;
791 Hybrid<float> sp_angl_equilibria;
792 Hybrid<float> sp_dihe_amplitudes;
793 Hybrid<float> sp_dihe_periodicities;
794 Hybrid<float> sp_dihe_phase_angles;
795 Hybrid<float> sp_attn14_elec_factors;
797 Hybrid<float> sp_attn14_vdw_factors;
799 Hybrid<double> valparam_double_data;
800 Hybrid<float> valparam_float_data;
801 Hybrid<int> valparam_int_data;
803 Hybrid<int2> valparam_int2_data;
806
807 // Valence term indexing arrays, all of them indexing atoms in the synthesis list, updated from
808 // the indexing in their original topologies and parameters in the condensed tables of the
809 // synthesis. All of the following arrays are POINTER-kind objects targeting valence_int_data.
810 Hybrid<int> ubrd_i_atoms;
811 Hybrid<int> ubrd_k_atoms;
812 Hybrid<int> ubrd_param_idx;
813 Hybrid<int> cimp_i_atoms;
814 Hybrid<int> cimp_j_atoms;
815 Hybrid<int> cimp_k_atoms;
816 Hybrid<int> cimp_l_atoms;
817 Hybrid<int> cimp_param_idx;
818 Hybrid<int> cmap_i_atoms;
819 Hybrid<int> cmap_j_atoms;
820 Hybrid<int> cmap_k_atoms;
821 Hybrid<int> cmap_l_atoms;
822 Hybrid<int> cmap_m_atoms;
823 Hybrid<int> cmap_param_idx;
824 Hybrid<int> bond_i_atoms;
825 Hybrid<int> bond_j_atoms;
826 Hybrid<int> bond_param_idx;
827 Hybrid<int> angl_i_atoms;
828 Hybrid<int> angl_j_atoms;
829 Hybrid<int> angl_k_atoms;
830 Hybrid<int> angl_param_idx;
831 Hybrid<int> dihe_i_atoms;
832 Hybrid<int> dihe_j_atoms;
833 Hybrid<int> dihe_k_atoms;
834 Hybrid<int> dihe_l_atoms;
835 Hybrid<int> dihe_param_idx;
836 Hybrid<int> valence_int_data;
838
839 // Non-bonded parameter indexing and van-der Waals tables
840 Hybrid<int> charge_indices;
844 Hybrid<int> lennard_jones_indices;
847 Hybrid<double> charge_parameters;
851 Hybrid<double2> lennard_jones_ab_coeff;
859 Hybrid<double> lennard_jones_c_coeff;
861 Hybrid<double> lennard_jones_14_a_coeff;
863 Hybrid<double> lennard_jones_14_b_coeff;
865 Hybrid<double> lennard_jones_14_c_coeff;
867 Hybrid<double> lennard_jones_sigma;
869 Hybrid<double> lennard_jones_14_sigma;
871 Hybrid<float> sp_charge_parameters;
873 Hybrid<float2> sp_lennard_jones_ab_coeff;
876 Hybrid<float> sp_lennard_jones_c_coeff;
879 Hybrid<float> sp_lennard_jones_14_a_coeff;
882 Hybrid<float> sp_lennard_jones_14_b_coeff;
885 Hybrid<float> sp_lennard_jones_14_c_coeff;
888 Hybrid<float> sp_lennard_jones_sigma;
890 Hybrid<float> sp_lennard_jones_14_sigma;
893
894 // Implicit solvent model parameters. Like the atomic partial charges arrays, most of these will
895 // target chem_double_data, chem_float_data, and (in the case of the indices) chem_int_data,
896 // for convenience. The neck tables are their own ARRAY-type objects to expedite applying a
897 // new implicit solvent model to an existing synthesis.
898 int neck_table_size;
902 Hybrid<int> neck_gb_indices;
904 Hybrid<double> atomic_pb_radii;
905 Hybrid<double> gb_screening_factors;
906 Hybrid<double> gb_alpha_parameters;
909 Hybrid<double> gb_beta_parameters;
910 Hybrid<double> gb_gamma_parameters;
911 Hybrid<double2> neck_limit_tables;
918 Hybrid<float> sp_atomic_pb_radii;
919 Hybrid<float> sp_gb_screening_factors;
920 Hybrid<float> sp_gb_alpha_parameters;
921 Hybrid<float> sp_gb_beta_parameters;
922 Hybrid<float> sp_gb_gamma_parameters;
923 Hybrid<float2> sp_neck_limit_tables;
924
925 // NMR restraint terms and details: these function exactly like other parameter sets and are
926 // indexed by lists of atoms in the bond work units arrays. They can be included in the
927 // synthesis of AtomGraphs due to their nature as potential terms, whereas the original
928 // topologies had to be read from files that did not contain such terms.
929 Hybrid<int2> rposn_step_bounds;
931 Hybrid<int2> rbond_step_bounds;
933 Hybrid<int2> rangl_step_bounds;
935 Hybrid<int2> rdihe_step_bounds;
937 Hybrid<double2> rposn_init_k;
940 Hybrid<double2> rposn_final_k;
942 Hybrid<double4> rposn_init_r;
945 Hybrid<double4> rposn_final_r;
947 Hybrid<double2> rposn_init_xy;
950 Hybrid<double> rposn_init_z;
953 Hybrid<double2> rposn_final_xy;
956 Hybrid<double> rposn_final_z;
959 Hybrid<double2> rbond_init_k;
962 Hybrid<double2> rbond_final_k;
964 Hybrid<double4> rbond_init_r;
967 Hybrid<double4> rbond_final_r;
969 Hybrid<double2> rangl_init_k;
971 Hybrid<double2> rangl_final_k;
973 Hybrid<double4> rangl_init_r;
976 Hybrid<double4> rangl_final_r;
978 Hybrid<double2> rdihe_init_k;
981 Hybrid<double2> rdihe_final_k;
983 Hybrid<double4> rdihe_init_r;
986 Hybrid<double4> rdihe_final_r;
988 Hybrid<float2> sp_rposn_init_k;
991 Hybrid<float2> sp_rposn_final_k;
993 Hybrid<float4> sp_rposn_init_r;
996 Hybrid<float4> sp_rposn_final_r;
998 Hybrid<float2> sp_rposn_init_xy;
1001 Hybrid<float> sp_rposn_init_z;
1004 Hybrid<float2> sp_rposn_final_xy;
1007 Hybrid<float> sp_rposn_final_z;
1010 Hybrid<float2> sp_rbond_init_k;
1013 Hybrid<float2> sp_rbond_final_k;
1015 Hybrid<float4> sp_rbond_init_r;
1018 Hybrid<float4> sp_rbond_final_r;
1020 Hybrid<float2> sp_rangl_init_k;
1022 Hybrid<float2> sp_rangl_final_k;
1024 Hybrid<float4> sp_rangl_init_r;
1027 Hybrid<float4> sp_rangl_final_r;
1029 Hybrid<float2> sp_rdihe_init_k;
1032 Hybrid<float2> sp_rdihe_final_k;
1034 Hybrid<float4> sp_rdihe_init_r;
1037 Hybrid<float4> sp_rdihe_final_r;
1039 Hybrid<double> nmr_double_data;
1043 Hybrid<double2> nmr_double2_data;
1047 Hybrid<double4> nmr_double4_data;
1051 Hybrid<float> nmr_float_data;
1055 Hybrid<float2> nmr_float2_data;
1059 Hybrid<float4> nmr_float4_data;
1063
1064 // Restraint indexing arrays, offering atom indices in the concatenated atom and restraint
1065 // parameter arrays of the synthesis. The Hybrid objects in this section are POINTER-kind
1066 // objects targeting the nmr_int_data array, much like the parameters in the preceding section
1067 // target nmr_[type]_data.
1068 Hybrid<int> rposn_atoms;
1069 Hybrid<int> rposn_kr_param_idx;
1071 Hybrid<int> rposn_xyz_param_idx;
1073 Hybrid<int> rbond_i_atoms;
1074 Hybrid<int> rbond_j_atoms;
1075 Hybrid<int> rbond_param_idx;
1077 Hybrid<int> rangl_i_atoms;
1078 Hybrid<int> rangl_j_atoms;
1079 Hybrid<int> rangl_k_atoms;
1080 Hybrid<int> rangl_param_idx;
1082 Hybrid<int> rdihe_i_atoms;
1083 Hybrid<int> rdihe_j_atoms;
1084 Hybrid<int> rdihe_k_atoms;
1085 Hybrid<int> rdihe_l_atoms;
1086 Hybrid<int> rdihe_param_idx;
1088 Hybrid<int> rposn_kr_param_map;
1089 Hybrid<int> rposn_xyz_param_map;
1090 Hybrid<int2> rposn_param_map_bounds;
1094 Hybrid<int> rbond_param_map;
1095 Hybrid<int2> rbond_param_map_bounds;
1096 Hybrid<int> rangl_param_map;
1097 Hybrid<int2> rangl_param_map_bounds;
1098 Hybrid<int> rdihe_param_map;
1099 Hybrid<int2> rdihe_param_map_bounds;
1100 Hybrid<int> nmr_int_data;
1101 Hybrid<int2> nmr_int2_data;
1102
1103 // Virtual site details: virtual sites are considered parameters, which like valence terms apply
1104 // to a small group of atoms and index into a table of parameters including the frame type and
1105 // up to three dimensional measurements.
1106 Hybrid<double4> virtual_site_parameters;
1114 Hybrid<float4> sp_virtual_site_parameters;
1116 Hybrid<int> virtual_site_atoms;
1125 Hybrid<int> virtual_site_frame1_atoms;
1127 Hybrid<int> virtual_site_frame2_atoms;
1128 Hybrid<int> virtual_site_frame3_atoms;
1129 Hybrid<int> virtual_site_frame4_atoms;
1130 Hybrid<int> virtual_site_parameter_indices;
1133 Hybrid<int> vsite_int_data;
1136
1137 // SETTLE constraint group parameter sets and atom indices
1138 Hybrid<double4> settle_group_geometry;
1142 Hybrid<double4> settle_group_masses;
1147 Hybrid<float4> sp_settle_group_geometry;
1150 Hybrid<float4> sp_settle_group_masses;
1155 Hybrid<int4> settle_group_indexing;
1158
1159 // Constraint group parameter sets and atom indices
1160 Hybrid<int> constraint_group_indices;
1163 Hybrid<int2> constraint_group_bounds;
1171 Hybrid<int> constraint_group_param_idx;
1174 Hybrid<int> constraint_param_bounds;
1176 Hybrid<double2> constraint_group_params;
1180 Hybrid<float2> sp_constraint_group_params;
1183
1184 // Valence work units (VWUs): work units providing instruction sets for the GPU to operate on a
1185 // continuous, non-rearranged list of atoms and implement all valence terms. Each VWU pertains
1186 // to one and only one individual topology from the list and each topology will have at least one
1187 // VWU, or more depending on what is needed to cover all of its atoms. The VWU imports a list of
1188 // atoms, then computes a series of bond, angle, dihedral, and other force field terms based on
1189 // the cached atoms. Forces are accumulated on all atoms in __shared__ memory and then
1190 // contributed back to the global arrays.
1191 int total_valence_work_units;
1192 int2 valence_work_unit_size;
1196
1199 ValenceKernelSize valence_thread_block_size;
1200
1205 Hybrid<int2> vwu_instruction_sets;
1206
1210 Hybrid<int> vwu_import_lists;
1211
1217 Hybrid<uint2> vwu_manipulation_masks;
1218
1224 Hybrid<uint2> cbnd_instructions;
1225
1229 Hybrid<uint2> angl_instructions;
1230
1250 Hybrid<uint2> cdhe_instructions;
1251
1256 Hybrid<uint> cdhe_overtones;
1257
1262 Hybrid<uint2> cmap_instructions;
1263
1268 Hybrid<uint> infr14_instructions;
1269
1276 Hybrid<uint2> rposn_instructions;
1277
1282 Hybrid<uint2> rbond_instructions;
1283
1288 Hybrid<uint2> rangl_instructions;
1289
1296 Hybrid<uint2> rdihe_instructions;
1297
1305 Hybrid<uint2> vste_instructions;
1306
1311 Hybrid<uint2> sett_instructions;
1312
1321 Hybrid<uint2> cnst_instructions;
1322
1323 // Instructions for evaluating energy terms--these are stored as extra arrays, as they will only
1324 // be called when evaluating energies, and offer a compact format in which each interaction in
1325 // the list is a bit in a bitstring, (1) indicating that the work unit should contribute the
1326 // interaction to its sum, (0) indicating that the work unit should omit the interaction.
1327 // All of these are POINTER-kind objects targeting insr_uint_data.
1328 Hybrid<uint> accumulate_cbnd_energy;
1330 Hybrid<uint> accumulate_angl_energy;
1331 Hybrid<uint> accumulate_cdhe_energy;
1333 Hybrid<uint> accumulate_cmap_energy;
1334 Hybrid<uint> accumulate_infr14_energy;
1335 Hybrid<uint> accumulate_rposn_energy;
1339 Hybrid<uint> accumulate_rbond_energy;
1340 Hybrid<uint> accumulate_rangl_energy;
1341 Hybrid<uint> accumulate_rdihe_energy;
1342
1344 Hybrid<uint> insr_uint_data;
1345
1347 Hybrid<uint2> insr_uint2_data;
1348
1349 // Non-bonded work units (NBWUs): the synthesis also stores information on how to carry out
1350 // non-bonded tiles (or domain decomposition workloads, depending on the nature of the boundary
1351 // conditions). These are less detailed, overall, than the valence work units.
1352 int total_nonbonded_work_units;
1355 NbwuKind nonbonded_work_type;
1358
1380 Hybrid<int> nonbonded_abstracts;
1381
1386 Hybrid<uint2> nbwu_instructions;
1387
1388 // Reduction work units. Various processes in molecular mechanics require system-wide
1389 // accumulations of particular quantities, and a collection of any number of systems of any size
1390 // requires a sophisticated system to make that happen efficiently. The reduction work units
1391 // are an ecumenical solution to all such reductions, but they may need to combine with different
1392 // substrates to do a specific process.
1393 int total_reduction_work_units;
1394 RdwuPerSystem rdwu_per_system;
1397
1401 Hybrid<int> reduction_abstracts;
1402
1404 StopWatch* timer;
1405
1407 void rebasePointers();
1408
1417 std::vector<int> createDummyRestraints(const std::vector<int> &restraint_indices_in,
1418 const std::vector<int> &topology_indices_in);
1419
1427 std::vector<int> checkTopologyList(const std::vector<int> &topology_indices_in);
1428
1439 std::vector<int> checkRestraintList(const std::vector<int> &restraint_indices_in,
1440 const std::vector<int> &topology_indices_in,
1441 const std::vector<int> &topology_index_rebase);
1442
1444 void checkCommonSettings();
1445
1461 void buildAtomAndTermArrays(const std::vector<int> &topology_indices_in,
1462 const std::vector<int> &topology_index_rebase,
1463 const std::vector<int> &new_restraint_indices,
1464 const std::vector<int> &restraint_index_rebase);
1465
1470 void condenseParameterTables();
1471
1477 void extendLJMatrices();
1478
1504 int mapUniqueRestraintKRSeries(int order, const std::vector<int> &network_table_offsets,
1505 std::vector<int> *synthesis_index,
1506 std::vector<int2> *filtered_step_bounds,
1507 std::vector<double2> *filtered_init_keq,
1508 std::vector<double2> *filtered_finl_keq,
1509 std::vector<double4> *filtered_init_r,
1510 std::vector<double4> *filtered_finl_r);
1511
1515 void condenseRestraintNetworks();
1516
1527 int setVwuAbstractLimits(int item_counter, int vwu_counter, VwuAbstractMap slot,
1528 int item_quantity);
1529
1538 void loadValenceWorkUnits(int2 vwu_atom_limits = { maximum_valence_work_unit_atoms,
1539 maximum_valence_work_unit_atoms });
1540
1544 void loadReductionWorkUnits(const GpuDetails &gpu = null_gpu);
1545
1551 void importImplicitSolventAtomParameters(int system_index);
1552
1561 void setImplicitSolventNeckParameters(const NeckGeneralizedBornTable &ngb_tab);
1562 void setImplicitSolventNeckParameters();
1564};
1565
1566} // namespace synthesis
1567} // namespace stormm
1568
1569#include "atomgraph_synthesis.tpp"
1570
1571#endif
AtomGraphSynthesis(const std::vector< AtomGraph * > &topologies_in, const std::vector< RestraintApparatus * > &restraints_in, const std::vector< int > &topology_indices_in, const std::vector< int > &restraint_indices_in, ExceptionResponse policy_in=ExceptionResponse::WARN, const GpuDetails &gpu=null_gpu, StopWatch *timer_in=nullptr)
The constructor takes a series of topologies and NMR restraints. The NMR restraints point to specific...
Definition ag_synthesis_constructors.cpp:11
const AtomGraph * getSystemTopologyPointer(int system_index) const
Get a topology pointer for a specific system contained within the synthesis.
Definition ag_synthesis_mechanics.cpp:2768
Pertinent aspects of one particular GPU. Condensing the data for each GPU in this manner helps to ens...
Definition gpu_details.h:27
An evolution of GpuBuffer in pmemd.cuda, the Composite array has elements that are accessible from ei...
Definition hybrid.h:202
Object to hold a complex array of constants referenced by various GB calculations using the "neck" fo...
Definition generalized_born.h:151
int getBondParameterCount() const
Get the number of unique harmonic bond parameter sets.
Definition ag_synthesis_mechanics.cpp:2889
int getBondTermCount() const
Get the total number of bond terms across all systems, including replicas.
Definition ag_synthesis_mechanics.cpp:2849
std::vector< int > getTopologyIndices() const
Get a const reference to the list of all topology indices, indicating the topologies describing each ...
Definition ag_synthesis_mechanics.cpp:2793
int getLJTypeCount() const
Get the number of unique atom types (a parameter, not an extensive quantity dependent on the number o...
Definition ag_synthesis_mechanics.cpp:2879
std::vector< AtomicRadiusSet > getPBRadiiSet() const
Get the name of the PB radii set for one or more systems.
Definition ag_synthesis_mechanics.cpp:2964
SyRestraintKit< double, double2, double4 > getDoublePrecisionRestraintKit(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get a minimal kit with double-precision parameter detail for computing valence interactions for all s...
Definition ag_synthesis_mechanics.cpp:3109
UnitCellType getUnitCellType() const
Get the unit cell type that will be taken for all systems (TRICLINIC subsumes ORTHORHOMBIC in a sort ...
Definition ag_synthesis_mechanics.cpp:2939
int getRandomCacheDepth() const
Get the depth of the random numbers cache that the non-bonded work units of this topology synthesis w...
Definition ag_synthesis_mechanics.cpp:3024
SyNonbondedKit< double, double2 > getDoublePrecisionNonbondedKit(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get a minimal kit with double-precision real numbers for computing non-bonded interactions for all sy...
Definition ag_synthesis_mechanics.cpp:3159
double getSaltConcentration() const
Get the salt concentration (supporting the implicit solvent model) for all systems.
Definition ag_synthesis_mechanics.cpp:2954
int getAtomCount() const
Get the number of atoms in one or more systems.
Definition ag_synthesis_mechanics.cpp:2813
const std::vector< AtomGraph * > & getUniqueTopologies() const
Get a const reference to the list of all pointers for unique topologies making up this synthesis.
Definition ag_synthesis_mechanics.cpp:2788
int getValenceWorkUnitCount() const
Get the overall number of valence work units needed to account for interactions in all systems.
Definition ag_synthesis_mechanics.cpp:2994
const Hybrid< int > & getSystemAtomCounts() const
Get the sizes of all individual systems as a const reference to the Hybrid array member variable.
Definition ag_synthesis_mechanics.cpp:2919
int getDihedralParameterCount() const
Get the number of unique cosine-based dihedral parameter sets.
Definition ag_synthesis_mechanics.cpp:2899
void setImplicitSolventModel()
Apply an implicit solvent model to the synthesis. Any mode of operation other than taking the origina...
Definition ag_synthesis_mechanics.cpp:3419
AtomGraphSynthesis(const std::vector< AtomGraph * > &topologies_in, const std::vector< RestraintApparatus * > &restraints_in, const std::vector< int > &topology_indices_in, const std::vector< int > &restraint_indices_in, ExceptionResponse policy_in=ExceptionResponse::WARN, const GpuDetails &gpu=null_gpu, StopWatch *timer_in=nullptr)
The constructor takes a series of topologies and NMR restraints. The NMR restraints point to specific...
Definition ag_synthesis_constructors.cpp:11
const Hybrid< int > & getConstrainedDegreesOfFreedom() const
Get the numbers of constrained degrees of freedom for all individual systems in the.
Definition ag_synthesis_mechanics.cpp:2934
int getUniqueTopologyCount() const
Get the number of unique topologies described by the synthesis.
Definition ag_synthesis_mechanics.cpp:2763
ImplicitSolventModel getImplicitSolventModel() const
Get the implicit solvent model in use across all systems.
Definition ag_synthesis_mechanics.cpp:2944
int getAtomOffset(int system_index) const
Get the starting point for atoms of a specific system in the lineup of all topologies.
Definition ag_synthesis_mechanics.cpp:2834
NbwuKind getNonbondedWorkType() const
Get the type of non-bonded work required by systems in this synthesis.
Definition ag_synthesis_mechanics.cpp:3014
SyValenceKit< float > getSinglePrecisionValenceKit(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get a minimal kit with single-precision parameter detail for computing valence interactions for all s...
Definition ag_synthesis_mechanics.cpp:3082
const AtomGraphSynthesis * getSelfPointer() const
Get a const pointer to the object itself, useful for retrieving a valid pointer when the object is av...
Definition ag_synthesis_mechanics.cpp:3230
int getDihedralTermCount() const
Get the total number of cosine-based dihedral terms across all systems and replicas.
Definition ag_synthesis_mechanics.cpp:2859
int getCmapSurfaceCount() const
Get the number of unique CMAP surfaces.
Definition ag_synthesis_mechanics.cpp:2914
SyAtomUpdateKit< double, double2, double4 > getDoublePrecisionAtomUpdateKit(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get a minimal kit with double-precision real numbers for updating atom and virtual site positions bas...
Definition ag_synthesis_mechanics.cpp:3201
int getChargeTypeCount() const
Get the number of unique charge parameters.
Definition ag_synthesis_mechanics.cpp:2884
double getCoulombConstant() const
Get the fundamental Coulomb constant defining the electrostatics of all systems.
Definition ag_synthesis_mechanics.cpp:2959
const Hybrid< int > & getSystemAtomOffsets() const
Get the starting locations of all individual systems as a const reference to the Hybrid array member ...
Definition ag_synthesis_mechanics.cpp:2924
std::vector< T > getAtomicMasses(HybridTargetLevel tier) const
Get the masses (or inverse masses) of atoms in the synthesis.
std::vector< T > getPartialCharges(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get partial charges stored within the synthesis.
int getUreyBradleyTermCount() const
Get the total number of Urey-Bradley terms across all systems and replicas.
Definition ag_synthesis_mechanics.cpp:2864
int getNonbondedWorkUnitCount() const
Get the number of non-bonded work units serving systems in this synthesis.
Definition ag_synthesis_mechanics.cpp:3019
SyNonbondedKit< float, float2 > getSinglePrecisionNonbondedKit(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get a minimal kit with single-precision real numbers for computing non-bonded interactions for all sy...
Definition ag_synthesis_mechanics.cpp:3179
int getPaddedAtomCount() const
Get the entire padded number of atoms covering all systems.
Definition ag_synthesis_mechanics.cpp:2828
int getSystemCount() const
Get the number of systems that this synthesis describes.
Definition ag_synthesis_mechanics.cpp:2808
SyAtomUpdateKit< float, float2, float4 > getSinglePrecisionAtomUpdateKit(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get a minimal kit with single-precision real numbers for updating atom and virtual site positions bas...
Definition ag_synthesis_mechanics.cpp:3215
RestraintApparatus * getSystemRestraintPointer(int system_index) const
Get a restraint apparatus pointer for a sepcific system contained within the synthesis.
Definition ag_synthesis_mechanics.cpp:2798
double getDielectricConstant() const
Get the dielectric constant (supporting the implicit solvent model) for all systems.
Definition ag_synthesis_mechanics.cpp:2949
const Hybrid< int > & getDegreesOfFreedom() const
Get the numbers of unconstrained degrees of freedom for all individual systems in the.
Definition ag_synthesis_mechanics.cpp:2929
int getAngleParameterCount() const
Get the number of unique harmonic bond angle parameter sets.
Definition ag_synthesis_mechanics.cpp:2894
int getCharmmImproperTermCount() const
Get the total number of CHARMM improper terms across all systems and replicas.
Definition ag_synthesis_mechanics.cpp:2869
ValenceKernelSize getValenceThreadBlockSize() const
Get the necessary thread block size for evaluating the valence work units.
Definition ag_synthesis_mechanics.cpp:3004
int getCharmmImproperParameterCount() const
Get the number of unique CHARMM improper parameter sets.
Definition ag_synthesis_mechanics.cpp:2909
int getUreyBradleyParameterCount() const
Get the number of unique Urey-Bradley harmonic angle parameter sets.
Definition ag_synthesis_mechanics.cpp:2904
const Hybrid< int > & getReductionWorkUnitAbstracts() const
Get the reduction work unit abstracts (this is not required for valence or non-bonded work units as t...
Definition ag_synthesis_mechanics.cpp:3049
SyValenceKit< double > getDoublePrecisionValenceKit(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get a minimal kit with double-precision parameter detail for computing valence interactions for all s...
Definition ag_synthesis_mechanics.cpp:3055
SyRestraintKit< float, float2, float4 > getSinglePrecisionRestraintKit(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get a minimal kit with single-precision parameter detail for computing valence interactions for all s...
Definition ag_synthesis_mechanics.cpp:3134
void loadNonbondedWorkUnits(const StaticExclusionMaskSynthesis &poly_se, InitializationTask init_request=InitializationTask::NONE, int random_cache_depth=0, const GpuDetails &gpu=null_gpu)
Construct non-bonded work units for all unique topologies (there are no restraints for non-bonded int...
Definition ag_synthesis_mechanics.cpp:2668
int getCmapTermCount() const
Get the total number of CMAP terms across all systems and replicas.
Definition ag_synthesis_mechanics.cpp:2874
RdwuPerSystem getRdwuPerSystem() const
Get a qualitative assessment of the number of reduction work units assigned to any one system.
Definition ag_synthesis_mechanics.cpp:3044
int2 getValenceWorkUnitSize() const
Get the maximum size of valence work units, whether for organic molecules or solvent (rigid water) re...
Definition ag_synthesis_mechanics.cpp:2999
int getVirtualSiteCount() const
Get the total number of virtual sites across all systems, including replicas.
Definition ag_synthesis_mechanics.cpp:2844
int getReductionWorkUnitCount() const
Get the number of reduction work units spanning all systems.
Definition ag_synthesis_mechanics.cpp:3039
const Hybrid< int2 > & getValenceWorkUnitAbstracts() const
Get the abstracts (condensed lists of import and instruction set limits) for the valence work units s...
Definition ag_synthesis_mechanics.cpp:3009
int getAngleTermCount() const
Get the total number of bond angle terms across all systems, including replicas.
Definition ag_synthesis_mechanics.cpp:2854
A struct to hold information relating to an Amber topology. This struct's member functions are limite...
Definition atomgraph.h:50
A collection of all restraints pertaining to a specific topology for the purposes of one simulation,...
Definition restraint_apparatus.h:109
An exclusion mask object for a compilation of systems. All systems are represented in full detail,...
Definition static_mask_synthesis.h:54
Object for managing calls to the C-standard function gettimeofday(), calculating deltas and categoriz...
Definition stopwatch.h:23
Definition stormm_vector_types.h:141
Definition stormm_vector_types.h:22
Abstract for the NeckGeneralizedBornTable object, in single- or double-precision.
Definition generalized_born.h:130
Collect the virtual site details and constraint parameters of the topology synthesis into a single ab...
Definition synthesis_abstracts.h:257
Collect the critical non-bonded parameters and masking information for work unit-based evaluation of ...
Definition synthesis_abstracts.h:181
Collect the critical restraint parameters and masking information for work unit-based evaluation of t...
Definition synthesis_abstracts.h:100
Collect the critical valence parameters and indexing information for work unit-based evaluation of th...
Definition synthesis_abstracts.h:19