ls1-MarDyn
ls1-MarDyn molecular dynamics code
MmpldWriter.h
1#ifndef MMPLDWRITER_H_
2#define MMPLDWRITER_H_
3
4#include <string>
5#include <vector>
6#include <array>
7
8#include "plugins/PluginBase.h"
9#include "molecules/MoleculeForwardDeclaration.h"
10
11#ifdef ENABLE_MPI
12#include "utils/MPI_Info_object.h"
13#endif
14
15enum InitSphereData : uint8_t
16{
17 ISD_USE_DEFAULT = 1,
18 ISD_READ_FROM_FILE = 2,
19 ISD_READ_FROM_XML = 3,
20};
21
22enum MMPLD_Vertex_type : uint8_t {
23 MMPLD_VERTEX_NONE = 0,
24 MMPLD_VERTEX_FLOAT_XYZ = 1,
25 MMPLD_VERTEX_FLOAT_XYZR = 2,
26 MMPLD_VERTEX_SHORT_XYZ = 3,
27};
28
29enum MMPLD_Color_type : uint8_t {
30 MMPLD_COLOR_NONE = 0,
31 MMPLD_COLOR_UINT8_RGB = 1,
32 MMPLD_COLOR_UINT8_RGBA = 2,
33 MMPLD_COLOR_FLOAT_I = 3,
34 MMPLD_COLOR_FLOAT_RGB = 4,
35 MMPLD_COLOR_FLOAT_RGBA = 5,
36};
37
40class MmpldWriter : public PluginBase
41{
42protected:
55 MmpldWriter(uint64_t startTimestep, uint64_t writeFrequency, uint64_t stopTimestep, uint64_t numFramesPerFile,
56 std::string outputPrefix);
57 virtual ~~MmpldWriter() {}
58
59 virtual void SetNumSphereTypes() {}
60 virtual void CalcNumSpheresPerType(ParticleContainer* particleContainer, uint64_t* numSpheresPerType) {}
61 virtual bool GetSpherePos(float *spherePos, Molecule* mol, uint8_t& nSphereTypeIndex) { return false; }
62
63 void InitSphereData();
64
65public:
66 void readXML(XMLfileUnits& xmlconfig);
67
68 void init(ParticleContainer *particleContainer,
69 DomainDecompBase *domainDecomp, Domain *domain);
70 void endStep(
71 ParticleContainer *particleContainer,
72 DomainDecompBase *domainDecomp, Domain *domain,
73 unsigned long simstep
74 );
75 void finish(ParticleContainer *particleContainer,
76 DomainDecompBase *domainDecomp, Domain *domain);
77
79 return std::string("MmpldWriter");
80 }
81 static PluginBase* createInstance() { return new MmpldWriter(); }
82
83 void SetInitSphereDataParameters(const uint8_t &bInitSphereData, const std::string &strSphereDataFilename) {
84 _bInitSphereData = bInitSphereData; _strSphereDataFilename = strSphereDataFilename;
85 }
86
87protected:
88 std::string getOutputFilename();
89 void MultiFileApproachReset(ParticleContainer* particleContainer,
90 DomainDecompBase* domainDecomp, Domain* domain);
91 void PrepareWriteControl();
92 long get_data_frame_header_size();
93 long get_seekTable_size();
94 void writeSeekTableEntry(int id, uint64_t offset);
95 long get_data_list_header_size();
96 long get_particle_data_size();
97 long get_data_list_size(uint64_t particle_count);
98 void write_frame_header(uint32_t num_data_lists);
99 void write_particle_list_header(uint64_t particle_count, int sphereId);
100 void write_frame(ParticleContainer* particleContainer, DomainDecompBase* domainDecomp);
101
102protected:
109 long _writeBufferSize;
110 std::string _outputPrefix;
111 std::string _timestampString;
112 uint32_t _frameCount;
113 uint8_t _numComponents;
114 uint8_t _numSitesTotal;
115 uint8_t _numSphereTypes;
116 std::vector<uint64_t> _seekTable;
117 std::vector<uint8_t> _numSitesPerComp;
118 std::vector<uint8_t> _nCompSitesOffset;
119 std::string _strSphereDataFilename;
120 uint8_t _bInitSphereData;
121 bool _bWriteControlPrepared;
122
123 long _fileCount;
124 uint32_t _numFramesPerFile;
125
126 uint16_t _mmpldversion;
127 uint32_t _numSeekEntries;
128 MMPLD_Vertex_type _vertex_type;
129 MMPLD_Color_type _color_type;
130 std::vector<float> _global_radius;
131 std::vector< std::array<uint8_t, 4> > _global_rgba;
132 std::vector< std::array<float, 2> > _global_intensity_range;
133
134
135#ifdef ENABLE_MPI
136 MPI_File _mpifh;
137 MPI_Info_object _mpiinfo;
138#endif
139};
140
142{
143public:
145 MmpldWriterSimpleSphere(uint64_t startTimestep, uint64_t writeFrequency, uint64_t stopTimestep, uint64_t numFramesPerFile,
146 std::string outputPrefix)
147 : MmpldWriter(startTimestep, writeFrequency, stopTimestep, numFramesPerFile, outputPrefix)
148 {
149// MmpldWriter::_numSphereTypes = &(MmpldWriter::_numComponents);
150 }
151 virtual ~~MmpldWriterSimpleSphere() {}
152
153 virtual void SetNumSphereTypes() {_numSphereTypes = _numComponents;}
154 virtual void CalcNumSpheresPerType(ParticleContainer* particleContainer, uint64_t* numSpheresPerType);
155 virtual bool GetSpherePos(float *spherePos, Molecule* mol, uint8_t& nSphereTypeIndex);
156};
157
159{
160public:
162 MmpldWriterMultiSphere(uint64_t startTimestep, uint64_t writeFrequency, uint64_t stopTimestep, uint64_t numFramesPerFile,
163 std::string outputPrefix)
164 : MmpldWriter(startTimestep, writeFrequency, stopTimestep, numFramesPerFile, outputPrefix)
165 {
166// MmpldWriter::_numSphereTypes = &(MmpldWriter::_numSitesTotal);
167 }
168 virtual ~~MmpldWriterMultiSphere() {}
169
170 virtual void SetNumSphereTypes() {_numSphereTypes = _numSitesTotal;}
171 virtual void CalcNumSpheresPerType(ParticleContainer* particleContainer, uint64_t* numSpheresPerType);
172 virtual bool GetSpherePos(float *spherePos, Molecule* mol, uint8_t& nSphereTypeIndex);
173};
174
175#endif /* MMPLDWRITER_H_ */
handle boundary region and multiple processes
Definition: DomainDecompBase.h:51
This class is used to read in the phasespace and to handle macroscopic values.
Definition: Domain.h:47
FullMolecule modeled as LJ sphere with point polarities.
Definition: FullMolecule.h:18
MPI Info object implements functionalities to handle MPI Info objects, e.g. initialize from a XML fil...
Definition: MPI_Info_object.h:14
Definition: MmpldWriter.h:159
Definition: MmpldWriter.h:142
Output plugin to generate a MegaMolâ„¢ Particle List Data file (*.mmpld).
Definition: MmpldWriter.h:41
void finish(ParticleContainer *particleContainer, DomainDecompBase *domainDecomp, Domain *domain)
Method finish will be called at the end of the simulation.
Definition: MmpldWriter.cpp:324
uint64_t _startTimestep
Definition: MmpldWriter.h:104
void init(ParticleContainer *particleContainer, DomainDecompBase *domainDecomp, Domain *domain)
Method init will be called at the begin of the simulation.
Definition: MmpldWriter.cpp:146
void readXML(XMLfileUnits &xmlconfig)
Method readXML will be called once for each plugin section in the input file.
Definition: MmpldWriter.cpp:66
uint64_t _writeFrequency
Definition: MmpldWriter.h:106
uint64_t _stopTimestep
Definition: MmpldWriter.h:108
void endStep(ParticleContainer *particleContainer, DomainDecompBase *domainDecomp, Domain *domain, unsigned long simstep)
Method endStep will be called at the end of each time step.
Definition: MmpldWriter.cpp:308
std::string getPluginName()
return the name of the plugin
Definition: MmpldWriter.h:78
This Interface is used to get access to particles and pairs of particles.
Definition: ParticleContainer.h:69
The PluginBase class provides the interface for any kind of output/plugin classes - called "(output) ...
Definition: PluginBase.h:47
XML file with unit attributes abstraction.
Definition: xmlfileUnits.h:25
::xsd::cxx::tree::string< char, simple_type > string
C++ type corresponding to the string XML Schema built-in type.
Definition: vtk-punstructured.h:270