OpenPFC  0.1.4
Phase Field Crystal simulation framework
Loading...
Searching...
No Matches
app.hpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2025 VTT Technical Research Centre of Finland Ltd
2// SPDX-License-Identifier: AGPL-3.0-or-later
3
22#ifndef PFC_UI_APP_HPP
23#define PFC_UI_APP_HPP
24
26#include "from_json.hpp"
27#include "json_helpers.hpp"
28#include "openpfc/openpfc.hpp"
32#include <filesystem>
33#include <fstream>
34#include <iostream>
35#include <nlohmann/json.hpp>
36#include <toml++/toml.hpp>
37
38namespace pfc {
39namespace ui {
40
45template <class ConcreteModel> class App {
46private:
47 MPI_Comm m_comm;
48 MPI_Worker m_worker;
49 bool rank0;
50 json m_settings;
51
52 double m_total_steptime = 0.0;
53 double m_total_fft_time = 0.0;
54 double m_steptime = 0.0;
55 double m_fft_time = 0.0;
56 double m_avg_steptime = 0.0;
57 int m_steps_done = 0;
58
59 // save detailed timing information for each mpi rank and step?
60 bool m_detailed_timing = false;
61 bool m_detailed_timing_print = false;
62 bool m_detailed_timing_write = false;
63 std::string m_detailed_timing_filename = "timing.bin";
64
65 // read settings from file (JSON or TOML format)
66 json read_settings(int argc, char *argv[]) {
67 if (argc <= 1) {
68 if (rank0) {
69 std::cerr << "Error: Configuration file required.\n";
70 std::cerr << "Usage: " << argv[0] << " <config.json|config.toml>\n";
71 }
73 }
74
75 std::filesystem::path file(argv[1]);
76 if (!std::filesystem::exists(file)) {
77 if (rank0) std::cerr << "Error: File " << file << " does not exist!\n";
79 }
80
81 auto ext = file.extension().string();
82 json settings;
83
84 if (ext == ".toml") {
85 if (rank0) std::cout << "Reading TOML configuration from " << file << "\n\n";
86 try {
87 auto toml_data = toml::parse_file(file.string());
88 settings = utils::toml_to_json(toml_data);
89 } catch (const toml::parse_error &err) {
90 if (rank0) {
91 std::cerr << "Error parsing TOML file: " << err.description() << "\n";
92 std::cerr << " at line " << err.source().begin.line << ", column "
93 << err.source().begin.column << "\n";
94 }
96 }
97 } else if (ext == ".json") {
98 if (rank0) std::cout << "Reading JSON configuration from " << file << "\n\n";
99 std::ifstream input_file(file);
100 try {
102 } catch (const nlohmann::json::parse_error &err) {
103 if (rank0) {
104 std::cerr << "Error parsing JSON file: " << err.what() << "\n";
105 std::cerr << " at byte position " << err.byte << "\n";
106 }
108 }
109 } else {
110 if (rank0) {
111 std::cerr << "Error: Unsupported file format: " << ext << "\n";
112 std::cerr << "Supported formats: .json, .toml\n";
113 }
115 }
116
117 return settings;
118 }
119
120public:
121 App(int argc, char *argv[], MPI_Comm comm = MPI_COMM_WORLD)
122 : m_comm(comm), m_worker(MPI_Worker(argc, argv, comm)),
123 rank0(m_worker.get_rank() == 0), m_settings(read_settings(argc, argv)) {}
124
125 App(const json &settings, MPI_Comm comm = MPI_COMM_WORLD)
126 : m_comm(comm), m_worker(MPI_Worker(0, nullptr, comm)),
127 rank0(m_worker.get_rank() == 0), m_settings(settings) {}
128
129 bool create_results_dir(const std::string &output) {
130 std::filesystem::path results_dir(output);
131 if (results_dir.has_filename()) results_dir = results_dir.parent_path();
132 if (!std::filesystem::exists(results_dir)) {
133 std::cout << "Results dir " << results_dir << " does not exist, creating\n";
134 std::filesystem::create_directories(results_dir);
135 return true;
136 } else {
137 std::cout << "Warning: results dir " << results_dir << " already exists\n";
138 return false;
139 }
140 }
141
142 void read_detailed_timing_configuration() {
143 if (m_settings.contains("detailed_timing")) {
144 auto timing = m_settings["detailed_timing"];
145 if (timing.contains("enabled")) m_detailed_timing = timing["enabled"];
146 if (timing.contains("print")) m_detailed_timing_print = timing["print"];
147 if (timing.contains("write")) m_detailed_timing_write = timing["write"];
148 if (timing.contains("filename"))
149 m_detailed_timing_filename = timing["filename"];
150 }
151 }
152
153 void add_result_writers(Simulator &sim) {
154 std::cout << "Adding results writers" << std::endl;
155 if (m_settings.contains("saveat") && m_settings.contains("fields") &&
156 m_settings["saveat"] > 0) {
157 for (const auto &field : m_settings["fields"]) {
158 std::string name = field["name"];
159 std::string data = field["data"];
160 if (rank0) create_results_dir(data);
161 std::cout << "Writing field " << name << " to " << data << std::endl;
162 sim.add_results_writer(name, std::make_unique<BinaryWriter>(data));
163 }
164 } else {
165 std::cout << "Warning: not writing results to anywhere." << std::endl;
166 std::cout << "To write results, add ResultsWriter to model." << std::endl;
167 }
168 }
169
170 void add_initial_conditions(Simulator &sim) {
171 if (!m_settings.contains("initial_conditions")) {
172 std::cout << "WARNING: no initial conditions are set!" << std::endl;
173 return;
174 }
175 std::cout << "Adding initial conditions" << std::endl;
176 for (const json &params : m_settings["initial_conditions"]) {
177 std::cout << "Creating initial condition from data " << params << std::endl;
178 if (!params.contains("type")) {
179 std::cout << "Warning: no type is set for initial condition!" << std::endl;
180 continue;
181 }
182 std::string type = params["type"];
184 if (!params.contains("target")) {
185 std::cout << "Warning: no target is set for initial condition! Using "
186 "target 'default'"
187 << std::endl;
188 } else {
189 std::string target = params["target"];
190 std::cout << "Setting initial condition target to " << target << std::endl;
191 field_modifier->set_field_name(target);
192 }
193 sim.add_initial_conditions(std::move(field_modifier));
194 }
195 }
196
197 void add_boundary_conditions(Simulator &sim) {
198 if (!m_settings.contains("boundary_conditions")) {
199 std::cout << "Warning: no boundary conditions are set!" << std::endl;
200 return;
201 }
202 std::cout << "Adding boundary conditions" << std::endl;
203 for (const json &params : m_settings["boundary_conditions"]) {
204 std::cout << "Creating boundary condition from data " << params << std::endl;
205 if (!params.contains("type")) {
206 std::cout << "Warning: no type is set for initial condition!" << std::endl;
207 continue;
208 }
209 std::string type = params["type"];
211 if (!params.contains("target")) {
212 std::cout << "Warning: no target is set for boundary condition! Using "
213 "target 'default'"
214 << std::endl;
215 } else {
216 std::string target = params["target"];
217 std::cout << "Setting boundary condition target to " << target << std::endl;
218 field_modifier->set_field_name(target);
219 }
220 sim.add_boundary_conditions(std::move(field_modifier));
221 }
222 }
223
224 int main() {
225 std::cout << "Reading configuration from json file:" << std::endl;
226 std::cout << m_settings.dump(4) << "\n\n";
227
228 World world(ui::from_json<World>(m_settings));
229 std::cout << "World: " << world << std::endl;
230
231 int num_ranks = m_worker.get_num_ranks();
232 int rank_id = m_worker.get_rank();
233 auto decomp = decomposition::create(world, num_ranks);
234
235 // Create FFT with default FFTW backend for now
236 // Note: Runtime backend selection via create_with_backend() can be added when
237 // needed
238 auto options =
239 m_settings.contains("plan_options")
240 ? ui::from_json<heffte::plan_options>(m_settings["plan_options"])
242
243 auto fft_layout = fft::layout::create(decomp, 0);
244 auto fft = fft::create(fft_layout, rank_id, options);
245 Time time(ui::from_json<Time>(m_settings));
246 ConcreteModel model(fft, world);
248
249 if (m_settings.contains("model") && m_settings["model"].contains("params")) {
250 from_json(m_settings["model"]["params"], model);
251 }
252 read_detailed_timing_configuration();
253
254 std::cout << "Initializing model... " << std::endl;
255 model.initialize(time.get_dt());
256
257 // Report memory usage
258 {
259 size_t model_mem = model.get_allocated_memory_bytes();
260 size_t fft_mem = fft.get_allocated_memory_bytes();
262 pfc::Logger logger{pfc::LogLevel::Info, rank_id};
263 pfc::utils::report_memory_usage(usage, world, logger, m_comm);
264 }
265
266 add_result_writers(simulator);
267 add_initial_conditions(simulator);
268 add_boundary_conditions(simulator);
269
270 if (m_settings.contains("simulator")) {
271 const json &j = m_settings["simulator"];
272 if (j.contains("result_counter")) {
273 if (!j["result_counter"].is_number_integer()) {
274 throw std::invalid_argument(
275 "Invalid JSON input: missing or invalid 'result_counter' field.");
276 }
277 int result_counter = (int)j["result_counter"] + 1;
278 simulator.set_result_counter(result_counter);
279 }
280 if (j.contains("increment")) {
281 if (!j["increment"].is_number_integer()) {
282 throw std::invalid_argument(
283 "Invalid JSON input: missing or invalid 'increment' field.");
284 }
285 int increment = j["increment"];
286 time.set_increment(increment);
287 }
288 }
289
290 std::cout << "Applying initial conditions" << std::endl;
291 simulator.apply_initial_conditions();
292 if (time.get_increment() == 0) {
293 std::cout << "First increment: apply boundary conditions" << std::endl;
294 simulator.apply_boundary_conditions();
295 simulator.write_results();
296 }
297
298 while (!time.done()) {
299 time.next(); // increase increment counter by 1
300 simulator.apply_boundary_conditions();
301
302 double l_steptime = 0.0; // l = local for this mpi process
303 double l_fft_time = 0.0;
304 MPI_Barrier(m_comm);
306 step(simulator, model);
307 MPI_Barrier(m_comm);
309 l_fft_time = fft.get_fft_time();
310
311 if (m_detailed_timing) {
312 double timing[2] = {l_steptime, l_fft_time};
313 MPI_Send(timing, 2, MPI_DOUBLE, 0, 42, m_comm);
314 if (m_worker.get_rank() == 0) {
315 // Use the num_ranks from outer scope (line 230) to avoid shadowing
316 double all_timing[num_ranks][2];
317 for (int rank = 0; rank < num_ranks; rank++) {
318 MPI_Recv(all_timing[rank], 2, MPI_DOUBLE, rank, 42, m_comm,
320 }
321 auto inc = time.get_increment();
322 if (m_detailed_timing_print) {
323 auto old_precision = std::cout.precision(6);
324 std::cout << "Timing information for all processes:" << std::endl;
325 std::cout << "step;rank;step_time;fft_time" << std::endl;
326 for (int rank = 0; rank < num_ranks; rank++) {
327 std::cout << inc << ";" << rank << ";" << all_timing[rank][0] << ";"
328 << all_timing[rank][1] << std::endl;
329 }
330 std::cout.precision(old_precision);
331 }
332 if (m_detailed_timing_write) {
333 // so we end up to a binary file, and opening with e.g. Python
334 // np.fromfile("timing.bin").reshape(n_steps, n_procs, 2)
335 std::ofstream outfile(m_detailed_timing_filename, std::ios::app);
336 outfile.write((const char *)timing, sizeof(double) * 2 * num_ranks);
337 outfile.close();
338 }
339 }
340 }
341
342 // max reduction over all mpi processes
343 MPI_Reduce(&l_steptime, &m_steptime, 1, MPI_DOUBLE, MPI_MAX, 0, m_comm);
344 MPI_Reduce(&l_fft_time, &m_fft_time, 1, MPI_DOUBLE, MPI_MAX, 0, m_comm);
345
346 if (time.do_save()) {
347 simulator.apply_boundary_conditions();
348 simulator.write_results();
349 }
350
351 // Calculate eta from average step time.
352 // Use exponential moving average when steps > 3.
353 m_avg_steptime = m_steptime;
354 if (m_steps_done > 3) {
355 m_avg_steptime = 0.01 * m_steptime + 0.99 * m_avg_steptime;
356 }
357 int increment = time.get_increment();
358 double t = time.get_current(), t1 = time.get_t1();
359 double eta_i = (t1 - t) / time.get_dt();
360 double eta_t = eta_i * m_avg_steptime;
361 double other_time = m_steptime - m_fft_time;
362 std::cout << "Step " << increment << " done in " << m_steptime << " s ";
363 std::cout << "(" << m_fft_time << " s FFT, " << other_time << " s other). ";
364 std::cout << "Simulation time: " << t << " / " << t1;
365 std::cout << " (" << (t / t1 * 100) << " % done). ";
366 std::cout << "ETA: " << pfc::utils::TimeLeft(eta_t) << std::endl;
367
368 m_total_steptime += m_steptime;
369 m_total_fft_time += m_fft_time;
370 m_steps_done += 1;
371 }
372
373 double avg_steptime = m_total_steptime / m_steps_done;
374 double avg_fft_time = m_total_fft_time / m_steps_done;
376 double p_fft = avg_fft_time / avg_steptime * 100.0;
377 double p_oth = avg_oth_time / avg_steptime * 100.0;
378 std::cout << "\nSimulated " << m_steps_done
379 << " steps. Average times:" << std::endl;
380 std::cout << "Step time: " << avg_steptime << " s" << std::endl;
381 std::cout << "FFT time: " << avg_fft_time << " s / " << p_fft << " %"
382 << std::endl;
383 std::cout << "Other time: " << avg_oth_time << " s / " << p_oth << " %"
384 << std::endl;
385
386 return 0;
387 }
388};
389
390} // namespace ui
391} // namespace pfc
392
393#endif // PFC_UI_APP_HPP
An MPI worker class that wraps MPI_Init and MPI_Finalize.
Definition worker.hpp:55
int get_num_ranks() const
Returns the number of processes in the MPI communicator.
Definition worker.hpp:122
int get_rank() const
Returns the rank of this worker process in the MPI communicator.
Definition worker.hpp:115
The Simulator class is responsible for running the simulation of the model.
Definition simulator.hpp:63
Definition time.hpp:233
The main json-based application.
Definition app.hpp:45
Definition timeleft.hpp:39
Registry for field modifiers (initial conditions and boundary conditions)
std::unique_ptr< FieldModifier > create_field_modifier(const std::string &type, const json &params)
Create an instance of a field modifier based on its type.
Definition field_modifier_registry.hpp:128
JSON deserialization functions for OpenPFC types.
JSON utility functions for configuration parsing.
Memory usage reporting utilities for parallel simulations.
Main convenience header that includes all OpenPFC public API components.
Lightweight logger configuration.
Definition logging.hpp:29
Memory usage statistics for a single MPI rank.
Definition memory_reporter.hpp:65
Represents the global simulation domain (the "world").
Definition world.hpp:91
World(const Int3 &lower, const Int3 &upper, const CoordinateSystem< T > &cs)
Constructs a World object.
Estimated time remaining display.
Convert TOML data to nlohmann::json format.