ALL 0.9.3
A Loadbalacing Library
Loading...
Searching...
No Matches
ALL.hpp
Go to the documentation of this file.
1/*
2Copyright 2018-2020 Rene Halver, Forschungszentrum Juelich GmbH, Germany
3Copyright 2018-2020 Godehard Sutmann, Forschungszentrum Juelich GmbH, Germany
4
5Redistribution and use in source and binary forms, with or without modification,
6are permitted provided that the following conditions are met:
7
81. Redistributions of source code must retain the above copyright notice, this
9 list of conditions and the following disclaimer.
10
112. Redistributions in binary form must reproduce the above copyright notice,
12this list of conditions and the following disclaimer in the documentation and/or
13 other materials provided with the distribution.
14
153. Neither the name of the copyright holder nor the names of its contributors
16may be used to endorse or promote products derived from this software without
17 specific prior written permission.
18
19THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
20ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
21WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
22DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR
23ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
24(INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
25LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
26ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
27(INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
28SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
29*/
30
31#ifndef ALL_MAIN_HEADER_INC
32#define ALL_MAIN_HEADER_INC
33
35#include "ALL_LB.hpp"
36#include "ALL_Point.hpp"
37#include "ALL_Tensor.hpp"
38#ifdef ALL_VORONOI_ACTIVE
39#include "ALL_Voronoi.hpp"
40#endif
41#include "ALL_ForceBased.hpp"
42#include "ALL_Histogram.hpp"
43#include "ALL_Staggered.hpp"
44#include <algorithm>
45#include <iomanip>
46#include <memory>
47#include <numeric>
48#include <sstream>
49#include <string>
50#include <vector>
51
52#include <cerrno>
53#include <sys/stat.h>
54
55#ifdef ALL_VTK_OUTPUT
56#include <limits>
57#include <vtkCellArray.h>
58#include <vtkCellData.h>
59#include <vtkFloatArray.h>
60#include <vtkInformation.h>
61#include <vtkIntArray.h>
62#include <vtkMPIController.h>
63#include <vtkPolyhedron.h>
64#include <vtkProgrammableFilter.h>
65#include <vtkSmartPointer.h>
66#include <vtkUnstructuredGrid.h>
67#include <vtkVoxel.h>
68#include <vtkXMLPUnstructuredGridWriter.h>
69#include <vtkXMLUnstructuredGridWriter.h>
70#include <vtkUnstructuredGridWriter.h> //For vtkXMLUnstructuredGridWriter debug
71
72#ifdef VTK_CELL_ARRAY_V2
73#include <vtkNew.h>
74#endif
75#endif
76
77namespace ALL {
78
79#define ALL_ESTIMATION_LIMIT_DEFAULT 0
80
82enum LB_t : int {
86 TENSOR = 1,
89#ifdef ALL_VORONOI_ACTIVE
92#endif
99};
100
103
104template <class T, class W> class ALL {
105public:
109 : loadbalancing_step(0)
110#ifdef ALL_VTK_OUTPUT
111 ,
112 vtk_init(false)
113#endif
114 {
115 // determine correct MPI data type for template T
116 if (std::is_same<T, double>::value)
117 MPIDataTypeT = MPI_DOUBLE;
118 else if (std::is_same<T, float>::value)
119 MPIDataTypeT = MPI_FLOAT;
120 else if (std::is_same<T, int>::value)
121 MPIDataTypeT = MPI_INT;
122 else if (std::is_same<T, long>::value)
123 MPIDataTypeT = MPI_LONG;
124
125 // determine correct MPI data type for template W
126 if (std::is_same<W, double>::value)
127 MPIDataTypeW = MPI_DOUBLE;
128 else if (std::is_same<W, float>::value)
129 MPIDataTypeW = MPI_FLOAT;
130 else if (std::is_same<W, int>::value)
131 MPIDataTypeW = MPI_INT;
132 else if (std::is_same<W, long>::value)
133 MPIDataTypeW = MPI_LONG;
134 }
135
144 ALL(const LB_t m, const int d, const T g) : ALL() {
145 method = m;
146 switch (method) {
147 case LB_t::TENSOR:
148 balancer.reset(new Tensor_LB<T, W>(d, (W)0, g));
149 break;
150 case LB_t::STAGGERED:
151 balancer.reset(new Staggered_LB<T, W>(d, (W)0, g));
152 break;
153 case LB_t::FORCEBASED:
154 balancer.reset(new ForceBased_LB<T, W>(d, (W)0, g));
155 break;
156#ifdef ALL_VORONOI_ACTIVE
157 case LB_t::VORONOI:
158 balancer.reset(new Voronoi_LB<T, W>(d, (W)0, g));
159 break;
160#endif
161 case LB_t::HISTOGRAM:
162 balancer.reset(new Histogram_LB<T, W>(d, std::vector<W>(10), g));
163 break;
164 case LB_t::TENSOR_MAX:
165 balancer.reset(new TensorMax_LB<T, W>(d, (W)0, g));
166 break;
167 default:
168 throw InvalidArgumentException(__FILE__, __func__, __LINE__,
169 "Unknown type of loadbalancing passed.");
170 }
171 balancer->setDimension(d);
172 balancer->setGamma(g);
173 estimation_limit = ALL_ESTIMATION_LIMIT_DEFAULT;
174 }
175
185 ALL(const LB_t m, const int d, const std::vector<Point<T>> &inp, const T g)
186 : ALL(m, d, g) {
187 balancer->setVertices(inp);
188 calculate_outline();
189 }
190
192 ~ALL() = default;
193
197 void setVertices(const std::vector<Point<T>> &inp) {
198 int dim = inp.at(0).getDimension();
199 if (dim != balancer->getDimension())
201 __FILE__, __func__, __LINE__,
202 "Dimension of ALL::Points in input vector do not match dimension of "
203 "ALL "
204 "object.");
205 balancer->setVertices(inp);
206 calculate_outline();
207 }
208
213 void setCommunicator(const MPI_Comm comm) {
214 int comm_type;
215 MPI_Topo_test(comm, &comm_type);
216 if (comm_type == MPI_CART) {
217 communicator = comm;
218 balancer->setCommunicator(communicator);
219 // set procGridLoc and procGridDim
220 const int dimension = balancer->getDimension();
221 int location[dimension];
222 int size[dimension];
223 int periodicity[dimension];
224 MPI_Cart_get(communicator, dimension, size, periodicity, location);
225 procGridLoc.assign(location, location + dimension);
226 procGridDim.assign(size, size + dimension);
227 procGridSet = true;
228 } else {
229 int rank;
230 MPI_Comm_rank(comm, &rank);
231 if (procGridSet) {
232 // compute 1D index from the location in the process grid (using z-y-x
233 // ordering, as MPI_Dims_create)
234 int idx;
235 switch (balancer->getDimension()) {
236 case 3:
237 idx = procGridLoc.at(2) + procGridLoc.at(1) * procGridDim.at(2) +
238 procGridLoc.at(0) * procGridDim.at(2) * procGridDim.at(1);
239 break;
240 case 2:
241 idx = procGridLoc.at(1) + procGridLoc.at(0) * procGridDim.at(1);
242 break;
243 case 1:
244 idx = procGridLoc.at(0);
245 break;
246 default:
248 __FILE__, __func__, __LINE__,
249 "ALL cannot deal with more then three dimensions (or less then "
250 "one).");
251 break;
252 }
253
254 // create new temporary communicator with correct ranks
255 MPI_Comm temp_comm;
256 MPI_Comm_split(comm, 0, idx, &temp_comm);
257
258 std::vector<int> periods(3, 1);
259
260 // transform temporary communicator to cartesian communicator
261 MPI_Cart_create(temp_comm, balancer->getDimension(), procGridDim.data(),
262 periods.data(), 0, &communicator);
263 balancer->setCommunicator(communicator);
264 } else {
266 __FILE__, __func__, __LINE__,
267 "When using a non-cartesian communicator process grid parameters "
268 "required (not set).");
269 }
270 }
271 }
272
275 T getGamma() { return balancer->getGamma(); };
276
279 void setGamma(const T g) { balancer->setGamma(g); };
280
283 void setWork(const W work) {
284 // set new value for work (single value for whole domain)
285 balancer->setWork(work);
286 }
287
288 // method to set a multi-dimensional work for the local domain
290 void setWork(const std::vector<W> &work) { balancer->setWork(work); }
291
294 void setup() { balancer->setup(); }
295
301 std::vector<Point<T>> &balance() {
302 nVertices = balancer->getNVertices();
303 switch (method) {
304 case LB_t::TENSOR:
305 balancer->balance(loadbalancing_step);
306 calculate_outline();
307 break;
308 case LB_t::STAGGERED:
309
310 /* ToDo: Check if repetition is useful at all and
311 * change it to be included for all methods,
312 * not only STAGGERED */
313 /*
314 // estimation to improve speed of convergence
315 // changing vertices during estimation
316 std::vector<Point<T>> tmp_vertices = balancer->getVertices();
317 // saved original vertices
318 std::vector<Point<T>> old_vertices = balancer->getVertices();
319 // old temporary vertices
320 std::vector<Point<T>> old_tmp_vertices = balancer->getVertices();
321 // result temporary vertices
322 std::vector<Point<T>> result_vertices = balancer->getVertices();
323 // temporary work
324 W tmp_work = balancer->getWork().at(0);
325 // old temporary work
326 W old_tmp_work = balancer->getWork().at(0);
327
328 W max_work;
329 W sum_work;
330 double d_max;
331 int n_ranks;
332
333 // compute work density on local domain
334 double V = 1.0;
335 for (int i = 0; i < balancer->getDimension(); ++i)
336 V *= ( outline.at(1)[i] - outline.at(0)[i] );
337 double rho = workArray.at(0) / V;
338
339 // collect maximum work in system
340 MPI_Allreduce(workArray.data(), &max_work, 1, MPIDataTypeW, MPI_MAX,
341 communicator); MPI_Allreduce(workArray.data(), &sum_work, 1,
342 MPIDataTypeW, MPI_SUM, communicator);
343 MPI_Comm_size(communicator,&n_ranks);
344 d_max = (double)(max_work) * (double)n_ranks / (double)sum_work - 1.0;
345
346
347 // internal needs to be readded to argument list
348 const bool internal = false
349 for (int i = 0; i < estimation_limit && d_max > 0.1 && !internal; ++i)
350 {
351 balancer->setWork(tmp_work);
352 balancer->setup();
353 balancer->balance(true);
354 bool sane = true;
355 // check if the estimated boundaries are not too deformed
356 for (int j = 0; j < balancer->getDimension(); ++j)
357 sane = sane && (result_vertices.at(0)[j] <=
358 old_vertices.at(1)[j]);
359 MPI_Allreduce(MPI_IN_PLACE,&sane,1,MPI_CXX_BOOL,MPI_LAND,communicator);
360 if (sane)
361 {
362 old_tmp_vertices = tmp_vertices;
363 tmp_vertices = result_vertices;
364 V = 1.0;
365 for (int i = 0; i < balancer->getDimension(); ++i)
366 V *= ( tmp_vertices.at(1)[i] - tmp_vertices.at(0)[i] );
367 old_tmp_work = tmp_work;
368 tmp_work = rho * V;
369 }
370 else if (!sane || i == estimation_limit - 1)
371 {
372 balancer->getVertices() = old_tmp_vertices;
373 workArray->at(0) = old_tmp_work;
374 calculate_outline();
375 i = estimation_limit;
376 }
377 }
378
379 ((Staggered_LB<T,W>*)balancer.get())->setVertices(outline->data());
380 */
381 balancer->balance(loadbalancing_step);
382 calculate_outline();
383 break;
384 case LB_t::FORCEBASED:
385 balancer->balance(loadbalancing_step);
386 calculate_outline();
387 break;
388#ifdef ALL_VORONOI_ACTIVE
389 case LB_t::VORONOI:
390 balancer->balance(loadbalancing_step);
391 break;
392#endif
393 case LB_t::HISTOGRAM:
394 balancer->balance(loadbalancing_step);
395 calculate_outline();
396 break;
397 case LB_t::TENSOR_MAX:
398 balancer->balance(loadbalancing_step);
399 calculate_outline();
400 break;
401
402 default:
403 throw InvalidArgumentException(__FILE__, __func__, __LINE__,
404 "Unknown type of loadbalancing passed.");
405 }
406 loadbalancing_step++;
407 return balancer->getVertices();
408 }
409
410 /**********************/
411 /* getter functions */
412 /**********************/
413
417 std::vector<Point<T>> &getPrevVertices() {
418 return balancer->getPrevVertices();
419 };
420
425 std::vector<Point<T>> &getVertices() { return balancer->getVertices(); };
426
430 int getDimension() { return balancer->getDimension(); }
431
435 void getWork(std::vector<W> &result) { result = balancer->getWork(); }
436
439 W getWork() { return balancer->getWork().at(0); }
440
444 std::vector<int> &getNeighbors() { return balancer->getNeighbors(); }
445
450 std::vector<T> &getNeighborVertices() {
451 return balancer->getNeighborVertices();
452 }
453
457 return balancer->getEfficiency();
458 }
459
463 // currently only implemented for HISTOGRAM
464 return balancer->getEstimatedEfficiency();
465 }
466
467 /**********************/
468 /* setter functions */
469 /**********************/
470
476 void setSysSize(const std::vector<T> &sysSize) {
477 balancer.get()->setSysSize(sysSize);
478 }
479
488 void setMethodData(const void *data) {
489 balancer.get()->setAdditionalData(data);
490 }
491
496 void setProcGridParams(const std::vector<int> &loc,
497 const std::vector<int> &size) {
498 // todo(s.schulz): We should probably throw if the dimension does not match
499 procGridLoc.insert(procGridLoc.begin(), loc.begin(), loc.end());
500 procGridDim.insert(procGridDim.begin(), size.begin(), size.end());
501 procGridSet = true;
502 }
503
506 void setProcTag(int tag) { procTag = tag; }
507
512 void setMinDomainSize(const std::vector<T> &minSize) {
513 (balancer.get())->setMinDomainSize(minSize);
514 }
515
516#ifdef ALL_VTK_OUTPUT
521 void printVTKoutlines(const int step) {
522 // define grid points, i.e. vertices of local domain
523 auto points = vtkSmartPointer<vtkPoints>::New();
524 // allocate space for eight points (describing the cuboid around the domain)
525 points->Allocate(8 * balancer->getDimension());
526
527 int n_ranks;
528 int local_rank;
529
530 if (!vtk_init) {
531 controller = vtkMPIController::New();
532 controller->Initialize();
533 controller->SetGlobalController(controller);
534 vtk_init = true;
535 }
536
537 MPI_Comm_rank(communicator, &local_rank);
538 MPI_Comm_size(communicator, &n_ranks);
539
540 std::vector<Point<W>> tmp_outline(2);
541 std::vector<W> tmp_0(
542 {outline.at(0)[0], outline.at(0)[1], outline.at(0)[2]});
543 std::vector<W> tmp_1(
544 {outline.at(1)[0], outline.at(1)[1], outline.at(1)[2]});
545 tmp_outline.at(0) = Point<W>(tmp_0);
546 tmp_outline.at(1) = Point<W>(tmp_1);
547
548 for (auto z = 0; z <= 1; ++z)
549 for (auto y = 0; y <= 1; ++y)
550 for (auto x = 0; x <= 1; ++x) {
551 points->InsertPoint(x + 2 * y + 4 * z, tmp_outline.at(x)[0],
552 tmp_outline.at(y)[1], tmp_outline.at(z)[2]);
553 }
554
555 auto hexa = vtkSmartPointer<vtkVoxel>::New();
556 for (int i = 0; i < 8; ++i)
557 hexa->GetPointIds()->SetId(i, i);
558
559 auto cellArray = vtkSmartPointer<vtkCellArray>::New();
560 cellArray->InsertNextCell(hexa);
561
562 // define work array (length = 1)
563 auto work = vtkSmartPointer<vtkFloatArray>::New();
564 work->SetNumberOfComponents(1);
565 work->SetNumberOfTuples(1);
566 work->SetName("work");
567 W total_work = std::accumulate(balancer->getWork().begin(),
568 balancer->getWork().end(), (W)0);
569 work->SetValue(0, total_work);
570
571 // define rank array (length = 4, x,y,z, rank)
572 int rank = 0;
573 MPI_Comm_rank(communicator, &rank);
574 auto coords = vtkSmartPointer<vtkFloatArray>::New();
575 coords->SetNumberOfComponents(3);
576 coords->SetNumberOfTuples(1);
577 coords->SetName("coords");
578 coords->SetValue(0, procGridLoc.at(0));
579 coords->SetValue(1, procGridLoc.at(1));
580 coords->SetValue(2, procGridLoc.at(2));
581
582 auto rnk = vtkSmartPointer<vtkFloatArray>::New();
583 rnk->SetNumberOfComponents(1);
584 rnk->SetNumberOfTuples(1);
585 rnk->SetName("MPI rank");
586 rnk->SetValue(0, rank);
587
588 // define tag array (length = 1)
589 auto tag = vtkSmartPointer<vtkIntArray>::New();
590 tag->SetNumberOfComponents(1);
591 tag->SetNumberOfTuples(1);
592 tag->SetName("tag");
593 tag->SetValue(0, procTag);
594
595 // determine extent of system
596 W known_unused global_extent[6];
597 W local_min[3];
598 W local_max[3];
599 W global_min[3];
600 W global_max[3];
601 W width[3];
602 W volume = (W)1.0;
603
604 for (int i = 0; i < 3; ++i) {
605 local_min[i] = outline.at(0)[i];
606 local_max[i] = outline.at(1)[i];
607 width[i] = local_max[i] - local_min[i];
608 volume *= width[i];
609 }
610
611 W surface = (W)2.0 * width[0] * width[1] + (W)2.0 * width[1] * width[2] +
612 (W)2.0 * width[0] * width[2];
613
614 auto expanse = vtkSmartPointer<vtkFloatArray>::New();
615 expanse->SetNumberOfComponents(6);
616 expanse->SetNumberOfTuples(1);
617 expanse->SetName("expanse");
618 expanse->SetValue(0, width[0]);
619 expanse->SetValue(1, width[0]);
620 expanse->SetValue(2, width[0]);
621 expanse->SetValue(3, volume);
622 expanse->SetValue(4, surface);
623 expanse->SetValue(5, surface / volume);
624
625 MPI_Allreduce(&local_min, &global_min, 3, MPIDataTypeW, MPI_MIN,
626 communicator);
627 MPI_Allreduce(&local_max, &global_max, 3, MPIDataTypeW, MPI_MAX,
628 communicator);
629
630 for (int i = 0; i < 3; ++i) {
631 global_extent[2 * i] = global_min[i];
632 global_extent[2 * i + 1] = global_max[i];
633 }
634
635 // create a structured grid and assign data to it
636 auto unstructuredGrid = vtkSmartPointer<vtkUnstructuredGrid>::New();
637 unstructuredGrid->SetPoints(points);
638 unstructuredGrid->SetCells(VTK_VOXEL, cellArray);
639 unstructuredGrid->GetCellData()->AddArray(work);
640 unstructuredGrid->GetCellData()->AddArray(expanse);
641 unstructuredGrid->GetCellData()->AddArray(rnk);
642 unstructuredGrid->GetCellData()->AddArray(coords);
643 unstructuredGrid->GetCellData()->AddArray(tag);
644
645 //DEBUG
646 /*
647 createDirectory("vtk_outline_debug");
648 std::ostringstream ss_local_debug;
649 ss_local_debug << "vtk_outline_debug/ALL_vtk_outline_" << std::setw(7)
650 << std::setfill('0') << step << "_" << local_rank << ".vtk";
651
652 auto debug_writer = vtkSmartPointer<vtkUnstructuredGridWriter>::New();
653 debug_writer->SetFileName(ss_local_debug.str().c_str());
654 debug_writer->SetInputData(unstructuredGrid);
655 debug_writer->Write();
656 */
657 //ENDDEBUG
658
659 createDirectory("vtk_outline");
660 std::ostringstream ss_local;
661 ss_local << "vtk_outline/ALL_vtk_outline_" << std::setw(7)
662 << std::setfill('0') << step << "_" << local_rank << ".vtu";
663
664 auto writer = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
665 writer->SetInputData(unstructuredGrid);
666 writer->SetFileName(ss_local.str().c_str());
667 //writer->SetDataModeToAscii();
668 writer->SetDataModeToBinary();
669 writer->Write();
670
671 // if (local_rank == 0)
672 //{
673 std::ostringstream ss_para;
674 ss_para << "vtk_outline/ALL_vtk_outline_" << std::setw(7)
675 << std::setfill('0') << step << ".pvtu";
676 // create the parallel writer
677 auto parallel_writer =
678 vtkSmartPointer<vtkXMLPUnstructuredGridWriter>::New();
679 parallel_writer->SetFileName(ss_para.str().c_str());
680 parallel_writer->SetNumberOfPieces(n_ranks);
681 parallel_writer->SetStartPiece(local_rank);
682 parallel_writer->SetEndPiece(local_rank);
683 parallel_writer->SetInputData(unstructuredGrid);
684 //parallel_writer->SetDataModeToAscii();
685 parallel_writer->SetDataModeToBinary();
686 parallel_writer->Write();
687 //}
688 }
689
694 void printVTKvertices(const int step) {
695 int n_ranks;
696 int local_rank;
697
698 MPI_Comm_rank(communicator, &local_rank);
699 MPI_Comm_size(communicator, &n_ranks);
700
701 // local vertices
702 // (vertices + work)
703 T local_vertices[nVertices * balancer->getDimension() + 1];
704
705 for (int v = 0; v < nVertices; ++v) {
706 for (int d = 0; d < balancer->getDimension(); ++d) {
707 local_vertices[v * balancer->getDimension() + d] =
708 balancer->getVertices().at(v)[d];
709 }
710 }
711 local_vertices[nVertices * balancer->getDimension()] =
712 (T)balancer->getWork().at(0);
713
714 /*
715 T *global_vertices;
716 if (local_rank == 0) {
717 global_vertices =
718 new T[n_ranks * (nVertices * balancer->getDimension() + 1)];
719 }
720 */
721 T global_vertices[n_ranks * (nVertices * balancer->getDimension() + 1)];
722
723 // collect all works and vertices on a single process
724 MPI_Gather(local_vertices, nVertices * balancer->getDimension() + 1,
725 MPIDataTypeT, global_vertices,
726 nVertices * balancer->getDimension() + 1, MPIDataTypeT, 0,
727 communicator);
728
729 if (local_rank == 0) {
730 auto vtkpoints = vtkSmartPointer<vtkPoints>::New();
731#ifdef VTK_CELL_ARRAY_V2
732 vtkNew<vtkUnstructuredGrid> unstructuredGrid;
733 unstructuredGrid->Allocate(n_ranks + 10);
734#else
735 auto unstructuredGrid = vtkSmartPointer<vtkUnstructuredGrid>::New();
736#endif
737
738 // enter vertices into unstructured grid
739 for (int i = 0; i < n_ranks; ++i) {
740 for (int v = 0; v < nVertices; ++v) {
741 vtkpoints->InsertNextPoint(
742 global_vertices[i * (nVertices * balancer->getDimension() + 1) +
743 v * balancer->getDimension() + 0],
744 global_vertices[i * (nVertices * balancer->getDimension() + 1) +
745 v * balancer->getDimension() + 1],
746 global_vertices[i * (nVertices * balancer->getDimension() + 1) +
747 v * balancer->getDimension() + 2]);
748 }
749 }
750 unstructuredGrid->SetPoints(vtkpoints);
751
752 // data arrays for work and cell id
753 auto work = vtkSmartPointer<vtkFloatArray>::New();
754 auto cell = vtkSmartPointer<vtkFloatArray>::New();
755 work->SetNumberOfComponents(1);
756 work->SetNumberOfTuples(n_ranks);
757 work->SetName("work");
758 cell->SetNumberOfComponents(1);
759 cell->SetNumberOfTuples(n_ranks);
760 cell->SetName("cell id");
761
762 for (int n = 0; n < n_ranks; ++n) {
763
764#ifdef VTK_CELL_ARRAY_V2
765 // define grid points, i.e. vertices of local domain
766 vtkIdType pointIds[8] = {8 * n + 0, 8 * n + 1, 8 * n + 2, 8 * n + 3,
767 8 * n + 4, 8 * n + 5, 8 * n + 6, 8 * n + 7};
768
769 vtkIdType faces[48] = {
770 3, 8 * n + 0, 8 * n + 2, 8 * n + 1,
771 3, 8 * n + 1, 8 * n + 2, 8 * n + 3,
772 3, 8 * n + 0, 8 * n + 4, 8 * n + 2,
773 3, 8 * n + 2, 8 * n + 4, 8 * n + 6,
774 3, 8 * n + 2, 8 * n + 6, 8 * n + 3,
775 3, 8 * n + 3, 8 * n + 6, 8 * n + 7,
776 3, 8 * n + 1, 8 * n + 5, 8 * n + 3,
777 3, 8 * n + 3, 8 * n + 5, 8 * n + 7,
778 3, 8 * n + 0, 8 * n + 4, 8 * n + 1,
779 3, 8 * n + 1, 8 * n + 4, 8 * n + 5,
780 3, 8 * n + 4, 8 * n + 6, 8 * n + 5,
781 3, 8 * n + 5, 8 * n + 6, 8 * n + 7};
782
783 unstructuredGrid->InsertNextCell(VTK_POLYHEDRON, 8, pointIds, 12,
784 faces);
785#else
786 // define grid points, i.e. vertices of local domain
787 vtkIdType pointIds[8] = {8 * n + 0, 8 * n + 1, 8 * n + 2, 8 * n + 3,
788 8 * n + 4, 8 * n + 5, 8 * n + 6, 8 * n + 7};
789
790 auto faces = vtkSmartPointer<vtkCellArray>::New();
791 // setup faces of polyhedron
792 vtkIdType f0[3] = {8 * n + 0, 8 * n + 2, 8 * n + 1};
793 vtkIdType f1[3] = {8 * n + 1, 8 * n + 2, 8 * n + 3};
794
795 vtkIdType f2[3] = {8 * n + 0, 8 * n + 4, 8 * n + 2};
796 vtkIdType f3[3] = {8 * n + 2, 8 * n + 4, 8 * n + 6};
797
798 vtkIdType f4[3] = {8 * n + 2, 8 * n + 6, 8 * n + 3};
799 vtkIdType f5[3] = {8 * n + 3, 8 * n + 6, 8 * n + 7};
800
801 vtkIdType f6[3] = {8 * n + 1, 8 * n + 5, 8 * n + 3};
802 vtkIdType f7[3] = {8 * n + 3, 8 * n + 5, 8 * n + 7};
803
804 vtkIdType f8[3] = {8 * n + 0, 8 * n + 4, 8 * n + 1};
805 vtkIdType f9[3] = {8 * n + 1, 8 * n + 4, 8 * n + 5};
806
807 vtkIdType fa[3] = {8 * n + 4, 8 * n + 6, 8 * n + 5};
808 vtkIdType fb[3] = {8 * n + 5, 8 * n + 6, 8 * n + 7};
809
810 faces->InsertNextCell(3, f0);
811 faces->InsertNextCell(3, f1);
812 faces->InsertNextCell(3, f2);
813 faces->InsertNextCell(3, f3);
814 faces->InsertNextCell(3, f4);
815 faces->InsertNextCell(3, f5);
816 faces->InsertNextCell(3, f6);
817 faces->InsertNextCell(3, f7);
818 faces->InsertNextCell(3, f8);
819 faces->InsertNextCell(3, f9);
820 faces->InsertNextCell(3, fa);
821 faces->InsertNextCell(3, fb);
822
823 unstructuredGrid->InsertNextCell(VTK_POLYHEDRON, 8, pointIds, 12,
824 faces->GetPointer());
825#endif
826 work->SetValue(
827 n, global_vertices[n * (nVertices * balancer->getDimension() + 1) +
828 8 * balancer->getDimension()]);
829 cell->SetValue(n, (T)n);
830 }
831 unstructuredGrid->GetCellData()->AddArray(work);
832 unstructuredGrid->GetCellData()->AddArray(cell);
833
834 createDirectory("vtk_vertices");
835 std::ostringstream filename;
836 filename << "vtk_vertices/ALL_vtk_vertices_" << std::setw(7)
837 << std::setfill('0') << step << ".vtu";
838 auto writer = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
839 writer->SetInputData(unstructuredGrid);
840 writer->SetFileName(filename.str().c_str());
841 writer->SetDataModeToAscii();
842 // writer->SetDataModeToBinary();
843 writer->Write();
844
845 // delete[] global_vertices;
846 }
847 }
848#endif
849
850private:
852 LB_t method;
853
857 std::vector<Point<T>> outline;
858
860 std::unique_ptr<LB<T, W>> balancer;
861
864 std::vector<W> workArray;
865
867 MPI_Comm communicator;
868
870 int nVertices;
871
873 int nNeighbors;
874
877 void calculate_outline() {
878 // calculation only possible if there are vertices defining the domain
879 if (balancer->getVertices().size() > 0) {
880 outline.resize(2);
881 // setup the outline with the first point
882 for (int i = 0; i < balancer->getDimension(); ++i) {
883 outline.at(0) = balancer->getVertices().at(0);
884 outline.at(1) = balancer->getVertices().at(0);
885 }
886 // compare value of each outline point with all vertices to find the
887 // maximum extend of the domain in each dimension
888 for (int i = 1; i < (int)balancer->getVertices().size(); ++i) {
889 Point<T> p = balancer->getVertices().at(i);
890 for (int j = 0; j < balancer->getDimension(); ++j) {
891 outline.at(0)[j] = std::min(outline.at(0)[j], p[j]);
892 outline.at(1)[j] = std::max(outline.at(1)[j], p[j]);
893 }
894 }
895 }
896 }
897
901 void createDirectory(const char *filename) {
902 errno = 0;
903 mode_t mode = S_IRWXU | S_IRWXG; // rwx for user and group
904 if (mkdir(filename, mode)) {
905 if (errno != EEXIST) {
906 throw FilesystemErrorException(__FILE__, __func__, __LINE__,
907 "Could not create output directory.");
908 }
909 }
910 // check permission in case directory already existed, but had wrong
911 // permissions
912 struct stat attr;
913 stat(filename, &attr);
914 if ((attr.st_mode & S_IRWXU) != S_IRWXU) {
916 __FILE__, __func__, __LINE__,
917 "Possibly already existing directory does not have correct "
918 "permissions (rwx) for user");
919 }
920 }
921
925 int estimation_limit;
926
928 std::vector<int> procGridLoc;
930 std::vector<int> procGridDim;
932 bool procGridSet;
934 int procTag;
935
937 MPI_Datatype MPIDataTypeT;
938
940 MPI_Datatype MPIDataTypeW;
941
943 int loadbalancing_step;
944
945#ifdef ALL_VTK_OUTPUT
947 bool vtk_init;
948
950 vtkMPIController *controller;
951
953 vtkMPIController *controller_vertices;
954#endif
955};
956
957} // namespace ALL
958
959#endif
960
961// VIM modeline for indentation
962// vim: sw=2 et
#define ALL_ESTIMATION_LIMIT_DEFAULT
Definition ALL.hpp:79
#define known_unused
Definition ALL_Defines.h:49
std::vector< T > & getNeighborVertices()
Definition ALL.hpp:450
T getGamma()
Definition ALL.hpp:275
void setSysSize(const std::vector< T > &sysSize)
Definition ALL.hpp:476
ALL(const LB_t m, const int d, const std::vector< Point< T > > &inp, const T g)
Definition ALL.hpp:185
void setProcGridParams(const std::vector< int > &loc, const std::vector< int > &size)
Definition ALL.hpp:496
std::vector< Point< T > > & getPrevVertices()
Definition ALL.hpp:417
void setMethodData(const void *data)
Definition ALL.hpp:488
void setCommunicator(const MPI_Comm comm)
Definition ALL.hpp:213
std::vector< int > & getNeighbors()
Definition ALL.hpp:444
void printVTKvertices(const int step)
Definition ALL.hpp:694
void getWork(std::vector< W > &result)
Definition ALL.hpp:435
ALL()
Definition ALL.hpp:108
std::vector< Point< T > > & getVertices()
Definition ALL.hpp:425
W getEstimatedEfficiency()
Definition ALL.hpp:462
W getEfficiency()
Definition ALL.hpp:456
W getWork()
Definition ALL.hpp:439
void setMinDomainSize(const std::vector< T > &minSize)
Definition ALL.hpp:512
void setWork(const W work)
Definition ALL.hpp:283
void setup()
Definition ALL.hpp:294
std::vector< Point< T > > & balance()
Definition ALL.hpp:301
void setVertices(const std::vector< Point< T > > &inp)
Definition ALL.hpp:197
void printVTKoutlines(const int step)
Definition ALL.hpp:521
int getDimension()
Definition ALL.hpp:430
ALL(const LB_t m, const int d, const T g)
Definition ALL.hpp:144
void setWork(const std::vector< W > &work)
Definition ALL.hpp:290
void setGamma(const T g)
Definition ALL.hpp:279
void setProcTag(int tag)
Definition ALL.hpp:506
~ALL()=default
destructor
Definition ALL.hpp:77
LB_t
enum type to describe the different balancing methods
Definition ALL.hpp:82
@ FORCEBASED
unstructured-mesh load balancing
Definition ALL.hpp:88
@ STAGGERED
staggered grid load balancing
Definition ALL.hpp:84
@ VORONOI
voronoi cell based load balancing
Definition ALL.hpp:91
@ UNIMPLEMENTED
Unimplemented load balancing NEVER SUPPORTED.
Definition ALL.hpp:98
@ HISTOGRAM
histogram based load balancing
Definition ALL.hpp:94
@ TENSOR_MAX
tensor based load balancing using maximum values
Definition ALL.hpp:96
@ TENSOR
tensor based load balancing
Definition ALL.hpp:86
Exception to be used for missmatches in dimension for ALL::Point class.