openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
SF_mesh_io.h
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 // openCARP is an open cardiac electrophysiology simulator.
3 //
4 // Copyright (C) 2020 openCARP project
5 //
6 // This program is licensed under the openCARP Academic Public License (APL)
7 // v1.0: You can use and redistribute it and/or modify it in non-commercial
8 // academic environments under the terms of APL as published by the openCARP
9 // project v1.0, or (at your option) any later version. Commercial use requires
10 // a commercial license (info@opencarp.org).
11 //
12 // This program is distributed without any warranty; see the openCARP APL for
13 // more details.
14 //
15 // You should have received a copy of the openCARP APL along with this program
16 // and can find it online: http://www.opencarp.org/license
17 // ----------------------------------------------------------------------------
18 
27 #ifndef _SF_MESH_IO
28 #define _SF_MESH_IO
29 
30 #include <cstddef>
31 #include <list>
32 #include <string>
33 #include <cstring>
34 #include <unistd.h>
35 
36 #include <mpi.h>
37 
38 #include "SF_container.h"
39 #include "SF_io_base.h"
40 #include "SF_network.h"
41 #include "SF_globals.h"
42 #include "SF_vector.h"
43 
44 #define HDR_SIZE 1024
45 
46 namespace SF {
47 
55 inline size_t read_num_pts(std::string basename)
56 {
57  size_t numpts;
58 
59  bool read_binary = fileExists(basename + ".bpts");
60  std::string pts_file = read_binary ? basename + ".bpts" : basename + ".pts";
61 
62  FILE* pts_fd = fopen(pts_file.c_str(), "r");
63  if(!pts_fd) {
64  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", pts_file.c_str());
65  exit(1);
66  }
67 
68  // read number of points
69  char buffer[2048];
70  if(read_binary) fread(buffer, sizeof(char), HDR_SIZE, pts_fd);
71  else fgets( buffer, 2048, pts_fd);
72  sscanf(buffer, "%lu", &numpts);
73  fclose(pts_fd);
74 
75  return numpts;
76 }
77 
90 inline void read_headers(FILE* ele_fd, FILE* fib_fd, bool read_binary, size_t & numelem, int & nFib)
91 {
92  char buffer[HDR_SIZE];
93  // first read number of elems
94  if(read_binary) {
95  fread(buffer, sizeof(char), HDR_SIZE, ele_fd);
96  }
97  else {
98  fgets(buffer, HDR_SIZE, ele_fd);
99  }
100  sscanf(buffer, "%lu", &numelem);
101 
102  // now read number of fiber axes
103  if(fib_fd) {
104  if(read_binary) {
105  fread(buffer, sizeof(char), HDR_SIZE, fib_fd);
106  sscanf(buffer, "%d", &nFib);
107  }
108  else {
109  // fiber files can have no header, thus we can potentially read 1, 3 or 6 values
110  fgets(buffer, HDR_SIZE, fib_fd);
111  float fbuff[6];
112  short nread = sscanf(buffer, "%f %f %f %f %f %f", fbuff, fbuff+1, fbuff+2, fbuff+3, fbuff+4, fbuff+5);
113  switch(nread) {
114  case 1: nFib = fbuff[0]; break;
115  case 3: nFib = 1; break;
116  case 6: nFib = 2; break;
117  default:
118  fprintf(stderr, "%s error: read %d values in line 1 of fiber file."
119  "This is unexpected! Setting number of fibers to 1.\n",
120  __func__, nread);
121  nFib = 1;
122  break;
123  }
124 
125  // if the fiber file had no header, we have to reset the FD
126  if(nread != 1)
127  fseek(fib_fd, 0, SEEK_SET);
128  }
129  }
130  else
131  nFib = 0;
132 }
145 inline void write_elem_headers(FILE* & ele_fd, FILE* & fib_fd, bool binary, size_t numelem, int nFib)
146 {
147 
148  // we need to write the CARP header
149  char header[HDR_SIZE];
150  int checksum = 666;
151 
152  // element header
153  if(binary) {
154  // first elements, second endianness, third checksum
155  snprintf(header, sizeof header, "%lu %d %d", numelem, MT_ENDIANNESS, checksum);
156  fwrite(header, HDR_SIZE, sizeof(char), ele_fd);
157  }
158  else
159  fprintf(ele_fd, "%lu\n", numelem);
160 
161  // fiber header
162  if(nFib >= 1) {
163  if(binary) {
164  // first numaxes, second numelems, third endianness, last checksum
165  snprintf(header, sizeof header, "%d %lu %d %d", nFib, (unsigned long int)numelem, MT_ENDIANNESS, checksum);
166  fwrite(header, sizeof(char), HDR_SIZE, fib_fd);
167  }
168  else
169  fprintf(fib_fd, "%d\n", nFib);
170  }
171 }
182 inline void write_pts_header(FILE* & pts_fd, bool binary, size_t numpts)
183 {
184 
185  // we need to write the CARP header
186  char header[HDR_SIZE];
187  int checksum = 666;
188 
189  // header
190  if(binary) {
191  snprintf(header, sizeof header, "%lu %d %d", numpts, MT_ENDIANNESS, checksum);
192  fwrite(header, HDR_SIZE, sizeof(char), pts_fd);
193  }
194  else
195  fprintf(pts_fd, "%lu\n", numpts);
196 }
197 
208 template<class T, class S>
209 inline void read_elem_block(FILE* & fd, bool read_binary, size_t bstart, size_t bsize, meshdata<T, S> & mesh)
210 {
211  const int bufsize = 2048;
212  const int max_elem_size = 9;
213 
214  char buffer[bufsize];
215  char etype_str[8];
216  int n[max_elem_size];
217 
218  vector<T> & ref_eidx = mesh.register_numbering(NBR_ELEM_REF);
219  ref_eidx.resize(bsize);
220 
221  mesh.dsp.resize(bsize+1);
222  mesh.tag.resize(bsize);
223  mesh.type.resize(bsize);
224  mesh.con.resize(bsize*max_elem_size);
225 
226  size_t ele_read = 0, con_size = 0;
227  mesh.dsp[0] = 0;
228  T* elem = mesh.con.data();
229 
230  for(size_t i=0; i<bsize; i++)
231  {
232  if(read_binary) {
233  int etbuff;
234  size_t r = fread(&etbuff, sizeof(int), 1, fd);
235  if(r != 1) break;
236  mesh.type[i] = (elem_t)etbuff;
237  }
238  else {
239  char* ptr = fgets( buffer, bufsize, fd);
240  if(ptr == NULL) break;
241  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);
242  mesh.type[i] = getElemTypeID(etype_str);
243  }
244 
245  T nodes;
246  switch(mesh.type[i])
247  {
248  case Line:
249  nodes = 2;
250  break;
251 
252  case Tri:
253  nodes = 3;
254  break;
255 
256  case Quad:
257  nodes = 4;
258  break;
259 
260  case Tetra:
261  nodes = 4;
262  break;
263 
264  case Pyramid:
265  nodes = 5;
266  break;
267 
268  case Prism:
269  nodes = 6;
270  break;
271 
272  case Hexa:
273  nodes = 8;
274  break;
275 
276  default:
277  fprintf(stderr, "Error: Unsupported element type!\n");
278  exit(1);
279  }
280  mesh.dsp[i+1] = mesh.dsp[i] + nodes;
281  con_size += nodes;
282 
283  // copy the element connectivity
284  if(read_binary) {
285  fread(n, sizeof(int), nodes+1, fd);
286  }
287  for(int j=0; j<nodes; j++) elem[j] = n[j];
288 
289  mesh.tag[i] = n[nodes];
290  ref_eidx[i] = bstart + i;
291  ele_read++;
292 
293  elem += nodes;
294  }
295  mesh.dsp.resize(ele_read+1);
296  mesh.tag.resize(ele_read);
297  mesh.type.resize(ele_read);
298  ref_eidx.resize(ele_read);
299  mesh.con.resize(con_size);
300 
301  mesh.l_numelem = ele_read;
302 }
303 
304 
312 template<class T, class S>
313 inline void write_elem_block(FILE* fd, bool write_binary, const meshdata<T, S> & mesh)
314 {
315  const T* con = mesh.con.data();
316 
317  if(write_binary)
318  {
319  int wbuff[9];
320 
321  for(size_t eidx=0; eidx < mesh.l_numelem; eidx++)
322  {
323  fwrite(&mesh.type[eidx], 1, sizeof(elem_t), fd);
324  T esize = mesh.dsp[eidx+1] - mesh.dsp[eidx];
325  vec_assign(wbuff, con, esize); // copy-convert to int
326  wbuff[esize] = mesh.tag[eidx]; // convert to int
327 
328  fwrite(wbuff, esize+1, sizeof(int), fd);
329  con += esize;
330  }
331  }
332  else
333  {
334  for(size_t eidx=0; eidx < mesh.l_numelem; eidx++)
335  {
336  switch(mesh.type[eidx])
337  {
338  case Line:
339  fprintf(fd, "Ln %d %d %d\n", con[0], con[1], mesh.tag[eidx]);
340  con += 2;
341  break;
342 
343  case Tri:
344  fprintf(fd, "Tr %d %d %d %d\n", con[0], con[1], con[2], mesh.tag[eidx]);
345  con += 3;
346  break;
347 
348  case Quad:
349  fprintf(fd, "Qd %d %d %d %d %d\n", con[0], con[1], con[2], con[3], mesh.tag[eidx]);
350  con += 4;
351  break;
352 
353  case Tetra:
354  fprintf(fd, "Tt %d %d %d %d %d\n", con[0], con[1], con[2], con[3], mesh.tag[eidx]);
355  con += 4;
356  break;
357 
358  case Pyramid:
359  fprintf(fd, "Py %d %d %d %d %d %d\n", con[0], con[1], con[2], con[3],
360  con[4], mesh.tag[eidx]);
361  con += 5;
362  break;
363 
364  case Prism:
365  fprintf(fd, "Pr %d %d %d %d %d %d %d\n", con[0], con[1], con[2], con[3],
366  con[4], con[5], mesh.tag[eidx]);
367  con += 6;
368  break;
369 
370  case Hexa:
371  fprintf(fd, "Hx %d %d %d %d %d %d %d %d %d\n", con[0], con[1], con[2], con[3],
372  con[4], con[5], con[6], con[7], mesh.tag[eidx]);
373  con += 8;
374  break;
375 
376  default:
377  fprintf(stderr, "Error: Unsupported element type!\n");
378  exit(1);
379  }
380  }
381  }
382 }
383 
384 
394 template<class T, class S>
395 inline void read_fib_block(FILE* & fd, bool read_binary, int nFib, size_t bsize, meshdata<T, S> & mesh)
396 {
397  const int bufsize = 2048;
398  char buffer[bufsize];
399 
400  if(nFib == 2) {
401  float fib[6];
402  size_t nr = 0;
403 
404  mesh.fib.resize(bsize*3);
405  mesh.she.resize(bsize*3);
406 
407  for(size_t i=0; i<bsize; i++) {
408  if(read_binary) {
409  size_t r = fread(fib, sizeof(float), 6, fd);
410  if(r != 6) break;
411  }
412  else {
413  char* ptr = fgets( buffer, bufsize, fd);
414  if(ptr == NULL) break;
415  sscanf(buffer, "%f %f %f %f %f %f", fib, fib+1, fib+2, fib+3, fib+4, fib+5);
416  }
417  mesh.fib[nr*3+0] = fib[0];
418  mesh.fib[nr*3+1] = fib[1];
419  mesh.fib[nr*3+2] = fib[2];
420  mesh.she[nr*3+0] = fib[3];
421  mesh.she[nr*3+1] = fib[4];
422  mesh.she[nr*3+2] = fib[5];
423 
424  nr++;
425  }
426  mesh.fib.resize(nr*3);
427  mesh.she.resize(nr*3);
428  }
429  else if (nFib == 1) {
430  float fib[3];
431  size_t nr = 0;
432 
433  mesh.fib.resize(bsize*3);
434 
435  for(size_t i=0; i<bsize; i++) {
436  if(read_binary) {
437  size_t r = fread(fib, sizeof(float), 3, fd);
438  if(r != 3) break;
439  }
440  else {
441  char* ptr = fgets( buffer, bufsize, fd);
442  if(ptr == NULL) break;
443  sscanf(buffer, "%f %f %f", fib, fib+1, fib+2);
444  }
445 
446  mesh.fib[nr*3+0] = fib[0];
447  mesh.fib[nr*3+1] = fib[1];
448  mesh.fib[nr*3+2] = fib[2];
449 
450  nr++;
451  }
452  mesh.fib.resize(nr*3);
453  }
454 }
455 
456 
464 template<class T, class S>
465 inline void write_fib_block(FILE* & fd, bool write_binary, const meshdata<T, S> & mesh)
466 {
467  float fib[6];
468  int nFib = 0;
469  if(mesh.fib.size() > 0) {
470  nFib = 1;
471  if(mesh.she.size() == mesh.fib.size())
472  nFib = 2;
473  }
474 
475  if(nFib == 2)
476  {
477  for(size_t i=0; i<mesh.l_numelem; i++)
478  {
479  fib[0] = mesh.fib[i*3+0];
480  fib[1] = mesh.fib[i*3+1];
481  fib[2] = mesh.fib[i*3+2];
482  fib[3] = mesh.she[i*3+0];
483  fib[4] = mesh.she[i*3+1];
484  fib[5] = mesh.she[i*3+2];
485 
486  if(write_binary)
487  fwrite(fib, 6, sizeof(float), fd);
488  else
489  fprintf(fd, "%f %f %f %f %f %f\n", fib[0], fib[1], fib[2], fib[3], fib[4], fib[5]);
490  }
491  }
492  else if(nFib == 1)
493  {
494  for(size_t i=0; i<mesh.l_numelem; i++)
495  {
496  fib[0] = mesh.fib[i*3+0];
497  fib[1] = mesh.fib[i*3+1];
498  fib[2] = mesh.fib[i*3+2];
499 
500  if(write_binary)
501  fwrite(fib, 3, sizeof(float), fd);
502  else
503  fprintf(fd, "%f %f %f\n", fib[0], fib[1], fib[2]);
504  }
505  }
506 }
507 
517 template<class T, class S>
518 inline void read_elements(meshdata<T, S> & mesh, std::string basename, bool require_fibers = true)
519 {
520  // we use the communicator provided by the mesh
521  MPI_Comm comm = mesh.comm;
522 
523  int size, rank;
524  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
525 
526  FILE* ele_fd = NULL, *fib_fd = NULL;
527  size_t gnumelems = 0;
528  int nFib;
529  bool read_binary = false;
530  int err = 0;
531 
532  // open the files on root and read number of elements
533  if(rank == 0) {
534  read_binary = fileExists(basename + ".belem");
535  std::string ele_file = read_binary ? basename + ".belem" : basename + ".elem";
536  std::string fib_file = read_binary ? basename + ".blon" : basename + ".lon";
537 
538  ele_fd = fopen(ele_file.c_str(), "r");
539  if(!ele_fd) {
540  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", ele_file.c_str());
541  err++;
542  }
543 
544  fib_fd = fopen(fib_file.c_str(), "r");
545  if((!fib_fd) && require_fibers) {
546  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", fib_file.c_str());
547  err++;
548  }
549  }
550 
551  MPI_Allreduce(MPI_IN_PLACE, &err, 1, MPI_INT, MPI_SUM, comm);
552  if(err) {
553  exit(EXIT_FAILURE);
554  }
555 
556  if(rank == 0)
557  read_headers(ele_fd, fib_fd, read_binary, gnumelems, nFib);
558 
559  MPI_Bcast(&gnumelems, sizeof(size_t), MPI_BYTE, 0, comm);
560  MPI_Bcast(&nFib, sizeof(int), MPI_BYTE, 0, comm);
561  mesh.g_numelem = gnumelems;
562 
563  // compute the size of the block stored on each process
564  size_t blocksize = (gnumelems + size - 1) / size;
565 
566  // rank 0 reads mesh and distributes
567  if(rank == 0) {
568  // read own block
569  read_elem_block(ele_fd, read_binary, 0, blocksize, mesh);
570  if(nFib>=1)
571  read_fib_block(fib_fd, read_binary, nFib, blocksize, mesh);
572 
573  // read blocks of other ranks and communicate
574  meshdata<T, S> meshbuff;
575  for(int pid=1; pid<size; pid++) {
576  // read block in buffer
577  read_elem_block(ele_fd, read_binary, pid*blocksize, blocksize, meshbuff);
578  if(nFib>=1)
579  read_fib_block(fib_fd, read_binary, nFib, blocksize, meshbuff);
580  vector<T> & ref_eidx = meshbuff.get_numbering(NBR_ELEM_REF);
581 
582  // communicate
583  MPI_Send(&meshbuff.l_numelem, sizeof(size_t), MPI_BYTE, pid, SF_MPITAG, comm);
584  MPI_Send(meshbuff.dsp.data(), meshbuff.dsp.size()*sizeof(T), MPI_BYTE, pid, SF_MPITAG, comm);
585  MPI_Send(meshbuff.tag.data(), meshbuff.tag.size()*sizeof(T), MPI_BYTE, pid, SF_MPITAG, comm);
586  MPI_Send(ref_eidx.data(), ref_eidx.size()*sizeof(T), MPI_BYTE, pid, SF_MPITAG, comm);
587  MPI_Send(meshbuff.type.data(), meshbuff.type.size()*sizeof(elem_t),MPI_BYTE, pid, SF_MPITAG, comm);
588  if(nFib>=1)
589  MPI_Send(meshbuff.fib.data(), meshbuff.fib.size()*sizeof(S), MPI_BYTE, pid, SF_MPITAG, comm);
590  if(nFib==2)
591  MPI_Send(meshbuff.she.data(), meshbuff.she.size()*sizeof(S), MPI_BYTE, pid, SF_MPITAG, comm);
592 
593  size_t con_size = meshbuff.con.size();
594  MPI_Send(&con_size, sizeof(size_t), MPI_BYTE, pid, SF_MPITAG, comm);
595  MPI_Send(meshbuff.con.data(), con_size*sizeof(T), MPI_BYTE, pid, SF_MPITAG, comm);
596  }
597 
598  fclose(ele_fd);
599  if(nFib>=1) fclose(fib_fd);
600  }
601  else {
602  // the rest receives
603  MPI_Status stat;
604  MPI_Recv(&mesh.l_numelem, sizeof(size_t), MPI_BYTE, 0, SF_MPITAG, comm, &stat);
605 
606  vector<T> & ref_eidx = mesh.register_numbering(NBR_ELEM_REF);
607  ref_eidx.resize(mesh.l_numelem);
608 
609  mesh.dsp.resize(mesh.l_numelem+1);
610  mesh.tag.resize(mesh.l_numelem);
611  mesh.type.resize(mesh.l_numelem);
612 
613  if(nFib>=1) mesh.fib.resize(mesh.l_numelem*3);
614  if(nFib==2) mesh.she.resize(mesh.l_numelem*3);
615 
616  MPI_Recv(mesh.dsp.data(), mesh.dsp.size()*sizeof(T), MPI_BYTE, 0, SF_MPITAG, comm, &stat);
617  MPI_Recv(mesh.tag.data(), mesh.tag.size()*sizeof(T), MPI_BYTE, 0, SF_MPITAG, comm, &stat);
618  MPI_Recv(ref_eidx.data(), ref_eidx.size()*sizeof(T), MPI_BYTE, 0, SF_MPITAG, comm, &stat);
619  MPI_Recv(mesh.type.data(), mesh.type.size()*sizeof(elem_t),MPI_BYTE, 0, SF_MPITAG, comm, &stat);
620  if(nFib>=1)
621  MPI_Recv(mesh.fib.data(), mesh.fib.size()*sizeof(S), MPI_BYTE, 0, SF_MPITAG, comm, &stat);
622  if(nFib==2)
623  MPI_Recv(mesh.she.data(), mesh.she.size()*sizeof(S), MPI_BYTE, 0, SF_MPITAG, comm, &stat);
624 
625  size_t con_size;
626  MPI_Recv(&con_size, sizeof(size_t), MPI_BYTE, 0, SF_MPITAG, comm, &stat);
627  mesh.con.resize(con_size);
628  MPI_Recv(mesh.con.data(), con_size*sizeof(T), MPI_BYTE, 0, SF_MPITAG, comm, &stat);
629  }
630 
631  // localize connectivity (w.r.t. reference numbering)
632  mesh.localize(NBR_REF);
633 }
634 
635 
645 template<class T, class S>
646 inline void write_elements(const meshdata<T, S> & mesh, bool binary, std::string basename)
647 {
648  const MPI_Comm comm = mesh.comm;
649 
650  int size, rank;
651  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
652 
653  FILE* ele_fd = NULL, *fib_fd = NULL;
654  int nFib = 0;
655  if(mesh.fib.size() > 0) {
656  nFib = 1;
657  if(mesh.she.size() == mesh.fib.size())
658  nFib = 2;
659  }
660  std::string ele_file = binary ? basename + ".belem" : basename + ".elem";
661  std::string fib_file = binary ? basename + ".blon" : basename + ".lon";
662 
663  // open the files on root and write headers
664  if(rank == 0) {
665  ele_fd = fopen(ele_file.c_str(), "w");
666  if(!ele_fd) {
667  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", ele_file.c_str());
668  exit(1);
669  }
670  if(nFib>=1) {
671  fib_fd = fopen(fib_file.c_str(), "w");
672  if(!fib_fd) {
673  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", fib_file.c_str());
674  fclose(ele_fd);
675  exit(1);
676  }
677  }
678  write_elem_headers(ele_fd, fib_fd, binary, mesh.g_numelem, nFib);
679  fclose(ele_fd);
680  if (nFib>=1) fclose(fib_fd);
681  }
682 
683  // write mesh sequentially rank after rank
684  for(int pid=0; pid < size; pid++)
685  {
686  if(pid == rank)
687  {
688  ele_fd = fopen(ele_file.c_str(), "a");
689  if(!ele_fd) {
690  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", ele_file.c_str());
691  exit(1);
692  }
693  if(nFib >= 1) {
694  fib_fd = fopen(fib_file.c_str(), "a");
695  if(!fib_fd) {
696  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", fib_file.c_str());
697  fclose(ele_fd);
698  exit(1);
699  }
700  }
701  write_elem_block(ele_fd, binary, mesh);
702  if(nFib>=1) write_fib_block(fib_fd, binary, mesh);
703 
704  fclose(ele_fd);
705  if(nFib>=1) fclose(fib_fd);
706  }
707  MPI_Barrier(comm);
708  }
709 }
710 
711 template<class T, class S>
712 inline void write_surface(const meshdata<T, S> & surfmesh, std::string surffile)
713 {
714  const MPI_Comm comm = surfmesh.comm;
715  const bool binary = false;
716 
717  int size, rank;
718  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
719 
720  FILE* ele_fd = NULL;
721 
722  // open the files on root and write headers
723  if(rank == 0) {
724  ele_fd = fopen(surffile.c_str(), "w");
725  if(!ele_fd) {
726  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", surffile.c_str());
727  exit(1);
728  }
729  fprintf(ele_fd, "%lu\n", surfmesh.g_numelem);
730  fclose(ele_fd);
731  }
732 
733  // write surface sequentially rank after rank
734  for(int pid=0; pid < size; pid++)
735  {
736  if(pid == rank) {
737  ele_fd = fopen(surffile.c_str(), "a");
738  if(!ele_fd) {
739  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", surffile.c_str());
740  exit(1);
741  }
742 
743  write_elem_block(ele_fd, binary, surfmesh);
744  fclose(ele_fd);
745  }
746  MPI_Barrier(comm);
747  }
748 }
749 
750 
759 template<class S>
760 inline void read_pts_block(FILE* & fd, bool read_binary, size_t bsize, vector<S> & xyz)
761 {
762  const int bufsize = 2048;
763  char buffer[bufsize];
764  float pts[3];
765  xyz.resize(bsize*3);
766 
767  size_t nr = 0;
768 
769  for(size_t i=0; i<bsize; i++) {
770  if(read_binary) {
771  size_t r = fread(pts, sizeof(float), 3, fd);
772  if(r != 3) break;
773  }
774  else {
775  char* ptr = fgets( buffer, bufsize, fd);
776  if(ptr == NULL) break;
777  sscanf(buffer, "%f %f %f", pts, pts+1, pts+2);
778  }
779  xyz[nr*3+0] = pts[0];
780  xyz[nr*3+1] = pts[1];
781  xyz[nr*3+2] = pts[2];
782 
783  nr++;
784  }
785  xyz.resize(nr*3);
786 }
787 
788 
796 template<class S>
797 inline void write_pts_block(FILE* & fd, bool write_binary, const vector<S> & xyz)
798 {
799  size_t nnodes = xyz.size() / 3;
800  float pt[3];
801 
802  for(size_t i=0; i<nnodes; i++)
803  {
804  pt[0] = xyz[i*3+0];
805  pt[1] = xyz[i*3+1];
806  pt[2] = xyz[i*3+2];
807 
808  if(write_binary)
809  fwrite(pt, 3, sizeof(float), fd);
810  else
811  fprintf(fd, "%f %f %f\n", pt[0], pt[1], pt[2]);
812  }
813 }
814 
815 
816 
823 template<class T, class S>
824 inline void read_points(const std::string basename, const MPI_Comm comm, vector<S> & pts, vector<T> & ptsidx)
825 {
826  int size, rank;
827  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
828 
829  FILE* pts_fd = NULL;
830  size_t gnumpts = 0;
831  bool read_binary = false;
832 
833  // open the file on root and read number of points
834  if(rank == 0) {
835  read_binary = fileExists(basename + ".bpts");
836  std::string pts_file = read_binary ? basename + ".bpts" : basename + ".pts";
837 
838  pts_fd = fopen(pts_file.c_str(), "r");
839  if(!pts_fd) {
840  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", pts_file.c_str());
841  exit(1);
842  }
843 
844  // read number of points
845  char buffer[2048];
846  if(read_binary) fread(buffer, sizeof(char), HDR_SIZE, pts_fd);
847  else fgets( buffer, 2048, pts_fd);
848  sscanf(buffer, "%lu", &gnumpts);
849  }
850  MPI_Bcast(&gnumpts, sizeof(size_t), MPI_BYTE, 0, comm);
851 
852  // compute the size of the block to read .. this could also be set to a constant
853  size_t blocksize = (gnumpts + size - 1) / size;
854 
855  if(rank == 0) {
856  // first read own block
857  read_pts_block(pts_fd, read_binary, blocksize, pts);
858 
859  // now read remaining blocks
860  vector<S> buff;
861  for(int pid = 1; pid < size; pid++) {
862  read_pts_block(pts_fd, read_binary, blocksize, buff);
863  long int numsend = buff.size();
864 
865  MPI_Send(&numsend, 1, MPI_LONG, pid, SF_MPITAG, comm);
866  MPI_Send(buff.data(), numsend*sizeof(S), MPI_BYTE, pid, SF_MPITAG, comm);
867  }
868  }
869  else {
870  MPI_Status stat;
871  long int numrecv = 0;
872  MPI_Recv(&numrecv, 1, MPI_LONG, 0, SF_MPITAG, comm, &stat);
873 
874  pts.resize(numrecv);
875  MPI_Recv(pts.data(), numrecv*sizeof(S), MPI_BYTE, 0, SF_MPITAG, comm, &stat);
876  }
877 
878  long int mysize = pts.size() / 3;
879  vector<long int> layout;
880  layout_from_count(mysize, layout, comm);
881 
882  interval(ptsidx, layout[rank], layout[rank+1]);
883 
884  if(rank == 0) fclose(pts_fd);
885 }
886 
893 template<class T, class S>
894 inline void insert_points(const vector<S> & pts, const vector<T> & ptsidx,
895  std::list<meshdata<T, S>*> & meshlist)
896 {
897  // list must not be empty
898  assert(meshlist.size() > 0);
899  assert(pts.size() == (ptsidx.size() * 3));
900 
901  // we use the communicator provided by the first mesh in the list
902  // the code does not support a meshlist with different communicators
903  MPI_Comm comm = (*meshlist.begin())->comm;
904 
905  int size, rank;
906  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
907 
908  vector<S> ptsbuff;
909  vector<T> idxbuff;
910 
911  for(int pid = 0; pid < size; pid++) {
912  size_t numsend = pts.size();
913  MPI_Bcast(&numsend, sizeof(size_t), MPI_BYTE, pid, comm);
914 
915  if(rank == pid) {
916  ptsbuff = pts;
917  idxbuff = ptsidx;
918  }
919  else {
920  ptsbuff.resize(numsend);
921  idxbuff.resize(numsend / 3);
922  }
923 
924  MPI_Bcast(ptsbuff.data(), ptsbuff.size()*sizeof(S), MPI_BYTE, pid, comm);
925  MPI_Bcast(idxbuff.data(), idxbuff.size()*sizeof(T), MPI_BYTE, pid, comm);
926 
927  // iterate over all meshes in meshlist. if the interval of read indices contains indices of the mesh,
928  // we copy the respective nodes from the point buffer into the mesh
929  for(auto it = meshlist.begin(); it != meshlist.end(); ++it)
930  {
931  meshdata<T, S> & mesh = *(*it); // reference to current mesh in meshlist
932  const vector<T> & rnod = mesh.get_numbering(NBR_REF);
933 
934  // this function is used at the earliest stages of setting up a mesh,
935  // thus we cannot yet use the parallel layout functionalities
937  for(size_t i=0; i<rnod.size(); i++)
938  g2l[rnod[i]] = i;
939 
940  mesh.xyz.resize(mesh.l_numpts*3);
941 
942  vector<T> ridx, widx;
943  ridx.reserve(rnod.size()), widx.reserve(rnod.size());
944 
945  for(size_t j=0; j<idxbuff.size(); j++) {
946  if(g2l.count(idxbuff[j])) {
947  ridx.push_back(j);
948  widx.push_back(g2l[idxbuff[j]]);
949  }
950  }
951 
952  // insert those nodes which are in the local range
953  for(size_t j=0; j<ridx.size(); j++) {
954  T w = widx[j], r = ridx[j];
955  mesh.xyz[w*3+0] = ptsbuff[r*3+0];
956  mesh.xyz[w*3+1] = ptsbuff[r*3+1];
957  mesh.xyz[w*3+2] = ptsbuff[r*3+2];
958  }
959  }
960  }
961 }
962 
969 template<class T, class S>
970 inline void writeVTKmesh_binary(const meshdata<T, S> & mesh, std::string file)
971 {
972  FILE* vtk_file = fopen(file.c_str(), "w");
973  if(vtk_file == NULL) return;
974 
975 
976  fprintf (vtk_file, "# vtk DataFile Version 3.0\n");
977  fprintf (vtk_file, "vtk output\n");
978  fprintf (vtk_file, "binary\n");
979  fprintf (vtk_file, "DATASET UNSTRUCTURED_GRID\n\n");
980  fprintf (vtk_file, "POINTS %lu float\n", mesh.l_numpts);
981 
982  float pts[3];
983  const S* p = mesh.xyz.data();
984 
985  for (unsigned long int i=0; i<mesh.l_numpts; i++ ) {
986  pts[0] = htobe(p[0]);
987  pts[1] = htobe(p[1]);
988  pts[2] = htobe(p[2]);
989  fwrite(pts, sizeof(float), 3, vtk_file);
990  p += 3;
991  }
992  fprintf(vtk_file, "\n");
993 
994  fprintf (vtk_file, "CELL_TYPES %lu\n", mesh.l_numelem);
995  unsigned long int valcount = 0;
996  int vtk_type;
997  for(unsigned long int i=0; i< mesh.l_numelem; i++) {
998  elem_t etype = mesh.type[i];
999  switch(etype) {
1000  // Line for Purkinje
1001  case Line:
1002  vtk_type = 3; // Lines are encoded as index 3
1003  break;
1004  case Tri:
1005  vtk_type = 5; // Triangles are encoded as index 5
1006  break;
1007  case Tetra:
1008  vtk_type = 10; // Tetras are encoded as index 10
1009  break;
1010  case Quad:
1011  vtk_type = 9; // Quads are encoded as index 9
1012  break;
1013  case Pyramid:
1014  vtk_type = 14; // Pyramids are encoded as index 14
1015  break;
1016  case Prism:
1017  vtk_type = 13; // Prisms are encoded as index 13
1018  break;
1019  case Hexa:
1020  vtk_type = 12; // Hexahedras are encoded as index 12
1021  break;
1022  default: break;
1023  }
1024  vtk_type = htobe(vtk_type);
1025  fwrite(&vtk_type, sizeof(int), 1, vtk_file);
1026  valcount += mesh.dsp[i+1] - mesh.dsp[i] + 1;
1027  }
1028  fprintf(vtk_file, "\n");
1029 
1030  fprintf(vtk_file, "CELLS %lu %lu\n", mesh.l_numelem, valcount);
1031  const T* elem = mesh.con.data();
1032  for(unsigned long int i=0; i<mesh.l_numelem; i++)
1033  {
1034  int nodes = mesh.dsp[i+1] - mesh.dsp[i], n;
1035  int be_nodes = htobe(nodes);
1036  fwrite(&be_nodes, sizeof(int), 1, vtk_file);
1037 
1038  for(int j=0; j<nodes; j++) {
1039  n = elem[j]; n = htobe(n);
1040  fwrite(&n, sizeof(int), 1, vtk_file);
1041  }
1042  elem += nodes;
1043  }
1044  fprintf(vtk_file, "\n");
1045 
1046  fprintf (vtk_file, "CELL_DATA %lu \n", mesh.l_numelem);
1047  fprintf (vtk_file, "SCALARS elemTag int 1\n");
1048  fprintf (vtk_file, "LOOKUP_TABLE default\n");
1049  for (unsigned long int i=0; i<mesh.l_numelem; i++ ) {
1050  int t = htobe(int(mesh.tag[i]));
1051  fwrite(&t, sizeof(int), 1, vtk_file);
1052  }
1053 
1054  // write fiber data
1055  int nFib = 0;
1056  if(mesh.fib.size() > 0) {
1057  nFib = 1;
1058  if(mesh.she.size() == mesh.fib.size())
1059  nFib = 2;
1060  }
1061 
1062  fprintf (vtk_file, "VECTORS fiber float\n");
1063  for (unsigned long int i=0; i<mesh.l_numelem; i++ ) {
1064  pts[0] = htobe(float(mesh.fib[i*3+0])),
1065  pts[1] = htobe(float(mesh.fib[i*3+1])),
1066  pts[2] = htobe(float(mesh.fib[i*3+2]));
1067  fwrite(pts, sizeof(float), 3, vtk_file);
1068  }
1069  fprintf(vtk_file, "\n");
1070 
1071  if(nFib == 2) {
1072  fprintf (vtk_file, "VECTORS sheet float\n");
1073  for (unsigned long int i=0; i<mesh.l_numelem; i++ ) {
1074  pts[0] = htobe(float(mesh.she[i*3+0])),
1075  pts[1] = htobe(float(mesh.she[i*3+1])),
1076  pts[2] = htobe(float(mesh.she[i*3+2]));
1077  fwrite(pts, sizeof(float), 3, vtk_file);
1078  }
1079  fprintf(vtk_file, "\n");
1080  }
1081 
1082  fclose(vtk_file);
1083 }
1084 
1085 }
1086 #endif
Basic containers.
#define SF_MPITAG
the MPI tag when communicating
Definition: SF_globals.h:29
#define MT_ENDIANNESS
Definition: SF_io_base.h:9
#define htobe(x)
Definition: SF_io_base.h:5
#define HDR_SIZE
Definition: SF_mesh_io.h:44
Functions related to network communication.
The vector class and related algorithms.
The mesh storage class. It contains both element and vertex data.
Definition: SF_container.h:396
vector< T > dsp
connectivity starting index of each element
Definition: SF_container.h:416
vector< S > she
sheet direction
Definition: SF_container.h:421
vector< S > fib
fiber direction
Definition: SF_container.h:420
size_t l_numelem
local number of elements
Definition: SF_container.h:399
vector< elem_t > type
element type
Definition: SF_container.h:418
vector< T > & register_numbering(SF_nbr nbr_type)
Register a new numbering to the mesh and return the associated index vector.
Definition: SF_container.h:444
vector< T > con
Definition: SF_container.h:412
vector< S > xyz
node cooridnates
Definition: SF_container.h:427
size_t l_numpts
local number of points
Definition: SF_container.h:401
size_t g_numelem
global number of elements
Definition: SF_container.h:398
void localize(SF_nbr nbr_type)
Localize the connectivity data w.r.t. a given numbering.
Definition: SF_container.h:496
MPI_Comm comm
the parallel mesh is defined on a MPI world
Definition: SF_container.h:404
vector< T > & get_numbering(SF_nbr nbr_type)
Get the vector defining a certain numbering.
Definition: SF_container.h:464
vector< T > tag
element tag
Definition: SF_container.h:417
A vector storing arbitrary data.
Definition: SF_vector.h:43
size_t size() const
The current size of the vector.
Definition: SF_vector.h:104
void resize(size_t n)
Resize a vector.
Definition: SF_vector.h:209
void reserve(size_t n)
Definition: SF_vector.h:241
T * data()
Pointer to the vector's start.
Definition: SF_vector.h:91
T & push_back(T val)
Definition: SF_vector.h:283
hm_int count(const K &key) const
Check if key exists.
Definition: hashmap.hpp:579
Definition: dense_mat.hpp:34
void read_pts_block(FILE *&fd, bool read_binary, size_t bsize, vector< S > &xyz)
Read a chunk of points from a file descriptor.
Definition: SF_mesh_io.h:760
void read_headers(FILE *ele_fd, FILE *fib_fd, bool read_binary, size_t &numelem, int &nFib)
Read the header from the element and fiber files.
Definition: SF_mesh_io.h:90
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.
Definition: SF_mesh_io.h:824
void read_fib_block(FILE *&fd, bool read_binary, int nFib, size_t bsize, meshdata< T, S > &mesh)
Read a chunk of fibers from a file descriptor.
Definition: SF_mesh_io.h:395
void write_pts_block(FILE *&fd, bool write_binary, const vector< S > &xyz)
Write a chunk of points to a file.
Definition: SF_mesh_io.h:797
void interval(vector< T > &vec, size_t start, size_t end)
Create an integer interval between start and end.
Definition: SF_vector.h:350
void write_surface(const meshdata< T, S > &surfmesh, std::string surffile)
Definition: SF_mesh_io.h:712
void write_elements(const meshdata< T, S > &mesh, bool binary, std::string basename)
Read the element data (elements and fibers) of a CARP mesh.
Definition: SF_mesh_io.h:646
void writeVTKmesh_binary(const meshdata< T, S > &mesh, std::string file)
Write a mesh in binary vtk format.
Definition: SF_mesh_io.h:970
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.
Definition: SF_mesh_io.h:894
void write_pts_header(FILE *&pts_fd, bool binary, size_t numpts)
Write the header of the points file.
Definition: SF_mesh_io.h:182
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.
Definition: SF_mesh_io.h:209
void write_elem_headers(FILE *&ele_fd, FILE *&fib_fd, bool binary, size_t numelem, int nFib)
Write the header of the element and fiber files.
Definition: SF_mesh_io.h:145
void write_fib_block(FILE *&fd, bool write_binary, const meshdata< T, S > &mesh)
Write the local chunk of fibers to a file.
Definition: SF_mesh_io.h:465
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.
Definition: SF_vector.h:371
bool fileExists(std::string filename)
Function which checks if a given file exists.
Definition: SF_io_base.h:66
elem_t getElemTypeID(char *eletype)
Generate element type enum from string.
Definition: SF_container.h:167
void layout_from_count(const T count, vector< T > &layout, MPI_Comm comm)
Definition: SF_network.h:201
size_t read_num_pts(std::string basename)
Function returns the number of points in a CARP points file.
Definition: SF_mesh_io.h:55
elem_t
element type enum
Definition: SF_container.h:53
@ Line
Definition: SF_container.h:61
@ Tri
Definition: SF_container.h:60
@ Prism
Definition: SF_container.h:58
@ Pyramid
Definition: SF_container.h:57
@ Tetra
Definition: SF_container.h:54
@ Quad
Definition: SF_container.h:59
@ Hexa
Definition: SF_container.h:55
void read_elements(meshdata< T, S > &mesh, std::string basename, bool require_fibers=true)
Read the element data (elements and fibers) of a CARP mesh.
Definition: SF_mesh_io.h:518
@ NBR_ELEM_REF
The element numbering of the reference mesh (the one stored on HD).
Definition: SF_container.h:204
@ NBR_REF
The nodal numbering of the reference mesh (the one stored on HD).
Definition: SF_container.h:201
void write_elem_block(FILE *fd, bool write_binary, const meshdata< T, S > &mesh)
Write the local element block to a file.
Definition: SF_mesh_io.h:313