Appendix

Mint Application Code Example

Below is the complete Mint Application Code Example presented in the Getting Started with Mint section. The code can be found in the Axom source code under src/axom/mint/examples/user_guide/mint_getting_started.cpp.

  1// sphinx_tutorial_walkthrough_includes_start
  2
  3#include "axom/config.hpp"                          // compile time definitions
  4#include "axom/core/execution/execution_space.hpp"  // for execution_space traits
  5
  6#include "axom/mint.hpp"                  // for mint classes and functions
  7#include "axom/core/numerics/Matrix.hpp"  // for numerics::Matrix
  8
  9// sphinx_tutorial_walkthrough_includes_end
 10
 11// namespace aliases
 12namespace mint = axom::mint;
 13namespace numerics = axom::numerics;
 14namespace xargs = mint::xargs;
 15
 16using IndexType = axom::IndexType;
 17
 18// compile-time switch for execution policy
 19#if defined(AXOM_USE_RAJA) && defined(AXOM_USE_CUDA) && defined(AXOM_USE_UMPIRE)
 20constexpr int NUM_BLOCKS = 512;
 21using ExecPolicy = axom::CUDA_EXEC<NUM_BLOCKS>;
 22#elif defined(AXOM_USE_RAJA) && defined(AXOM_USE_HIP) && defined(AXOM_USE_UMPIRE)
 23constexpr int NUM_BLOCKS = 512;
 24using ExecPolicy = axom::HIP_EXEC<NUM_BLOCKS>;
 25#elif defined(AXOM_USE_RAJA) && defined(AXOM_USE_OPENMP)
 26using ExecPolicy = axom::OMP_EXEC;
 27#else
 28using ExecPolicy = axom::SEQ_EXEC;
 29#endif
 30
 31constexpr IndexType NUM_NODES_PER_CELL = 4;
 32constexpr double ONE_OVER_4 = 1. / static_cast<double>(NUM_NODES_PER_CELL);
 33
 34/*!
 35 * \brief Holds command-line arguments
 36 */
 37static struct
 38{
 39  int res;
 40  bool useUnstructured;
 41} Arguments;
 42
 43//------------------------------------------------------------------------------
 44// FUNCTION PROTOTYPES
 45//------------------------------------------------------------------------------
 46void parse_args(int argc, char** argv);
 47mint::Mesh* getUniformMesh();
 48mint::Mesh* getUnstructuredMesh();
 49
 50//------------------------------------------------------------------------------
 51// PROGRAM MAIN
 52//------------------------------------------------------------------------------
 53int main(int argc, char** argv)
 54{
 55  parse_args(argc, argv);
 56
 57  // sphinx_tutorial_walkthrough_set_memory_start
 58  // NOTE: use unified memory if we are using CUDA or HIP
 59  const int allocID = axom::execution_space<ExecPolicy>::allocatorID();
 60  axom::setDefaultAllocator(allocID);
 61  // sphinx_tutorial_walkthrough_set_memory_end
 62
 63  // sphinx_tutorial_walkthrough_construct_mesh_start
 64
 65  mint::Mesh* mesh = (Arguments.useUnstructured) ? getUnstructuredMesh() : getUniformMesh();
 66
 67  // sphinx_tutorial_walkthrough_construct_mesh_end
 68
 69  // sphinx_tutorial_walkthrough_add_fields_start
 70
 71  // add a cell-centered and a node-centered field
 72  double* phi = mesh->createField<double>("phi", mint::NODE_CENTERED);
 73  double* hc = mesh->createField<double>("hc", mint::CELL_CENTERED);
 74
 75  constexpr int NUM_COMPONENTS = 2;
 76  double* xc = mesh->createField<double>("xc", mint::CELL_CENTERED, NUM_COMPONENTS);
 77
 78  // sphinx_tutorial_walkthrough_add_fields_end
 79
 80  // sphinx_tutorial_walkthrough_compute_hf_start
 81
 82  // loop over the nodes and evaluate Himmelblaus Function
 83  mint::for_all_nodes<ExecPolicy, xargs::xy>(
 84    mesh,
 85    AXOM_LAMBDA(IndexType nodeIdx, double x, double y) {
 86      const double x_2 = x * x;
 87      const double y_2 = y * y;
 88      const double A = x_2 + y - 11.0;
 89      const double B = x + y_2 - 7.0;
 90
 91      phi[nodeIdx] = A * A + B * B;
 92    });
 93
 94  // sphinx_tutorial_walkthrough_compute_hf_end
 95
 96  // sphinx_tutorial_walkthrough_cell_centers_start
 97
 98  // loop over cells and compute cell centers
 99  mint::for_all_cells<ExecPolicy, xargs::coords>(
100    mesh,
101    AXOM_LAMBDA(IndexType cellIdx, const numerics::Matrix<double>& coords, const IndexType* nodeIds) {
102      // NOTE: A column vector of the coords matrix corresponds to a nodes coords
103
104      // Sum the cell's nodal coordinates
105      double xsum = 0.0;
106      double ysum = 0.0;
107      double hsum = 0.0;
108
109      const IndexType numNodes = coords.getNumColumns();
110      for(IndexType inode = 0; inode < numNodes; ++inode)
111      {
112        const double* node = coords.getColumn(inode);
113        xsum += node[mint::X_COORDINATE];
114        ysum += node[mint::Y_COORDINATE];
115
116        hsum += phi[nodeIds[inode]];
117      }  // END for all cell nodes
118
119      // compute cell centroid by averaging the nodal coordinate sums
120      const IndexType offset = cellIdx * NUM_COMPONENTS;
121      const double invnnodes = 1.f / static_cast<double>(numNodes);
122      xc[offset] = xsum * invnnodes;
123      xc[offset + 1] = ysum * invnnodes;
124
125      hc[cellIdx] = hsum * invnnodes;
126    });
127
128  // sphinx_tutorial_walkthrough_cell_centers_end
129
130  // sphinx_tutorial_walkthrough_vtk_start
131
132  // write the mesh in a VTK file for visualization
133  std::string vtkfile = (Arguments.useUnstructured) ? "unstructured_mesh.vtk" : "uniform_mesh.vtk";
134  mint::write_vtk(mesh, vtkfile);
135
136  // sphinx_tutorial_walkthrough_vtk_end
137
138  delete mesh;
139  mesh = nullptr;
140
141  return 0;
142}
143
144//------------------------------------------------------------------------------
145//  FUNCTION PROTOTYPE IMPLEMENTATION
146//------------------------------------------------------------------------------
147void parse_args(int argc, char** argv)
148{
149  Arguments.res = 25;
150  Arguments.useUnstructured = false;
151
152  for(int i = 1; i < argc; ++i)
153  {
154    if(strcmp(argv[i], "--unstructured") == 0)
155    {
156      Arguments.useUnstructured = true;
157    }
158
159    else if(strcmp(argv[i], "--resolution") == 0)
160    {
161      Arguments.res = std::atoi(argv[++i]);
162    }
163
164  }  // END for all arguments
165
166  SLIC_ERROR_IF(Arguments.res < 2, "invalid mesh resolution! Please, pick a value greater than 2.");
167}
168
169//------------------------------------------------------------------------------
170// sphinx_tutorial_walkthrough_construct_umesh_start
171mint::Mesh* getUniformMesh()
172{
173  // construct a N x N grid within a domain defined in [-5.0, 5.0]
174  const double lo[] = {-5.0, -5.0};
175  const double hi[] = {5.0, 5.0};
176  mint::Mesh* m = new mint::UniformMesh(lo, hi, Arguments.res, Arguments.res);
177  return m;
178}
179// sphinx_tutorial_walkthrough_construct_umesh_end
180
181//------------------------------------------------------------------------------
182mint::Mesh* getUnstructuredMesh()
183{
184  mint::Mesh* umesh = getUniformMesh();
185  const IndexType umesh_ncells = umesh->getNumberOfCells();
186  const IndexType umesh_nnodes = umesh->getNumberOfNodes();
187
188  const IndexType ncells = umesh_ncells * 4;  // split each quad into 4 triangles
189  const IndexType nnodes = umesh_nnodes + umesh_ncells;
190
191  constexpr int DIMENSION = 2;
192  using MeshType = mint::UnstructuredMesh<mint::SINGLE_SHAPE>;
193  MeshType* m = new MeshType(DIMENSION, mint::TRIANGLE, nnodes, ncells);
194  m->resize(nnodes, ncells);
195
196  double* x = m->getCoordinateArray(mint::X_COORDINATE);
197  double* y = m->getCoordinateArray(mint::Y_COORDINATE);
198  IndexType* cells = m->getCellNodesArray();
199
200  // fill coordinates from uniform mesh
201  mint::for_all_nodes<ExecPolicy, xargs::xy>(
202    umesh,
203    AXOM_LAMBDA(IndexType nodeIdx, double nx, double ny) {
204      x[nodeIdx] = nx;
205      y[nodeIdx] = ny;
206    });
207
208  // loop over cells, compute cell centers and fill connectivity
209  mint::for_all_cells<ExecPolicy, xargs::coords>(
210    umesh,
211    AXOM_LAMBDA(IndexType cellIdx, const numerics::Matrix<double>& coords, const IndexType* nodeIds) {
212      // NOTE: A column vector of the coords matrix corresponds to a nodes coords
213
214      // Sum the cell's nodal coordinates
215      double xsum = 0.0;
216      double ysum = 0.0;
217      for(IndexType inode = 0; inode < NUM_NODES_PER_CELL; ++inode)
218      {
219        const double* node = coords.getColumn(inode);
220        xsum += node[mint::X_COORDINATE];
221        ysum += node[mint::Y_COORDINATE];
222      }  // END for all cell nodes
223
224      // compute cell centroid by averaging the nodal coordinate sums
225      const IndexType nc = umesh_nnodes + cellIdx; /* centroid index */
226      x[nc] = xsum * ONE_OVER_4;
227      y[nc] = ysum * ONE_OVER_4;
228
229      // triangulate
230      const IndexType& n0 = nodeIds[0];
231      const IndexType& n1 = nodeIds[1];
232      const IndexType& n2 = nodeIds[2];
233      const IndexType& n3 = nodeIds[3];
234
235      const IndexType offset = cellIdx * 12;
236
237      cells[offset] = n0;
238      cells[offset + 1] = nc;
239      cells[offset + 2] = n3;
240
241      cells[offset + 3] = n0;
242      cells[offset + 4] = n1;
243      cells[offset + 5] = nc;
244
245      cells[offset + 6] = n1;
246      cells[offset + 7] = n2;
247      cells[offset + 8] = nc;
248
249      cells[offset + 9] = n2;
250      cells[offset + 10] = n3;
251      cells[offset + 11] = nc;
252    });
253
254  // delete uniform mesh
255  delete umesh;
256  umesh = nullptr;
257
258  return m;
259}

