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