File: dump_atom_vtk.cpp

package info (click to toggle)
liggghts 3.5.0%2Brepack1-10~bpo8%2B1
  • links: PTS, VCS
  • area: main
  • in suites: jessie-backports
  • size: 28,712 kB
  • sloc: cpp: 136,869; ansic: 1,446; python: 676; sh: 603; makefile: 358
file content (366 lines) | stat: -rw-r--r-- 11,606 bytes parent folder | download | duplicates (2)
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
/* ----------------------------------------------------------------------
    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:
    Anton Gladky(TU Bergakademie Freiberg), gladky.anton@gmail.com
------------------------------------------------------------------------- */

#ifdef LAMMPS_VTK
#include <string.h>
#include "dump_atom_vtk.h"
#include "atom.h"
#include "group.h"
#include "error.h"
#include "memory.h"
#include "comm.h"
#include <sstream>
#include <vtkVersion.h>
#ifndef VTK_MAJOR_VERSION
#include <vtkConfigure.h>
#endif
#include<vtkCellArray.h>
#include<vtkDoubleArray.h>
#include<vtkIntArray.h>
#include<vtkPoints.h>
#include<vtkPointData.h>
#include<vtkCellData.h>
#include<vtkSmartPointer.h>
#include<vtkUnstructuredGrid.h>
#include<vtkXMLUnstructuredGridWriter.h>

using namespace LAMMPS_NS;

/* ---------------------------------------------------------------------- */

DumpATOMVTK::DumpATOMVTK(LAMMPS *lmp, int narg, char **arg) : Dump(lmp, narg, arg)
{
  if (narg != 5) error->all(FLERR,"Illegal dump command");
  if (binary || multiproc) error->all(FLERR,"Invalid dump filename");

  sort_flag = 1;
  sortcol = 0;

  size_one = 17;

  char *str = (char *) "%d %g %g %g";
  int n = strlen(str) + 1;
  format_default = new char[n];
  strcpy(format_default,str);
}

/* ---------------------------------------------------------------------- */

void DumpATOMVTK::init_style()
{

}

/* ---------------------------------------------------------------------- */

void DumpATOMVTK::write_header(bigint /* n */)
{
}

/* ---------------------------------------------------------------------- */

int DumpATOMVTK::count()
{
  n_calls_ = 0;

  return Dump::count();
}

/* ---------------------------------------------------------------------- */

void DumpATOMVTK::pack(int *ids)
{
  int n = 0;
  int m = 0;

  int *tag = atom->tag;
  int *type = atom->type;
  int *mask = atom->mask;
  double *rmass = atom->rmass;
  double *mass = atom->mass;
  double **x = atom->x;
  double **v = atom->v;
  double **f = atom->f;
  double **o = atom->omega;
  int nlocal = atom->nlocal;

  for (int i = 0; i < nlocal; i++) {
    if (mask[i] & groupbit) {

      if (ids) ids[n++] = tag[i];

      double massTemp;

      if (rmass) {
        massTemp=rmass[i];
      } else {
        massTemp=mass[type[i]];
      }

      int me = comm->me;

      buf[m++] = x[i][0];
      buf[m++] = x[i][1];
      buf[m++] = x[i][2];

      buf[m++] = atom->radius[i];
      buf[m++] = massTemp;
      buf[m++] = static_cast<double>(tag[i]);
      buf[m++] = static_cast<double>(type[i]);

      buf[m++] = v[i][0];
      buf[m++] = v[i][1];
      buf[m++] = v[i][2];

      buf[m++] = o[i][0];
      buf[m++] = o[i][1];
      buf[m++] = o[i][2];

      buf[m++] = f[i][0];
      buf[m++] = f[i][1];
      buf[m++] = f[i][2];

      buf[m++] = static_cast<double>(me);
    }
  }

  setFileCurrent();
  tmpEXP.setFileName(filecurrent);
  return;
}

/* ---------------------------------------------------------------------- */

