PolyFEM
Loading...
Searching...
No Matches
JSONUtils.cpp
Go to the documentation of this file.
1#include "JSONUtils.hpp"
2
5
6#include <fstream>
7
8#include <Eigen/Geometry>
9
10namespace polyfem
11{
12 namespace utils
13 {
15 {
16 if (!args.contains("common"))
17 return;
18
19 const std::string common_params_path = resolve_path(args["common"], args["root_path"]);
20
21 if (common_params_path.empty())
22 return;
23
24 std::ifstream file(common_params_path);
25 if (!file.is_open())
26 log_and_throw_error("Unable to open common params {} file", common_params_path);
27
28 json common_params;
29 file >> common_params;
30 file.close();
31
32 // Recursively apply common params
33 const bool has_root_path = common_params.contains("root_path");
34 if (has_root_path)
35 common_params["root_path"] = resolve_path(common_params["root_path"], common_params_path);
36 else
37 common_params["root_path"] = common_params_path;
38 apply_common_params(common_params);
39
40 // If there is a root path in the common params, it overrides the one in the current params.
41 // This is somewhat backwards as normally current params override common params, but this is
42 // an easy way to make sure that the relative paths in common are correct.
43 if (has_root_path)
44 args["root_path"] = common_params["root_path"];
45
46 json patch;
47 if (args.contains("patch"))
48 {
49 patch = args["patch"];
50 args.erase("patch");
51 }
52
53 common_params.merge_patch(args);
54 if (!patch.empty())
55 common_params = common_params.patch(patch);
56 args = common_params;
57
58 args.erase("common"); // Remove common params from the final json
59 }
60
61 Eigen::Matrix3d to_rotation_matrix(const json &jr, std::string mode)
62 {
63 std::transform(mode.begin(), mode.end(), mode.begin(), ::tolower);
64
65 if (jr.is_array() && jr.empty())
66 {
67 return Eigen::Matrix3d::Identity(3, 3);
68 }
69
70 Eigen::VectorXd r;
71 if (jr.is_number())
72 {
73 r.setZero(3);
74 assert(mode.size() == 1); // must be either "x", "y", or "z"
75 int i = mode[0] - 'x';
76 assert(i >= 0 && i < 3);
77 r[i] = jr.get<double>();
78 }
79 else
80 {
81 assert(jr.is_array());
82 r = jr;
83 }
84
85 if (mode == "axis_angle")
86 {
87 assert(r.size() == 4);
88 double angle = deg2rad(r[0]); // NOTE: assumes input angle is in degrees
89 Eigen::Vector3d axis = r.tail<3>().normalized();
90 return Eigen::AngleAxisd(angle, axis).toRotationMatrix();
91 }
92
93 if (mode == "quaternion")
94 {
95 assert(r.size() == 4);
96 Eigen::Vector4d q = r.normalized();
97 return Eigen::Quaterniond(q).toRotationMatrix();
98 }
99
100 // The following expect the input is given in degrees
101 r = deg2rad(r);
102
103 if (mode == "rotation_vector")
104 {
105 assert(r.size() == 3);
106 double angle = r.norm();
107 if (angle != 0)
108 {
109 return Eigen::AngleAxisd(angle, r / angle).toRotationMatrix();
110 }
111 else
112 {
113 return Eigen::Matrix3d::Identity();
114 }
115 }
116
117 Eigen::Matrix3d R = Eigen::Matrix3d::Identity();
118
119 for (int i = 0; i < mode.size(); i++)
120 {
121 int j = mode[i] - 'x';
122 assert(j >= 0 && j < 3);
123 Eigen::Vector3d axis = Eigen::Vector3d::Zero();
124 axis[j] = 1;
125 R = Eigen::AngleAxisd(r[j], axis).toRotationMatrix() * R;
126 }
127
128 return R;
129 }
130
131 bool is_param_valid(const json &params, const std::string &key)
132 {
133 return params.contains(key) && !params[key].is_null();
134 }
135 } // namespace utils
136} // namespace polyfem
std::string resolve_path(const std::string &path, const std::string &input_file_path, const bool only_if_exists=false)
Eigen::Matrix3d to_rotation_matrix(const json &jr, std::string mode)
Definition JSONUtils.cpp:61
T deg2rad(T deg)
Definition JSONUtils.hpp:18
bool is_param_valid(const json &params, const std::string &key)
Determine if a key exists and is non-null in a json object.
void apply_common_params(json &args)
Definition JSONUtils.cpp:14
nlohmann::json json
Definition Common.hpp:9
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71