Shamrock 2025.10.0
Astrophysical Code
Loading...
Searching...
No Matches
AVConfig.hpp
Go to the documentation of this file.
1// -------------------------------------------------------//
2//
3// SHAMROCK code for hydrodynamics
4// Copyright (c) 2021-2026 Timothée David--Cléris <tim.shamrock@proton.me>
5// SPDX-License-Identifier: CeCILL Free Software License Agreement v2.1
6// Shamrock is licensed under the CeCILL 2.1 License, see LICENSE for more information
7//
8// -------------------------------------------------------//
9
10#pragma once
11
21#include "shambackends/vec.hpp"
24#include <nlohmann/json.hpp>
25#include <variant>
26
27namespace shammodels::sph {
28
32 template<class Tscal>
34 static constexpr std::string_view variant_type_name = "none";
35 };
36
37 template<class Tscal>
38 inline void to_json(nlohmann::json &j, const AVConfig_None<Tscal> &p) {}
39
40 template<class Tscal>
41 inline void from_json(const nlohmann::json &j, AVConfig_None<Tscal> &p) {
42 p = {};
43 }
44
48 template<class Tscal>
50 static constexpr std::string_view variant_type_name = "constant";
51
52 Tscal alpha_u = 1.0;
53 Tscal alpha_AV = 1.0;
54 Tscal beta_AV = 2.0;
55 };
56
57 template<class Tscal>
58 inline void to_json(nlohmann::json &j, const AVConfig_Constant<Tscal> &p) {
59 j = {
60 {"alpha_u", p.alpha_u},
61 {"alpha_AV", p.alpha_AV},
62 {"beta_AV", p.beta_AV},
63 };
64 }
65
66 template<class Tscal>
67 inline void from_json(const nlohmann::json &j, AVConfig_Constant<Tscal> &p) {
68 j.at("alpha_u").get_to(p.alpha_u);
69 j.at("alpha_AV").get_to(p.alpha_AV);
70 j.at("beta_AV").get_to(p.beta_AV);
71 }
72
77 template<class Tscal>
79 static constexpr std::string_view variant_type_name = "varying_mm97";
80
81 Tscal alpha_min = 0.1;
82 Tscal alpha_max = 1.0;
83 Tscal sigma_decay = 0.1;
84 Tscal alpha_u = 1.0;
85 Tscal beta_AV = 2.0;
86 };
87
88 template<class Tscal>
89 inline void to_json(nlohmann::json &j, const AVConfig_VaryingMM97<Tscal> &p) {
90 j = {
91 {"alpha_min", p.alpha_min},
92 {"alpha_max", p.alpha_max},
93 {"sigma_decay", p.sigma_decay},
94 {"alpha_u", p.alpha_u},
95 {"beta_AV", p.beta_AV},
96 };
97 }
98
99 template<class Tscal>
100 inline void from_json(const nlohmann::json &j, AVConfig_VaryingMM97<Tscal> &p) {
101 j.at("alpha_min").get_to(p.alpha_min);
102 j.at("alpha_max").get_to(p.alpha_max);
103 j.at("sigma_decay").get_to(p.sigma_decay);
104 j.at("alpha_u").get_to(p.alpha_u);
105 j.at("beta_AV").get_to(p.beta_AV);
106 }
107
112 template<class Tscal>
114 static constexpr std::string_view variant_type_name = "varying_cd10";
115
116 Tscal alpha_min = 0.1;
117 Tscal alpha_max = 1.0;
118 Tscal sigma_decay = 0.1;
119 Tscal alpha_u = 1.0;
120 Tscal beta_AV = 2.0;
121 };
122
123 template<class Tscal>
124 inline void to_json(nlohmann::json &j, const AVConfig_VaryingCD10<Tscal> &p) {
125 j = {
126 {"alpha_min", p.alpha_min},
127 {"alpha_max", p.alpha_max},
128 {"sigma_decay", p.sigma_decay},
129 {"alpha_u", p.alpha_u},
130 {"beta_AV", p.beta_AV},
131 };
132 }
133
134 template<class Tscal>
135 inline void from_json(const nlohmann::json &j, AVConfig_VaryingCD10<Tscal> &p) {
136 j.at("alpha_min").get_to(p.alpha_min);
137 j.at("alpha_max").get_to(p.alpha_max);
138 j.at("sigma_decay").get_to(p.sigma_decay);
139 j.at("alpha_u").get_to(p.alpha_u);
140 j.at("beta_AV").get_to(p.beta_AV);
141 }
142
146 template<class Tscal>
148 static constexpr std::string_view variant_type_name = "constant_disc";
149
150 Tscal alpha_AV = 1.0;
151 Tscal alpha_u = 1.0;
152 Tscal beta_AV = 2.0;
153 };
154
155 template<class Tscal>
156 inline void to_json(nlohmann::json &j, const AVConfig_ConstantDisc<Tscal> &p) {
157 j = {
158 {"alpha_AV", p.alpha_AV},
159 {"alpha_u", p.alpha_u},
160 {"beta_AV", p.beta_AV},
161 };
162 }
163
164 template<class Tscal>
165 inline void from_json(const nlohmann::json &j, AVConfig_ConstantDisc<Tscal> &p) {
166 j.at("alpha_AV").get_to(p.alpha_AV);
167 j.at("alpha_u").get_to(p.alpha_u);
168 j.at("beta_AV").get_to(p.beta_AV);
169 }
170
184 template<class Tvec>
185 struct AVConfig;
186} // namespace shammodels::sph
187
188template<class Tvec>
190
192 using Tscal = shambase::VecComponent<Tvec>;
195
201
203 using Variant = std::variant<None, Constant, VaryingMM97, VaryingCD10, ConstantDisc>;
204
207
209 void set(Variant v) { config = v; }
210
221 Tscal alpha_min, Tscal alpha_max, Tscal sigma_decay, Tscal alpha_u, Tscal beta_AV) {
222 set(VaryingCD10{alpha_min, alpha_max, sigma_decay, alpha_u, beta_AV});
223 }
224
230 inline bool has_alphaAV_field() {
231 bool is_varying_alpha
232 = bool(std::get_if<VaryingMM97>(&config)) || bool(std::get_if<VaryingCD10>(&config));
233 return is_varying_alpha;
234 }
235
242 inline bool has_divv_field() {
243 bool is_varying_alpha
244 = bool(std::get_if<VaryingMM97>(&config)) || bool(std::get_if<VaryingCD10>(&config));
245 return is_varying_alpha;
246 }
247
254 inline bool has_curlv_field() {
255 bool is_varying_alpha = bool(std::get_if<VaryingCD10>(&config));
256 return is_varying_alpha;
257 }
258
264 inline bool has_dtdivv_field() {
265 bool is_varying_alpha = bool(std::get_if<VaryingCD10>(&config));
266 return is_varying_alpha;
267 }
268
274 inline bool has_field_soundspeed() {
275
276 // this should not be needed idealy, but we need the pressure on the ghosts and
277 // we don't want to communicate it as it can be recomputed from the other fields
278 // hence we copy the soundspeed at the end of the step to a field in the patchdata
279 // cf eos module there is another soundspeed field available as a Compute field
280 // unifying the patchdata and the ghosts is really needed ...
281
282 bool is_varying_alpha
283 = bool(std::get_if<VaryingMM97>(&config)) || bool(std::get_if<VaryingCD10>(&config));
284 return is_varying_alpha;
285 }
286
290 inline void print_status() {
291 logger::raw_ln("--- artificial viscosity config");
292
293 if (None *v = std::get_if<None>(&config)) {
294 logger::raw_ln(" Config Type : None (No artificial viscosity)");
295 } else if (Constant *v = std::get_if<Constant>(&config)) {
296 logger::raw_ln(" Config Type : Constant (Constant artificial viscosity)");
297 logger::raw_ln(" alpha_u =", v->alpha_u);
298 logger::raw_ln(" alpha_AV =", v->alpha_AV);
299 logger::raw_ln(" beta_AV =", v->beta_AV);
300 } else if (VaryingMM97 *v = std::get_if<VaryingMM97>(&config)) {
301 logger::raw_ln(" Config Type : VaryingMM97 (Morris & Monaghan 1997)");
302 logger::raw_ln(" alpha_min =", v->alpha_min);
303 logger::raw_ln(" alpha_max =", v->alpha_max);
304 logger::raw_ln(" sigma_decay =", v->sigma_decay);
305 logger::raw_ln(" alpha_u =", v->alpha_u);
306 logger::raw_ln(" beta_AV =", v->beta_AV);
307 } else if (VaryingCD10 *v = std::get_if<VaryingCD10>(&config)) {
308 logger::raw_ln(" Config Type : VaryingCD10 (Cullen & Dehnen 2010)");
309 logger::raw_ln(" alpha_min =", v->alpha_min);
310 logger::raw_ln(" alpha_max =", v->alpha_max);
311 logger::raw_ln(" sigma_decay =", v->sigma_decay);
312 logger::raw_ln(" alpha_u =", v->alpha_u);
313 logger::raw_ln(" beta_AV =", v->beta_AV);
314 } else if (ConstantDisc *v = std::get_if<ConstantDisc>(&config)) {
315 logger::raw_ln(" Config Type : constant disc");
316 logger::raw_ln(" alpha_AV =", v->alpha_AV);
317 logger::raw_ln(" alpha_u =", v->alpha_u);
318 logger::raw_ln(" beta_AV =", v->beta_AV);
319 } else {
321 }
322
323 logger::raw_ln("--- artificial viscosity config (deduced)");
324
325 logger::raw_ln("-------------");
326 }
327
328 inline std::optional<Tscal> get_alpha_u() {
329 if (const Constant *v = std::get_if<Constant>(&config)) {
330 return v->alpha_u;
331 } else if (const VaryingMM97 *v = std::get_if<VaryingMM97>(&config)) {
332 return v->alpha_u;
333 } else if (const VaryingCD10 *v = std::get_if<VaryingCD10>(&config)) {
334 return v->alpha_u;
335 } else if (const ConstantDisc *v = std::get_if<ConstantDisc>(&config)) {
336 return v->alpha_u;
337 } else if (const None *v = std::get_if<None>(&config)) {
338 return std::nullopt;
339 } else {
341 return std::nullopt;
342 }
343 }
344};
345
346namespace shammodels::sph {
347
354 template<class Tvec>
355 inline void to_json(nlohmann::json &j, const AVConfig<Tvec> &p) {
356 std::visit(
357 [&](const auto &value) {
358 j = value;
359 j["type"] = value.variant_type_name;
360 },
361 p.config);
362 }
369 template<class Tvec>
370 inline void from_json(const nlohmann::json &j, AVConfig<Tvec> &p) {
371 using T = AVConfig<Tvec>;
372
373 using Tscal = shambase::VecComponent<Tvec>;
374
375 if (!j.contains("type") && !j.contains("av_type")) {
377 "neither \"type\" nor \"av_type\" in this json, can not infer type json=\n"
378 + j.dump(4));
379 }
380
381 std::string av_type;
382 if (j.contains("type")) {
383 j.at("type").get_to(av_type);
384 } else {
385 j.at("av_type").get_to(av_type);
386 }
387
388 shamrock::json_deserialize_variant(j, av_type, p.config);
389 }
390
391} // namespace shammodels::sph
std::uint32_t u32
32 bit unsigned integer
void throw_with_loc(std::string message, SourceLocation loc=SourceLocation{})
Throw an exception and append the source location to it.
void throw_unimplemented(SourceLocation loc=SourceLocation{})
Throw a std::runtime_error saying that the function is unimplemented.
namespace for the sph model
void to_json(nlohmann::json &j, const EOSConfig< Tvec > &p)
Serialize EOSConfig to json.
Definition EOSConfig.cpp:43
void from_json(const nlohmann::json &j, EOSConfig< Tvec > &p)
Deserializes an EOSConfig<Tvec> from a JSON object.
Contains traits and utilities for backend related types.
Constant artificial viscosity for alpha disc viscosity.
Definition AVConfig.hpp:147
Constant artificial viscosity: .
Definition AVConfig.hpp:49
No artificial viscosity: .
Definition AVConfig.hpp:33
Configuration for the Artificial Viscosity (AV)
Definition AVConfig.hpp:189
void set_varying_cd10(Tscal alpha_min, Tscal alpha_max, Tscal sigma_decay, Tscal alpha_u, Tscal beta_AV)
Sets the configuration to use a varying Cullen & Dehnen 2010 artificial viscosity.
Definition AVConfig.hpp:220
static constexpr u32 dim
Number of dimensions of the problem.
Definition AVConfig.hpp:194
std::variant< None, Constant, VaryingMM97, VaryingCD10, ConstantDisc > Variant
Variant of all types of artificial viscosity possible.
Definition AVConfig.hpp:203
Variant config
The actual configuration (default to constant viscosity)
Definition AVConfig.hpp:206
void set(Variant v)
Set the configuration.
Definition AVConfig.hpp:209
shambase::VecComponent< Tvec > Tscal
Type of the components of the vector of coordinates.
Definition AVConfig.hpp:192
void print_status()
Prints the status of the artificial viscosity configuration.
Definition AVConfig.hpp:290
bool has_alphaAV_field()
Checks if the current configuration has an alpha artificial viscosity field.
Definition AVConfig.hpp:230
bool has_divv_field()
Checks if the current configuration need the divergence of the velocity field.
Definition AVConfig.hpp:242
bool has_field_soundspeed()
Checks if the current configuration requires the soundspeed field.
Definition AVConfig.hpp:274
bool has_curlv_field()
Checks if the current configuration need the curl of the velocity field.
Definition AVConfig.hpp:254
bool has_dtdivv_field()
Checks if the current configuration has a dtdivv field.
Definition AVConfig.hpp:264