PRISMS-PF Manual
Loading...
Searching...
No Matches
periodic_distance.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#pragma once
5
6#include <deal.II/base/point.h>
7
10
12
13#include <prismspf/config.h>
14
16
22template <unsigned int dim, typename real>
23real
24distance(const dealii::Point<dim, real> &point1,
25 const dealii::Point<dim, real> &point2,
26 const RectangularMesh<dim> &rectangular_mesh)
27{
28 if (rectangular_mesh.periodic_directions.empty())
29 {
30 return point1.distance(point2);
31 }
32 // else
33 real dist = real(0.0);
34 for (unsigned int d = 0; d < dim; ++d)
35 {
36 using std::min;
37 real delta = point2[d] - point1[d];
38 if (rectangular_mesh.periodic_directions.contains(d))
39 {
40 const real length = rectangular_mesh.size[d];
41 const real half_length = length / 2.0;
43 }
44 dist += delta * delta;
45 }
46 using std::sqrt;
47 return sqrt(dist);
48}
49
@ Value
Use value of the variable as a criterion for refinement.
Definition grid_refiner_criterion.h:31
Definition conditional_ostreams.cc:20
PRISMS_PF_BEGIN_NAMESPACE real distance(const dealii::Point< dim, real > &point1, const dealii::Point< dim, real > &point2, const RectangularMesh< dim > &rectangular_mesh)
Calculate the distance between two points considering periodic boundaries. Note: templates must be pr...
Definition periodic_distance.h:24
Class for rectangular mesh parameters.
Definition spatial_discretization.h:43
dealii::Tensor< 1, dim, double > size
Domain extents in each cartesian direction.
Definition spatial_discretization.h:84
std::set< unsigned int > periodic_directions
Which directions have periodic conditions.
Definition spatial_discretization.h:99
Real pmod(const Real &value, const OtherReal &modulus)
Positive moldulo (remainder) returns the normal remainder. (c++ fmod is defined abnormally for negati...
Definition utilities.h:119