PolyFEM
Loading...
Searching...
No Matches
OptState.cpp
Go to the documentation of this file.
1#include "OptState.hpp"
2
7
10
11#include <polysolve/nonlinear/Solver.hpp>
12
13#include <spdlog/sinks/stdout_color_sinks.h>
14#include <spdlog/sinks/basic_file_sink.h>
15#include <spdlog/sinks/ostream_sink.h>
16
18{
19 NLOHMANN_JSON_SERIALIZE_ENUM(
20 spdlog::level::level_enum,
21 {{spdlog::level::level_enum::trace, "trace"},
22 {spdlog::level::level_enum::debug, "debug"},
23 {spdlog::level::level_enum::info, "info"},
24 {spdlog::level::level_enum::warn, "warning"},
25 {spdlog::level::level_enum::err, "error"},
26 {spdlog::level::level_enum::critical, "critical"},
27 {spdlog::level::level_enum::off, "off"},
28 {spdlog::level::level_enum::trace, 0},
29 {spdlog::level::level_enum::debug, 1},
30 {spdlog::level::level_enum::info, 2},
31 {spdlog::level::level_enum::warn, 3},
32 {spdlog::level::level_enum::err, 3},
33 {spdlog::level::level_enum::critical, 4},
34 {spdlog::level::level_enum::off, 5}})
35}
36
37namespace polyfem
38{
40 {
41
42 }
43
48
50 const std::string &log_file,
51 const spdlog::level::level_enum log_level,
52 const spdlog::level::level_enum file_log_level,
53 const bool is_quiet)
54 {
55 std::vector<spdlog::sink_ptr> sinks;
56
57 if (!is_quiet)
58 {
59 console_sink_ = std::make_shared<spdlog::sinks::stdout_color_sink_mt>();
60 sinks.emplace_back(console_sink_);
61 }
62
63 if (!log_file.empty())
64 {
65 file_sink_ = std::make_shared<spdlog::sinks::basic_file_sink_mt>(log_file, /*truncate=*/true);
66 // Set the file sink separately from the console so it can save all messages
67 file_sink_->set_level(file_log_level);
68 sinks.push_back(file_sink_);
69 }
70
71 init_logger(sinks, log_level);
72 spdlog::flush_every(std::chrono::seconds(3));
73 }
74
75 void OptState::init_logger(std::ostream &os, const spdlog::level::level_enum log_level)
76 {
77 std::vector<spdlog::sink_ptr> sinks;
78 sinks.emplace_back(std::make_shared<spdlog::sinks::ostream_sink_mt>(os, false));
79 init_logger(sinks, log_level);
80 }
81
83 const std::vector<spdlog::sink_ptr> &sinks,
84 const spdlog::level::level_enum log_level)
85 {
86 set_adjoint_logger(std::make_shared<spdlog::logger>("adjoint-polyfem", sinks.begin(), sinks.end()));
87
88 // Set the logger at the lowest level, so all messages are passed to the sinks
89 adjoint_logger().set_level(spdlog::level::trace);
90 set_log_level(log_level);
91 }
92
93 void OptState::set_log_level(const spdlog::level::level_enum log_level)
94 {
95 adjoint_logger().set_level(log_level);
96 if (console_sink_)
97 console_sink_->set_level(log_level); // Shared by all loggers
98 }
99
100 void OptState::init(const json &p_args_in, const bool strict_validation)
101 {
102 json args_in = p_args_in; // mutable copy
103 args = solver::AdjointOptUtils::apply_opt_json_spec(args_in, strict_validation);
104
105 // Save output directory and resolve output paths dynamically
106 const std::string output_dir = utils::resolve_path(this->args["output"]["directory"], root_path(), false);
107 if (!output_dir.empty())
108 {
109 std::filesystem::create_directories(output_dir);
110 }
111 this->output_dir = output_dir;
112
113 std::string out_path_log = this->args["output"]["log"]["path"];
114 if (!out_path_log.empty())
115 {
116 out_path_log = utils::resolve_path(out_path_log, root_path(), false);
117 }
118
120 out_path_log,
121 this->args["output"]["log"]["level"],
122 this->args["output"]["log"]["file_level"],
123 this->args["output"]["log"]["quiet"]);
124
125 adjoint_logger().info("Saving adjoint output to {}", output_dir);
126
127 const int thread_in = this->args["solver"]["max_threads"];
129 }
130
131 void OptState::create_states(const polyfem::solver::CacheLevel level, const int max_threads)
132 {
134 args["states"],
135 level,
136 max_threads <= 0 ? std::numeric_limits<unsigned int>::max() : max_threads);
137
139 }
140
142 {
143 /* DOFS */
144 ndof = 0;
145 for (const auto &arg : args["parameters"])
146 {
148 ndof += size;
149 variable_sizes.push_back(size);
150 }
151
152 /* variable to simulations */
153 variable_to_simulations.init(args["variable_to_simulation"], states, variable_sizes);
154 }
155
157 {
158 /* forms */
159 std::shared_ptr<solver::AdjointForm> obj = solver::AdjointOptUtils::create_form(
160 args["functionals"], variable_to_simulations, states);
161
162 /* stopping conditions */
163 std::vector<std::shared_ptr<solver::AdjointForm>> stopping_conditions;
164 for (const auto &arg : args["stopping_conditions"])
165 stopping_conditions.push_back(
167
168 nl_problem = std::make_unique<solver::AdjointNLProblem>(
169 obj, stopping_conditions, variable_to_simulations, states, args);
170 }
171
178
179 double OptState::eval(Eigen::VectorXd &x) const
180 {
181 nl_problem->solution_changed(x);
182 return nl_problem->value(x);
183 }
184
185 void OptState::solve(Eigen::VectorXd &x)
186 {
188 args["solver"]["nonlinear"],
189 args["solver"]["linear"],
190 args["solver"]["advanced"]["characteristic_length"]);
191 nl_solver->minimize(*nl_problem, x);
192 }
193} // namespace polyfem
int x
std::vector< std::shared_ptr< State > > states
State used in the opt.
Definition OptState.hpp:94
void solve(Eigen::VectorXd &x)
Definition OptState.cpp:185
void initial_guess(Eigen::VectorXd &x)
Definition OptState.cpp:172
void set_log_level(const spdlog::level::level_enum log_level)
change log level
Definition OptState.cpp:93
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
spdlog::sink_ptr file_sink_
Definition OptState.hpp:87
std::string root_path() const
Definition OptState.hpp:75
double eval(Eigen::VectorXd &x) const
Definition OptState.cpp:179
std::vector< int > variable_sizes
variables
Definition OptState.hpp:97
std::unique_ptr< solver::AdjointNLProblem > nl_problem
Definition OptState.hpp:102
std::string output_dir
Directory for output files.
Definition OptState.hpp:106
spdlog::sink_ptr console_sink_
logger sink to stdout
Definition OptState.hpp:86
void init_logger(const std::string &log_file, const spdlog::level::level_enum log_level, const spdlog::level::level_enum file_log_level, const bool is_quiet)
initializing the logger
Definition OptState.cpp:49
void create_states(const polyfem::solver::CacheLevel level, const int max_threads=-1)
create the opt states
Definition OptState.cpp:131
solver::VariableToSimulationGroup variable_to_simulations
Definition OptState.hpp:100
void init_variables()
init variables
Definition OptState.cpp:141
OptState()
Constructor.
Definition OptState.cpp:44
void init(const json &args, const std::vector< std::shared_ptr< State > > &states, const std::vector< int > &variable_sizes)
void update(const Eigen::VectorXd &x)
Update parameters in simulators.
void set_logger(spdlog::logger &logger)
static GeogramUtils & instance()
static NThread & get()
Definition par_for.hpp:19
void set_num_threads(const int max_threads)
Definition par_for.hpp:27
std::string resolve_path(const std::string &path, const std::string &input_file_path, const bool only_if_exists=false)
spdlog::logger & adjoint_logger()
Retrieves the current logger for adjoint.
Definition Logger.cpp:28
nlohmann::json json
Definition Common.hpp:9
void set_adjoint_logger(std::shared_ptr< spdlog::logger > p_logger)
Setup a logger object to be used by adjoint Polyfem.
Definition Logger.cpp:66
static std::shared_ptr< AdjointForm > create_form(const json &args, const VariableToSimulationGroup &var2sim, const std::vector< std::shared_ptr< State > > &states)
static std::shared_ptr< polysolve::nonlinear::Solver > make_nl_solver(const json &solver_params, const json &linear_solver_params, const double characteristic_length)
static json apply_opt_json_spec(const json &input_args, bool strict_validation)
static Eigen::VectorXd inverse_evaluation(const json &args, const int ndof, const std::vector< int > &variable_sizes, VariableToSimulationGroup &var2sim)
static int compute_variable_size(const json &args, const std::vector< std::shared_ptr< State > > &states)
static std::vector< std::shared_ptr< State > > create_states(const json &state_args, const CacheLevel &level, const size_t max_threads)