SERiF 0.0.1a
3+1D Stellar Structure and Evolution
Loading...
Searching...
No Matches
probe.cpp
Go to the documentation of this file.
1/* ***********************************************************************
2//
3// Copyright (C) 2025 -- The 4D-STAR Collaboration
4// File Author: Emily Boudreaux
5// Last Modified: March 18, 2025
6//
7// 4DSSE is free software; you can use it and/or modify
8// it under the terms and restrictions the GNU General Library Public
9// License version 3 (GPLv3) as published by the Free Software Foundation.
10//
11// 4DSSE is distributed in the hope that it will be useful,
12// but WITHOUT ANY WARRANTY; without even the implied warranty of
13// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14// See the GNU Library General Public License for more details.
15//
16// You should have received a copy of the GNU Library General Public License
17// along with this software; if not, write to the Free Software
18// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19//
20// *********************************************************************** */
21#include "quill/Backend.h"
22#include "quill/Frontend.h"
23#include "quill/Logger.h"
24#include "quill/sinks/ConsoleSink.h"
25#include "quill/sinks/FileSink.h"
26#include "quill/LogMacros.h"
27
28#include <stdexcept>
29#include <string>
30#include <iostream>
31#include <chrono>
32#include <cmath>
33#include <vector>
34#include <utility>
35#include <filesystem>
36
37#include "mfem.hpp"
38
39#include "config.h"
40#include "probe.h"
41
42#include "warning_control.h"
43
44
45namespace serif::probe {
46
47void pause() {
48 std::cout << "Execution paused. Please press enter to continue..."
49 << std::endl; // Use endl to flush
50 std::cin.get();
51}
52
53void wait(int seconds) {
54 std::this_thread::sleep_for(std::chrono::seconds(seconds));
55}
56
57void glVisView(mfem::GridFunction& u, mfem::Mesh& mesh,
58 const std::string& windowTitle, const std::string& keyset) {
59 serif::config::Config& config = serif::config::Config::getInstance();
60 quill::Logger* logger = LogManager::getInstance().getLogger("log");
61 if (config.get<bool>("Probe:GLVis:Visualization", true)) {
62 std::string usedKeyset;
63 LOG_INFO(logger, "Visualizing solution using GLVis...");
64 LOG_INFO(logger, "Window title: {}", windowTitle);
65 if (keyset.empty()) {
66 usedKeyset = config.get<std::string>("Probe:GLVis:DefaultKeyset", "");
67 } else {
68 usedKeyset = keyset;
69 }
70 LOG_INFO(logger, "Keyset: {}", usedKeyset);
71 const auto vishost = config.get<std::string>("Probe:GLVis:Host", "localhost");
72 const auto visport = config.get<int>("Probe:GLVis:Port", 19916);
73 mfem::socketstream sol_sock(vishost.c_str(), visport);
74 sol_sock.precision(8);
75 sol_sock << "solution\n" << mesh << u
76 << "window_title '" << windowTitle <<
77 "'\n" << "keys " << usedKeyset << "\n";
78 sol_sock << std::flush;
79 }
80}
81
82void glVisView(mfem::Vector &vec, mfem::FiniteElementSpace &fes, const std::string &windowTitle, const std::string& keyset) {
83 mfem::GridFunction gf(&fes);
84
86 gf.SetData(vec);
88 glVisView(gf, *fes.GetMesh(), windowTitle, keyset);
89}
90
91double getMeshRadius(mfem::Mesh& mesh) {
92 mesh.EnsureNodes();
93 const mfem::GridFunction *nodes = mesh.GetNodes();
94 double *node_data = nodes->GetData(); // THIS IS KEY
95
96 int data_size = nodes->Size();
97 double max_radius = 0.0;
98
99 for (int i = 0; i < data_size; ++i) {
100 if (node_data[i] > max_radius) {
101 max_radius = node_data[i];
102 }
103 }
104 return max_radius;
105
106}
107
108std::pair<std::vector<double>, std::vector<double>> getRaySolution(mfem::GridFunction& u, mfem::Mesh& mesh,
109 const std::vector<double>& rayDirection,
110 int numSamples, std::string filename) {
111 serif::config::Config& config = serif::config::Config::getInstance();
113 quill::Logger* logger = logManager.getLogger("log");
114 LOG_INFO(logger, "Getting ray solution...");
115 // Check if the directory to write to exists
116 // If it does not exist and MakeDir is true create it
117 // Otherwise throw an exception
118 bool makeDir = config.get<bool>("Probe:GetRaySolution:MakeDir", true);
119 std::filesystem::path path = filename;
120
121 if (makeDir) {
122 std::filesystem::path dir = path.parent_path();
123 if (!std::filesystem::exists(dir)) {
124 LOG_INFO(logger, "Creating directory {}", dir.string());
125 std::filesystem::create_directories(dir);
126 }
127 } else {
128 if (!std::filesystem::exists(path.parent_path())) {
129 throw(std::runtime_error("Directory " + path.parent_path().string() + " does not exist"));
130 }
131 }
132
133 std::vector<double> samples;
134 samples.reserve(numSamples);
135 double x, y, z, r, sampleValue;
136 double radius = getMeshRadius(mesh);
137 mfem::Vector rayOrigin(3); rayOrigin = 0.0;
138 mfem::DenseMatrix rayPoints(3, numSamples);
139 std::vector<double> radialPoints;
140 radialPoints.reserve(numSamples);
141 mfem::ElementTransformation* Trans = nullptr;
142 mfem::Vector physicalCoords;
143
144 for (int i = 0; i < numSamples; i++) {
145 r = i * radius / numSamples;
146 // Let rayDirection = (theta, phi) that the ray will be cast too
147 x = r * std::sin(rayDirection[0]) * std::cos(rayDirection[1]);
148 y = r * std::sin(rayDirection[0]) * std::sin(rayDirection[1]);
149 z = r * std::cos(rayDirection[0]);
150 rayPoints(0, i) = x;
151 rayPoints(1, i) = y;
152 rayPoints(2, i) = z;
153 }
154
155 mfem::Array<int> elementIds;
156 mfem::Array<mfem::IntegrationPoint> ips;
157 mesh.FindPoints(rayPoints, elementIds, ips);
158 for (int i = 0; i < elementIds.Size(); i++) {
159 Trans = mesh.GetElementTransformation(elementIds[i]);
160 Trans->Transform(ips[i], physicalCoords);
161 double r = std::sqrt(physicalCoords[0] * physicalCoords[0] +
162 physicalCoords[1] * physicalCoords[1] +
163 physicalCoords[2] * physicalCoords[2]);
164 radialPoints.push_back(r);
165
166 int elementId = elementIds[i];
167 mfem::IntegrationPoint ip = ips[i];
168 if (elementId >= 0) { // Check if the point was found in an element
169 sampleValue = u.GetValue(elementId, ip);
170 LOG_DEBUG(logger, "Probe::getRaySolution() : Ray point {} found in element {} with r={:0.2f} and theta={:0.2f}", i, elementId, r, sampleValue);
171 samples.push_back(sampleValue);
172 } else { // If the point was not found in an element
173 LOG_INFO(logger, "Probe::getRaySolution() : Ray point {} not found", i);
174 samples.push_back(0.0);
175 }
176 }
177 std::pair<std::vector<double>, std::vector<double>> samplesPair(radialPoints, samples);
178
179 if (!filename.empty()) {
180 std::ofstream file(filename);
181 if (file.is_open()) {
182 file << "r,u\n";
183 for (int i = 0; i < numSamples; i++) {
184 file << samplesPair.first[i] << "," << samplesPair.second[i] << "\n";
185 }
186 file << std::endl;
187 file.close();
188 } else {
189 throw(std::runtime_error("Could not open file " + filename));
190 }
191 }
192
193 return samplesPair;
194}
195
196std::pair<std::vector<double>, std::vector<double>> getRaySolution(mfem::Vector &vec, mfem::FiniteElementSpace &fes,
197 const std::vector<double>& rayDirection,
198 int numSamples, std::string filename) {
199 mfem::GridFunction gf(&fes);
201 gf.SetData(vec);
203 return getRaySolution(gf, *fes.GetMesh(), rayDirection, numSamples, filename);
204}
205
207 serif::config::Config& config = serif::config::Config::getInstance();
208 quill::Backend::start();
209 auto CLILogger = quill::Frontend::create_or_get_logger(
210 "root",
211 quill::Frontend::create_or_get_sink<quill::ConsoleSink>("sink_id_1"));
212
213 newFileLogger(config.get<std::string>("Probe:LogManager:DefaultLogName", "4DSSE.log"), "log");
214 loggerMap.emplace("stdout", CLILogger);
215}
216
217LogManager::~LogManager() = default;
218
219quill::Logger* LogManager::getLogger(const std::string& loggerName) {
220 auto it = loggerMap.find(loggerName); // Find *once*
221 if (it == loggerMap.end()) {
222 throw std::runtime_error("Cannot find logger " + loggerName);
223 }
224 return it->second; // Return the raw pointer from the shared_ptr
225}
226
227std::vector<std::string> LogManager::getLoggerNames() {
228 std::vector<std::string> loggerNames;
229 loggerNames.reserve(loggerMap.size());
230 for (const auto& pair : loggerMap) { // Use range-based for loop and const auto&
231 loggerNames.push_back(pair.first);
232 }
233 return loggerNames;
234}
235
236std::vector<quill::Logger*> LogManager::getLoggers() {
237 std::vector<quill::Logger*> loggers;
238 loggers.reserve(loggerMap.size());
239 for (const auto& pair : loggerMap) {
240 loggers.push_back(pair.second); // Get the raw pointer
241 }
242 return loggers;
243}
244
245quill::Logger* LogManager::newFileLogger(const std::string& filename,
246 const std::string& loggerName) {
247 auto file_sink = quill::Frontend::create_or_get_sink<quill::FileSink>(
248 filename,
249 []() {
250 quill::FileSinkConfig cfg;
251 cfg.set_open_mode('w');
252 return cfg;
253 }(),
254 quill::FileEventNotifier{});
255 // Get the raw pointer.
256 quill::Logger* rawLogger = quill::Frontend::create_or_get_logger(loggerName, std::move(file_sink));
257
258 // Create a shared_ptr from the raw pointer.
259 loggerMap.emplace(loggerName, rawLogger);
260 return rawLogger;
261}
262
263} // namespace Probe
Class to manage logging operations.
Definition probe.h:79
LogManager()
Private constructor for singleton pattern.
Definition probe.cpp:206
quill::Logger * getLogger(const std::string &loggerName)
Get a logger by name.
Definition probe.cpp:219
std::vector< quill::Logger * > getLoggers()
Get all loggers.
Definition probe.cpp:236
quill::Logger * newFileLogger(const std::string &filename, const std::string &loggerName)
Create a new file logger.
Definition probe.cpp:245
static LogManager & getInstance()
Get the singleton instance of LogManager.
Definition probe.h:103
std::map< std::string, quill::Logger * > loggerMap
Definition probe.h:92
std::vector< std::string > getLoggerNames()
Get the names of all loggers.
Definition probe.cpp:227
The Probe namespace contains utility functions for debugging and logging.
Definition probe.cpp:45
double getMeshRadius(mfem::Mesh &mesh)
Definition probe.cpp:91
void wait(int seconds)
Wait for a specified number of seconds.
Definition probe.cpp:53
void glVisView(mfem::GridFunction &u, mfem::Mesh &mesh, const std::string &windowTitle, const std::string &keyset)
Visualize a solution using GLVis.
Definition probe.cpp:57
std::pair< std::vector< double >, std::vector< double > > getRaySolution(mfem::GridFunction &u, mfem::Mesh &mesh, const std::vector< double > &rayDirection, int numSamples, std::string filename)
Definition probe.cpp:108
void pause()
Pause the execution and wait for user input.
Definition probe.cpp:47
#define DEPRECATION_WARNING_OFF
#define DEPRECATION_WARNING_ON