OpenLB 1.8.1
Loading...
Searching...
No Matches
vtkWriter.hh
Go to the documentation of this file.
1/* This file is part of the OpenLB library
2 *
3 * Copyright (C) 2013-2016 Thomas Henn, Mathias J. Krause
4 * 2016-2017 Albert Mink, Maximilian Gaedtke, Markus Morhard, Mathias J. Krause
5 * 2021 Nicolas Hafen, Mathias J. Krause
6 * E-mail contact: info@openlb.net
7 * The most recent release of OpenLB can be downloaded at
8 * <http://www.openlb.net/>
9 *
10 * This program is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU General Public License
12 * as published by the Free Software Foundation; either version 2
13 * of the License, or (at your option) any later version.
14 *
15 * This program is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU General Public License for more details.
19 *
20 * You should have received a copy of the GNU General Public
21 * License along with this program; if not, write to the Free
22 * Software Foundation, Inc., 51 Franklin Street, Fifth Floor,
23 * Boston, MA 02110-1301, USA.
24 */
25
26#ifndef VTK_WRITER_HH
27#define VTK_WRITER_HH
28
29
30
31namespace olb {
32
33
34template<typename T, typename FUNCTOR, vtkType VTKTYPE>
35VTKwriter<T,FUNCTOR,VTKTYPE>::VTKwriter( const std::string & name, bool binary, bool compress )
36 : clout( std::cout, "VTKwriter" ), _createFile(false), _name(name), _binary(binary), _compress(compress)
37{}
38
39
40template<typename T, typename FUNCTOR, vtkType VTKTYPE>
42{
43 _pointerVec.push_back(&f);
44}
45
46template<typename T, typename FUNCTOR, vtkType VTKTYPE>
47void VTKwriter<T,FUNCTOR,VTKTYPE>::addFunctor(FUNCTOR& f, const std::string& functorName)
48{
49 f.getName() = functorName;
50 _pointerVec.push_back(&f);
51}
52
53
54template<typename T, typename FUNCTOR, vtkType VTKTYPE>
55template<bool vtmAsMaster>
56void VTKwriter<T,FUNCTOR,VTKTYPE>::writeVTM( int iT, int rankSize, std::string fileExtension, std::string nameCollection)
57{
58 //Evaluate vtmAsMaster
59 std::string nameVTM;
60 std::string pathPrefixFileData;
61 if constexpr(!vtmAsMaster){
62 //Name of vtm
63 nameVTM = "data/" + createFileName(nameCollection, iT ) + ".vtm";
64 //path prefix
65 pathPrefixFileData = "";
66 //Create PVD
67 std::string fullNamePVD = singleton::directories().getVtkOutDir()
68 + createFileName(nameCollection) + ".pvd";
69 dataPVD(iT, fullNamePVD, nameVTM);
70 } else {
71 nameVTM = createFileName(nameCollection, iT) + ".vtm";
72 pathPrefixFileData = "data/";
73 }
74 //Create full name VTM
75 std::string fullNameVTM = singleton::directories().getVtkOutDir() + nameVTM;
76 //Write preamble VTM
77 preambleVTM(fullNameVTM);
78 //Loop over all cuboids (globIcs)
79 for (int iC = 0; iC < rankSize; iC++) {
80 std::string nameFileData = pathPrefixFileData + createFileName(nameCollection, iT, iC) + fileExtension;
81 // puts name of .vti piece to a .pvd file [fullNameVTM]
82 dataVTM(iC, fullNameVTM, nameFileData);
83 // adds a namePiece to master.pvd file.
84 }
85 // To do so we overwrite closePVD() and add new entry.
86 closeVTM(fullNameVTM); // timestep
87}
88
89
90template<typename T, typename FUNCTOR, vtkType VTKTYPE>
91template<bool vtmAsMaster, typename F>
92void VTKwriter<T,FUNCTOR,VTKTYPE>::write(int iT, std::string nameCollection, F applyFunctors)
93{
94 std::string fileExtension;
95
96 int rank = 0;
97 int rankSize;
98#ifdef PARALLEL_MODE_MPI
99 rank = singleton::mpi().getRank();
100#endif
101
102 if ( _pointerVec.empty() ) {
103 clout << "Error: Did you add a Functor ?" << std::endl;
104 } else {
105
106 //PARALLEL FUNCTORTYPE
107 if constexpr (parallel) {
108
109 //VTI
110 if constexpr(VTKTYPE==VTI){
111 fileExtension = ".vti";
112
113 //Retrieve cuboidDecomposition from first functor (TODO: include into common functortype)
114 auto& superStructure = dynamic_cast<SuperF<D,T>*>(_pointerVec[0])->getSuperStructure();
115 auto const& cGeometry = superStructure.getCuboidDecomposition();
116 LoadBalancer<T> const& loadBalancer = superStructure.getLoadBalancer();
117 rankSize = loadBalancer.getRankSize();
118 //Ensure no gaps between vti files (TODO: still necessary?)
119 superStructure.communicate();
120 //Retrieve delta
121 const T delta = cGeometry.getMotherCuboid().getDeltaR();
122
123 //Cor each cuboid
124 for (int iCloc = 0; iCloc < loadBalancer.size(); iCloc++) {
125
126 //Retrieve global IC
127 int globIC = loadBalancer.glob(iCloc);
128 //Retrieve BlockGeometry
129 auto blockGeometry = cGeometry.get(globIC);
130 //Retrieve extent
131 Vector<int,D> extent0(-1);
132 Vector<int,D> extent1( blockGeometry.getExtent() );
133
134 //Name of Data file
135 std::string fullNameFileData = singleton::directories().getVtkOutDir()
136 + "data/" + createFileName(nameCollection, iT, globIC ) + fileExtension;
137
138 // get dimension/extent for each cuboid
139 LatticeR<D> originCuboid(0.);
140 auto originPhysR = cGeometry.getPhysR(originCuboid.withPrefix(globIC));
141 //Write preamble VTU
142 preambleVTI(fullNameFileData, extent0, extent1, originPhysR.data(), delta);
143 //Loop over functors or call individual one
144 applyFunctors(fullNameFileData, extent1, globIC);
145 //Close piece and vti
146 closePiece(fullNameFileData);
147 closeVTI(fullNameFileData);
148 }
149
150 //VTU
151 } else if constexpr(VTKTYPE==VTU){
152 fileExtension = ".vtu";
153
154 //Retrive loadBalancer from first functor (TODO: include into common functortype)
155 LoadBalancer<T> const& loadBalancer = _pointerVec[0]->getLoadBalancer();
156 rankSize = loadBalancer.getRankSize();
157
158 //Cor each cuboid
159 for (int iCloc = 0; iCloc < loadBalancer.size(); iCloc++) {
160 //Retrieve global IC
161 int globIC = loadBalancer.glob(iCloc);
162 // Retrieve container size (number of particles) (TODO: include into common functortype)
163 Vector<int,1> extent1 = _pointerVec[0]->getContainerF(iCloc).getContainerSize();
164 //Name of VTU
165 std::string fullNameFileData = singleton::directories().getVtkOutDir()
166 + "data/" + createFileName(nameCollection, iT, globIC ) + fileExtension;
167 //Write preamble VTU
168 preambleVTU(fullNameFileData,extent1);
169 //Loop over functors or call individual one
170 applyFunctors(fullNameFileData, extent1, globIC);
171 // Write celldata (connectivity, offset, types)
172 this->cellDataVTU(fullNameFileData, extent1);
173 // Write points
174 this->dataArrayPoints(fullNameFileData, extent1, globIC);
175 // CloseVTU
176 closeVTU(fullNameFileData);
177 }
178
179 //VTP
180 } else {
181 fileExtension = ".vtp";
182 clout << "Error: VTP type not implemented yet" << std::endl;
183 }
184
185 //Write VTM (only master)
186 if (rank == 0) { writeVTM<vtmAsMaster>(iT,rankSize,fileExtension,nameCollection); }
187
188 //UNPARRALLIZED FUNCTORTYPE (for now only relevant for VTU)
189 } else {
190 fileExtension = ".vtu";
191
192 //if master
193 if (rank == 0) {
194 // Full name pvd
195 std::string fullNamePVD = singleton::directories().getVtkOutDir()
196 + createFileName(nameCollection) + ".pvd";
197 // Retrieve container size (number of particles)(TODO: include into common functortype)
198 Vector<int,1> extent1(_pointerVec[0]->getContainerSize());
199 // Name VTU
200 std::string nameFileData = "data/" + createFileName(nameCollection, iT ) + fileExtension;
201 // Create PVD
202 dataPVD(iT, fullNamePVD, nameFileData);
203 // Full Name VTU
204 std::string fullNameFileData = singleton::directories().getVtkOutDir() + nameFileData;
205 // Write preample VTU
206 preambleVTU(fullNameFileData, extent1);
207 //Loop over functors or call individual one
208 applyFunctors(fullNameFileData, extent1);
209 // Write celldata (connectivity, offset, types)
210 this->cellDataVTU(fullNameFileData, extent1);
211 // Write points
212 this->dataArrayPoints(fullNameFileData, extent1);
213 // CloseVTU
214 closeVTU(fullNameFileData);
215 } // master only
216 } // if parralized type
217 } //if ( _pointerVec.empty() )
218}
219
220
221template<typename T, typename FUNCTOR, vtkType VTKTYPE>
223{
224 //VTI
225 if constexpr(VTKTYPE==VTI){
226 //Define functors to be applied
227 auto applyFunctors = [&](std::string fullNameFileData, Vector<int,D> extent1, int globIC=0 )
228 {
229 for (int iF=0; iF<_pointerVec.size(); ++iF) {
230 this->dataArraySingleFunctor(fullNameFileData, *_pointerVec[iF], extent1, globIC);
231 }
232 };
233 //Call write
234 write<false>(iT,_name,applyFunctors);
235
236 //VTU
237 } else if constexpr(VTKTYPE==VTU){
238 //Define functors to be applied
239 auto applyFunctors = [&](std::string fullNameFileData, Vector<int,1> extent1, int globIC=0 )
240 {
241 for (std::size_t iF=1; iF<_pointerVec.size(); ++iF) {
242 this->dataArraySingleFunctor(fullNameFileData, *_pointerVec[iF], extent1, globIC);
243 }
244 };
245 //Call write
246 write<false>(iT,_name, applyFunctors);
247 }
248}
249
250
251template<typename T, typename FUNCTOR, vtkType VTKTYPE>
253{
254 //VTI
255 if constexpr(VTKTYPE==VTI){
256 //Define functor
257 auto applyFunctors = [&](std::string fullNameFileData, Vector<int,D> extent1, int globIC=0 ){
258 this->dataArraySingleFunctor(fullNameFileData, f, extent1, globIC); };
259 //Call write
260 write<true>(iT,f.getName(),applyFunctors);
261
262 //VTU
263 } else if constexpr(VTKTYPE==VTU){
264 //Define functor
265 auto applyFunctors = [&](std::string fullNameFileData, Vector<int,1> extent1, int globIC=0 ){
266 this->dataArraySingleFunctor(fullNameFileData, f, extent1, globIC); };
267 //Call write
268 write<true>(iT,f.getName(),applyFunctors);
269 }
270
271}
272
273template<typename T, typename FUNCTOR, vtkType VTKTYPE>
274void VTKwriter<T,FUNCTOR,VTKTYPE>::write(std::shared_ptr<FUNCTOR> ptr_f, int iT)
275{
276 write(*ptr_f, iT);
277}
278
279
280template<typename T, typename FUNCTOR, vtkType VTKTYPE>
282{
283 int rank = 0;
284#ifdef PARALLEL_MODE_MPI
285 rank = singleton::mpi().getRank();
286#endif
287 if (rank == 0) {
288 std::string fullNamePVD = singleton::directories().getVtkOutDir()
289 + createFileName(_name) + ".pvd";
290 preamblePVD(fullNamePVD);
291 closePVD(fullNamePVD);
292 _createFile = true;
293 }
294}
295
296template<typename T, typename FUNCTOR, vtkType VTKTYPE>
298 const std::string& fullName, Vector<int,1> extent1)
299{
300 std::ofstream fout(fullName, std::ios::trunc);
301 if (!fout) {
302 clout << "Error: could not open " << fullName << std::endl;
303 }
304 fout << "<?xml version=\"1.0\"?>\n";
305 fout << "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" ";
306 if (_compress) {
307 fout << "byte_order=\"LittleEndian\" compressor=\"vtkZLibDataCompressor\">\n";
308 }
309 else {
310 fout << "byte_order=\"LittleEndian\">\n";
311 }
312 fout << "<UnstructuredGrid>" << std::endl;
313 fout << "<Piece NumberOfPoints=\"" << extent1[0]
314 << "\" NumberOfCells=\"" << extent1[0] << "\">"
315 << std::endl;
316 fout << "<PointData Vectors=\"Particles\">" << std::endl;
317 fout.close();
318}
319
320template<typename T, typename FUNCTOR, vtkType VTKTYPE>
322 const std::string& fullNamePiece)
323{
324 std::ofstream fout(fullNamePiece.c_str(), std::ios::app);
325 if (!fout) {
326 clout << "Error: could not open " << fullNamePiece << std::endl;
327 }
328 fout << "</Piece>" << std::endl;
329 fout << "</UnstructuredGrid>\n";
330 fout << "</VTKFile>\n";
331 fout.close();
332}
333
334
335
336
337//TODO: make dimension insensitive
338template<typename T, typename FUNCTOR, vtkType VTKTYPE>
340 const std::string& fullName,
341 const LatticeR<D> extent0, const LatticeR<D> extent1,
342 PhysR<T,D> origin, T delta)
343{
344 std::ofstream fout(fullName, std::ios::trunc);
345 if (!fout) {
346 clout << "Error: could not open " << fullName << std::endl;
347 }
348 fout << "<?xml version=\"1.0\"?>\n";
349 fout << "<VTKFile type=\"ImageData\" version=\"0.1\" ";
350 if (_compress) {
351 fout << "byte_order=\"LittleEndian\" compressor=\"vtkZLibDataCompressor\">\n";
352 }
353 else {
354 fout << "byte_order=\"LittleEndian\">\n";
355 }
356
357 fout << "<ImageData WholeExtent=\""
358 << extent0[0] <<" "<< extent1[0];
359 for (unsigned iDim=1; iDim<D; ++iDim){
360 fout << " " << extent0[iDim] << " " << extent1[iDim];
361 }
362 fout << "\" Origin=\"" << origin[0];
363 for (unsigned iDim=1; iDim<D; ++iDim){
364 fout << " " << origin[iDim];
365 }
366 fout << "\" Spacing=\"" << delta << " " << delta << " " << delta << "\">\n";
367 fout << "<Piece Extent=\""
368 << extent0[0] <<" "<< extent1[0];
369 for (unsigned iDim=1; iDim<D; ++iDim){
370 fout << " " << extent0[iDim] <<" "<< extent1[iDim];
371 }
372 fout <<"\">\n";
373 fout << "<PointData>\n";
374 fout.close();
375}
376
377
378template<typename T, typename FUNCTOR, vtkType VTKTYPE>
380 const std::string& fullNamePiece)
381{
382 std::ofstream fout(fullNamePiece, std::ios::app );
383 if (!fout) {
384 clout << "Error: could not open " << fullNamePiece << std::endl;
385 }
386 fout << "</ImageData>\n";
387 fout << "</VTKFile>\n";
388 fout.close();
389}
390
391template<typename T, typename FUNCTOR, vtkType VTKTYPE>
392void VTKwriter<T,FUNCTOR,VTKTYPE>::closePiece(const std::string& fullNamePiece)
393{
394 std::ofstream fout(fullNamePiece, std::ios::app );
395 if (!fout) {
396 clout << "Error: could not open " << fullNamePiece << std::endl;
397 }
398 fout << "</PointData>\n";
399 fout << "</Piece>\n";
400 fout.close();
401}
402
403
404template<typename T, typename FUNCTOR, vtkType VTKTYPE>
406 const std::string& fullNamePVD)
407{
408 std::ofstream fout(fullNamePVD.c_str(), std::ios::trunc);
409 if (!fout) {
410 clout << "Error: could not open " << fullNamePVD << std::endl;
411 }
412
413 fout << "<?xml version=\"1.0\"?>\n";
414 fout << "<VTKFile type=\"Collection\" version=\"0.1\" "
415 << "byte_order=\"LittleEndian\">\n" << "<Collection>\n";
416 fout.close();
417}
418
419template<typename T, typename FUNCTOR, vtkType VTKTYPE>
421 const std::string& fullNamePVD)
422{
423 std::ofstream fout(fullNamePVD.c_str(), std::ios::app);
424 if (!fout) {
425 clout << "Error: could not open " << fullNamePVD << std::endl;
426 }
427 fout << "</Collection>\n";
428 fout << "</VTKFile>\n";
429 fout.close();
430}
431
432
433
434
435template<typename T, typename FUNCTOR, vtkType VTKTYPE>
436void VTKwriter<T,FUNCTOR,VTKTYPE>::preambleVTM(const std::string& fullNameVTM)
437{
438 std::ofstream fout(fullNameVTM, std::ios::trunc);
439 if (!fout) {
440 clout << "Error: could not open " << fullNameVTM << std::endl;
441 }
442 fout << "<?xml version=\"1.0\"?>\n";
443 fout << "<VTKFile type=\"vtkMultiBlockDataSet\" version=\"1.0\" "
444 << "byte_order=\"LittleEndian\">\n"
445 << "<vtkMultiBlockDataSet>\n" ;
446 fout.close();
447}
448
449template<typename T, typename FUNCTOR, vtkType VTKTYPE>
450void VTKwriter<T,FUNCTOR,VTKTYPE>::closeVTM(const std::string& fullNameVTM)
451{
452 std::ofstream fout(fullNameVTM, std::ios::app );
453 if (!fout) {
454 clout << "Error: could not open " << fullNameVTM << std::endl;
455 }
456 fout << "</vtkMultiBlockDataSet>\n";
457 fout << "</VTKFile>\n";
458 fout.close();
459}
460
461
462template<typename T, typename FUNCTOR, vtkType VTKTYPE>
463void VTKwriter<T,FUNCTOR,VTKTYPE>::dataVTM(int iC, const std::string& fullNameVTM,
464 const std::string& namePiece)
465{
466 std::ofstream fout(fullNameVTM, std::ios::app);
467 if (!fout) {
468 clout << "Error: could not open " << fullNameVTM << std::endl;
469 }
470 fout << "<Block index=\"" << iC << "\" >\n";
471 fout << "<DataSet index= \"0\" " << "file=\"" << namePiece << "\">\n"
472 << "</DataSet>\n";
473 fout << "</Block>\n";
474 fout.close();
475}
476
477template<typename T, typename FUNCTOR, vtkType VTKTYPE>
479 const std::string& fullNamePVD, const std::string& namePiece)
480{
481 std::ofstream fout(fullNamePVD.c_str(),
482 std::ios::in | std::ios::out | std::ios::ate);
483 if (fout) {
484 fout.seekp(-25, std::ios::end); // jump -25 form the end of file to overwrite closePVD
485
486 fout << "<DataSet timestep=\"" << iT << "\" "
487 << "group=\"\" part=\"\" "
488 << "file=\"" << namePiece << "\"/>\n";
489 fout.close();
490 closePVD(fullNamePVD);
491 }
492 else {
493 clout << "Error: could not open " << fullNamePVD << std::endl;
494 }
495}
496
497template<typename T, typename FUNCTOR, vtkType VTKTYPE>
498template<unsigned sourceDim>
500 const std::string& fullName,
501 Vector<int,sourceDim> extent1, int iC )
502
503{
504 std::ofstream fout(fullName.c_str(), std::ios::app);
505 if (!fout) {
506 clout << "Error: could not open " << fullName << std::endl;
507 }
508 fout << "<Points>" << std::endl;
509 // Call dataArray with first functor (which should contain points)
510 this->dataArraySingleFunctor(fullName, *_pointerVec[0], extent1, iC);
511 fout << "</Points>" << std::endl;
512}
513
514
515
516template<typename T, typename FUNCTOR, vtkType VTKTYPE>
517void VTKwriter<T,FUNCTOR,VTKTYPE>::cellDataVTU( const std::string& fullName, Vector<int, 1> extent1 ){
518 std::ofstream fout(fullName.c_str(), std::ios::app);
519 if (!fout) {
520 clout << "Error: could not open " << fullName << std::endl;
521 }
522 fout << "</PointData>" << std::endl;
523 fout << "<CellData /> " << std::endl; //TODO: open tag missing?
524 fout << "<Cells>" << std::endl;
525 fout << "<DataArray type=\"Int32\" Name=\"connectivity\" format=\"ascii\">"
526 << std::endl;
527 int32_t i32 = 0;
528 for (int iTmp=0; iTmp<extent1[0]; ++iTmp){ fout << i32++ << " "; }
529 fout << "</DataArray>" << std::endl;
530 fout << "<DataArray type=\"Int32\" Name=\"offsets\" format=\"ascii\">"
531 << std::endl;
532 i32 = 1;
533 for (int iTmp=0; iTmp<extent1[0]; ++iTmp){ fout << i32++ << " "; }
534 fout << "</DataArray>" << std::endl;
535 fout << "<DataArray type=\"UInt8\" Name=\"types\" format=\"ascii\">"
536 << std::endl;
537 for (int iTmp=0; iTmp<extent1[0]; ++iTmp){ fout << 1 << " "; }
538 fout << "</DataArray>" << std::endl;
539 fout << "</Cells>" << std::endl;
540}
541
542
543
544template<typename T, typename FUNCTOR, vtkType VTKTYPE>
545template<unsigned sourceDim>
547 const std::string& fullName,
548 FUNCTOR& f,
549 Vector<int,sourceDim> extent1, int iC)
550{
551 std::ofstream fout(fullName.c_str(), std::ios::app);
552 if (!fout) {
553 clout << "Error: could not open " << fullName << std::endl;
554 }
555 fout << "<DataArray type=\"Float32\" Name=\"" << f.getName() << "\" NumberOfComponents=\"" << f.getTargetDim() << "\" ";
556 if (_compress || _binary) {
557 fout << "format=\"binary\" encoding=\"base64\">\n";
558 }
559 else {
560 fout << ">\n";
561 }
562
563 // Create functor input and output
564 int i[4] = {0}; //4 as maximum dimension
565 std::vector<T> evaluated(f.getTargetDim());
566
567 //Set ic if parallel
568 if constexpr(parallel){ i[0] = iC;}
569
570 //Create float array
571 size_t numberOfFloats;
572 std::unique_ptr<float[]> streamFloat;
573 int itter = 0;
574
575 //VTI
576 if constexpr(VTKTYPE==VTI){
577 //2D
578 if constexpr(D==2){
579 numberOfFloats = f.getTargetDim() * (extent1[0]+2) * (extent1[1]+2);
580 streamFloat = std::make_unique<float[]>(numberOfFloats);
581 for (i[2]=-1; i[2]<extent1[1]+1; ++i[2]) {
582 for (i[1]=-1; i[1]<extent1[0]+1; ++i[1]) {
583 f(evaluated.data(),i);
584 for (int iDim=0; iDim < f.getTargetDim(); ++iDim){
585 streamFloat[itter] = float( evaluated[iDim] );
586 ++itter;
587 }
588 }
589 }
590 //3D
591 } else if constexpr(D==3){
592 numberOfFloats = f.getTargetDim() * (extent1[0]+2) * (extent1[1]+2) * (extent1[2]+2);
593 streamFloat = std::make_unique<float[]>(numberOfFloats);
594 for (i[3]=-1; i[3]<extent1[2]+1; ++i[3]) {
595 for (i[2]=-1; i[2]<extent1[1]+1; ++i[2]) {
596 for (i[1]=-1; i[1]<extent1[0]+1; ++i[1]) {
597 f(evaluated.data(),i);
598 for (int iDim=0; iDim < f.getTargetDim(); ++iDim){
599 streamFloat[itter] = float( evaluated[iDim] );
600 ++itter;
601 }
602 }
603 }
604 }
605 } else {
606 clout << "Error: only 2D and 3D supportet" << std::endl;
607 }
608
609 //VTU
610 } else if constexpr(VTKTYPE==VTU){
611 numberOfFloats = f.getTargetDim() * extent1[0];
612 streamFloat = std::make_unique<float[]>(numberOfFloats);
613 //PARALLEL
614 if constexpr(parallel){
615 for (i[1]=0; i[1]<extent1[0]; ++i[1]) {
616 f(evaluated.data(),i);
617 for (int iDim=0; iDim < f.getTargetDim(); ++iDim){
618 streamFloat[itter] = float( evaluated[iDim] );
619 ++itter;
620 }
621 }
622 //NONPARALLEL
623 } else {
624 for (i[0]=0; i[0]<extent1[0]; ++i[0]) {
625 f(evaluated.data(),i);
626 for (int iDim=0; iDim < f.getTargetDim(); ++iDim){
627 streamFloat[itter] = float( evaluated[iDim] );
628 ++itter;
629 }
630 }
631 }
632 } else {
633 clout << "Error: VTP format not implemented yet" << std::endl;
634 }
635
636 //Define binarSize from number of floats
637 uint32_t binarySize = static_cast<uint32_t>( numberOfFloats*sizeof(float) );
638
639 //Write data (TODO: repair for non parallised vtu version)
640 if (_compress) {
641 // char buffer for functor data
642 const unsigned char* charData = reinterpret_cast<unsigned char*>(streamFloat.get());
643 // buffer for compression
644 std::unique_ptr<unsigned char[]> comprData(new unsigned char[ binarySize ]); // stack may be too small
645
646 // compress data (not yet decoded as base64) by zlib
647 uLongf sizeCompr = compressBound(binarySize);
648 compress2( comprData.get(), &sizeCompr, charData, binarySize, -1);
649
650 // encode prefix to base64 documented in http://www.earthmodels.org/software/vtk-and-paraview/vtk-file-formats
651 Base64Encoder<uint32_t> prefixEncoder(fout, 4);
652 uint32_t prefix[4] = {1,binarySize,binarySize,static_cast<uint32_t>(sizeCompr)};
653 prefixEncoder.encode(prefix, 4);
654
655 // encode compressed data to base64
656 Base64Encoder<unsigned char> dataEncoder( fout, sizeCompr );
657 dataEncoder.encode(comprData.get(), sizeCompr);
658 } else if (_binary) {
659 // encode prefix to base64 documented in http://www.earthmodels.org/software/vtk-and-paraview/vtk-file-formats
660 Base64Encoder<uint32_t> prefixEncoder(fout, 1);
661 prefixEncoder.encode(&binarySize, 1);
662 // write numbers from functor
663 Base64Encoder<float> dataEncoder(fout, numberOfFloats);
664 dataEncoder.encode(streamFloat.get(),numberOfFloats);
665 } else {
666 for ( size_t iOut = 0; iOut < numberOfFloats; ++iOut ) {
667 fout << streamFloat[iOut] << " ";
668 }
669 }
670 fout.close();
671
672 std::ofstream ffout( fullName, std::ios::out | std::ios::app );
673 if (!ffout) {
674 clout << "Error: could not open " << fullName << std::endl;
675 }
676 ffout << "\n</DataArray>\n";
677 ffout.close();
678}
679
680
681} // namespace OLB
682
683#endif
void encode(const T *data, size_t length)
Definition base64.hh:59
Base class for all LoadBalancer.
Definition vtiWriter.h:42
int getRankSize() const
int glob(int loc) const
void closeVTI(const std::string &fullNamePiece)
performes </ImageData> and </VTKFile>
Definition vtkWriter.hh:379
void dataArraySingleFunctor(const std::string &fullName, FUNCTOR &f, Vector< int, sourceDim > extent1, int iC=0)
writes functors stored at pointerVec
Definition vtkWriter.hh:546
void closeVTU(const std::string &fullNamePiece)
performes </ImageData> and </VTKFile>
Definition vtkWriter.hh:321
void preambleVTM(const std::string &fullNameVTM)
performes <VTKFile ...> and <Collection>
Definition vtkWriter.hh:436
VTKwriter(const std::string &name, bool binary=true, bool compress=true)
constructor
Definition vtkWriter.hh:35
void preamblePVD(const std::string &fullNamePVD)
performes <VTKFile ...> and <Collection>
Definition vtkWriter.hh:405
void createMasterFile()
have to be called before calling write(int iT=0), since it creates
Definition vtkWriter.hh:281
void closeVTM(const std::string &fullNameVTM)
performes </Collection> and </VTKFile>
Definition vtkWriter.hh:450
void dataVTM(int iC, const std::string &fullNameVTM, const std::string &namePiece)
performes <DataSet timestep= ... file=namePiece > used for linking vti into pvd files
Definition vtkWriter.hh:463
void dataArrayPoints(const std::string &fullName, Vector< int, sourceDim > extent1, int iC=0)
writes points necessary for VTU
Definition vtkWriter.hh:499
void write(int iT, std::string nameCollection, F applyFunctors)
to clear stored functors, not yet used due to lack of necessity
Definition vtkWriter.hh:92
void closePiece(const std::string &fullNamePiece)
performes </PointData> and </Piece>
Definition vtkWriter.hh:392
void preambleVTU(const std::string &fullName, Vector< int, 1 > extent1)
performes <VTKFile ...>, <ImageData ...> and <PieceExtent ...>
Definition vtkWriter.hh:297
void preambleVTI(const std::string &fullName, const LatticeR< D > extent0, const LatticeR< D > extent1, PhysR< T, D > origin, T delta)
performes <VTKFile ...>, <ImageData ...>, <PieceExtent ...> and <PointData ...>
Definition vtkWriter.hh:339
void closePVD(const std::string &fullNamePVD)
performes </Collection> and </VTKFile>
Definition vtkWriter.hh:420
void cellDataVTU(const std::string &fullName, Vector< int, 1 > extent1)
TODO: add description: connectivity, offsete, type of unscructured nodes.
Definition vtkWriter.hh:517
void addFunctor(FUNCTOR &f)
put functor to _pointerVec to simplify writing process of several functors
Definition vtkWriter.hh:41
void writeVTM(int iT, int rankSize, std::string fileExtension, std::string nameCollection)
wrapper for VTM file creation
Definition vtkWriter.hh:56
void dataPVD(int iT, const std::string &fullNamePVD, const std::string &namePiece)
performes <DataSet timestep= ... file=namePiece >
Definition vtkWriter.hh:478
Plain old scalar vector.
constexpr const T * data() const any_platform
Definition vector.h:172
constexpr Vector< T, D+1 > withPrefix(T prefix) const any_platform
Return new (prefix,this...) vector.
Definition vector.h:198
std::string getVtkOutDir() const
Definition singleton.h:103
int getRank() const
Returns the process ID.
MpiManager & mpi()
Directories & directories()
Definition singleton.h:162
Top level namespace for all of OpenLB.
std::string createFileName(std::string name)
for .pvd masterFile
Definition fileName.hh:34
std::conditional_t< D==2, SuperF2D< T, U >, SuperF3D< T, U > > SuperF
Definition aliases.h:187
@ VTU
Definition vtkWriter.h:52
@ VTI
Definition vtkWriter.h:52