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