1 /* ----------------------------------------------------------------------
2     This is the
3 
4     ██╗     ██╗ ██████╗  ██████╗  ██████╗ ██╗  ██╗████████╗███████╗
5     ██║     ██║██╔════╝ ██╔════╝ ██╔════╝ ██║  ██║╚══██╔══╝██╔════╝
6     ██║     ██║██║  ███╗██║  ███╗██║  ███╗███████║   ██║   ███████╗
7     ██║     ██║██║   ██║██║   ██║██║   ██║██╔══██║   ██║   ╚════██║
8     ███████╗██║╚██████╔╝╚██████╔╝╚██████╔╝██║  ██║   ██║   ███████║
9     ╚══════╝╚═╝ ╚═════╝  ╚═════╝  ╚═════╝ ╚═╝  ╚═╝   ╚═╝   ╚══════╝®
10 
11     DEM simulation engine, released by
12     DCS Computing Gmbh, Linz, Austria
13     http://www.dcs-computing.com, office@dcs-computing.com
14 
15     LIGGGHTS® is part of CFDEM®project:
16     http://www.liggghts.com | http://www.cfdem.com
17 
18     Core developer and main author:
19     Christoph Kloss, christoph.kloss@dcs-computing.com
20 
21     LIGGGHTS® is open-source, distributed under the terms of the GNU Public
22     License, version 2 or later. It is distributed in the hope that it will
23     be useful, but WITHOUT ANY WARRANTY; without even the implied warranty
24     of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. You should have
25     received a copy of the GNU General Public License along with LIGGGHTS®.
26     If not, see http://www.gnu.org/licenses . See also top-level README
27     and LICENSE files.
28 
29     LIGGGHTS® and CFDEM® are registered trade marks of DCS Computing GmbH,
30     the producer of the LIGGGHTS® software and the CFDEM®coupling software
31     See http://www.cfdem.com/terms-trademark-policy for details.
32 
33 -------------------------------------------------------------------------
34     Contributing author and copyright for this file:
35     (if not contributing author is listed, this file has been contributed
36     by the core developer)
37     Arno Mayrhofer (DCS Computing GmbH, Linz)
38 
39     Copyright 2017-     DCS Computing GmbH, Linz
40 ------------------------------------------------------------------------- */
41 
42 #ifdef LAMMPS_VTK
43 
44 #include "dump_vtk.h"
45 #include "error.h"
46 #include "comm.h"
47 #include "update.h"
48 #include "universe.h"
49 #include <vtkUnstructuredGridWriter.h>
50 #include <vtkXMLUnstructuredGridWriter.h>
51 #include <vtkXMLPUnstructuredGridWriter.h>
52 #include <vtkPolyDataWriter.h>
53 #include <vtkXMLPolyDataWriter.h>
54 #include <vtkXMLPPolyDataWriter.h>
55 #include <vtkRectilinearGridWriter.h>
56 #include <vtkXMLRectilinearGridWriter.h>
57 #include <vtkXMLPRectilinearGridWriter.h>
58 #include <vtkXMLImageDataWriter.h>
59 #include <vtkXMLPImageDataWriter.h>
60 #include <sstream>
61 
62 namespace LAMMPS_NS
63 {
64 
DumpVTK(LAMMPS * lmp)65 DumpVTK::DumpVTK(LAMMPS *lmp) :
66     lmp_(lmp),
67     vtk_compressor_(VTK_COMP_NONE),
68     binary_(false)
69 {
70     filesuffixes[0] = (char*) ".vtk";
71     filesuffixes[1] = (char*) ".vtp";
72     filesuffixes[2] = (char*) ".vtu";
73     filesuffixes[3] = (char*) ".vti";
74     filesuffixes[4] = (char*) ".vtr";
75     filesuffixes[5] = (char*) ".vtm";
76     filesuffixes[6] = (char*) ".pvtp";
77     filesuffixes[7] = (char*) ".pvtu";
78     filesuffixes[8] = (char*) ".pvti";
79     filesuffixes[9] = (char*) ".pvtr";
80 }
81 
modify_param(int narg,char ** arg)82 int DumpVTK::modify_param(int narg, char **arg)
83 {
84     if (strcmp(arg[0],"binary") == 0) {
85         if (narg < 2)
86             lmp_->error->all(FLERR,"Illegal dump_modify command [binary]");
87         if (strcmp(arg[1],"yes") == 0)
88             binary_ = true;
89         else if (strcmp(arg[1],"no") == 0)
90             binary_ = false;
91         else
92             lmp_->error->all(FLERR,"Illegal dump_modify command [binary]");
93         return 2;
94     }
95 
96     if (strcmp(arg[0],"compressor") == 0)
97     {
98         if (narg < 2)
99             lmp_->error->all(FLERR,"Illegal dump_modify command [compressor]");
100 
101         if      (strcmp(arg[1],"zlib") == 0)
102             vtk_compressor_ = VTK_COMP_ZLIB;
103         else if (strcmp(arg[1],"lz4") == 0)
104 #if VTK_MAJOR_VERSION >= 8
105             vtk_compressor_ = VTK_COMP_LZ4;
106 #else
107             lmp_->error->all(FLERR, "Lz4 compressor is only available for VTK >= 8");
108 #endif
109         else if (strcmp(arg[1],"none") == 0)
110             vtk_compressor_ = VTK_COMP_NONE;
111         else
112             lmp_->error->all(FLERR,"Illegal dump_modify command [compressor]");
113 
114         // set binary on if compressor is used
115         if (vtk_compressor_ != VTK_COMP_NONE && !binary_)
116         {
117             lmp_->error->warning(FLERR, "Vtk dump will switch to binary writing as compressor is used");
118             binary_ = true;
119         }
120         return 2;
121     }
122 
123     return 0;
124 }
125 
setVtkWriterOptions(vtkSmartPointer<vtkDataWriter> writer)126 void DumpVTK::setVtkWriterOptions(vtkSmartPointer<vtkDataWriter> writer)
127 {
128     if (vtk_compressor_ != VTK_COMP_NONE && lmp_->comm->me == 0)
129         lmp_->error->warning(FLERR, "Vtk compressor enabled but data format does not support compression. To avoid this message do not use the *.vtk file ending");
130 
131     if (binary_)
132         writer->SetFileTypeToBinary();
133     else
134         writer->SetFileTypeToASCII();
135 }
136 
setVtkWriterOptions(vtkSmartPointer<vtkXMLWriter> writer)137 void DumpVTK::setVtkWriterOptions(vtkSmartPointer<vtkXMLWriter> writer)
138 {
139     if (binary_)
140         writer->SetDataModeToBinary();
141     else
142         writer->SetDataModeToAscii();
143 
144     switch (vtk_compressor_)
145     {
146     case VTK_COMP_ZLIB:
147         writer->SetCompressorTypeToZLib();
148         break;
149 #if VTK_MAJOR_VERSION >= 8
150     case VTK_COMP_LZ4:
151         writer->SetCompressorTypeToLZ4();
152         break;
153 #endif
154     default:
155         writer->SetCompressorTypeToNone();
156         break;
157     }
158 }
159 
write_vtp(vtkSmartPointer<vtkDataObject> data,const int vtk_file_format,const char * const filename)160 void DumpVTK::write_vtp(vtkSmartPointer<vtkDataObject> data, const int vtk_file_format, const char * const filename)
161 {
162     if (vtk_file_format == VTK_FILE_FORMATS::PVTP)
163     {
164         vtkSmartPointer<vtkXMLPPolyDataWriter> pwriter = vtkSmartPointer<vtkXMLPPolyDataWriter>::New();
165         pwriter->SetFileName(filename);
166 
167         setVtkWriterOptions(vtkXMLWriter::SafeDownCast(pwriter));
168 
169 #if VTK_MAJOR_VERSION < 6
170         pwriter->SetInput(data);
171 #else
172         pwriter->SetInputData(data);
173 #endif
174 
175         pwriter->SetNumberOfPieces(lmp_->comm->nprocs);
176         pwriter->SetStartPiece(lmp_->comm->me);
177         pwriter->SetEndPiece(lmp_->comm->me);
178         pwriter->Write();
179     }
180     else if (vtk_file_format == VTK_FILE_FORMATS::VTP)
181     {
182         vtkSmartPointer<vtkXMLPolyDataWriter> writer = vtkSmartPointer<vtkXMLPolyDataWriter>::New();
183 
184         setVtkWriterOptions(vtkXMLWriter::SafeDownCast(writer));
185 
186 #if VTK_MAJOR_VERSION < 6
187         writer->SetInput(data);
188 #else
189         writer->SetInputData(data);
190 #endif
191 
192         writer->SetFileName(filename);
193         writer->Write();
194     }
195     else
196         lmp_->error->all(FLERR, "Internal error");
197 }
198 
write_vtu(vtkSmartPointer<vtkDataObject> data,const int vtk_file_format,const char * const filename)199 void DumpVTK::write_vtu(vtkSmartPointer<vtkDataObject> data, const int vtk_file_format, const char * const filename)
200 {
201     if (vtk_file_format == VTK_FILE_FORMATS::PVTU)
202     {
203         vtkSmartPointer<vtkXMLPUnstructuredGridWriter> pwriter = vtkSmartPointer<vtkXMLPUnstructuredGridWriter>::New();
204         pwriter->SetFileName(filename);
205 
206         setVtkWriterOptions(vtkXMLWriter::SafeDownCast(pwriter));
207 
208 #if VTK_MAJOR_VERSION < 6
209         pwriter->SetInput(data);
210 #else
211         pwriter->SetInputData(data);
212 #endif
213 
214         pwriter->SetNumberOfPieces(lmp_->comm->nprocs);
215         pwriter->SetStartPiece(lmp_->comm->me);
216         pwriter->SetEndPiece(lmp_->comm->me);
217         pwriter->Write();
218     }
219     else if (vtk_file_format == VTK_FILE_FORMATS::VTU)
220     {
221         vtkSmartPointer<vtkXMLUnstructuredGridWriter> writer = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
222 
223         setVtkWriterOptions(vtkXMLWriter::SafeDownCast(writer));
224 
225 #if VTK_MAJOR_VERSION < 6
226         writer->SetInput(data);
227 #else
228         writer->SetInputData(data);
229 #endif
230 
231         writer->SetFileName(filename);
232         writer->Write();
233     }
234     else
235         lmp_->error->all(FLERR, "Internal error");
236 }
237 
write_vti(vtkSmartPointer<vtkAlgorithmOutput> data,const int vtk_file_format,const char * const filename)238 void DumpVTK::write_vti(vtkSmartPointer<vtkAlgorithmOutput> data, const int vtk_file_format, const char * const filename)
239 {
240     if (vtk_file_format == VTK_FILE_FORMATS::PVTI)
241     {
242         vtkSmartPointer<vtkXMLPImageDataWriter> pwriter = vtkSmartPointer<vtkXMLPImageDataWriter>::New();
243         pwriter->SetFileName(filename);
244 
245         setVtkWriterOptions(vtkXMLWriter::SafeDownCast(pwriter));
246 
247         pwriter->SetInputConnection(data);
248 
249         pwriter->SetNumberOfPieces(lmp_->comm->nprocs);
250         pwriter->SetStartPiece(lmp_->comm->me);
251         pwriter->SetEndPiece(lmp_->comm->me);
252         pwriter->SetWriteSummaryFile(lmp_->comm->me == 0 ? 1 : 0);
253         pwriter->Write();
254     }
255     else if (vtk_file_format == VTK_FILE_FORMATS::VTI)
256     {
257         vtkSmartPointer<vtkXMLImageDataWriter> writer = vtkSmartPointer<vtkXMLImageDataWriter>::New();
258 
259         setVtkWriterOptions(vtkXMLWriter::SafeDownCast(writer));
260 
261         writer->SetInputConnection(data);
262 
263         writer->SetFileName(filename);
264         writer->Write();
265     }
266     else
267         lmp_->error->all(FLERR, "Internal error");
268 }
269 
write_vtr(vtkSmartPointer<vtkDataObject> data,const int vtk_file_format,const char * const filename)270 void DumpVTK::write_vtr(vtkSmartPointer<vtkDataObject> data, const int vtk_file_format, const char * const filename)
271 {
272     if (vtk_file_format == VTK_FILE_FORMATS::PVTR)
273     {
274         vtkSmartPointer<vtkXMLPRectilinearGridWriter> pwriter = vtkSmartPointer<vtkXMLPRectilinearGridWriter>::New();
275         pwriter->SetFileName(filename);
276 
277         setVtkWriterOptions(vtkXMLWriter::SafeDownCast(pwriter));
278 
279         #if VTK_MAJOR_VERSION <= 5
280         pwriter->SetInputConnection(data->GetProducerPort());
281         #else
282         pwriter->SetInputData(data);
283         #endif
284 
285         pwriter->SetNumberOfPieces(lmp_->comm->nprocs);
286         pwriter->SetStartPiece(lmp_->comm->me);
287         pwriter->SetEndPiece(lmp_->comm->me);
288         pwriter->Write();
289     }
290     else if (vtk_file_format == VTK_FILE_FORMATS::VTR)
291     {
292         vtkSmartPointer<vtkXMLRectilinearGridWriter> writer = vtkSmartPointer<vtkXMLRectilinearGridWriter>::New();
293 
294         setVtkWriterOptions(vtkXMLWriter::SafeDownCast(writer));
295 
296         #if VTK_MAJOR_VERSION <= 5
297         writer->SetInputConnection(data->GetProducerPort());
298         #else
299         writer->SetInputData(data);
300         #endif
301 
302         writer->SetFileName(filename);
303         writer->Write();
304     }
305     else
306         lmp_->error->all(FLERR, "Internal error");
307 }
308 
write_vtk_unstructured_grid(vtkSmartPointer<vtkDataObject> data,const int vtk_file_format,const char * const filename,char * const label)309 void DumpVTK::write_vtk_unstructured_grid(vtkSmartPointer<vtkDataObject> data, const int vtk_file_format, const char * const filename, char * const label)
310 {
311     if (vtk_file_format == VTK_FILE_FORMATS::VTK)
312     {
313         vtkSmartPointer<vtkUnstructuredGridWriter> writer = vtkSmartPointer<vtkUnstructuredGridWriter>::New();
314 
315         if(label)
316             writer->SetHeader(label);
317         else
318             writer->SetHeader("Generated by LIGGGHTS");
319         setVtkWriterOptions(vtkDataWriter::SafeDownCast(writer));
320 
321         #if VTK_MAJOR_VERSION <= 5
322         writer->SetInputConnection(data->GetProducerPort());
323         #else
324         writer->SetInputData(data);
325         #endif
326 
327         writer->SetFileName(filename);
328         writer->Write();
329     }
330     else
331         lmp_->error->all(FLERR, "Internal error");
332 }
333 
write_vtk_rectilinear_grid(vtkSmartPointer<vtkDataObject> data,const int vtk_file_format,const char * const filename,char * const label)334 void DumpVTK::write_vtk_rectilinear_grid(vtkSmartPointer<vtkDataObject> data, const int vtk_file_format, const char * const filename, char * const label)
335 {
336     if (vtk_file_format == VTK_FILE_FORMATS::VTK)
337     {
338         vtkSmartPointer<vtkRectilinearGridWriter> writer = vtkSmartPointer<vtkRectilinearGridWriter>::New();
339 
340         if(label)
341             writer->SetHeader(label);
342         else
343             writer->SetHeader("Generated by LIGGGHTS");
344         setVtkWriterOptions(vtkDataWriter::SafeDownCast(writer));
345 
346         #if VTK_MAJOR_VERSION <= 5
347         writer->SetInputConnection(data->GetProducerPort());
348         #else
349         writer->SetInputData(data);
350         #endif
351 
352         writer->SetFileName(filename);
353         writer->Write();
354     }
355     else
356         lmp_->error->all(FLERR, "Internal error");
357 }
358 
write_vtk_poly(vtkSmartPointer<vtkDataObject> data,const int vtk_file_format,const char * const filename,char * const label)359 void DumpVTK::write_vtk_poly(vtkSmartPointer<vtkDataObject> data, const int vtk_file_format, const char * const filename, char * const label)
360 {
361     if (vtk_file_format == VTK_FILE_FORMATS::VTK)
362     {
363         vtkSmartPointer<vtkPolyDataWriter> writer = vtkSmartPointer<vtkPolyDataWriter>::New();
364 
365         if(label)
366             writer->SetHeader(label);
367         else
368             writer->SetHeader("Generated by LIGGGHTS");
369         setVtkWriterOptions(vtkDataWriter::SafeDownCast(writer));
370 
371         #if VTK_MAJOR_VERSION <= 5
372         writer->SetInputConnection(data->GetProducerPort());
373         #else
374         writer->SetInputData(data);
375         #endif
376 
377         writer->SetFileName(filename);
378         writer->Write();
379     }
380     else
381         lmp_->error->all(FLERR, "Internal error");
382 }
383 
getLocalController()384 vtkMPIController * DumpVTK::getLocalController()
385 {
386     vtkMPIController *vtkGlobalController = static_cast<vtkMPIController*>(vtkMultiProcessController::GetGlobalController());
387     if (!vtkGlobalController)
388         lmp_->error->all(FLERR, "Global VTK MPI Controller not found");
389 
390     if (lmp_->universe->existflag == 0)
391         return vtkGlobalController;
392     else
393     {
394         vtkMPIController *vtkLocalController = vtkGlobalController->PartitionController(lmp_->universe->iworld, 0);
395         if (!vtkLocalController)
396             lmp_->error->all(FLERR, "Local VTK MPI Controller not found");
397         return vtkLocalController;
398     }
399 }
400 
setFileCurrent(char * & filecurrent,char * const filename,const int multifile,const int padflag)401 void DumpVTK::setFileCurrent(char * &filecurrent, char * const filename, const int multifile, const int padflag)
402 {
403     delete [] filecurrent;
404     filecurrent = NULL;
405 
406     if (multifile == 0)
407     {
408         // contains no '*' -> simply copy
409         filecurrent = new char[strlen(filename) + 1];
410         strcpy(filecurrent, filename);
411     }
412     else
413     {
414         // contains '*' -> replace with time step
415         filecurrent = new char[strlen(filename) + 16];
416         char *ptr = strchr(filename,'*');
417         *ptr = '\0';
418         if (padflag == 0)
419         {
420             sprintf(filecurrent,"%s" BIGINT_FORMAT "%s",
421                     filename,lmp_->update->ntimestep,ptr+1);
422         }
423         else
424         {
425             char bif[8],pad[16];
426             strcpy(bif,BIGINT_FORMAT);
427             sprintf(pad,"%%s%%0%d%s%%s",padflag,&bif[1]);
428             sprintf(filecurrent,pad,filename,lmp_->update->ntimestep,ptr+1);
429         }
430         *ptr = '*';
431     }
432 }
433 
identify_file_type(char * const filename,std::list<int> & allowed_extensions,char * const style,int & multiproc,int & nclusterprocs,int & filewriter,int & fileproc,MPI_Comm & world,MPI_Comm & clustercomm)434 int DumpVTK::identify_file_type(char * const filename, std::list<int> &allowed_extensions, char * const style, int &multiproc, int &nclusterprocs, int &filewriter, int &fileproc, MPI_Comm &world, MPI_Comm &clustercomm)
435 {
436     // ensure no old % format is used
437     // this is set in dump.cpp
438     if (multiproc)
439         type_error("It is no longer allowed to enable parallel writing by setting the \% character, please see the documentation for help.", style, allowed_extensions);
440     // find last dot
441     char *suffix = strrchr(filename, '.');
442     if (strlen(suffix) == 5)
443     {
444         multiproc = 1;
445         nclusterprocs = 1;
446         filewriter = 1;
447         fileproc = lmp_->comm->me;
448         MPI_Comm_split(world,lmp_->comm->me,0,&clustercomm);
449         std::list<int>::iterator it = allowed_extensions.begin();
450         for (; it != allowed_extensions.end(); it++)
451         {
452             if (*it >= VTK_FILE_FORMATS::vtk_serial_file_types && strcmp(suffix, filesuffixes[*it]) == 0)
453                 return *it;
454         }
455         type_error("Could not find allowed filetype for parallel writing.", style, allowed_extensions);
456     }
457     else if (strlen(suffix) == 4)
458     {
459         std::list<int>::iterator it = allowed_extensions.begin();
460         for (; it != allowed_extensions.end(); it++)
461         {
462             if (*it < VTK_FILE_FORMATS::vtk_serial_file_types && strcmp(suffix, filesuffixes[*it]) == 0)
463                 return *it;
464         }
465         type_error("Could not find allowed filetype for serial writing.", style, allowed_extensions);
466     }
467     else
468         type_error("Could not find allowed filetype for writing of VTK file.", style, allowed_extensions);
469     return -1;
470 }
471 
type_error(std::string msg,char * const style,std::list<int> & allowed_extensions)472 void DumpVTK::type_error(std::string msg, char * const style, std::list<int> &allowed_extensions)
473 {
474     std::stringstream ss;
475     ss << "dump " << std::string(style) << ": " << msg << " Allowed file extensions for this dump style are:";
476     std::list<int>::iterator it = allowed_extensions.begin();
477     for (; it != allowed_extensions.end(); it++)
478         ss << " " << std::string(filesuffixes[*it]);
479     lmp_->error->all(FLERR, ss.str().c_str());
480 }
481 
482 }
483 
484 #endif // LAMMPS_VTK
485