AXOM_LAMBDA Macro

The AXOM_LAMBDA convenience macro expands to:

  • [=] capture by value when the Axom Toolkit is compiled without CUDA.

  • [=] __host__ __device__ when the Axom Toolkit is compiled with CUDA

Raw Sidre Data

  1{
  2  "name": "",
  3  "groups": 
  4  {
  5    "state": 
  6    {
  7      "name": "state",
  8      "groups": 
  9      {
 10        "t1": 
 11        {
 12          "name": "t1",
 13          "views": 
 14          {
 15            "block_id": 
 16            {
 17              "name": "block_id",
 18              "schema": "{\"dtype\":\"int32\", \"number_of_elements\": 1, \"offset\": 0, \"stride\": 4, \"element_bytes\": 4, \"endianness\": \"little\"}",
 19              "value": "-1",
 20              "state": "SCALAR",
 21              "is_applied": 1
 22            },
 23            "partition_id": 
 24            {
 25              "name": "partition_id",
 26              "schema": "{\"dtype\":\"int32\", \"number_of_elements\": 1, \"offset\": 0, \"stride\": 4, \"element_bytes\": 4, \"endianness\": \"little\"}",
 27              "value": "-1",
 28              "state": "SCALAR",
 29              "is_applied": 1
 30            }
 31          }
 32        }
 33      }
 34    },
 35    "coordsets": 
 36    {
 37      "name": "coordsets",
 38      "groups": 
 39      {
 40        "c1": 
 41        {
 42          "name": "c1",
 43          "views": 
 44          {
 45            "type": 
 46            {
 47              "name": "type",
 48              "schema": "{\"dtype\":\"char8_str\", \"number_of_elements\": 9, \"offset\": 0, \"stride\": 1, \"element_bytes\": 1, \"endianness\": \"little\"}",
 49              "value": "\"explicit\"",
 50              "state": "STRING",
 51              "is_applied": 1
 52            }
 53          },
 54          "groups": 
 55          {
 56            "values": 
 57            {
 58              "name": "values",
 59              "views": 
 60              {
 61                "x": 
 62                {
 63                  "name": "x",
 64                  "schema": "{\"dtype\":\"float64\", \"number_of_elements\": 6, \"offset\": 0, \"stride\": 8, \"element_bytes\": 8, \"endianness\": \"little\"}",
 65                  "value": "[0.0, 2.0, 1.0, 3.5, 2.5, 5.0]",
 66                  "state": "BUFFER",
 67                  "is_applied": 1
 68                },
 69                "y": 
 70                {
 71                  "name": "y",
 72                  "schema": "{\"dtype\":\"float64\", \"number_of_elements\": 6, \"offset\": 0, \"stride\": 8, \"element_bytes\": 8, \"endianness\": \"little\"}",
 73                  "value": "[0.0, 0.0, 1.0, 1.0, 2.0, 0.0]",
 74                  "state": "BUFFER",
 75                  "is_applied": 1
 76                }
 77              }
 78            }
 79          }
 80        }
 81      }
 82    },
 83    "topologies": 
 84    {
 85      "name": "topologies",
 86      "groups": 
 87      {
 88        "t1": 
 89        {
 90          "name": "t1",
 91          "views": 
 92          {
 93            "coordset": 
 94            {
 95              "name": "coordset",
 96              "schema": "{\"dtype\":\"char8_str\", \"number_of_elements\": 3, \"offset\": 0, \"stride\": 1, \"element_bytes\": 1, \"endianness\": \"little\"}",
 97              "value": "\"c1\"",
 98              "state": "STRING",
 99              "is_applied": 1
100            },
101            "type": 
102            {
103              "name": "type",
104              "schema": "{\"dtype\":\"char8_str\", \"number_of_elements\": 13, \"offset\": 0, \"stride\": 1, \"element_bytes\": 1, \"endianness\": \"little\"}",
105              "value": "\"unstructured\"",
106              "state": "STRING",
107              "is_applied": 1
108            }
109          },
110          "groups": 
111          {
112            "elements": 
113            {
114              "name": "elements",
115              "views": 
116              {
117                "shape": 
118                {
119                  "name": "shape",
120                  "schema": "{\"dtype\":\"char8_str\", \"number_of_elements\": 4, \"offset\": 0, \"stride\": 1, \"element_bytes\": 1, \"endianness\": \"little\"}",
121                  "value": "\"tri\"",
122                  "state": "STRING",
123                  "is_applied": 1
124                },
125                "connectivity": 
126                {
127                  "name": "connectivity",
128                  "schema": "{\"dtype\":\"int32\", \"number_of_elements\": 12, \"offset\": 0, \"stride\": 4, \"element_bytes\": 4, \"endianness\": \"little\"}",
129                  "value": "[1, 3, 2, 2, 0, 1, 3, 4, 2, 1, 5, 3]",
130                  "state": "BUFFER",
131                  "is_applied": 1
132                },
133                "stride": 
134                {
135                  "name": "stride",
136                  "schema": "{\"dtype\":\"int32\", \"number_of_elements\": 1, \"offset\": 0, \"stride\": 4, \"element_bytes\": 4, \"endianness\": \"little\"}",
137                  "value": "3",
138                  "state": "SCALAR",
139                  "is_applied": 1
140                }
141              }
142            }
143          }
144        }
145      }
146    },
147    "fields": 
148    {
149      "name": "fields",
150      "groups": 
151      {
152        "den": 
153        {
154          "name": "den",
155          "views": 
156          {
157            "association": 
158            {
159              "name": "association",
160              "schema": "{\"dtype\":\"char8_str\", \"number_of_elements\": 8, \"offset\": 0, \"stride\": 1, \"element_bytes\": 1, \"endianness\": \"little\"}",
161              "value": "\"element\"",
162              "state": "STRING",
163              "is_applied": 1
164            },
165            "volume_dependent": 
166            {
167              "name": "volume_dependent",
168              "schema": "{\"dtype\":\"char8_str\", \"number_of_elements\": 5, \"offset\": 0, \"stride\": 1, \"element_bytes\": 1, \"endianness\": \"little\"}",
169              "value": "\"true\"",
170              "state": "STRING",
171              "is_applied": 1
172            },
173            "topology": 
174            {
175              "name": "topology",
176              "schema": "{\"dtype\":\"char8_str\", \"number_of_elements\": 3, \"offset\": 0, \"stride\": 1, \"element_bytes\": 1, \"endianness\": \"little\"}",
177              "value": "\"t1\"",
178              "state": "STRING",
179              "is_applied": 1
180            },
181            "values": 
182            {
183              "name": "values",
184              "schema": "{\"dtype\":\"float64\", \"number_of_elements\": 4, \"offset\": 0, \"stride\": 8, \"element_bytes\": 8, \"endianness\": \"little\"}",
185              "value": "[0.5, 1.2, 2.5, 0.9]",
186              "state": "BUFFER",
187              "is_applied": 1
188            }
189          }
190        }
191      }
192    }
193  }
194}