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