diff --git a/CMakeLists.txt b/CMakeLists.txt index 99ff8b35..1d65cd14 100644 --- a/CMakeLists.txt +++ b/CMakeLists.txt @@ -172,6 +172,8 @@ meshfields_add_exe(OmegahTriTests test/testOmegahTri.cpp) meshfields_add_exe(ExceptionTest test/testExceptions.cpp) meshfields_add_exe(PointMapping test/testPointMapping.cpp) meshfields_add_exe(OmegahTetTest test/testOmegahTet.cpp) +meshfields_add_exe(IntegratorTest test/testIntegrator.cpp) +meshfields_add_exe(BoxIntegratorTest test/testBoxIntegrator.cpp) if(MeshFields_USE_Cabana) meshfields_add_exe(ControllerPerformance test/testControllerPerformance.cpp) @@ -190,6 +192,9 @@ test_func(CountIntegrator ./CountIntegrator) test_func(OmegahTriTests ./OmegahTriTests) test_func(PointMapping ./PointMapping) test_func(OmegahTetTest, ./OmegahTetTest) +test_func(IntegratorTest ./IntegratorTest) +test_func(BoxIntegratorTest ./BoxIntegratorTest) + if(MeshFields_USE_EXCEPTIONS) # exception caught - no error test_func(ExceptionTest ./ExceptionTest) diff --git a/docs/Ordering.pdf b/docs/Ordering.pdf new file mode 100644 index 00000000..941cf6a7 Binary files /dev/null and b/docs/Ordering.pdf differ diff --git a/docs/Ordering.tex b/docs/Ordering.tex new file mode 100644 index 00000000..bf6002b9 --- /dev/null +++ b/docs/Ordering.tex @@ -0,0 +1,76 @@ +\documentclass{article} +\usepackage{graphicx} +\usepackage[letterpaper,textwidth=5in,right=2.5in,textheight=9in]{geometry} + +\usepackage{amsmath,amssymb, changepage} +\usepackage{enumerate} +\usepackage{upgreek} +\pagestyle{empty} +\usepackage{listings} +\usepackage{tikz} +\begin{document} +\noindent +Joshua Kloepfer\\ +9/30/2025\\ +\textbf{Meshfields DOF Holder and Integration Point Ordering}\\\\ +Linear Triangle:\\ +\begin{tikzpicture} + \draw[gray, thick] (0, 0) -- (2.5, 2.5); + \draw[gray, thick] (0, 0) -- (5, 0); + \draw[gray, thick] (2.5, 2.5) -- (5, 0); + \filldraw[black] (0, 0) circle (2pt) node[anchor=east]{V0}; + \filldraw[black] (5, 0) circle (2pt) node[anchor=west]{V1}; + \filldraw[black] (2.5, 2.5) circle (2pt) node[anchor=south]{V2}; +\end{tikzpicture} +\\ +Quadratic Triangle\\ +\begin{tikzpicture} + \draw[gray, thick] (0, 0) -- (2.5, 2.5); + \draw[gray, thick] (0, 0) -- (5, 0); + \draw[gray, thick] (2.5, 2.5) -- (5, 0); + \filldraw[black] (0, 0) circle (2pt) node[anchor=east]{V0}; + \filldraw[black] (5, 0) circle (2pt) node[anchor=west]{V1}; + \filldraw[black] (2.5, 2.5) circle (2pt) node[anchor=south]{V2}; + \filldraw[black] (2.5, 0) circle (2pt) node[anchor=north]{E0}; + \filldraw[black] (3.75, 1.25) circle (2pt) node[anchor=west]{E1}; + \filldraw[black] (1.25, 1.25) circle (2pt) node[anchor=east]{E2}; +\end{tikzpicture} +\\ +Linear Tetrahedron:\\ +\begin{tikzpicture} + \draw[gray, thick] (0, 0, 0) -- (5, 0, 0); + \draw[gray, thick] (0, 0, 0) -- (0, 5, 0); + \draw[gray, thick] (5, 0, 0) -- (0, 5, 0); + \draw[gray, thick] (0, 0, 0) -- (0, 0, 5); + \draw[gray, thick] (5, 0, 0) -- (0, 0, 5); + \draw[gray, thick] (0, 5, 0) -- (0, 0, 5); + \filldraw[black] (0, 0, 5) circle (2pt) node[anchor=east]{V0}; + \filldraw[black] (5, 0, 0) circle (2pt) node[anchor=west]{V1}; + \filldraw[black] (0, 0, 0) circle (2pt) node[anchor=east]{V2}; + \filldraw[black] (0, 5, 0) circle (2pt) node[anchor=south]{V3}; +\end{tikzpicture} +\\ +Vertex ordering must follow right hand rule.\\ +Quadratic Tetrahedron:\\ +\begin{tikzpicture} + \draw[gray, thick] (0, 0, 0) -- (5, 0, 0); + \draw[gray, thick] (0, 0, 0) -- (0, 5, 0); + \draw[gray, thick] (5, 0, 0) -- (0, 5, 0); + \draw[gray, thick] (0, 0, 0) -- (0, 0, 5); + \draw[gray, thick] (5, 0, 0) -- (0, 0, 5); + \draw[gray, thick] (0, 5, 0) -- (0, 0, 5); + \filldraw[black] (0, 0, 5) circle (2pt) node[anchor=east]{V0}; + \filldraw[black] (5, 0, 0) circle (2pt) node[anchor=west]{V1}; + \filldraw[black] (0, 0, 0) circle (2pt) node[anchor=east]{V2}; + \filldraw[black] (0, 5, 0) circle (2pt) node[anchor=south]{V3}; + + \filldraw[black] (2.5, 0, 2.5) circle (2pt) node[anchor=north]{E0}; + \filldraw[black] (0, 0, 2.5) circle (2pt) node[anchor=east]{E1}; + \filldraw[black] (0, 2.5, 2.5) circle (2pt) node[anchor=east]{E2}; + \filldraw[black] (2.5, 0, 0) circle (2pt) node[anchor=south]{E3}; + \filldraw[black] (2.5, 2.5, 0) circle (2pt) node[anchor=east]{E4}; + \filldraw[black] (0, 2.5, 0) circle (2pt) node[anchor=west]{E5}; +\end{tikzpicture} +\\ +Vertex ordering must follow right hand rule and edge ordering must follow diagram. +\end{document} diff --git a/mfemTest/IntegratorMfemTest.cpp b/mfemTest/IntegratorMfemTest.cpp new file mode 100644 index 00000000..bd1c4950 --- /dev/null +++ b/mfemTest/IntegratorMfemTest.cpp @@ -0,0 +1,38 @@ +#include "mfem.hpp" +#include +#include + +using namespace std; +using namespace mfem; + +double intFunc(const Vector &x) { + return x[0]; +} + +int main(int argc, char *argv[]) { + + OptionsParser args(argc, argv); + int order = 1; + int size = 1; + int numberOfElements = 1; + int type = 1; + args.AddOption(&order, "-o", "--order", "Finite element polynomial degree"); + args.AddOption(&size, "-s", "--size", "Size of mesh"); + args.AddOption(&numberOfElements, "-n", "--numElem", "number of elements in each direction"); + args.AddOption(&type, "-t", "--type", "type of element(1 for tet, 0 for triangle"); + args.ParseCheck(); + + Mesh mesh = Mesh::MakeCartesian3D(numberOfElements, numberOfElements, type == 1 ? numberOfElements : 1, type == 1 ? Element::TETRAHEDRON : Element::TRIANGLE, size, size, type == 1 ? size : 1, false); + H1_FECollection fec(order, mesh.Dimension()); + FiniteElementSpace fespace(&mesh, &fec); + + FunctionCoefficient custom(intFunc); + + LinearForm b(&fespace); + + b.AddDomainIntegrator(new DomainLFIntegrator(custom)); + + b.Assemble(); + cout << b.Sum() << endl; + return 0; +} diff --git a/mfemTest/README.md b/mfemTest/README.md new file mode 100644 index 00000000..a745ad3c --- /dev/null +++ b/mfemTest/README.md @@ -0,0 +1,7 @@ +Run using g++ IntegratorMfemTest.cpp -I -L -l mfem +in order to run, use -m to specify the mesh file, -o to specify order, -s to specify size, -n to specify number of elements in each direction, and -t to specify type of element with 1 for tetrahedron and 0 for triangle. +mfem version: https://github.com/mfem/mfem/tree/2caa75e35a54df93d19f23655170254556dfc081 +For the build: +mkdir ; cd +cmake +make -j 4 diff --git a/src/MeshField_Element.hpp b/src/MeshField_Element.hpp index 2197a0a6..34d4756c 100644 --- a/src/MeshField_Element.hpp +++ b/src/MeshField_Element.hpp @@ -303,6 +303,30 @@ struct FieldElement { return Kokkos::View("foo", J.extent(0)); } + template + KOKKOS_INLINE_FUNCTION auto getGradients(Kokkos::View lc, + size_t pt) const { + if constexpr (ShapeType::Order == 1) { + return shapeFn.getLocalGradients(); + } else { + Kokkos::Array coord; + for (int i = 0; i < MeshEntDim + 1; ++i) { + coord[i] = lc(pt, i); + } + return shapeFn.getLocalGradients(coord); + } + } + + KOKKOS_INLINE_FUNCTION auto setNodalGradients( + Kokkos::View nodalGradients, + auto grad, size_t pt, size_t node, size_t d) const { + if constexpr (ShapeType::Order == 1) { + nodalGradients(pt, node, d) = grad[node * MeshEntDim + d]; + } else { + nodalGradients(pt, node, d) = grad[node][d]; + } + } + /** * @brief * Given an array of parametric coordinates 'localCoords', one per mesh @@ -380,22 +404,22 @@ struct FieldElement { // one matrix per point Kokkos::View res("result", numPts, MeshEntDim, MeshEntDim); Kokkos::deep_copy(res, 0.0); // initialize all entries to zero - // fill the views of node coordinates and node gradients Kokkos::View nodeCoords( "nodeCoords", numPts); Kokkos::View nodalGradients( "nodalGradients", numPts); - const auto grad = shapeFn.getLocalGradients(); Kokkos::parallel_for( numMeshEnts, KOKKOS_CLASS_LAMBDA(const int ent) { const auto vals = getNodeValues(ent); assert(vals.size() == MeshEntDim * ShapeType::numNodes); for (auto pt = offsets(ent); pt < offsets(ent + 1); pt++) { + auto ptCoords = Kokkos::subview(localCoords, pt, Kokkos::ALL()); + const auto grad = getGradients(localCoords, pt); for (size_t node = 0; node < ShapeType::numNodes; node++) { for (size_t d = 0; d < MeshEntDim; d++) { nodeCoords(pt, node, d) = vals[node * MeshEntDim + d]; - nodalGradients(pt, node, d) = grad[node * MeshEntDim + d]; + setNodalGradients(nodalGradients, grad, pt, node, d); } } } diff --git a/src/MeshField_Integrate.hpp b/src/MeshField_Integrate.hpp index 0c16d68c..849638b7 100644 --- a/src/MeshField_Integrate.hpp +++ b/src/MeshField_Integrate.hpp @@ -92,16 +92,16 @@ class TetrahedronIntegration : public EntityIntegration<4> { virtual std::vector> getPoints() const { return {IntegrationPoint(Vector4{0.138196601125011, 0.138196601125011, - 0.138196601125011, 0.585410196624969}, + 0.138196601125011, 0.585410196624967}, 0.25 / 6.0), - IntegrationPoint(Vector4{0.585410196624969, 0.138196601125011, + IntegrationPoint(Vector4{0.585410196624967, 0.138196601125011, 0.138196601125011, 0.138196601125011}, 0.25 / 6.0), - IntegrationPoint(Vector4{0.138196601125011, 0.585410196624969, + IntegrationPoint(Vector4{0.138196601125011, 0.585410196624967, 0.138196601125011, 0.138196601125011}, 0.25 / 6.0), IntegrationPoint(Vector4{0.138196601125011, 0.138196601125011, - 0.585410196624969, 0.138196601125011}, + 0.585410196624967, 0.138196601125011}, 0.25 / 6.0)}; } virtual int getAccuracy() const { return 2; } diff --git a/test/testBoxIntegrator.cpp b/test/testBoxIntegrator.cpp new file mode 100644 index 00000000..fa05b4ac --- /dev/null +++ b/test/testBoxIntegrator.cpp @@ -0,0 +1,119 @@ +#include "KokkosController.hpp" +#include "MeshField.hpp" +#include "MeshField_Element.hpp" //remove? +#include "MeshField_Fail.hpp" //remove? +#include "MeshField_For.hpp" //remove? +#include "MeshField_ShapeField.hpp" //remove? +#include "Omega_h_build.hpp" +#include "Omega_h_file.hpp" +#include "Omega_h_simplex.hpp" +#include +#include +#include +#include + +using ExecutionSpace = Kokkos::DefaultExecutionSpace; +using MemorySpace = Kokkos::DefaultExecutionSpace::memory_space; + +template +class testIntegrator : public MeshField::Integrator { +private: + testIntegrator(){}; + +protected: + MeshField::Real integral; + FieldElement &fes; + +public: + MeshField::Real getIntegral() { return integral; } + testIntegrator(FieldElement &fes_in, int order) + : Integrator(order), integral(0), fes(fes_in){}; + void atPoints(Kokkos::View p, + Kokkos::View w, + Kokkos::View dV) { + const auto globalCoords = MeshField::evaluate(fes, p, p.extent(0) / fes.numMeshEnts); + Kokkos::parallel_reduce( + "integrate", p.extent(0), + KOKKOS_LAMBDA(const int &ent, MeshField::Real &integ) { + const auto x = globalCoords(ent, 0); + const auto y = globalCoords(ent, 1); + const auto z = globalCoords.extent(1) == 3 ? globalCoords(ent, 2) : 1; + integ += w(ent) * x * dV(ent); + }, + integral); + } +}; + +template