Skip to content

Commit

Permalink
Renaming and style fixes
Browse files Browse the repository at this point in the history
  • Loading branch information
hughcars committed Aug 29, 2023
1 parent 55910f1 commit 7485eb3
Show file tree
Hide file tree
Showing 3 changed files with 17 additions and 19 deletions.
16 changes: 8 additions & 8 deletions palace/fem/lumpedelement.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -70,10 +70,10 @@ class UniformElementData : public LumpedElementData
{
// Check that the bounding box discovered matches the area. This validates that the
// boundary elements form a right angled quadrilateral port.
constexpr double rel_tol = 1.0e-6;
constexpr double rel_tol = 1.0e-3;
double A = GetArea(fespace);
MFEM_VERIFY((!bounding_box.planar || std::abs(A - bounding_box.Area()) / A < rel_tol),
"Assumed bounding box area "
MFEM_VERIFY((!bounding_box.planar || (std::abs(A - bounding_box.Area()) / A < rel_tol)),
"Discovered bounding box area "
<< bounding_box.Area() << " and integrated area " << A
<< " do not match: Planar port geometry is not a quadrilateral!");

Expand Down Expand Up @@ -103,12 +103,13 @@ class UniformElementData : public LumpedElementData
Mpi::Warning(
"User specified direction {:.3e} does not align with either bounding box "
"axis up to {:.3e} degrees.\n"
"Axis 1: {:.3e} ({:.3e} degrees)\nAxis 2: {:.3e} ({:.3e} degrees)\nAxis 3: {:.3e} ({:.3e} degrees)!\n",
"Axis 1: {:.3e} ({:.3e} degrees)\nAxis 2: {:.3e} ({:.3e} degrees)\nAxis 3: "
"{:.3e} ({:.3e} degrees)!\n",
input_dir, angle_warning_deg, normal_0, deviation_deg[0], normal_1,
deviation_deg[1], normal_2, deviation_deg[2]);
}
MFEM_VERIFY(std::any_of(deviation_deg.begin(), deviation_deg.end(),
[](double x) { return x < angle_error_deg; }),
[angle_error_deg](double x) { return x < angle_error_deg; }),
"Specified direction does not align sufficiently with bounding box axes: "
<< deviation_deg[0] << ' ' << deviation_deg[1] << ' '
<< deviation_deg[2] << " tolerance " << angle_error_deg << '!');
Expand All @@ -121,9 +122,8 @@ class UniformElementData : public LumpedElementData
std::min_element(deviation_deg.begin(), deviation_deg.end()))];

MFEM_ASSERT(
(l - mesh::GetDirectionalExtent(*fespace.GetParMesh(), marker, true, input_dir)) /
l <
1e-3,
(l - mesh::GetProjectedLength(*fespace.GetParMesh(), marker, true, input_dir)) / l <
rel_tol,
"Bounding box discovered length should match projected length");
w = A / l;
}
Expand Down
12 changes: 5 additions & 7 deletions palace/utils/geodata.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -763,8 +763,6 @@ BoundingBox BoundingBoxFromPointCloud(MPI_Comm comm,
// two opposite edges of the cuboid. Now look for a component that maximizes distance
// from the planar system: complete the axes with a cross, then use a dot product to
// pick the greatest deviation.
const Eigen::Vector3d n_3 = n_1.cross(n_2).normalized();

auto OutOfPlaneDistance = [&n_1, &n_2, &origin](const Eigen::Vector3d &v)
{
return ((v - origin) - (v - origin).dot(n_1) * n_1 - (v - origin).dot(n_2) * n_2)
Expand Down Expand Up @@ -936,21 +934,21 @@ double LengthFromPointCloud(MPI_Comm comm, const std::vector<Eigen::Vector3d> &v

} // namespace

double GetDirectionalExtent(mfem::ParMesh &mesh, const mfem::Array<int> &marker, bool bdr,
const std::array<double, 3> &dir)
double GetProjectedLength(mfem::ParMesh &mesh, const mfem::Array<int> &marker, bool bdr,
const std::array<double, 3> &dir)
{
std::vector<Eigen::Vector3d> vertices;
int dominant_rank = CollectPointCloudOnRoot(mesh, marker, bdr, vertices);
return LengthFromPointCloud(mesh.GetComm(), vertices, dominant_rank, dir);
}

double GetDirectionalExtent(mfem::ParMesh &mesh, int attr, bool bdr,
const std::array<double, 3> &dir)
double GetProjectedLength(mfem::ParMesh &mesh, int attr, bool bdr,
const std::array<double, 3> &dir)
{
mfem::Array<int> marker(bdr ? mesh.bdr_attributes.Max() : mesh.attributes.Max());
marker = 0;
marker[attr - 1] = 1;
return GetDirectionalExtent(mesh, marker, bdr, dir);
return GetProjectedLength(mesh, marker, bdr, dir);
}

BoundingBox GetBoundingBox(mfem::ParMesh &mesh, const mfem::Array<int> &marker, bool bdr)
Expand Down
8 changes: 4 additions & 4 deletions palace/utils/geodata.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -113,10 +113,10 @@ BoundingBox GetBoundingBox(mfem::ParMesh &mesh, const mfem::Array<int> &marker,
BoundingBox GetBoundingBox(mfem::ParMesh &mesh, int attr, bool bdr);

// Helper function for computing the direction aligned length of a marked group.
double GetDirectionalExtent(mfem::ParMesh &mesh, const mfem::Array<int> &marker, bool bdr,
const std::array<double, 3> &dir);
double GetDirectionalExtent(mfem::ParMesh &mesh, int attr, bool bdr,
const std::array<double, 3> &dir);
double GetProjectedLength(mfem::ParMesh &mesh, const mfem::Array<int> &marker, bool bdr,
const std::array<double, 3> &dir);
double GetProjectedLength(mfem::ParMesh &mesh, int attr, bool bdr,
const std::array<double, 3> &dir);

// Given a mesh and a marker, compute the diameter of a bounding circle/sphere, assuming
// that the extrema points are in the marked group.
Expand Down

0 comments on commit 7485eb3

Please sign in to comment.