LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
write_matplotlib.cc
1
9#include "write_matplotlib.h"
10
11#include <fstream>
12#include <iostream>
13
14namespace lf::io {
15
16void writeMatplotlib(const lf::mesh::Mesh &mesh, std::string filename) {
18
19 // append .csv to filename if necessary
20 if (filename.compare(filename.size() - 4, 4, ".csv") != 0) {
21 filename += ".csv";
22 }
23
24 std::ofstream file(filename);
25
26 if (file.is_open()) {
27 const dim_t dim_mesh = mesh.DimMesh();
28 LF_VERIFY_MSG(dim_mesh == 2,
29 "write_matplotlib() only available for 2D meshes");
30
31 // loop through all elements of every codimension
32 for (int codim = 0; codim <= dim_mesh; ++codim) {
33 for (const lf::mesh::Entity *obj : mesh.Entities(codim)) {
34 const size_t obj_idx = mesh.Index(*obj);
35 const lf::base::RefEl obj_ref_el = obj->RefEl();
36 const lf::geometry::Geometry *obj_geo_ptr = obj->Geometry();
37 const Eigen::MatrixXd vertices =
38 obj_geo_ptr->Global(obj_ref_el.NodeCoords());
39
40 switch (obj_ref_el) {
42 file << codim << ',' << obj_idx << ',' << vertices(0, 0) << ','
43 << vertices(1, 0) << std::endl;
44 break;
45 }
47 file << codim << ',' << obj_idx << ',';
48 // to access points of segment use SubEntities(1)
49 for (const auto *sub : obj->SubEntities(codim)) {
50 file << mesh.Index(*sub) << ',';
51 }
52 file << std::endl;
53
54 break;
55 }
58 file << codim << ',' << obj_idx << ',';
59 // to access points of cell use SubEntities(1)
60 for (const auto *sub : obj->SubEntities(codim + 1)) {
61 file << mesh.Index(*sub) << ',';
62 }
63 file << std::endl;
64
65 break;
66 }
67 default: {
68 LF_VERIFY_MSG(false, "Error for object " + std::to_string(obj_idx) +
69 " in codim " + std::to_string(codim) +
70 " of type " + obj_ref_el.ToString());
71 break;
72 }
73 }
74 }
75 }
76 }
77}
78
79} // namespace lf::io
Represents a reference element with all its properties.
Definition: ref_el.h:106
static constexpr RefEl kSegment()
Returns the (1-dimensional) reference segment.
Definition: ref_el.h:150
constexpr RefEl(RefElType type) noexcept
Create a RefEl from a lf::base::RefElType enum.
Definition: ref_el.h:172
unsigned int dim_t
Definition: ref_el.h:129
const Eigen::MatrixXd & NodeCoords() const
Get the coordinates of the nodes of this reference element.
Definition: ref_el.h:238
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
Definition: ref_el.h:141
static constexpr RefEl kTria()
Returns the reference triangle.
Definition: ref_el.h:158
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
Definition: ref_el.h:166
std::string ToString() const
Return a string representation of this Reference element.
Definition: ref_el.h:455
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:39
Abstract interface for objects representing a single mesh.
virtual unsigned DimMesh() const =0
The dimension of the manifold described by the mesh, or equivalently the maximum dimension of the ref...
virtual nonstd::span< const Entity *const > Entities(unsigned codim) const =0
All entities of a given codimension.
virtual size_type Index(const Entity &e) const =0
Acess to the index of a mesh entity of any co-dimension.
unsigned int dim_t
type for dimensions and co-dimensions and numbers derived from them
Definition: base.h:36
Mesh input (from file) and output (in various formats) facilities.
Definition: gmsh_file_v2.cc:35
void writeMatplotlib(const lf::mesh::Mesh &mesh, std::string filename)
Write affine triangulation data to file in matplotlib format.