STORMM Source Documentation
Loading...
Searching...
No Matches
background_mesh.h
1// -*-c++-*-
2#ifndef STORMM_PUREMESH_H
3#define STORMM_PUREMESH_H
4
5#include <algorithm>
6#include <string>
7#include <vector>
8#include "copyright.h"
9#include "Accelerator/gpu_details.h"
10#include "Accelerator/hybrid.h"
11#include "Accelerator/mesh_kernel_manager.h"
12#include "Constants/behavior.h"
13#include "Constants/fixed_precision.h"
14#include "Constants/scaling.h"
15#include "Constants/symbol_values.h"
16#include "DataTypes/common_types.h"
17#include "DataTypes/stormm_vector_types.h"
18#include "Math/matrix_ops.h"
19#include "Math/one_dimensional_splines.h"
20#include "Math/radial_derivatives.h"
21#include "Math/rounding.h"
22#include "Math/tricubic_cell.h"
23#include "Math/vector_ops.h"
24#include "Namelists/nml_mesh.h"
25#include "Numerics/split_fixed_precision.h"
26#include "Parsing/polynumeric.h"
27#include "Potential/energy_enumerators.h"
28#include "Reporting/reporting_enumerators.h"
29#include "Topology/atomgraph.h"
30#include "Trajectory/coordinateframe.h"
31#include "Trajectory/coordinate_copy.h"
32#include "Trajectory/coordinate_series.h"
33#include "Trajectory/phasespace.h"
34#include "local_arrangement.h"
35#include "mesh_forcefield.h"
36#include "mesh_foundation.h"
37#include "mesh_parameters.h"
38#include "mesh_rulers.h"
39#include "structure_enumerators.h"
40
41namespace stormm {
42namespace structure {
43
44using card::GpuDetails;
45using card::Hybrid;
46using card::HybridKind;
47using card::HybridTargetLevel;
48using card::MeshKlManager;
49using constants::CartesianDimension;
50using constants::PrecisionModel;
51using constants::UnitCellAxis;
52using data_types::getStormmScalarTypeName;
53using data_types::isScalarType;
54using data_types::isFloatingPointScalarType;
55using data_types::isSignedIntegralScalarType;
56using energy::NonbondedPotential;
57using energy::cubicSoftCore;
58using energy::quinticSoftCore;
59using energy::VdwCombiningRule;
60using parse::NumberFormat;
61using namelist::default_mesh_density_averaging_order;
62using namelist::default_mesh_elec_damping_range;
63using namelist::default_mesh_vdw_damping_ratio;
64using numerics::default_globalpos_scale_bits;
65using numerics::fixedPrecisionGrid;
66using review::GridFileSyntax;
67using review::getEnumerationName;
68using stmath::addScalarToVector;
69using stmath::convertData;
70using stmath::elementwiseMultiply;
71using stmath::evaluateCubicSpline;
72using stmath::evaluateQuinticSpline;
73using stmath::fvStencilCoordinates;
74using stmath::hessianNormalWidths;
75using stmath::invertSquareMatrix;
76using stmath::maxAbsoluteDifference;
77using stmath::radialFirstDerivative;
78using stmath::radialSecondDerivative;
79using stmath::radialThirdDerivative;
80using stmath::roundUp;
81using stmath::TricubicCell;
82using stmath::TricubicStencil;
83using topology::AtomGraph;
84using topology::NonbondedKit;
85using trajectory::coordCopy;
86using trajectory::CoordinateFrame;
87using trajectory::CoordinateFrameReader;
88using trajectory::CoordinateSeries;
89using trajectory::CoordinateSeriesReader;
90using trajectory::PhaseSpace;
91
93template <typename Tdata> struct BackgroundMeshWriter {
94
96 BackgroundMeshWriter(const MeshParamKit &measurements, GridDetail kind_in,
97 NonbondedPotential field_in, const MeshRulerKit &rulers,
98 Tdata* coeffs_in, double coeff_scale_in, double probe_radius_in,
99 double well_depth_in, double occ_cost_in, const MeshBasicsKit &mbss_in);
100
109
111 const GridDetail kind;
112 const NonbondedPotential field;
114
115 // The content of the mesh is what is actually writeable.
116 Tdata* coeffs;
120 const double coeff_scale;
126 const float coeff_scale_f;
127
128 // The following data is critical to producing the mesh, and can be necessary for interpreting
129 // the potentials read from it.
130 const double probe_radius;
133 const double well_depth;
135 const double occ_cost;
139};
140
142template <typename Tdata> struct BackgroundMeshReader {
143
145 BackgroundMeshReader(const MeshParamKit &dims_in, GridDetail kind_in,
146 NonbondedPotential field_in, const MeshRulerKit &rulers,
147 const Tdata* coeffs_in, double coeff_scale_in, double probe_radius_in,
148 double well_depth_in, double occ_cost_in, const MeshBasicsKit &mbss_in);
149
158
160 const GridDetail kind;
161 const NonbondedPotential field;
163 const Tdata* coeffs;
168 const double coeff_scale;
175 const float coeff_scale_f;
176
177 // The following data can be necessary for interpreting potentials read from the mesh.
178 const double probe_radius;
181 const double well_depth;
183 const double occ_cost;
187};
188
195template <typename T> class BackgroundMesh {
196public:
197
238 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in,
239 const MeshParameters &measurements_in = MeshParameters(),
240 const PrecisionModel build_precision_in = PrecisionModel::SINGLE);
241
242 // General meshes with all available details
243 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in,
244 const AtomGraph *ag_in, const CoordinateFrame *cf_in,
245 double probe_radius_in, double well_depth_in, VdwCombiningRule mixing_protocol_in,
246 const MeshParameters &measurements_in = MeshParameters(),
247 const std::vector<double> &probe_sigma = {},
248 const std::vector<double> &probe_epsilon = {},
249 double clash_distance_in = default_mesh_elec_damping_range,
250 double clash_ratio_in = default_mesh_vdw_damping_ratio,
251 PrecisionModel prec = PrecisionModel::SINGLE,
252 const MeshKlManager &launcher = MeshKlManager(),
253 HybridTargetLevel availability = HybridTargetLevel::HOST);
254
255 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in,
256 const AtomGraph &ag_in, const CoordinateFrame &cf_in,
257 double probe_radius_in, double well_depth_in, VdwCombiningRule mixing_protocol_in,
258 const MeshParameters &measurements_in = MeshParameters(),
259 const std::vector<double> &probe_sigma = {},
260 const std::vector<double> &probe_epsilon = {},
261 double clash_distance_in = default_mesh_elec_damping_range,
262 double clash_ratio_in = default_mesh_vdw_damping_ratio,
263 PrecisionModel prec = PrecisionModel::SINGLE,
264 const MeshKlManager &launcher = MeshKlManager(),
265 HybridTargetLevel availability = HybridTargetLevel::HOST);
266
267 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, double probe_radius_in,
268 double well_depth_in, VdwCombiningRule mixing_protocol_in, const AtomGraph *ag_in,
269 const CoordinateFrame *cf_in, double buffer, double spacing,
270 int scale_bits_in = default_globalpos_scale_bits,
271 const std::vector<double> &probe_sigma = {},
272 const std::vector<double> &probe_epsilon = {},
273 double clash_distance_in = default_mesh_elec_damping_range,
274 double clash_ratio_in = default_mesh_vdw_damping_ratio,
275 PrecisionModel prec = PrecisionModel::SINGLE,
276 const MeshKlManager &launcher = MeshKlManager(),
277 HybridTargetLevel availability = HybridTargetLevel::HOST);
278
279 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, double probe_radius_in,
280 double well_depth_in, VdwCombiningRule mixing_protocol_in, const AtomGraph &ag_in,
281 const CoordinateFrame &cf_in, double buffer, double spacing,
282 int scale_bits_in = default_globalpos_scale_bits,
283 const std::vector<double> &probe_sigma = {},
284 const std::vector<double> &probe_epsilon = {},
285 double clash_distance_in = default_mesh_elec_damping_range,
286 double clash_ratio_in = default_mesh_vdw_damping_ratio,
287 PrecisionModel prec = PrecisionModel::SINGLE,
288 const MeshKlManager &launcher = MeshKlManager(),
289 HybridTargetLevel availability = HybridTargetLevel::HOST);
290
291 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, double probe_radius_in,
292 double well_depth_in, VdwCombiningRule mixing_protocol_in, const AtomGraph *ag_in,
293 const CoordinateFrame *cf_in, double buffer, const std::vector<double> &spacing,
294 int scale_bits_in = default_globalpos_scale_bits,
295 const std::vector<double> &probe_sigma = {},
296 const std::vector<double> &probe_epsilon = {},
297 double clash_distance_in = default_mesh_elec_damping_range,
298 double clash_ratio_in = default_mesh_vdw_damping_ratio,
299 PrecisionModel prec = PrecisionModel::SINGLE,
300 const MeshKlManager &launcher = MeshKlManager(),
301 HybridTargetLevel availability = HybridTargetLevel::HOST);
302
303 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, double probe_radius_in,
304 double well_depth_in, VdwCombiningRule mixing_protocol_in, const AtomGraph &ag_in,
305 const CoordinateFrame &cf_in, double buffer, const std::vector<double> &spacing,
306 int scale_bits_in = default_globalpos_scale_bits,
307 const std::vector<double> &probe_sigma = {},
308 const std::vector<double> &probe_epsilon = {},
309 double clash_distance_in = default_mesh_elec_damping_range,
310 double clash_ratio_in = default_mesh_vdw_damping_ratio,
311 PrecisionModel prec = PrecisionModel::SINGLE,
312 const MeshKlManager &launcher = MeshKlManager(),
313 HybridTargetLevel availability = HybridTargetLevel::HOST);
314
315 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, double probe_radius_in,
316 double well_depth_in, VdwCombiningRule mixing_protocol_in, const AtomGraph *ag_in,
317 const CoordinateFrame *cf_in, const std::vector<double> &mesh_bounds,
318 double spacing, int scale_bits_in = default_globalpos_scale_bits,
319 const std::vector<double> &probe_sigma = {},
320 const std::vector<double> &probe_epsilon = {},
321 double clash_distance_in = default_mesh_elec_damping_range,
322 double clash_ratio_in = default_mesh_vdw_damping_ratio,
323 PrecisionModel prec = PrecisionModel::SINGLE,
324 const MeshKlManager &launcher = MeshKlManager(),
325 HybridTargetLevel availability = HybridTargetLevel::HOST);
326
327 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, double probe_radius_in,
328 double well_depth_in, VdwCombiningRule mixing_protocol_in, const AtomGraph &ag_in,
329 const CoordinateFrame &cf_in, const std::vector<double> &mesh_bounds,
330 double spacing, int scale_bits_in = default_globalpos_scale_bits,
331 const std::vector<double> &probe_sigma = {},
332 const std::vector<double> &probe_epsilon = {},
333 double clash_distance_in = default_mesh_elec_damping_range,
334 double clash_ratio_in = default_mesh_vdw_damping_ratio,
335 PrecisionModel prec = PrecisionModel::SINGLE,
336 const MeshKlManager &launcher = MeshKlManager(),
337 HybridTargetLevel availability = HybridTargetLevel::HOST);
338
339 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, double probe_radius_in,
340 double well_depth_in, VdwCombiningRule mixing_protocol_in, const AtomGraph *ag_in,
341 const CoordinateFrame *cf_in, const std::vector<double> &mesh_bounds,
342 const std::vector<double> &spacing,
343 int scale_bits_in = default_globalpos_scale_bits,
344 const std::vector<double> &probe_sigma = {},
345 const std::vector<double> &probe_epsilon = {},
346 double clash_distance_in = default_mesh_elec_damping_range,
347 double clash_ratio_in = default_mesh_vdw_damping_ratio,
348 PrecisionModel prec = PrecisionModel::SINGLE,
349 const MeshKlManager &launcher = MeshKlManager(),
350 HybridTargetLevel availability = HybridTargetLevel::HOST);
351
352 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, double probe_radius_in,
353 double well_depth_in, VdwCombiningRule mixing_protocol_in, const AtomGraph &ag_in,
354 const CoordinateFrame &cf_in, const std::vector<double> &mesh_bounds,
355 const std::vector<double> &spacing,
356 int scale_bits_in = default_globalpos_scale_bits,
357 const std::vector<double> &probe_sigma = {},
358 const std::vector<double> &probe_epsilon = {},
359 double clash_distance_in = default_mesh_elec_damping_range,
360 double clash_ratio_in = default_mesh_vdw_damping_ratio,
361 PrecisionModel prec = PrecisionModel::SINGLE,
362 const MeshKlManager &launcher = MeshKlManager(),
363 HybridTargetLevel availability = HybridTargetLevel::HOST);
364
365 // Electrostatics-based meshes
366 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const AtomGraph *ag_in,
367 const CoordinateFrame *cf_in, double buffer, double spacing,
368 int scale_bits_in = default_globalpos_scale_bits,
369 double clash_distance_in = default_mesh_elec_damping_range,
370 PrecisionModel prec = PrecisionModel::SINGLE,
371 const MeshKlManager &launcher = MeshKlManager(),
372 HybridTargetLevel availability = HybridTargetLevel::HOST);
373
374 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const AtomGraph &ag_in,
375 const CoordinateFrame &cf_in, double buffer, double spacing,
376 int scale_bits_in = default_globalpos_scale_bits,
377 double clash_distance_in = default_mesh_elec_damping_range,
378 PrecisionModel prec = PrecisionModel::SINGLE,
379 const MeshKlManager &launcher = MeshKlManager(),
380 HybridTargetLevel availability = HybridTargetLevel::HOST);
381
382 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const AtomGraph *ag_in,
383 const CoordinateFrame *cf_in, double buffer,
384 const std::vector<double> &spacing,
385 int scale_bits_in = default_globalpos_scale_bits,
386 double clash_distance_in = default_mesh_elec_damping_range,
387 PrecisionModel prec = PrecisionModel::SINGLE,
388 const MeshKlManager &launcher = MeshKlManager(),
389 HybridTargetLevel availability = HybridTargetLevel::HOST);
390
391 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const AtomGraph &ag_in,
392 const CoordinateFrame &cf_in, double buffer,
393 const std::vector<double> &spacing,
394 int scale_bits_in = default_globalpos_scale_bits,
395 double clash_distance_in = default_mesh_elec_damping_range,
396 PrecisionModel prec = PrecisionModel::SINGLE,
397 const MeshKlManager &launcher = MeshKlManager(),
398 HybridTargetLevel availability = HybridTargetLevel::HOST);
399
400 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const AtomGraph *ag_in,
401 const CoordinateFrame *cf_in, const std::vector<double> &mesh_bounds,
402 double spacing, int scale_bits_in = default_globalpos_scale_bits,
403 double clash_distance_in = default_mesh_elec_damping_range,
404 PrecisionModel prec = PrecisionModel::SINGLE,
405 const MeshKlManager &launcher = MeshKlManager(),
406 HybridTargetLevel availability = HybridTargetLevel::HOST);
407
408 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const AtomGraph &ag_in,
409 const CoordinateFrame &cf_in, const std::vector<double> &mesh_bounds,
410 double spacing, int scale_bits_in = default_globalpos_scale_bits,
411 double clash_distance_in = default_mesh_elec_damping_range,
412 PrecisionModel prec = PrecisionModel::SINGLE,
413 const MeshKlManager &launcher = MeshKlManager(),
414 HybridTargetLevel availability = HybridTargetLevel::HOST);
415
416 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const AtomGraph *ag_in,
417 const CoordinateFrame *cf_in, const std::vector<double> &mesh_bounds,
418 const std::vector<double> &spacing,
419 int scale_bits_in = default_globalpos_scale_bits,
420 double clash_distance_in = default_mesh_elec_damping_range,
421 PrecisionModel prec = PrecisionModel::SINGLE,
422 const MeshKlManager &launcher = MeshKlManager(),
423 HybridTargetLevel availability = HybridTargetLevel::HOST);
424
425 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const AtomGraph &ag_in,
426 const CoordinateFrame &cf_in, const std::vector<double> &mesh_bounds,
427 const std::vector<double> &spacing,
428 int scale_bits_in = default_globalpos_scale_bits,
429 double clash_distance_in = default_mesh_elec_damping_range,
430 PrecisionModel prec = PrecisionModel::SINGLE,
431 const MeshKlManager &launcher = MeshKlManager(),
432 HybridTargetLevel availability = HybridTargetLevel::HOST);
433
434 // Clash-based meshes
435 BackgroundMesh(GridDetail kind_in, double probe_radius_in, VdwCombiningRule mixing_protocol_in,
436 const AtomGraph *ag_in, const CoordinateFrame *cf_in, double buffer,
437 double spacing, int scale_bits_in = default_globalpos_scale_bits,
438 const std::vector<double> &probe_sigma = {},
439 PrecisionModel prec = PrecisionModel::SINGLE,
440 const MeshKlManager &launcher = MeshKlManager(),
441 HybridTargetLevel availability = HybridTargetLevel::HOST);
442
443 BackgroundMesh(GridDetail kind_in, double probe_radius_in, VdwCombiningRule mixing_protocol_in,
444 const AtomGraph &ag_in, const CoordinateFrame &cf_in, double buffer,
445 double spacing, int scale_bits_in = default_globalpos_scale_bits,
446 const std::vector<double> &probe_sigma = {},
447 PrecisionModel prec = PrecisionModel::SINGLE,
448 const MeshKlManager &launcher = MeshKlManager(),
449 HybridTargetLevel availability = HybridTargetLevel::HOST);
450
451 BackgroundMesh(GridDetail kind_in, double probe_radius_in, VdwCombiningRule mixing_protocol_in,
452 const AtomGraph *ag_in, const CoordinateFrame *cf_in, double buffer,
453 const std::vector<double> &spacing,
454 int scale_bits_in = default_globalpos_scale_bits,
455 const std::vector<double> &probe_sigma = {},
456 PrecisionModel prec = PrecisionModel::SINGLE,
457 const MeshKlManager &launcher = MeshKlManager(),
458 HybridTargetLevel availability = HybridTargetLevel::HOST);
459
460 BackgroundMesh(GridDetail kind_in, double probe_radius_in, VdwCombiningRule mixing_protocol_in,
461 const AtomGraph &ag_in, const CoordinateFrame &cf_in, double buffer,
462 const std::vector<double> &spacing,
463 int scale_bits_in = default_globalpos_scale_bits,
464 const std::vector<double> &probe_sigma = {},
465 PrecisionModel prec = PrecisionModel::SINGLE,
466 const MeshKlManager &launcher = MeshKlManager(),
467 HybridTargetLevel availability = HybridTargetLevel::HOST);
468
469 BackgroundMesh(GridDetail kind_in, double probe_radius_in, VdwCombiningRule mixing_protocol_in,
470 const AtomGraph *ag_in, const CoordinateFrame *cf_in,
471 const std::vector<double> &mesh_bounds, double spacing,
472 int scale_bits_in = default_globalpos_scale_bits,
473 const std::vector<double> &probe_sigma = {},
474 PrecisionModel prec = PrecisionModel::SINGLE,
475 const MeshKlManager &launcher = MeshKlManager(),
476 HybridTargetLevel availability = HybridTargetLevel::HOST);
477
478 BackgroundMesh(GridDetail kind_in, double probe_radius_in, VdwCombiningRule mixing_protocol_in,
479 const AtomGraph &ag_in, const CoordinateFrame &cf_in,
480 const std::vector<double> &mesh_bounds, double spacing,
481 int scale_bits_in = default_globalpos_scale_bits,
482 const std::vector<double> &probe_sigma = {},
483 PrecisionModel prec = PrecisionModel::SINGLE,
484 const MeshKlManager &launcher = MeshKlManager(),
485 HybridTargetLevel availability = HybridTargetLevel::HOST);
486
487 BackgroundMesh(GridDetail kind_in, double probe_radius_in, VdwCombiningRule mixing_protocol_in,
488 const AtomGraph *ag_in, const CoordinateFrame *cf_in,
489 const std::vector<double> &mesh_bounds, const std::vector<double> &spacing,
490 int scale_bits_in = default_globalpos_scale_bits,
491 const std::vector<double> &probe_sigma = {},
492 PrecisionModel prec = PrecisionModel::SINGLE,
493 const MeshKlManager &launcher = MeshKlManager(),
494 HybridTargetLevel availability = HybridTargetLevel::HOST);
495
496 BackgroundMesh(GridDetail kind_in, double probe_radius_in, VdwCombiningRule mixing_protocol_in,
497 const AtomGraph &ag_in, const CoordinateFrame &cf_in,
498 const std::vector<double> &mesh_bounds, const std::vector<double> &spacing,
499 int scale_bits_in = default_globalpos_scale_bits,
500 const std::vector<double> &probe_sigma = {},
501 PrecisionModel prec = PrecisionModel::SINGLE,
502 const MeshKlManager &launcher = MeshKlManager(),
503 HybridTargetLevel availability = HybridTargetLevel::HOST);
504
505 // General meshes with all available details
506 template <typename Tcoord>
507 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in,
508 const AtomGraph &ag_in, const CoordinateSeries<Tcoord> &cs_in,
509 double probe_radius_in, double well_depth_in, VdwCombiningRule mixing_protocol_in,
510 const MeshParameters &measurements_in = MeshParameters(),
511 int averaging_order = default_mesh_density_averaging_order,
512 const std::vector<double> &probe_sigma = {},
513 const std::vector<double> &probe_epsilon = {},
514 double clash_distance_in = default_mesh_elec_damping_range,
515 double clash_ratio_in = default_mesh_vdw_damping_ratio,
516 PrecisionModel prec = PrecisionModel::SINGLE,
517 const MeshKlManager &launcher = MeshKlManager(),
518 HybridTargetLevel availability = HybridTargetLevel::HOST);
519
520 template <typename Tcoord>
521 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, double probe_radius_in,
522 double well_depth_in, VdwCombiningRule mixing_protocol_in, const AtomGraph &ag_in,
523 const CoordinateSeries<Tcoord> &cs_in, double buffer, double spacing,
524 int scale_bits_in = default_globalpos_scale_bits,
525 int averaging_order = default_mesh_density_averaging_order,
526 const std::vector<double> &probe_sigma = {},
527 const std::vector<double> &probe_epsilon = {},
528 double clash_distance_in = default_mesh_elec_damping_range,
529 double clash_ratio_in = default_mesh_vdw_damping_ratio,
530 PrecisionModel prec = PrecisionModel::SINGLE,
531 const MeshKlManager &launcher = MeshKlManager(),
532 HybridTargetLevel availability = HybridTargetLevel::HOST);
533
534 template <typename Tcoord>
535 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, double probe_radius_in,
536 double well_depth_in, VdwCombiningRule mixing_protocol_in, const AtomGraph &ag_in,
537 const CoordinateSeries<Tcoord> &cs_in, double buffer,
538 const std::vector<double> &spacing,
539 int scale_bits_in = default_globalpos_scale_bits,
540 int averaging_order = default_mesh_density_averaging_order,
541 const std::vector<double> &probe_sigma = {},
542 const std::vector<double> &probe_epsilon = {},
543 double clash_distance_in = default_mesh_elec_damping_range,
544 double clash_ratio_in = default_mesh_vdw_damping_ratio,
545 PrecisionModel prec = PrecisionModel::SINGLE,
546 const MeshKlManager &launcher = MeshKlManager(),
547 HybridTargetLevel availability = HybridTargetLevel::HOST);
548
549 template <typename Tcoord>
550 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, double probe_radius_in,
551 double well_depth_in, VdwCombiningRule mixing_protocol_in, const AtomGraph &ag_in,
552 const CoordinateSeries<Tcoord> &cs_in, const std::vector<double> &mesh_bounds,
553 double spacing, int scale_bits_in = default_globalpos_scale_bits,
554 int averaging_order = default_mesh_density_averaging_order,
555 const std::vector<double> &probe_sigma = {},
556 const std::vector<double> &probe_epsilon = {},
557 double clash_distance_in = default_mesh_elec_damping_range,
558 double clash_ratio_in = default_mesh_vdw_damping_ratio,
559 PrecisionModel prec = PrecisionModel::SINGLE,
560 const MeshKlManager &launcher = MeshKlManager(),
561 HybridTargetLevel availability = HybridTargetLevel::HOST);
562
563 template <typename Tcoord>
564 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, double probe_radius_in,
565 double well_depth_in, VdwCombiningRule mixing_protocol_in, const AtomGraph &ag_in,
566 const CoordinateSeries<Tcoord> &cs_in, const std::vector<double> &mesh_bounds,
567 const std::vector<double> &spacing,
568 int scale_bits_in = default_globalpos_scale_bits,
569 int averaging_order = default_mesh_density_averaging_order,
570 const std::vector<double> &probe_sigma = {},
571 const std::vector<double> &probe_epsilon = {},
572 double clash_distance_in = default_mesh_elec_damping_range,
573 double clash_ratio_in = default_mesh_vdw_damping_ratio,
574 PrecisionModel prec = PrecisionModel::SINGLE,
575 const MeshKlManager &launcher = MeshKlManager(),
576 HybridTargetLevel availability = HybridTargetLevel::HOST);
577
578 // Electrostatics-based meshes
579 template <typename Tcoord>
580 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const AtomGraph &ag_in,
581 const CoordinateSeries<Tcoord> &cs_in, double buffer, double spacing,
582 int scale_bits_in = default_globalpos_scale_bits,
583 double clash_distance_in = default_mesh_elec_damping_range,
584 PrecisionModel prec = PrecisionModel::SINGLE,
585 const MeshKlManager &launcher = MeshKlManager(),
586 HybridTargetLevel availability = HybridTargetLevel::HOST);
587
588 template <typename Tcoord>
589 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const AtomGraph &ag_in,
590 const CoordinateSeries<Tcoord> &cs_in, double buffer,
591 const std::vector<double> &spacing,
592 int scale_bits_in = default_globalpos_scale_bits,
593 double clash_distance_in = default_mesh_elec_damping_range,
594 PrecisionModel prec = PrecisionModel::SINGLE,
595 const MeshKlManager &launcher = MeshKlManager(),
596 HybridTargetLevel availability = HybridTargetLevel::HOST);
597
598 template <typename Tcoord>
599 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const AtomGraph &ag_in,
600 const CoordinateSeries<Tcoord> &cs_in, const std::vector<double> &mesh_bounds,
601 double spacing, int scale_bits_in = default_globalpos_scale_bits,
602 double clash_distance_in = default_mesh_elec_damping_range,
603 PrecisionModel prec = PrecisionModel::SINGLE,
604 const MeshKlManager &launcher = MeshKlManager(),
605 HybridTargetLevel availability = HybridTargetLevel::HOST);
606
607 template <typename Tcoord>
608 BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const AtomGraph &ag_in,
609 const CoordinateSeries<Tcoord> &cs_in, const std::vector<double> &mesh_bounds,
610 const std::vector<double> &spacing,
611 int scale_bits_in = default_globalpos_scale_bits,
612 double clash_distance_in = default_mesh_elec_damping_range,
613 PrecisionModel prec = PrecisionModel::SINGLE,
614 const MeshKlManager &launcher = MeshKlManager(),
615 HybridTargetLevel availability = HybridTargetLevel::HOST);
616
617 // Clash-based meshes
618 template <typename Tcoord>
619 BackgroundMesh(GridDetail kind_in, double probe_radius_in, const AtomGraph &ag_in,
620 const CoordinateSeries<Tcoord> &cs_in, double buffer, double spacing,
621 int scale_bits_in = default_globalpos_scale_bits,
622 int averaging_order = default_mesh_density_averaging_order,
623 const std::vector<double> &probe_sigma = {},
624 PrecisionModel prec = PrecisionModel::SINGLE,
625 const MeshKlManager &launcher = MeshKlManager(),
626 HybridTargetLevel availability = HybridTargetLevel::HOST);
627
628 template <typename Tcoord>
629 BackgroundMesh(GridDetail kind_in, double probe_radius_in, const AtomGraph &ag_in,
630 const CoordinateSeries<Tcoord> &cs_in, double buffer,
631 const std::vector<double> &spacing,
632 int scale_bits_in = default_globalpos_scale_bits,
633 int averaging_order = default_mesh_density_averaging_order,
634 const std::vector<double> &probe_sigma = {},
635 PrecisionModel prec = PrecisionModel::SINGLE,
636 const MeshKlManager &launcher = MeshKlManager(),
637 HybridTargetLevel availability = HybridTargetLevel::HOST);
638
639 template <typename Tcoord>
640 BackgroundMesh(GridDetail kind_in, double probe_radius_in, const AtomGraph &ag_in,
641 const CoordinateSeries<Tcoord> &cs_in, const std::vector<double> &mesh_bounds,
642 double spacing, int scale_bits_in = default_globalpos_scale_bits,
643 int averaging_order = default_mesh_density_averaging_order,
644 const std::vector<double> &probe_sigma = {},
645 PrecisionModel prec = PrecisionModel::SINGLE,
646 const MeshKlManager &launcher = MeshKlManager(),
647 HybridTargetLevel availability = HybridTargetLevel::HOST);
648
649 template <typename Tcoord>
650 BackgroundMesh(GridDetail kind_in, double probe_radius_in, const AtomGraph &ag_in,
651 const CoordinateSeries<Tcoord> &cs_in, const std::vector<double> &mesh_bounds,
652 const std::vector<double> &spacing,
653 int scale_bits_in = default_globalpos_scale_bits,
654 int averaging_order = default_mesh_density_averaging_order,
655 const std::vector<double> &probe_sigma = {},
656 PrecisionModel prec = PrecisionModel::SINGLE,
657 const MeshKlManager &launcher = MeshKlManager(),
658 HybridTargetLevel availability = HybridTargetLevel::HOST);
660
668 BackgroundMesh(const BackgroundMesh<T> &original) = default;
669 BackgroundMesh(BackgroundMesh<T> &&original) = default;
670 BackgroundMesh& operator=(const BackgroundMesh<T> &original) = default;
671 BackgroundMesh& operator=(BackgroundMesh<T> &&original) = default;
673
676
687
690
693
696 template <typename Tcoord>
698
700 size_t getEnsembleTypeCode() const;
701
703 GridDetail getMeshKind() const;
704
707 NonbondedPotential getNonbondedPotential() const;
708
710 double getProbeRadius() const;
711
714 double getWellDepth() const;
715
719 double getClashDistance() const;
720
725 double getClashRatio() const;
726
728 double getOcclusionPenalty() const;
729
732 VdwCombiningRule getCombiningRule() const;
733
735 PrecisionModel getBuildPrecision() const;
736
739
743
752 const BackgroundMeshReader<T> data(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
753 BackgroundMeshWriter<T> data(HybridTargetLevel tier = HybridTargetLevel::HOST);
755
765 templateFreeData(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
766
768 templateFreeData(HybridTargetLevel tier = HybridTargetLevel::HOST);
770
771#ifdef STORMM_USE_HPC
772# ifdef STORMM_USE_CUDA
780 BackgroundMeshWriter<void> deviceViewToTemplateFreeHostData();
781 const BackgroundMeshReader<void> deviceViewToTemplateFreeHostData() const;
783# endif
784#endif
785
789 const MeshFFKit<T>
790 getNonbondedKit(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
791
796 getReferenceNonbondedKit(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
797
802 const MeshFFKit<void>
803 templateFreeNonbondedKit(HybridTargetLevel tier = HybridTargetLevel::HOST) const;
804
805#ifdef STORMM_USE_HPC
806# ifdef STORMM_USE_CUDA
809 const MeshFFKit<double> deviceViewToReferenceNonbondedKit() const;
810
813 const MeshFFKit<void> deviceViewToTemplateFreeNonbondedKit() const;
814# endif
816 void upload();
817
819 void download();
820
823 void uploadFrame();
824
827 void downloadFrame();
828
830 void uploadForceField();
831
833 void downloadForceField();
834#endif
835
866 void setMeshParameters(const AtomGraph *ag_in, const CoordinateFrame *cf_in, double padding,
867 double spacing, int scale_bits_in = -100);
868
869 void setMeshParameters(const AtomGraph *ag_in, const CoordinateFrame *cf_in, double padding,
870 const std::vector<double> &spacing, int scale_bits_in = -100);
871
872 void setMeshParameters(double padding, double spacing, int scale_bits_in = -100);
873
874 void setMeshParameters(double padding, const std::vector<double> &spacing,
875 int scale_bits_in = -100);
876
877 void setMeshParameters(const std::vector<double> &mesh_bounds, double spacing,
878 int scale_bits_in = -100);
879
880 void setMeshParameters(const std::vector<double> &mesh_bounds,
881 const std::vector<double> &spacing, int scale_bits_in = -100);
883
888 void setProbeRadius(double probe_radius_in);
889
894 void setWellDepth(double well_depth_in);
895
901 void setOcclusionPenalty(double occlusion_penalty_in);
902
918 void validateCombiningRule(VdwCombiningRule mixing_protocol_in,
919 const std::vector<double> &probe_sigma,
920 const std::vector<double> &probe_epsilon);
921
926 void setCoefficientScalingBits(int scaling_bits_in);
927
934 void computeField(const MeshKlManager &launcher = MeshKlManager(),
935 PrecisionModel prec = PrecisionModel::SINGLE,
936 HybridTargetLevel availability = HybridTargetLevel::HOST,
937 const std::vector<double> &probe_sigma = {},
938 const std::vector<double> &probe_epsilon = {});
939
940private:
941
944 MeshParameters measurements;
945
956 GridDetail kind;
957
961 NonbondedPotential field;
962
964 double probe_radius;
965
970 double well_depth;
971
974 double occlusion_penalty;
975
979 PrecisionModel build_precision;
980
982 MeshRulers tick_marks;
983
986 MeshForceField<T> nonbonded_model;
987
999 Hybrid<T> coefficients;
1000
1001 // Mesh coefficients may come in fixed-precision format. These coefficients mirror contents of
1002 // other templated classes to track the conversion.
1003 int coefficient_scale_bits;
1005 double coefficient_scale;
1007 double inverse_coefficient_scale;
1013
1015 MeshFoundation basis;
1016
1020 void allocate();
1021
1043 void setNonbondedModel(VdwCombiningRule mixing_protocol_in, double clash_ratio_in,
1044 double clash_distance_in, const std::vector<double> &probe_sigma,
1045 const std::vector<double> &probe_epsilon);
1046
1049 void validateMeshKind() const;
1050
1054 void validateScalingBits() const;
1055
1066 void colorOcclusionMesh(const MeshKlManager &launcher,
1067 PrecisionModel prec = PrecisionModel::SINGLE,
1068 HybridTargetLevel availability = HybridTargetLevel::HOST,
1069 const std::vector<double> &probe_sigma = {});
1070
1075 void mapOccupancyField(const MeshKlManager &launcher,
1076 PrecisionModel prec = PrecisionModel::SINGLE,
1077 HybridTargetLevel availability = HybridTargetLevel::HOST,
1078 const std::vector<double> &probe_sigma = {});
1079
1105 void computeCoefficients(const std::vector<double> &u_grid,
1106 const std::vector<double> &dudx_grid,
1107 const std::vector<double> &dudy_grid,
1108 const std::vector<double> &dudz_grid,
1109 const std::vector<double> &dudxx_grid,
1110 const std::vector<double> &dudxy_grid,
1111 const std::vector<double> &dudxz_grid,
1112 const std::vector<double> &dudyy_grid,
1113 const std::vector<double> &dudyz_grid,
1114 const std::vector<double> &dudxxx_grid,
1115 const std::vector<double> &dudxxy_grid,
1116 const std::vector<double> &dudxxz_grid,
1117 const std::vector<double> &dudxyy_grid,
1118 const std::vector<double> &dudxyz_grid,
1119 const TricubicStencil &tc_weights);
1120
1121 void computeCoefficients(const std::vector<double> &u_grid,
1122 const std::vector<double> &dudx_grid,
1123 const std::vector<double> &dudy_grid,
1124 const std::vector<double> &dudz_grid,
1125 const std::vector<double> &dudxy_grid,
1126 const std::vector<double> &dudxz_grid,
1127 const std::vector<double> &dudyz_grid,
1128 const std::vector<double> &dudxyz_grid,
1129 const TricubicStencil &tc_weights);
1131
1140 template <typename Tcalc, typename Tcalc2, typename Tcalc3, typename Tcalc4>
1141 void colorNonbondedField(const NonbondedKit<Tcalc> &nbk, const TricubicStencil &tc_weights,
1142 const std::vector<double> &sigma_table = {},
1143 const std::vector<double> &eps_table = {});
1144
1153 void mapPureNonbondedPotential(const MeshKlManager &launcher,
1154 PrecisionModel prec = PrecisionModel::SINGLE,
1155 HybridTargetLevel availability = HybridTargetLevel::HOST,
1156 const std::vector<double> &sigma_table = {},
1157 const std::vector<double> &eps_table = {});
1158};
1159
1173template <typename T> BackgroundMeshWriter<T> restoreType(BackgroundMeshWriter<void> *rasa);
1174template <typename T> BackgroundMeshWriter<T> restoreType(const BackgroundMeshWriter<void> &rasa);
1175template <typename T> BackgroundMeshReader<T> restoreType(const BackgroundMeshReader<void> *rasa);
1176template <typename T> BackgroundMeshReader<T> restoreType(const BackgroundMeshReader<void> &rasa);
1178
1185template <typename Ttarget, typename Tcontrib>
1186void checkMeshCompatibility(const BackgroundMeshWriter<Ttarget> &target,
1187 const BackgroundMeshReader<Tcontrib> &contrib);
1188
1197template <typename Tcalc4, typename Tcalc2, typename Tcalc>
1198Tcalc2 evalElecCSC(const Tcalc4 abcd, Tcalc r, Tcalc rswitch, Tcalc q);
1199
1203template <typename Tcalc4, typename Tcalc2, typename Tcalc>
1204Tcalc evalElecCSCND(const Tcalc4 abcd, Tcalc r, Tcalc rswitch, Tcalc q);
1205
1215template <typename Tcalc4, typename Tcalc2, typename Tcalc>
1216Tcalc2 evalLennardJonesCSC(const Tcalc4 abcd, Tcalc r, Tcalc rswitch, Tcalc lja, Tcalc ljb);
1217
1221template <typename Tcalc4, typename Tcalc2, typename Tcalc>
1222Tcalc evalLennardJonesCSCND(const Tcalc4 abcd, Tcalc r, Tcalc rswitch, Tcalc lja, Tcalc ljb);
1223
1224#ifdef STORMM_USE_HPC
1243void launchColorOcclusionMesh(BackgroundMeshWriter<void> *bgmw, const AtomGraph *ag,
1244 const CoordinateFrameReader &cfr, PrecisionModel prec,
1245 const MeshKlManager &launcher);
1246#endif
1247
1264void accumulateOcclusionMesh(BackgroundMeshWriter<llint> *target,
1265 const BackgroundMeshReader<llint> &contribution);
1266
1267void accumulateOcclusionMesh(BackgroundMesh<llint> *target,
1268 const BackgroundMesh<llint> &contribution,
1269 const GpuDetails &gpu = null_gpu,
1270 HybridTargetLevel availability = HybridTargetLevel::HOST);
1272
1273#ifdef STORMM_USE_HPC
1276void launchAccOcclusionMesh(BackgroundMeshWriter<llint> *target,
1277 const BackgroundMeshReader<llint> &contribution,
1278 const GpuDetails &gpu);
1279#endif
1280
1291template <typename T>
1292void accumulateNonbondedFieldMesh(BackgroundMeshWriter<llint> *target,
1293 const BackgroundMeshReader<void> &contribution,
1294 size_t ct_contrib, const GpuDetails &gpu);
1295
1296template <typename T>
1297void accumulateNonbondedFieldMesh(BackgroundMesh<llint> *target,
1298 const BackgroundMesh<T> &contribution,
1299 const GpuDetails &gpu, HybridTargetLevel availability);
1301
1302#ifdef STORMM_USE_HPC
1320void launchColorNonbondedFieldMesh(BackgroundMeshWriter<void> *target, size_t ct_targ,
1321 const MeshFFKit<double> &mnbk, PrecisionModel prec,
1322 const double* stn_xfrm, const AtomGraph *ag,
1323 const CoordinateFrameReader &cfr,
1324 const MeshKlManager &launcher,
1325 const HybridTargetLevel availability = HybridTargetLevel::HOST);
1326
1335void launchAccNonbondedFieldMesh(BackgroundMeshWriter<llint> *target,
1336 const BackgroundMeshReader<void> &contribution,
1337 const size_t ct_contrib, const GpuDetails &gpu);
1338#endif
1339
1356template <typename T>
1357void inferPatchBoundary(std::vector<T> *u_patch, int dim_a, int dim_b, int dim_c, int start_a,
1358 int start_b, int start_c, int end_a, int end_b, int end_c, int inc_a,
1359 int inc_b, int inc_c);
1360
1387template <typename T>
1388void transferPatchResult(const std::vector<T> &patch, int patch_dim_a, int patch_dim_b,
1389 int patch_dim_c, BackgroundMeshWriter<T> *occfieldw, int a_index,
1390 int b_index, int c_index, int readout_a_start, int readout_b_start,
1391 int readout_c_start, int element_offset);
1392
1404template <typename T>
1405void patchNumericDeriv(BackgroundMeshWriter<T> *occfieldw, int a_index, int b_index, int c_index,
1406 T max_occlusion);
1407
1419template <typename T>
1420void occlusionFieldDerivatives(BackgroundMesh<T> *occfield, const GpuDetails &gpu,
1421 HybridTargetLevel availability);
1422
1423#ifdef STORMM_USE_HPC
1430void launchOccFieldDerivativeCalc(BackgroundMeshWriter<void> *bgmw, size_t ct_mesh);
1431#endif
1432
1433} // namespace structure
1434} // namespace stormm
1435
1436#include "background_mesh.tpp"
1437#include "background_mesh_freefunc.tpp"
1438
1439#endif
Pertinent aspects of one particular GPU. Condensing the data for each GPU in this manner helps to ens...
Definition gpu_details.h:27
Collect and dispense the launch parameters for various mesh-based kernels. The launches for these ker...
Definition mesh_kernel_manager.h:29
BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const MeshParameters &measurements_in=MeshParameters(), const PrecisionModel build_precision_in=PrecisionModel::SINGLE)
The constructor takes all dimension parameters plus an indication of what type of potential,...
A workspace for constructing a pure potential mesh based on the frozen atoms of a large molecule....
Definition background_mesh.h:195
const BackgroundMeshReader< void > templateFreeData(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get an abstract of the mesh with any templating removed.
void setCoefficientScalingBits(int scaling_bits_in)
Set the number of bits after the decimal to be used in fixed-precision representations of coordinates...
void setProbeRadius(double probe_radius_in)
Set the probe radius, meaning either the Lennard-Jones potential sigma radius or the clash probe radi...
const MeshFFKit< double > getReferenceNonbondedKit(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get the abstract of the mesh in double precision.
const AtomGraph * getTopologyPointer() const
Get a const pointer to the topology responsible for creating this mesh.
size_t getEnsembleTypeCode() const
Get the codified data type of the ensemble.
const MeshFFKit< T > getNonbondedKit(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get the abstract of the mesh in the precision of the mesh's coefficient data.
NonbondedPotential getNonbondedPotential() const
Get the non-bonded potential expressed on the mesh. This will return an error if the mesh is not asso...
const CoordinateSeries< Tcoord > * getEnsemblePointer() const
Get a const reference to the array of structures (each stored as a CoordinateFrame) reponsible for cr...
VdwCombiningRule getCombiningRule() const
Get the mixing protocol used for van-der Waals (as expressed by the Lennard-Jones model) interactions...
double getProbeRadius() const
Get the probe radius. Valid for all types of meshes.
void setWellDepth(double well_depth_in)
Set the Lennard-Jones well depth for the probe that will generate the potential. This includes a vali...
double getOcclusionPenalty() const
Get the penalty associated with a collision on an occlusion mesh.
BackgroundMesh(GridDetail kind_in, NonbondedPotential field_in, const MeshParameters &measurements_in=MeshParameters(), const PrecisionModel build_precision_in=PrecisionModel::SINGLE)
The constructor takes all dimension parameters plus an indication of what type of potential,...
void setMeshParameters(const AtomGraph *ag_in, const CoordinateFrame *cf_in, double padding, double spacing, int scale_bits_in=-100)
Set the dimensions of the mesh. Variants of this member function that do not include a coordinate poi...
const MeshFFKit< void > templateFreeNonbondedKit(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get a template-free form of the of the mesh nonbonded abstract in the mesh coefficient data type,...
int getCoefficientScalingBits() const
Get the number of bits after the decimal in fixed-precision mesh coefficient representations.
double getClashDistance() const
Get the absolute distance at which softcore electrostatic interactions take over. This is provided fo...
double getCoefficientScalingFactor() const
Get the scaling factor for mesh coefficients in fixed-precision representations.
const MeshParameters & getDimensions() const
Get an object describing the mesh dimensions.
void setOcclusionPenalty(double occlusion_penalty_in)
Set the penalty associated with striking occupied volume on an occlusion mesh. This is intended to be...
const BackgroundMeshReader< T > data(HybridTargetLevel tier=HybridTargetLevel::HOST) const
Get an abstract of the mesh.
double getClashRatio() const
Get the ratio of the van-der Waals (Lennard-Jones) sigma parameter for interacting pairs of particles...
const CoordinateFrame * getCoordinatePointer() const
Get a const pointer to the coordinates responsible for creating this mesh.
GridDetail getMeshKind() const
Get the type of mesh.
BackgroundMesh(const BackgroundMesh< T > &original)=default
Copy and move constructors as well as assignment operators can be set to their defaults because the o...
void computeField(const MeshKlManager &launcher=MeshKlManager(), PrecisionModel prec=PrecisionModel::SINGLE, HybridTargetLevel availability=HybridTargetLevel::HOST, const std::vector< double > &probe_sigma={}, const std::vector< double > &probe_epsilon={})
Compute the appropriate field for the mesh. This is called automatically by the constructor if enough...
double getWellDepth() const
Get the non-bonded probe well depth. This is valid only for non-bonded fields, in particular those as...
PrecisionModel getBuildPrecision() const
Get the precisionmodel under which the mesh was calculated.
void validateCombiningRule(VdwCombiningRule mixing_protocol_in, const std::vector< double > &probe_sigma, const std::vector< double > &probe_epsilon)
Check the combining rule that will be used to make the probe interact with the receptor on any mesh (...
const MeshFoundation & getMolecularBasis() const
Get the collection of molecular details underlying the mesh, including the topology and associated co...
A small object to hold the essential descriptors of the molecule.
Definition mesh_foundation.h:50
Collect and dispense the launch parameters for various mesh-based kernels. The launches for these ker...
Definition mesh_kernel_manager.h:29
Encode the critical dimensions of a regular, rectilinear mesh. The locations of mesh points as well a...
Definition mesh_parameters.h:121
A struct to hold information relating to an Amber topology. This struct's member functions are limite...
Definition atomgraph.h:50
Store the coordinates and box information for a frame, only. This abridged struct can serve when the ...
Definition coordinateframe.h:111
Store the coordinates and box information for a series of frames, in one of several levels of precisi...
Definition coordinate_series.h:137
BackgroundMeshReader(const MeshParamKit &dims_in, GridDetail kind_in, NonbondedPotential field_in, const MeshRulerKit &rulers, const Tdata *coeffs_in, double coeff_scale_in, double probe_radius_in, double well_depth_in, double occ_cost_in, const MeshBasicsKit &mbss_in)
The constructor takes arguments for all member variables.
The templated, writeable abstract of a BackgroundMesh object.
Definition background_mesh.h:142
const Tdata * coeffs
Definition background_mesh.h:163
const double well_depth
Definition background_mesh.h:181
const MeshParamKit dims
Dimensions of the mesh. These are pre-established.
Definition background_mesh.h:159
const NonbondedPotential field
The field described by the mesh, also pre-established.
Definition background_mesh.h:161
const MeshBasicsKit mbss
Definition background_mesh.h:184
const double probe_radius
Definition background_mesh.h:178
const float coeff_scale_f
Single-precision variant of coeff_scale.
Definition background_mesh.h:175
const double coeff_scale
Definition background_mesh.h:168
const double occ_cost
Energetic penalty of an occlusion interaction.
Definition background_mesh.h:183
const GridDetail kind
The type of mesh, also pre-established.
Definition background_mesh.h:160
const MeshRulerKit rlrs
Coordinate series for tick marks along each mesh axis.
Definition background_mesh.h:162
BackgroundMeshReader(const MeshParamKit &dims_in, GridDetail kind_in, NonbondedPotential field_in, const MeshRulerKit &rulers, const Tdata *coeffs_in, double coeff_scale_in, double probe_radius_in, double well_depth_in, double occ_cost_in, const MeshBasicsKit &mbss_in)
The constructor takes arguments for all member variables.
BackgroundMeshReader(const BackgroundMeshReader< Tdata > &original)=default
The default copy and move constructors will be valid for this object. Const members negate the use of...
The templated, writeable abstract of a BackgroundMesh object.
Definition background_mesh.h:93
BackgroundMeshWriter(const BackgroundMeshWriter< Tdata > &original)=default
The default copy and move constructors will be valid for this object. Const members negate the use of...
Tdata * coeffs
Definition background_mesh.h:116
const double occ_cost
Energetic penalty of an occlusion interaction.
Definition background_mesh.h:135
const double probe_radius
Definition background_mesh.h:130
const GridDetail kind
The type of mesh, also pre-established.
Definition background_mesh.h:111
const MeshBasicsKit mbss
Definition background_mesh.h:136
const MeshParamKit dims
Dimensions of the mesh. These are pre-established.
Definition background_mesh.h:110
const double well_depth
Definition background_mesh.h:133
const NonbondedPotential field
The field described by the mesh, also pre-established.
Definition background_mesh.h:112
const float coeff_scale_f
Single-precision variant of coeff_scale.
Definition background_mesh.h:126
const double coeff_scale
Definition background_mesh.h:120
BackgroundMeshWriter(const MeshParamKit &measurements, GridDetail kind_in, NonbondedPotential field_in, const MeshRulerKit &rulers, Tdata *coeffs_in, double coeff_scale_in, double probe_radius_in, double well_depth_in, double occ_cost_in, const MeshBasicsKit &mbss_in)
The constructor takes arguments for all member variables.
const MeshRulerKit rlrs
Coordinate series for tick marks along each mesh axis.
Definition background_mesh.h:113
Abstract for the MeshFoundation object to encapsulate pointers to its data on the host or the device.
Definition mesh_foundation.h:29
An object to store the probe's specific non-bonded parameters and softcore function coefficients used...
Definition mesh_forcefield.h:35
The abstract of a MeshParameters object is read-only (modify the original to get a new abstract if th...
Definition mesh_parameters.h:61
Definition mesh_rulers.h:21
Collect C-style pointers for the elements of a read-only CoordinateFrame object.
Definition coordinateframe.h:65