8#ifdef POLYFEM_WITH_OPTIMIZATION
17using namespace solver;
19bool has_arg(
const CLI::App &command_line,
const std::string &value)
21 const auto *opt = command_line.get_option_no_throw(value.size() == 1 ? (
"-" + value) : (
"--" + value));
25 return opt->count() > 0;
30 std::ifstream file(json_file);
37 if (!out.contains(
"root_path"))
38 out[
"root_path"] = json_file;
48 if (!out.contains(
"root_path"))
49 out[
"root_path"] = yaml_file;
59 const std::string &hdf5_file,
60 const std::string output_dir,
61 const unsigned max_threads,
63 const bool fallback_solver,
64 const spdlog::level::level_enum &log_level,
67#ifdef POLYFEM_WITH_OPTIMIZATION
68int optimization_simulation(
const CLI::App &command_line,
69 const unsigned max_threads,
71 const spdlog::level::level_enum &log_level,
75int main(
int argc,
char **argv)
79 CLI::App command_line{
"polyfem"};
81 command_line.ignore_case();
82 command_line.ignore_underscore();
85 unsigned max_threads = std::numeric_limits<unsigned>::max();
86 command_line.add_option(
"--max_threads", max_threads,
"Maximum number of threads");
88 auto input = command_line.add_option_group(
"input");
90 std::string json_file =
"";
91 input->add_option(
"-j,--json", json_file,
"Simulation JSON file")->check(CLI::ExistingFile);
93 std::string yaml_file =
"";
94 input->add_option(
"-y,--yaml", yaml_file,
"Simulation YAML file")->check(CLI::ExistingFile);
96 std::string hdf5_file =
"";
97 input->add_option(
"--hdf5", hdf5_file,
"Simulation HDF5 file")->check(CLI::ExistingFile);
99 input->require_option(1);
101 std::string output_dir =
"";
102 command_line.add_option(
"-o,--output_dir", output_dir,
"Directory for output files")->check(CLI::ExistingDirectory | CLI::NonexistentPath);
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");
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.");
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));
123 CLI11_PARSE(command_line, argc, argv);
127 if (!json_file.empty() || !yaml_file.empty())
129 const bool ok = !json_file.empty() ?
load_json(json_file, in_args) :
load_yaml(yaml_file, in_args);
134 if (in_args.contains(
"states"))
136#ifndef POLYFEM_WITH_OPTIMIZATION
139 return optimization_simulation(command_line, max_threads, is_strict, log_level, in_args);
144 is_strict, fallback_solver, log_level, in_args);
148 is_strict, fallback_solver, log_level, in_args);
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,
160 std::vector<std::string> names;
161 std::vector<Eigen::MatrixXi> cells;
162 std::vector<Eigen::MatrixXd> vertices;
164 if (in_args.empty() && hdf5_file.empty())
166 logger().error(
"No input file specified!");
167 return command_line.exit(CLI::RequiredError(
"--json or --hdf5"));
170 if (in_args.empty() && !hdf5_file.empty())
172 using MatrixXl = Eigen::Matrix<int64_t, Eigen::Dynamic, Eigen::Dynamic>;
174 h5pp::File file(hdf5_file, h5pp::FileAccess::READONLY);
175 std::string json_string = file.readDataset<std::string>(
"json");
177 in_args = json::parse(json_string);
178 in_args[
"root_path"] = hdf5_file;
180 names = file.findGroups(
"",
"/meshes");
181 cells.resize(names.size());
182 vertices.resize(names.size());
184 for (
int i = 0; i < names.size(); ++i)
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");
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);
205 state.
init(in_args, is_strict);
206 state.
load_mesh(
false, names, cells, vertices);
209 if (state.
mesh ==
nullptr)
223 Eigen::MatrixXd pressure;
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,
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);
252 opt_state.
init(opt_args, is_strict);
261 if (opt_state.
args[
"compute_objective"].get<
bool>())
263 logger().info(
"Objective is {}", opt_state.
eval(
x));
main class that contains the polyfem adjoint solver and all its state
void solve(Eigen::VectorXd &x)
void initial_guess(Eigen::VectorXd &x)
void create_states(const int max_threads=-1)
create the opt states
json args
main input arguments containing all defaults
void init(const json &args, const bool strict_validation)
initialize the polyfem solver with a json settings
double eval(Eigen::VectorXd &x) const
void init_variables()
init variables
main class that contains the polyfem solver and all its state
io::OutRuntimeData timings
runtime statistics
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...
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
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
void build_basis()
builds the bases step 2 of solve modifies bases, pressure_bases, geom_bases_, boundary_nodes,...
io::OutStatsData stats
Other statistics.
void assemble_rhs()
compute rhs, step 3 of solve build rhs vector based on defined basis and given rhs of the problem mod...
void solve_problem(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step={}, const InitialConditionOverride *ic_override=nullptr)
solves the problems
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
void compute_mesh_stats(const polyfem::mesh::Mesh &mesh)
compute stats (counts els type, mesh lenght, etc), step 1 of solve
bool load_json(const std::string &json_file, json &out)
int main(int argc, char **argv)
bool load_yaml(const std::string &yaml_file, json &out)
bool has_arg(const CLI::App &command_line, const std::string &value)
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)
json yaml_file_to_json(const std::string &yaml_filepath)
Load a YAML file to JSON.
spdlog::logger & logger()
Retrieves the current logger.
void log_and_throw_error(const std::string &msg)