PolyFEM
Loading...
Searching...
No Matches
main.cpp
Go to the documentation of this file.
1#include <filesystem>
2
3#include <CLI/CLI.hpp>
4
5#include <h5pp/h5pp.h>
6
7#include <polyfem/State.hpp>
9
13
14using namespace polyfem;
15using namespace solver;
16
17bool has_arg(const CLI::App &command_line, const std::string &value)
18{
19 const auto *opt = command_line.get_option_no_throw(value.size() == 1 ? ("-" + value) : ("--" + value));
20 if (!opt)
21 return false;
22
23 return opt->count() > 0;
24}
25
26bool load_json(const std::string &json_file, json &out)
27{
28 std::ifstream file(json_file);
29
30 if (!file.is_open())
31 return false;
32
33 file >> out;
34
35 if (!out.contains("root_path"))
36 out["root_path"] = json_file;
37
38 return true;
39}
40
41bool load_yaml(const std::string &yaml_file, json &out)
42{
43 try
44 {
45 out = io::yaml_file_to_json(yaml_file);
46 if (!out.contains("root_path"))
47 out["root_path"] = yaml_file;
48 }
49 catch (...)
50 {
51 return false;
52 }
53 return true;
54}
55
56int forward_simulation(const CLI::App &command_line,
57 const std::string &hdf5_file,
58 const std::string output_dir,
59 const unsigned max_threads,
60 const bool is_strict,
61 const bool fallback_solver,
62 const spdlog::level::level_enum &log_level,
63 json &in_args);
64
65int optimization_simulation(const CLI::App &command_line,
66 const unsigned max_threads,
67 const bool is_strict,
68 const spdlog::level::level_enum &log_level,
69 json &opt_args);
70
71int main(int argc, char **argv)
72{
73 using namespace polyfem;
74
75 CLI::App command_line{"polyfem"};
76
77 command_line.ignore_case();
78 command_line.ignore_underscore();
79
80 // Eigen::setNbThreads(1);
81 unsigned max_threads = std::numeric_limits<unsigned>::max();
82 command_line.add_option("--max_threads", max_threads, "Maximum number of threads");
83
84 auto input = command_line.add_option_group("input");
85
86 std::string json_file = "";
87 input->add_option("-j,--json", json_file, "Simulation JSON file")->check(CLI::ExistingFile);
88
89 std::string yaml_file = "";
90 input->add_option("-y,--yaml", yaml_file, "Simulation YAML file")->check(CLI::ExistingFile);
91
92 std::string hdf5_file = "";
93 input->add_option("--hdf5", hdf5_file, "Simulation HDF5 file")->check(CLI::ExistingFile);
94
95 input->require_option(1);
96
97 std::string output_dir = "";
98 command_line.add_option("-o,--output_dir", output_dir, "Directory for output files")->check(CLI::ExistingDirectory | CLI::NonexistentPath);
99
100 bool is_strict = true;
101 command_line.add_flag("-s,--strict_validation,!--ns,!--no_strict_validation", is_strict, "Disables strict validation of input JSON");
102
103 bool fallback_solver = false;
104 command_line.add_flag("--enable_overwrite_solver", fallback_solver, "If solver in input is not present, falls back to default.");
105
106 const std::vector<std::pair<std::string, spdlog::level::level_enum>>
107 SPDLOG_LEVEL_NAMES_TO_LEVELS = {
108 {"trace", spdlog::level::trace},
109 {"debug", spdlog::level::debug},
110 {"info", spdlog::level::info},
111 {"warning", spdlog::level::warn},
112 {"error", spdlog::level::err},
113 {"critical", spdlog::level::critical},
114 {"off", spdlog::level::off}};
115 spdlog::level::level_enum log_level = spdlog::level::debug;
116 command_line.add_option("--log_level", log_level, "Log level")
117 ->transform(CLI::CheckedTransformer(SPDLOG_LEVEL_NAMES_TO_LEVELS, CLI::ignore_case));
118
119 CLI11_PARSE(command_line, argc, argv);
120
121 json in_args = json({});
122
123 if (!json_file.empty() || !yaml_file.empty())
124 {
125 const bool ok = !json_file.empty() ? load_json(json_file, in_args) : load_yaml(yaml_file, in_args);
126
127 if (!ok)
128 log_and_throw_error(fmt::format("unable to open {} file", json_file));
129
130 if (in_args.contains("states"))
131 return optimization_simulation(command_line, max_threads, is_strict, log_level, in_args);
132 else
133 return forward_simulation(command_line, "", output_dir, max_threads,
134 is_strict, fallback_solver, log_level, in_args);
135 }
136 else
137 return forward_simulation(command_line, hdf5_file, output_dir, max_threads,
138 is_strict, fallback_solver, log_level, in_args);
139}
140
141int forward_simulation(const CLI::App &command_line,
142 const std::string &hdf5_file,
143 const std::string output_dir,
144 const unsigned max_threads,
145 const bool is_strict,
146 const bool fallback_solver,
147 const spdlog::level::level_enum &log_level,
148 json &in_args)
149{
150 std::vector<std::string> names;
151 std::vector<Eigen::MatrixXi> cells;
152 std::vector<Eigen::MatrixXd> vertices;
153
154 if (in_args.empty() && hdf5_file.empty())
155 {
156 logger().error("No input file specified!");
157 return command_line.exit(CLI::RequiredError("--json or --hdf5"));
158 }
159
160 if (in_args.empty() && !hdf5_file.empty())
161 {
162 using MatrixXl = Eigen::Matrix<int64_t, Eigen::Dynamic, Eigen::Dynamic>;
163
164 h5pp::File file(hdf5_file, h5pp::FileAccess::READONLY);
165 std::string json_string = file.readDataset<std::string>("json");
166
167 in_args = json::parse(json_string);
168 in_args["root_path"] = hdf5_file;
169
170 names = file.findGroups("", "/meshes");
171 cells.resize(names.size());
172 vertices.resize(names.size());
173
174 for (int i = 0; i < names.size(); ++i)
175 {
176 const std::string &name = names[i];
177 cells[i] = file.readDataset<MatrixXl>("/meshes/" + name + "/c").cast<int>();
178 vertices[i] = file.readDataset<Eigen::MatrixXd>("/meshes/" + name + "/v");
179 }
180 }
181
182 json tmp = json::object();
183 if (has_arg(command_line, "log_level"))
184 tmp["/output/log/level"_json_pointer] = int(log_level);
185 if (has_arg(command_line, "max_threads"))
186 tmp["/solver/max_threads"_json_pointer] = max_threads;
187 if (has_arg(command_line, "output_dir"))
188 tmp["/output/directory"_json_pointer] = std::filesystem::absolute(output_dir);
189 if (has_arg(command_line, "enable_overwrite_solver"))
190 tmp["/solver/linear/enable_overwrite_solver"_json_pointer] = fallback_solver;
191 assert(tmp.is_object());
192 in_args.merge_patch(tmp);
193
194 State state;
195 state.init(in_args, is_strict);
196 state.load_mesh(/*non_conforming=*/false, names, cells, vertices);
197
198 // Mesh was not loaded successfully; load_mesh() logged the error.
199 if (state.mesh == nullptr)
200 {
201 // Cannot proceed without a mesh.
202 return EXIT_FAILURE;
203 }
204
205 state.stats.compute_mesh_stats(*state.mesh);
206
207 state.build_basis();
208
209 state.assemble_rhs();
210 state.assemble_mass_mat();
211
212 Eigen::MatrixXd sol;
213 Eigen::MatrixXd pressure;
214
215 state.solve_problem(sol, pressure);
216
217 state.compute_errors(sol);
218
219 logger().info("total time: {}s", state.timings.total_time());
220
221 state.save_json(sol);
222 state.export_data(sol, pressure);
223
224 return EXIT_SUCCESS;
225}
226
227int optimization_simulation(const CLI::App &command_line,
228 const unsigned max_threads,
229 const bool is_strict,
230 const spdlog::level::level_enum &log_level,
231 json &opt_args)
232{
233 json tmp = json::object();
234 if (has_arg(command_line, "log_level"))
235 tmp["/output/log/level"_json_pointer] = int(log_level);
236 if (has_arg(command_line, "max_threads"))
237 tmp["/solver/max_threads"_json_pointer] = max_threads;
238 opt_args.merge_patch(tmp);
239
240 OptState opt_state;
241 opt_state.init(opt_args, is_strict);
242
243 opt_state.create_states(opt_state.args["compute_objective"].get<bool>() ? polyfem::solver::CacheLevel::Solution : polyfem::solver::CacheLevel::Derivatives, opt_state.args["solver"]["max_threads"].get<int>());
244 opt_state.init_variables();
245 opt_state.create_problem();
246
247 Eigen::VectorXd x;
248 opt_state.initial_guess(x);
249
250 if (opt_state.args["compute_objective"].get<bool>())
251 {
252 logger().info("Objective is {}", opt_state.eval(x));
253 return EXIT_SUCCESS;
254 }
255
256 opt_state.solve(x);
257 return EXIT_SUCCESS;
258}
int x
main class that contains the polyfem adjoint solver and all its state
Definition OptState.hpp:22
void solve(Eigen::VectorXd &x)
Definition OptState.cpp:185
void initial_guess(Eigen::VectorXd &x)
Definition OptState.cpp:172
json args
main input arguments containing all defaults
Definition OptState.hpp:38
void init(const json &args, const bool strict_validation)
initialize the polyfem solver with a json settings
Definition OptState.cpp:100
double eval(Eigen::VectorXd &x) const
Definition OptState.cpp:179
void create_states(const polyfem::solver::CacheLevel level, const int max_threads=-1)
create the opt states
Definition OptState.cpp:131
void init_variables()
init variables
Definition OptState.cpp:141
main class that contains the polyfem solver and all its state
Definition State.hpp:79
io::OutRuntimeData timings
runtime statistics
Definition State.hpp:603
void init(const json &args, const bool strict_validation)
initialize the polyfem solver with a json settings
void assemble_mass_mat()
assemble mass, step 4 of solve build mass matrix based on defined basis modifies mass (and maybe more...
Definition State.cpp:1530
void compute_errors(const Eigen::MatrixXd &sol)
computes all errors
void load_mesh(bool non_conforming=false, const std::vector< std::string > &names=std::vector< std::string >(), const std::vector< Eigen::MatrixXi > &cells=std::vector< Eigen::MatrixXi >(), const std::vector< Eigen::MatrixXd > &vertices=std::vector< Eigen::MatrixXd >())
loads the mesh from the json arguments
Definition StateLoad.cpp:91
void save_json(const Eigen::MatrixXd &sol, std::ostream &out)
saves the output statistic to a stream
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:466
void solve_problem(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves the problems
Definition State.cpp:1726
void build_basis()
builds the bases step 2 of solve modifies bases, pressure_bases, geom_bases_, boundary_nodes,...
Definition State.cpp:577
io::OutStatsData stats
Other statistics.
Definition State.hpp:605
void assemble_rhs()
compute rhs, step 3 of solve build rhs vector based on defined basis and given rhs of the problem mod...
Definition State.cpp:1654
void export_data(const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure)
saves all data on the disk according to the input params
double total_time()
computes total time
Definition OutData.hpp:383
void compute_mesh_stats(const polyfem::mesh::Mesh &mesh)
compute stats (counts els type, mesh lenght, etc), step 1 of solve
Definition OutData.cpp:2705
bool load_json(const std::string &json_file, json &out)
Definition main.cpp:26
int main(int argc, char **argv)
Definition main.cpp:71
int optimization_simulation(const CLI::App &command_line, const unsigned max_threads, const bool is_strict, const spdlog::level::level_enum &log_level, json &opt_args)
Definition main.cpp:227
bool load_yaml(const std::string &yaml_file, json &out)
Definition main.cpp:41
bool has_arg(const CLI::App &command_line, const std::string &value)
Definition main.cpp:17
int forward_simulation(const CLI::App &command_line, const std::string &hdf5_file, const std::string output_dir, const unsigned max_threads, const bool is_strict, const bool fallback_solver, const spdlog::level::level_enum &log_level, json &in_args)
Definition main.cpp:141
json yaml_file_to_json(const std::string &yaml_filepath)
Load a YAML file to JSON.
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
nlohmann::json json
Definition Common.hpp:9
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71