/home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/common/eom.hh Source File#

DiFfRG: /home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/common/eom.hh Source File
DiFfRG
eom.hh
Go to the documentation of this file.
1#pragma once
2
3// external libraries
4#include <deal.II/base/point.h>
5#include <deal.II/dofs/dof_handler.h>
6#include <deal.II/grid/grid_tools.h>
7#include <deal.II/numerics/fe_field_function.h>
8#include <gsl/gsl_multiroots.h>
9#include <gsl/gsl_vector.h>
10#include <gsl/gsl_vector_double.h>
11
12// standard library
13#include <algorithm>
14#include <iostream>
15#include <iterator>
16
17namespace DiFfRG
18{
19 namespace internal
20 {
21 template <int dim> int gsl_unwrap(const gsl_vector *gsl_x, void *params, gsl_vector *gsl_f)
22 {
23 const int subdim = gsl_x->size;
24
25 dealii::Point<dim> x{};
26 for (int i = 0; i < subdim; ++i)
27 x[i] = gsl_vector_get(gsl_x, i);
28
29 auto fp = static_cast<std::function<std::array<double, dim>(const dealii::Point<dim> &)> *>(params);
30 const auto f = (*fp)(x);
31
32 for (int i = 0; i < subdim; ++i)
33 gsl_vector_set(gsl_f, i, f[i]);
34
35 return GSL_SUCCESS;
36 }
37
38 using namespace dealii;
39 template <int dim, typename Cell>
40 bool is_in_cell(const Cell &cell, const Point<dim> &point, const Mapping<dim> &mapping)
41 {
42 try {
43 Point<dim> qp = mapping.transform_real_to_unit_cell(cell, point);
44 if (GeometryInfo<dim>::is_inside_unit_cell(qp))
45 return true;
46 else
47 return false;
48 } catch (const typename Mapping<dim>::ExcTransformationFailed &) {
49 // transformation failed, so assume the point is outside
50 return false;
51 }
52 }
53
54 template <int dim>
55 dealii::Point<dim> get_origin(const dealii::DoFHandler<dim> &dof_handler,
56 typename dealii::DoFHandler<dim>::cell_iterator &EoM_cell)
57 {
58 using namespace dealii;
59 using CellIterator = typename dealii::DoFHandler<dim>::cell_iterator;
60
61 auto l1_norm = [](const Point<dim> &p) {
62 double norm = 0.;
63 for (uint d = 0; d < dim; ++d)
64 norm += std::abs(p[d]);
65 return norm;
66 };
67
68 std::vector<Point<dim>> candidates;
69
70 auto iterate_cell = [&](const CellIterator &cell) {
71 std::array<Point<dim>, GeometryInfo<dim>::vertices_per_cell> vertices;
72 for (uint i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i) {
73 vertices[i] = cell->vertex(i);
74 }
75 // choose the one with the smallest l1 norm
76 auto min_it = std::min_element(vertices.begin(), vertices.end(),
77 [&](const auto &p1, const auto &p2) { return l1_norm(p1) < l1_norm(p2); });
78 // add the point to the candidates
79 candidates.push_back(*min_it);
80 };
81
82 for (const auto &cell : dof_handler.active_cell_iterators())
83 iterate_cell(cell);
84
85 // find the candidate with the smallest l1 norm
86 auto min_it = std::min_element(candidates.begin(), candidates.end(),
87 [&](const auto &p1, const auto &p2) { return l1_norm(p1) < l1_norm(p2); });
88
89 Point<dim> origin = *min_it;
90 EoM_cell = GridTools::find_active_cell_around_point(dof_handler, origin);
91
92 return origin;
93 }
94
95 } // namespace internal
96
114 template <typename VectorType, typename EoMFUN, typename EoMPFUN>
115 dealii::Point<1> get_EoM_point_1D(
116 typename dealii::DoFHandler<1>::cell_iterator &EoM_cell, const VectorType &sol,
117 const dealii::DoFHandler<1> &dof_handler, const dealii::Mapping<1> &mapping, const EoMFUN &get_EoM,
118 const EoMPFUN &EoM_postprocess = [](const auto &p, const auto &values) { return p; },
119 const double EoM_abs_tol = 1e-8, const uint max_iter = 100)
120 {
121 constexpr uint dim = 1;
122
123 using namespace dealii;
124 using CellIterator = typename dealii::DoFHandler<dim>::cell_iterator;
125 using EoMType = std::array<double, dim>;
126
127 auto EoM = Point<dim>();
128 Vector<typename VectorType::value_type> values(dof_handler.get_fe().n_components());
129 Functions::FEFieldFunction<dim, VectorType> fe_function(dof_handler, sol, mapping);
130
131 // We start by investigating the origin.
132 const auto origin = internal::get_origin(dof_handler, EoM_cell);
133 fe_function.set_active_cell(EoM_cell);
134 fe_function.vector_value(origin, values);
135 const auto origin_val = get_EoM(origin, values);
136
137 const auto secondary_point = (origin + EoM_cell->center()) / 2.;
138 fe_function.vector_value(secondary_point, values);
139 const auto secondary_val = get_EoM(secondary_point, values);
140
141 bool origin_EoM = true;
142 bool secondary_EoM = true;
143 for (uint d = 0; d < dim; ++d) {
144 origin_EoM = origin_EoM && (origin_val[d] > EoM_abs_tol);
145 secondary_EoM = secondary_EoM && (secondary_val[d] > 0);
146 }
147
148 if (origin_EoM && secondary_EoM) {
149 EoM = origin;
150 return EoM;
151 }
152
153 auto check_cell = [&](const CellIterator &cell) {
154 // Obtain the values at the vertices of the cell.
155 std::array<double, GeometryInfo<dim>::vertices_per_cell> EoM_vals;
156 std::array<Point<dim>, GeometryInfo<dim>::vertices_per_cell> vertices;
157 // Vector<typename VectorType::value_type> values(dof_handler.get_fe().n_components());
158 // Functions::FEFieldFunction<dim, VectorType> fe_function(dof_handler, sol, mapping);
159
160 fe_function.set_active_cell(cell);
161 for (uint i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i) {
162 vertices[i] = cell->vertex(i);
163 fe_function.vector_value(vertices[i], values);
164 EoM_vals[i] = get_EoM(vertices[i], values)[0];
165 }
166
167 bool cell_has_EoM = false;
168 // Find if the cell has an EoM point, i.e. if the values at some vertices have different signs.
169 for (uint i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
170 for (uint j = 0; j < i; ++j)
171 if (EoM_vals[i] * EoM_vals[j] < 0.) {
172 cell_has_EoM = true;
173 i = GeometryInfo<dim>::vertices_per_cell;
174 break;
175 }
176
177 return cell_has_EoM;
178 };
179
180 auto find_EoM = [&](const CellIterator &cell, CellIterator &m_EoM_cell, Point<dim> &m_EoM,
181 double &EoM_val) -> bool {
182 std::array<Point<dim>, GeometryInfo<dim>::vertices_per_cell> vertices;
183 for (uint i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
184 vertices[i] = cell->vertex(i);
185 fe_function.set_active_cell(cell);
186
187 auto p1 = vertices[0];
188 auto p2 = vertices[1];
189 if (p2[0] < p1[0]) std::swap(p1, p2);
190 auto p = m_EoM[0] <= p2[0] && m_EoM[0] >= p1[0] ? m_EoM : (p1 + p2) / 2.;
191
192 fe_function.vector_value(p, values);
193 EoM_val = get_EoM(p, values)[0];
194 double err = 1e100;
195 uint iter = 0;
196
197 while (err > EoM_abs_tol) {
198 if (EoM_val < 0.)
199 p1 = p;
200 else
201 p2 = p;
202 p = (p1 + p2) / 2.;
203
204 fe_function.vector_value(p, values);
205 EoM_val = get_EoM(p, values)[0];
206 err = std::abs(EoM_val);
207
208 if (iter > max_iter) {
209 m_EoM = EoM_postprocess(p, values);
210 m_EoM_cell = cell;
211 return false;
212 }
213 iter++;
214 }
215
216 m_EoM = EoM_postprocess(p, values);
217 m_EoM_cell = cell;
218
219 return true;
220 };
221
222 std::vector<typename dealii::DoFHandler<dim>::cell_iterator> cell_candidates;
223 std::vector<double> value_candidates;
224 std::vector<Point<dim>> EoM_candidates;
225
226 CellIterator t_EoM_cell;
227 Point<dim> t_EoM;
228 double t_EoM_value = 0.;
229
230 // Check the cell from the previous iteration first
231 if (EoM_cell != dof_handler.active_cell_iterators().end())
232 if (check_cell(EoM_cell)) {
233 if (find_EoM(EoM_cell, t_EoM_cell, t_EoM, t_EoM_value)) {
234 EoM_cell = t_EoM_cell;
235 EoM = t_EoM;
236 return EoM;
237 }
238 cell_candidates.push_back(t_EoM_cell);
239 value_candidates.push_back(std::abs(t_EoM_value));
240 EoM_candidates.push_back(t_EoM);
241 }
242
243 for (const auto &cell : dof_handler.active_cell_iterators())
244 if (check_cell(cell)) {
245 if (find_EoM(cell, t_EoM_cell, t_EoM, t_EoM_value)) {
246 EoM_cell = t_EoM_cell;
247 EoM = t_EoM;
248 return EoM;
249 }
250 cell_candidates.push_back(t_EoM_cell);
251 value_candidates.push_back(std::abs(t_EoM_value));
252 EoM_candidates.push_back(t_EoM);
253 }
254
255 if (cell_candidates.size() == 1) {
256 EoM_cell = cell_candidates[0];
257 EoM = EoM_candidates[0];
258 return EoM;
259 } else if (cell_candidates.size() == 0) {
260 EoM_cell = GridTools::find_active_cell_around_point(dof_handler, origin);
261 return origin;
262 } else if (cell_candidates.size() > 1) {
263 // If we have more than one candidate, we choose the one with the smallest EoM value.
264 auto min_it = std::min_element(value_candidates.begin(), value_candidates.end());
265 auto min_idx = std::distance(value_candidates.begin(), min_it);
266 EoM_cell = cell_candidates[min_idx];
267 EoM = EoM_candidates[min_idx];
268 return EoM;
269 }
270
271 find_EoM(cell_candidates[0], EoM_cell, EoM, t_EoM_value);
272
273 return EoM;
274 }
275
293 template <int dim, typename VectorType, typename EoMFUN, typename EoMPFUN>
294 dealii::Point<dim> get_EoM_point_ND(
295 typename dealii::DoFHandler<dim>::cell_iterator &EoM_cell, const VectorType &sol,
296 const dealii::DoFHandler<dim> &dof_handler, const dealii::Mapping<dim> &mapping, const EoMFUN &get_EoM,
297 const EoMPFUN &EoM_postprocess = [](const auto &p, const auto &values) { return p; },
298 const double EoM_abs_tol = 1e-8, const uint max_iter = 100)
299 {
300 using namespace dealii;
301 using CellIterator = typename dealii::DoFHandler<dim>::cell_iterator;
302 using EoMType = std::array<double, dim>;
303
304 Vector<typename VectorType::value_type> values(dof_handler.get_fe().n_components());
305 Functions::FEFieldFunction<dim, VectorType> fe_function(dof_handler, sol, mapping);
306
307 // We start by investigating the origin.
308 const auto origin = internal::get_origin(dof_handler, EoM_cell);
309 fe_function.set_active_cell(EoM_cell);
310 fe_function.vector_value(origin, values);
311 const auto origin_val = get_EoM(origin, values);
312
313 const auto secondary_point = (origin + EoM_cell->center()) / 2.;
314 fe_function.vector_value(secondary_point, values);
315 const auto secondary_val = get_EoM(secondary_point, values);
316
317 // this actually constrains the EoM to one axis, possibly.
318 // in a future version, we should exploit this to reduce the dimensionality of the problem.
319 std::array<bool, dim> axis_restrictions{{}};
320 for (uint d = 0; d < dim; ++d)
321 axis_restrictions[d] = (origin_val[d] >= 0) && (secondary_val[d] >= 0);
322
323 const bool any_axis_restricted =
324 std::any_of(std::begin(axis_restrictions), std::end(axis_restrictions), [](bool i) { return i; });
325 const bool all_axis_restricted =
326 std::all_of(std::begin(axis_restrictions), std::end(axis_restrictions), [](bool i) { return i; });
327
328 auto EoM = origin;
329 EoM_cell = GridTools::find_active_cell_around_point(dof_handler, EoM);
330 if (all_axis_restricted) {
331 EoM = origin;
332 std::cout << "All axis restricted" << std::endl;
333 return EoM;
334 }
335
336 std::cout << "Restrictions: ";
337 for (uint d = 0; d < dim; ++d)
338 std::cout << axis_restrictions[d] << " ";
339 std::cout << std::endl;
340
341 auto check_cell = [&](const CellIterator &cell) -> bool {
342 if (any_axis_restricted && !cell->has_boundary_lines()) return false;
343
344 // Obtain the values at the vertices of the cell.
345 std::array<EoMType, GeometryInfo<dim>::vertices_per_cell> EoM_vals;
346 std::array<Point<dim>, GeometryInfo<dim>::vertices_per_cell> vertices;
347 Vector<typename VectorType::value_type> values(dof_handler.get_fe().n_components());
348 Functions::FEFieldFunction<dim, VectorType> fe_function(dof_handler, sol, mapping);
349
350 fe_function.set_active_cell(cell);
351 for (uint i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i) {
352 vertices[i] = cell->vertex(i);
353 fe_function.vector_value(vertices[i], values);
354 EoM_vals[i] = get_EoM(vertices[i], values);
355 }
356
357 if (any_axis_restricted)
358 for (uint d = 0; d < dim; ++d)
359 if (axis_restrictions[d]) {
360 bool has_zero_boundary = false;
361 for (uint i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
362 if (is_close(origin[d], vertices[i][d], 1e-12)) {
363 has_zero_boundary = true;
364 break;
365 }
366 if (!has_zero_boundary) {
367 std::cout << "No zero boundary for cell with vertices: ";
368 for (uint i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i) {
369 std::cout << "vertex " << i << ": ";
370 for (uint d = 0; d < dim; ++d)
371 std::cout << vertices[i][d] << " ";
372 std::cout << std::endl;
373 }
374
375 return false;
376 }
377 }
378
379 std::array<bool, dim> cell_has_EoM{{}};
380 // Find if the cell has an EoM point, i.e. if the values at some vertices have different signs.
381 if (!any_axis_restricted) {
382 for (uint i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
383 for (uint j = 0; j < i; ++j)
384 for (uint d = 0; d < dim; ++d)
385 cell_has_EoM[d] = cell_has_EoM[d] || (EoM_vals[i][d] * EoM_vals[j][d] < 0.);
386 } else {
387 std::vector<bool> valid_vertices(GeometryInfo<dim>::vertices_per_cell, true);
388
389 for (uint d = 0; d < dim; ++d)
390 if (axis_restrictions[d])
391 for (uint i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
392 // if there is a vertex with smaller coordinate in direction d, we set the vertex to invalid
393 for (uint j = 0; j < i; ++j) {
394 if (vertices[i][d] > vertices[j][d]) {
395 valid_vertices[i] = false;
396 break;
397 }
398 }
399
400 for (uint i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
401 if (valid_vertices[i])
402 for (uint j = 0; j < i; ++j)
403 if (valid_vertices[j])
404 for (uint d = 0; d < dim; ++d)
405 cell_has_EoM[d] = cell_has_EoM[d] || (EoM_vals[i][d] * EoM_vals[j][d] < 0.);
406 }
407
408 bool has_EoM = true;
409 for (uint d = 0; d < dim; ++d)
410 if (!axis_restrictions[d]) has_EoM = has_EoM && cell_has_EoM[d];
411 // if all components have a potential crossing, we return true
412 return has_EoM;
413 };
414
415 gsl_set_error_handler_off();
416
417 auto find_EoM = [&](const CellIterator &cell, CellIterator &m_EoM_cell, Point<dim> &m_EoM,
418 double &EoM_val) -> bool {
419 Functions::FEFieldFunction<dim, VectorType> fe_function(dof_handler, sol, mapping);
420 Vector<typename VectorType::value_type> values(dof_handler.get_fe().n_components());
421 fe_function.set_active_cell(cell);
422
423 // find the cell boundaries
424 std::array<std::array<double, 2>, dim> cell_boundaries{{}};
425 for (uint d = 0; d < dim; ++d) {
426 cell_boundaries[d][0] = cell->center()[d];
427 cell_boundaries[d][1] = cell->center()[d];
428 }
429 for (uint d = 0; d < dim; ++d) {
430 for (uint i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i) {
431 cell_boundaries[d][0] = std::min(cell_boundaries[d][0], cell->vertex(i)[d]);
432 cell_boundaries[d][1] = std::max(cell_boundaries[d][1], cell->vertex(i)[d]);
433 }
434 }
435
436 std::cout << "Cell boundaries: ";
437 for (uint d = 0; d < dim; ++d)
438 std::cout << "[" << cell_boundaries[d][0] << ", " << cell_boundaries[d][1] << "] ";
439 std::cout << std::endl;
440
441 int subdim = dim;
442 for (uint d = 0; d < dim; ++d)
443 if (axis_restrictions[d]) subdim--;
444
445 if (subdim == 1) {
446 uint dir = 0;
447 for (uint d = 0; d < dim; ++d)
448 if (!axis_restrictions[d]) dir = d;
449
450 // utlize the fact that we have a 1D problem and we can just do a bisection
451 Point<dim> p1{};
452 Point<dim> p2{};
453 for (uint d = 0; d < dim; ++d) {
454 p1[d] = cell_boundaries[d][0];
455 p2[d] = cell_boundaries[d][0];
456 }
457 p2[dir] = cell_boundaries[dir][1];
458 auto p = (p1 + p2) / 2.;
459
460 std::cout << "SUBDIM TRYING AT p = " << p << std::endl;
461
462 fe_function.vector_value(p, values);
463 EoM_val = get_EoM(p, values)[dir];
464 double err = 1e100;
465 uint iter = 0;
466
467 while (err > EoM_abs_tol) {
468 if (EoM_val < 0.)
469 p1 = p;
470 else
471 p2 = p;
472 p = (p1 + p2) / 2.;
473
474 fe_function.vector_value(p, values);
475 EoM_val = get_EoM(p, values)[dir];
476 err = std::abs(EoM_val);
477
478 if (iter > max_iter) {
479 m_EoM = EoM_postprocess(p, values);
480 m_EoM_cell = cell;
481 return false;
482 }
483 iter++;
484 }
485
486 m_EoM = EoM_postprocess(p, values);
487 m_EoM_cell = cell;
488
489 return true;
490 }
491
492 std::function<std::array<double, dim>(const dealii::Point<dim> &)> eval_on_point =
493 [&](const Point<dim> &p) -> std::array<double, dim> {
494 Point<dim> p_proj = p;
495 for (uint d = 0, i = 0; d < dim; ++d)
496 p_proj[d] = axis_restrictions[d] ? origin[d] : p[i++];
497
498 // check if the point is inside the cell
499 std::array<double, dim> out_distance{{}};
500 for (uint d = 0; d < dim; ++d) {
501 if (p_proj[d] < cell_boundaries[d][0]) {
502 out_distance[d] = std::abs(p_proj[d] - cell_boundaries[d][0]);
503 p_proj[d] = cell_boundaries[d][0];
504 std::cout << "Point outside cell" << std::endl;
505 std::cout << "Point: " << p << std::endl;
506 } else if (p_proj[d] > cell_boundaries[d][1]) {
507 out_distance[d] = std::abs(p_proj[d] - cell_boundaries[d][1]);
508 p_proj[d] = cell_boundaries[d][1];
509 std::cout << "Point outside cell" << std::endl;
510 std::cout << "Point: " << p << std::endl;
511 }
512 }
513
514 try {
515 fe_function.vector_value(p_proj, values);
516 } catch (const std::exception &e) {
517 // if p is outside the triangulation, give a default value
518 std::cerr << "Warning: EoM evaluation failed at point " << p_proj << ": " << e.what() << std::endl;
519 return std::array<double, dim>{{
520 std::numeric_limits<double>::quiet_NaN(),
521 }};
522 } catch (...) {
523 std::cerr << "Warning: EoM evaluation failed at point " << p_proj << " with unknown exception" << std::endl;
524 return std::array<double, dim>{{
525 std::numeric_limits<double>::quiet_NaN(),
526 }};
527 }
528
529 auto EoM = get_EoM(p, values);
530
531 // if (out_distance[d] > 0), linearly extrapolate the value
532 for (uint d = 0; d < dim; ++d)
533 if (out_distance[d] > 0) {
534 std::cout << "Extrapolating value from " << EoM[d];
535 EoM[d] = is_close(EoM[d], 0.) ? out_distance[d] : EoM[d] * (1 + out_distance[d]);
536 std::cout << " to " << EoM[d] << std::endl;
537 }
538
539 // reshuflle the values to the correct order
540 std::array<double, dim> EoM_out{{}};
541 for (uint d = 0, i = 0; d < dim; ++d)
542 if (!axis_restrictions[d]) EoM_out[d] = EoM[i++];
543
544 return EoM_out;
545 };
546
547 const auto cell_center = cell->center();
548 std::vector<double> subdim_center(subdim);
549 for (uint d = 0, i = 0; d < dim; ++d)
550 if (!axis_restrictions[d]) subdim_center[i++] = cell_center[d];
551
552 // Create GSL multiroot solver
553 const gsl_multiroot_fsolver_type *T = gsl_multiroot_fsolver_hybrids;
554 gsl_multiroot_fsolver *s = gsl_multiroot_fsolver_alloc(T, subdim);
555
556 // Create GSL function
557 gsl_multiroot_function f = {&internal::gsl_unwrap<dim>, (size_t)subdim, &eval_on_point};
558
559 // Create initial guess
560 gsl_vector *x = gsl_vector_alloc(subdim);
561 for (int d = 0; d < subdim; ++d)
562 gsl_vector_set(x, d, subdim_center[d]);
563
564 // Set the solver with the function and initial guess
565 gsl_multiroot_fsolver_set(s, &f, x);
566
567 // start the iteration
568 uint iter = 0;
569 int status;
570 do {
571 iter++;
572 status = gsl_multiroot_fsolver_iterate(s);
573
574 if (status) break;
575
576 status = gsl_multiroot_test_residual(s->f, EoM_abs_tol);
577
578 } while (status == GSL_CONTINUE && iter < max_iter);
579
580 EoM_val = 0.;
581 for (uint d = 0, sd = 0; d < dim; ++d)
582 m_EoM[d] = axis_restrictions[d] ? origin[d] : gsl_vector_get(s->x, sd++);
583
584 // don't leak memory. Stupid C
585 gsl_multiroot_fsolver_free(s);
586 gsl_vector_free(x);
587
588 if (status != GSL_SUCCESS) return false;
589
590 m_EoM = EoM_postprocess(m_EoM, values);
591 m_EoM_cell = internal::is_in_cell(cell, m_EoM, mapping)
592 ? cell
593 : GridTools::find_active_cell_around_point(dof_handler, m_EoM);
594
595 return true;
596 };
597
598 std::vector<typename dealii::DoFHandler<dim>::cell_iterator> cell_candidates;
599 std::vector<double> value_candidates;
600 std::vector<Point<dim>> EoM_candidates;
601
602 // mutex for the EoM cell
603 // std::mutex EoM_cell_mutex;
604
605 const uint n_active_cells = dof_handler.get_triangulation().n_active_cells();
606
607 // this needs to be done smarter, but for now we just iterate over all cells
608 // preferrably, we would divide the full domain into smaller subdomains and then traverse a kind of tree.
609 for (uint index = 0; index < n_active_cells; ++index) {
610 auto cell = dof_handler.begin_active();
611 std::advance(cell, index);
612 if (check_cell(cell)) {
613 Point<dim> t_EoM = cell->center();
614 CellIterator t_EoM_cell = cell;
615 double t_EoM_value = 0.;
616
617 // lock the mutex
618 // std::lock_guard<std::mutex> lock(EoM_cell_mutex);
619
620 if (find_EoM(cell, t_EoM_cell, t_EoM, t_EoM_value)) {
621 std::cout << "Found EoM point in cell " << index << std::endl;
622 std::cout << "EoM: " << t_EoM << std::endl;
623 for (uint d = 0; d < dim; ++d)
624 if (t_EoM[d] < origin[d]) t_EoM[d] = origin[d];
625 EoM = t_EoM;
626 EoM_cell = dealii::GridTools::find_active_cell_around_point(dof_handler, EoM);
627 return EoM;
628 }
629
630 cell_candidates.push_back(t_EoM_cell);
631 EoM_candidates.push_back(t_EoM);
632 value_candidates.push_back(std::abs(t_EoM_value));
633 }
634 }
635
636 std::cout << "Found " << cell_candidates.size() << " candidates." << std::endl;
637
638 double EoM_value = 0.;
639
640 if (cell_candidates.size() == 1) {
641 if (find_EoM(cell_candidates[0], EoM_cell, EoM, EoM_value)) return EoM;
642 } else if (cell_candidates.size() == 0) {
643 EoM_cell = GridTools::find_active_cell_around_point(dof_handler, origin);
644 return origin;
645 } else if (cell_candidates.size() > 1) {
646 // If we have more than one candidate, we choose the one with the smallest EoM value.
647 auto min_it = std::min_element(value_candidates.begin(), value_candidates.end());
648 auto min_idx = std::distance(value_candidates.begin(), min_it);
649 EoM_cell = cell_candidates[min_idx];
650 EoM = EoM_candidates[min_idx];
651 return EoM;
652 }
653
654 return EoM;
655 }
656
673 template <int dim, typename VectorType, typename EoMFUN, typename EoMPFUN>
674 dealii::Point<dim> get_EoM_point(
675 typename dealii::DoFHandler<dim>::cell_iterator &EoM_cell, const VectorType &sol,
676 const dealii::DoFHandler<dim> &dof_handler, const dealii::Mapping<dim> &mapping, const EoMFUN &get_EoM,
677 const EoMPFUN &EoM_postprocess = [](const auto &p, const auto &values) { return p; },
678 const double EoM_abs_tol = 1e-5, const uint max_iter = 100)
679 {
680 if (max_iter == 0) return internal::get_origin(dof_handler, EoM_cell);
681
682 if constexpr (dim == 1) {
683 return get_EoM_point_1D(EoM_cell, sol, dof_handler, mapping, get_EoM, EoM_postprocess, EoM_abs_tol, max_iter);
684 } else if constexpr (dim > 1) {
685 return get_EoM_point_ND(EoM_cell, sol, dof_handler, mapping, get_EoM, EoM_postprocess, EoM_abs_tol, max_iter);
686 } else
687 throw std::runtime_error("get_EoM_point is not yet implemented for dim > 1. Please set EoM_max_iter = 0.");
688 }
689} // namespace DiFfRG
dealii::Point< dim > get_origin(const dealii::DoFHandler< dim > &dof_handler, typename dealii::DoFHandler< dim >::cell_iterator &EoM_cell)
Definition eom.hh:55
int gsl_unwrap(const gsl_vector *gsl_x, void *params, gsl_vector *gsl_f)
Definition eom.hh:21
bool is_in_cell(const Cell &cell, const Point< dim > &point, const Mapping< dim > &mapping)
Definition eom.hh:40
Definition complex_math.hh:10
dealii::Point< 1 > get_EoM_point_1D(typename dealii::DoFHandler< 1 >::cell_iterator &EoM_cell, const VectorType &sol, const dealii::DoFHandler< 1 > &dof_handler, const dealii::Mapping< 1 > &mapping, const EoMFUN &get_EoM, const EoMPFUN &EoM_postprocess=[](const auto &p, const auto &values) { return p;}, const double EoM_abs_tol=1e-8, const uint max_iter=100)
Get the EoM point for a given solution and model in 1D. This is done by first checking the origin,...
Definition eom.hh:115
dealii::Point< dim > get_EoM_point(typename dealii::DoFHandler< dim >::cell_iterator &EoM_cell, const VectorType &sol, const dealii::DoFHandler< dim > &dof_handler, const dealii::Mapping< dim > &mapping, const EoMFUN &get_EoM, const EoMPFUN &EoM_postprocess=[](const auto &p, const auto &values) { return p;}, const double EoM_abs_tol=1e-5, const uint max_iter=100)
Get the EoM point for a given solution and model.
Definition eom.hh:674
bool KOKKOS_INLINE_FUNCTION is_close(T1 a, T2 b, T3 eps_)
Function to evaluate whether two floats are equal to numerical precision. Tests for both relative and...
Definition math.hh:168
unsigned int uint
Definition utils.hh:24
dealii::Point< dim > get_EoM_point_ND(typename dealii::DoFHandler< dim >::cell_iterator &EoM_cell, const VectorType &sol, const dealii::DoFHandler< dim > &dof_handler, const dealii::Mapping< dim > &mapping, const EoMFUN &get_EoM, const EoMPFUN &EoM_postprocess=[](const auto &p, const auto &values) { return p;}, const double EoM_abs_tol=1e-8, const uint max_iter=100)
Get the EoM point for a given solution and model in 2D. This is done by first checking the origin,...
Definition eom.hh:294