void DumpATOMVTK::write_data(int n, double *mybuf)
{
  if (comm->me != 0) return;
  n_calls_++;
  int m = 0;
  for (int i = 0; i < n; i++) {
    DumpATOMVTK::DataVTK tmpVTKDat(
      V3(mybuf[m+0], mybuf[m+1], mybuf[m+2]),
      mybuf[m+3], mybuf[m+4], static_cast<int> (mybuf[m+5]), static_cast<int> (mybuf[m+6]),
      V3(mybuf[m+7], mybuf[m+8], mybuf[m+9]),
      V3(mybuf[m+10], mybuf[m+11], mybuf[m+12]),
      V3(mybuf[m+13], mybuf[m+14], mybuf[m+15]),
      static_cast<int> (mybuf[m+16]));
    tmpEXP.add(tmpVTKDat);
    m += size_one;
  }

  if(n_calls_ == comm->nprocs) {
    tmpEXP.writeSER();
    tmpEXP.clear();
    delete [] filecurrent;
  }
}

/* ---------------------------------------------------------------------- */

DumpATOMVTK::DataVTK::DataVTK(V3 Pos, double Rad, double Mass, int Id, int Type,
  V3 VelL, V3 VelA, V3 Force, int proc) {
  _Pos = Pos;
  _Rad = Rad;
  _Mass = Mass;
  _Id = Id;
  _Type = Type;
  _VelL = VelL;
  _VelA = VelA;
  _Force = Force;
  _proc = proc;
}
/* ---------------------------------------------------------------------- */

std::string  DumpATOMVTK::DataVTK::serialize() {
  std::string tmp;
  std::ostringstream stringStream;

  stringStream <<_Pos[0]<<' '<<_Pos[1]<<' '<<_Pos[2]<<' '<<_Rad<<' '<<_Mass<<' '<<_Id
   <<' '<<_Type<<' '
   <<_VelL[0]<<' '<<_VelL[1]<<' '<<_VelL[2]<<' '
   <<_VelA[0]<<' '<<_VelA[1]<<' '<<_VelA[2]<<' '
   <<_Force[0]<<' '<<_Force[1]<<' '<<_Force[2]<<' '<<_proc<<'\n';

  tmp  = stringStream.str();
  return tmp;
}

/* ---------------------------------------------------------------------- */

void DumpATOMVTK::vtkExportData::add(DumpATOMVTK::DataVTK & d) {
  vtkData.push_back(d);
}

/* ---------------------------------------------------------------------- */

void DumpATOMVTK::vtkExportData::setFileName(const char * fileName) {
  _fileName = fileName;
  _setFileName = true;
}

/* ---------------------------------------------------------------------- */

DumpATOMVTK::vtkExportData::vtkExportData() {
  _setFileName=false;
}
/* ---------------------------------------------------------------------- */

int DumpATOMVTK::vtkExportData::size() {
  return vtkData.size();
}

/* ---------------------------------------------------------------------- */

