SPH
main.cpp
Go to the documentation of this file.
1 
5 #include "io/FileManager.h"
6 #include "io/FileSystem.h"
7 #include "io/Logger.h"
8 #include "io/Output.h"
9 #include "objects/Exceptions.h"
12 #include "physics/Functions.h"
13 #include "post/Analysis.h"
14 #include "post/StatisticTests.h"
15 #include "quantities/Quantity.h"
16 #include "quantities/Storage.h"
17 #include "sph/initial/Initial.h"
18 #include "system/Factory.h"
19 #include "system/Process.h"
20 #include "system/Statistics.h"
21 #include "thread/Pool.h"
22 #include <fstream>
23 #include <iostream>
24 #include <mutex>
25 
26 using namespace Sph;
27 
28 static Expected<Storage> parsePkdgravOutput(const Path& path) {
29  Storage storage;
30  Statistics stats;
31  PkdgravInput io;
32  Outcome result = io.load(path, storage, stats);
33  if (!result) {
34  return makeUnexpected<Storage>(result.error());
35  } else {
36  return storage;
37  }
38 }
39 
40 int pkdgravToSfd(const Path& filePath, const Path& sfdPath) {
41  std::cout << "Processing pkdgrav file ... " << std::endl;
42  Expected<Storage> storage = parsePkdgravOutput(filePath);
43  if (!storage) {
44  std::cout << "Invalid file: " << storage.error() << std::endl;
45  return 0;
46  }
47  Post::HistogramParams params;
50  FileLogger logRadiiSfd(sfdPath, FileLogger::Options::KEEP_OPENED);
51  for (Post::HistPoint& p : sfd) {
52  logRadiiSfd.write(p.value, " ", p.count);
53  }
54  return 0;
55 }
56 
57 int pkdgravToOmega(const Path& filePath, const Path& omegaPath) {
58  std::cout << "Processing pkdgrav file ... " << std::endl;
59  Expected<Storage> storage = parsePkdgravOutput(filePath);
60  if (!storage) {
61  std::cout << "Invalid file: " << storage.error() << std::endl;
62  return 0;
63  }
64  /*Post::HistogramParams params;
65  params.source = Post::HistogramParams::Source::PARTICLES;
66  params.id = Post::HistogramId::ANGULAR_VELOCITIES;
67  params.binCnt = 50;
68  params.validator = [](const Float value) { return value > 0._f; };*/
69 
70  // Array<Post::SfdPoint> sfd = Post::getDifferentialSfd(storage.value(), params);
71  FileLogger logOmegaSfd(omegaPath, FileLogger::Options::KEEP_OPENED);
72  /*for (Post::SfdPoint& p : sfd) {
73  logOmegaSfd.write(p.value, " ", p.count);
74  }*/
75  ArrayView<Vector> omega = storage->getValue<Vector>(QuantityId::ANGULAR_FREQUENCY);
76  std::sort(
77  omega.begin(), omega.end(), [](Vector& v1, Vector& v2) { return getLength(v1) > getLength(v2); });
78  for (Vector v : omega) {
79  if (getLength(v) != 0._f) {
80  logOmegaSfd.write(getLength(v));
81  }
82  }
83  return 0;
84 }
85 
86 int pkdgravToMoons(const Path& filePath, const Float limit) {
87  std::cout << "Processing pkdgrav file ... " << std::endl;
88  Expected<Storage> storage = parsePkdgravOutput(filePath);
89  if (!storage) {
90  std::cout << "Invalid file: " << storage.error() << std::endl;
91  return 0;
92  }
94  Array<Post::MoonEnum> moons = Post::findMoons(storage.value(), 1.2_f, limit);
95  Size moonCnt = std::count(moons.begin(), moons.end(), Post::MoonEnum::MOON);
96  std::cout << "Moon count = " << moonCnt << std::endl;
97  return 0;
98 }
99 
100 int ssfToSfd(const Post::HistogramSource source, const Path& filePath, const Path& sfdPath) {
101  std::cout << "Processing SPH file ... " << std::endl;
102  AutoPtr<IInput> input = Factory::getInput(filePath);
103  Storage storage;
104  Statistics stats;
105  Outcome outcome = input->load(filePath, storage, stats);
106  if (!outcome) {
107  std::cout << "Cannot load particle data, " << outcome.error() << std::endl;
108  return 0;
109  }
110 
111  Post::HistogramParams params;
115  for (Post::HistPoint& p : sfd) {
116  logSfd.write(p.value, " ", p.count);
117  }
118  return 0;
119 }
120 
121 int ssfToOmega(const Path& filePath,
122  const Path& omegaPath,
123  const Path& omegaDPath,
124  const Path& omegaDirPath) {
125  std::cout << "Processing SPH file ... " << std::endl;
126  BinaryInput input;
127  Storage storage;
128  Statistics stats;
129  Outcome outcome = input.load(filePath, storage, stats);
130  if (!outcome) {
131  std::cout << "Cannot load particle data, " << outcome.error() << std::endl;
132  return 0;
133  }
134 
135  Post::HistogramParams params;
136  params.range = Interval(0._f, 13._f);
137  params.binCnt = 12;
138 
141  // const Float massCutoff = 1._f / 300000._f;
142  const Float m_total = std::accumulate(m.begin(), m.end(), 0._f);
143  params.validator = [w](const Size i) { //
144  return getSqrLength(w[i]) > 0._f; //= m_total * massCutoff;
145  };
146 
147  params.centerBins = false;
148 
151 
152  FileLogger logOmega(omegaPath, FileLogger::Options::KEEP_OPENED);
153  for (Post::HistPoint& p : sfd) {
154  logOmega.write(p.value, " ", p.count); // / sum);
155  }
156 
157  params.range = Interval();
160 
161  FileLogger logOmegaDir(omegaDirPath, FileLogger::Options::KEEP_OPENED);
162  for (Post::HistPoint& p : dirs) {
163  logOmegaDir.write(p.value, " ", p.count);
164  }
165 
166  // print omega-D relation
167  Array<Float> h(storage.getParticleCnt());
170  for (Size i = 0; i < m.size(); ++i) {
171  h[i] = r[i][H]; // root<3>(3.f * m[i] / (rho[i] * 4.f * PI));
172  }
173 
174 
175  FileLogger logOmegaD(omegaDPath, FileLogger::Options::KEEP_OPENED);
176  for (Size i = 0; i < m.size(); ++i) {
177  if (m[i] > 3._f * params.massCutoff * m_total) {
178  // in m vs. rev/day
179  logOmegaD.write(2._f * h[i], " ", getLength(w[i]) * 3600._f * 24._f / (2._f * PI));
180  }
181  }
182 
183  return 0;
184 }
185 
186 
187 int ssfToVelocity(const Path& filePath, const Path& outPath) {
188  std::cout << "Processing SPH file ... " << std::endl;
189  AutoPtr<IInput> input = Factory::getInput(filePath);
190  Storage storage;
191  Statistics stats;
192  Outcome outcome = input->load(filePath, storage, stats);
193  if (!outcome) {
194  std::cout << "Cannot load particle data, " << outcome.error() << std::endl;
195  return -1;
196  }
197 
198  // convert to system with center at LR
199  Array<Size> idxs = Post::findLargestComponent(storage, 2._f, EMPTY_FLAGS);
201  ArrayView<Vector> r, v, dv;
202  tie(r, v, dv) = storage.getAll<Vector>(QuantityId::POSITION);
203  Vector r0(0._f);
204  Vector v0(0._f);
205  Float m0 = 0._f;
206  for (Size i : idxs) {
207  m0 += m[i];
208  r0 += m[i] * r[i];
209  v0 += m[i] * v[i];
210  }
211  r0 /= m0;
212  v0 /= m0;
213 
214  for (Size i = 0; i < r.size(); ++i) {
215  r[i] -= r0;
216  v[i] -= v0;
217  }
218 
219  Post::HistogramParams params;
220  params.binCnt = 2000;
223 
225  for (Post::HistPoint& p : hist) {
226  logSfd.write(p.value, " ", p.count);
227  }
228 
229  return 0;
230 }
231 
232 void ssfToVelDir(const Path& filePath, const Path& outPath) {
233  std::cout << "Processing SPH file ... " << std::endl;
234  BinaryInput input;
235  Storage storage;
236  Statistics stats;
237  Outcome outcome = input.load(filePath, storage, stats);
238  if (!outcome) {
239  std::cout << "Cannot load particle data, " << outcome.error() << std::endl;
240  return;
241  }
242 
243  // convert to system with center at LR
244  Array<Size> idxs = Post::findLargestComponent(storage, 2._f, EMPTY_FLAGS);
246  ArrayView<Vector> r, v, dv;
247  tie(r, v, dv) = storage.getAll<Vector>(QuantityId::POSITION);
248  Vector r0(0._f);
249  Vector v0(0._f);
250  Float m0 = 0._f;
251  for (Size i : idxs) {
252  m0 += m[i];
253  r0 += m[i] * r[i];
254  v0 += m[i] * v[i];
255  }
256  r0 /= m0;
257  v0 /= m0;
258 
259  for (Size i = 0; i < r.size(); ++i) {
260  r[i] -= r0;
261  v[i] -= v0;
262  }
263 
264  // compute velocity directions (in x-y plane)
265  Array<Float> dirs;
266  for (Size i = 0; i < v.size(); ++i) {
267  Float dir = atan2(v[i][Y], v[i][X]);
268  if (dir < 0._f) {
269  dir += 2._f * PI;
270  }
271  dirs.push(dir * RAD_TO_DEG);
272  }
273 
274  Post::HistogramParams params;
275  params.range = Interval(0._f, 360._f);
276  params.binCnt = 72; // 5 deg bins
279  for (Post::HistPoint& p : hist) {
280  logSfd.write(p.value, " ", p.count);
281  }
282 }
283 
286  std::string name;
289 };
290 
291 static Array<HarrisAsteroid> loadHarris(std::ifstream& ifs) {
292  Array<HarrisAsteroid> harris;
293  while (ifs) {
294  std::string dummy;
295  std::string number;
296  ifs >> number >> dummy;
297 
298  std::string name;
299  ifs >> name;
300  for (int i = 0; i < 6; ++i) {
301  ifs >> dummy;
302  }
303 
304  std::string radius;
305  ifs >> radius;
306  for (int i = 0; i < 5; ++i) {
307  ifs >> dummy;
308  }
309 
310  std::string period;
311  ifs >> period;
312  for (int i = 0; i < 10; ++i) {
313  ifs >> dummy;
314  }
315 
316  harris.push(HarrisAsteroid{
317  fromString<Size>(number),
318  name,
319  fromString<Float>(radius),
320  fromString<Float>(period),
321  });
322  }
323  return harris;
324 }
325 
329 };
330 
331 static Array<FamilyAsteroid> loadFamilies(std::ifstream& ifs) {
332  SPH_ASSERT(ifs);
333  Array<FamilyAsteroid> asteroids;
334  std::string line;
335  bool firstLine = true;
336  int format = 1;
337  while (std::getline(ifs, line)) {
338  if (line.empty() || line[0] == '#') {
339  // comment line
340  if (firstLine) {
341  // this is the other format of the file, with asteroid names, etc.
342  format = 2;
343  }
344  firstLine = false;
345  continue;
346  }
347  firstLine = false;
348  std::stringstream ss(line);
349  // both formats start with asteroid number
350  std::string number;
351  ss >> number;
352  if (format == 2) {
353  std::string name;
354  ss >> name;
355  asteroids.push(FamilyAsteroid{ fromString<Size>(number), name });
356  } else {
357  asteroids.push(FamilyAsteroid{ fromString<Size>(number), NOTHING });
358  }
359  // check that we have at least one information
360  SPH_ASSERT(asteroids.back().name || asteroids.back().number);
361  }
362  return asteroids;
363 }
364 
365 static Optional<HarrisAsteroid> findInHarris(const FamilyAsteroid& ast1,
367  auto iter = std::find_if(catalog.begin(), catalog.end(), [&ast1](const HarrisAsteroid& ast2) { //
368  // search primarily by number
369  if (ast1.number && ast2.number && ast1.number.value() == ast2.number.value()) {
370  return true;
371  }
372  // if we don't have the number, search by name
373  if (ast1.name && ast1.name.value() == ast2.name) {
374  return true;
375  }
376  // either don't match or we don't have the information
377  return false;
378  });
379  if (iter == catalog.end()) {
380  return NOTHING;
381  } else {
382  if (iter->period && iter->radius) {
383  return *iter;
384  } else {
385  return NOTHING;
386  }
387  }
388 }
389 
390 static Optional<Post::KsResult> printDvsOmega(const Path& familyData,
391  const Path& outputPath,
393  Array<PlotPoint>& outPoints) {
394 
395  std::ifstream ifs(familyData.native());
396  SPH_ASSERT(ifs);
397  Array<FamilyAsteroid> family = loadFamilies(ifs);
398  Array<HarrisAsteroid> found;
399  Interval rangePeriod;
400  Interval rangeR;
401  for (FamilyAsteroid& famAst : family) {
402  if (auto harAst = findInHarris(famAst, catalog)) {
403  found.push(harAst.value());
404  rangePeriod.extend(harAst->period.value());
405  rangeR.extend(harAst->radius.value());
406  }
407  }
408  if (found.size() < 10) {
409  // too few data to do any statistics
410  return NOTHING;
411  }
412 
414  std::ofstream ofs(outputPath.native());
415  SPH_ASSERT(ofs);
416 
417  auto compare = [](HarrisAsteroid& ast1, HarrisAsteroid& ast2) {
418  return ast1.radius.value() < ast2.radius.value();
419  };
420  Iterator<HarrisAsteroid> largestRemnant = std::max_element(found.begin(), found.end(), compare);
421 
422  auto periodToOmega = [](const Float P) { return 1._f / (P / 24._f); };
423 
424  Array<PlotPoint> points;
425  for (HarrisAsteroid& ast : found) {
426  const Float omega = periodToOmega(ast.period.value());
427  if (&ast != &*largestRemnant) {
428  std::string printedName = ast.name;
429  if (ast.number) {
430  printedName = "(" + std::to_string(ast.number.value()) + ") " + printedName;
431  }
432  ofs << ast.radius.value() << " " << omega << " " << printedName << std::endl;
433  }
434  points.push(PlotPoint(ast.radius.value(), omega));
435  }
436  Interval rangeOmega;
437  rangeOmega.extend(periodToOmega(rangePeriod.lower()));
438  rangeOmega.extend(periodToOmega(rangePeriod.upper()));
439 
440  Post::KsFunction uniformCdf = Post::getUniformKsFunction(rangeR, rangeOmega);
441  Post::KsResult result = Post::kolmogorovSmirnovTest(points, uniformCdf);
442 
443 
444  if (points.size() > 36) {
445  const Path histPath = Path("histogram") / outputPath.fileName();
447  std::ofstream histofs(histPath.native());
448  Array<Float> values;
449  for (PlotPoint& p : points) {
450  values.push(p.y);
451  }
452  Post::HistogramParams params;
453  params.range = Interval(0._f, 13._f);
454  Array<Post::HistPoint> histogram = Post::getDifferentialHistogram(values, params);
455  for (Post::HistPoint p : histogram) {
456  histofs << p.value << " " << p.count << std::endl;
457  }
458  }
459 
460  outPoints = std::move(points);
461  return result;
462 }
463 
464 static std::string elliptize(const std::string& s, const Size maxSize) {
465  SPH_ASSERT(maxSize >= 5);
466  if (s.size() < maxSize) {
467  return s;
468  }
469  const Size cutSize = (maxSize - 3) / 2;
470  return s.substr(0, cutSize) + "..." + s.substr(s.size() - cutSize);
471 }
472 
474  Path harrisPath("/home/pavel/projects/astro/asteroids/grant3/harris.out");
475  std::ifstream ifs(harrisPath.native());
476  Array<HarrisAsteroid> harris = loadHarris(ifs);
477  ifs.close();
478 
479  Size below3 = 0, below7 = 0, below12 = 0, total = 0;
480  Float avgPeriod = 0._f;
481  Size count = 0;
482  for (auto& h : harris) {
483  if (!h.period) {
484  continue;
485  }
486 
487  if (h.period.value() < 3) {
488  below3++;
489  }
490  if (h.period.value() < 7) {
491  below7++;
492  }
493  if (h.period.value() < 12) {
494  below12++;
495  }
496 
497  if (h.radius.valueOr(0._f) > 100._f) {
498  avgPeriod += h.period.value();
499  count++;
500  }
501  total++;
502  }
503  std::cout << "Below 3h: " << (100.f * below3) / total << "%" << std::endl;
504  std::cout << "Below 7h: " << (100.f * below7) / total << "%" << std::endl;
505  std::cout << "Below 12h: " << (100.f * below12) / total << "%" << std::endl;
506 
507  if (count > 0) {
508  std::cout << "Average period for R>100km: " << avgPeriod / count << std::endl;
509  } else {
510  SPH_ASSERT(false);
511  }
512 
513  Array<PlotPoint> points;
514  printDvsOmega(
515  Path("/home/pavel/projects/astro/asteroids/families.txt"), Path("LR_D_omega.txt"), harris, points);
516 
517  const Path parentPath("/home/pavel/projects/astro/asteroids/families");
518  Array<Path> paths = FileSystem::getFilesInDirectory(parentPath);
519  std::mutex printMutex;
520  std::mutex plotMutex;
521  UniquePathManager uniquePaths;
522  std::ofstream kss("KS.txt");
523  kss << "# name D_ks probability" << std::endl;
524 
525  std::ofstream alls("D_omega_all.txt");
526 
528  parallelFor(pool, 0, paths.size(), [&](const Size index) {
529  {
530  std::unique_lock<std::mutex> lock(printMutex);
531  std::cout << paths[index].native() << std::endl;
532  }
533  const Path name = paths[index].fileName();
534  const Path targetPath = Path("D_omega") / Path(name.native() + ".txt");
535  const Path ksPath = Path("KS") / Path(name.native() + ".txt");
536  Array<PlotPoint> points;
537  Optional<Post::KsResult> ks = printDvsOmega(parentPath / paths[index], targetPath, harris, points);
538 
539  std::unique_lock<std::mutex> lock(plotMutex);
540  if (!FileSystem::pathExists(targetPath) || FileSystem::fileSize(targetPath) == 0) {
541  return;
542  }
543  if (ks) {
544  kss << std::left << std::setw(35) << elliptize(Path(name).removeExtension().native(), 30)
545  << std::setw(15) << ks->D << std::setw(15) << ks->prob << std::endl;
546  }
547  for (PlotPoint p : points) {
548  alls << p.x << " " << p.y << " " << index << std::endl;
549  }
550 
551  FileSystem::copyFile(targetPath, Path("family.txt"));
552  // make plot
553  Process gnuplot(Path("/bin/gnuplot"), { "doplot.plt" });
554  gnuplot.wait();
555  SPH_ASSERT(FileSystem::pathExists(Path("plot.png")));
556  FileSystem::copyFile(Path("plot.png"), uniquePaths.getPath(Path(targetPath).replaceExtension("png")));
557 
558  if (points.size() > 25) {
559  const Path histPath = Path("histogram") / targetPath.fileName();
560  FileSystem::copyFile(histPath, Path("hist.txt"));
561  Process gnuplot2(Path("/bin/gnuplot"), { "dohistogram.plt" });
562  gnuplot2.wait();
563  SPH_ASSERT(FileSystem::pathExists(Path("hist.png")));
565  Path("hist.png"), uniquePaths.getPath(Path(histPath).replaceExtension("png")));
566  }
567  });
568 
569  FileSystem::copyFile(Path("D_omega_all.txt"), Path("family.txt"));
570  Process gnuplot(Path("/bin/gnuplot"), { "doplot_all.plt" });
571  gnuplot.wait();
572 }
573 
574 static Float maxwellBoltzmann(const Float x, const Float a) {
575  return sqrt(2._f / PI) * sqr(x) * exp(-sqr(x) / (2._f * sqr(a))) / pow<3>(a);
576 }
577 
578 static Float sampleMaxwellBoltzmann(UniformRng& rng, const Float a) {
579  while (true) {
580  const Float x = rng() * a * 10._f;
581  const Float y = rng() / a;
582 
583  if (maxwellBoltzmann(x, a) > y) {
584  return x;
585  }
586  }
587 }
588 
589 void makeSwift(const Path& filePath) {
590  // for Hygiea
591  /*const Float a = 3.14178_f;
592  const Float e = 0.135631_f;
593  const Float I = asin(0.0889622);
594  const Float W = 64.621768_f * DEG_TO_RAD;
595  const Float w = 128.543611_f * DEG_TO_RAD;
596  const Float u = 0._f;
597 
598  const Float X = (cos(W) * cos(w) - sin(W) * cos(I) * sin(w)) * a * (cos(u) - e) -
599  (cos(W) * sin(w) + sin(W) * cos(I) * cos(w)) * a * sqrt(1 - sqr(e)) * sin(u);
600  const Float Y = (sin(W) * cos(w) + cos(W) * cos(I) * sin(w)) * a * (cos(u) - e) +
601  (-sin(W) * sin(w) + cos(W) * cos(I) * cos(w)) * a * sqrt(1 - sqr(e)) * sin(u);
602  const Float Z = sin(I) * sin(w) * a * (cos(u) - e) + sin(I) * cos(w) * a * sqrt(1 - sqr(e)) * sin(u);
603 
604  const Vector r = Vector(X, Y, Z);
605 
606  const Float npart = 1500;
607 
608  BinaryInput input;
609  Storage storage;
610  Statistics stats;
611  if (!input.load(filePath, storage, stats)) {
612  std::cout << "Cannot parse ssf file" << std::endl;
613  }
614 
615  ArrayView<const Vector> v = storage.getDt<Vector>(QuantityId::POSITION);
616 
617  FileLogger logger(Path("tp.in"));
618  UniformRng rng;
619  logger.write(npart);
620  for (Size i = 0; i < npart; ++i) {
621  logger.write(r / Constants::au);
622  const Size idx = clamp<Size>(Size(rng() * v.size()), 0, storage.getParticleCnt() - 1);
623  logger.write(v[idx] / Constants::au * Constants::year);
624  logger.write("0");
625  logger.write("0.0");
626  }*/
627  std::ifstream ifs(filePath.native());
628  std::string line;
629 
630  Array<Float> rs;
631  while (std::getline(ifs, line)) {
632  const Float d = std::stof(line);
633  rs.push(d / 2 * 1000);
634  }
635 
636  std::ofstream yarko("yarko.in");
637  yarko << rs.size() << std::endl;
638  for (Float r : rs) {
639  yarko << r << " 2860.0 1500.0 0.0010 680.0 0.10 0.90" << std::endl;
640  }
641 
642  std::ofstream spin("spin.in");
643  spin << rs.size() << "\n-1\n1\n";
644  UniformRng rng;
645  for (Size i = 0; i < rs.size(); ++i) {
646  const Float phi = rng() * 2._f * PI;
647  const Float cosTheta = rng() * 2._f - 1._f;
648  const Float theta = acos(cosTheta);
649  spin << sphericalToCartesian(1._f, theta, phi) << " " << sampleMaxwellBoltzmann(rng, 0.0001_f)
650  << std::endl;
651  }
652 
653  std::ofstream yorp("yorp.in");
654  for (Size i = 1; i <= rs.size(); ++i) {
655  yorp << i << " " << clamp(int(rng() * 200), 0, 199) << std::endl;
656  }
657 }
658 
659 void origComponents(const Path& lastDumpPath, const Path& firstDumpPath, const Path& colorizedDumpPath) {
660  BinaryInput input;
661  Storage lastDump, firstDump;
662  Statistics stats;
663  Outcome res1 = input.load(lastDumpPath, lastDump, stats);
664  Outcome res2 = input.load(firstDumpPath, firstDump, stats);
665  if (!res1 || !res2) {
666  throw IoError((res1 || res2).error());
667  }
668 
669  // use last dump to find components
670  Array<Size> components;
673 
674  // "colorize" the flag quantity using the components
675  SPH_ASSERT(firstDump.getParticleCnt() == components.size());
676  firstDump.getValue<Size>(QuantityId::FLAG) = components.clone();
677 
678  // save as new file
679  BinaryOutput output(colorizedDumpPath);
680  output.dump(firstDump, stats);
681 }
682 
683 void extractLr(const Path& inputPath, const Path& outputPath) {
684  BinaryInput input;
685  Storage storage;
686  Statistics stats;
687  Outcome res = input.load(inputPath, storage, stats);
688  if (!res) {
689  throw IoError(res.error());
690  }
691 
692  // allow using this for storage without masses --> add ad hoc mass if it's missing
693  if (!storage.has(QuantityId::MASS)) {
694  storage.insert<Float>(QuantityId::MASS, OrderEnum::ZERO, 1._f);
695  }
696 
697  Array<Size> components;
698  const Size componentCnt =
699  Post::findComponents(storage, 1.5_f, Post::ComponentFlag::SORT_BY_MASS, components);
700  std::cout << "Component cnt = " << componentCnt << std::endl;
701 
702  Array<Size> toRemove;
703  for (Size i = 0; i < components.size(); ++i) {
704  if (components[i] != 0) {
705  // not LR
706  toRemove.push(i);
707  }
708  }
709  storage.remove(toRemove, Storage::IndicesFlag::INDICES_SORTED);
710 
715 
716  BinaryOutput output(outputPath);
717  output.dump(storage, stats);
718 
719  if (storage.has(QuantityId::DENSITY)) {
721 
722  Float volume = 0._f;
723  for (Size i = 0; i < m.size(); ++i) {
724  volume += m[i] / rho[i];
725  }
726 
727  std::cout << "eq. diameter = " << Sph::cbrt(3._f * volume / (4._f * PI)) * 2._f / 1000._f << "km"
728  << std::endl;
729  }
730 
731  const Float omega = getLength(Post::getAngularFrequency(m, r, v));
732  std::cout << "period = " << 2._f * PI / omega / 3600._f << "h" << std::endl;
733 
735  SPH_ASSERT(isReal(I), I);
736  Eigen e = eigenDecomposition(I);
737  /*std::cout << "I = " << I << std::endl;
738  std::cout << "matrix = " << e.vectors << std::endl;
739  std::cout << "values = " << e.values << std::endl;*/
740  const Float A = e.values[2];
741  const Float B = e.values[1];
742  const Float C = e.values[0];
743  const Float a = sqrt(B + C - A);
744  const Float b = sqrt(A + C - B);
745  const Float c = sqrt(A + B - C);
746  SPH_ASSERT(a > 0._f && b > 0._f && c > 0._f, a, b, c);
747  std::cout << "a/b = " << a / b << std::endl;
748  std::cout << "b/c = " << b / c << std::endl;
749 }
750 
751 void printHelp() {
752  std::cout << "Expected usage: post mode [parameters]" << std::endl
753  << " where 'mode' is one of:" << std::endl
754  << " - pkdgravToSfd - computes the cumulative SFD from pkdgrav output file" << std::endl
755  << " - pkdgravToOmega - computes the spin rate distribution from pkdgrav output file"
756  << std::endl
757  << " - pkdgravToMoons - finds satellites of the largest remnant (fragment) from pkdgrav "
758  "output file"
759  << std::endl
760  << "- ssfToSfd - computes the cumulative SFD from SPH output file" << std::endl
761  << "- ssfToVelocity - computes the velocity distribution from SPH output file" << std::endl
762  << "- harris - TODO" << std::endl
763  << "- stats - prints ejected mass and the period of the largest remnant" << std::endl
764  << "- swift - makes yarko.in, yorp.in and spin.in input file for swift" << std::endl;
765 }
766 
767 int main(int argc, char** argv) {
768  if (argc < 2) {
769  printHelp();
770  return 0;
771  }
772  try {
773 
774  std::string mode(argv[1]);
775  if (mode == "pkdgravToSfd") {
776  if (argc < 4) {
777  std::cout << "Expected parameters: post pkdgravToSfd ss.50000.bt sfd.txt";
778  return 0;
779  }
780  return pkdgravToSfd(Path(argv[2]), Path(argv[3]));
781  } else if (mode == "pkdgravToOmega") {
782  if (argc < 4) {
783  std::cout << "Expected parameters: post pkdgravToOmega ss.50000.bt omega.txt";
784  return 0;
785  }
786  return pkdgravToOmega(Path(argv[2]), Path(argv[3]));
787  } else if (mode == "pkdgravToMoons") {
788  if (argc < 4) {
789  std::cout << "Expected parameters: post pkdgravToMoons ss.50000.bt 0.1";
790  return 0;
791  }
792  const Float limit = std::atof(argv[3]);
793  return pkdgravToMoons(Path(argv[2]), limit);
794  } else if (mode == "ssfToSfd") {
795  if (argc < 4) {
796  std::cout << "Expected parameters: post ssfToSfd [--components] output.ssf sfd.txt";
797  return 0;
798  }
799  if (std::string(argv[2]) == "--components") {
800  return ssfToSfd(Post::HistogramSource::COMPONENTS, Path(argv[3]), Path(argv[4]));
801  } else {
802  return ssfToSfd(Post::HistogramSource::PARTICLES, Path(argv[2]), Path(argv[3]));
803  }
804  } else if (mode == "ssfToVelocity") {
805  return ssfToVelocity(Path(argv[2]), Path(argv[3]));
806  } else if (mode == "ssfToOmega") {
807  if (argc < 6) {
808  std::cout << "Expected parameters: post ssfToOmega output.ssf omega.txt omega_D.txt "
809  "omega_dir.txt";
810  return 0;
811  }
812  return ssfToOmega(Path(argv[2]), Path(argv[3]), Path(argv[4]), Path(argv[5]));
813  } else if (mode == "ssfToVelDir") {
814  if (argc < 4) {
815  std::cout << "Expected parameters: post ssfToVelDir output.ssf veldir.txt" << std::endl;
816  return 0;
817  }
818  ssfToVelDir(Path(argv[2]), Path(argv[3]));
819  } else if (mode == "harris") {
821  } else if (mode == "swift") {
822  if (argc < 3) {
823  std::cout << "Expected parameters: post maketp D.dat";
824  return 0;
825  }
826  makeSwift(Path(argv[2]));
827  } else if (mode == "origComponents") {
828  if (argc < 5) {
829  std::cout << "Expected parameters: post origComponents lastDump.ssf firstDump.ssf "
830  "colorizedDump.ssf"
831  << std::endl;
832  }
833  origComponents(Path(argv[2]), Path(argv[3]), Path(argv[4]));
834  } else if (mode == "extractLr") {
835  if (argc < 4) {
836  std::cout << "Expected parameters: post extractLr input.ssf lr.ssf" << std::endl;
837  }
838  extractLr(Path(argv[2]), Path(argv[3]));
839  } else {
840  printHelp();
841  return 0;
842  }
843 
844  } catch (const std::exception& e) {
845  std::cout << "ERROR: " << e.what() << std::endl;
846  return -1;
847  }
848 }
Various function for interpretation of the results of a simulation.
INLINE bool isReal(const AntisymmetricTensor &t)
#define SPH_ASSERT(x,...)
Definition: Assert.h:94
const float radius
Definition: CurveDialog.cpp:18
const EmptyFlags EMPTY_FLAGS
Definition: Flags.h:16
uint32_t Size
Integral type used to index arrays (by default).
Definition: Globals.h:16
double Float
Precision used withing the code. Use Float instead of float or double where precision is important.
Definition: Globals.h:13
Vector moveToCenterOfMassSystem(ArrayView< const Float > m, ArrayView< Vector > r)
Modifies particle positions so that their center of mass lies at the origin.
Definition: Initial.cpp:363
Generating initial conditions of SPH particles.
Logging routines of the run.
constexpr INLINE T clamp(const T &f, const T &f1, const T &f2)
Definition: MathBasic.h:35
constexpr INLINE T sqr(const T &f) noexcept
Return a squared value.
Definition: MathUtils.h:67
INLINE Float cbrt(const Float f)
Returns a cubed root of a value.
Definition: MathUtils.h:84
INLINE T sqrt(const T f)
Return a squared root of a value.
Definition: MathUtils.h:78
constexpr INLINE Float pow< 3 >(const Float v)
Definition: MathUtils.h:132
INLINE T exp(const T f)
Definition: MathUtils.h:269
INLINE T acos(const T f)
Definition: MathUtils.h:306
INLINE T atan2(const T y, const T x)
Definition: MathUtils.h:321
constexpr Float RAD_TO_DEG
Definition: MathUtils.h:371
constexpr Float PI
Mathematical constants.
Definition: MathUtils.h:361
const NothingType NOTHING
Definition: Optional.h:16
Saving and loading particle data.
Simple thread pool with fixed number of threads.
Creating and managing processes.
@ FLAG
ID of original body, used to implement discontinuities between bodies in SPH.
@ POSITION
Positions (velocities, accelerations) of particles, always a vector quantity,.
@ DENSITY
Density, always a scalar quantity.
@ MASS
Paricles masses, always a scalar quantity.
Holder of quantity values and their temporal derivatives.
@ ZERO
Quantity without derivatives, or "zero order" of quantity.
INLINE void parallelFor(IScheduler &scheduler, const Size from, const Size to, TFunctor &&functor)
Executes a functor concurrently from all available threads.
Definition: Scheduler.h:98
StaticArray< T0 &, sizeof...(TArgs)+1 > tie(T0 &t0, TArgs &... rest)
Creates a static array from a list of l-value references.
Definition: StaticArray.h:281
Statistics gathered and periodically displayed during the run.
Container for storing particle quantities and materials.
Eigen eigenDecomposition(const SymmetricTensor &t)
Computes eigenvectors and corresponding eigenvalues of symmetric matrix.
INLINE Float getLength(const Vector &v)
Returns the length of the vector. Enabled only for vectors of floating-point precision.
Definition: Vector.h:579
INLINE Float getSqrLength(const Vector &v)
Definition: Vector.h:574
INLINE Vector sphericalToCartesian(const Float r, const Float theta, const Float phi)
Constructs a vector from spherical coordinates.
Definition: Vector.h:788
@ H
Definition: Vector.h:25
@ Y
Definition: Vector.h:23
@ X
Definition: Vector.h:22
int main(int argc, char *argv[])
Definition: main.cpp:7
Object providing safe access to continuous memory of data.
Definition: ArrayView.h:17
INLINE Iterator< StorageType > begin()
Definition: ArrayView.h:55
INLINE Iterator< StorageType > end()
Definition: ArrayView.h:63
INLINE Iterator< StorageType > end() noexcept
Definition: Array.h:462
INLINE void push(U &&u)
Adds new element to the end of the array, resizing the array if necessary.
Definition: Array.h:306
INLINE T & back() noexcept
Definition: Array.h:176
INLINE TCounter size() const noexcept
Definition: Array.h:193
INLINE Iterator< StorageType > begin() noexcept
Definition: Array.h:450
Array clone() const
Performs a deep copy of all elements of the array.
Definition: Array.h:147
Wrapper of pointer that deletes the resource from destructor.
Definition: AutoPtr.h:15
INLINE const TError & error() const
Returns the error message.
Definition: Outcome.h:88
Input for the binary file, generated by BinaryOutput.
Definition: Output.h:352
virtual Outcome load(const Path &path, Storage &storage, Statistics &stats) override
Loads data from the file into the storage.
Definition: Output.cpp:718
Output saving data to binary data without loss of precision.
Definition: Output.h:335
virtual Expected< Path > dump(const Storage &storage, const Statistics &stats) override
Saves data from particle storage into the file.
Definition: Output.cpp:512
Wrapper of type that either contains a value of given type, or an error message.
Definition: Expected.h:25
const Error & error() const
Returns the error message.
Definition: Expected.h:94
Type & value()
Returns the reference to expected value.
Definition: Expected.h:69
File output logger.
Definition: Logger.h:160
void write(TArgs &&... args)
Creates and logs a message by concatenating arguments.
Definition: Logger.h:37
Object representing a 1D interval of real numbers.
Definition: Interval.h:17
INLINE Float lower() const
Returns lower bound of the interval.
Definition: Interval.h:74
INLINE void extend(const Float &value)
Extends the interval to contain given value.
Definition: Interval.h:41
INLINE Float upper() const
Returns upper bound of the interval.
Definition: Interval.h:79
Exception thrown when file cannot be read, it has invalid format, etc.
Definition: Exceptions.h:31
Simple (forward) iterator over continuous array of objects of type T.
Definition: Iterator.h:18
INLINE Type & value()
Returns the reference to the stored value.
Definition: Optional.h:172
Object representing a path on a filesystem.
Definition: Path.h:17
std::string native() const
Returns the native version of the path.
Definition: Path.cpp:71
Path fileName() const
Returns the filename of the path.
Definition: Path.cpp:47
Path parentPath() const
Returns the parent directory. If the path is empty or root, return empty path.
Definition: Path.cpp:35
Input for pkdgrav output files (ss.xxxxx.bt).
Definition: Output.h:561
virtual Outcome load(const Path &path, Storage &storage, Statistics &stats) override
Loads data from the file into the storage.
Definition: Output.cpp:1391
Holds a handle to a created process.
Definition: Process.h:19
void wait()
Blocks the calling thread until the managed process exits. The function may block indefinitely.
Definition: Process.cpp:50
Object holding various statistics about current run.
Definition: Statistics.h:22
Container storing all quantities used within the simulations.
Definition: Storage.h:230
Quantity & insert(const QuantityId key, const OrderEnum order, const TValue &defaultValue)
Creates a quantity in the storage, given its key, value type and order.
Definition: Storage.cpp:270
StaticArray< Array< TValue > &, 3 > getAll(const QuantityId key)
Retrieves quantity buffers from the storage, given its key and value type.
Definition: Storage.cpp:163
Array< TValue > & getDt(const QuantityId key)
Retrieves a quantity derivative from the storage, given its key and value type.
Definition: Storage.cpp:217
Size getParticleCnt() const
Returns the number of particles.
Definition: Storage.cpp:445
@ INDICES_SORTED
Use if the given array is already sorted (optimization)
void remove(ArrayView< const Size > idxs, const Flags< IndicesFlag > flags=EMPTY_FLAGS)
Removes specified particles from the storage.
Definition: Storage.cpp:756
bool has(const QuantityId key) const
Checks if the storage contains quantity with given key.
Definition: Storage.cpp:130
Array< TValue > & getValue(const QuantityId key)
Retrieves a quantity values from the storage, given its key and value type.
Definition: Storage.cpp:191
Symmetric tensor of 2nd order.
Thread pool capable of executing tasks concurrently.
Definition: Pool.h:70
static SharedPtr< ThreadPool > getGlobalInstance()
Returns the global instance of the thread pool.
Definition: Pool.cpp:211
Random number generator with uniform distribution.
Definition: Rng.h:16
Object generating unique paths.
Definition: FileManager.h:13
Path getPath(const Path &expected)
Generates a unique path, based on given input path.
Definition: FileManager.cpp:15
Creating code components based on values from settings.
AutoPtr< IInput > getInput(const Path &path)
Definition: Factory.cpp:601
bool pathExists(const Path &path)
Checks if a file or directory exists (or more precisely, if a file or directory is accessible).
Definition: FileSystem.cpp:21
Array< Path > getFilesInDirectory(const Path &directory)
Alternatitve to iterateDirectory, returning all files in directory in an array.
Definition: FileSystem.cpp:342
Size fileSize(const Path &path)
Returns the size of a file.
Definition: FileSystem.cpp:29
Outcome createDirectory(const Path &path, const Flags< CreateDirectoryFlag > flags=CreateDirectoryFlag::ALLOW_EXISTING)
Creates a directory with given path. Creates all parent directories as well.
Definition: FileSystem.cpp:119
Outcome copyFile(const Path &from, const Path &to)
Copies a file on given path to a different path.
Definition: FileSystem.cpp:219
Array< MoonEnum > findMoons(const Storage &storage, const Float radius=1._f, const Float limit=0._f)
Find a potential satellites of the largest body.
Definition: Analysis.cpp:367
Array< HistPoint > getDifferentialHistogram(ArrayView< const Float > values, const HistogramParams &params)
Computes the differential histogram from given values.
Definition: Analysis.cpp:861
HistogramSource
Source data used to construct the histogram.
Definition: Analysis.h:217
@ COMPONENTS
Equivalent radii of connected chunks of particles (SPH framework)
@ PARTICLES
Radii of individual particles, considering particles as spheres (N-body framework)
KsFunction getUniformKsFunction(Interval rangeX, Interval rangeY)
Array< Size > findLargestComponent(const Storage &storage, const Float particleRadius, const Flags< ComponentFlag > flags)
Returns the indices of particles belonging to the largest remnant.
Definition: Analysis.cpp:216
Vector getAngularFrequency(ArrayView< const Float > m, ArrayView< const Vector > r, ArrayView< const Vector > v, const Vector &r0, const Vector &v0, ArrayView< const Size > idxs=nullptr)
Computes the immediate vector of angular frequency of a rigid body.
Definition: Analysis.cpp:514
KsResult kolmogorovSmirnovTest(ArrayView< const Float > data, const Function< Float(Float)> &expectedCdf)
One-dimensional Kolmogorov-Smirnov test with given CDF of expected probability distribution.
Size findComponents(const Storage &storage, const Float particleRadius, const Flags< ComponentFlag > flags, Array< Size > &indices)
Finds and marks connected components (a.k.a. separated bodies) in the array of vertices.
Definition: Analysis.cpp:86
SymmetricTensor getInertiaTensor(ArrayView< const Float > m, ArrayView< const Vector > r, const Vector &r0, ArrayView< const Size > idxs=nullptr)
Computes the total inertia tensor of particles with respect to given center.
Definition: Analysis.cpp:484
@ MOON
Body is on elliptical trajectory, it is a potential sattelite.
Array< HistPoint > getCumulativeHistogram(const Storage &storage, const ExtHistogramId id, const HistogramSource source, const HistogramParams &params)
Computes cumulative (integral) histogram of particles in the storage.
Definition: Analysis.cpp:814
@ VELOCITIES
Particle velocities.
@ ROTATIONAL_AXIS
Distribution of axis directions, from -pi to pi.
@ ROTATIONAL_FREQUENCY
Rotational frequency in revs/day.
Definition: MemoryPool.h:5
void printHelp()
Definition: main.cpp:751
void origComponents(const Path &lastDumpPath, const Path &firstDumpPath, const Path &colorizedDumpPath)
Definition: main.cpp:659
int pkdgravToOmega(const Path &filePath, const Path &omegaPath)
Definition: main.cpp:57
int ssfToSfd(const Post::HistogramSource source, const Path &filePath, const Path &sfdPath)
Definition: main.cpp:100
void ssfToVelDir(const Path &filePath, const Path &outPath)
Definition: main.cpp:232
void processHarrisFile()
Definition: main.cpp:473
int ssfToVelocity(const Path &filePath, const Path &outPath)
Definition: main.cpp:187
int pkdgravToSfd(const Path &filePath, const Path &sfdPath)
Definition: main.cpp:40
void extractLr(const Path &inputPath, const Path &outputPath)
Definition: main.cpp:683
int ssfToOmega(const Path &filePath, const Path &omegaPath, const Path &omegaDPath, const Path &omegaDirPath)
Definition: main.cpp:121
int pkdgravToMoons(const Path &filePath, const Float limit)
Definition: main.cpp:86
void makeSwift(const Path &filePath)
Definition: main.cpp:589
Vector values
Eigenvalues.
Optional< Size > number
Definition: main.cpp:327
Optional< std::string > name
Definition: main.cpp:328
Optional< Float > radius
Definition: main.cpp:287
std::string name
Definition: main.cpp:286
Optional< Size > number
Definition: main.cpp:285
Optional< Float > period
Definition: main.cpp:288
Point in 2D plot.
Definition: Point.h:16
Point in the histogram.
Definition: Analysis.h:277
Parameters of the histogram.
Definition: Analysis.h:227
Function< bool(Size index)> validator
Function used for inclusiong/exclusion of values in the histogram.
Definition: Analysis.h:273
Float massCutoff
Cutoff value (lower bound) of particle mass for inclusion in the histogram.
Definition: Analysis.h:246
Interval range
Range of values from which the histogram is constructed.
Definition: Analysis.h:232
Size binCnt
Number of histogram bins.
Definition: Analysis.h:237