openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
fem_utils.cc
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 #include "SF_init.h"
28 #include "fem.h"
29 
30 namespace opencarp {
31 
32 void parse_comment_line(char* buff, const int buffsize, std::map<std::string,std::string> & metadata)
33 {
34  // we skip comment lines with no key value pair
35  if(! has_char(buff, buffsize, '=')) return;
36 
37  // there might be multiple comment characters and spaces preceding a key-value pair
38  remove_preceding_char(buff, buffsize, COMMENT_CHAR);
39  remove_preceding_char(buff, buffsize, ' ');
40 
41  const int wsize = 1024;
42  char lword[wsize], cword[wsize], rword[wsize];
43 
44  int numread = sscanf(buff, "%s %s %s", lword, cword, rword);
45 
46  if(numread == 3 && cword[0] == '=') {
47  // input was perfect
48  std::string key = lword, value = rword;
49  metadata[key] = value;
50  }
51  else if(numread == 1) {
52  // we have key=value without spaces
53  int ridx=0, widx=0;
54  while(ridx < wsize && lword[ridx] != '=') cword[widx++] = lword[ridx++];
55  cword[widx] = '\0';
56 
57  std::string key = cword;
58 
59  widx=0, ridx += 1;
60  while(ridx < wsize && lword[ridx] != '\0') cword[widx++] = lword[ridx++];
61  cword[widx] = '\0';
62 
63  std::string value = cword;
64 
65  metadata[key] = value;
66  }
67  else {
68  fprintf(stderr, "%s warning: Malformed key-value pair:\n%s\n", __func__, buff);
69  }
70 }
71 
72 void read_metadata(const std::string filename, std::map<std::string,std::string> & metadata, MPI_Comm comm)
73 {
74  // There are 3 ways to add metadata to a vtx file:
75  // 1) the new way:
76  // # key1 = value1
77  // # key2 = value2
78  // # key3 = value3
79  // 2) nvtx # mesh [intra|extra] nelem npts
80  // 3) number of nodes
81  // [intra|extra]
82  int size, rank;
83  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
84 
85  int err = 0;
86  FILE_SPEC stream = NULL;
87 
88  if(rank == 0) {
89  stream = f_open(filename.c_str(), "r");
90  if(stream == NULL) err++;
91  }
92 
93  // return on file open error
94  if(get_global(err, MPI_SUM)) {
95  log_msg(0,5,0, "%s error: Could not read vtx data. Aborting!", __func__);
96  EXIT(1);
97  }
98 
99  const size_t str_size = 5000;
100  char readbuff[str_size];
101 
102  // parse type 1 metadata lines. for simplicity we parse on each rank.
103  char* ptr = f_gets_par(readbuff, str_size, stream, comm);
104  remove_preceding_char(readbuff, str_size, ' ');
105 
106  while(ptr != NULL && readbuff[0] == COMMENT_CHAR) {
107  parse_comment_line(readbuff, str_size, metadata);
108  ptr = f_gets_par(readbuff, str_size, stream, comm);
109  remove_preceding_char(readbuff, str_size, ' ');
110  }
111 
112  // parse type 2 metadata
113  if(has_char(readbuff, str_size, '#')) {
114  int ibuff, ne, nn;
115  char iestr[128], nestr[128], nnstr[128];
116  int nread = sscanf(readbuff, "%d # %*s %31s %d %d", &ibuff, iestr, &ne, &nn);
117 
118  if(nread == 4) {
119  metadata["grid"] = iestr;
120  metadata["nelem"] = std::to_string(ne);
121  metadata["nnode"] = std::to_string(nn);
122  }
123  }
124 
125  // parse type 3 metadata
126  ptr = f_gets_par(readbuff, str_size, stream, comm);
127 
128  if(ptr) {
129  remove_preceding_char(readbuff, str_size, ' ');
130  char grid[128];
131  sscanf(readbuff, "%s", grid);
132 
133  if(strcmp(grid, "intra") == 0) {
134  metadata["grid"] = "intra";
135  } else if(strcmp(grid, "extra") == 0) {
136  metadata["grid"] = "extra";
137  }
138  }
139 
140  f_close(stream);
141 }
142 
143 char* skip_comments(FILE_SPEC stream, char* readbuff, size_t buffsize, MPI_Comm comm)
144 {
145  char* ptr = NULL;
146  do {
147  ptr = f_gets_par(readbuff, buffsize, stream, comm);
148  remove_preceding_char(readbuff, buffsize, ' ');
149  }
150  while(ptr != NULL && readbuff[0] == COMMENT_CHAR);
151 
152  return ptr;
153 }
154 
155 void indices_from_region_tag(SF::vector<mesh_int_t> & idx, const sf_mesh & mesh, const int tag)
156 {
157  std::set<mesh_int_t> idx_set;
158 
159  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
160  if(mesh.tag[eidx] == tag) {
161  for(mesh_int_t i = mesh.dsp[eidx]; i<mesh.dsp[eidx + 1]; i++)
162  idx_set.insert(mesh.con[i]);
163  }
164  }
165 
166  idx.assign(idx_set.begin(), idx_set.end());
167 }
168 
170 {
172  idx_set.reserve(mesh.l_numpts); // optional heuristic
173 
174  for (size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
175  if (tags.count(mesh.tag[eidx])) {
176  for (mesh_int_t i = mesh.dsp[eidx]; i < mesh.dsp[eidx + 1]; i++)
177  idx_set.insert(mesh.con[i]);
178  }
179  }
180 
181  idx.assign(idx_set.begin(), idx_set.end());
182 }
183 
184 void indices_from_geom_shape(SF::vector<mesh_int_t> & idx, const sf_mesh & mesh, const geom_shape shape, const bool nodal)
185 {
186  Point p;
187 
188  if(nodal) {
189  idx.reserve(mesh.l_numpts);
190 
191  for(size_t i=0; i<mesh.l_numpts; i++) {
192  p.x = mesh.xyz[i*3 + 0], p.y = mesh.xyz[i*3 + 1], p.z = mesh.xyz[i*3 + 2];
193  if(point_in_shape(p, shape))
194  idx.push_back(i);
195  }
196  }
197  else {
198  idx.reserve(mesh.l_numelem);
199 
200  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
201  p.x = p.y = p.z = 0.0;
202 
203  for(mesh_int_t i = mesh.dsp[eidx]; i<mesh.dsp[eidx + 1]; i++) {
204  mesh_int_t c = mesh.con[i];
205  p.x += mesh.xyz[c*3 + 0];
206  p.y += mesh.xyz[c*3 + 1];
207  p.z += mesh.xyz[c*3 + 2];
208  }
209  p /= double(mesh.dsp[eidx + 1] - mesh.dsp[eidx]);
210 
211  if(point_in_shape(p, shape))
212  idx.push_back(eidx);
213  }
214  }
215 }
216 
218 {
219  assert(mass.mesh_ptr() != NULL);
220 
221  sf_vec* weights; SF::init_vector(&weights, *mass.mesh_ptr(), 1, sf_vec::algebraic);
222  sf_vec* vols; SF::init_vector(&vols, *mass.mesh_ptr(), 1, sf_vec::algebraic);
223 
224  const SF::vector<mesh_int_t> & petsc_nod = mass.mesh_ptr()->get_numbering(SF::NBR_PETSC);
225 
226  // generate global petsc indices
227  SF::vector<SF_int> petsc_idx(local_idx.size());
228  for(size_t i=0; i<local_idx.size(); i++) petsc_idx[i] = petsc_nod[local_idx[i]];
229 
230  weights->set(0);
231  weights->set(petsc_idx, 1.0);
232 
233  mass.mult(*weights, *vols);
234  SF_real V = vols->sum();
235 
236  return V;
237 }
238 
239 void warn_when_passing_intra_vtx(const std::string filename)
240 {
241  std::map<std::string,std::string> metadata;
242  read_metadata(filename, metadata, PETSC_COMM_WORLD);
243 
244  if(metadata.count("grid")) {
245  if(metadata["grid"].compare("intra") == 0) {
246  log_msg(0,3,0, "Warning: openCARP requires input vtx indices to be of from the input mesh, i.e. of type \"extra\".");
247  }
248  } else {
249  log_msg(0,3,0, "%s warning: Could not derive grid type info from file %s", __func__, filename.c_str());
250  }
251 }
252 
253 } // namespace opencarp
254 
int mesh_int_t
Definition: SF_container.h:46
double SF_real
Use the general double as real type.
Definition: SF_globals.h:38
const meshdata< mesh_int_t, mesh_real_t > * mesh_ptr() const
virtual void mult(const abstract_vector< T, S > &x, abstract_vector< T, S > &b) const =0
virtual S sum() const =0
virtual void set(const vector< T > &idx, const vector< S > &vals, const bool additive=false, const bool local=false)=0
vector< T > dsp
connectivity starting index of each element
Definition: SF_container.h:416
size_t l_numelem
local number of elements
Definition: SF_container.h:399
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
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
size_t size() const
The current size of the vector.
Definition: SF_vector.h:104
void assign(InputIterator s, InputIterator e)
Assign a memory range.
Definition: SF_vector.h:161
void reserve(size_t n)
Definition: SF_vector.h:241
T & push_back(T val)
Definition: SF_vector.h:283
void reserve(size_t n)
Definition: hashmap.hpp:1065
hm_int count(const K &key) const
Definition: hashmap.hpp:992
void insert(InputIterator first, InputIterator last)
Definition: hashmap.hpp:962
class to store shape definitions
Definition: basics.h:381
Top-level header of FEM module.
#define COMMENT_CHAR
Definition: fem_utils.h:37
void init_vector(SF::abstract_vector< T, S > **vec)
Definition: SF_init.h:99
@ NBR_PETSC
PETSc numbering of nodes.
Definition: SF_container.h:203
bool has_char(char *buff, const int buffsize, const char c)
Definition: basics.h:368
char * f_gets_par(char *s, int size, FILE_SPEC stream, MPI_Comm comm)
Definition: basics.cc:199
void parse_comment_line(char *buff, const int buffsize, std::map< std::string, std::string > &metadata)
Definition: fem_utils.cc:32
bool point_in_shape(const Point &p, const geom_shape &shape)
test if a point is inside a simple geometric shape
Definition: basics.cc:250
SF_real get_volume_from_nodes(sf_mat &mass, SF::vector< mesh_int_t > &local_idx)
Definition: fem_utils.cc:217
char * skip_comments(FILE_SPEC stream, char *readbuff, size_t buffsize, MPI_Comm comm)
Definition: fem_utils.cc:143
void read_metadata(const std::string filename, std::map< std::string, std::string > &metadata, MPI_Comm comm)
Read metadata from the header.
Definition: fem_utils.cc:72
T get_global(T in, MPI_Op OP, MPI_Comm comm=PETSC_COMM_WORLD)
Do a global reduction on a variable.
Definition: basics.h:233
void indices_from_region_tags(SF::vector< mesh_int_t > &idx, const sf_mesh &mesh, const hashmap::unordered_set< int > &tags)
Populate vertex data with the vertices of multiple tag regions.
Definition: fem_utils.cc:169
FILE_SPEC f_open(const char *fname, const char *mode)
Open a FILE_SPEC.
Definition: basics.cc:135
void warn_when_passing_intra_vtx(const std::string filename)
Definition: fem_utils.cc:239
void indices_from_region_tag(SF::vector< mesh_int_t > &idx, const sf_mesh &mesh, const int tag)
Populate vertex data with the vertices of a given tag region.
Definition: fem_utils.cc:155
void indices_from_geom_shape(SF::vector< mesh_int_t > &idx, const sf_mesh &mesh, const geom_shape shape, const bool nodal)
Populate vertex data with the vertices inside a defined box shape.
Definition: fem_utils.cc:184
void log_msg(FILE_SPEC out, int level, unsigned char flag, const char *fmt,...)
Definition: basics.cc:72
void f_close(FILE_SPEC &f)
Close a FILE_SPEC.
Definition: basics.cc:162
void remove_preceding_char(char *buff, const int buffsize, const char c)
Definition: basics.h:345
File descriptor struct.
Definition: basics.h:133