MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fmsconvert.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "../config/config.hpp"
13 
14 #ifdef MFEM_USE_FMS
15 
16 #include "fmsconvert.hpp"
17 #include <climits>
18 
19 using std::endl;
20 
21 // #define DEBUG_FMS_MFEM 1
22 // #define DEBUG_MFEM_FMS 1
23 
24 namespace mfem
25 {
26 
27 static inline int
28 FmsBasisTypeToMfemBasis(FmsBasisType b)
29 {
30  int retval = -1;
31  switch (b)
32  {
33  case FMS_NODAL_GAUSS_OPEN:
35  break;
36  case FMS_NODAL_GAUSS_CLOSED:
38  break;
39  case FMS_POSITIVE:
40  retval = mfem::BasisType::Positive;;
41  break;
42  case FMS_NODAL_UNIFORM_OPEN:
44  break;
45  case FMS_NODAL_UNIFORM_CLOSED:
47  break;
48  case FMS_NODAL_CHEBYSHEV_OPEN:
49  case FMS_NODAL_CHEBYSHEV_CLOSED:
50  mfem::out <<
51  "FMS_NODAL_CHEBYSHEV_OPEN, FMS_NODAL_CHEBYSHEV_CLOSED need conversion to MFEM types."
52  << endl;
53  break;
54  }
55  return retval;
56 }
57 
58 // The following function is unused (for now), so it is commented out to
59 // suppress compilation warning.
60 #if 0
61 /**
62 @brief Get the order and layout of the field.
63 */
64 static int
65 FmsFieldGetOrderAndLayout(FmsField f, FmsInt *f_order, FmsLayoutType *f_layout)
66 {
67  int err = 0;
68  FmsFieldDescriptor fd;
69  FmsLayoutType layout;
70  FmsScalarType data_type;
71  const void *data = nullptr;
72  FmsInt order = 0;
73 
74  FmsFieldGet(f, &fd, NULL, &layout, &data_type,
75  &data);
76 
77  FmsFieldDescriptorType f_fd_type;
78  FmsFieldDescriptorGetType(fd, &f_fd_type);
79  if (f_fd_type != FMS_FIXED_ORDER)
80  {
81  err = 1;
82  }
83  else
84  {
85  FmsFieldType field_type;
86  FmsBasisType basis_type;
87  FmsFieldDescriptorGetFixedOrder(fd, &field_type,
88  &basis_type, &order);
89  }
90 
91  *f_order = order;
92  *f_layout = layout;
93 
94  return err;
95 }
96 #endif
97 
98 /**
99 @brief This function converts an FmsField to an MFEM GridFunction.
100 @note I took some of the Pumi example code from the mesh conversion function
101  that converted coordinates and am trying to make it more general.
102  Coordinates are just another field so it seems like a good starting
103  point. We still have to support a bunch of other function types, etc.
104 */
105 template <typename DataType>
106 int
107 FmsFieldToGridFunction(FmsMesh fms_mesh, FmsField f, Mesh *mesh,
108  GridFunction &func, bool setFE)
109 {
110  int err = 0;
111 
112  // NOTE: transplanted from the FmsMeshToMesh function
113  // We should do this work once and save it.
114  //--------------------------------------------------
115  FmsInt dim, n_vert, n_elem, space_dim;
116 
117  // Find the first component that has coordinates - that will be the new mfem
118  // mesh.
119  FmsInt num_comp;
120  FmsMeshGetNumComponents(fms_mesh, &num_comp);
121  FmsComponent main_comp = NULL;
122  FmsField coords = NULL;
123  for (FmsInt comp_id = 0; comp_id < num_comp; comp_id++)
124  {
125  FmsComponent comp;
126  FmsMeshGetComponent(fms_mesh, comp_id, &comp);
127  FmsComponentGetCoordinates(comp, &coords);
128  if (coords) { main_comp = comp; break; }
129  }
130  if (!main_comp) { return 1; }
131  FmsComponentGetDimension(main_comp, &dim);
132  FmsComponentGetNumEntities(main_comp, &n_elem);
133  FmsInt n_ents[FMS_NUM_ENTITY_TYPES];
134  FmsInt n_main_parts;
135  FmsComponentGetNumParts(main_comp, &n_main_parts);
136  for (FmsInt et = FMS_VERTEX; et < FMS_NUM_ENTITY_TYPES; et++)
137  {
138  n_ents[et] = 0;
139  for (FmsInt part_id = 0; part_id < n_main_parts; part_id++)
140  {
141  FmsInt num_ents;
142  FmsComponentGetPart(main_comp, part_id, (FmsEntityType)et, NULL, NULL,
143  NULL, NULL, &num_ents);
144  n_ents[et] += num_ents;
145  }
146  }
147  n_vert = n_ents[FMS_VERTEX];
148  //--------------------------------------------------
149 
150  // Interrogate the field.
151  FmsFieldDescriptor f_fd;
152  FmsLayoutType f_layout;
153  FmsScalarType f_data_type;
154  const void *f_data;
155  FmsFieldGet(f, &f_fd, &space_dim, &f_layout, &f_data_type,
156  &f_data);
157  // FmsFieldGet(coords, NULL, &space_dim, NULL, NULL, NULL);
158 
159  FmsInt f_num_dofs;
160  FmsFieldDescriptorGetNumDofs(f_fd, &f_num_dofs);
161 
162  // Access FMS data through this typed pointer.
163  auto src_data = reinterpret_cast<const DataType *>(f_data);
164 
165  FmsFieldDescriptorType f_fd_type;
166  FmsFieldDescriptorGetType(f_fd, &f_fd_type);
167  if (f_fd_type != FMS_FIXED_ORDER)
168  {
169  return 9;
170  }
171  FmsFieldType f_field_type;
172  FmsBasisType f_basis_type;
173  FmsInt f_order;
174  FmsFieldDescriptorGetFixedOrder(f_fd, &f_field_type,
175  &f_basis_type, &f_order);
176 
177  if (f_field_type != FMS_CONTINUOUS && f_field_type != FMS_DISCONTINUOUS &&
178  f_field_type != FMS_HDIV)
179  {
180  return 10;
181  }
182 
183  int btype = FmsBasisTypeToMfemBasis(f_basis_type);
184  if (btype < 0)
185  {
186  mfem::out << "\tInvalid BasisType: " << BasisType::Name(btype) << std::endl;
187  return 11;
188  }
189 
190  //------------------------------------------------------------------
191  if (setFE)
192  {
193  // We could assemble a name based on fe_coll.hpp rules and pass to
194  // FiniteElementCollection::New()
195 
196  mfem::FiniteElementCollection *fec = nullptr;
197  switch (f_field_type)
198  {
199  case FMS_DISCONTINUOUS:
200  fec = new L2_FECollection(f_order, dim, btype);
201  break;
202  case FMS_CONTINUOUS:
203  fec = new H1_FECollection(f_order, dim, btype);
204  break;
205  case FMS_HDIV:
206  fec = new RT_FECollection(f_order, dim);
207  break;
208  case FMS_HCURL:
209  case FMS_DISCONTINUOUS_WEIGHTED:
210  MFEM_ABORT("Field types FMS_HCURL and FMS_DISCONTINUOUS_WEIGHTED"
211  " are not supported yet.");
212  break;
213  }
214 
215  int ordering = (f_layout == FMS_BY_VDIM) ? Ordering::byVDIM : Ordering::byNODES;
216  auto fes = new FiniteElementSpace(mesh, fec, space_dim, ordering);
217  func.SetSpace(fes);
218  }
219  //------------------------------------------------------------------
220  const FmsInt nstride = (f_layout == FMS_BY_VDIM) ? space_dim : 1;
221  const FmsInt vstride = (f_layout == FMS_BY_VDIM) ? 1 : f_num_dofs;
222 
223  // Data reordering to store the data into func.
224  if ((FmsInt)(func.Size()) != f_num_dofs*space_dim)
225  {
226  return 12;
227  }
228 
229  mfem::FiniteElementSpace *fes = func.FESpace();
230  const int vdim = fes->GetVDim();
231  const mfem::FiniteElementCollection *fec = fes->FEColl();
232  const int vert_dofs = fec->DofForGeometry(mfem::Geometry::POINT);
233  const int edge_dofs = fec->DofForGeometry(mfem::Geometry::SEGMENT);
234  const int tri_dofs = fec->DofForGeometry(mfem::Geometry::TRIANGLE);
235  const int quad_dofs = fec->DofForGeometry(mfem::Geometry::SQUARE);
236  const int tet_dofs = fec->DofForGeometry(mfem::Geometry::TETRAHEDRON);
237  const int hex_dofs = fec->DofForGeometry(mfem::Geometry::CUBE);
238  int ent_dofs[FMS_NUM_ENTITY_TYPES];
239  ent_dofs[FMS_VERTEX] = vert_dofs;
240  ent_dofs[FMS_EDGE] = edge_dofs;
241  ent_dofs[FMS_TRIANGLE] = tri_dofs;
242  ent_dofs[FMS_QUADRILATERAL] = quad_dofs;
243  ent_dofs[FMS_TETRAHEDRON] = tet_dofs;
244  ent_dofs[FMS_HEXAHEDRON] = hex_dofs;
245  FmsInt fms_dof_offset = 0;
246  int mfem_ent_cnt[4] = {0,0,0,0}; // mfem entity counters, by dimension
247  int mfem_last_vert_cnt = 0;
250  if (dim >= 2 && edge_dofs > 0)
251  {
252  mfem::Array<int> ev;
253  for (int i = 0; i < mesh->GetNEdges(); i++)
254  {
255  mesh->GetEdgeVertices(i, ev);
256  int id = mfem_edge.GetId(ev[0], ev[1]);
257  if (id != i) { return 13; }
258  }
259  }
260  if (dim >= 3 &&
261  ((n_ents[FMS_TRIANGLE] > 0 && tri_dofs > 0) ||
262  (n_ents[FMS_QUADRILATERAL] > 0 && quad_dofs > 0)))
263  {
264  mfem::Array<int> fv;
265  for (int i = 0; i < mesh->GetNFaces(); i++)
266  {
267  mesh->GetFaceVertices(i, fv);
268  if (fv.Size() == 3) { fv.Append(INT_MAX); }
269  // HashTable uses the smallest 3 of the 4 indices to hash Hashed4
270  int id = mfem_face.GetId(fv[0], fv[1], fv[2], fv[3]);
271  if (id != i) { return 14; }
272  }
273  }
274 
275  // Loop over all parts of the main component
276  for (FmsInt part_id = 0; part_id < n_main_parts; part_id++)
277  {
278  // Loop over all entity types in the part
279  for (FmsInt et = FMS_VERTEX; et < FMS_NUM_ENTITY_TYPES; et++)
280  {
281  FmsDomain domain;
282  FmsIntType ent_id_type;
283  const void *ents;
284  const FmsOrientation *ents_ori;
285  FmsInt num_ents;
286  FmsComponentGetPart(main_comp, part_id, (FmsEntityType)et, &domain,
287  &ent_id_type, &ents, &ents_ori, &num_ents);
288  if (num_ents == 0) { continue; }
289  if (ent_dofs[et] == 0)
290  {
291  if (et == FMS_VERTEX) { mfem_last_vert_cnt = mfem_ent_cnt[et]; }
292  mfem_ent_cnt[FmsEntityDim[et]] += num_ents;
293  continue;
294  }
295 
296  if (ents != NULL &&
297  (ent_id_type != FMS_INT32 && ent_id_type != FMS_UINT32))
298  {
299  return 15;
300  }
301  if (ents_ori != NULL)
302  {
303  return 16;
304  }
305 
306  if (et == FMS_VERTEX)
307  {
308  const int mfem_dof_offset = mfem_ent_cnt[0]*vert_dofs;
309  for (FmsInt i = 0; i < num_ents*vert_dofs; i++)
310  {
311  for (int j = 0; j < vdim; j++)
312  {
313  const int idx = i*nstride+j*vstride;
314  func(mfem_dof_offset*nstride+idx) =
315  static_cast<double>(src_data[fms_dof_offset*nstride+idx]);
316  }
317  }
318  fms_dof_offset += num_ents*vert_dofs;
319  mfem_last_vert_cnt = mfem_ent_cnt[et];
320  mfem_ent_cnt[0] += num_ents;
321  continue;
322  }
323  mfem::Array<int> dofs;
324  if (FmsEntityDim[et] == dim)
325  {
326  for (FmsInt e = 0; e < num_ents; e++)
327  {
328  fes->GetElementInteriorDofs(mfem_ent_cnt[dim]+e, dofs);
329  for (int i = 0; i < ent_dofs[et]; i++, fms_dof_offset++)
330  {
331  for (int j = 0; j < vdim; j++)
332  {
333  func(fes->DofToVDof(dofs[i],j)) =
334  static_cast<double>(src_data[fms_dof_offset*nstride+j*vstride]);
335  }
336  }
337  }
338  mfem_ent_cnt[dim] += num_ents;
339  continue;
340  }
341  const FmsInt nv = FmsEntityNumVerts[et];
342  mfem::Array<int> ents_verts(num_ents*nv), m_ev;
343  const int *ei = (const int *)ents;
344  if (ents == NULL)
345  {
346  FmsDomainGetEntitiesVerts(domain, (FmsEntityType)et, NULL, FMS_INT32,
347  0, ents_verts.GetData(), num_ents);
348  }
349  else
350  {
351  for (FmsInt i = 0; i < num_ents; i++)
352  {
353  FmsDomainGetEntitiesVerts(domain, (FmsEntityType)et, NULL,
354  FMS_INT32, ei[i], &ents_verts[i*nv], 1);
355  }
356  }
357  for (int i = 0; i < ents_verts.Size(); i++)
358  {
359  ents_verts[i] += mfem_last_vert_cnt;
360  }
361  const int *perm;
362  switch ((FmsEntityType)et)
363  {
364  case FMS_EDGE:
365  {
366  for (FmsInt part_ent_id = 0; part_ent_id < num_ents; part_ent_id++)
367  {
368  const int *ev = &ents_verts[2*part_ent_id];
369  int mfem_edge_id = mfem_edge.FindId(ev[0], ev[1]);
370  if (mfem_edge_id < 0)
371  {
372  return 17;
373  }
374  mesh->GetEdgeVertices(mfem_edge_id, m_ev);
375  int ori = (ev[0] == m_ev[0]) ? 0 : 1;
377  fes->GetEdgeInteriorDofs(mfem_edge_id, dofs);
378  for (int i = 0; i < edge_dofs; i++)
379  {
380  for (int j = 0; j < vdim; j++)
381  {
382  func(fes->DofToVDof(dofs[i],j)) =
383  static_cast<double>(src_data[(fms_dof_offset+perm[i])*nstride+j*vstride]);
384  }
385  }
386  fms_dof_offset += edge_dofs;
387  }
388  break;
389  }
390  case FMS_TRIANGLE:
391  {
392  for (FmsInt part_ent_id = 0; part_ent_id < num_ents; part_ent_id++)
393  {
394  const int *tv = &ents_verts[3*part_ent_id];
395  int mfem_face_id = mfem_face.FindId(tv[0], tv[1], tv[2], INT_MAX);
396  if (mfem_face_id < 0)
397  {
398  return 18;
399  }
400  mesh->GetFaceVertices(mfem_face_id, m_ev);
401  int ori = 0;
402  while (tv[ori] != m_ev[0]) { ori++; }
403  ori = (tv[(ori+1)%3] == m_ev[1]) ? 2*ori : 2*ori+1;
405  fes->GetFaceInteriorDofs(mfem_face_id, dofs);
406  for (int i = 0; i < tri_dofs; i++)
407  {
408  for (int j = 0; j < vdim; j++)
409  {
410  func(fes->DofToVDof(dofs[i],j)) =
411  static_cast<double>(src_data[(fms_dof_offset+perm[i])*nstride+j*vstride]);
412  }
413  }
414  fms_dof_offset += tri_dofs;
415  }
416  break;
417  }
418  case FMS_QUADRILATERAL:
419  {
420  for (FmsInt part_ent_id = 0; part_ent_id < num_ents; part_ent_id++)
421  {
422  const int *qv = &ents_verts[4*part_ent_id];
423  int mfem_face_id = mfem_face.FindId(qv[0], qv[1], qv[2], qv[3]);
424  if (mfem_face_id < 0) { return 19; }
425  mesh->GetFaceVertices(mfem_face_id, m_ev);
426  int ori = 0;
427  while (qv[ori] != m_ev[0]) { ori++; }
428  ori = (qv[(ori+1)%4] == m_ev[1]) ? 2*ori : 2*ori+1;
430  fes->GetFaceInteriorDofs(mfem_face_id, dofs);
431  for (int i = 0; i < quad_dofs; i++)
432  {
433  for (int j = 0; j < vdim; j++)
434  {
435  func(fes->DofToVDof(dofs[i],j)) =
436  static_cast<double>(src_data[(fms_dof_offset+perm[i])*nstride+j*vstride]);
437  }
438  }
439  fms_dof_offset += quad_dofs;
440  }
441  break;
442  }
443  default: break;
444  }
445  mfem_ent_cnt[FmsEntityDim[et]] += num_ents;
446  }
447  }
448 
449  return err;
450 }
451 
452 int
453 FmsMeshToMesh(FmsMesh fms_mesh, Mesh **mfem_mesh)
454 {
455  FmsInt dim, n_vert, n_elem, n_bdr_elem, space_dim;
456 
457  // Find the first component that has coordinates - that will be the new mfem
458  // mesh.
459  FmsInt num_comp;
460  FmsMeshGetNumComponents(fms_mesh, &num_comp);
461  FmsComponent main_comp = NULL;
462  FmsField coords = NULL;
463  for (FmsInt comp_id = 0; comp_id < num_comp; comp_id++)
464  {
465  FmsComponent comp;
466  FmsMeshGetComponent(fms_mesh, comp_id, &comp);
467  FmsComponentGetCoordinates(comp, &coords);
468  if (coords) { main_comp = comp; break; }
469  }
470  if (!main_comp) { return 1; }
471  FmsComponentGetDimension(main_comp, &dim);
472  FmsComponentGetNumEntities(main_comp, &n_elem);
473  FmsInt n_ents[FMS_NUM_ENTITY_TYPES];
474  FmsInt n_main_parts;
475  FmsComponentGetNumParts(main_comp, &n_main_parts);
476 
477 #define RENUMBER_ENTITIES
478 #ifdef RENUMBER_ENTITIES
479  // I noticed that to get domains working right, since they appear to be
480  // defined in a local vertex numbering scheme, we have to offset the vertex
481  // ids that MFEM makes for shapes to move them to the coordinates in the
482  // current domain.
483 
484  // However, parts would just be a set of element ids in the current domain
485  // and it does not seem appropriate to offset the points in that case.
486  // Should domains be treated specially?
487  int *verts_per_part = new int[n_main_parts];
488 #endif
489 
490  // Sum the counts for each entity type across parts.
491  for (FmsInt et = FMS_VERTEX; et < FMS_NUM_ENTITY_TYPES; et++)
492  {
493  n_ents[et] = 0;
494  for (FmsInt part_id = 0; part_id < n_main_parts; part_id++)
495  {
496  FmsInt num_ents;
497  FmsComponentGetPart(main_comp, part_id, (FmsEntityType)et, NULL, NULL,
498  NULL, NULL, &num_ents);
499  n_ents[et] += num_ents;
500 #ifdef RENUMBER_ENTITIES
501  if (et == FMS_VERTEX)
502  {
503  verts_per_part[part_id] = num_ents;
504  }
505 #endif
506  }
507  }
508  n_vert = n_ents[FMS_VERTEX];
509 
510 #ifdef RENUMBER_ENTITIES
511  int *verts_start = new int[n_main_parts];
512  verts_start[0] = 0;
513  for (FmsInt i = 1; i < n_main_parts; ++i)
514  {
515  verts_start[i] = verts_start[i-1] + verts_per_part[i-1];
516  }
517 #endif
518 
519  // The first related component of dimension dim-1 will be the boundary of the
520  // new mfem mesh.
521  FmsComponent bdr_comp = NULL;
522  FmsInt num_rel_comps;
523  const FmsInt *rel_comp_ids;
524  FmsComponentGetRelations(main_comp, &rel_comp_ids, &num_rel_comps);
525  for (FmsInt i = 0; i < num_rel_comps; i++)
526  {
527  FmsComponent comp;
528  FmsMeshGetComponent(fms_mesh, rel_comp_ids[i], &comp);
529  FmsInt comp_dim;
530  FmsComponentGetDimension(comp, &comp_dim);
531  if (comp_dim == dim-1) { bdr_comp = comp; break; }
532  }
533  if (bdr_comp)
534  {
535  FmsComponentGetNumEntities(bdr_comp, &n_bdr_elem);
536  }
537  else
538  {
539  n_bdr_elem = 0;
540  }
541 
542  FmsFieldGet(coords, NULL, &space_dim, NULL, NULL, NULL);
543  int err = 0;
544  Mesh *mesh = nullptr;
545  mesh = new Mesh(dim, n_vert, n_elem, n_bdr_elem, space_dim);
546 
547  // Element tags
548  FmsInt num_tags;
549  FmsMeshGetNumTags(fms_mesh, &num_tags);
550  FmsTag elem_tag = NULL, bdr_tag = NULL;
551  for (FmsInt tag_id = 0; tag_id < num_tags; tag_id++)
552  {
553  FmsTag tag;
554  FmsMeshGetTag(fms_mesh, tag_id, &tag);
555  FmsComponent comp;
556  FmsTagGetComponent(tag, &comp);
557  if (!elem_tag && comp == main_comp)
558  {
559 #if DEBUG_FMS_MFEM
560  const char *tn = NULL;
561  FmsTagGetName(tag, &tn);
562  mfem::out << "Found element tag " << tn << std::endl;
563 #endif
564  elem_tag = tag;
565  }
566  else if (!bdr_tag && comp == bdr_comp)
567  {
568 #if DEBUG_FMS_MFEM
569  const char *tn = NULL;
570  FmsTagGetName(tag, &tn);
571  mfem::out << "Found boundary tag " << tn << std::endl;
572 #endif
573  bdr_tag = tag;
574  }
575  }
576  FmsIntType attr_type;
577  const void *v_attr, *v_bdr_attr;
578  mfem::Array<int> attr, bdr_attr;
579  FmsInt num_attr;
580  // Element attributes
581  if (elem_tag)
582  {
583  FmsTagGet(elem_tag, &attr_type, &v_attr, &num_attr);
584  if (attr_type == FMS_UINT8)
585  {
586  mfem::Array<uint8_t> at((uint8_t*)v_attr, num_attr);
587  attr = at;
588  }
589  else if (attr_type == FMS_INT32 || attr_type == FMS_UINT32)
590  {
591  attr.MakeRef((int*)v_attr, num_attr);
592  }
593  else
594  {
595  err = 1; // "attribute type not supported!"
596  goto func_exit;
597  }
598  }
599  // Boundary attributes
600  if (bdr_tag)
601  {
602  FmsTagGet(bdr_tag, &attr_type, &v_bdr_attr, &num_attr);
603  if (attr_type == FMS_UINT8)
604  {
605  mfem::Array<uint8_t> at((uint8_t*)v_bdr_attr, num_attr);
606  bdr_attr = at;
607  }
608  else if (attr_type == FMS_INT32 || attr_type == FMS_UINT32)
609  {
610  bdr_attr.MakeRef((int*)v_bdr_attr, num_attr);
611  }
612  else
613  {
614  err = 2; // "bdr attribute type not supported!"
615  goto func_exit;
616  }
617  }
618 
619  // Add elements
620  for (FmsInt part_id = 0; part_id < n_main_parts; part_id++)
621  {
622  for (int et = FMS_VERTEX; et < FMS_NUM_ENTITY_TYPES; et++)
623  {
624  if (FmsEntityDim[et] != dim) { continue; }
625 
626  FmsDomain domain;
627  FmsIntType elem_id_type;
628  const void *elem_ids;
629  const FmsOrientation *elem_ori;
630  FmsInt num_elems;
631  FmsComponentGetPart(main_comp, part_id, (FmsEntityType)et, &domain,
632  &elem_id_type, &elem_ids, &elem_ori, &num_elems);
633 
634  if (num_elems == 0) { continue; }
635 
636  if (elem_ids != NULL &&
637  (elem_id_type != FMS_INT32 && elem_id_type != FMS_UINT32))
638  {
639  err = 3; goto func_exit;
640  }
641  if (elem_ori != NULL)
642  {
643  err = 4; goto func_exit;
644  }
645 
646  const FmsInt nv = FmsEntityNumVerts[et];
647  mfem::Array<int> ents_verts(num_elems*nv);
648  if (elem_ids == NULL)
649  {
650  FmsDomainGetEntitiesVerts(domain, (FmsEntityType)et, NULL, FMS_INT32,
651  0, ents_verts.GetData(), num_elems);
652  }
653  else
654  {
655  const int *ei = (const int *)elem_ids;
656  for (FmsInt i = 0; i < num_elems; i++)
657  {
658  FmsDomainGetEntitiesVerts(domain, (FmsEntityType)et, NULL, FMS_INT32,
659  ei[i], &ents_verts[i*nv], 1);
660  }
661  }
662  const int elem_offset = mesh->GetNE();
663  switch ((FmsEntityType)et)
664  {
665  case FMS_EDGE:
666  err = 5;
667  goto func_exit;
668  break;
669  case FMS_TRIANGLE:
670 #ifdef RENUMBER_ENTITIES
671  // The domain vertices/edges were defined in local ordering. We
672  // now have a set of triangle vertices defined in terms of local
673  // vertex numbers. Renumber them to a global numbering.
674  for (FmsInt i = 0; i < num_elems*3; i++)
675  {
676  ents_verts[i] += verts_start[part_id];
677  }
678 #endif
679 
680  for (FmsInt i = 0; i < num_elems; i++)
681  {
682  mesh->AddTriangle(
683  &ents_verts[3*i], elem_tag ? attr[elem_offset+i] : 1);
684  }
685  break;
686  case FMS_QUADRILATERAL:
687 #ifdef RENUMBER_ENTITIES
688  for (FmsInt i = 0; i < num_elems*4; i++)
689  {
690  ents_verts[i] += verts_start[part_id];
691  }
692 #endif
693  for (FmsInt i = 0; i < num_elems; i++)
694  {
695  mesh->AddQuad(&ents_verts[4*i], elem_tag ? attr[elem_offset+i] : 1);
696  }
697  break;
698  case FMS_TETRAHEDRON:
699 #ifdef RENUMBER_ENTITIES
700  for (FmsInt i = 0; i < num_elems*4; i++)
701  {
702  ents_verts[i] += verts_start[part_id];
703  }
704 #endif
705  for (FmsInt i = 0; i < num_elems; i++)
706  {
707  mesh->AddTet(&ents_verts[4*i], elem_tag ? attr[elem_offset+i] : 1);
708  }
709  break;
710 
711  // TODO: What about wedges and pyramids?
712 
713 
714  case FMS_HEXAHEDRON:
715 #ifdef RENUMBER_ENTITIES
716  for (FmsInt i = 0; i < num_elems*8; i++)
717  {
718  ents_verts[i] += verts_start[part_id];
719  }
720 
721 #endif
722  for (FmsInt i = 0; i < num_elems; i++)
723  {
724  const int *hex_verts = &ents_verts[8*i];
725 #if 0
726  const int reorder[8] = {0, 1, 2, 3, 5, 4, 6, 7};
727  const int new_verts[8] = {hex_verts[reorder[0]],
728  hex_verts[reorder[1]],
729  hex_verts[reorder[2]],
730  hex_verts[reorder[3]],
731  hex_verts[reorder[4]],
732  hex_verts[reorder[5]],
733  hex_verts[reorder[6]],
734  hex_verts[reorder[7]]
735  };
736  hex_verts = new_verts;
737 #endif
738  mesh->AddHex(hex_verts, elem_tag ? attr[elem_offset+i] : 1);
739  }
740  break;
741  default:
742  break;
743  }
744  }
745  }
746 
747  // Add boundary elements
748  if (bdr_comp && n_bdr_elem > 0)
749  {
750  FmsInt n_bdr_parts;
751  FmsComponentGetNumParts(bdr_comp, &n_bdr_parts);
752 
753  for (FmsInt part_id = 0; part_id < n_bdr_parts; part_id++)
754  {
755  for (int et = FMS_VERTEX; et < FMS_NUM_ENTITY_TYPES; et++)
756  {
757  if (FmsEntityDim[et] != dim-1) { continue; }
758 
759  FmsDomain domain;
760  FmsIntType elem_id_type;
761  const void *elem_ids;
762  const FmsOrientation *elem_ori;
763  FmsInt num_elems;
764  FmsComponentGetPart(bdr_comp, part_id, (FmsEntityType)et, &domain,
765  &elem_id_type, &elem_ids, &elem_ori, &num_elems);
766  if (num_elems == 0) { continue; }
767 
768  if (elem_ids != NULL &&
769  (elem_id_type != FMS_INT32 && elem_id_type != FMS_UINT32))
770  {
771  err = 6; goto func_exit;
772  }
773  if (elem_ori != NULL)
774  {
775  err = 7; goto func_exit;
776  }
777 
778  const FmsInt nv = FmsEntityNumVerts[et];
779  mfem::Array<int> ents_verts(num_elems*nv);
780  if (elem_ids == NULL)
781  {
782  FmsDomainGetEntitiesVerts(domain, (FmsEntityType)et, NULL, FMS_INT32,
783  0, ents_verts.GetData(), num_elems);
784  }
785  else
786  {
787  const int *ei = (const int *)elem_ids;
788  for (FmsInt i = 0; i < num_elems; i++)
789  {
790  FmsDomainGetEntitiesVerts(domain, (FmsEntityType)et, NULL,
791  FMS_INT32, ei[i], &ents_verts[i*nv], 1);
792  }
793  }
794  const int elem_offset = mesh->GetNBE();
795  switch ((FmsEntityType)et)
796  {
797  case FMS_EDGE:
798  for (FmsInt i = 0; i < num_elems; i++)
799  {
800  mesh->AddBdrSegment(
801  &ents_verts[2*i], bdr_tag ? bdr_attr[elem_offset+i] : 1);
802  }
803  break;
804  case FMS_TRIANGLE:
805  for (FmsInt i = 0; i < num_elems; i++)
806  {
807  mesh->AddBdrTriangle(
808  &ents_verts[3*i], bdr_tag ? bdr_attr[elem_offset+i] : 1);
809  }
810  break;
811  case FMS_QUADRILATERAL:
812  for (FmsInt i = 0; i < num_elems; i++)
813  {
814  mesh->AddBdrQuad(
815  &ents_verts[4*i], bdr_tag ? bdr_attr[elem_offset+i] : 1);
816  }
817  break;
818  default:
819  break;
820  }
821  }
822  }
823  }
824 
825 #ifdef RENUMBER_ENTITIES
826  delete [] verts_per_part;
827  delete [] verts_start;
828 #endif
829 
830  // Transfer coordinates
831  {
832  // Set the vertex coordinates to zero
833  const double origin[3] = {0.,0.,0.};
834  for (FmsInt vi = 0; vi < n_vert; vi++)
835  {
836  mesh->AddVertex(origin);
837  }
838 
839  // Finalize the mesh topology
840  mesh->FinalizeTopology();
841 
842  FmsFieldDescriptor coords_fd = NULL;
843  FmsLayoutType coords_layout;
844  FmsFieldGet(coords, &coords_fd, NULL, &coords_layout, NULL, NULL);
845  if (!coords_fd)
846  {
847  mfem::err << "Error reading the FMS mesh coords' FieldDescriptor." << std::endl;
848  err = 8;
849  goto func_exit;
850  }
851  FmsInt coords_order = 0;
852  FmsBasisType coords_btype = FMS_NODAL_GAUSS_CLOSED;
853  FmsFieldType coords_ftype = FMS_CONTINUOUS;
854  FmsFieldDescriptorGetFixedOrder(coords_fd, &coords_ftype, &coords_btype,
855  &coords_order);
856  // Maybe this is extra but it seems mesh->SetCurvature assumes
857  // btype=1. Maybe protects us against corrupt data.
858  if (coords_btype != FMS_NODAL_GAUSS_CLOSED)
859  {
860  mfem::err << "Error reading FMS mesh coords." << std::endl;
861  err = 9;
862  goto func_exit;
863  }
864 
865  // Switch to mfem::Mesh with nodes (interpolates the linear coordinates)
866  const bool discont = (coords_ftype == FMS_DISCONTINUOUS);
867  mesh->SetCurvature(coords_order, discont, space_dim,
868  (coords_layout == FMS_BY_VDIM) ?
870 
871  // Finalize mesh construction
872  mesh->Finalize();
873 
874  // Set the high-order mesh nodes
875  mfem::GridFunction &nodes = *mesh->GetNodes();
876  FmsFieldToGridFunction<double>(fms_mesh, coords, mesh, nodes, false);
877  }
878 
879 func_exit:
880 
881  if (err)
882  {
883  delete mesh;
884  }
885  else
886  {
887  *mfem_mesh = mesh;
888  }
889  return err;
890 }
891 
892 bool
893 BasisTypeToFmsBasisType(int bt, FmsBasisType &btype)
894 {
895  bool retval = false;
896  switch (bt)
897  {
899  // mfem::out << "mfem::BasisType::GaussLegendre -> FMS_NODAL_GAUSS_OPEN" << endl;
900  btype = FMS_NODAL_GAUSS_OPEN;
901  retval = true;
902  break;
904  // mfem::out << "mfem::BasisType::GaussLobato -> FMS_NODAL_GAUSS_CLOSED" << endl;
905  btype = FMS_NODAL_GAUSS_CLOSED;
906  retval = true;
907  break;
909  // mfem::out << "mfem::BasisType::Positive -> FMS_POSITIVE" << endl;
910  btype = FMS_POSITIVE;
911  retval = true;
912  break;
914  // mfem::out << "mfem::BasisType::OpenUniform -> ?" << endl;
915  btype = FMS_NODAL_UNIFORM_OPEN;
916  retval = true;
917  break;
919  // mfem::out << "mfem::BasisType::ClosedUniform -> ?" << endl;
920  btype = FMS_NODAL_UNIFORM_CLOSED;
921  retval = true;
922  break;
924  // mfem::out << "mfem::BasisType::OpenHalfUniform -> ?" << endl;
925  break;
927  // mfem::out << "mfem::BasisType::Serendipity -> ?" << endl;
928  break;
930  // mfem::out << "mfem::BasisType::ClosedGL -> ?" << endl;
931  break;
932 
933  }
934  /*
935  Which MFEM types map to:?
936  FMS_NODAL_CHEBYSHEV_OPEN,
937  FMS_NODAL_CHEBYSHEV_CLOSED,
938  */
939 
940  return retval;
941 }
942 
943 /**
944 @note We add the FMS field descriptor and field in here so we can only do it
945  after successfully validating the inputs (handling certain grid function
946  types, etc.)
947 */
948 int
949 GridFunctionToFmsField(FmsDataCollection dc,
950  FmsComponent comp,
951  const std::string &fd_name,
952  const std::string &field_name,
953  const Mesh *mesh,
954  const GridFunction *gf,
955  FmsField *outfield)
956 {
957  if (!dc) { return 1; }
958  if (!comp) { return 2; }
959  if (!mesh) { return 3; }
960  if (!gf) { return 4; }
961  if (!outfield) { return 5; }
962 
963  double *c = gf->GetData();
964 
965  const mfem::FiniteElementSpace *fespace = gf->FESpace();
966  const mfem::FiniteElementCollection *fecoll = fespace->FEColl();
967 
968 #ifdef DEBUG_MFEM_FMS
969  mfem::out << "Adding FMS field for " << field_name << "..." << endl;
970 #endif
971 
972  /* Q: No getter for the basis, do different kinds of FECollection have
973  implied basis? There are two subclasses that actually have the getter,
974  maybe those aren't implied? */
975  FmsInt order = 1;
976  int vdim = 1;
977  FmsFieldType ftype = FMS_CONTINUOUS;
978  FmsBasisType btype = FMS_NODAL_GAUSS_CLOSED;
979  switch (fecoll->GetContType())
980  {
982  {
983  ftype = FMS_CONTINUOUS;
984  order = static_cast<FmsInt>(fespace->GetOrder(0));
985  vdim = gf->VectorDim();
986  auto fec = dynamic_cast<const mfem::H1_FECollection *>(fecoll);
987  if (fec != nullptr)
988  {
989  if (!BasisTypeToFmsBasisType(fec->GetBasisType(), btype))
990  {
991  mfem::err << "Error converting MFEM basis type to FMS for"
992  " FMS_CONTINUOUS." << std::endl;
993  return 6;
994  }
995  }
996  break;
997  }
999  {
1000  ftype = FMS_DISCONTINUOUS;
1001  order = static_cast<FmsInt>(fespace->GetOrder(0));
1002  vdim = gf->VectorDim();
1003  auto fec = dynamic_cast<const mfem::L2_FECollection *>(fecoll);
1004  if (fec != nullptr)
1005  {
1006  if (!BasisTypeToFmsBasisType(fec->GetBasisType(), btype))
1007  {
1008  mfem::err << "Error converting MFEM basis type to FMS for"
1009  " FMS_DISCONTINUOUS." << std::endl;
1010  return 7;
1011  }
1012  }
1013  break;
1014  }
1016  {
1017  mfem::out << "Warning, unsupported ContType (TANGENTIAL) for "
1018  << field_name << ". Using FMS_CONTINUOUS." << std::endl;
1019  break;
1020  }
1022  {
1023  ftype = FMS_HDIV;
1024  // This RT_FECollection type seems to arise from "RT" fields such as "RT_3D_P1".
1025  // Checking fe_coll.hpp, this contains verbiage about H_DIV so we assign it the
1026  // FMS type FMS_HDIV.
1027 
1028  // I've seen RT_3D_P1 return the wrong order so get it from the name.
1029  int idim, iorder;
1030  if (sscanf(fecoll->Name(), "RT_%dD_P%d", &idim, &iorder) == 2)
1031  {
1032  order = (FmsInt)iorder;
1033  }
1034  else
1035  {
1036  order = static_cast<FmsInt>(fespace->GetOrder(0));
1037  }
1038 
1039  // Get the vdim from the fespace since the grid function is returning
1040  // 3 but we need it to be what was read from the file so we can pass the
1041  // right vdim to the FMS field descriptor to compute the expected number of dofs.
1042  vdim = fespace->GetVDim();
1043  break;
1044  }
1045  default:
1046  mfem::out << "Warning, unsupported ContType for field " << field_name
1047  << ". Using FMS_CONTINUOUS." << std::endl;
1048  ftype = FMS_CONTINUOUS;
1049  break;
1050  }
1051 
1052  // Now that we're not failing, create the fd and field.
1053  FmsFieldDescriptor fd = NULL;
1054  FmsField f = NULL;
1055  FmsDataCollectionAddFieldDescriptor(dc, fd_name.c_str(), &fd);
1056  FmsDataCollectionAddField(dc, field_name.c_str(), &f);
1057  *outfield = f;
1058 
1059  /* Q: Why is order defined on a per element basis? */
1060  FmsFieldDescriptorSetComponent(fd, comp);
1061  FmsFieldDescriptorSetFixedOrder(fd, ftype, btype, order);
1062 
1063  FmsInt ndofs;
1064  FmsFieldDescriptorGetNumDofs(fd, &ndofs);
1065 
1066  const char *name = NULL;
1067  FmsFieldGetName(f, &name);
1068  FmsLayoutType layout = fespace->GetOrdering() == mfem::Ordering::byVDIM ?
1069  FMS_BY_VDIM : FMS_BY_NODES;
1070 
1071 #ifdef DEBUG_MFEM_FMS
1072  switch (ftype)
1073  {
1074  case FMS_CONTINUOUS:
1075  mfem::out << "\tFMS_CONTINUOUS" << std::endl;
1076  break;
1077  case FMS_DISCONTINUOUS:
1078  mfem::out << "\tFMS_DISCONTINUOUS" << std::endl;
1079  break;
1080  case FMS_HDIV:
1081  mfem::out << "\tFMS_HDIV" << std::endl;
1082  break;
1083  }
1084  mfem::out << "\tField is order " << order << " with vdim " << vdim <<
1085  " and nDoFs " << ndofs << std::endl;
1086  mfem::out << "\tgf->size() " << gf->Size() << " ndofs * vdim " << ndofs * vdim
1087  << std::endl;
1088  mfem::out << "\tlayout " << layout << " (0 = BY_NODES, 1 = BY_VDIM)" <<
1089  std::endl;
1090 #endif
1091 
1092  if (FmsFieldSet(f, fd, vdim, layout, FMS_DOUBLE, c))
1093  {
1094  mfem::err << "Error setting field " << field_name << " in FMS." << std::endl;
1095  return 8;
1096  }
1097  return 0;
1098 }
1099 
1100 bool
1101 MfemMetaDataToFmsMetaData(DataCollection *mdc, FmsDataCollection fdc)
1102 {
1103  if (!mdc) { return false; }
1104  if (!fdc) { return false; }
1105 
1106  int *cycle = NULL;
1107  double *time = NULL, *timestep = NULL;
1108  FmsMetaData top_level = NULL;
1109  FmsMetaData *cycle_time_timestep = NULL;
1110  int mdata_err = 0;
1111  mdata_err = FmsDataCollectionAttachMetaData(fdc, &top_level);
1112  if (!top_level || mdata_err)
1113  {
1114  mfem::err << "Failed to attach metadata to the FmsDataCollection" << std::endl;
1115  return 40;
1116  }
1117 
1118  mdata_err = FmsMetaDataSetMetaData(top_level, "MetaData", 3,
1119  &cycle_time_timestep);
1120  if (!cycle_time_timestep || mdata_err)
1121  {
1122  mfem::err << "Failed to acquire FmsMetaData array" << std::endl;
1123  return false;
1124  }
1125 
1126  if (!cycle_time_timestep[0])
1127  {
1128  mfem::err << "The MetaData pointer for cycle is NULL" << std::endl;
1129  return false;
1130  }
1131  mdata_err = FmsMetaDataSetIntegers(cycle_time_timestep[0], "cycle",
1132  FMS_INT32, 1, (void**)&cycle);
1133  if (!cycle || mdata_err)
1134  {
1135  mfem::err << "The data pointer for cycle is NULL" << std::endl;
1136  return false;
1137  }
1138  *cycle = mdc->GetCycle();
1139 
1140  if (!cycle_time_timestep[1])
1141  {
1142  mfem::err << "The FmsMetaData pointer for time is NULL" << std::endl;
1143  return false;
1144  }
1145  mdata_err = FmsMetaDataSetScalars(cycle_time_timestep[1], "time", FMS_DOUBLE,
1146  1, (void**)&time);
1147  if (!time || mdata_err)
1148  {
1149  mfem::err << "The data pointer for time is NULL." << std::endl;
1150  return false;
1151  }
1152  *time = mdc->GetTime();
1153 
1154  if (!cycle_time_timestep[2])
1155  {
1156  mfem::err << "The FmsMetData pointer for timestep is NULL" << std::endl;
1157  return false;
1158  }
1159  mdata_err = FmsMetaDataSetScalars(cycle_time_timestep[2], "timestep",
1160  FMS_DOUBLE, 1, (void**)&timestep);
1161  if (!timestep || mdata_err)
1162  {
1163  mfem::err << "The data pointer for timestep is NULL" << std::endl;
1164  return false;
1165  }
1166  *timestep = mdc->GetTimeStep();
1167 
1168  return true;
1169 }
1170 
1171 //---------------------------------------------------------------------------
1172 bool
1173 FmsMetaDataGetInteger(FmsMetaData mdata, const std::string &key,
1174  std::vector<int> &values)
1175 {
1176  if (!mdata) { return false; }
1177 
1178  bool retval = false;
1179  FmsMetaDataType type;
1180  FmsIntType int_type;
1181  FmsInt i, size;
1182  FmsMetaData *children = nullptr;
1183  const void *data = nullptr;
1184  const char *mdata_name = nullptr;
1185  if (FmsMetaDataGetType(mdata, &type) == 0)
1186  {
1187  switch (type)
1188  {
1189  case FMS_INTEGER:
1190  if (FmsMetaDataGetIntegers(mdata, &mdata_name, &int_type, &size, &data) == 0)
1191  {
1192  if (strcasecmp(key.c_str(), mdata_name) == 0)
1193  {
1194  retval = true;
1195 
1196  // Interpret the integers and store them in the std::vector<int>
1197  switch (int_type)
1198  {
1199  case FMS_INT8:
1200  for (i = 0; i < size; i++)
1201  {
1202  values.push_back(static_cast<int>(reinterpret_cast<const int8_t*>(data)[i]));
1203  }
1204  break;
1205  case FMS_INT16:
1206  for (i = 0; i < size; i++)
1207  {
1208  values.push_back(static_cast<int>(reinterpret_cast<const int16_t*>(data)[i]));
1209  }
1210  break;
1211  case FMS_INT32:
1212  for (i = 0; i < size; i++)
1213  {
1214  values.push_back(static_cast<int>(reinterpret_cast<const int32_t*>(data)[i]));
1215  }
1216  break;
1217  case FMS_INT64:
1218  for (i = 0; i < size; i++)
1219  {
1220  values.push_back(static_cast<int>(reinterpret_cast<const int64_t*>(data)[i]));
1221  }
1222  break;
1223  case FMS_UINT8:
1224  for (i = 0; i < size; i++)
1225  {
1226  values.push_back(static_cast<int>(reinterpret_cast<const uint8_t*>(data)[i]));
1227  }
1228  break;
1229  case FMS_UINT16:
1230  for (i = 0; i < size; i++)
1231  {
1232  values.push_back(static_cast<int>(reinterpret_cast<const uint16_t*>(data)[i]));
1233  }
1234  break;
1235  case FMS_UINT32:
1236  for (i = 0; i < size; i++)
1237  {
1238  values.push_back(static_cast<int>(reinterpret_cast<const uint32_t*>(data)[i]));
1239  }
1240  break;
1241  case FMS_UINT64:
1242  for (i = 0; i < size; i++)
1243  {
1244  values.push_back(static_cast<int>(reinterpret_cast<const uint64_t*>(data)[i]));
1245  }
1246  break;
1247  default:
1248  retval = false;
1249  break;
1250  }
1251  }
1252  }
1253  break;
1254  case FMS_META_DATA:
1255  if (FmsMetaDataGetMetaData(mdata, &mdata_name, &size, &children) == 0)
1256  {
1257  // Recurse to look for the key we want.
1258  for (i = 0; i < size && !retval; i++)
1259  {
1260  retval = FmsMetaDataGetInteger(children[i], key, values);
1261  }
1262  }
1263  break;
1264  default:
1265  break;
1266  }
1267  }
1268 
1269  return retval;
1270 }
1271 
1272 //---------------------------------------------------------------------------
1273 bool
1274 FmsMetaDataGetScalar(FmsMetaData mdata, const std::string &key,
1275  std::vector<double> &values)
1276 {
1277  if (!mdata) { return false; }
1278 
1279  bool retval = false;
1280  FmsMetaDataType type;
1281  FmsScalarType scal_type;
1282  FmsInt i, size;
1283  FmsMetaData *children = nullptr;
1284  const void *data = nullptr;
1285  const char *mdata_name = nullptr;
1286  if (FmsMetaDataGetType(mdata, &type) == 0)
1287  {
1288  switch (type)
1289  {
1290  case FMS_SCALAR:
1291  if (FmsMetaDataGetScalars(mdata, &mdata_name, &scal_type, &size, &data) == 0)
1292  {
1293  if (strcasecmp(key.c_str(), mdata_name) == 0)
1294  {
1295  retval = true;
1296 
1297  // Interpret the integers and store them in the std::vector<int>
1298  switch (scal_type)
1299  {
1300  case FMS_FLOAT:
1301  for (i = 0; i < size; i++)
1302  {
1303  values.push_back(static_cast<double>(reinterpret_cast<const float*>(data)[i]));
1304  }
1305  break;
1306  case FMS_DOUBLE:
1307  for (i = 0; i < size; i++)
1308  {
1309  values.push_back(reinterpret_cast<const double*>(data)[i]);
1310  }
1311  break;
1312  default:
1313  retval = false;
1314  break;
1315  }
1316  }
1317  }
1318  break;
1319  case FMS_META_DATA:
1320  if (FmsMetaDataGetMetaData(mdata, &mdata_name, &size, &children) == 0)
1321  {
1322  // Recurse to look for the key we want.
1323  for (i = 0; i < size && !retval; i++)
1324  {
1325  retval = FmsMetaDataGetScalar(children[i], key, values);
1326  }
1327  }
1328  break;
1329  default:
1330  break;
1331  }
1332  }
1333 
1334  return retval;
1335 }
1336 
1337 //---------------------------------------------------------------------------
1338 bool
1339 FmsMetaDataGetString(FmsMetaData mdata, const std::string &key,
1340  std::string &value)
1341 {
1342  if (!mdata) { return false; }
1343 
1344  bool retval = false;
1345  FmsMetaDataType type;
1346  FmsInt i, size;
1347  FmsMetaData *children = nullptr;
1348  const char *mdata_name = nullptr;
1349  const char *str_value = nullptr;
1350 
1351  if (FmsMetaDataGetType(mdata, &type) == 0)
1352  {
1353  switch (type)
1354  {
1355  case FMS_STRING:
1356  if (FmsMetaDataGetString(mdata, &mdata_name, &str_value) == 0)
1357  {
1358  if (strcasecmp(key.c_str(), mdata_name) == 0)
1359  {
1360  retval = true;
1361  value = str_value;
1362  }
1363  }
1364  break;
1365  case FMS_META_DATA:
1366  if (FmsMetaDataGetMetaData(mdata, &mdata_name, &size, &children) == 0)
1367  {
1368  // Recurse to look for the key we want.
1369  for (i = 0; i < size && !retval; i++)
1370  {
1371  retval = FmsMetaDataGetString(children[i], key, value);
1372  }
1373  }
1374  break;
1375  default:
1376  break;
1377  }
1378  }
1379 
1380  return retval;
1381 }
1382 
1383 /* -------------------------------------------------------------------------- */
1384 /* FMS to MFEM conversion function */
1385 /* -------------------------------------------------------------------------- */
1386 
1387 int FmsDataCollectionToDataCollection(FmsDataCollection dc,
1388  DataCollection **mfem_dc)
1389 {
1390  int retval = 0;
1391  FmsMesh fms_mesh;
1392  FmsDataCollectionGetMesh(dc, &fms_mesh);
1393 
1394  // NOTE: The MFEM data collection has a single Mesh. Mesh has a constructor
1395  // to take multiple Mesh objects but it appears to glue them together.
1396  Mesh *mesh = nullptr;
1397  int err = FmsMeshToMesh(fms_mesh, &mesh);
1398  if (err == 0)
1399  {
1400  std::string collection_name("collection");
1401  const char *cn = nullptr;
1402  FmsDataCollectionGetName(dc, &cn);
1403  if (cn != nullptr)
1404  {
1405  collection_name = cn;
1406  }
1407 
1408  // Make a data collection that contains the mesh.
1409  DataCollection *mdc = new DataCollection(collection_name, mesh);
1410  mdc->SetOwnData(true);
1411 
1412  // Now do fields, etc. and add them to mdc.
1413  FmsField *fields = nullptr;
1414  FmsInt num_fields = 0;
1415  if (FmsDataCollectionGetFields(dc, &fields, &num_fields) == 0)
1416  {
1417  for (FmsInt i = 0; i < num_fields; ++i)
1418  {
1419  const char *name = nullptr;
1420  FmsFieldGetName(fields[i], &name);
1421 
1422  GridFunction *gf = new GridFunction;
1423 
1424  // Get the data type.
1425  FmsFieldDescriptor f_fd;
1426  FmsLayoutType f_layout;
1427  FmsScalarType f_data_type;
1428  const void *f_data;
1429  FmsFieldGet(fields[i], &f_fd, NULL, &f_layout, &f_data_type,
1430  &f_data);
1431 
1432  // Interpret the field according to its data type.
1433  int err = 1;
1434  switch (f_data_type)
1435  {
1436  case FMS_FLOAT:
1437  err = FmsFieldToGridFunction<float>(fms_mesh, fields[i], mesh, *gf, true);
1438  break;
1439  case FMS_DOUBLE:
1440  err = FmsFieldToGridFunction<double>(fms_mesh, fields[i], mesh, *gf, true);
1441  break;
1442  case FMS_COMPLEX_FLOAT:
1443  case FMS_COMPLEX_DOUBLE:
1444  // Does MFEM support complex?
1445  break;
1446  default:
1447  break;
1448  }
1449 
1450  if (err == 0)
1451  {
1452  mdc->RegisterField(name, gf);
1453  }
1454  else
1455  {
1456  mfem::out << "There was an error converting " << name << " code: " << err <<
1457  std::endl;
1458  delete gf;
1459  }
1460 
1461  const char *fname = NULL;
1462  FmsFieldGetName(fields[i], &fname);
1463  }
1464  }
1465 
1466  // If we have metadata in FMS, pass what we can through to MFEM.
1467  FmsMetaData mdata = NULL;
1468  FmsDataCollectionGetMetaData(dc, &mdata);
1469  if (mdata)
1470  {
1471  std::vector<int> ivalues;
1472  std::vector<double> dvalues;
1473  std::string svalue;
1474  if (FmsMetaDataGetInteger(mdata, "cycle", ivalues))
1475  {
1476  if (!ivalues.empty())
1477  {
1478  mdc->SetCycle(ivalues[0]);
1479  }
1480  }
1481  if (FmsMetaDataGetScalar(mdata, "time", dvalues))
1482  {
1483  if (!dvalues.empty())
1484  {
1485  mdc->SetTime(dvalues[0]);
1486  }
1487  }
1488  if (FmsMetaDataGetScalar(mdata, "timestep", dvalues))
1489  {
1490  if (!dvalues.empty())
1491  {
1492  mdc->SetTimeStep(dvalues[0]);
1493  }
1494  }
1495  }
1496 
1497  *mfem_dc = mdc;
1498  }
1499  else
1500  {
1501  mfem::out << "FmsDataCollectionToDataCollection: mesh failed to convert. err="
1502  << err
1503  << endl;
1504 
1505  retval = 1;
1506  }
1507 
1508 #if DEBUG_FMS_MFEM
1509  if (*mfem_dc)
1510  {
1511  VisItDataCollection visit_dc(std::string("DEBUG_DC"), mesh);
1512  visit_dc.SetOwnData(false);
1513  const auto &fields = (*mfem_dc)->GetFieldMap();
1514  for (const auto &field : fields)
1515  {
1516  visit_dc.RegisterField(field.first, field.second);
1517  }
1518  visit_dc.Save();
1519  }
1520 #endif
1521 
1522  return retval;
1523 }
1524 
1525 
1526 
1527 /* -------------------------------------------------------------------------- */
1528 /* MFEM to FMS conversion function */
1529 /* -------------------------------------------------------------------------- */
1530 
1531 int
1532 MeshToFmsMesh(const Mesh *mmesh, FmsMesh *fmesh, FmsComponent *volume)
1533 {
1534  if (!mmesh) { return 1; }
1535  if (!fmesh) { return 2; }
1536  if (!volume) { return 3; }
1537 
1538 
1539  int err = 0;
1540  const int num_verticies = mmesh->GetNV();
1541 #ifdef DEBUG_MFEM_FMS
1542  const int num_edges = mmesh->GetNEdges();
1543 #endif
1544  const int num_faces = mmesh->GetNFaces();
1545  const int num_elements = mmesh->GetNE();
1546 
1547 #ifdef DEBUG_MFEM_FMS
1548  mfem::out << "nverts: " << num_verticies << std::endl;
1549  mfem::out << "nedges: " << num_edges << std::endl;
1550  mfem::out << "nfaces: " << num_faces << std::endl;
1551  mfem::out << "nele: " << num_elements << std::endl;
1552 #endif
1553 
1554  FmsMeshConstruct(fmesh);
1555  FmsMeshSetPartitionId(*fmesh, 0, 1);
1556 
1557  FmsDomain *domains = NULL;
1558  FmsMeshAddDomains(*fmesh, "Domain", 1, &domains);
1559  FmsDomainSetNumVertices(domains[0], num_verticies);
1560 
1561  const int edge_reorder[2] = {1, 0};
1562  const int quad_reorder[4] = {0,1,2,3};
1563  const int tet_reorder[4] = {3,2,1,0};
1564  const int hex_reorder[6] = {0,5,1,3,4,2};
1565  const int *reorder[8] = {NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL};
1566  reorder[FMS_EDGE] = edge_reorder;
1567  reorder[FMS_QUADRILATERAL] = quad_reorder;
1568  reorder[FMS_TETRAHEDRON] = tet_reorder;
1569  reorder[FMS_HEXAHEDRON] = hex_reorder;
1570 
1571  const mfem::Table *edges = mmesh->GetEdgeVertexTable();
1572  if (!edges)
1573  {
1574  mfem::err << "Error, mesh has no edges." << std::endl;
1575  return 1;
1576  }
1577  mfem::Table *faces = mmesh->GetFaceEdgeTable();
1578  if (!faces && num_faces > 0)
1579  {
1580  mfem::err <<
1581  "Error, mesh contains faces but the \"GetFaceEdgeTable\" returned NULL." <<
1582  std::endl;
1583  return 2;
1584  }
1585 
1586  // Build edges
1587  std::vector<int> edge_verts(edges->Size() * 2);
1588  for (int i = 0; i < edges->Size(); i++)
1589  {
1590  mfem::Array<int> nids;
1591  edges->GetRow(i, nids);
1592  for (int j = 0; j < 2; j++)
1593  {
1594  edge_verts[i*2 + j] = nids[j];
1595  }
1596  }
1597 
1598  // TODO: Move this code to after the for loop so edges can be added at top level entities
1599  FmsDomainSetNumEntities(domains[0], FMS_EDGE, FMS_INT32, edge_verts.size() / 2);
1600  FmsDomainAddEntities(domains[0], FMS_EDGE, reorder, FMS_INT32,
1601  edge_verts.data(), edge_verts.size() / 2);
1602 #ifdef DEBUG_MFEM_FMS
1603  mfem::out << "EDGES: ";
1604  for (int i = 0; i < edge_verts.size(); i++)
1605  {
1606  if (i % 2 == 0) { mfem::out << std::endl << "\t" << i/2 << ": "; }
1607  mfem::out << edge_verts[i] << " ";
1608  }
1609  mfem::out << std::endl;
1610 #endif
1611 
1612  // Build faces
1613  if (faces)
1614  {
1615  // TODO: Support Triangles and Quads, and move this code after the for
1616  // loop so these can be added as top level entities
1617  int rowsize = faces->RowSize(0);
1618  std::vector<int> face_edges(faces->Size() * rowsize);
1619  for (int i = 0; i < faces->Size(); i++)
1620  {
1621  mfem::Array<int> eids;
1622  faces->GetRow(i, eids);
1623  for (int j = 0; j < rowsize; j++)
1624  {
1625  face_edges[i*rowsize + j] = eids[j];
1626  }
1627  }
1628  FmsEntityType ent_type = (rowsize == 3) ? FMS_TRIANGLE : FMS_QUADRILATERAL;
1629  FmsDomainSetNumEntities(domains[0], ent_type, FMS_INT32,
1630  face_edges.size() / rowsize);
1631  FmsDomainAddEntities(domains[0], ent_type, NULL, FMS_INT32, face_edges.data(),
1632  face_edges.size() / rowsize);
1633 #ifdef DEBUG_MFEM_FMS
1634  mfem::out << "FACES: ";
1635  for (int i = 0; i < face_edges.size(); i++)
1636  {
1637  if (i % rowsize == 0) { mfem::out << std::endl << "\t" << i/rowsize << ": "; }
1638  mfem::out << "(" << edge_verts[face_edges[i]*2] << ", " <<
1639  edge_verts[face_edges[i]*2+1] << ") ";
1640  }
1641  mfem::out << std::endl;
1642 #endif
1643  }
1644 
1645  // Add top level elements
1646  std::vector<int> tags;
1647  std::vector<int> tris;
1648  std::vector<int> quads;
1649  std::vector<int> tets;
1650  std::vector<int> hexes;
1651  for (int i = 0; i < num_elements; i++)
1652  {
1653  auto etype = mmesh->GetElementType(i);
1654  tags.push_back(mmesh->GetAttribute(i));
1655  switch (etype)
1656  {
1657  case mfem::Element::POINT:
1658  {
1659  // TODO: ?
1660  break;
1661  }
1663  {
1664  // TODO: ?
1665  break;
1666  }
1668  {
1669  mfem::Array<int> eids, oris;
1670  mmesh->GetElementEdges(i, eids, oris);
1671  for (int e = 0; e < 3; e++)
1672  {
1673  tris.push_back(eids[e]);
1674  }
1675  break;
1676  }
1678  {
1679  mfem::Array<int> eids, oris;
1680  mmesh->GetElementEdges(i, eids, oris);
1681  for (int e = 0; e < 4; e++)
1682  {
1683  quads.push_back(eids[e]);
1684  }
1685  break;
1686  }
1688  {
1689  mfem::Array<int> fids, oris;
1690  mmesh->GetElementFaces(i, fids, oris);
1691  for (int f = 0; f < 4; f++)
1692  {
1693  tets.push_back(fids[f]);
1694  }
1695  break;
1696  }
1698  {
1699  mfem::Array<int> fids, oris;
1700  mmesh->GetElementFaces(i, fids, oris);
1701  for (int f = 0; f < 6; f++)
1702  {
1703  hexes.push_back(fids[f]);
1704  }
1705  break;
1706  }
1707  default:
1708  mfem::err << "Error, element not implemented." << std::endl;
1709  return 3;
1710  }
1711  }
1712 
1713  if (tris.size())
1714  {
1715  FmsDomainSetNumEntities(domains[0], FMS_TRIANGLE, FMS_INT32, tris.size() / 3);
1716  FmsDomainAddEntities(domains[0], FMS_TRIANGLE, reorder, FMS_INT32, tris.data(),
1717  tris.size() / 3);
1718 #ifdef DEBUG_MFEM_FMS
1719  mfem::out << "TRIS: ";
1720  for (int i = 0; i < tris.size(); i++)
1721  {
1722  if (i % 3 == 0) { mfem::out << std::endl << "\t" << i/3 << ": "; }
1723  mfem::out << tris[i] << " ";
1724  }
1725  mfem::out << std::endl;
1726 #endif
1727  }
1728 
1729  if (quads.size())
1730  {
1731  // TODO: Not quite right either, if there are hexes and quads then this
1732  // will overwrite the faces that made up the hexes
1733  FmsDomainSetNumEntities(domains[0], FMS_QUADRILATERAL, FMS_INT32,
1734  quads.size() / 4);
1735  FmsDomainAddEntities(domains[0], FMS_QUADRILATERAL, reorder, FMS_INT32,
1736  quads.data(), quads.size() / 4);
1737 #ifdef DEBUG_MFEM_FMS
1738  mfem::out << "QUADS: ";
1739  for (int i = 0; i < quads.size(); i++)
1740  {
1741  if (i % 4 == 0) { mfem::out << std::endl << "\t" << i/4 << ": "; }
1742  mfem::out << quads[i] << " ";
1743  }
1744  mfem::out << std::endl;
1745 #endif
1746  }
1747 
1748  if (tets.size())
1749  {
1750  FmsDomainSetNumEntities(domains[0], FMS_TETRAHEDRON, FMS_INT32,
1751  tets.size() / 4);
1752  FmsDomainAddEntities(domains[0], FMS_TETRAHEDRON, reorder, FMS_INT32,
1753  tets.data(), tets.size() / 4);
1754 #ifdef DEBUG_MFEM_FMS
1755  mfem::out << "TETS: ";
1756  for (int i = 0; i < tets.size(); i++)
1757  {
1758  if (i % 4 == 0) { mfem::out << std::endl << "\t" << i/4 << ": "; }
1759  mfem::out << tets[i] << " ";
1760  }
1761  mfem::out << std::endl;
1762 #endif
1763  }
1764 
1765  if (hexes.size())
1766  {
1767  FmsDomainSetNumEntities(domains[0], FMS_HEXAHEDRON, FMS_INT32,
1768  hexes.size() / 6);
1769  FmsDomainAddEntities(domains[0], FMS_HEXAHEDRON, reorder, FMS_INT32,
1770  hexes.data(), hexes.size() / 6);
1771 #ifdef DEBUG_MFEM_FMS
1772  mfem::out << "HEXES: ";
1773  for (int i = 0; i < hexes.size(); i++)
1774  {
1775  if (i % 6 == 0) { mfem::out << std::endl << "\t" << i/6 << ": "; }
1776  mfem::out << hexes[i] << " ";
1777  }
1778  mfem::out << std::endl;
1779 #endif
1780  }
1781 
1782  err = FmsMeshFinalize(*fmesh);
1783  if (err)
1784  {
1785  mfem::err << "FmsMeshFinalize returned error code " << err << std::endl;
1786  return 4;
1787  }
1788 
1789  err = FmsMeshValidate(*fmesh);
1790  if (err)
1791  {
1792  mfem::err << "FmsMeshValidate returned error code " << err << std::endl;
1793  return 5;
1794  }
1795 
1796  FmsComponent v = NULL;
1797  FmsMeshAddComponent(*fmesh, "volume", &v);
1798  FmsComponentAddDomain(v, domains[0]);
1799 
1800  FmsTag tag;
1801  FmsMeshAddTag(*fmesh, "element_attribute", &tag);
1802  FmsTagSetComponent(tag, v);
1803  FmsTagSet(tag, FMS_INT32, FMS_INT32, tags.data(), tags.size());
1804 
1805  // Add boundary component
1806  std::vector<int> bdr_eles[FMS_NUM_ENTITY_TYPES];
1807  std::vector<int> bdr_attributes;
1808  const int NBE = mmesh->GetNBE();
1809  for (int i = 0; i < NBE; i++)
1810  {
1811  const Element::Type betype = mmesh->GetBdrElementType(i);
1812  bdr_attributes.push_back(mmesh->GetBdrAttribute(i));
1813  switch (betype)
1814  {
1815  case Element::POINT:
1816  bdr_eles[FMS_VERTEX].push_back(mmesh->GetBdrElementEdgeIndex(i));
1817  break;
1818  case Element::SEGMENT:
1819  bdr_eles[FMS_EDGE].push_back(mmesh->GetBdrElementEdgeIndex(i));
1820  break;
1821  case Element::TRIANGLE:
1822  bdr_eles[FMS_TRIANGLE].push_back(mmesh->GetBdrElementEdgeIndex(i));
1823  break;
1825  bdr_eles[FMS_QUADRILATERAL].push_back(mmesh->GetBdrElementEdgeIndex(i));
1826  break;
1827  case Element::TETRAHEDRON:
1828  bdr_eles[FMS_TETRAHEDRON].push_back(mmesh->GetBdrElementEdgeIndex(i));
1829  break;
1830  case Element::HEXAHEDRON:
1831  bdr_eles[FMS_HEXAHEDRON].push_back(mmesh->GetBdrElementEdgeIndex(i));
1832  break;
1833  default:
1834  MFEM_WARNING("Unsupported boundary element " << betype << " at boundary index "
1835  << i);
1836  break;
1837  }
1838  }
1839 
1840  if (NBE)
1841  {
1842  FmsComponent boundary = NULL;
1843  FmsMeshAddComponent(*fmesh, "boundary", &boundary);
1844  FmsInt part_id;
1845  FmsComponentAddPart(boundary, domains[0], &part_id);
1846  for (int i = FMS_NUM_ENTITY_TYPES - 1; i > 0; i--)
1847  {
1848  if (bdr_eles[i].size())
1849  {
1850  FmsComponentAddPartEntities(boundary, part_id, (FmsEntityType)i,
1851  FMS_INT32, FMS_INT32, FMS_INT32, NULL,
1852  bdr_eles[i].data(),
1853  NULL, bdr_eles[i].size());
1854  break;
1855  }
1856  }
1857  FmsComponentAddRelation(v, 1);
1858  FmsTag boundary_tag = NULL;
1859  FmsMeshAddTag(*fmesh, "boundary_attribute", &boundary_tag);
1860  FmsTagSetComponent(boundary_tag, boundary);
1861  FmsTagSet(boundary_tag, FMS_INT32, FMS_INT32, bdr_attributes.data(),
1862  bdr_attributes.size());
1863  }
1864  *volume = v;
1865  return 0;
1866 }
1867 
1868 int
1870  FmsDataCollection *dc)
1871 {
1872  // TODO: Write me
1873  int err = 0;
1874  const Mesh *mmesh = mfem_dc->GetMesh();
1875 
1876  FmsMesh fmesh = NULL;
1877  FmsComponent volume = NULL;
1878  err = MeshToFmsMesh(mmesh, &fmesh, &volume);
1879  if (!fmesh || !volume || err)
1880  {
1881  mfem::err << "Error converting mesh topology from MFEM to FMS" << std::endl;
1882  if (fmesh)
1883  {
1884  FmsMeshDestroy(&fmesh);
1885  }
1886  return 1;
1887  }
1888 
1889  err = FmsDataCollectionCreate(fmesh, mfem_dc->GetCollectionName().c_str(), dc);
1890  if (!*dc || err)
1891  {
1892  mfem::err << "There was an error creating the FMS data collection." <<
1893  std::endl;
1894  FmsMeshDestroy(&fmesh);
1895  if (*dc)
1896  {
1897  FmsDataCollectionDestroy(dc);
1898  }
1899  return 3;
1900  }
1901 
1902  // Add the coordinates field to the data collection
1903  const mfem::GridFunction *mcoords = mmesh->GetNodes();
1904  if (mcoords)
1905  {
1906  FmsField fcoords = NULL;
1907  err = GridFunctionToFmsField(*dc, volume, "CoordsDescriptor", "Coords", mmesh,
1908  mcoords, &fcoords);
1909  err |= FmsComponentSetCoordinates(volume, fcoords);
1910  }
1911  else
1912  {
1913  // Sometimes the nodes are stored as just a vector of vertex coordinates
1914  mfem::Vector mverts;
1915  mmesh->GetVertices(mverts);
1916  FmsFieldDescriptor fdcoords = NULL;
1917  FmsField fcoords = NULL;
1918  err = FmsDataCollectionAddFieldDescriptor(*dc, "CoordsDescriptor", &fdcoords);
1919  err |= FmsFieldDescriptorSetComponent(fdcoords, volume);
1920  err |= FmsFieldDescriptorSetFixedOrder(fdcoords, FMS_CONTINUOUS,
1921  FMS_NODAL_GAUSS_CLOSED, 1);
1922  err |= FmsDataCollectionAddField(*dc, "Coords", &fcoords);
1923  err |= FmsFieldSet(fcoords, fdcoords, mmesh->SpaceDimension(), FMS_BY_NODES,
1924  FMS_DOUBLE, mverts);
1925  err |= FmsComponentSetCoordinates(volume, fcoords);
1926  }
1927 
1928  if (err)
1929  {
1930  mfem::err << "There was an error setting the mesh coordinates." << std::endl;
1931  FmsMeshDestroy(&fmesh);
1932  FmsDataCollectionDestroy(dc);
1933  return 4;
1934  }
1935 
1936  const auto &fields = mfem_dc->GetFieldMap();
1937  for (const auto &pair : fields)
1938  {
1939  std::string fd_name(pair.first + "Descriptor");
1940  FmsField field;
1941  err = GridFunctionToFmsField(*dc, volume, fd_name, pair.first.c_str(), mmesh,
1942  pair.second, &field);
1943  if (err)
1944  {
1945  mfem::err << "WARNING: There was an error adding the " << pair.first <<
1946  " field. Continuing..." << std::endl;
1947  }
1948  }
1949 
1950  // /* TODO:
1951  // const auto &qfields = mfem_dc->GetQFieldMap();
1952  // for(const auto &pair : qfields) {
1953  // FmsFieldDescriptor fd = NULL;
1954  // FmsField f = NULL;
1955  // std::string fd_name(pair.first + "Collection");
1956  // FmsDataCollectionAddFieldDescriptor(*dc, fd_name.c_str(), &fd);
1957  // FmsDataCollectionAddField(*dc, pair.first.c_str(), &f);
1958  // GridFunctionToFmsField(*dc, fd, f, volume, pair.second); // TODO: Volume isn't always going to be correct
1959  // } */
1960 
1961  MfemMetaDataToFmsMetaData(mfem_dc, *dc);
1962 
1963  return 0;
1964 }
1965 
1966 } // namespace mfem
1967 #endif
Normal component of vector field.
Definition: fe_coll.hpp:47
Field is discontinuous across element interfaces.
Definition: fe_coll.hpp:48
bool FmsMetaDataGetInteger(FmsMetaData mdata, const std::string &key, std::vector< int > &values)
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:552
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:1203
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:233
Table * GetEdgeVertexTable() const
Returns the edge-to-vertex Table (3D)
Definition: mesh.cpp:5573
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1324
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void GetEdgeVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of edge i.
Definition: mesh.cpp:5536
double GetTime() const
Get physical time (for time-dependent simulations)
virtual int GetContType() const =0
Tangential components of vector field.
Definition: fe_coll.hpp:46
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:849
bool FmsMetaDataGetScalar(FmsMetaData mdata, const std::string &key, std::vector< double > &values)
bool MfemMetaDataToFmsMetaData(DataCollection *mdc, FmsDataCollection fdc)
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:5748
const std::string & GetCollectionName() const
Get the name of the collection.
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: mesh.cpp:5753
T * GetData()
Returns the data.
Definition: array.hpp:108
int AddTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:1310
int VectorDim() const
Definition: gridfunc.cpp:311
void GetFaceVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of face i.
Definition: mesh.hpp:1029
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition: table.cpp:191
int GetBdrElementEdgeIndex(int i) const
Definition: mesh.cpp:5714
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
virtual void Save()
Save the collection and a VisIt root file.
void GetVertices(Vector &vert_coord) const
Definition: mesh.cpp:7272
static const char * Name(int b_type)
Check and convert a BasisType constant to a string identifier.
Definition: fe.hpp:94
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:199
int FmsMeshToMesh(FmsMesh fms_mesh, Mesh **mfem_mesh)
Definition: fmsconvert.cpp:453
Nodes: x_i = i/(n-1), i=0,...,n-1.
Definition: fe.hpp:37
bool FmsMetaDataGetString(FmsMetaData mdata, const std::string &key, std::string &value)
double GetTimeStep() const
Get the simulation time step (for time-dependent simulations)
Nodes: x_i = (i+1/2)/n, i=0,...,n-1.
Definition: fe.hpp:38
void GetFaceInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:2662
virtual int DofForGeometry(Geometry::Type GeomType) const =0
int AddBdrTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:1454
int AddVertex(double x, double y=0.0, double z=0.0)
Definition: mesh.cpp:1263
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
Bernstein polynomials.
Definition: fe.hpp:35
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:746
double b
Definition: lissajous.cpp:42
int AddTet(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1338
const FieldMapType & GetFieldMap() const
Get a const reference to the internal field map.
Closed GaussLegendre.
Definition: fe.hpp:40
int GetCycle() const
Get time cycle (for time-dependent simulations)
Data collection with VisIt I/O routines.
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:4882
Type
Constants for the classes derived from Element.
Definition: element.hpp:41
int AddBdrSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1440
int DataCollectionToFmsDataCollection(DataCollection *mfem_dc, FmsDataCollection *dc)
int FmsDataCollectionToDataCollection(FmsDataCollection dc, DataCollection **mfem_dc)
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:92
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:534
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:629
void SetTime(double t)
Set physical time (for time-dependent simulations)
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:341
void SetOwnData(bool o)
Set the ownership of collection data.
void GetElementInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:2629
int SpaceDimension() const
Definition: mesh.hpp:912
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const =0
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Definition: globals.hpp:71
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
int AddHex(int v1, int v2, int v3, int v4, int v5, int v6, int v7, int v8, int attr=1)
Definition: mesh.cpp:1373
void FinalizeTopology(bool generate_bdr=true)
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
Definition: mesh.cpp:2507
int AddBdrQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1468
int GetOrder(int i) const
Returns the polynomial degree of the i&#39;th finite element.
Definition: fespace.hpp:524
int GridFunctionToFmsField(FmsDataCollection dc, FmsComponent comp, const std::string &fd_name, const std::string &field_name, const Mesh *mesh, const GridFunction *gf, FmsField *outfield)
Definition: fmsconvert.cpp:949
Field is continuous across element interfaces.
Definition: fe_coll.hpp:45
bool BasisTypeToFmsBasisType(int bt, FmsBasisType &btype)
Definition: fmsconvert.cpp:893
Nodes: x_i = (i+1)/(n+1), i=0,...,n-1.
Definition: fe.hpp:36
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:843
void SetTimeStep(double ts)
Set the simulation time step (for time-dependent simulations)
virtual const char * Name() const
Definition: fe_coll.hpp:61
int MeshToFmsMesh(const Mesh *mmesh, FmsMesh *fmesh, FmsComponent *volume)
void GetElementEdges(int i, Array< int > &edges, Array< int > &cor) const
Return the indices and the orientations of all edges of element i.
Definition: mesh.cpp:5452
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:2600
void GetEdgeInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:2650
int dim
Definition: ex24.cpp:53
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:193
Mesh * GetMesh()
Get a pointer to the mesh in the collection.
int GetNEdges() const
Return the number of edges.
Definition: mesh.hpp:852
void MakeRef(T *, int)
Make this Array a reference to a pointer.
Definition: array.hpp:859
int GetNFaces() const
Return the number of faces in a 3D mesh.
Definition: mesh.hpp:855
Vector data type.
Definition: vector.hpp:60
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:554
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7343
int GetId(int p1, int p2)
Get id of item whose parents are p1, p2... Create it if it doesn&#39;t exist.
Definition: hash.hpp:371
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
Table * GetFaceEdgeTable() const
Returns the face-to-edge Table (3D)
Definition: mesh.cpp:5545
int RowSize(int i) const
Definition: table.hpp:108
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
Serendipity basis (squares / cubes)
Definition: fe.hpp:39
int FindId(int p1, int p2) const
Find id of item whose parents are p1, p2... Return -1 if it doesn&#39;t exist.
Definition: hash.hpp:462
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:1197
int FmsFieldToGridFunction(FmsMesh fms_mesh, FmsField f, Mesh *mesh, GridFunction &func, bool setFE)
This function converts an FmsField to an MFEM GridFunction.
Definition: fmsconvert.cpp:107
void GetElementFaces(int i, Array< int > &faces, Array< int > &ori) const
Return the indices and the orientations of all faces of element i.
Definition: mesh.cpp:5668
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:285