MPS-Basic
Loading...
Searching...
No Matches
loader.cpp
Go to the documentation of this file.
1#include "loader.hpp"
2
5
6#include <cmath>
7#include <iostream>
8#include <yaml-cpp/yaml.h>
9
10using std::cerr;
11using std::cout;
12using std::endl;
13
14namespace fs = std::filesystem;
15
16Input Loader::load(const fs::path& settingPath, const fs::path& outputDirectory) {
17 Input input;
18
19 input.settings = loadSettingYaml(settingPath);
20
21 auto particlesPath = input.settings.particlesPath;
22 this->particlesLoader = getParticlesLoader(particlesPath);
23 auto [startTime, particles] = this->particlesLoader->load(particlesPath, input.settings.defaultDensity);
24 input.startTime = startTime;
25 input.particles = particles;
26
27 copyInputFileToOutputDirectory(settingPath, outputDirectory);
28 copyInputFileToOutputDirectory(particlesPath, outputDirectory);
29
30 return input;
31}
32
39std::unique_ptr<ParticlesLoader::Interface> Loader::getParticlesLoader(const fs::path& particlesPath) {
40 auto extension = particlesPath.extension();
41 if (extension == ".csv") {
42 return std::make_unique<ParticlesLoader::Csv>();
43 } else if (extension == ".prof") {
44 return std::make_unique<ParticlesLoader::Prof>();
45 } else {
46 cerr << "unsupported file format: " << particlesPath.extension() << endl;
47 std::exit(-1);
48 }
49}
50
57void Loader::copyInputFileToOutputDirectory(const fs::path& inputFilePath, const fs::path& outputDirectory) {
58 auto outputFilePath = outputDirectory / inputFilePath.filename();
59 if (fs::exists(outputFilePath)) {
60 cerr << "file " << outputFilePath << " already exists in the output directory" << endl;
61 std::exit(-1);
62 } else {
63 fs::copy_file(inputFilePath, outputFilePath);
64 }
65}
66
67Settings Loader::loadSettingYaml(const fs::path& settingPath) {
68 YAML::Node yaml = YAML::LoadFile(settingPath.string());
69
70 Settings s;
71
72 // computational conditions
73 s.dim = yaml["dim"].as<int>();
74 s.particleDistance = yaml["particleDistance"].as<double>();
75 s.dt = yaml["dt"].as<double>();
76 s.endTime = yaml["endTime"].as<double>();
77 s.outputPeriod = yaml["outputPeriod"].as<double>();
78 s.cflCondition = yaml["cflCondition"].as<double>();
79 s.numPhysicalCores = yaml["numPhysicalCores"].as<int>();
80
81 // physical properties
82 s.defaultDensity = yaml["defaultDensity"].as<double>();
83 s.kinematicViscosity = yaml["kinematicViscosity"].as<double>();
84
85 // gravity
86 s.gravity[0] = yaml["gravity"][0].as<double>();
87 s.gravity[1] = yaml["gravity"][1].as<double>();
88 s.gravity[2] = yaml["gravity"][2].as<double>();
89
90 // free surface detection
91 s.surfaceDetection_numberDensity_threshold = yaml["surfaceDetection-numberDensity-threshold"].as<double>();
92 s.surfaceDetection_particleDistribution = yaml["surfaceDetection-particleDistribution"].as<bool>();
94 yaml["surfaceDetection-particleDistribution-threshold"].as<double>();
95
96 // pressure calculation method
97 s.pressureCalculationMethod = yaml["pressureCalculationMethod"].as<std::string>();
98 // for Implicit
99 s.compressibility = yaml["compressibility"].as<double>();
100 s.relaxationCoefficientForPressure = yaml["relaxationCoefficientForPressure"].as<double>();
101 // for Explicit
102 s.soundSpeed = yaml["soundSpeed"].as<double>();
103
104 // collision
105 s.collisionDistance = yaml["collisionDistanceRatio"].as<double>() * s.particleDistance;
106 s.coefficientOfRestitution = yaml["coefficientOfRestitution"].as<double>();
107
108 // effective radius
109 s.re_forNumberDensity = yaml["radiusRatioForNumberDensity"].as<double>() * s.particleDistance;
110 s.re_forGradient = yaml["radiusRatioForGradient"].as<double>() * s.particleDistance;
111 s.re_forLaplacian = yaml["radiusRatioForLaplacian"].as<double>() * s.particleDistance;
113
114 // domain
115 s.domain.xMin = yaml["domainMin"][0].as<double>();
116 s.domain.xMax = yaml["domainMax"][0].as<double>();
117 s.domain.yMin = yaml["domainMin"][1].as<double>();
118 s.domain.yMax = yaml["domainMax"][1].as<double>();
119 s.domain.zMin = yaml["domainMin"][2].as<double>();
120 s.domain.zMax = yaml["domainMax"][2].as<double>();
124
125 // particlesPath
126 auto yamlDir = settingPath.parent_path();
127 auto relativeProfPath = yaml["particlesPath"].as<std::string>();
128 s.particlesPath = fs::weakly_canonical(yamlDir / relativeProfPath);
129
130 // outputVtkFormatInBinary
131 // check if outputVtkFormat is defined in the yaml file since it is optional
132 if (yaml["outputVtkInBinary"]) {
133 s.outputVtkInBinary = yaml["outputVtkInBinary"].as<bool>();
134 } else {
135 s.outputVtkInBinary = false;
136 }
137 return s;
138}
double xMax
maximum x coordinate of the domain
Definition domain.hpp:12
double zMin
minimum z coordinate of the domain
Definition domain.hpp:15
double zMax
maximum z coordinate of the domain
Definition domain.hpp:16
double xMin
minimum x coordinate of the domain
Definition domain.hpp:11
double yMin
minimum y coordinate of the domain
Definition domain.hpp:13
double zLength
Definition domain.hpp:17
double yLength
Definition domain.hpp:17
double yMax
maximum y coordinate of the domain
Definition domain.hpp:14
double xLength
Definition domain.hpp:17
std::unique_ptr< ParticlesLoader::Interface > particlesLoader
Definition loader.hpp:37
std::unique_ptr< ParticlesLoader::Interface > getParticlesLoader(const fs::path &particlesPath)
Get the Particles Loader object according to the input file extension.
Definition loader.cpp:39
void copyInputFileToOutputDirectory(const fs::path &inputFilePath, const fs::path &outputDirectory)
Copy the input file to the output directory.
Definition loader.cpp:57
Input load(const fs::path &settingPath, const fs::path &outputDirectory)
Load the setting file and the particle file.
Definition loader.cpp:16
Settings loadSettingYaml(const fs::path &settingPath)
Definition loader.cpp:67
Represents the input data for MPS simulation.
Definition input.hpp:12
Particles particles
Initial particles arrangement in the simulation.
Definition input.hpp:14
Settings settings
Settings for the simulation.
Definition input.hpp:13
double startTime
Start time of the simulation.
Definition input.hpp:15
Struct for settings of calculation.
Definition settings.hpp:15
Eigen::Vector3d gravity
Gravity.
Definition settings.hpp:32
double reMax
Maximum of effective radius.
Definition settings.hpp:57
double outputPeriod
Output period of the simulation.
Definition settings.hpp:21
double soundSpeed
Speed of sound for Explicit method.
Definition settings.hpp:47
double compressibility
Compressibility of the fluid for Implicit method.
Definition settings.hpp:44
double collisionDistance
Distance for collision detection.
Definition settings.hpp:50
int numPhysicalCores
Number of cores to calculate.
Definition settings.hpp:23
double particleDistance
Initial distance between particles.
Definition settings.hpp:18
std::filesystem::path particlesPath
Path for input particle file.
Definition settings.hpp:60
bool outputVtkInBinary
Flag for saving VTK file in binary format.
Definition settings.hpp:61
double endTime
End time of the simulation.
Definition settings.hpp:20
int dim
Dimension of the simulation.
Definition settings.hpp:17
double relaxationCoefficientForPressure
Relaxation coefficient for pressure for Implicit method.
Definition settings.hpp:45
double re_forNumberDensity
Effective radius for number density.
Definition settings.hpp:54
bool surfaceDetection_particleDistribution
flag for free surface detection based on particle distribution
Definition settings.hpp:38
double defaultDensity
default density for fluid and wall particles.
Definition settings.hpp:29
std::string pressureCalculationMethod
Method for pressure calculation.
Definition settings.hpp:42
double surfaceDetection_particleDistribution_threshold
Definition settings.hpp:39
double kinematicViscosity
Kinematic viscosity.
Definition settings.hpp:28
double cflCondition
CFL condition.
Definition settings.hpp:22
double coefficientOfRestitution
Coefficient of restitution.
Definition settings.hpp:51
double surfaceDetection_numberDensity_threshold
threshold ratio of number density for free surface detection
Definition settings.hpp:36
double re_forGradient
Effective radius for gradient.
Definition settings.hpp:55
double dt
Time step.
Definition settings.hpp:19
Domain domain
domain of the simulation
Definition settings.hpp:25
double re_forLaplacian
Effective radius for Laplacian.
Definition settings.hpp:56