Peano
KernelBenchmarks-main.cpp
Go to the documentation of this file.
1 // This file is part of the ExaHyPE2 project. For conditions of distribution and
2 // use, please see the copyright notice at www.peano-framework.org
3 #include "KernelBenchmarks-main.h"
4 
5 #include "AbstractFVSolver.h"
6 #include "Constants.h"
7 #include "DataRepository.h"
8 #include "FVSolver.h"
9 #include "celldata/FVSolverQ.h"
10 #include "celldata/FVSolverQReconstructed.h"
11 #include "exahype2/CellData.h"
12 #include "exahype2/EnclaveBookkeeping.h"
13 #include "exahype2/fv/BoundaryConditions.h"
14 #include "exahype2/fv/PatchUtils.h"
15 #include "facedata/FVSolverQNew.h"
16 #include "kernels/FVSolver/ReconstructPatch.h"
17 #include "peano4/datamanagement/CellMarker.h"
18 #include "peano4/datamanagement/FaceEnumerator.h"
19 #include "peano4/peano4.h"
20 #include "peano4/utils/Globals.h"
21 #include "peano4/utils/Loop.h"
22 #include "tarch/logging/Log.h"
23 #include "tarch/plotter/griddata/blockstructured/PeanoTextPatchFileWriter.h"
24 #include "tarch/timing/Measurement.h"
25 #include "tarch/timing/Watch.h"
26 #include "tasks/FVSolverEnclaveTask.h"
27 #include "toolbox/blockstructured/Projection.h"
28 
30 
31 namespace {
32  tarch::logging::Log _log("");
33  tarch::timing::Measurement timeStepMeasurement;
34  double printNextTimeStep = 0.0;
35  bool haveReceivedNonCriticalAssertion = false;
36 
37  std::vector<peano4::datamanagement::CellMarker> _cells;
38  const tarch::la::Vector<DIMENSIONS, double> CellSize = DomainSize / static_cast<double>(GridLength);
39 } // namespace
40 
41 void reportRuntime(const tarch::timing::Measurement& measurement) {
42  const std::string context = "reportRuntime()";
43  constexpr int LabelWidth = 25; // Set a fixed width for all labels
44 
45  if (measurement.getNumberOfMeasurements() == 0) {
46  logInfo(context, "Timing Measurement Report: No data recorded.");
47  return;
48  }
49 
50  logInfo(context, "---------- Timing Measurement Report ----------");
51  double avg = measurement.getValue();
52  double std_dev = measurement.getStandardDeviation();
53  double min_val = measurement.min();
54  double max_val = measurement.max();
55  int min_idx = measurement.getMinMeasurement();
56  int max_idx = measurement.getMaxMeasurement();
57  double min_dev_percent = (avg == 0) ? 0.0 : (min_val / avg - 1.0) * 100.0;
58  double max_dev_percent = (avg == 0) ? 0.0 : (max_val / avg - 1.0) * 100.0;
59 
60  std::stringstream msg;
61 
62  msg.str("");
63  msg << std::left << std::setw(LabelWidth) << "Average Value:" << std::fixed << std::setprecision(8) << avg << "s";
64  logInfo(context, msg.str());
65 
66  msg.str("");
67  msg << std::left << std::setw(LabelWidth) << "Average Cell Updates:" << std::fixed << std::setprecision(0)
68  << _cells.size() / avg;
69  logInfo(context, msg.str());
70 
71  msg.str("");
72  msg << std::left << std::setw(LabelWidth) << "Standard Deviation:" << std::scientific << std::setprecision(8)
73  << std_dev;
74  logInfo(context, msg.str());
75 
76  msg.str("");
77  msg << std::left << std::setw(LabelWidth) << "# Measurements:" << measurement.getNumberOfMeasurements();
78  logInfo(context, msg.str());
79 
80  msg.str("");
81  msg
82  << std::left << std::setw(LabelWidth) << "Min Value:" << std::fixed << std::setprecision(8) << min_val << "s (at #"
83  << min_idx << ") "
84  << "[" << std::fixed << std::setprecision(2) << min_dev_percent << "%]";
85  logInfo(context, msg.str());
86 
87  msg.str("");
88  msg
89  << std::left << std::setw(LabelWidth) << "Max Value:" << std::fixed << std::setprecision(8) << max_val << "s (at #"
90  << max_idx << ") " << std::showpos // Show the '+' sign
91  << "[" << std::fixed << std::setprecision(2) << max_dev_percent << "%]";
92  logInfo(context, msg.str());
93 
94  msg.str("");
95  msg << std::left << std::setw(LabelWidth) << "Most Recent Value:" << std::fixed << std::setprecision(8)
96  << measurement.getMostRecentValue() << "s";
97  logInfo(context, msg.str());
98 
99  msg.str("");
100  msg << std::left << std::setw(LabelWidth) << "Accumulated Time:" << std::fixed << std::setprecision(8)
101  << measurement.getAccumulatedValue() << "s";
102  logInfo(context, msg.str());
103 
104  logInfo(context, "-----------------------------------------------");
105 }
106 
107 void plotGrid() {
108  repositories::startPlottingStep();
109  static int counter = -1;
110  counter++;
111 
112  std::ostringstream snapshotFileName;
113  snapshotFileName << "solutions/solution-FVSolver-" << counter;
114 
115  auto writer = tarch::plotter::griddata::blockstructured::PeanoTextPatchFileWriter(
116  DIMENSIONS,
117  snapshotFileName.str(),
118  "solutions/solution-FVSolver",
119  tarch::plotter::griddata::blockstructured::PeanoTextPatchFileWriter::IndexFileMode::AppendNewData,
120  0.5 * (repositories::getMinTimeStamp() + repositories::getMaxTimeStamp())
121  );
122 
123  auto dataWriter = writer.createCellDataWriter(
124  "FVSolverQ",
125  AbstractFVSolver::NumberOfFiniteVolumesPerAxisPerPatch,
126  AbstractFVSolver::NumberOfUnknowns,
127  PlotDescription,
128  "",
129  nullptr
130  );
131  dataWriter->setPrecision(PlotterPrecision);
132 
133  for (auto& marker : _cells) {
134  const int patchIndex = writer.plotPatch(marker.x() - marker.h() * 0.5, marker.h());
135  int cellIndex = dataWriter->getFirstCellWithinPatch(patchIndex);
136  int currentDoF = 0;
137  double* QOut = DataRepository::getInstance().getCellQ(marker.x(), marker.h());
138  dfor(k, AbstractFVSolver::NumberOfFiniteVolumesPerAxisPerPatch) {
139  double* data = QOut + currentDoF;
140 
141  dataWriter->plotCell(cellIndex, data);
142  cellIndex++;
143  currentDoF += AbstractFVSolver::NumberOfUnknowns;
144  }
145  }
146  dataWriter->close();
147  writer.writeToFile();
148  delete dataWriter;
149  repositories::finishPlottingStep();
150 }
151 
153  const double* NOALIAS Qinside,
154  double* NOALIAS Qoutside,
155  const tarch::la::Vector<DIMENSIONS, double>& x,
156  const tarch::la::Vector<DIMENSIONS, double>& h,
157  double t,
158  double dt,
159  int normal
160 ) {
161  repositories::instanceOfFVSolver.boundaryConditions(Qinside, Qoutside, x, h, t, normal);
162 }
163 
164 void applyBoundaryConditionsToAxis(int axis, double t, double dt) {
165  const double ratio = DomainSize[0] / CellSize[0];
166  const int depth = static_cast<int>(std::round(std::log(ratio) / std::log(3.0)));
167  dfore(volume, GridLength, axis, 0) {
168  tarch::la::Vector<DIMENSIONS, double> faceOffset = CellSize / 2.0;
169  // --- Low Boundary ---
170  tarch::la::Vector<DIMENSIONS, double>
171  faceCenterLow = tarch::la::multiplyComponents(tarch::la::convertScalar<double>(volume), CellSize) + faceOffset
172  + DomainOffset;
173  int faceIndexLow = DataRepository::Indexing::getFaceIndex(faceCenterLow, CellSize, axis);
174  faceCenterLow[axis] = DomainOffset[axis];
175 
176  // --- High Boundary ---
177  tarch::la::Vector<DIMENSIONS, double> faceCenterHigh = faceCenterLow;
178  faceCenterHigh[axis] = DomainSize[axis] - faceOffset[axis] + DomainOffset[axis]; // high cell Centre
179  int faceIndexHigh = DataRepository::Indexing::getFaceIndex(faceCenterHigh, CellSize, axis);
180  faceCenterHigh[axis] = DomainSize[axis] + DomainOffset[axis]; // High face centre
181 
182  // Apply BC to low face
183  double* lowFacePtr = DataRepository::getInstance().getFaceQNew(faceIndexLow, axis);
184  ::exahype2::fv::applyBoundaryConditions<double>(
186  faceCenterLow,
187  CellSize,
188  t,
189  dt,
190  AbstractFVSolver::NumberOfFiniteVolumesPerAxisPerPatch,
191  1,
192  AbstractFVSolver::NumberOfUnknowns + AbstractFVSolver::NumberOfAuxiliaryVariables,
193  axis,
194  lowFacePtr
195  );
196 
197  // Apply BC to high face
198  double* highFacePtr = DataRepository::getInstance().getFaceQNew(faceIndexHigh, axis + DIMENSIONS);
199  ::exahype2::fv::applyBoundaryConditions<double>(
201  faceCenterHigh,
202  CellSize,
203  t,
204  dt,
205  AbstractFVSolver::NumberOfFiniteVolumesPerAxisPerPatch,
206  1,
207  AbstractFVSolver::NumberOfUnknowns + AbstractFVSolver::NumberOfAuxiliaryVariables,
208  axis + DIMENSIONS,
209  highFacePtr
210  );
211  }
212 }
213 
215  const double t = 0.0;
216  const double dt = 0.0;
217 
218  for (int axis = 0; axis < DIMENSIONS; ++axis) {
220  }
221 }
222 
224  for (auto& marker : _cells) {
225  double* QOut = DataRepository::getInstance().getCellQ(marker.x(), marker.h());
226  double* faces[TWO_TIMES_D];
227  for (int i = 0; i < TWO_TIMES_D; i++) {
228  faces[i] = DataRepository::getInstance().getFaceQNew(marker.x(), marker.h(), i);
229  }
230  toolbox::blockstructured::projectPatchSolutionOntoFaces(
231  AbstractFVSolver::NumberOfFiniteVolumesPerAxisPerPatch,
232  1,
233  AbstractFVSolver::NumberOfUnknowns,
234  AbstractFVSolver::NumberOfAuxiliaryVariables,
235  QOut,
236  faces
237  );
238  }
239 }
240 
242  for (auto& marker : _cells) {
243  double* QOut = DataRepository::getInstance().getCellQ(marker.x(), marker.h());
244  // Initialise cell data with initial conditions.
245  dfor(volume, AbstractFVSolver::NumberOfFiniteVolumesPerAxisPerPatch) {
246  int linearIndex = peano4::utils::dLinearised(volume, AbstractFVSolver::NumberOfFiniteVolumesPerAxisPerPatch)
247  * (AbstractFVSolver::NumberOfUnknowns + AbstractFVSolver::NumberOfAuxiliaryVariables);
248  FVSolver::initialCondition(
249  QOut + linearIndex,
250  ::exahype2::fv::
251  getVolumeCentre(marker.x(), marker.h(), AbstractFVSolver::NumberOfFiniteVolumesPerAxisPerPatch, volume),
252  ::exahype2::fv::getVolumeSize(marker.h(), AbstractFVSolver::NumberOfFiniteVolumesPerAxisPerPatch),
253  true
254  );
255  }
256  }
257 }
258 
259 peano4::datamanagement::CellMarker createMarkerFromIndex(const tarch ::la ::Vector<DIMENSIONS, int>& index) {
260  tarch::la::Vector<DIMENSIONS, double> cellIndex = tarch::la::convertScalar<double>(index);
261  tarch::la::Vector<DIMENSIONS, double> cellPositionBase = DomainOffset + CellSize / 2.0;
262  tarch::la::Vector<DIMENSIONS, double>
263  cellPosition = cellPositionBase + tarch::la::multiplyComponents(cellIndex, CellSize);
264  peano4::grid::GridTraversalEvent event;
265  event.setX(cellPosition);
266  event.setH(CellSize);
267  return peano4::datamanagement::CellMarker(event);
268 }
269 
270 void createGrid() {
271  repositories::startGridConstructionStep();
272  _cells.reserve(std::pow(GridLength, DIMENSIONS));
273  dfor(volume, GridLength) { _cells.push_back(createMarkerFromIndex(volume)); }
274  repositories::finishGridConstructionStep();
275 }
276 
277 void initGrid() {
278  repositories::startGridInitialisationStep();
281  repositories::finishGridInitialisationStep();
282 }
283 
284 int main(int argc, char** argv) {
285  int returnCode = EXIT_SUCCESS;
286 
287  peano4::init(
288  &argc,
289  &argv,
290  benchmarks::exahype2::kernelbenchmarks::DomainOffset,
291  benchmarks::exahype2::kernelbenchmarks::DomainSize
292  );
293 
294  repositories::initLogFilters();
295  repositories::startSimulation();
296 
297  createGrid();
298  logInfo("main()", "Created grid with " << _cells.size() << " cells");
299  initGrid();
300  logInfo("main()", "Initialised grid");
301 
303 
304  ::exahype2::CellData<double, double, 3> patchData(_cells.size(), tarch::MemoryLocation::Pinned);
305 
306  double* QInFlattened = tarch::allocateMemory<double>(
307  _cells.size() * grid.getCellQReconstructedCardinality(),
308  tarch::MemoryLocation::Pinned
309  );
310 
311  double* QOutFlattened = tarch::allocateMemory<double>(
312  _cells.size() * grid.getCellQCardinality(),
313  tarch::MemoryLocation::Pinned
314  );
315  double* QOutFlattenedHost = grid.getCellQ(0);
316 
317  tarch::la::Vector<DIMENSIONS, double>* cellCentresHost = new tarch::la::Vector<DIMENSIONS, double>[_cells.size()];
318  tarch::la::Vector<DIMENSIONS, double>* cellSizeHost = new tarch::la::Vector<DIMENSIONS, double>[_cells.size()];
319 
320  for (int patchIndex = 0; patchIndex < _cells.size(); patchIndex++) {
321  cellCentresHost[patchIndex] = _cells[patchIndex].x();
322  cellSizeHost[patchIndex] = _cells[patchIndex].h();
323  for (int type = 0; type < 3; type++) {
324  for (int face = 0; face < TWO_TIMES_D; face++) {
325  patchData.faces[patchIndex][type][face] = grid.getFaceQNewDevice(
326  _cells[patchIndex].x(),
327  _cells[patchIndex].h(),
328  face
329  );
330  }
331  }
332  }
333  tarch::copyTo<tarch::la::Vector<DIMENSIONS, double>, tarch::la::Vector<DIMENSIONS, double>>(
334  patchData.cellCentre,
335  cellCentresHost,
336  _cells.size(),
337  tarch::MemoryLocation::Pinned
338  );
339  tarch::copyTo<tarch::la::Vector<DIMENSIONS, double>, tarch::la::Vector<DIMENSIONS, double>>(
340  patchData.cellSize,
341  cellSizeHost,
342  _cells.size(),
343  tarch::MemoryLocation::Pinned
344  );
345 
346  tarch::copyTo<double, double>(
347  QOutFlattened,
348  QOutFlattenedHost,
349  _cells.size() * grid.getCellQCardinality(),
350  tarch::MemoryLocation::Pinned
351  );
352 
353  patchData.QOutFlattened = QOutFlattened;
354  patchData.QInFlattened = QInFlattened;
355 
356  tarch::timing::Measurement measurement(1.0);
357  tarch::timing::Watch timeStepWatch("::", "main()", false);
358  double timeStamp = 0.0;
359  double nextMinPlotTimeStamp = FirstPlotTimeStamp;
360  double nextMaxPlotTimeStamp = FirstPlotTimeStamp;
361  bool continueToSolve = true;
362 
363  while (continueToSolve) {
364  if (repositories::getMinTimeStamp() >= nextMinPlotTimeStamp
365  or repositories::getMaxTimeStamp() >= nextMaxPlotTimeStamp) {
366  if (repositories::getMinTimeStamp() >= nextMinPlotTimeStamp) {
367  nextMinPlotTimeStamp += TimeInBetweenPlots;
368  }
369  if (repositories::getMaxTimeStamp() >= nextMaxPlotTimeStamp) {
370  nextMaxPlotTimeStamp += TimeInBetweenPlots;
371  }
372 
373  if constexpr (ShouldPlot) {
374  plotGrid();
375  }
376  }
377 
378  timeStepWatch.start();
379 
381  grid.copyFacesToDevice();
382 
383  patchData.t = timeStamp;
384  patchData.dt = repositories::instanceOfFVSolver.getAdmissibleTimeStepSize();
385  repositories::startTimeStep();
386 
387  tarch::timing::Watch watch("", "", false);
388  benchmarks::exahype2::kernelbenchmarks::tasks::FVSolverEnclaveTask::fuse(patchData);
389  watch.stop();
390  measurement.setValue(watch.getCalendarTime());
391 
392  timeStamp += repositories::instanceOfFVSolver.getAdmissibleTimeStepSize();
393  repositories::finishTimeStep();
394 
395  tarch::copyFrom<double, double>(
396  QOutFlattened,
397  QOutFlattenedHost,
398  _cells.size() * grid.getCellQCardinality(),
399  tarch::MemoryLocation::Pinned
400  );
401 
402  grid.copyFacesFromDevice();
403  timeStepWatch.stop();
404  timeStepMeasurement.setValue(timeStepWatch.getCalendarTime());
405 
406  if (timeStepMeasurement.getAccumulatedValue() >= printNextTimeStep) {
407  repositories::printTimeStep(false);
408  if (timeStepMeasurement.getNumberOfMeasurements() > 0) {
409  logInfo("main()", "Average cell updates per second: " << _cells.size() / measurement.getValue());
410  }
411  printNextTimeStep += repositories::PrintTimeStepIntervalInSeconds;
412  }
413 
414  if (repositories::getMinTimeStamp() >= MinTerminalTime and repositories::getMaxTimeStamp() >= MaxTerminalTime) {
415  continueToSolve = false;
416  }
417 
418  if (tarch::hasNonCriticalAssertionBeenViolated() and not haveReceivedNonCriticalAssertion) {
419  continueToSolve = false;
420  haveReceivedNonCriticalAssertion = true;
421  logError("main()", "Noncritical assertion has been triggered in code. Dump final state and terminate.");
422  }
423  }
424 
425  // Plot last timestep
426  if constexpr (ShouldPlot) {
427  plotGrid();
428  }
429  repositories::finishSimulation();
430 
431  if (haveReceivedNonCriticalAssertion) {
432  logWarning("main()", "Terminated gracefully after noncritical assertion");
433  returnCode = EXIT_FAILURE;
434  } else {
435  // Print final timestep information
436  repositories::printTimeStep(false);
437  logInfo("main()", "Terminated successfully");
438  returnCode = EXIT_SUCCESS;
439  }
440 
441  reportRuntime(measurement);
442 
443  delete[] cellCentresHost;
444  delete[] cellSizeHost;
445  tarch::freeMemory(QOutFlattened, tarch::MemoryLocation::Pinned);
446  tarch::freeMemory(QInFlattened, tarch::MemoryLocation::Pinned);
447 
448  peano4::shutdown();
449 
450  return returnCode;
451 }
constexpr double timeStamp
int main(int argc, char **argv)
tarch::logging::Log _log("::")
const tarch::la::Vector< DIMENSIONS, double > CellSize
tarch::timing::Measurement timeStepMeasurement
Definition: ccz4-main.cpp:58
A singleton repository to manage and provide access to grid cell and face data.
double * getFaceQNewDevice(const tarch::la::Vector< DIMENSIONS, double > &cellCentre, const tarch::la::Vector< DIMENSIONS, double > &cellSize, int axis)
peano4::datamanagement::CellMarker createMarkerFromIndex(const tarch ::la ::Vector< DIMENSIONS, int > &index)
void applyBoundaryConditionsToAxis(int axis, double t, double dt)
void applyBoundaryConditions()
void applyInitialConditions()
void initGrid()
void plotGrid()
void invokeFVSolverBoundaryConditions(const double *NOALIAS Qinside, double *NOALIAS Qoutside, const tarch::la::Vector< DIMENSIONS, double > &x, const tarch::la::Vector< DIMENSIONS, double > &h, double t, double dt, int normal)
void createGrid()
void projectPatchOntoFaces()
int depth
Definition: acoustic.py:29
void reportRuntime(const std::string &kernelIdentificator, const tarch::timing::Measurement &timingComputeKernel, const tarch::timing::Measurement &timingKernelLaunch, int numberOfCells, int numberOfThreads, tarch::logging::Log _log)
Reports the runtime and throughput of the benchmarks.
Definition: Utils.h:68
string msg
parameter setting according to scenarios
h
Definition: swe.py:79
static int getFaceIndex(const tarch::la::Vector< DIMENSIONS, double > &cellCentre, const tarch::la::Vector< DIMENSIONS, double > &cellSize, int axis)