59 bool read_binary =
fileExists(basename +
".bpts");
60 std::string pts_file = read_binary ? basename +
".bpts" : basename +
".pts";
62 FILE* pts_fd = fopen(pts_file.c_str(),
"r");
64 fprintf(stderr,
"Error: could not open file: %s. Aborting!\n", pts_file.c_str());
70 if(read_binary) fread(buffer,
sizeof(
char),
HDR_SIZE, pts_fd);
71 else fgets( buffer, 2048, pts_fd);
72 sscanf(buffer,
"%lu", &numpts);
90 inline void read_headers(FILE* ele_fd, FILE* fib_fd,
bool read_binary,
size_t & numelem,
bool & twoFib)
97 fread(buffer,
sizeof(
char),
HDR_SIZE, ele_fd);
102 sscanf(buffer,
"%lu", &numelem);
107 fread(buffer,
sizeof(
char),
HDR_SIZE, fib_fd);
108 sscanf(buffer,
"%d", &numaxes);
114 short nread = sscanf(buffer,
"%f %f %f %f %f %f", fbuff, fbuff+1, fbuff+2, fbuff+3, fbuff+4, fbuff+5);
116 case 1: numaxes = fbuff[0];
break;
117 case 3: numaxes = 1;
break;
118 case 6: numaxes = 2;
break;
120 fprintf(stderr,
"%s error: read %d values in line 1 of fiber file." 121 "This is unexpected! Setting number of fibers to 1.\n",
129 fseek(fib_fd, 0, SEEK_SET);
132 twoFib = numaxes == 2;
147 inline void write_elem_headers(FILE* & ele_fd, FILE* & fib_fd,
bool binary,
size_t numelem,
bool twoFib)
153 int numfibres = twoFib ? 2 : 1;
158 snprintf(header,
sizeof header,
"%lu %d %d", numelem,
MT_ENDIANNESS, checksum);
159 fwrite(header,
HDR_SIZE,
sizeof(
char), ele_fd);
162 fprintf(ele_fd,
"%lu\n", numelem);
167 snprintf(header,
sizeof header,
"%d %lu %d %d", numfibres, (
unsigned long int)numelem,
MT_ENDIANNESS, checksum);
168 fwrite(header,
sizeof(
char),
HDR_SIZE, fib_fd);
171 fprintf(fib_fd,
"%d\n", numfibres);
192 snprintf(header,
sizeof header,
"%lu %d %d", numpts,
MT_ENDIANNESS, checksum);
193 fwrite(header,
HDR_SIZE,
sizeof(
char), pts_fd);
196 fprintf(pts_fd,
"%lu\n", numpts);
209 template<
class T,
class S>
212 const int bufsize = 2048;
213 const int max_elem_size = 9;
215 char buffer[bufsize];
217 int n[max_elem_size];
222 mesh.
dsp.resize(bsize+1);
223 mesh.
tag.resize(bsize);
225 mesh.
con.resize(bsize*max_elem_size);
227 size_t ele_read = 0, con_size = 0;
229 T* elem = mesh.
con.data();
231 for(
size_t i=0; i<bsize; i++)
235 size_t r = fread(&etbuff,
sizeof(
int), 1, fd);
240 char* ptr = fgets( buffer, bufsize, fd);
241 if(ptr == NULL)
break;
242 sscanf(buffer,
"%s %d %d %d %d %d %d %d %d %d", etype_str, n, n+1, n+2, n+3, n+4, n+5, n+6, n+7, n+8);
278 fprintf(stderr,
"Error: Unsupported element type!\n");
281 mesh.
dsp[i+1] = mesh.
dsp[i] + nodes;
286 fread(n,
sizeof(
int), nodes+1, fd);
288 for(
int j=0; j<nodes; j++) elem[j] = n[j];
290 mesh.
tag[i] = n[nodes];
291 ref_eidx[i] = bstart + i;
296 mesh.
dsp.resize(ele_read+1);
297 mesh.
tag.resize(ele_read);
299 ref_eidx.
resize(ele_read);
300 mesh.
con.resize(con_size);
304 if(ele_read < bsize) {
305 fprintf(stderr,
"\n\n%s error: Corrupted mesh!\nExpected to read %zu elements, but read only %zu! Aborting!\n\n", __func__, bsize, ele_read);
318 template<
class T,
class S>
321 const T* con = mesh.
con.data();
327 for(
size_t eidx=0; eidx < mesh.
l_numelem; eidx++)
329 fwrite(&mesh.
type[eidx], 1,
sizeof(
elem_t), fd);
330 T esize = mesh.
dsp[eidx+1] - mesh.
dsp[eidx];
332 wbuff[esize] = mesh.
tag[eidx];
334 fwrite(wbuff, esize+1,
sizeof(
int), fd);
340 for(
size_t eidx=0; eidx < mesh.
l_numelem; eidx++)
342 switch(mesh.
type[eidx])
345 fprintf(fd,
"Ln %d %d %d\n", con[0], con[1], mesh.
tag[eidx]);
350 fprintf(fd,
"Tr %d %d %d %d\n", con[0], con[1], con[2], mesh.
tag[eidx]);
355 fprintf(fd,
"Qd %d %d %d %d %d\n", con[0], con[1], con[2], con[3], mesh.
tag[eidx]);
360 fprintf(fd,
"Tt %d %d %d %d %d\n", con[0], con[1], con[2], con[3], mesh.
tag[eidx]);
365 fprintf(fd,
"Py %d %d %d %d %d %d\n", con[0], con[1], con[2], con[3],
366 con[4], mesh.
tag[eidx]);
371 fprintf(fd,
"Pr %d %d %d %d %d %d %d\n", con[0], con[1], con[2], con[3],
372 con[4], con[5], mesh.
tag[eidx]);
377 fprintf(fd,
"Hx %d %d %d %d %d %d %d %d %d\n", con[0], con[1], con[2], con[3],
378 con[4], con[5], con[6], con[7], mesh.
tag[eidx]);
383 fprintf(stderr,
"Error: Unsupported element type!\n");
400 template<
class T,
class S>
403 const int bufsize = 2048;
404 char buffer[bufsize];
414 for(
size_t i=0; i<bsize; i++) {
416 size_t r = fread(fib,
sizeof(
float), 6, fd);
418 fprintf(stderr,
"%s error: Orthotropy blon file does not contain two fiber data for all elements. Aborting!\n", __func__);
424 char* ptr = fgets(buffer, bufsize, fd);
425 if(ptr == NULL)
break;
427 size_t r = sscanf(buffer,
"%f %f %f %f %f %f", fib, fib+1, fib+2, fib+3, fib+4, fib+5);
429 fprintf(stderr,
"%s error: Orthotropy lon file does not contain data two fiber for all elements. Aborting!\n", __func__);
434 mesh.
fib[nr*3+0] = fib[0];
435 mesh.
fib[nr*3+1] = fib[1];
436 mesh.
fib[nr*3+2] = fib[2];
437 mesh.
she[nr*3+0] = fib[3];
438 mesh.
she[nr*3+1] = fib[4];
439 mesh.
she[nr*3+2] = fib[5];
450 for(
size_t i=0; i<bsize; i++) {
452 size_t r = fread(fib,
sizeof(
float), 3, fd);
454 fprintf(stderr,
"Error: Orthotropy blon file does not contain data for all elements. i=%i, bsize=%i Aborting!\n", i, bsize);
460 char* ptr = fgets( buffer, bufsize, fd);
461 if(ptr == NULL)
break;
463 size_t r = sscanf(buffer,
"%f %f %f", fib, fib+1, fib+2);
465 fprintf(stderr,
"Error: Orthotropy lon file does not contain data for all elements. Aborting!\n");
470 mesh.
fib[nr*3+0] = fib[0];
471 mesh.
fib[nr*3+1] = fib[1];
472 mesh.
fib[nr*3+2] = fib[2];
479 if(!err && nr != bsize) {
480 fprintf(stderr,
"\n\n%s error: Corrupted mesh!\nExpected to read %zu fibers, but read only %zu! Aborting!\n\n", __func__, bsize, nr);
497 template<
class T,
class S>
507 fib[0] = mesh.
fib[i*3+0];
508 fib[1] = mesh.
fib[i*3+1];
509 fib[2] = mesh.
fib[i*3+2];
510 fib[3] = mesh.
she[i*3+0];
511 fib[4] = mesh.
she[i*3+1];
512 fib[5] = mesh.
she[i*3+2];
515 fwrite(fib, 6,
sizeof(
float), fd);
517 fprintf(fd,
"%f %f %f %f %f %f\n", fib[0], fib[1], fib[2], fib[3], fib[4], fib[5]);
524 fib[0] = mesh.
fib[i*3+0];
525 fib[1] = mesh.
fib[i*3+1];
526 fib[2] = mesh.
fib[i*3+2];
529 fwrite(fib, 3,
sizeof(
float), fd);
531 fprintf(fd,
"%f %f %f\n", fib[0], fib[1], fib[2]);
545 template<
class T,
class S>
549 MPI_Comm comm = mesh.
comm;
552 MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
554 FILE* ele_fd = NULL, *fib_fd = NULL;
555 size_t gnumelems = 0;
557 bool read_binary =
false;
562 read_binary =
fileExists(basename +
".belem");
563 std::string ele_file = read_binary ? basename +
".belem" : basename +
".elem";
564 std::string fib_file = read_binary ? basename +
".blon" : basename +
".lon";
566 ele_fd = fopen(ele_file.c_str(),
"r");
568 fprintf(stderr,
"Error: could not open file: %s. Aborting!\n", ele_file.c_str());
572 fib_fd = fopen(fib_file.c_str(),
"r");
574 fprintf(stderr,
"Error: could not open file: %s. Aborting!\n", fib_file.c_str());
579 MPI_Allreduce(MPI_IN_PLACE, &err, 1, MPI_INT, MPI_SUM, comm);
585 read_headers(ele_fd, fib_fd, read_binary, gnumelems, twoFib);
587 MPI_Bcast(&gnumelems,
sizeof(
size_t), MPI_BYTE, 0, comm);
588 MPI_Bcast(&twoFib,
sizeof(
bool), MPI_BYTE, 0, comm);
593 divide(gnumelems, size, block_cnt);
604 for(
int pid=1; pid<size; pid++) {
606 read_elem_block(ele_fd, read_binary, block_dsp[pid], block_cnt[pid], meshbuff);
607 read_fib_block (fib_fd, read_binary, twoFib, block_cnt[pid], meshbuff);
613 MPI_Send(meshbuff.
dsp.data(), meshbuff.
dsp.size()*
sizeof(T), MPI_BYTE, pid,
SF_MPITAG, comm);
614 MPI_Send(meshbuff.
tag.data(), meshbuff.
tag.size()*
sizeof(T), MPI_BYTE, pid,
SF_MPITAG, comm);
615 MPI_Send(ref_eidx.
data(), ref_eidx.
size()*
sizeof(T), MPI_BYTE, pid,
SF_MPITAG, comm);
622 size_t con_size = meshbuff.
con.size();
623 MPI_Send(&con_size,
sizeof(
size_t), MPI_BYTE, pid,
SF_MPITAG, comm);
624 MPI_Send(meshbuff.
con.data(), con_size*
sizeof(T), MPI_BYTE, pid,
SF_MPITAG, comm);
645 MPI_Recv(mesh.
dsp.data(), mesh.
dsp.size()*
sizeof(T), MPI_BYTE, 0,
SF_MPITAG, comm, &stat);
646 MPI_Recv(mesh.
tag.data(), mesh.
tag.size()*
sizeof(T), MPI_BYTE, 0,
SF_MPITAG, comm, &stat);
647 MPI_Recv(ref_eidx.
data(), ref_eidx.
size()*
sizeof(T), MPI_BYTE, 0,
SF_MPITAG, comm, &stat);
654 MPI_Recv(&con_size,
sizeof(
size_t), MPI_BYTE, 0,
SF_MPITAG, comm, &stat);
655 mesh.
con.resize(con_size);
656 MPI_Recv(mesh.
con.data(), con_size*
sizeof(T), MPI_BYTE, 0,
SF_MPITAG, comm, &stat);
673 template<
class T,
class S>
676 const MPI_Comm comm = mesh.
comm;
679 MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
681 FILE* ele_fd = NULL, *fib_fd = NULL;
683 std::string ele_file = binary ? basename +
".belem" : basename +
".elem";
684 std::string fib_file = binary ? basename +
".blon" : basename +
".lon";
688 ele_fd = fopen(ele_file.c_str(),
"w");
690 fprintf(stderr,
"Error: could not open file: %s. Aborting!\n", ele_file.c_str());
693 fib_fd = fopen(fib_file.c_str(),
"w");
695 fprintf(stderr,
"Error: could not open file: %s. Aborting!\n", fib_file.c_str());
705 for(
int pid=0; pid < size; pid++)
709 ele_fd = fopen(ele_file.c_str(),
"a");
711 fprintf(stderr,
"Error: could not open file: %s. Aborting!\n", ele_file.c_str());
714 fib_fd = fopen(fib_file.c_str(),
"a");
716 fprintf(stderr,
"Error: could not open file: %s. Aborting!\n", fib_file.c_str());
731 template<
class T,
class S>
734 const MPI_Comm comm = surfmesh.
comm;
735 const bool binary =
false;
738 MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
744 ele_fd = fopen(surffile.c_str(),
"w");
746 fprintf(stderr,
"Error: could not open file: %s. Aborting!\n", surffile.c_str());
749 fprintf(ele_fd,
"%lu\n", surfmesh.
g_numelem);
754 for(
int pid=0; pid < size; pid++)
757 ele_fd = fopen(surffile.c_str(),
"a");
759 fprintf(stderr,
"Error: could not open file: %s. Aborting!\n", surffile.c_str());
782 const int bufsize = 2048;
783 char buffer[bufsize];
789 for(
size_t i=0; i<bsize; i++) {
791 size_t r = fread(pts,
sizeof(
float), 3, fd);
795 char* ptr = fgets( buffer, bufsize, fd);
796 if(ptr == NULL)
break;
797 sscanf(buffer,
"%f %f %f", pts, pts+1, pts+2);
799 xyz[nr*3+0] = pts[0];
800 xyz[nr*3+1] = pts[1];
801 xyz[nr*3+2] = pts[2];
819 size_t nnodes = xyz.
size() / 3;
822 for(
size_t i=0; i<nnodes; i++)
829 fwrite(pt, 3,
sizeof(
float), fd);
831 fprintf(fd,
"%f %f %f\n", pt[0], pt[1], pt[2]);
843 template<
class T,
class S>
847 MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
851 bool read_binary =
false;
856 std::string pts_file = read_binary ? basename +
".bpts" : basename +
".pts";
858 pts_fd = fopen(pts_file.c_str(),
"r");
860 fprintf(stderr,
"Error: could not open file: %s. Aborting!\n", pts_file.c_str());
866 if(read_binary) fread(buffer,
sizeof(
char),
HDR_SIZE, pts_fd);
867 else fgets( buffer, 2048, pts_fd);
868 sscanf(buffer,
"%lu", &gnumpts);
870 MPI_Bcast(&gnumpts,
sizeof(
size_t), MPI_BYTE, 0, comm);
873 size_t blocksize = (gnumpts + size - 1) / size;
881 for(
int pid = 1; pid < size; pid++) {
883 long int numsend = buff.
size();
885 MPI_Send(&numsend, 1, MPI_LONG, pid,
SF_MPITAG, comm);
886 MPI_Send(buff.
data(), numsend*
sizeof(S), MPI_BYTE, pid,
SF_MPITAG, comm);
891 long int numrecv = 0;
892 MPI_Recv(&numrecv, 1, MPI_LONG, 0,
SF_MPITAG, comm, &stat);
895 MPI_Recv(pts.
data(), numrecv*
sizeof(S), MPI_BYTE, 0,
SF_MPITAG, comm, &stat);
898 long int mysize = pts.
size() / 3;
902 interval(ptsidx, layout[rank], layout[rank+1]);
904 if(rank == 0) fclose(pts_fd);
914 template<
class T,
class S>
919 assert(meshlist.size() > 0);
920 assert(pts.
size() == (ptsidx.
size() * 3));
924 MPI_Comm comm = (*meshlist.begin())->comm;
927 MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
932 for(
int pid = 0; pid < size; pid++) {
933 size_t numsend = pts.
size();
934 MPI_Bcast(&numsend,
sizeof(
size_t), MPI_BYTE, pid, comm);
942 idxbuff.
resize(numsend / 3);
945 MPI_Bcast(ptsbuff.
data(), ptsbuff.
size()*
sizeof(S), MPI_BYTE, pid, comm);
946 MPI_Bcast(idxbuff.
data(), idxbuff.
size()*
sizeof(T), MPI_BYTE, pid, comm);
950 for(
auto it = meshlist.begin(); it != meshlist.end(); ++it)
958 for(
size_t i=0; i<rnod.
size(); i++)
966 for(
size_t j=0; j<idxbuff.
size(); j++) {
967 if(g2l.
count(idxbuff[j])) {
974 for(
size_t j=0; j<ridx.
size(); j++) {
975 T w = widx[j], r = ridx[j];
976 mesh.
xyz[w*3+0] = ptsbuff[r*3+0];
977 mesh.
xyz[w*3+1] = ptsbuff[r*3+1];
978 mesh.
xyz[w*3+2] = ptsbuff[r*3+2];
990 template<
class T,
class S>
993 FILE* vtk_file = fopen(file.c_str(),
"w");
994 if(vtk_file == NULL)
return;
997 fprintf (vtk_file,
"# vtk DataFile Version 3.0\n");
998 fprintf (vtk_file,
"vtk output\n");
999 fprintf (vtk_file,
"binary\n");
1000 fprintf (vtk_file,
"DATASET UNSTRUCTURED_GRID\n\n");
1001 fprintf (vtk_file,
"POINTS %lu float\n", mesh.
l_numpts);
1006 for (
unsigned long int i=0; i<mesh.
l_numpts; i++ ) {
1007 pts[0] =
htobe(p[0]);
1008 pts[1] =
htobe(p[1]);
1009 pts[2] =
htobe(p[2]);
1010 fwrite(pts,
sizeof(
float), 3, vtk_file);
1013 fprintf(vtk_file,
"\n");
1015 fprintf (vtk_file,
"CELL_TYPES %lu\n", mesh.
l_numelem);
1016 unsigned long int valcount = 0;
1018 for(
unsigned long int i=0; i< mesh.
l_numelem; i++) {
1045 vtk_type =
htobe(vtk_type);
1046 fwrite(&vtk_type,
sizeof(
int), 1, vtk_file);
1047 valcount += mesh.
dsp[i+1] - mesh.
dsp[i] + 1;
1049 fprintf(vtk_file,
"\n");
1051 fprintf(vtk_file,
"CELLS %lu %lu\n", mesh.
l_numelem, valcount);
1052 const T* elem = mesh.
con.data();
1053 for(
unsigned long int i=0; i<mesh.
l_numelem; i++)
1055 int nodes = mesh.
dsp[i+1] - mesh.
dsp[i], n;
1056 int be_nodes =
htobe(nodes);
1057 fwrite(&be_nodes,
sizeof(
int), 1, vtk_file);
1059 for(
int j=0; j<nodes; j++) {
1060 n = elem[j]; n =
htobe(n);
1061 fwrite(&n,
sizeof(
int), 1, vtk_file);
1065 fprintf(vtk_file,
"\n");
1067 fprintf (vtk_file,
"CELL_DATA %lu \n", mesh.
l_numelem);
1068 fprintf (vtk_file,
"SCALARS elemTag int 1\n");
1069 fprintf (vtk_file,
"LOOKUP_TABLE default\n");
1070 for (
unsigned long int i=0; i<mesh.
l_numelem; i++ ) {
1072 fwrite(&t,
sizeof(
int), 1, vtk_file);
1078 fprintf (vtk_file,
"VECTORS fiber float\n");
1079 for (
unsigned long int i=0; i<mesh.
l_numelem; i++ ) {
1080 pts[0] =
htobe(
float(mesh.
fib[i*3+0])),
1081 pts[1] =
htobe(
float(mesh.
fib[i*3+1])),
1082 pts[2] =
htobe(
float(mesh.
fib[i*3+2]));
1083 fwrite(pts,
sizeof(
float), 3, vtk_file);
1085 fprintf(vtk_file,
"\n");
1088 fprintf (vtk_file,
"VECTORS sheet float\n");
1089 for (
unsigned long int i=0; i<mesh.
l_numelem; i++ ) {
1090 pts[0] =
htobe(
float(mesh.
she[i*3+0])),
1091 pts[1] =
htobe(
float(mesh.
she[i*3+1])),
1092 pts[2] =
htobe(
float(mesh.
she[i*3+2]));
1093 fwrite(pts,
sizeof(
float), 3, vtk_file);
1095 fprintf(vtk_file,
"\n");
The mesh storage class. It contains both element and vertex data.
The vector class and related algorithms.
void write_elements(const meshdata< T, S > &mesh, bool binary, std::string basename)
Write the element data (elements and fibers) of a CARP mesh.
void writeVTKmesh_binary(const meshdata< T, S > &mesh, std::string file)
Write a mesh in binary vtk format.
vector< T > & get_numbering(SF_nbr nbr_type)
Get the vector defining a certain numbering.
void write_elem_headers(FILE *&ele_fd, FILE *&fib_fd, bool binary, size_t numelem, bool twoFib)
Write the header of the element and fiber files.
void write_elem_block(FILE *fd, bool write_binary, const meshdata< T, S > &mesh)
Write the local element block to a file.
MPI_Comm comm
the parallel mesh is defined on a MPI world
void layout_from_count(const T count, vector< T > &layout, MPI_Comm comm)
void write_pts_block(FILE *&fd, bool write_binary, const vector< S > &xyz)
Write a chunk of points to a file.
vector< T > tag
element tag
vector< T > dsp
connectivity starting index of each element
void interval(vector< T > &vec, size_t start, size_t end)
Create an integer interval between start and end.
vector< S > fib
fiber direction
void divide(const size_t gsize, const size_t num_parts, vector< T > &loc_sizes)
divide gsize into num_parts local parts with even distribution of the remainder
The nodal numbering of the reference mesh (the one stored on HD).
#define SF_MPITAG
the MPI tag when communicating
vector< S > xyz
node cooridnates
T * data()
Pointer to the vector's start.
void read_elements(meshdata< T, S > &mesh, std::string basename)
Read the element data (elements and fibers) of a CARP mesh.
void write_pts_header(FILE *&pts_fd, bool binary, size_t numpts)
Write the header of the points file.
Functions related to network communication.
void write_surface(const meshdata< T, S > &surfmesh, std::string surffile)
void write_fib_block(FILE *&fd, bool write_binary, const meshdata< T, S > &mesh)
Write the local chunk of fibers to a file.
void read_headers(FILE *ele_fd, FILE *fib_fd, bool read_binary, size_t &numelem, bool &twoFib)
Read the header from the element and fiber files.
void insert_points(const vector< S > &pts, const vector< T > &ptsidx, std::list< meshdata< T, S > *> &meshlist)
Insert the points from the read-in buffers into a list of distributed meshes.
void read_pts_block(FILE *&fd, bool read_binary, size_t bsize, vector< S > &xyz)
Read a chunk of points from a file descriptor.
The element numbering of the reference mesh (the one stored on HD).
bool fileExists(std::string filename)
Function which checks if a given file exists.
void read_elem_block(FILE *&fd, bool read_binary, size_t bstart, size_t bsize, meshdata< T, S > &mesh)
Read a block of size bsize from an CARP element file.
size_t read_num_pts(std::string basename)
Function returns the number of points in a CARP points file.
void read_fib_block(FILE *&fd, bool read_binary, bool twoFib, size_t bsize, meshdata< T, S > &mesh)
Read a chunk of fibers from a file descriptor.
elem_t getElemTypeID(char *eletype)
Generate element type enum from string.
void localize(SF_nbr nbr_type)
Localize the connectivity data w.r.t. a given numbering.
hm_int count(const K &key) const
Check if key exists.
size_t size() const
The current size of the vector.
size_t l_numpts
local number of points
A vector storing arbitrary data.
size_t l_numelem
local number of elements
void vec_assign(S *lhs, const V *rhs, size_t size)
Assign the values in rhs to lhs. The data-type of rhs is cast to the type of lhs. ...
vector< S > she
sheet direction
size_t g_numelem
global number of elements
vector< T > & register_numbering(SF_nbr nbr_type)
Register a new numbering to the mesh and return the associated index vector.
void resize(size_t n)
Resize a vector.
vector< elem_t > type
element type
void dsp_from_cnt(const vector< T > &cnt, vector< T > &dsp)
Compute displacements from counts.
void read_points(const std::string basename, const MPI_Comm comm, vector< S > &pts, vector< T > &ptsidx)
Read the points and insert them into a list of meshes.