/usr/include/liggghts/dump_particle.h is in libliggghts-dev 3.7.0+repack1-1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 | /* ----------------------------------------------------------------------
This is the
██╗ ██╗ ██████╗ ██████╗ ██████╗ ██╗ ██╗████████╗███████╗
██║ ██║██╔════╝ ██╔════╝ ██╔════╝ ██║ ██║╚══██╔══╝██╔════╝
██║ ██║██║ ███╗██║ ███╗██║ ███╗███████║ ██║ ███████╗
██║ ██║██║ ██║██║ ██║██║ ██║██╔══██║ ██║ ╚════██║
███████╗██║╚██████╔╝╚██████╔╝╚██████╔╝██║ ██║ ██║ ███████║
╚══════╝╚═╝ ╚═════╝ ╚═════╝ ╚═════╝ ╚═╝ ╚═╝ ╚═╝ ╚══════╝®
DEM simulation engine, released by
DCS Computing Gmbh, Linz, Austria
http://www.dcs-computing.com, office@dcs-computing.com
LIGGGHTS® is part of CFDEM®project:
http://www.liggghts.com | http://www.cfdem.com
Core developer and main author:
Christoph Kloss, christoph.kloss@dcs-computing.com
LIGGGHTS® is open-source, distributed under the terms of the GNU Public
License, version 2 or later. It is distributed in the hope that it will
be useful, but WITHOUT ANY WARRANTY; without even the implied warranty
of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. You should have
received a copy of the GNU General Public License along with LIGGGHTS®.
If not, see http://www.gnu.org/licenses . See also top-level README
and LICENSE files.
LIGGGHTS® and CFDEM® are registered trade marks of DCS Computing GmbH,
the producer of the LIGGGHTS® software and the CFDEM®coupling software
See http://www.cfdem.com/terms-trademark-policy for details.
-------------------------------------------------------------------------
Contributing author and copyright for this file:
(if no contributing author is listed, this file has been contributed
by the core developer)
Arno Mayrhofer (DCS Computing GmbH)
Copyright 2016- DCS Computing GmbH, Linz
------------------------------------------------------------------------- */
#if defined(LAMMPS_VTK)
#ifndef LMP_DUMP_PARTICLE_H
#define LMP_DUMP_PARTICLE_H
#include "pointers.h"
#include "sort_buffer.h"
#include <map>
#include <set>
#include <string>
#include <list>
#include <vtkSmartPointer.h>
#include <vtkPoints.h>
#include <vtkCellArray.h>
#include <vtkPolyVertex.h>
#include <vtkMultiBlockDataSet.h>
class vtkAbstractArray;
class vtkRectilinearGrid;
class vtkUnstructuredGrid;
namespace LAMMPS_NS {
/**
* @brief DumpParticle class
* write atom data to vtk files.
*
* Similar to the DumpCustom class but uses the vtk library to write data to vtk simple
* legacy or xml format depending on the filename extension specified. (Since this
* conflicts with the way binary output is specified, dump_modify allows to set the
* binary flag for this dump command explicitly).
* In contrast to DumpCustom class the attributes to be packed are stored in a std::map
* to avoid duplicate entries and enforce correct ordering of vector components (except
* for computes and fixes - these have to be given in the right order in the input script).
* (Note: std::map elements are sorted by their keys.)
* This dump command does not support compressed files, buffering or custom format strings,
* multiproc is only supported by the xml formats, multifile option has to be used.
*/
class DumpParticle : public Pointers {
public:
DumpParticle(class LAMMPS *, int, int, int, int, int, int);
virtual ~DumpParticle();
int parse_parameters(int narg, const char *const *const arg, const bool pp_keyword_optional = false, std::list<std::string> keyword_list = std::list<std::string>());
void prepare_mbSet(vtkSmartPointer<vtkMultiBlockDataSet> mbSet, bool usePolyData = false);
bigint memory_usage();
virtual int modify_param(int, char **);
virtual void init_style();
int count();
protected:
int nevery; // dump frequency for output
int nclusterprocs; // number of procs that write to one file
int multiproc; // number of procs writing files
int filewriter; // 1 if this proc writes a file, else 0
int fileproc; // ID of proc in my cluster who writes to file
int iregion; // -1 if no region, else which region
char *idregion; // region ID
int igroup; // group id
int groupbit; // group mask
int nthresh; // # of defined thresholds
int *thresh_array; // array to threshold on for each nthresh
int *thresh_op; // threshold operation for each nthresh
double *thresh_value; // threshold value for each nthresh
int nchoose; // # of selected atoms
int maxlocal; // size of atom selection and variable arrays
int *choose; // local indices of selected atoms
double *dchoose; // value for each atom to threshhold against
int *clist; // compressed list of indices of selected atoms
int nfield; // # of keywords listed by user
int size_one; // number of doubles used per particle
std::map<int, int> field2index; // which compute,fix,variable calcs this field
std::map<int, int> argindex; // index into compute,fix scalar_atom,vector_atom
// 0 for scalar_atom, 1-N for vector_atom values
int ncompute; // # of Compute objects used by dump
char **id_compute; // their IDs
class Compute **compute; // list of ptrs to the Compute objects
int nfix; // # of Fix objects used by dump
char **id_fix; // their IDs
class Fix **fix; // list of ptrs to the Fix objects
int nvariable; // # of Variables used by dump
char **id_variable; // their names
int *variable; // list of indices for the Variables
double **vbuf; // local storage for variable evaluation
int ntypes; // # of atom types
char **typenames; // array of element names for each type
int maxbuf; // max size of buffer
double *buf; // array buffer
// private methods
void pack(int *);
virtual void write_data(int, double *, vtkSmartPointer<vtkMultiBlockDataSet>, bool usePolyData);
void identify_vectors();
void identify_tensor();
int add_compute(char *);
int add_fix(char *);
int add_variable(char *);
void prepare_domain_data(vtkRectilinearGrid *);
void prepare_domain_data_triclinic(vtkUnstructuredGrid *);
void write_domain_vtk();
void write_domain_vtk_triclinic();
void write_domain_vtr();
void write_domain_vtu_triclinic();
typedef void (DumpParticle::*FnPtrPack)(int);
std::map<int, FnPtrPack> pack_choice; // ptrs to pack functions
std::map<int, int> vtype; // data type (INT, DOUBLE,...) for each type of entry for each atommyarrays
std::map<int, std::string> name; // attribute labels (e.g. "x", "v", ...) for each type of entry for each atommyarrays
std::set<int> vector_set; // set of vector attributes for each type of entry for each atommyarrays
int current_pack_choice_key;
// vtk data containers
vtkSmartPointer<vtkPoints> points; // list of points, one point for each point cell
vtkSmartPointer<vtkCellArray> pointsCells; // list of point cells
std::map<int, vtkSmartPointer<vtkAbstractArray> > myarrays; // list of a list of arrays that is presents data for each atom (x, v,...)
// is then added to the point cells upon writing
int n_calls_;
double (*boxcorners)[3]; // corners of triclinic domain box
bool convex_hull_detected;
int convex_hull_max_n_tri;
bool tensor_detected;
double boxxlo,boxxhi; // local copies of domain values
double boxylo,boxyhi; // lo/hi are bounding box for triclinic
double boxzlo,boxzhi;
SortBuffer *sortBuffer;
void setFileCurrent();
void buf2arrays(int, double *); // transfer data from buf array to vtk arrays
void reset_vtk_data_containers();
// customize by adding a method prototype
void pack_compute(int);
void pack_fix(int);
void pack_variable(int);
void pack_id(int);
void pack_molecule(int);
void pack_type(int);
void pack_mass(int);
void pack_x(int);
void pack_y(int);
void pack_z(int);
void pack_points_convexhull(int);
void pack_xs(int);
void pack_ys(int);
void pack_zs(int);
void pack_xs_triclinic(int);
void pack_ys_triclinic(int);
void pack_zs_triclinic(int);
void pack_xu(int);
void pack_yu(int);
void pack_zu(int);
void pack_xu_triclinic(int);
void pack_yu_triclinic(int);
void pack_zu_triclinic(int);
void pack_xsu(int);
void pack_ysu(int);
void pack_zsu(int);
void pack_xsu_triclinic(int);
void pack_ysu_triclinic(int);
void pack_zsu_triclinic(int);
void pack_ix(int);
void pack_iy(int);
void pack_iz(int);
void pack_vx(int);
void pack_vy(int);
void pack_vz(int);
void pack_fx(int);
void pack_fy(int);
void pack_fz(int);
void pack_q(int);
void pack_density(int);
void pack_p(int);
void pack_rho(int);
void pack_mux(int);
void pack_muy(int);
void pack_muz(int);
void pack_mu(int);
void pack_radius(int);
void pack_diameter(int);
void pack_omegax(int);
void pack_omegay(int);
void pack_omegaz(int);
void pack_angmomx(int);
void pack_angmomy(int);
void pack_angmomz(int);
void pack_tqx(int);
void pack_tqy(int);
void pack_tqz(int);
void pack_spin(int);
void pack_eradius(int);
void pack_ervel(int);
void pack_erforce(int);
void pack_shapex(int);
void pack_shapey(int);
void pack_shapez(int);
void pack_roundness1(int);
void pack_roundness2(int);
void pack_quat1(int);
void pack_quat2(int);
void pack_quat3(int);
void pack_quat4(int);
void pack_tensor(int n);
union ubuf {
double d;
int64_t i;
ubuf(double arg) : d(arg) {}
ubuf(int64_t arg) : i(arg) {}
ubuf(int arg) : i(arg) {}
};
};
}
#endif
#endif
/* ERROR/WARNING messages:
E: No dump custom arguments specified
The dump custom command requires that atom quantities be specified to
output to dump file.
E: Invalid attribute in dump custom command
Self-explantory.
E: Could not find dump custom compute ID
The compute ID needed by dump custom to compute a per-atom quantity
does not exist.
E: Could not find dump custom fix ID
Self-explanatory.
E: Dump custom and fix not computed at compatible times
The fix must produce per-atom quantities on timesteps that dump custom
needs them.
E: Could not find dump custom variable name
Self-explanatory.
E: Region ID for dump custom does not exist
Self-explanatory.
E: Threshhold for an atom property that isn't allocated
A dump threshhold has been requested on a quantity that is
not defined by the atom style used in this simulation.
E: Dumping an atom property that isn't allocated
The chosen atom style does not define the per-atom quantity being
dumped.
E: Dumping an atom quantity that isn't allocated
Only per-atom quantities that are defined for the atom style being
used are allowed.
E: Dump custom compute does not compute per-atom info
Self-explanatory.
E: Dump custom compute does not calculate per-atom vector
Self-explanatory.
E: Dump custom compute does not calculate per-atom array
Self-explanatory.
E: Dump custom compute vector is accessed out-of-range
Self-explanatory.
E: Dump custom fix does not compute per-atom info
Self-explanatory.
E: Dump custom fix does not compute per-atom vector
Self-explanatory.
E: Dump custom fix does not compute per-atom array
Self-explanatory.
E: Dump custom fix vector is accessed out-of-range
Self-explanatory.
E: Dump custom variable is not atom-style variable
Only atom-style variables generate per-atom quantities, needed for
dump output.
E: Illegal ... command
Self-explanatory. Check the input script syntax and compare to the
documentation for the command. You can use -echo screen as a
command-line option when running LAMMPS to see the offending line.
E: Dump_modify region ID does not exist
Self-explanatory.
E: Dump modify element names do not match atom types
Number of element names must equal number of atom types.
E: Invalid attribute in dump modify command
Self-explantory.
E: Could not find dump modify compute ID
Self-explanatory.
E: Dump modify compute ID does not compute per-atom info
Self-explanatory.
E: Dump modify compute ID does not compute per-atom vector
Self-explanatory.
E: Dump modify compute ID does not compute per-atom array
Self-explanatory.
E: Dump modify compute ID vector is not large enough
Self-explanatory.
E: Could not find dump modify fix ID
Self-explanatory.
E: Dump modify fix ID does not compute per-atom info
Self-explanatory.
E: Dump modify fix ID does not compute per-atom vector
Self-explanatory.
E: Dump modify fix ID does not compute per-atom array
Self-explanatory.
E: Dump modify fix ID vector is not large enough
Self-explanatory.
E: Could not find dump modify variable name
Self-explanatory.
E: Dump modify variable is not atom-style variable
Self-explanatory.
E: Invalid dump_modify threshhold operator
Operator keyword used for threshold specification in not recognized.
*/
|