void DumpATOMVTK::vtkExportData::writeSER() {

  vtkSmartPointer<vtkPoints>  spheresPos = vtkSmartPointer<vtkPoints>::New();
  vtkSmartPointer<vtkCellArray> spheresCells = vtkSmartPointer<vtkCellArray>::New();

  vtkSmartPointer<vtkDoubleArray> radii = vtkSmartPointer<vtkDoubleArray>::New();
  radii->SetNumberOfComponents(1);
  radii->SetName("radii");

  vtkSmartPointer<vtkDoubleArray> spheresMass = vtkSmartPointer<vtkDoubleArray>::New();
  spheresMass->SetNumberOfComponents(1);
  spheresMass->SetName("mass");

  vtkSmartPointer<vtkIntArray> spheresId = vtkSmartPointer<vtkIntArray>::New();
  spheresId->SetNumberOfComponents(1);
  spheresId->SetName("id");

  vtkSmartPointer<vtkIntArray> spheresType = vtkSmartPointer<vtkIntArray>::New();
  spheresType->SetNumberOfComponents(1);
  spheresType->SetName("type");

  vtkSmartPointer<vtkIntArray> spheresProc = vtkSmartPointer<vtkIntArray>::New();
  spheresProc->SetNumberOfComponents(1);
  spheresProc->SetName("proc");

  vtkSmartPointer<vtkDoubleArray> spheresVelL = vtkSmartPointer<vtkDoubleArray>::New();
  spheresVelL->SetNumberOfComponents(3);
  spheresVelL->SetName("velocity_lin");

  vtkSmartPointer<vtkDoubleArray> spheresVelA = vtkSmartPointer<vtkDoubleArray>::New();
  spheresVelA->SetNumberOfComponents(3);
  spheresVelA->SetName("velocity_ang");

  vtkSmartPointer<vtkDoubleArray> spheresForce = vtkSmartPointer<vtkDoubleArray>::New();
  spheresForce->SetNumberOfComponents(3);
  spheresForce->SetName("force");

  for (unsigned int i=0; i < vtkData.size(); i++) {
    vtkIdType pid[1];
    pid[0] = spheresPos->InsertNextPoint(vtkData[i]._Pos[0], vtkData[i]._Pos[1], vtkData[i]._Pos[2]);
    radii->InsertNextValue(vtkData[i]._Rad);

    double vv[3] = {vtkData[i]._VelL[0], vtkData[i]._VelL[1], vtkData[i]._VelL[2]};
    spheresVelL->InsertNextTupleValue(vv);

    double oo[3] = {vtkData[i]._VelA[0], vtkData[i]._VelA[1], vtkData[i]._VelA[2]};
    spheresVelA->InsertNextTupleValue(oo);

    double ff[3] = {vtkData[i]._Force[0], vtkData[i]._Force[1], vtkData[i]._Force[2]};
    spheresForce->InsertNextTupleValue(ff);

    spheresMass->InsertNextValue(vtkData[i]._Mass);

    spheresId->InsertNextValue(vtkData[i]._Id);
    spheresType->InsertNextValue(vtkData[i]._Type);
    spheresProc->InsertNextValue(vtkData[i]._proc);

    spheresCells->InsertNextCell(1,pid);
  }

  vtkSmartPointer<vtkUnstructuredGrid> spheresUg = vtkSmartPointer<vtkUnstructuredGrid>::New();

  spheresUg->SetPoints(spheresPos);
  spheresUg->SetCells(VTK_VERTEX, spheresCells);
  spheresUg->GetPointData()->AddArray(radii);
  spheresUg->GetPointData()->AddArray(spheresId);
  spheresUg->GetPointData()->AddArray(spheresType);
  spheresUg->GetPointData()->AddArray(spheresProc);
  spheresUg->GetPointData()->AddArray(spheresMass);
  spheresUg->GetPointData()->AddArray(spheresVelL);
  spheresUg->GetPointData()->AddArray(spheresVelA);
  spheresUg->GetPointData()->AddArray(spheresForce);

  vtkSmartPointer<vtkXMLUnstructuredGridWriter> writer = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
  writer->SetDataModeToAscii();
  #ifdef LAMMPS_VTK6
    writer->SetInputData(spheresUg);
  #else
    writer->SetInput(spheresUg);
  #endif
  writer->SetFileName(_fileName);
  writer->Write();
}

/* ---------------------------------------------------------------------- */

void DumpATOMVTK::vtkExportData::show() {
  for (unsigned int i=0; i < vtkData.size(); i++) {
    std::cerr << vtkData[i].serialize();
  }
}

/* ---------------------------------------------------------------------- */

void DumpATOMVTK::vtkExportData::clear() {
  vtkData.clear();
}

/* ---------------------------------------------------------------------- */

void DumpATOMVTK::setFileCurrent() {
  if (multifile == 0) filecurrent = filename;
  else {
    filecurrent = new char[strlen(filename) + 16];
    char *ptr = strchr(filename,'*');
    *ptr = '\0';
    if (padflag == 0)
      sprintf(filecurrent,"%s" BIGINT_FORMAT "%s",
              filename,update->ntimestep,ptr+1);
    else {
      char bif[8],pad[16];
      strcpy(bif,BIGINT_FORMAT);
      sprintf(pad,"%%s%%0%d_%%d%s%%s",padflag,&bif[1]);
      sprintf(filecurrent,pad,filename,comm->me,update->ntimestep,ptr+1);
    }
    *ptr = '*';
  }
}
#endif