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