LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
write_matplotlib.cc
1
9#include "write_matplotlib.h"
10
11#include <fstream>
12#include <iostream>
13
14namespace lf::io {
15
16// NOLINTNEXTLINE(readability-function-cognitive-complexity)
17void writeMatplotlib(const lf::mesh::Mesh &mesh, std::string filename) {
18 using dim_t = lf::base::RefEl::dim_t;
19
20 // append .csv to filename if necessary
21 if (filename.compare(filename.size() - 4, 4, ".csv") != 0) {
22 filename += ".csv";
23 }
24
25 std::ofstream file(filename);
26
27 if (file.is_open()) {
28 const dim_t dim_mesh = mesh.DimMesh();
29 LF_VERIFY_MSG(dim_mesh == 2,
30 "write_matplotlib() only available for 2D meshes");
31
32 // loop through all elements of every codimension
33 for (int codim = 0; codim <= dim_mesh; ++codim) {
34 for (const lf::mesh::Entity *obj : mesh.Entities(codim)) {
35 const size_t obj_idx = mesh.Index(*obj);
36 const lf::base::RefEl obj_ref_el = obj->RefEl();
37 const lf::geometry::Geometry *obj_geo_ptr = obj->Geometry();
38 const Eigen::MatrixXd vertices =
39 obj_geo_ptr->Global(obj_ref_el.NodeCoords());
40
41 switch (obj_ref_el) {
43 file << codim << ',' << obj_idx << ',' << vertices(0, 0) << ','
44 << vertices(1, 0) << '\n';
45 break;
46 }
48 file << codim << ',' << obj_idx << ',';
49 // to access points of segment use SubEntities(1)
50 for (const auto *sub : obj->SubEntities(codim)) {
51 file << mesh.Index(*sub) << ',';
52 }
53 file << '\n';
54
55 break;
56 }
59 file << codim << ',' << obj_idx << ',';
60 // to access points of cell use SubEntities(1)
61 for (const auto *sub : obj->SubEntities(codim + 1)) {
62 file << mesh.Index(*sub) << ',';
63 }
64 file << '\n';
65
66 break;
67 }
68 default: {
69 LF_VERIFY_MSG(false, "Error for object " + std::to_string(obj_idx) +
70 " in codim " + std::to_string(codim) +
71 " of type " + obj_ref_el.ToString());
72 break;
73 }
74 }
75 }
76 }
77 }
78}
79
80} // namespace lf::io
Represents a reference element with all its properties.
Definition ref_el.h:109
static constexpr RefEl kSegment()
Returns the (1-dimensional) reference segment.
Definition ref_el.h:153
constexpr RefEl(RefElType type) noexcept
Create a RefEl from a lf::base::RefElType enum.
Definition ref_el.h:175
unsigned int dim_t
Definition ref_el.h:132
const Eigen::MatrixXd & NodeCoords() const
Get the coordinates of the nodes of this reference element.
Definition ref_el.h:241
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
Definition ref_el.h:144
static constexpr RefEl kTria()
Returns the reference triangle.
Definition ref_el.h:161
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
Definition ref_el.h:169
std::string ToString() const
Return a string representation of this Reference element.
Definition ref_el.h:458
Interface class for shape information on a mesh cell in the spirit of parametric finite element metho...
virtual Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const =0
Map a number of points in local coordinates into the global coordinate system.
Interface class representing a topological entity in a cellular complex
Definition entity.h:42
Abstract interface for objects representing a single mesh.
virtual std::span< const Entity *const > Entities(unsigned codim) const =0
All entities of a given codimension.
virtual unsigned DimMesh() const =0
The dimension of the manifold described by the mesh, or equivalently the maximum dimension of the ref...
virtual size_type Index(const Entity &e) const =0
Acess to the index of a mesh entity of any co-dimension.
Mesh input (from file) and output (in various formats) facilities.
void writeMatplotlib(const lf::mesh::Mesh &mesh, std::string filename)
Write affine triangulation data to file in matplotlib format.