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  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);
306  exit(EXIT_FAILURE);
307  }
308 }
309 
310 
318 template<class T, class S>
319 inline void write_elem_block(FILE* fd, bool write_binary, const meshdata<T, S> & mesh)
320 {
321  const T* con = mesh.con.data();
322 
323  if(write_binary)
324  {
325  int wbuff[9];
326 
327  for(size_t eidx=0; eidx < mesh.l_numelem; eidx++)
328  {
329  fwrite(&mesh.type[eidx], 1, sizeof(elem_t), fd);
330  T esize = mesh.dsp[eidx+1] - mesh.dsp[eidx];
331  vec_assign(wbuff, con, esize); // copy-convert to int
332  wbuff[esize] = mesh.tag[eidx]; // convert to int
333 
334  fwrite(wbuff, esize+1, sizeof(int), fd);
335  con += esize;
336  }
337  }
338  else
339  {
340  for(size_t eidx=0; eidx < mesh.l_numelem; eidx++)
341  {
342  switch(mesh.type[eidx])
343  {
344  case Line:
345  fprintf(fd, "Ln %d %d %d\n", con[0], con[1], mesh.tag[eidx]);
346  con += 2;
347  break;
348 
349  case Tri:
350  fprintf(fd, "Tr %d %d %d %d\n", con[0], con[1], con[2], mesh.tag[eidx]);
351  con += 3;
352  break;
353 
354  case Quad:
355  fprintf(fd, "Qd %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 Tetra:
360  fprintf(fd, "Tt %d %d %d %d %d\n", con[0], con[1], con[2], con[3], mesh.tag[eidx]);
361  con += 4;
362  break;
363 
364  case Pyramid:
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]);
367  con += 5;
368  break;
369 
370  case Prism:
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]);
373  con += 6;
374  break;
375 
376  case Hexa:
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]);
379  con += 8;
380  break;
381 
382  default:
383  fprintf(stderr, "Error: Unsupported element type!\n");
384  exit(1);
385  }
386  }
387  }
388 }
389 
390 
400 template<class T, class S>
401 inline void read_fib_block(FILE* & fd, bool read_binary, bool twoFib, size_t bsize, meshdata<T, S> & mesh)
402 {
403  const int bufsize = 2048;
404  char buffer[bufsize];
405  int err = 0;
406  size_t nr = 0;
407 
408  if(twoFib) {
409  float fib[6];
410 
411  mesh.fib.resize(bsize*3);
412  mesh.she.resize(bsize*3);
413 
414  for(size_t i=0; i<bsize; i++) {
415  if(read_binary) {
416  size_t r = fread(fib, sizeof(float), 6, fd);
417  if(r != 6) {
418  fprintf(stderr, "%s error: Orthotropy blon file does not contain two fiber data for all elements. Aborting!\n", __func__);
419  err++;
420  break;
421  }
422  }
423  else {
424  char* ptr = fgets(buffer, bufsize, fd);
425  if(ptr == NULL) break;
426 
427  size_t r = sscanf(buffer, "%f %f %f %f %f %f", fib, fib+1, fib+2, fib+3, fib+4, fib+5);
428  if(r != 6) {
429  fprintf(stderr, "%s error: Orthotropy lon file does not contain data two fiber for all elements. Aborting!\n", __func__);
430  err++;
431  break;
432  }
433  }
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];
440 
441  nr++;
442  }
443  mesh.fib.resize(nr*3);
444  mesh.she.resize(nr*3);
445  } else {
446  float fib[3];
447 
448  mesh.fib.resize(bsize*3);
449 
450  for(size_t i=0; i<bsize; i++) {
451  if(read_binary) {
452  size_t r = fread(fib, sizeof(float), 3, fd);
453  if(r != 3) {
454  fprintf(stderr, "Error: Orthotropy blon file does not contain data for all elements. i=%i, bsize=%i Aborting!\n", i, bsize);
455  err++;
456  break;
457  }
458  }
459  else {
460  char* ptr = fgets( buffer, bufsize, fd);
461  if(ptr == NULL) break;
462 
463  size_t r = sscanf(buffer, "%f %f %f", fib, fib+1, fib+2);
464  if(r != 3) {
465  fprintf(stderr, "Error: Orthotropy lon file does not contain data for all elements. Aborting!\n");
466  err++;
467  break;
468  }
469  }
470  mesh.fib[nr*3+0] = fib[0];
471  mesh.fib[nr*3+1] = fib[1];
472  mesh.fib[nr*3+2] = fib[2];
473 
474  nr++;
475  }
476  mesh.fib.resize(nr*3);
477  }
478 
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);
481  err++;
482  }
483 
484  if(err) {
485  exit(EXIT_FAILURE);
486  }
487 }
488 
489 
497 template<class T, class S>
498 inline void write_fib_block(FILE* & fd, bool write_binary, const meshdata<T, S> & mesh)
499 {
500  float fib[6];
501  bool twoFib = mesh.she.size() == mesh.fib.size();
502 
503  if(twoFib)
504  {
505  for(size_t i=0; i<mesh.l_numelem; i++)
506  {
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];
513 
514  if(write_binary)
515  fwrite(fib, 6, sizeof(float), fd);
516  else
517  fprintf(fd, "%f %f %f %f %f %f\n", fib[0], fib[1], fib[2], fib[3], fib[4], fib[5]);
518  }
519  }
520  else
521  {
522  for(size_t i=0; i<mesh.l_numelem; i++)
523  {
524  fib[0] = mesh.fib[i*3+0];
525  fib[1] = mesh.fib[i*3+1];
526  fib[2] = mesh.fib[i*3+2];
527 
528  if(write_binary)
529  fwrite(fib, 3, sizeof(float), fd);
530  else
531  fprintf(fd, "%f %f %f\n", fib[0], fib[1], fib[2]);
532  }
533  }
534 }
535 
545 template<class T, class S>
546 inline void read_elements(meshdata<T, S> & mesh, std::string basename)
547 {
548  // we use the communicator provided by the mesh
549  MPI_Comm comm = mesh.comm;
550 
551  int size, rank;
552  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
553 
554  FILE* ele_fd = NULL, *fib_fd = NULL;
555  size_t gnumelems = 0;
556  bool twoFib = false;
557  bool read_binary = false;
558  int err = 0;
559 
560  // open the files on root and read number of elements
561  if(rank == 0) {
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";
565 
566  ele_fd = fopen(ele_file.c_str(), "r");
567  if(!ele_fd) {
568  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", ele_file.c_str());
569  err++;
570  }
571 
572  fib_fd = fopen(fib_file.c_str(), "r");
573  if(!fib_fd) {
574  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", fib_file.c_str());
575  err++;
576  }
577  }
578 
579  MPI_Allreduce(MPI_IN_PLACE, &err, 1, MPI_INT, MPI_SUM, comm);
580  if(err) {
581  exit(EXIT_FAILURE);
582  }
583 
584  if(rank == 0)
585  read_headers(ele_fd, fib_fd, read_binary, gnumelems, twoFib);
586 
587  MPI_Bcast(&gnumelems, sizeof(size_t), MPI_BYTE, 0, comm);
588  MPI_Bcast(&twoFib, sizeof(bool), MPI_BYTE, 0, comm);
589  mesh.g_numelem = gnumelems;
590 
591  // compute the size of the block stored on each process
592  SF::vector<size_t> block_cnt, block_dsp;
593  divide(gnumelems, size, block_cnt);
594  dsp_from_cnt(block_cnt, block_dsp);
595 
596  // rank 0 reads mesh and distributes
597  if(rank == 0) {
598  // read own block
599  read_elem_block(ele_fd, read_binary, 0, block_cnt[0], mesh);
600  read_fib_block (fib_fd, read_binary, twoFib, block_cnt[0], mesh);
601 
602  // read blocks of other ranks and communicate
603  meshdata<T, S> meshbuff;
604  for(int pid=1; pid<size; pid++) {
605  // read block in buffer
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);
608 
609  vector<T> & ref_eidx = meshbuff.get_numbering(NBR_ELEM_REF);
610 
611  // communicate
612  MPI_Send(&meshbuff.l_numelem, sizeof(size_t), MPI_BYTE, pid, SF_MPITAG, comm);
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);
616  MPI_Send(meshbuff.type.data(), meshbuff.type.size()*sizeof(elem_t),MPI_BYTE, pid, SF_MPITAG, comm);
617  MPI_Send(meshbuff.fib.data(), meshbuff.fib.size()*sizeof(S), MPI_BYTE, pid, SF_MPITAG, comm);
618 
619  if(twoFib)
620  MPI_Send(meshbuff.she.data(), meshbuff.she.size()*sizeof(S), MPI_BYTE, pid, SF_MPITAG, comm);
621 
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);
625  }
626 
627  fclose(ele_fd);
628  fclose(fib_fd);
629  }
630  else {
631  // the rest receives
632  MPI_Status stat;
633  MPI_Recv(&mesh.l_numelem, sizeof(size_t), MPI_BYTE, 0, SF_MPITAG, comm, &stat);
634 
635  vector<T> & ref_eidx = mesh.register_numbering(NBR_ELEM_REF);
636  ref_eidx.resize(mesh.l_numelem);
637 
638  mesh.dsp.resize(mesh.l_numelem+1);
639  mesh.tag.resize(mesh.l_numelem);
640  mesh.type.resize(mesh.l_numelem);
641 
642  mesh.fib.resize(mesh.l_numelem*3);
643  if(twoFib) mesh.she.resize(mesh.l_numelem*3);
644 
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);
648  MPI_Recv(mesh.type.data(), mesh.type.size()*sizeof(elem_t),MPI_BYTE, 0, SF_MPITAG, comm, &stat);
649  MPI_Recv(mesh.fib.data(), mesh.fib.size()*sizeof(S), MPI_BYTE, 0, SF_MPITAG, comm, &stat);
650  if(twoFib)
651  MPI_Recv(mesh.she.data(), mesh.she.size()*sizeof(S), MPI_BYTE, 0, SF_MPITAG, comm, &stat);
652 
653  size_t con_size;
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);
657  }
658 
659  // localize connectivity (w.r.t. reference numbering)
660  mesh.localize(NBR_REF);
661 }
662 
663 
673 template<class T, class S>
674 inline void write_elements(const meshdata<T, S> & mesh, bool binary, std::string basename)
675 {
676  const MPI_Comm comm = mesh.comm;
677 
678  int size, rank;
679  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
680 
681  FILE* ele_fd = NULL, *fib_fd = NULL;
682  bool twoFib = mesh.she.size() == mesh.fib.size();
683  std::string ele_file = binary ? basename + ".belem" : basename + ".elem";
684  std::string fib_file = binary ? basename + ".blon" : basename + ".lon";
685 
686  // open the files on root and write headers
687  if(rank == 0) {
688  ele_fd = fopen(ele_file.c_str(), "w");
689  if(!ele_fd) {
690  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", ele_file.c_str());
691  exit(1);
692  }
693  fib_fd = fopen(fib_file.c_str(), "w");
694  if(!fib_fd) {
695  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", fib_file.c_str());
696  fclose(ele_fd);
697  exit(1);
698  }
699  write_elem_headers(ele_fd, fib_fd, binary, mesh.g_numelem, twoFib);
700  fclose(ele_fd);
701  fclose(fib_fd);
702  }
703 
704  // write mesh sequentially rank after rank
705  for(int pid=0; pid < size; pid++)
706  {
707  if(pid == rank)
708  {
709  ele_fd = fopen(ele_file.c_str(), "a");
710  if(!ele_fd) {
711  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", ele_file.c_str());
712  exit(1);
713  }
714  fib_fd = fopen(fib_file.c_str(), "a");
715  if(!fib_fd) {
716  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", fib_file.c_str());
717  fclose(ele_fd);
718  exit(1);
719  }
720 
721  write_elem_block(ele_fd, binary, mesh);
722  write_fib_block(fib_fd, binary, mesh);
723 
724  fclose(ele_fd);
725  fclose(fib_fd);
726  }
727  MPI_Barrier(comm);
728  }
729 }
730 
731 template<class T, class S>
732 inline void write_surface(const meshdata<T, S> & surfmesh, std::string surffile)
733 {
734  const MPI_Comm comm = surfmesh.comm;
735  const bool binary = false;
736 
737  int size, rank;
738  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
739 
740  FILE* ele_fd = NULL;
741 
742  // open the files on root and write headers
743  if(rank == 0) {
744  ele_fd = fopen(surffile.c_str(), "w");
745  if(!ele_fd) {
746  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", surffile.c_str());
747  exit(1);
748  }
749  fprintf(ele_fd, "%lu\n", surfmesh.g_numelem);
750  fclose(ele_fd);
751  }
752 
753  // write surface sequentially rank after rank
754  for(int pid=0; pid < size; pid++)
755  {
756  if(pid == rank) {
757  ele_fd = fopen(surffile.c_str(), "a");
758  if(!ele_fd) {
759  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", surffile.c_str());
760  exit(1);
761  }
762 
763  write_elem_block(ele_fd, binary, surfmesh);
764  fclose(ele_fd);
765  }
766  MPI_Barrier(comm);
767  }
768 }
769 
770 
779 template<class S>
780 inline void read_pts_block(FILE* & fd, bool read_binary, size_t bsize, vector<S> & xyz)
781 {
782  const int bufsize = 2048;
783  char buffer[bufsize];
784  float pts[3];
785  xyz.resize(bsize*3);
786 
787  size_t nr = 0;
788 
789  for(size_t i=0; i<bsize; i++) {
790  if(read_binary) {
791  size_t r = fread(pts, sizeof(float), 3, fd);
792  if(r != 3) break;
793  }
794  else {
795  char* ptr = fgets( buffer, bufsize, fd);
796  if(ptr == NULL) break;
797  sscanf(buffer, "%f %f %f", pts, pts+1, pts+2);
798  }
799  xyz[nr*3+0] = pts[0];
800  xyz[nr*3+1] = pts[1];
801  xyz[nr*3+2] = pts[2];
802 
803  nr++;
804  }
805  xyz.resize(nr*3);
806 }
807 
808 
816 template<class S>
817 inline void write_pts_block(FILE* & fd, bool write_binary, const vector<S> & xyz)
818 {
819  size_t nnodes = xyz.size() / 3;
820  float pt[3];
821 
822  for(size_t i=0; i<nnodes; i++)
823  {
824  pt[0] = xyz[i*3+0];
825  pt[1] = xyz[i*3+1];
826  pt[2] = xyz[i*3+2];
827 
828  if(write_binary)
829  fwrite(pt, 3, sizeof(float), fd);
830  else
831  fprintf(fd, "%f %f %f\n", pt[0], pt[1], pt[2]);
832  }
833 }
834 
835 
836 
843 template<class T, class S>
844 inline void read_points(const std::string basename, const MPI_Comm comm, vector<S> & pts, vector<T> & ptsidx)
845 {
846  int size, rank;
847  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
848 
849  FILE* pts_fd = NULL;
850  size_t gnumpts = 0;
851  bool read_binary = false;
852 
853  // open the file on root and read number of points
854  if(rank == 0) {
855  read_binary = fileExists(basename + ".bpts");
856  std::string pts_file = read_binary ? basename + ".bpts" : basename + ".pts";
857 
858  pts_fd = fopen(pts_file.c_str(), "r");
859  if(!pts_fd) {
860  fprintf(stderr, "Error: could not open file: %s. Aborting!\n", pts_file.c_str());
861  exit(1);
862  }
863 
864  // read number of points
865  char buffer[2048];
866  if(read_binary) fread(buffer, sizeof(char), HDR_SIZE, pts_fd);
867  else fgets( buffer, 2048, pts_fd);
868  sscanf(buffer, "%lu", &gnumpts);
869  }
870  MPI_Bcast(&gnumpts, sizeof(size_t), MPI_BYTE, 0, comm);
871 
872  // compute the size of the block to read .. this could also be set to a constant
873  size_t blocksize = (gnumpts + size - 1) / size;
874 
875  if(rank == 0) {
876  // first read own block
877  read_pts_block(pts_fd, read_binary, blocksize, pts);
878 
879  // now read remaining blocks
880  vector<S> buff;
881  for(int pid = 1; pid < size; pid++) {
882  read_pts_block(pts_fd, read_binary, blocksize, buff);
883  long int numsend = buff.size();
884 
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);
887  }
888  }
889  else {
890  MPI_Status stat;
891  long int numrecv = 0;
892  MPI_Recv(&numrecv, 1, MPI_LONG, 0, SF_MPITAG, comm, &stat);
893 
894  pts.resize(numrecv);
895  MPI_Recv(pts.data(), numrecv*sizeof(S), MPI_BYTE, 0, SF_MPITAG, comm, &stat);
896  }
897 
898  long int mysize = pts.size() / 3;
899  vector<long int> layout;
900  layout_from_count(mysize, layout, comm);
901 
902  interval(ptsidx, layout[rank], layout[rank+1]);
903 
904  if(rank == 0) fclose(pts_fd);
905 }
906 
907 
914 template<class T, class S>
915 inline void insert_points(const vector<S> & pts, const vector<T> & ptsidx,
916  std::list<meshdata<T, S>*> & meshlist)
917 {
918  // list must not be empty
919  assert(meshlist.size() > 0);
920  assert(pts.size() == (ptsidx.size() * 3));
921 
922  // we use the communicator provided by the first mesh in the list
923  // the code does not support a meshlist with different communicators
924  MPI_Comm comm = (*meshlist.begin())->comm;
925 
926  int size, rank;
927  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
928 
929  vector<S> ptsbuff;
930  vector<T> idxbuff;
931 
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);
935 
936  if(rank == pid) {
937  ptsbuff = pts;
938  idxbuff = ptsidx;
939  }
940  else {
941  ptsbuff.resize(numsend);
942  idxbuff.resize(numsend / 3);
943  }
944 
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);
947 
948  // iterate over all meshes in meshlist. if the interval of read indices contains indices of the mesh,
949  // we copy the respective nodes from the point buffer into the mesh
950  for(auto it = meshlist.begin(); it != meshlist.end(); ++it)
951  {
952  meshdata<T, S> & mesh = *(*it); // reference to current mesh in meshlist
953  const vector<T> & rnod = mesh.get_numbering(NBR_REF);
954 
955  // this function is used at the earliest stages of setting up a mesh,
956  // thus we cannot yet use the parallel layout functionalities
958  for(size_t i=0; i<rnod.size(); i++)
959  g2l[rnod[i]] = i;
960 
961  mesh.xyz.resize(mesh.l_numpts*3);
962 
963  vector<T> ridx, widx;
964  ridx.reserve(rnod.size()), widx.reserve(rnod.size());
965 
966  for(size_t j=0; j<idxbuff.size(); j++) {
967  if(g2l.count(idxbuff[j])) {
968  ridx.push_back(j);
969  widx.push_back(g2l[idxbuff[j]]);
970  }
971  }
972 
973  // insert those nodes which are in the local range
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];
979  }
980  }
981  }
982 }
983 
990 template<class T, class S>
991 inline void writeVTKmesh_binary(const meshdata<T, S> & mesh, std::string file)
992 {
993  FILE* vtk_file = fopen(file.c_str(), "w");
994  if(vtk_file == NULL) return;
995 
996 
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);
1002 
1003  float pts[3];
1004  const S* p = mesh.xyz.data();
1005 
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);
1011  p += 3;
1012  }
1013  fprintf(vtk_file, "\n");
1014 
1015  fprintf (vtk_file, "CELL_TYPES %lu\n", mesh.l_numelem);
1016  unsigned long int valcount = 0;
1017  int vtk_type;
1018  for(unsigned long int i=0; i< mesh.l_numelem; i++) {
1019  elem_t etype = mesh.type[i];
1020  switch(etype) {
1021  // Line for Purkinje
1022  case Line:
1023  vtk_type = 3; // Lines are encoded as index 3
1024  break;
1025  case Tri:
1026  vtk_type = 5; // Triangles are encoded as index 5
1027  break;
1028  case Tetra:
1029  vtk_type = 10; // Tetras are encoded as index 10
1030  break;
1031  case Quad:
1032  vtk_type = 9; // Quads are encoded as index 9
1033  break;
1034  case Pyramid:
1035  vtk_type = 14; // Pyramids are encoded as index 14
1036  break;
1037  case Prism:
1038  vtk_type = 13; // Prisms are encoded as index 13
1039  break;
1040  case Hexa:
1041  vtk_type = 12; // Hexahedras are encoded as index 12
1042  break;
1043  default: break;
1044  }
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;
1048  }
1049  fprintf(vtk_file, "\n");
1050 
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++)
1054  {
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);
1058 
1059  for(int j=0; j<nodes; j++) {
1060  n = elem[j]; n = htobe(n);
1061  fwrite(&n, sizeof(int), 1, vtk_file);
1062  }
1063  elem += nodes;
1064  }
1065  fprintf(vtk_file, "\n");
1066 
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++ ) {
1071  int t = htobe(int(mesh.tag[i]));
1072  fwrite(&t, sizeof(int), 1, vtk_file);
1073  }
1074 
1075  // write fiber data
1076  bool twofibers = mesh.she.size() == (mesh.l_numelem*3);
1077 
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);
1084  }
1085  fprintf(vtk_file, "\n");
1086 
1087  if(twofibers) {
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);
1094  }
1095  fprintf(vtk_file, "\n");
1096  }
1097 
1098  fclose(vtk_file);
1099 }
1100 
1101 }
1102 #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)
Write the element data (elements and fibers) of a CARP mesh.
Definition: SF_mesh_io.h:674
void writeVTKmesh_binary(const meshdata< T, S > &mesh, std::string file)
Write a mesh in binary vtk format.
Definition: SF_mesh_io.h:991
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:319
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:817
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
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
Definition: SF_vector.h:358
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:546
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:732
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:498
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:915
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:780
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:401
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 dsp_from_cnt(const vector< T > &cnt, vector< T > &dsp)
Compute displacements from counts.
Definition: SF_vector.h:310
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:844