XGC1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
grid_files.hpp
Go to the documentation of this file.
1 #ifndef GRID_FILES_HPP
2 #define GRID_FILES_HPP
3 
4 #ifdef USE_MPI
5 #include "my_mpi.hpp"
6 #endif
7 #include "file_reader.hpp"
8 #include "magnetic_field.hpp"
9 #include "grid_structs.hpp"
10 #include "grid_setup.hpp"
11 #include "space_settings.hpp"
12 #include "NamelistReader.hpp"
13 #include "check_aif_flx_format.hpp"
14 
15 extern "C" void set_fortran_flx_aif_format(int is_flx_aif_format_int);
16 
17 struct GridFiles{
18  // Contents of .node
19  int nnodes;
20  View<RZPair*, HostType> x;
21  View<int*, HostType> rgn;
22 
23  // Contents of .ele
25  View<Vertex*, HostType> nodes;
26 
27  // Contents of .flx.aif
28  View<int*, HostType> i_x;
29 
30  int nsurfaces;
32  View<int*, HostType> i_surf_separatrices;
33  View<int*, HostType> nnodes_on_surface;
36  View<int*, HostType> non_aligned_vert;
37  View<int*, HostType> non_aligned_nsurf;
38  View<int**, CLayout, HostType> non_aligned_surf_idx;
39  View<int**, CLayout, HostType> surf_idx;
40 
41  void convert_to_x_rgn_format(const View<double**, CLayout, HostType>& tmp_nodes, View<RZPair*, HostType>& x, View<int*, HostType>& rgn){
42  Kokkos::parallel_for("transpose_x_and_rgn", Kokkos::RangePolicy<HostExSpace>(0,tmp_nodes.extent(0)), KOKKOS_LAMBDA( const int i ){
43  // tmp_nodes[i,0]; <-- first column is ignored
44  x(i).r = tmp_nodes(i,1);
45  x(i).z = tmp_nodes(i,2);
46  rgn(i) = tmp_nodes(i,3);
47  });
48  }
49 
50  void convert_to_node_format(const View<int**, CLayout, HostType>& tmp_ele, View<Vertex*, HostType>& nodes){
51  Kokkos::parallel_for("transpose_x_and_rgn", Kokkos::RangePolicy<HostExSpace>(0,tmp_ele.extent(0)), KOKKOS_LAMBDA( const int i ){
52  // tmp_ele[i,0]; <-- first column is ignored
53  nodes(i).vertex[0] = tmp_ele(i,1);
54  nodes(i).vertex[1] = tmp_ele(i,2);
55  nodes(i).vertex[2] = tmp_ele(i,3);
56  });
57  }
58 
59  // Default constructor, to use to deallocate the file arrays after use
61 
62  // Read the files
63  GridFiles(const std::string& node_file, const std::string& element_file, const std::string& flx_aif_file){
64  int max_nodes_on_surface;
65  int n_xpts;
66  int n_separatrices;
67  // Only one MPI rank reads the files, then broadcasts to the others
68  if(is_rank_zero()){
69 
70  FileReader file_reader;
71 
73  file_reader.open(node_file);
74  file_reader.read(&nnodes, 1); // Get nnodes from file
75  file_reader.skip(3); // Skip 3 unused values in file
76 
77  // Use dummy array since reformatting is needed
78  View<double**, CLayout, HostType> tmp_nodes(NoInit("tmp_nodes"),nnodes,4); // 4 entries per node: index, r, z, region
79  file_reader.read(tmp_nodes.data(), tmp_nodes.size());
80 
81  // Reformat to x and rgn
82  x = View<RZPair*, HostType>(NoInit("x"),nnodes);
83  rgn = View<int*, HostType>(NoInit("rgn"),nnodes);
84  convert_to_x_rgn_format(tmp_nodes, x, rgn);
85 
87  file_reader.open(element_file);
88  file_reader.read(&ntriangles, 1); // Get nnodes from file
89  file_reader.skip(2); // Skip 2 unused values in file
90 
91  // Use dummy array since reformatting is needed
92  View<int**, CLayout, HostType> tmp_ele(NoInit("tmp_ele"),ntriangles,4); // 4 entries per node: index, r, z, region
93  file_reader.read(tmp_ele.data(), tmp_ele.size());
94 
95  // Reformat since first column is ignored
96  nodes = View<Vertex*, HostType>("nodes",ntriangles);
97  convert_to_node_format(tmp_ele, nodes);
98 
100  if(nnodes==1) exit_XGC("\nGrid must have more than 1 nnode!\n"); // This is assumed in the flx format check
101  bool old_flx_aif_format = is_old_flx_format(flx_aif_file);
102  set_fortran_flx_aif_format(old_flx_aif_format ? 1 : 0);
103 
104  file_reader.open(flx_aif_file);
105 
106  if(old_flx_aif_format){
107  file_reader.read(&nsurfaces, 1); // Get number of surfaces
108  n_xpts = 0;
109  n_separatrices = 0;
110  i_x = View<int*, HostType>("i_x",0); // No xpoints specified
111  i_surf_separatrices = View<int*, HostType>(NoInit("i_surf_separatrices"), 0); // No separatrices specified
112 
113  // These don't get read or used by the old format
114  nsurfaces1 = 0;
115  nsurfaces2 = 0;
116  nsurfaces3a = 0;
117  nsurfaces3b = 0;
118  }else{
119  file_reader.read(&n_xpts, 1); // Get number of x-points
120  i_x = View<int*, HostType>("i_x", n_xpts);
121 
122  // Current flx.aif files list 2 x-point indices, regardless of n_xpts
123  int i_xpt1, i_xpt2;
124  file_reader.read(&i_xpt1, 1); // Get index of 1st x-point
125  file_reader.read(&i_xpt2, 1); // Get index of 2nd x-point
126 
127  // Check for consistency
128  if(n_xpts==1 && i_xpt2 != -1) exit_XGC("Error in flx.aif file: n_xpts=1, but i_xpt2 is specified.\n");
129  if(n_xpts==2 && i_xpt2 <= 0) exit_XGC("Error in flx.aif file: n_xpts=2, but i_xpt2 is not specified.\n");
130 
131  if(n_xpts>=1) i_x(0) = i_xpt1 - 1; // switch to zero-index
132  if(n_xpts>=2) i_x(1) = i_xpt2 - 1;
133 
134  // Next are the nsurface variables
135  file_reader.read(&nsurfaces1, 1); // Get number of surfaces
136  file_reader.read(&nsurfaces2, 1); // Get number of surfaces
137  file_reader.read(&nsurfaces3a, 1); // Get number of surfaces
138  file_reader.read(&nsurfaces3b, 1); // Get number of surfaces
139 
140  // nsurfaces is actually the sum of these 4 inputs
142 
143  // Get surface index of separatrices
144  // Hard-coded to 2, but generalize the read just in case
145  n_separatrices = 2;
146  i_surf_separatrices = View<int*, HostType>(NoInit("i_surf_separatrices"), n_separatrices);
147  file_reader.read(&i_surf_separatrices(0), 1);
148  file_reader.read(&i_surf_separatrices(1), 1);
149  }
150 
151  // Read in nnodes on each surface
152  nnodes_on_surface = View<int*, HostType>(NoInit("nnodes_on_surface"), nsurfaces);
153  file_reader.read(nnodes_on_surface.data(), nnodes_on_surface.size());
154 
155  // Loop through surfaces to find max number of nodes in a surface
156  max_nodes_on_surface = 0;
157  for(int i=0; i<nnodes_on_surface.size(); i++) max_nodes_on_surface = std::max(nnodes_on_surface(i), max_nodes_on_surface);
158 
159  // Allocate surf_idx and initialize to 0
160  surf_idx = View<int**, CLayout, HostType>("surf_idx",nsurfaces, max_nodes_on_surface);
161 
162  // Read in values
163  for(int i=0; i<nnodes_on_surface.size(); i++){
164  file_reader.read(&surf_idx(i,0), nnodes_on_surface(i));
165  }
166 
167  // There is a "-1" in the file here to mark the end of the section
168  int should_be_negative_one;
169  file_reader.read(&should_be_negative_one, 1);
170  if(should_be_negative_one != -1){
171  printf("\nExpected -1 at end of section in surface file, encountered %d instead\n", should_be_negative_one);
172  exit_XGC("Error reading surface file.\n");
173  }
174 
175  file_reader.read(&num_non_aligned, 1);
176  printf(" - reading num_non_aligned: %d",num_non_aligned);
177 
178  // Read the non-aligned vertex data
179  non_aligned_vert = View<int*, HostType>(NoInit("non_aligned_vert"), num_non_aligned);
180  non_aligned_nsurf = View<int*, HostType>(NoInit("non_aligned_nsurf"), num_non_aligned);
181  // Initialize to 0
182  non_aligned_surf_idx = View<int**, CLayout, HostType>("non_aligned_surf_idx", num_non_aligned, 4);
183  for(int i=0; i<num_non_aligned; i++){
184  file_reader.read(&non_aligned_vert(i), 1);
185  if(old_flx_aif_format){
186  file_reader.read(&non_aligned_nsurf(i), 1);
187  }else{
188  // New format does not contain non_aligned_nsurf
189  non_aligned_nsurf(i)=2;
190  }
191  file_reader.read(&non_aligned_surf_idx(i,0), non_aligned_nsurf(i));
192  }
193 
194  // Another "-1" at the end of the file
195  file_reader.read(&should_be_negative_one, 1);
196  if(should_be_negative_one != -1){
197  printf("\nExpected -1 at end of the surface file, encountered %d instead\n", should_be_negative_one);
198  exit_XGC("Error reading surface file.\n");
199  }
200  }
201 
202  /*** Broadcast grid contents to other MPI ranks ***/
203 
204  // First broadcast scalars so we know the sizes of the arrays
205  int n_scalars = 11;
206  View<int*, HostType> scalars("scalars", n_scalars);
207  scalars(0) = nnodes;
208  scalars(1) = ntriangles;
209  scalars(2) = nsurfaces;
210  scalars(3) = max_nodes_on_surface;
211  scalars(4) = num_non_aligned;
212  scalars(5) = n_xpts;
213  scalars(6) = n_separatrices;
214  scalars(7) = nsurfaces1;
215  scalars(8) = nsurfaces2;
216  scalars(9) = nsurfaces3a;
217  scalars(10) = nsurfaces3b;
218 
219  // Broadcast file size to other processes
220 #ifdef USE_MPI
221  int root_rank = 0;
222  MPI_Bcast(scalars.data(), scalars.size(), MPI_INT, root_rank, SML_COMM_WORLD);
223 #endif
224 
225  if(!is_rank_zero()){
226  // Set scalars on non-root ranks
227  nnodes = scalars(0);
228  ntriangles = scalars(1);
229  nsurfaces = scalars(2);
230  max_nodes_on_surface = scalars(3);
231  num_non_aligned = scalars(4);
232  n_xpts = scalars(5);
233  n_separatrices = scalars(6);
234  nsurfaces1 = scalars(7);
235  nsurfaces2 = scalars(8);
236  nsurfaces3a = scalars(9);
237  nsurfaces3b = scalars(10);
238 
240  // node and ele views
241  x = View<RZPair*, HostType>(NoInit("x"),nnodes);
242  rgn = View<int*, HostType>(NoInit("rgn"),nnodes);
243  nodes = View<Vertex*, HostType>("nodes",ntriangles);
244 
245  // flx.aif views
246  i_x = View<int*, HostType>("i_x", n_xpts);
247  i_surf_separatrices = View<int*, HostType>(NoInit("i_surf_separatrices"), n_separatrices);
248  non_aligned_vert = View<int*, HostType>(NoInit("non_aligned_vert"), num_non_aligned);
249  non_aligned_nsurf = View<int*, HostType>(NoInit("non_aligned_nsurf"), num_non_aligned);
250  non_aligned_surf_idx = View<int**, CLayout, HostType>(NoInit("non_aligned_surf_idx"), num_non_aligned, 4);
251  nnodes_on_surface = View<int*, HostType>(NoInit("nnodes_on_surface"), nsurfaces);
252  surf_idx = View<int**, CLayout, HostType>(NoInit("surf_idx"),nsurfaces, max_nodes_on_surface);
253  }
254 
255  // Broadcast grid arrays
256 #ifdef USE_MPI
257  MPI_Bcast( x.data(), x.size()*2, MPI_DOUBLE, root_rank, SML_COMM_WORLD);
258  MPI_Bcast( rgn.data(), rgn.size(), MPI_INT, root_rank, SML_COMM_WORLD);
259  MPI_Bcast( nodes.data(), nodes.size()*3, MPI_INT, root_rank, SML_COMM_WORLD);
260 
261  if(i_x.size()>0) MPI_Bcast( i_x.data(), i_x.size(), MPI_INT, root_rank, SML_COMM_WORLD);
262  if(i_surf_separatrices.size()>0) MPI_Bcast( i_surf_separatrices.data(), i_surf_separatrices.size(), MPI_INT, root_rank, SML_COMM_WORLD);
263  MPI_Bcast( non_aligned_vert.data(), non_aligned_vert.size(), MPI_INT, root_rank, SML_COMM_WORLD);
264  MPI_Bcast( non_aligned_nsurf.data(), non_aligned_nsurf.size(), MPI_INT, root_rank, SML_COMM_WORLD);
265  MPI_Bcast( non_aligned_surf_idx.data(), non_aligned_surf_idx.size(), MPI_INT, root_rank, SML_COMM_WORLD);
266  MPI_Bcast( nnodes_on_surface.data(), nnodes_on_surface.size(), MPI_INT, root_rank, SML_COMM_WORLD);
267  MPI_Bcast( surf_idx.data(), surf_idx.size(), MPI_INT, root_rank, SML_COMM_WORLD);
268 #endif
269  }
270 
271  // Options: hex vs circular
275  };
276 
286  void construct_hex_grid(int nshells, double raxis, double zaxis, double rscale, double zscale) {
287  // First, build construction of what .node and .ele files would look like
288 
289  // Build hexagonal sample grid of concentric shells
290  // Minimum 1 shell, which produces 7 nodes and 6 triangles
291 
292  // Grid construction uses this enum for directions
293  enum { NE = 0, N, NW, SW, S, SE };
294 
295  // Node indices of each surface
296  int max_nodes_on_surface = 6*nshells;
297  surf_idx = View<int**, CLayout, HostType>(NoInit("surf_idx"),nsurfaces, max_nodes_on_surface);
298 
299  // Work from center outward
300  x(0).r = raxis;
301  x(0).z = zaxis;
302  rgn(0) = 0; // not wall
303  nnodes_on_surface(0) = 1;
304  surf_idx(0, 0) = 1; // Since surf_idx is 1-indexed
305  int inode = 1;
306  int itri = 0;
307  for(int ish = 1; ish<=nshells; ish++){
308  int nodes_in_shell = 6*ish;
309  int first_inner_shell_node = ish==1 ? 0 : inode - 6*(ish-1);
310  int first_node = inode;
311  int inner_shell_node = first_inner_shell_node;
312  // Starting position:
313  double r_track = raxis + rscale*ish;
314  double z_track = zaxis;
315  for (int inodesh = 0; inodesh < nodes_in_shell; inodesh++){
316  x(inode).r = r_track;
317  x(inode).z = z_track;
318  rgn(inode) = (ish==nshells? OnWall : Inside); // last shell is wall
319 
320  surf_idx(ish, inodesh) = inode+1; // Since surf_idx is 1-indexed
321 
322  if (inodesh>0){
323  nodes(itri).vertex[0] = inner_shell_node;
324  nodes(itri).vertex[1] = inode-1;
325  nodes(itri).vertex[2] = inode;
326  itri++;
327  if (inodesh%ish>0){ // corners only add one triangle
328  nodes(itri).vertex[0] = inner_shell_node;
329  nodes(itri).vertex[1] = (inodesh == nodes_in_shell - 1) ? first_inner_shell_node : (inner_shell_node+1);
330  nodes(itri).vertex[2] = inode;
331  inner_shell_node++;
332  itri++;
333  }
334  }
335 
336  // Determine position of next node
337  int hex_region = inodesh/ish;
338  double r_move, z_move;
339  if(hex_region == NE)
340  {r_move = -0.5; z_move = 1.0;}
341  else if (hex_region == N)
342  {r_move = -1.0; z_move = 0.0;}
343  else if (hex_region == NW)
344  {r_move = -0.5; z_move = -1.0;}
345  else if (hex_region == SW)
346  {r_move = 0.5; z_move = -1.0;}
347  else if (hex_region == S)
348  {r_move = 1.0; z_move = 0.0;}
349  else if (hex_region == SE)
350  {r_move = 0.5; z_move = 1.0;}
351  r_track += rscale*r_move;
352  z_track += zscale*z_move;
353  inode++;
354  }
355 
356  // Add last element to close ring
357  nodes(itri).vertex[0] = first_inner_shell_node;
358  nodes(itri).vertex[1] = inode-1;
359  nodes(itri).vertex[2] = first_node;
360  itri++;
361 
362  nnodes_on_surface(ish) = nodes_in_shell;
363  }
364 
365  // Shift node list to 1-indexing (input files are 1-indexed)
366  for(int i = 0; i<ntriangles; i++){
367  nodes(i).vertex[0]++;
368  nodes(i).vertex[1]++;
369  nodes(i).vertex[2]++;
370  }
371 
372  // Field geometry
373  int n_xpts = 0; // Assume no X-points for now
374  int n_separatrices = 0; // Assume no separatrices for now
375  i_x = View<int*, HostType>("i_x", n_xpts);
376  i_surf_separatrices = View<int*, HostType>(NoInit("i_surf_separatrices"), n_separatrices);
377 
378  // Assume no non-aligned surfaces for now
379  num_non_aligned = 0;
380  non_aligned_nsurf = View<int*, HostType>(NoInit("non_aligned_nsurf"), num_non_aligned);
381  non_aligned_surf_idx = View<int**, CLayout, HostType>(NoInit("non_aligned_surf_idx"), num_non_aligned, 4);
382  }
383 
390  return sqrt((a.r-b.r)*(a.r-b.r) + (a.z-b.z)*(a.z-b.z));
391  }
392 
402  void construct_circular_grid(int nshells, double raxis, double zaxis, double rscale, double zscale) {
403  // First, build construction of what .node and .ele files would look like
404 
405  // Build circular sample grid of concentric shells
406  // Minimum 1 shell, which produces 7 nodes and 6 triangles
407 
408  // Work from center outward
409  x(0).r = raxis;
410  x(0).z = zaxis;
411  rgn(0) = 0; // not wall
412 
413  // Node indices of each surface
414  int max_nodes_on_surface = 6*nshells;
415  surf_idx = View<int**, CLayout, HostType>(NoInit("surf_idx"),nsurfaces, max_nodes_on_surface);
416 
417  nnodes_on_surface(0) = 1;
418  surf_idx(0, 0) = 1; // Since surf_idx is 1-indexed
419  int inode = 1;
420  int itri = 0;
421  for(int ish = 1; ish<=nshells; ish++){
422  int nodes_in_shell = 6*ish;
423  int nodes_in_inner_shell = ish==1 ? 1 : 6*(ish-1);
424  int first_inner_shell_node = ish==1 ? 0 : inode - 6*(ish-1);
425  int first_node = inode;
426  int inner_shell_node = first_inner_shell_node;
427  // Starting position:
428  for (int inodesh = 0; inodesh <= nodes_in_shell; inodesh++){
429  // Create node
430  if(inodesh<nodes_in_shell){
431  double angle = TWOPI*double(inodesh)/double(nodes_in_shell);
432  x(inode).r = raxis + ish*rscale*cos(angle);
433  x(inode).z = zaxis + ish*zscale*sin(angle);
434  rgn(inode) = (ish==nshells? OnWall : Inside); // last shell is wall
435 
436  surf_idx(ish, inodesh) = inode+1; // Since surf_idx is 1-indexed
437  }
438 
439  // Create element(s)
440  if (inodesh>0){
441  int outer_node = (inodesh < nodes_in_shell) ? inode : first_node;
442  int prev_outer_node = inode-1;
443 
444  // First, check cross-distances
445  int next_inner_shell_node = inner_shell_node+1;
446  if(next_inner_shell_node-first_inner_shell_node==nodes_in_inner_shell) next_inner_shell_node = first_inner_shell_node;
447  double distance1 = node_to_node_dist(x(outer_node), x(inner_shell_node));
448  double distance2 = node_to_node_dist(x(prev_outer_node), x(next_inner_shell_node));
449 
450  // The next node is far away enough that we need to add two elements
451  if(distance1>distance2 && next_inner_shell_node!=inner_shell_node){
452  if(itri>=nodes.size()) {printf("\nMore elements than expected in circular grid generation. Dynamically resize?"); exit(1); }
453  nodes(itri).vertex[0] = inner_shell_node;
454  nodes(itri).vertex[1] = next_inner_shell_node;
455  nodes(itri).vertex[2] = prev_outer_node;
456  inner_shell_node = next_inner_shell_node;
457  itri++;
458  }
459 
460  if(itri>=nodes.size()) {printf("\nMore elements than expected in circular grid generation. Dynamically resize?"); exit(1); }
461  nodes(itri).vertex[0] = inner_shell_node;
462  nodes(itri).vertex[1] = prev_outer_node;
463  nodes(itri).vertex[2] = outer_node;
464  itri++;
465  }
466 
467  if(inodesh<nodes_in_shell){
468  inode++;
469  }
470  }
471 
472  nnodes_on_surface(ish) = nodes_in_shell;
473  }
474 
475  if(itri!=nodes.size()) {printf("\nFewer elements than expected in circular grid generation. Dynamically resize?"); exit(1); }
476 
477  // Shift node list to 1-indexing (input files are 1-indexed)
478  for(int i = 0; i<ntriangles; i++){
479  nodes(i).vertex[0]++;
480  nodes(i).vertex[1]++;
481  nodes(i).vertex[2]++;
482  }
483 
484  // Field geometry
485  int n_xpts = 0; // Assume no X-points for now
486  int n_separatrices = 0; // Assume no separatrices for now
487  i_x = View<int*, HostType>("i_x", n_xpts);
488  i_surf_separatrices = View<int*, HostType>(NoInit("i_surf_separatrices"), n_separatrices);
489 
491  nsurfaces2 = 0;
492  nsurfaces3a = 0;
493  nsurfaces3b = 0;
494 
495  // Assume no non-aligned surfaces for now
496  num_non_aligned = 0;
497  non_aligned_nsurf = View<int*, HostType>(NoInit("non_aligned_nsurf"), num_non_aligned);
498  non_aligned_surf_idx = View<int**, CLayout, HostType>(NoInit("non_aligned_surf_idx"), num_non_aligned, 4);
499  non_aligned_vert = View<int*, HostType>(NoInit("non_aligned_vert"), num_non_aligned);
500  }
501 
502  // Returns the number of shells needed to create a grid with at least nnodes nodes.
503  static int nshells_from_nnodes(int nnodes){
504  // nnodes = 1 + 3*(nshells+1)*nshells
505  // 0 = 3s^2 + 3s + 1 - n
506  // s = (-3 +/- sqrt(9 + 4*3*(n-1)))/6
507  // s = sqrt(9 + 12*(n-1))/6 - 0.5
508  if(nnodes<1) return 0;
509  else return std::ceil(sqrt(12*nnodes - 3)/6 - 0.5);
510  }
511 
522  GridFiles(int nshells, double raxis, double zaxis, double rscale, double zscale, GridCreateOpt grid_opt=Hexagonal)
523  : nnodes(1 + 3*(nshells+1)*nshells),
524  x(NoInit("x"),nnodes),
525  rgn(NoInit("rgn"),nnodes),
526  ntriangles(6*nshells*nshells),
527  nodes(NoInit("nodes"),ntriangles),
528  nsurfaces(nshells+1),
529  nnodes_on_surface(NoInit("nnodes_on_surface"),nsurfaces),
531  {
532  if(grid_opt==Hexagonal){
533  construct_hex_grid(nshells, raxis, zaxis, rscale, zscale);
534  }else{
535  construct_circular_grid(nshells, raxis, zaxis, rscale, zscale);
536  }
537  }
538 
547  : nnodes(1 + 3*(nshells+1)*nshells),
548  x(NoInit("x"),nnodes),
549  rgn(NoInit("rgn"),nnodes),
550  ntriangles(6*nshells*nshells),
551  nodes(NoInit("nodes"),ntriangles),
552  nsurfaces(nshells+1),
553  nnodes_on_surface(NoInit("nnodes_on_surface"),nsurfaces),
555  {
556  // Magnetic axis will be the center of the grid
557  double raxis = magnetic_field.equil.axis_r;
558  double zaxis = magnetic_field.equil.axis_z;
559 
560  // Since the axis is not necessarily centered in the field map, take the shorter
561  // distance to the edge to use as the dimensions of the grid
562  double width = 2*std::min(magnetic_field.bounds.max_r - magnetic_field.equil.axis_r,
563  magnetic_field.equil.axis_r - magnetic_field.bounds.min_r);
564 
565  double height = 2*std::min(magnetic_field.bounds.max_z - magnetic_field.equil.axis_z,
566  magnetic_field.equil.axis_z - magnetic_field.bounds.min_z);
567 
568  // Make sure the grid is within the psi bounds
569  width *= 0.99;
570  height *= 0.99;
571 
572  if(grid_opt==Hexagonal){
573  double triangle_width = width/(2*nshells);
574  double triangle_height = height/(2*nshells);
575  construct_hex_grid(nshells, raxis, zaxis, triangle_width, triangle_height);
576  }else{
577  // Circular grid; take the smaller of the width or height to make sure it fits
578  double scale = std::min(width,height)/(2*nshells);
579  construct_circular_grid(nshells, raxis, zaxis, scale, scale);
580  }
581  }
582 
584  return i_x.size();
585  }
586 
587  View<Equil::XPoint*, HostType> get_xpts(){
588  View<Equil::XPoint*, HostType> xpts(NoInit("xpts"), i_x.size());
589  for(int i = 0; i<i_x.size(); i++){
590  xpts(i).r = x(i_x(i)).r;
591  xpts(i).z = x(i_x(i)).z;
592  }
593  return xpts;
594  };
595 
596  /* Returns GridFiles object using the file-reader constructor. This may be updated so the input file
597  * can specify an analytic grid.
598  * */
600  // Read file names and path from input file
601  nlr.use_namelist("sml_param");
602 
603  if(nlr.get<bool>("sml_create_analytic_grid", false)==true){
604  // Keep analytic grid options simple: number of surfaces, and axis_r
605  int nshells = 30;
606  double axis_r = 1.5;
607  if(nlr.namelist_present("analytic_grid_param")){
608  nlr.use_namelist("analytic_grid_param");
609 
610  nshells = nlr.get<int>("nsurfaces", 30);
611  axis_r = nlr.get<double>("axis_r", 1.5);
612  }
613 
614  if(is_rank_zero()) printf("\nCreating a circular analytic grid with %d surfaces and %d nodes.\n", nshells, 1 + 3*(nshells+1)*nshells);
615 
616  double axis_z = 0.0;
617 
618  // Leave a margin so that all vertices are inside equilibrium
619  double scale = 0.999/nshells; // Height of each triangle
620  double rscale = scale; // rscale = zscale for circular
621  double zscale = scale;
622  return GridFiles(nshells, axis_r, axis_z, rscale, zscale, GridCreateOpt::Circular);
623  }
624 
625  std::string node_file = nlr.get<std::string>("sml_node_file","example_file.node");
626  std::string ele_file = nlr.get<std::string>("sml_ele_file","example_file.ele");
627  std::string surf_file = nlr.get<std::string>("sml_surf_file","example_file.flx.aif");
628  std::string input_file_dir = nlr.get<std::string>("sml_input_file_dir","./");
629 
630  // Add path to file names
631  node_file = input_file_dir + "/" + node_file;
632  ele_file = input_file_dir + "/" + ele_file;
633  surf_file = input_file_dir + "/" + surf_file;
634 
635  // Read in grid files and return result
636  return GridFiles(node_file, ele_file, surf_file);
637  }
638 
639 };
640 
641 #endif
int nsurfaces3b
Definition: grid_files.hpp:31
int nsurfaces2
Definition: grid_files.hpp:31
void convert_to_x_rgn_format(const View< double **, CLayout, HostType > &tmp_nodes, View< RZPair *, HostType > &x, View< int *, HostType > &rgn)
Definition: grid_files.hpp:41
bool is_rank_zero()
Definition: globals.hpp:27
View< int *, HostType > non_aligned_vert
Definition: grid_files.hpp:36
MPI_Comm SML_COMM_WORLD
Definition: my_mpi.cpp:4
int num_non_aligned
Definition: grid_files.hpp:35
T get(const string &param, const T default_val, int val_ind=0)
Definition: NamelistReader.hpp:373
View< Vertex *, HostType > nodes
Definition: grid_files.hpp:25
void set_fortran_flx_aif_format(int is_flx_aif_format_int)
int nsurfaces
Definition: grid_files.hpp:30
Definition: NamelistReader.hpp:193
int nsurfaces1
Definition: grid_files.hpp:31
Definition: magnetic_field.hpp:12
View< int **, CLayout, HostType > surf_idx
Definition: grid_files.hpp:39
RZBounds bounds
Simulation boundary.
Definition: magnetic_field.hpp:23
Definition: grid_files.hpp:17
Equilibrium equil
The object containing information about the magnetic equilibrium.
Definition: magnetic_field.hpp:32
View< int *, HostType > non_aligned_nsurf
Definition: grid_files.hpp:37
static GridFiles grid_files_from_namelist(NLReader::NamelistReader &nlr)
Definition: grid_files.hpp:599
void read(T *var, int n)
Definition: file_reader.hpp:28
View< int *, HostType > nnodes_on_surface
Definition: grid_files.hpp:33
void construct_circular_grid(int nshells, double raxis, double zaxis, double rscale, double zscale)
Definition: grid_files.hpp:402
View< int *, HostType > rgn
Definition: grid_files.hpp:21
View< int *, HostType > i_x
Definition: grid_files.hpp:28
GridFiles(int nshells, double raxis, double zaxis, double rscale, double zscale, GridCreateOpt grid_opt=Hexagonal)
Definition: grid_files.hpp:522
double z
Definition: grid_structs.hpp:30
GridFiles()
Definition: grid_files.hpp:60
double axis_z
z coordinate of axis
Definition: equil.hpp:95
Definition: grid_structs.hpp:28
int nsurfaces3a
Definition: grid_files.hpp:31
View< RZPair *, HostType > x
Definition: grid_files.hpp:20
Definition: grid_files.hpp:274
GridFiles(MagneticField< HostType > &magnetic_field, int nshells, GridCreateOpt grid_opt=Hexagonal)
Definition: grid_files.hpp:546
View< Equil::XPoint *, HostType > get_xpts()
Definition: grid_files.hpp:587
View< int *, HostType > i_surf_separatrices
Definition: grid_files.hpp:32
int ntriangles
Definition: grid_files.hpp:24
void use_namelist(const string &namelist)
Definition: NamelistReader.hpp:355
double max_z
Definition: rz_bounds.hpp:8
Definition: grid_files.hpp:273
void convert_to_node_format(const View< int **, CLayout, HostType > &tmp_ele, View< Vertex *, HostType > &nodes)
Definition: grid_files.hpp:50
double axis_r
r coordinate of axis
Definition: equil.hpp:94
double min_r
Definition: rz_bounds.hpp:5
Definition: col_grid.hpp:34
Definition: file_reader.hpp:6
bool namelist_present(const string &namelist)
Definition: NamelistReader.hpp:351
void construct_hex_grid(int nshells, double raxis, double zaxis, double rscale, double zscale)
Definition: grid_files.hpp:286
double r
Definition: grid_structs.hpp:29
void exit_XGC(std::string msg)
Definition: globals.hpp:37
Definition: magnetic_field.F90:1
int sum_nnodes_on_surfaces
Definition: grid_files.hpp:34
static int nshells_from_nnodes(int nnodes)
Definition: grid_files.hpp:503
GridCreateOpt
Definition: grid_files.hpp:272
GridFiles(const std::string &node_file, const std::string &element_file, const std::string &flx_aif_file)
Definition: grid_files.hpp:63
void parallel_for(const std::string name, int n_ptl, Function func, Option option, HostAoSoA aosoa_h, DeviceAoSoA aosoa_d)
Definition: streamed_parallel_for.hpp:252
double min_z
Definition: rz_bounds.hpp:7
int nnodes
Definition: grid_files.hpp:19
void skip(int n)
Definition: file_reader.hpp:45
View< int **, CLayout, HostType > non_aligned_surf_idx
Definition: grid_files.hpp:38
Kokkos::ViewAllocateWithoutInitializing NoInit
Definition: space_settings.hpp:68
constexpr double TWOPI
Definition: constants.hpp:9
bool is_old_flx_format(const std::string &filename)
Definition: check_aif_flx_format.hpp:13
double node_to_node_dist(RZPair a, RZPair b)
Definition: grid_files.hpp:389
void open(const std::string &input)
Definition: file_reader.hpp:13
int n_specified_xpts()
Definition: grid_files.hpp:583
double max_r
Definition: rz_bounds.hpp:6