// $Id$
// Version: $Name$
//
-// Copyright (C) 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009 by the deal.II authors
+// Copyright (C) 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
{
Assert (tria != 0, ExcNoTriangulationSelected());
AssertThrow (in, ExcIO());
-
+
// skip comments at start of file
skip_comment_lines (in, '#');
// in ucd-file (key) and in the
// vertices vector
std::map<int,int> vertex_indices;
-
- for (unsigned int vertex=0; vertex<n_vertices; ++vertex)
+
+ for (unsigned int vertex=0; vertex<n_vertices; ++vertex)
{
int vertex_number;
double x[3];
std::vector<CellData<dim> > cells;
SubCellData subcelldata;
- for (unsigned int cell=0; cell<n_cells; ++cell)
+ for (unsigned int cell=0; cell<n_cells; ++cell)
{
// note that since in the input
// file we found the number of
// should still be input here,
// so check this:
AssertThrow (in, ExcIO());
-
+
std::string cell_type;
int material_id;
-
+
in >> dummy // cell number
>> material_id;
in >> cell_type;
if (vertex_indices.find (cells.back().vertices[i]) != vertex_indices.end())
// vertex with this index exists
cells.back().vertices[i] = vertex_indices[cells.back().vertices[i]];
- else
+ else
{
// no such vertex index
AssertThrow (false, ExcInvalidVertexIndex(cell, cells.back().vertices[i]));
// vertex with this index exists
subcelldata.boundary_lines.back().vertices[i]
= vertex_indices[subcelldata.boundary_lines.back().vertices[i]];
- else
+ else
{
// no such vertex index
AssertThrow (false,
>> subcelldata.boundary_quads.back().vertices[1]
>> subcelldata.boundary_quads.back().vertices[2]
>> subcelldata.boundary_quads.back().vertices[3];
-
+
subcelldata.boundary_quads.back().material_id = material_id;
-
+
// transform from ucd to
// consecutive numbering
for (unsigned int i=0; i<4; ++i)
subcelldata.boundary_quads.back().vertices[i] =
numbers::invalid_unsigned_int;
};
-
+
}
else
// cannot read this
AssertThrow (false, ExcUnknownIdentifier(cell_type));
};
-
+
// check that no forbidden arrays are used
Assert (subcelldata.check_consistency(dim), ExcInternalError());
Assert (dim==2, ExcNotImplemented());
AssertThrow (in, ExcIO());
-
+
// skip comments at start of file
skip_comment_lines (in, '#');
// now read vertices
getline (in, line);
AssertThrow (line=="Vertices", ExcInvalidDBMESHInput(line));
-
+
unsigned int n_vertices;
double dummy;
-
+
in >> n_vertices;
std::vector<Point<spacedim> > vertices (n_vertices);
for (unsigned int vertex=0; vertex<n_vertices; ++vertex)
// discard the input
getline (in, line);
AssertThrow (line=="Edges", ExcInvalidDBMESHInput(line));
-
+
unsigned int n_edges;
in >> n_edges;
for (unsigned int edge=0; edge<n_edges; ++edge)
// discard the input
getline (in, line);
AssertThrow (line=="CrackedEdges", ExcInvalidDBMESHInput(line));
-
+
in >> n_edges;
for (unsigned int edge=0; edge<n_edges; ++edge)
{
SubCellData subcelldata;
unsigned int n_cells;
in >> n_cells;
- for (unsigned int cell=0; cell<n_cells; ++cell)
+ for (unsigned int cell=0; cell<n_cells; ++cell)
{
// read in vertex numbers. they
// are 1-based, so subtract one
for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
{
in >> cells.back().vertices[i];
-
+
AssertThrow ((cells.back().vertices[i] >= 1)
&&
(static_cast<unsigned int>(cells.back().vertices[i]) <= vertices.size()),
ExcInvalidVertexIndex(cell, cells.back().vertices[i]));
-
+
--cells.back().vertices[i];
};
AssertThrow (in, ExcInvalidDBMeshFormat());
skip_empty_lines(in);
-
+
// then there are again a whole lot
// of fields of which I have no
// ok, so we are not at the end of
// the file, that's it, mostly
-
+
// check that no forbidden arrays are used
Assert (subcelldata.check_consistency(dim), ExcInternalError());
std::vector<CellData<2> > cells (n_cells);
SubCellData subcelldata;
- for (unsigned int cell=0; cell<n_cells; ++cell)
+ for (unsigned int cell=0; cell<n_cells; ++cell)
{
// note that since in the input
// file we found the number of
AssertThrow (in, ExcIO());
Assert (GeometryInfo<2>::vertices_per_cell == 4,
ExcInternalError());
-
+
for (unsigned int i=0; i<4; ++i)
in >> cells[cell].vertices[i];
};
-
+
// set up array of vertices
std::vector<Point<2> > vertices (n_vertices);
- for (unsigned int vertex=0; vertex<n_vertices; ++vertex)
+ for (unsigned int vertex=0; vertex<n_vertices; ++vertex)
{
double x[3];
AssertThrow (in, ExcIO());
static const unsigned int xda_to_dealII_map[] = {0,1,5,4,3,2,6,7};
-
+
std::string line;
// skip comments at start of file
getline (in, line);
std::vector<CellData<3> > cells (n_cells);
SubCellData subcelldata;
- for (unsigned int cell=0; cell<n_cells; ++cell)
+ for (unsigned int cell=0; cell<n_cells; ++cell)
{
// note that since in the input
// file we found the number of
AssertThrow (in, ExcIO());
Assert(GeometryInfo<3>::vertices_per_cell == 8,
ExcInternalError());
-
+
unsigned int xda_ordered_nodes[8];
-
+
for (unsigned int i=0; i<8; ++i)
in >> xda_ordered_nodes[i];
};
-
+
// set up array of vertices
std::vector<Point<3> > vertices (n_vertices);
- for (unsigned int vertex=0; vertex<n_vertices; ++vertex)
+ for (unsigned int vertex=0; vertex<n_vertices; ++vertex)
{
double x[3];
unsigned int n_cells;
unsigned int dummy;
std::string line;
-
+
in >> line;
// first determine file format
ExcInvalidGMSHInput(line));
in >> line;
+ // if the next block is of kind
+ // $PhysicalNames, ignore it
+ if (line == "$PhysicalNames")
+ {
+ do
+ {
+ in >> line;
+ }
+ while (line != "$EndPhysicalNames");
+ in >> line;
+ }
+
+ // but the next thing should,
+ // in any case, be the list of
+ // nodes:
AssertThrow (line == "$Nodes",
ExcInvalidGMSHInput(line));
}
-
+
// now read the nodes list
in >> n_vertices;
std::vector<Point<spacedim> > vertices (n_vertices);
// in msh-file (nod) and in the
// vertices vector
std::map<int,int> vertex_indices;
-
- for (unsigned int vertex=0; vertex<n_vertices; ++vertex)
+
+ for (unsigned int vertex=0; vertex<n_vertices; ++vertex)
{
int vertex_number;
double x[3];
// read vertex
in >> vertex_number
>> x[0] >> x[1] >> x[2];
-
+
for (unsigned int d=0; d<spacedim; ++d)
vertices[vertex](d) = x[d];
// store mapping
static const std::string end_nodes_marker[] = {"$ENDNOD", "$EndNodes" };
AssertThrow (line==end_nodes_marker[gmsh_file_format-1],
ExcInvalidGMSHInput(line));
-
+
// Now read in next bit
in >> line;
static const std::string begin_elements_marker[] = {"$ELM", "$Elements" };
ExcInvalidGMSHInput(line));
in >> n_cells;
-
+
// set up array of cells
std::vector<CellData<dim> > cells;
SubCellData subcelldata;
- for (unsigned int cell=0; cell<n_cells; ++cell)
+ for (unsigned int cell=0; cell<n_cells; ++cell)
{
// note that since in the input
// file we found the number of
// should still be input here,
// so check this:
AssertThrow (in, ExcIO());
-
+
unsigned int cell_type;
unsigned int material_id;
unsigned int nod_num;
(version 1) or the first tag (version 2, if any tag is given at all) as
material id.
*/
-
+
in >> dummy // ELM-NUMBER
>> cell_type; // ELM-TYPE
>> nod_num;
break;
}
-
+
case 2:
{
// read the tags; ignore
in >> material_id;
else
material_id = 0;
-
+
for (unsigned int i=1; i<n_tags; ++i)
in >> dummy;
default:
AssertThrow (false, ExcNotImplemented());
}
-
-
+
+
/* `ELM-TYPE'
defines the geometrical type of the N-th element:
`1'
Line (2 nodes, 1 edge).
-
+
`3'
Quadrangle (4 nodes, 4 edges).
-
+
`5'
Hexahedron (8 nodes, 12 edges, 6 faces).
-
+
`15'
Point (1 node).
*/
-
+
if (((cell_type == 1) && (dim == 1)) ||
((cell_type == 3) && (dim == 2)) ||
((cell_type == 5) && (dim == 3)))
AssertThrow (nod_num == GeometryInfo<dim>::vertices_per_cell,
ExcMessage ("Number of nodes does not coincide with the "
"number required for this object"));
-
+
// allocate and read indices
cells.push_back (CellData<dim>());
for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
AssertThrow (vertex_indices.find (cells.back().vertices[i]) !=
vertex_indices.end(),
ExcInvalidVertexIndex(cell, cells.back().vertices[i]));
-
+
// vertex with this index exists
cells.back().vertices[i] = vertex_indices[cells.back().vertices[i]];
}
// vertex with this index exists
subcelldata.boundary_lines.back().vertices[i]
= vertex_indices[subcelldata.boundary_lines.back().vertices[i]];
- else
+ else
{
// no such vertex index
AssertThrow (false,
>> subcelldata.boundary_quads.back().vertices[1]
>> subcelldata.boundary_quads.back().vertices[2]
>> subcelldata.boundary_quads.back().vertices[3];
-
+
subcelldata.boundary_quads.back().material_id = material_id;
-
+
// transform from gmsh to
// consecutive numbering
for (unsigned int i=0; i<4; ++i)
subcelldata.boundary_quads.back().vertices[i] =
numbers::invalid_unsigned_int;
};
-
+
}
else
if (cell_type == 15)
static const std::string end_elements_marker[] = {"$ENDELM", "$EndElements" };
AssertThrow (line==end_elements_marker[gmsh_file_format-1],
ExcInvalidGMSHInput(line));
-
+
// check that no forbidden arrays are used
Assert (subcelldata.check_consistency(dim), ExcInternalError());
// check that we actually read some
// cells.
AssertThrow(cells.size() > 0, ExcGmshNoCellInformation());
-
+
// do some clean-up on
// vertices...
GridTools::delete_unused_vertices (vertices, cells, subcelldata);
// to find the right value for
// coord, and setting x2d and y2d
// accordingly.
-
+
// First, open the file
NcFile nc (filename.c_str());
AssertThrow(nc.is_valid(), ExcIO());
// use &* to convert
// vector<int>::iterator to int *
marker_var->get(&*marker.begin(), n_markers);
-
+
if (output)
{
std::cout << "n_cell=" << n_cells << std::endl;
// different
AssertThrow(n_bquads_per_bmarker.find(marker[i])==
n_bquads_per_bmarker.end(), ExcIO());
-
+
n_bquads_per_bmarker[marker[i]]=
count(bmarker.begin(), bmarker.end(), marker[i]);
}
AssertThrow(quad_vertices_dim->is_valid(), ExcIO());
const unsigned int vertices_per_quad=quad_vertices_dim->size();
AssertThrow(vertices_per_quad==GeometryInfo<dim>::vertices_per_cell, ExcIO());
-
+
NcVar *vertex_indices_var=nc.get_var("points_of_surfacequadrilaterals");
AssertThrow(vertex_indices_var->is_valid(), ExcIO());
AssertThrow(vertex_indices_var->num_dims()==2, ExcIO());
zero_plane=false;
break;
}
-
+
if (zero_plane)
zero_plane_markers[bmarker[quad]]=true;
}
zero_plane=true;
break;
}
-
+
if (zero_plane)
{
for (unsigned int i=0; i<vertices_per_quad; ++i)
SubCellData subcelldata;
GridTools::delete_unused_vertices(vertices, cells, subcelldata);
GridReordering<dim,spacedim>::reorder_cells (cells);
- tria->create_triangulation_compatibility (vertices, cells, subcelldata);
+ tria->create_triangulation_compatibility (vertices, cells, subcelldata);
#endif
}
Assert (tria != 0, ExcNoTriangulationSelected());
// this function assumes the TAU
// grid format.
-
+
// First, open the file
NcFile nc (filename.c_str());
AssertThrow(nc.is_valid(), ExcIO());
// then read n_cells
NcDim *elements_dim=nc.get_dim("no_of_elements");
AssertThrow(elements_dim->is_valid(), ExcIO());
- const unsigned int n_cells=elements_dim->size();
+ const unsigned int n_cells=elements_dim->size();
if (output)
std::cout << "n_cell=" << n_cells << std::endl;
// and n_hexes
const unsigned int n_hexes=hexes_dim->size();
AssertThrow(n_hexes==n_cells,
ExcMessage("deal.II can handle purely hexaedral grids, only."));
-
+
// next we read
// int points_of_hexaeders(
// no_of_hexaeders,
AssertThrow(hex_vertices_dim->is_valid(), ExcIO());
const unsigned int vertices_per_hex=hex_vertices_dim->size();
AssertThrow(vertices_per_hex==GeometryInfo<dim>::vertices_per_cell, ExcIO());
-
+
NcVar *vertex_indices_var=nc.get_var("points_of_hexaeders");
AssertThrow(vertex_indices_var->is_valid(), ExcIO());
AssertThrow(vertex_indices_var->num_dims()==2, ExcIO());
// we read the vertex indices of
// the boundary quadrilaterals and
// their boundary markers
-
+
// first we read
// int points_of_surfacequadrilaterals(
// no_of_surfacequadrilaterals,
AssertThrow(quad_vertices_dim->is_valid(), ExcIO());
const unsigned int vertices_per_quad=quad_vertices_dim->size();
AssertThrow(vertices_per_quad==GeometryInfo<dim>::vertices_per_face, ExcIO());
-
+
NcVar *bvertex_indices_var=nc.get_var("points_of_surfacequadrilaterals");
AssertThrow(bvertex_indices_var->is_valid(), ExcIO());
AssertThrow(bvertex_indices_var->num_dims()==2, ExcIO());
GridTools::delete_unused_vertices(vertices, cells, subcelldata);
GridReordering<dim,spacedim>::invert_all_cells_of_negative_grid (vertices, cells);
GridReordering<dim,spacedim>::reorder_cells (cells);
- tria->create_triangulation_compatibility (vertices, cells, subcelldata);
+ tria->create_triangulation_compatibility (vertices, cells, subcelldata);
#endif
}
}
structured=true;
blocked=false;
-
+
// convert the string to upper case
std::transform(header.begin(),header.end(),header.begin(),::toupper);
-
+
// replace all tabs, commas, newlines by
// whitespaces
std::replace(header.begin(),header.end(),'\t',' ');
std::replace(header.begin(),header.end(),',',' ');
std::replace(header.begin(),header.end(),'\n',' ');
-
+
// now remove whitespace in front of and
// after '='
std::string::size_type pos=header.find("=");
-
+
while(pos!=static_cast<std::string::size_type>(std::string::npos))
if(header[pos+1]==' ')
header.erase(pos+1,1);
}
else
pos=header.find("=",++pos);
-
+
// split the string into individual entries
std::vector<std::string> entries=Utilities::break_text_into_lines(header,1,' ');
// set i back, so that the next
// string is treated correctly
--i;
-
+
AssertThrow(n_vars>=dim,
ExcMessage("Tecplot file must contain at least one variable for each dimension"));
for (unsigned int i=1; i<dim; ++i)
else if (Utilities::match_at_string_start(entries[i],"E="))
n_cells=Utilities::get_integer_at_position(entries[i],2).first;
}
-
+
// now we have read all the fields we are
// interested in. do some checks and
// calculate the variables
const unsigned int spacedim=2;
Assert (tria != 0, ExcNoTriangulationSelected());
AssertThrow (in, ExcIO());
-
+
// skip comments at start of file
skip_comment_lines (in, '#');
// some strings for parsing the header
std::string line, header;
-
+
// first, concatenate all header lines
// create a searchstring with almost all
// letters. exclude e and E from the letters
header+=" "+line;
getline(in,line);
}
-
+
// now create some variables holding
// important information on the mesh, get
// this information from the header string
parse_tecplot_header(header,
tecplot2deal,n_vars,n_vertices,n_cells,IJK,
structured,blocked);
-
+
// reserve space for vertices. note, that in
// tecplot vertices are ordered beginning
// with 1, whereas in deal all indices start
// which is the first index to read in
// the loop (see below)
unsigned int next_index=0;
-
+
// note, that we have already read the
// first line containing the first variable
if (tecplot2deal[0]==0)
// variables)
if (next_index==dim && structured)
break;
-
+
if (i==tecplot2deal[next_index])
{
// we need this line, read it in
// next and so on. create a vector to
// hold these components
std::vector<double> vars(n_vars);
-
+
// now fill the first vertex. note, that we
// have already read the first line
// containing the first vertex
char *endptr;
for (unsigned int d=0; d<dim; ++d)
vertices[1](d) = std::strtod (first_vertex[tecplot2deal[d]].c_str(), &endptr);
-
+
// read the remaining vertices from the
// list
for (unsigned int v=2; v<n_vertices+1; ++v)
vertices[v](i)=vars[tecplot2deal[i]];
}
}
-
+
if (structured)
{
// this is the part of the code that only
// works in 2d
unsigned int I=IJK[0],
J=IJK[1];
-
+
unsigned int cell=0;
// set up array of cells
- for (unsigned int j=0; j<J-1; ++j)
+ for (unsigned int j=0; j<J-1; ++j)
for (unsigned int i=1; i<I; ++i)
{
cells[cell].vertices[0]=i+ j *I;
// set up array of cells, unstructured
// mode, so the connectivity is
// explicitly given
- for (unsigned int i=0; i<n_cells; ++i)
+ for (unsigned int i=0; i<n_cells; ++i)
{
// note that since in the input file
// we found the number of cells at
// the top, there should still be
// input here, so check this:
AssertThrow (in, ExcIO());
-
+
// get the connectivity from the
// input file. the vertices are
// ordered like in the ucd format
// do some clean-up on vertices
GridTools::delete_unused_vertices (vertices, cells, subcelldata);
}
-
+
// check that no forbidden arrays are
// used. as we do not read in any
// subcelldata, nothing should happen here.
template <int dim, int spacedim>
void GridIn<dim, spacedim>::skip_empty_lines (std::istream &in)
-{
+{
std::string line;
- while (in)
+ while (in)
{
// get line
getline (in, line);
// the comment starter
while (in.get() != '\n')
;
-
-
+
+
// put back first character of
// first non-comment line
in.putback (c);
max_x = vertices[cells[0].vertices[0]](0),
min_y = vertices[cells[0].vertices[0]](1),
max_y = vertices[cells[0].vertices[0]](1);
-
+
for (unsigned int i=0; i<cells.size(); ++i)
{
for (unsigned int v=0; v<4; ++v)
{
const Point<2> &p = vertices[cells[i].vertices[v]];
-
+
if (p(0) < min_x)
min_x = p(0);
if (p(0) > max_x)
// first two line right direction
for (unsigned int f=0; f<2; ++f)
out << "set arrow from "
- << vertices[cells[i].vertices[f]](0) << ','
+ << vertices[cells[i].vertices[f]](0) << ','
<< vertices[cells[i].vertices[f]](1)
<< " to "
- << vertices[cells[i].vertices[(f+1)%4]](0) << ','
+ << vertices[cells[i].vertices[(f+1)%4]](0) << ','
<< vertices[cells[i].vertices[(f+1)%4]](1)
<< std::endl;
// other two lines reverse direction
for (unsigned int f=2; f<4; ++f)
out << "set arrow from "
- << vertices[cells[i].vertices[(f+1)%4]](0) << ','
+ << vertices[cells[i].vertices[(f+1)%4]](0) << ','
<< vertices[cells[i].vertices[(f+1)%4]](1)
<< " to "
- << vertices[cells[i].vertices[f]](0) << ','
+ << vertices[cells[i].vertices[f]](0) << ','
<< vertices[cells[i].vertices[f]](1)
<< std::endl;
out << std::endl;
};
-
+
out << std::endl
<< "set nokey" << std::endl
name = search.find(filename, default_suffix(format));
std::ifstream in(name.c_str());
-
+
if (format == Default)
{
const std::string::size_type slashpos = name.find_last_of('/');
{
if (format == Default)
format = default_format;
-
+
switch (format)
{
case dbmesh:
read_dbmesh (in);
return;
-
+
case msh:
read_msh (in);
return;
-
+
case ucd:
read_ucd (in);
return;
-
+
case xda:
read_xda (in);
return;
case tecplot:
read_tecplot (in);
return;
-
+
case Default:
break;
}
template <int dim, int spacedim>
std::string
-GridIn<dim, spacedim>::default_suffix (const Format format)
+GridIn<dim, spacedim>::default_suffix (const Format format)
{
- switch (format)
+ switch (format)
{
case dbmesh:
return ".dbmesh";
- case msh:
+ case msh:
return ".msh";
- case ucd:
+ case ucd:
return ".inp";
case xda:
return ".xda";
return ".nc";
case tecplot:
return ".dat";
- default:
- Assert (false, ExcNotImplemented());
+ default:
+ Assert (false, ExcNotImplemented());
return ".unknown_format";
}
}
if (format_name == "msh")
return msh;
-
+
if (format_name == "inp")
return ucd;
template <int dim, int spacedim>
-std::string GridIn<dim, spacedim>::get_format_names ()
+std::string GridIn<dim, spacedim>::get_format_names ()
{
return "dbmesh|msh|ucd|xda|netcdf|tecplot";
}
--- /dev/null
+$MeshFormat
+2 0 8
+$EndMeshFormat
+$PhysicalNames
+2
+1 1 "1"
+2 2 "This surface"
+$EndPhysicalNames
+$Nodes
+100
+1 0 0 0
+2 20 0 0
+3 20 10 0
+4 0 10 0
+5 2.222222222217813 0 0
+6 4.444444444435626 0 0
+7 6.666666666655381 0 0
+8 8.888888888875783 0 0
+9 11.11111111109865 0 0
+10 13.33333333332399 0 0
+11 15.55555555554933 0 0
+12 17.77777777777466 0 0
+13 20 1.111111111108906 0
+14 20 2.222222222217813 0
+15 20 3.33333333332769 0
+16 20 4.444444444437892 0
+17 20 5.555555555549327 0
+18 20 6.666666666661994 0
+19 20 7.777777777774663 0
+20 20 8.888888888887331 0
+21 17.77777777777839 10 0
+22 15.55555555555679 10 0
+23 13.33333333333518 10 0
+24 11.11111111111358 10 0
+25 8.888888888891973 10 0
+26 6.666666666670369 10 0
+27 4.444444444448146 10 0
+28 2.222222222224073 10 0
+29 0 8.888888888889197 0
+30 0 7.777777777778395 0
+31 0 6.666666666667592 0
+32 0 5.555555555556789 0
+33 0 4.444444444445987 0
+34 0 3.333333333335184 0
+35 0 2.222222222224073 0
+36 0 1.111111111112036 0
+37 2.222222222218642 1.111111111111496 0
+38 2.222222222219417 2.222222222222914 0
+39 2.222222222220129 3.333333333334147 0
+40 2.222222222220792 4.444444444445237 0
+41 2.222222222221424 5.555555555556247 0
+42 2.222222222222049 6.666666666667213 0
+43 2.222222222222688 7.777777777778153 0
+44 2.222222222223361 8.888888888889079 0
+45 4.444444444437504 1.111111111111264 0
+46 4.444444444439081 2.222222222222499 0
+47 4.444444444440458 3.333333333333679 0
+48 4.444444444441706 4.444444444444804 0
+49 4.444444444442881 5.555555555555887 0
+50 4.444444444444038 6.666666666666938 0
+51 4.44444444444524 7.777777777777969 0
+52 4.444444444446567 8.888888888888987 0
+53 6.666666666657572 1.111111111111154 0
+54 6.666666666659543 2.2222222222223 0
+55 6.66666666666131 3.333333333333432 0
+56 6.666666666662916 4.444444444444548 0
+57 6.666666666664412 5.555555555555651 0
+58 6.66666666666585 6.666666666666745 0
+59 6.666666666667285 7.777777777777832 0
+60 6.666666666668773 8.888888888888916 0
+61 8.888888888878521 1.111111111111078 0
+62 8.888888888880823 2.22222222222216 0
+63 8.888888888882811 3.333333333333248 0
+64 8.888888888884571 4.444444444444345 0
+65 8.888888888886175 5.555555555555455 0
+66 8.888888888887676 6.666666666666577 0
+67 8.888888888889122 7.77777777777771 0
+68 8.888888888890547 8.888888888888852 0
+69 11.11111111110122 1.111111111110989 0
+70 11.11111111110336 2.222222222221994 0
+71 11.11111111110522 3.333333333333026 0
+72 11.11111111110686 4.444444444444096 0
+73 11.11111111110835 5.55555555555521 0
+74 11.11111111110973 6.666666666666364 0
+75 11.11111111111104 7.777777777777554 0
+76 11.11111111111232 8.88888888888877 0
+77 13.33333333332558 1.11111111111084 0
+78 13.33333333332708 2.222222222221715 0
+79 13.33333333332847 3.333333333332654 0
+80 13.33333333332975 4.444444444443679 0
+81 13.33333333333094 5.555555555554799 0
+82 13.33333333333206 6.666666666666009 0
+83 13.33333333333313 7.777777777777295 0
+84 13.33333333333417 8.888888888888633 0
+85 15.5555555555503 1.111111111110555 0
+86 15.55555555555125 2.22222222222118 0
+87 15.55555555555216 3.333333333331943 0
+88 15.55555555555302 4.444444444442883 0
+89 15.55555555555383 5.555555555554025 0
+90 15.55555555555461 6.666666666665354 0
+91 15.55555555555535 7.777777777776826 0
+92 15.55555555555608 8.88888888888839 0
+93 17.77777777777513 1.111111111109987 0
+94 17.77777777777559 2.222222222220101 0
+95 17.77777777777604 3.333333333330518 0
+96 17.77777777777647 4.444444444441285 0
+97 17.77777777777688 5.555555555552501 0
+98 17.77777777777728 6.66666666666413 0
+99 17.77777777777766 7.777777777775988 0
+100 17.77777777777803 8.888888888887967 0
+$EndNodes
+$Elements
+90
+1 1 3 1 2 0 2 13
+2 1 3 1 2 0 13 14
+3 1 3 1 2 0 14 15
+4 1 3 1 2 0 15 16
+5 1 3 1 2 0 16 17
+6 1 3 1 2 0 17 18
+7 1 3 1 2 0 18 19
+8 1 3 1 2 0 19 20
+9 1 3 1 2 0 20 3
+10 3 3 2 1 0 1 5 37 36
+11 3 3 2 1 0 36 37 38 35
+12 3 3 2 1 0 35 38 39 34
+13 3 3 2 1 0 34 39 40 33
+14 3 3 2 1 0 33 40 41 32
+15 3 3 2 1 0 32 41 42 31
+16 3 3 2 1 0 31 42 43 30
+17 3 3 2 1 0 30 43 44 29
+18 3 3 2 1 0 29 44 28 4
+19 3 3 2 1 0 5 6 45 37
+20 3 3 2 1 0 37 45 46 38
+21 3 3 2 1 0 38 46 47 39
+22 3 3 2 1 0 39 47 48 40
+23 3 3 2 1 0 40 48 49 41
+24 3 3 2 1 0 41 49 50 42
+25 3 3 2 1 0 42 50 51 43
+26 3 3 2 1 0 43 51 52 44
+27 3 3 2 1 0 44 52 27 28
+28 3 3 2 1 0 6 7 53 45
+29 3 3 2 1 0 45 53 54 46
+30 3 3 2 1 0 46 54 55 47
+31 3 3 2 1 0 47 55 56 48
+32 3 3 2 1 0 48 56 57 49
+33 3 3 2 1 0 49 57 58 50
+34 3 3 2 1 0 50 58 59 51
+35 3 3 2 1 0 51 59 60 52
+36 3 3 2 1 0 52 60 26 27
+37 3 3 2 1 0 7 8 61 53
+38 3 3 2 1 0 53 61 62 54
+39 3 3 2 1 0 54 62 63 55
+40 3 3 2 1 0 55 63 64 56
+41 3 3 2 1 0 56 64 65 57
+42 3 3 2 1 0 57 65 66 58
+43 3 3 2 1 0 58 66 67 59
+44 3 3 2 1 0 59 67 68 60
+45 3 3 2 1 0 60 68 25 26
+46 3 3 2 1 0 8 9 69 61
+47 3 3 2 1 0 61 69 70 62
+48 3 3 2 1 0 62 70 71 63
+49 3 3 2 1 0 63 71 72 64
+50 3 3 2 1 0 64 72 73 65
+51 3 3 2 1 0 65 73 74 66
+52 3 3 2 1 0 66 74 75 67
+53 3 3 2 1 0 67 75 76 68
+54 3 3 2 1 0 68 76 24 25
+55 3 3 2 1 0 9 10 77 69
+56 3 3 2 1 0 69 77 78 70
+57 3 3 2 1 0 70 78 79 71
+58 3 3 2 1 0 71 79 80 72
+59 3 3 2 1 0 72 80 81 73
+60 3 3 2 1 0 73 81 82 74
+61 3 3 2 1 0 74 82 83 75
+62 3 3 2 1 0 75 83 84 76
+63 3 3 2 1 0 76 84 23 24
+64 3 3 2 1 0 10 11 85 77
+65 3 3 2 1 0 77 85 86 78
+66 3 3 2 1 0 78 86 87 79
+67 3 3 2 1 0 79 87 88 80
+68 3 3 2 1 0 80 88 89 81
+69 3 3 2 1 0 81 89 90 82
+70 3 3 2 1 0 82 90 91 83
+71 3 3 2 1 0 83 91 92 84
+72 3 3 2 1 0 84 92 22 23
+73 3 3 2 1 0 11 12 93 85
+74 3 3 2 1 0 85 93 94 86
+75 3 3 2 1 0 86 94 95 87
+76 3 3 2 1 0 87 95 96 88
+77 3 3 2 1 0 88 96 97 89
+78 3 3 2 1 0 89 97 98 90
+79 3 3 2 1 0 90 98 99 91
+80 3 3 2 1 0 91 99 100 92
+81 3 3 2 1 0 92 100 21 22
+82 3 3 2 1 0 12 2 13 93
+83 3 3 2 1 0 93 13 14 94
+84 3 3 2 1 0 94 14 15 95
+85 3 3 2 1 0 95 15 16 96
+86 3 3 2 1 0 96 16 17 97
+87 3 3 2 1 0 97 17 18 98
+88 3 3 2 1 0 98 18 19 99
+89 3 3 2 1 0 99 19 20 100
+90 3 3 2 1 0 100 20 3 21
+$EndElements