PolyFEM
Loading...
Searching...
No Matches
OptState.cpp
Go to the documentation of this file.
1#include "OptState.hpp"
2
3#include <polyfem/Common.hpp>
4
9
15
16#include <polysolve/nonlinear/Solver.hpp>
17
18#include <spdlog/sinks/stdout_color_sinks.h>
19#include <spdlog/sinks/basic_file_sink.h>
20#include <spdlog/sinks/ostream_sink.h>
21
22#include <Eigen/Core>
23
24#include <memory>
25#include <string>
26#include <vector>
27
28namespace spdlog::level
29{
31 spdlog::level::level_enum,
32 {{spdlog::level::level_enum::trace, "trace"},
33 {spdlog::level::level_enum::debug, "debug"},
34 {spdlog::level::level_enum::info, "info"},
35 {spdlog::level::level_enum::warn, "warning"},
36 {spdlog::level::level_enum::err, "error"},
37 {spdlog::level::level_enum::critical, "critical"},
38 {spdlog::level::level_enum::off, "off"},
39 {spdlog::level::level_enum::trace, 0},
40 {spdlog::level::level_enum::debug, 1},
41 {spdlog::level::level_enum::info, 2},
42 {spdlog::level::level_enum::warn, 3},
43 {spdlog::level::level_enum::err, 3},
44 {spdlog::level::level_enum::critical, 4},
45 {spdlog::level::level_enum::off, 5}})
46}
47
48namespace polyfem
49{
53
58
60 const std::string &log_file,
61 const spdlog::level::level_enum log_level,
62 const spdlog::level::level_enum file_log_level,
63 const bool is_quiet)
64 {
65 std::vector<spdlog::sink_ptr> sinks;
66
67 if (!is_quiet)
68 {
69 console_sink_ = std::make_shared<spdlog::sinks::stdout_color_sink_mt>();
70 sinks.emplace_back(console_sink_);
71 }
72
73 if (!log_file.empty())
74 {
75 file_sink_ = std::make_shared<spdlog::sinks::basic_file_sink_mt>(log_file, /*truncate=*/true);
76 // Set the file sink separately from the console so it can save all messages
77 file_sink_->set_level(file_log_level);
78 sinks.push_back(file_sink_);
79 }
80
81 init_logger(sinks, log_level);
82 spdlog::flush_every(std::chrono::seconds(3));
83 }
84
85 void OptState::init_logger(std::ostream &os, const spdlog::level::level_enum log_level)
86 {
87 std::vector<spdlog::sink_ptr> sinks;
88 sinks.emplace_back(std::make_shared<spdlog::sinks::ostream_sink_mt>(os, false));
89 init_logger(sinks, log_level);
90 }
91
93 const std::vector<spdlog::sink_ptr> &sinks,
94 const spdlog::level::level_enum log_level)
95 {
96 set_adjoint_logger(std::make_shared<spdlog::logger>("adjoint-polyfem", sinks.begin(), sinks.end()));
97
98 // Set the logger at the lowest level, so all messages are passed to the sinks
99 adjoint_logger().set_level(spdlog::level::trace);
100 set_log_level(log_level);
101 }
102
103 void OptState::set_log_level(const spdlog::level::level_enum log_level)
104 {
105 adjoint_logger().set_level(log_level);
106 if (console_sink_)
107 console_sink_->set_level(log_level); // Shared by all loggers
108 }
109
110 void OptState::init(const json &p_args_in, const bool strict_validation)
111 {
112 json args_in = p_args_in; // mutable copy
113 args = solver::AdjointOptUtils::apply_opt_json_spec(args_in, strict_validation);
114
115 // Save output directory and resolve output paths dynamically
116 const std::string output_dir = utils::resolve_path(this->args["output"]["directory"], root_path(), false);
117 if (!output_dir.empty())
118 {
119 std::filesystem::create_directories(output_dir);
120 }
121 this->output_dir = output_dir;
122
123 std::string out_path_log = this->args["output"]["log"]["path"];
124 if (!out_path_log.empty())
125 {
126 out_path_log = utils::resolve_path(out_path_log, root_path(), false);
127 }
128
130 out_path_log,
131 this->args["output"]["log"]["level"],
132 this->args["output"]["log"]["file_level"],
133 this->args["output"]["log"]["quiet"]);
134
135 adjoint_logger().info("Saving adjoint output to {}", output_dir);
136
137 const int thread_in = this->args["solver"]["max_threads"];
139 }
140
141 void OptState::create_states(const int max_threads)
142 {
144 root_path(),
145 args["states"],
146 max_threads <= 0 ? std::numeric_limits<unsigned int>::max() : max_threads);
147
148 diff_caches.resize(states.size());
149 for (auto &diff_cache : diff_caches)
150 {
151 diff_cache = std::make_shared<DiffCache>();
152 }
153
155
157 }
158
160 {
161 for (int i = 0; i < states.size(); ++i)
162 {
163 const State &state = *states[i];
164
165 // No transient linear support.
166 if (state.problem->is_time_dependent() && state.is_problem_linear())
167 {
169 "State {}: transient linear problem is not supported in optimization.", i);
170 }
171
172 if (state.is_contact_enabled())
173 {
174 // No non-convergent contact formulation support.
175 if (!state.args["contact"]["use_gcp_formulation"].get<bool>()
176 && !state.args["contact"]["use_convergent_formulation"].get<bool>())
177 {
179 "State {}: non-convergent contact formulation is not supported in optimization.", i);
180 }
181
182 // No non-const barrier stiffness support.
183 if (state.args["/solver/contact/barrier_stiffness"_json_pointer].is_string())
184 {
186 "State {}: only constant barrier stiffness is supported in optimization.", i);
187 }
188 }
189
190 // No non-const boundary support.
191 if (state.args.contains("boundary_conditions") && state.args["boundary_conditions"].contains("rhs"))
192 {
193 const json &rhs = state.args["boundary_conditions"]["rhs"];
194 if (rhs.is_string() || (rhs.is_array() && rhs.size() > 0 && rhs[0].is_string()))
195 {
197 "State {}: only constant rhs over space is supported in optimization.", i);
198 }
199 }
200
201 // No high order geometric basis support.
202 for (const auto &element_bases : state.geom_bases())
203 {
204 for (const auto &basis : element_bases.bases)
205 {
206 if (basis.order() > 1)
207 {
209 "State {}: high-order geometry basis is not supported in optimization.", i);
210 }
211 }
212 }
213 }
214 }
215
217 {
218 /* DOFS */
219 ndof = 0;
220 for (const auto &arg : args["parameters"])
221 {
223 ndof += size;
224 variable_sizes.push_back(size);
225 }
226
227 /* variable to simulations */
229 args["variable_to_simulation"], states, diff_caches, variable_sizes);
230 }
231
233 {
234 /* forms */
235 std::shared_ptr<solver::AdjointForm> obj = from_json::build_form(
237
238 /* stopping conditions */
239 std::vector<std::shared_ptr<solver::AdjointForm>> stopping_conditions;
240 for (const auto &arg : args["stopping_conditions"])
241 stopping_conditions.push_back(
243
244 nl_problem = std::make_unique<solver::AdjointNLProblem>(
245 obj, stopping_conditions, variable_to_simulations, states, diff_caches, args);
246 }
247
254
255 double OptState::eval(Eigen::VectorXd &x) const
256 {
257 nl_problem->solution_changed(x);
258 return nl_problem->value(x);
259 }
260
261 void OptState::solve(Eigen::VectorXd &x)
262 {
264 args["solver"]["nonlinear"],
265 args["solver"]["linear"],
266 args["solver"]["advanced"]["characteristic_length"]);
267 nl_problem->normalize_forms();
268 nl_solver->minimize(*nl_problem, x);
269 }
270} // namespace polyfem
int x
std::vector< std::shared_ptr< State > > states
State used in the opt.
Definition OptState.hpp:103
void solve(Eigen::VectorXd &x)
Definition OptState.cpp:261
void initial_guess(Eigen::VectorXd &x)
Definition OptState.cpp:248
void set_log_level(const spdlog::level::level_enum log_level)
change log level
Definition OptState.cpp:103
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
std::vector< std::shared_ptr< DiffCache > > diff_caches
Definition OptState.hpp:104
spdlog::sink_ptr file_sink_
Definition OptState.hpp:96
std::string root_path() const
Definition OptState.hpp:81
void check_unsupported() const
Check and throw if any forward simulation State is not supported.
Definition OptState.cpp:159
double eval(Eigen::VectorXd &x) const
Definition OptState.cpp:255
std::vector< int > variable_sizes
variables
Definition OptState.hpp:107
std::unique_ptr< solver::AdjointNLProblem > nl_problem
Definition OptState.hpp:112
std::string output_dir
Directory for output files.
Definition OptState.hpp:116
spdlog::sink_ptr console_sink_
logger sink to stdout
Definition OptState.hpp:95
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:59
solver::VariableToSimulationGroup variable_to_simulations
Definition OptState.hpp:110
void init_variables()
init variables
Definition OptState.cpp:216
OptState()
Constructor.
Definition OptState.cpp:54
main class that contains the polyfem solver and all its state
Definition State.hpp:113
const std::vector< basis::ElementBases > & geom_bases() const
Get a constant reference to the geometry mapping bases.
Definition State.hpp:257
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:202
json args
main input arguments containing all defaults
Definition State.hpp:135
bool is_contact_enabled() const
does the simulation have contact
Definition State.hpp:699
bool is_problem_linear() const
Returns whether the system is linear. Collisions and pressure add nonlinearity to the problem.
Definition State.hpp:523
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::vector< std::shared_ptr< State > > build_states(const std::string &root_path, const json &args, const size_t max_threads)
solver::VariableToSimulationGroup build_variable_to_simulation_group(const json &args, const std::vector< std::shared_ptr< State > > &states, const std::vector< std::shared_ptr< DiffCache > > &diff_caches, const std::vector< int > &variable_sizes)
std::shared_ptr< solver::AdjointForm > build_form(const json &args, const solver::VariableToSimulationGroup &var2sim, const std::vector< std::shared_ptr< State > > &states, const std::vector< std::shared_ptr< DiffCache > > &diff_caches)
NLOHMANN_JSON_SERIALIZE_ENUM(CollisionProxyTessellation, {{CollisionProxyTessellation::REGULAR, "regular"}, {CollisionProxyTessellation::IRREGULAR, "irregular"}})
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:30
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:68
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:79
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)