PolyFEM
Loading...
Searching...
No Matches
StateOutput.cpp
Go to the documentation of this file.
1#include <polyfem/State.hpp>
2
5
6#include <filesystem>
7
8namespace polyfem
9{
10 void State::compute_errors(const Eigen::MatrixXd &sol)
11 {
12 if (!args["output"]["advanced"]["compute_error"])
13 return;
14
15 double tend = 0;
16
17 if (!args["time"].is_null())
18 {
19 tend = args["time"]["tend"];
20 }
21
23 }
24
25 std::string State::root_path() const
26 {
27 if (utils::is_param_valid(args, "root_path"))
28 return args["root_path"].get<std::string>();
29 return "";
30 }
31
32 std::string State::resolve_input_path(const std::string &path, const bool only_if_exists) const
33 {
34 return utils::resolve_path(path, root_path(), only_if_exists);
35 }
36
37 std::string State::resolve_output_path(const std::string &path) const
38 {
39 if (output_dir.empty() || path.empty() || std::filesystem::path(path).is_absolute())
40 {
41 return path;
42 }
43 return std::filesystem::weakly_canonical(std::filesystem::path(output_dir) / path).string();
44 }
45
46 void State::save_timestep(const double time, const int t, const double t0, const double dt, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure)
47 {
48 if (args["output"]["advanced"]["save_time_sequence"] && !(t % args["output"]["paraview"]["skip_frame"].get<int>()))
49 {
50 logger().trace("Saving VTU...");
51 POLYFEM_SCOPED_TIMER("Saving VTU");
52 const std::string step_name = args["output"]["advanced"]["timestep_prefix"];
53
55 solution_frames.emplace_back();
56
58 resolve_output_path(fmt::format(step_name + "{:d}.vtu", t)),
59 *this, sol, pressure, time, dt,
62
64 resolve_output_path(args["output"]["paraview"]["file_name"]),
65 [step_name](int i) { return fmt::format(step_name + "{:d}.vtm", i); },
66 t, t0, dt, args["output"]["paraview"]["skip_frame"].get<int>());
67 }
68 }
69
70 void State::save_json(const Eigen::MatrixXd &sol)
71 {
72 const std::string out_path = resolve_output_path(args["output"]["json"]);
73 if (!out_path.empty())
74 {
75 std::ofstream out(out_path);
76 if (!out.is_open())
77 {
78 logger().error("Unable to save simulation JSON to {}", out_path);
79 return;
80 }
81 save_json(sol, out);
82 out.close();
83 }
84 }
85
86 void State::save_json(const Eigen::MatrixXd &sol, std::ostream &out)
87 {
88 if (!mesh)
89 {
90 logger().error("Load the mesh first!");
91 return;
92 }
93 if (sol.size() <= 0)
94 {
95 logger().error("Solve the problem first!");
96 return;
97 }
98
99 logger().info("Saving json...");
100
101 using json = nlohmann::json;
102 json j;
105 assembler->name(), iso_parametric(), args["output"]["advanced"]["sol_at_node"],
106 j);
107 out << j.dump(4) << std::endl;
108 }
109
110 void State::save_subsolve(const int i, const int t, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure)
111 {
112 if (!args["output"]["advanced"]["save_solve_sequence_debug"].get<bool>())
113 return;
114
116 solution_frames.emplace_back();
117
118 double dt = 1;
119 if (!args["time"].is_null())
120 dt = args["time"]["dt"];
121
123 resolve_output_path(fmt::format("solve_{:d}.vtu", i)),
124 *this, sol, pressure, t, dt,
127 }
128
129 void State::export_data(const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure)
130 {
131 if (!mesh)
132 {
133 logger().error("Load the mesh first!");
134 return;
135 }
136 if (n_bases <= 0)
137 {
138 logger().error("Build the bases first!");
139 return;
140 }
141 // if (rhs.size() <= 0)
142 // {
143 // logger().error("Assemble the rhs first!");
144 // return;
145 // }
146 if (sol.size() <= 0)
147 {
148 logger().error("Solve the problem first!");
149 return;
150 }
151
152 // Export vtu mesh of solution + wire mesh of deformed input
153 // + mesh colored with the bases
154 const std::string vis_mesh_path = resolve_output_path(args["output"]["paraview"]["file_name"]);
155 const std::string nodes_path = resolve_output_path(args["output"]["data"]["nodes"]);
156 const std::string solution_path = resolve_output_path(args["output"]["data"]["solution"]);
157 const std::string stress_path = resolve_output_path(args["output"]["data"]["stress_mat"]);
158 const std::string mises_path = resolve_output_path(args["output"]["data"]["mises"]);
159
160 double tend = args.value("tend", 1.0);
161 double dt = 1;
162 if (!args["time"].is_null())
163 dt = args["time"]["dt"];
164
166 *this, sol, pressure,
167 !args["time"].is_null(),
168 tend, dt,
170 vis_mesh_path,
171 nodes_path,
172 solution_path,
173 stress_path,
174 mises_path,
176 }
177
178 void State::save_restart_json(const double t0, const double dt, const int t) const
179 {
180 const std::string restart_json_path = args["output"]["restart_json"];
181 if (restart_json_path.empty())
182 return;
183
184 json restart_json;
185 restart_json["root_path"] = root_path();
186 restart_json["common"] = root_path();
187 restart_json["time"] = {{"t0", t0 + dt * t}};
188
189 restart_json["space"] = R"({
190 "remesh": {
191 "collapse": {
192 "abs_max_edge_length": -1,
193 "rel_max_edge_length": -1
194 }
195 }
196 })"_json;
197 restart_json["space"]["remesh"]["collapse"]["abs_max_edge_length"] = std::min(
198 args["space"]["remesh"]["collapse"]["abs_max_edge_length"].get<double>(),
199 starting_min_edge_length * args["space"]["remesh"]["collapse"]["rel_max_edge_length"].get<double>());
200 restart_json["space"]["remesh"]["collapse"]["rel_max_edge_length"] = std::numeric_limits<float>::max();
201
202 std::string rest_mesh_path = args["output"]["data"]["rest_mesh"].get<std::string>();
203 if (!rest_mesh_path.empty())
204 {
205 rest_mesh_path = resolve_output_path(fmt::format(args["output"]["data"]["rest_mesh"], t));
206
207 std::vector<json> patch;
208 if (args["geometry"].is_array())
209 {
210 const std::vector<json> in_geometry = args["geometry"];
211 for (int i = 0; i < in_geometry.size(); ++i)
212 {
213 if (!in_geometry[i]["is_obstacle"].get<bool>())
214 {
215 patch.push_back({
216 {"op", "remove"},
217 {"path", fmt::format("/geometry/{}", i)},
218 });
219 }
220 }
221
222 const int remaining_geometry = in_geometry.size() - patch.size();
223 assert(remaining_geometry >= 0);
224
225 patch.push_back({
226 {"op", "add"},
227 {"path", fmt::format("/geometry/{}", remaining_geometry > 0 ? "0" : "-")},
228 {"value",
229 {
230 // TODO: this does not set the surface selections
231 {"mesh", rest_mesh_path},
232 }},
233 });
234 }
235 else
236 {
237 assert(args["geometry"].is_object());
238 patch.push_back({
239 {"op", "remove"},
240 {"path", "/geometry"},
241 });
242 patch.push_back({
243 {"op", "replace"},
244 {"path", "/geometry"},
245 {"value",
246 {
247 // TODO: this does not set the surface selections
248 {"mesh", rest_mesh_path},
249 }},
250 });
251 }
252
253 restart_json["patch"] = patch;
254 }
255
256 restart_json["input"] = {{
257 "data",
258 {
259 {"state", resolve_output_path(fmt::format(args["output"]["data"]["state"], t))},
260 },
261 }};
262
263 std::ofstream file(resolve_output_path(fmt::format(restart_json_path, t)));
264 file << restart_json;
265 }
266} // namespace polyfem
#define POLYFEM_SCOPED_TIMER(...)
Definition Timer.hpp:10
io::OutRuntimeData timings
runtime statistics
Definition State.hpp:603
int n_bases
number of bases
Definition State.hpp:178
int n_pressure_bases
number of pressure bases
Definition State.hpp:180
const std::vector< basis::ElementBases > & geom_bases() const
Get a constant reference to the geometry mapping bases.
Definition State.hpp:223
std::string root_path() const
Get the root path for the state (e.g., args["root_path"] or ".")
std::shared_ptr< assembler::Assembler > assembler
assemblers
Definition State.hpp:155
void compute_errors(const Eigen::MatrixXd &sol)
computes all errors
void save_subsolve(const int i, const int t, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure)
saves a subsolve when save_solve_sequence_debug is true
io::OutGeometryData out_geom
visualization stuff
Definition State.hpp:601
std::string resolve_output_path(const std::string &path) const
Resolve output path relative to output_dir if the path is not absolute.
std::string output_dir
Directory for output files.
Definition State.hpp:593
void save_timestep(const double time, const int t, const double t0, const double dt, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure)
saves a timestep
void save_json(const Eigen::MatrixXd &sol, std::ostream &out)
saves the output statistic to a stream
double starting_min_edge_length
Definition State.hpp:606
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:466
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:168
json args
main input arguments containing all defaults
Definition State.hpp:101
void save_restart_json(const double t0, const double dt, const int t) const
Save a JSON sim file for restarting the simulation at time t.
std::vector< io::SolutionFrame > solution_frames
saves the frames in a vector instead of VTU
Definition State.hpp:599
io::OutStatsData stats
Other statistics.
Definition State.hpp:605
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
Definition State.hpp:171
bool iso_parametric() const
check if using iso parametric bases
Definition State.cpp:536
void export_data(const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure)
saves all data on the disk according to the input params
std::string resolve_input_path(const std::string &path, const bool only_if_exists=false) const
Resolve input path relative to root_path() if the path is not absolute.
bool is_contact_enabled() const
does the simulation has contact
Definition State.hpp:574
Eigen::VectorXi disc_orders
vector of discretization orders, used when not all elements have the same degree, one per element
Definition State.hpp:190
bool solve_export_to_file
flag to decide if exporting the time dependent solution to files or save it in the solution_frames ar...
Definition State.hpp:597
void save_pvd(const std::string &name, const std::function< std::string(int)> &vtu_names, int time_steps, double t0, double dt, int skip_frame=1) const
save a PVD of a time dependent simulation
Definition OutData.cpp:2301
void export_data(const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure, const bool is_time_dependent, const double tend_in, const double dt, const ExportOptions &opts, const std::string &vis_mesh_path, const std::string &nodes_path, const std::string &solution_path, const std::string &stress_path, const std::string &mises_path, const bool is_contact_enabled, std::vector< SolutionFrame > &solution_frames) const
exports everytihng, txt, vtu, etc
Definition OutData.cpp:935
void save_vtu(const std::string &path, const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure, const double t, const double dt, const ExportOptions &opts, const bool is_contact_enabled, std::vector< SolutionFrame > &solution_frames) const
saves the vtu file for time t
Definition OutData.cpp:1108
void save_json(const nlohmann::json &args, const int n_bases, const int n_pressure_bases, const Eigen::MatrixXd &sol, const mesh::Mesh &mesh, const Eigen::VectorXi &disc_orders, const assembler::Problem &problem, const OutRuntimeData &runtime, const std::string &formulation, const bool isoparametric, const int sol_at_node_id, nlohmann::json &j)
saves the output statistic to a json object
Definition OutData.cpp:2777
void compute_errors(const int n_bases, const std::vector< polyfem::basis::ElementBases > &bases, const std::vector< polyfem::basis::ElementBases > &gbases, const polyfem::mesh::Mesh &mesh, const assembler::Problem &problem, const double tend, const Eigen::MatrixXd &sol)
compute errors
Definition OutData.cpp:2544
std::string resolve_path(const std::string &path, const std::string &input_file_path, const bool only_if_exists=false)
bool is_param_valid(const json &params, const std::string &key)
Determine if a key exists and is non-null in a json object.
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
nlohmann::json json
Definition Common.hpp:9