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 435 436 437 438 439 440 441 442 443 444 445 446 447 448
|
/* ----------------------------------------------------------------------
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 not contributing author is listed, this file has been contributed
by the core developer)
Arno Mayrhofer (DCS Computing GmbH, Linz)
Copyright 2017- DCS Computing GmbH, Linz
------------------------------------------------------------------------- */
#ifdef LAMMPS_VTK
#include "dump_vtk.h"
#include "error.h"
#include "comm.h"
#include "update.h"
#include "universe.h"
#include <vtkUnstructuredGridWriter.h>
#include <vtkXMLUnstructuredGridWriter.h>
#include <vtkXMLPUnstructuredGridWriter.h>
#include <vtkPolyDataWriter.h>
#include <vtkXMLPolyDataWriter.h>
#include <vtkXMLPPolyDataWriter.h>
#include <vtkRectilinearGridWriter.h>
#include <vtkXMLRectilinearGridWriter.h>
#include <vtkXMLPRectilinearGridWriter.h>
#include <vtkXMLImageDataWriter.h>
#include <vtkXMLPImageDataWriter.h>
#include <sstream>
namespace LAMMPS_NS
{
DumpVTK::DumpVTK(LAMMPS *lmp) :
lmp_(lmp),
vtk_compressor_(VTK_COMP_NONE),
binary_(false)
{
filesuffixes[0] = (char*) ".vtk";
filesuffixes[1] = (char*) ".vtp";
filesuffixes[2] = (char*) ".vtu";
filesuffixes[3] = (char*) ".vti";
filesuffixes[4] = (char*) ".vtr";
filesuffixes[5] = (char*) ".vtm";
filesuffixes[6] = (char*) ".pvtp";
filesuffixes[7] = (char*) ".pvtu";
filesuffixes[8] = (char*) ".pvti";
filesuffixes[9] = (char*) ".pvtr";
}
int DumpVTK::modify_param(int narg, char **arg)
{
if (strcmp(arg[0],"binary") == 0) {
if (narg < 2)
lmp_->error->all(FLERR,"Illegal dump_modify command [binary]");
if (strcmp(arg[1],"yes") == 0)
binary_ = true;
else if (strcmp(arg[1],"no") == 0)
binary_ = false;
else
lmp_->error->all(FLERR,"Illegal dump_modify command [binary]");
return 2;
}
if (strcmp(arg[0],"compressor") == 0)
{
if (narg < 2)
lmp_->error->all(FLERR,"Illegal dump_modify command [compressor]");
if (strcmp(arg[1],"zlib") == 0)
vtk_compressor_ = VTK_COMP_ZLIB;
else if (strcmp(arg[1],"lz4") == 0)
#if VTK_MAJOR_VERSION >= 8
vtk_compressor_ = VTK_COMP_LZ4;
#else
lmp_->error->all(FLERR, "Lz4 compressor is only available for VTK >= 8");
#endif
else if (strcmp(arg[1],"none") == 0)
vtk_compressor_ = VTK_COMP_NONE;
else
lmp_->error->all(FLERR,"Illegal dump_modify command [compressor]");
// set binary on if compressor is used
if (vtk_compressor_ != VTK_COMP_NONE && !binary_)
{
lmp_->error->warning(FLERR, "Vtk dump will switch to binary writing as compressor is used");
binary_ = true;
}
return 2;
}
return 0;
}
void DumpVTK::setVtkWriterOptions(vtkSmartPointer<vtkDataWriter> writer)
{
if (vtk_compressor_ != VTK_COMP_NONE && lmp_->comm->me == 0)
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");
if (binary_)
writer->SetFileTypeToBinary();
else
writer->SetFileTypeToASCII();
}
void DumpVTK::setVtkWriterOptions(vtkSmartPointer<vtkXMLWriter> writer)
{
if (binary_)
writer->SetDataModeToBinary();
else
writer->SetDataModeToAscii();
switch (vtk_compressor_)
{
case VTK_COMP_ZLIB:
writer->SetCompressorTypeToZLib();
break;
#if VTK_MAJOR_VERSION >= 8
case VTK_COMP_LZ4:
writer->SetCompressorTypeToLZ4();
break;
#endif
default:
writer->SetCompressorTypeToNone();
break;
}
}
void DumpVTK::write_vtp(vtkSmartPointer<vtkDataObject> data, const int vtk_file_format, const char * const filename)
{
if (vtk_file_format == VTK_FILE_FORMATS::PVTP)
{
vtkSmartPointer<vtkXMLPPolyDataWriter> pwriter = vtkSmartPointer<vtkXMLPPolyDataWriter>::New();
pwriter->SetFileName(filename);
setVtkWriterOptions(vtkXMLWriter::SafeDownCast(pwriter));
pwriter->SetInputData(data);
pwriter->SetNumberOfPieces(lmp_->comm->nprocs);
pwriter->SetStartPiece(lmp_->comm->me);
pwriter->SetEndPiece(lmp_->comm->me);
pwriter->Write();
}
else if (vtk_file_format == VTK_FILE_FORMATS::VTP)
{
vtkSmartPointer<vtkXMLPolyDataWriter> writer = vtkSmartPointer<vtkXMLPolyDataWriter>::New();
setVtkWriterOptions(vtkXMLWriter::SafeDownCast(writer));
writer->SetInputData(data);
writer->SetFileName(filename);
writer->Write();
}
else
lmp_->error->all(FLERR, "Internal error");
}
void DumpVTK::write_vtu(vtkSmartPointer<vtkDataObject> data, const int vtk_file_format, const char * const filename)
{
if (vtk_file_format == VTK_FILE_FORMATS::PVTU)
{
vtkSmartPointer<vtkXMLPUnstructuredGridWriter> pwriter = vtkSmartPointer<vtkXMLPUnstructuredGridWriter>::New();
pwriter->SetFileName(filename);
setVtkWriterOptions(vtkXMLWriter::SafeDownCast(pwriter));
pwriter->SetInputData(data);
pwriter->SetNumberOfPieces(lmp_->comm->nprocs);
pwriter->SetStartPiece(lmp_->comm->me);
pwriter->SetEndPiece(lmp_->comm->me);
pwriter->Write();
}
else if (vtk_file_format == VTK_FILE_FORMATS::VTU)
{
vtkSmartPointer<vtkXMLUnstructuredGridWriter> writer = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
setVtkWriterOptions(vtkXMLWriter::SafeDownCast(writer));
writer->SetInputData(data);
writer->SetFileName(filename);
writer->Write();
}
else
lmp_->error->all(FLERR, "Internal error");
}
void DumpVTK::write_vti(vtkSmartPointer<vtkAlgorithmOutput> data, const int vtk_file_format, const char * const filename)
{
if (vtk_file_format == VTK_FILE_FORMATS::PVTI)
{
vtkSmartPointer<vtkXMLPImageDataWriter> pwriter = vtkSmartPointer<vtkXMLPImageDataWriter>::New();
pwriter->SetFileName(filename);
setVtkWriterOptions(vtkXMLWriter::SafeDownCast(pwriter));
pwriter->SetInputConnection(data);
pwriter->SetNumberOfPieces(lmp_->comm->nprocs);
pwriter->SetStartPiece(lmp_->comm->me);
pwriter->SetEndPiece(lmp_->comm->me);
pwriter->SetWriteSummaryFile(lmp_->comm->me == 0 ? 1 : 0);
pwriter->Write();
}
else if (vtk_file_format == VTK_FILE_FORMATS::VTI)
{
vtkSmartPointer<vtkXMLImageDataWriter> writer = vtkSmartPointer<vtkXMLImageDataWriter>::New();
setVtkWriterOptions(vtkXMLWriter::SafeDownCast(writer));
writer->SetInputConnection(data);
writer->SetFileName(filename);
writer->Write();
}
else
lmp_->error->all(FLERR, "Internal error");
}
void DumpVTK::write_vtr(vtkSmartPointer<vtkDataObject> data, const int vtk_file_format, const char * const filename)
{
if (vtk_file_format == VTK_FILE_FORMATS::PVTR)
{
vtkSmartPointer<vtkXMLPRectilinearGridWriter> pwriter = vtkSmartPointer<vtkXMLPRectilinearGridWriter>::New();
pwriter->SetFileName(filename);
setVtkWriterOptions(vtkXMLWriter::SafeDownCast(pwriter));
pwriter->SetInputData(data);
pwriter->SetNumberOfPieces(lmp_->comm->nprocs);
pwriter->SetStartPiece(lmp_->comm->me);
pwriter->SetEndPiece(lmp_->comm->me);
pwriter->Write();
}
else if (vtk_file_format == VTK_FILE_FORMATS::VTR)
{
vtkSmartPointer<vtkXMLRectilinearGridWriter> writer = vtkSmartPointer<vtkXMLRectilinearGridWriter>::New();
setVtkWriterOptions(vtkXMLWriter::SafeDownCast(writer));
writer->SetInputData(data);
writer->SetFileName(filename);
writer->Write();
}
else
lmp_->error->all(FLERR, "Internal error");
}
void DumpVTK::write_vtk_unstructured_grid(vtkSmartPointer<vtkDataObject> data, const int vtk_file_format, const char * const filename, char * const label)
{
if (vtk_file_format == VTK_FILE_FORMATS::VTK)
{
vtkSmartPointer<vtkUnstructuredGridWriter> writer = vtkSmartPointer<vtkUnstructuredGridWriter>::New();
if(label)
writer->SetHeader(label);
else
writer->SetHeader("Generated by LIGGGHTS");
setVtkWriterOptions(vtkDataWriter::SafeDownCast(writer));
writer->SetInputData(data);
writer->SetFileName(filename);
writer->Write();
}
else
lmp_->error->all(FLERR, "Internal error");
}
void DumpVTK::write_vtk_rectilinear_grid(vtkSmartPointer<vtkDataObject> data, const int vtk_file_format, const char * const filename, char * const label)
{
if (vtk_file_format == VTK_FILE_FORMATS::VTK)
{
vtkSmartPointer<vtkRectilinearGridWriter> writer = vtkSmartPointer<vtkRectilinearGridWriter>::New();
if(label)
writer->SetHeader(label);
else
writer->SetHeader("Generated by LIGGGHTS");
setVtkWriterOptions(vtkDataWriter::SafeDownCast(writer));
writer->SetInputData(data);
writer->SetFileName(filename);
writer->Write();
}
else
lmp_->error->all(FLERR, "Internal error");
}
void DumpVTK::write_vtk_poly(vtkSmartPointer<vtkDataObject> data, const int vtk_file_format, const char * const filename, char * const label)
{
if (vtk_file_format == VTK_FILE_FORMATS::VTK)
{
vtkSmartPointer<vtkPolyDataWriter> writer = vtkSmartPointer<vtkPolyDataWriter>::New();
if(label)
writer->SetHeader(label);
else
writer->SetHeader("Generated by LIGGGHTS");
setVtkWriterOptions(vtkDataWriter::SafeDownCast(writer));
writer->SetInputData(data);
writer->SetFileName(filename);
writer->Write();
}
else
lmp_->error->all(FLERR, "Internal error");
}
vtkMPIController * DumpVTK::getLocalController()
{
vtkMPIController *vtkGlobalController = static_cast<vtkMPIController*>(vtkMultiProcessController::GetGlobalController());
if (!vtkGlobalController)
lmp_->error->all(FLERR, "Global VTK MPI Controller not found");
if (lmp_->universe->existflag == 0)
return vtkGlobalController;
else
{
vtkMPIController *vtkLocalController = vtkGlobalController->PartitionController(lmp_->universe->iworld, 0);
if (!vtkLocalController)
lmp_->error->all(FLERR, "Local VTK MPI Controller not found");
return vtkLocalController;
}
}
void DumpVTK::setFileCurrent(char * &filecurrent, char * const filename, const int multifile, const int padflag)
{
delete [] filecurrent;
filecurrent = NULL;
if (multifile == 0)
{
// contains no '*' -> simply copy
filecurrent = new char[strlen(filename) + 1];
strcpy(filecurrent, filename);
}
else
{
// contains '*' -> replace with time step
filecurrent = new char[strlen(filename) + 16];
char *ptr = strchr(filename,'*');
*ptr = '\0';
if (padflag == 0)
{
sprintf(filecurrent,"%s" BIGINT_FORMAT "%s",
filename,lmp_->update->ntimestep,ptr+1);
}
else
{
char bif[8],pad[16];
strcpy(bif,BIGINT_FORMAT);
sprintf(pad,"%%s%%0%d%s%%s",padflag,&bif[1]);
sprintf(filecurrent,pad,filename,lmp_->update->ntimestep,ptr+1);
}
*ptr = '*';
}
}
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)
{
// ensure no old % format is used
// this is set in dump.cpp
if (multiproc)
type_error("It is no longer allowed to enable parallel writing by setting the \% character, please see the documentation for help.", style, allowed_extensions);
// find last dot
char *suffix = strrchr(filename, '.');
if (strlen(suffix) == 5)
{
multiproc = 1;
nclusterprocs = 1;
filewriter = 1;
fileproc = lmp_->comm->me;
MPI_Comm_split(world,lmp_->comm->me,0,&clustercomm);
std::list<int>::iterator it = allowed_extensions.begin();
for (; it != allowed_extensions.end(); it++)
{
if (*it >= VTK_FILE_FORMATS::vtk_serial_file_types && strcmp(suffix, filesuffixes[*it]) == 0)
return *it;
}
type_error("Could not find allowed filetype for parallel writing.", style, allowed_extensions);
}
else if (strlen(suffix) == 4)
{
std::list<int>::iterator it = allowed_extensions.begin();
for (; it != allowed_extensions.end(); it++)
{
if (*it < VTK_FILE_FORMATS::vtk_serial_file_types && strcmp(suffix, filesuffixes[*it]) == 0)
return *it;
}
type_error("Could not find allowed filetype for serial writing.", style, allowed_extensions);
}
else
type_error("Could not find allowed filetype for writing of VTK file.", style, allowed_extensions);
return -1;
}
void DumpVTK::type_error(std::string msg, char * const style, std::list<int> &allowed_extensions)
{
std::stringstream ss;
ss << "dump " << std::string(style) << ": " << msg << " Allowed file extensions for this dump style are:";
std::list<int>::iterator it = allowed_extensions.begin();
for (; it != allowed_extensions.end(); it++)
ss << " " << std::string(filesuffixes[*it]);
lmp_->error->all(FLERR, ss.str().c_str());
}
}
#endif // LAMMPS_VTK
|