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 // Sets the maximum number of planes that is read for stellarator simulations, e.g. 4 means up to 9999 planes.
5 // (This will probably be removed later with a different mesh format.)
6 #define STELLARATOR_PLANE_LIMIT 4
7 
8 #ifdef USE_MPI
9 #include "my_mpi.hpp"
10 #endif
11 #include "file_reader.hpp"
12 #include "magnetic_field.hpp"
13 #include "grid_structs.hpp"
14 #include "space_settings.hpp"
15 #include "NamelistReader.hpp"
16 #include "check_aif_flx_format.hpp"
17 #include "send_recv_toroidal.hpp"
18 
19 extern "C" void set_fortran_flx_aif_format(int is_flx_aif_format_int);
20 
21 struct PlaneFiles{
22  // These are the regions designated in the .node file
23  enum FileRegion{
24  Inside=0,
26  };
27 
28  // Contents of .node
29  int nnodes;
30  View<RZPair*, HostType> x;
31  View<int*, HostType> rgn;
32 
33  // Contents of .ele
35  View<Vertex*, HostType> nodes;
36 
37  // Contents of .flx.aif
38  View<int*, HostType> i_x;
39 
40  int nsurfaces;
42  View<int*, HostType> i_surf_separatrices;
43  View<int*, HostType> nnodes_on_surface;
46  View<int*, HostType> non_aligned_vert;
47  View<int*, HostType> non_aligned_nsurf;
48  View<int**, CLayout, HostType> non_aligned_surf_idx;
49  View<int**, CLayout, HostType> surf_idx;
50 
51  View<double*, HostType> psiN;
52 
53  void convert_to_x_rgn_format(const View<double**, CLayout, HostType>& tmp_nodes, View<RZPair*, HostType>& x, View<int*, HostType>& rgn, bool is_stellarator){
54  Kokkos::parallel_for("transpose_x_and_rgn", Kokkos::RangePolicy<HostExSpace>(0,tmp_nodes.extent(0)), KOKKOS_LAMBDA( const int i ){
55  // tmp_nodes[i,0]; <-- first column is ignored
56  x(i).r = tmp_nodes(i,1);
57  x(i).z = tmp_nodes(i,2);
58  rgn(i) = is_stellarator ? tmp_nodes(i,5) : tmp_nodes(i,3); // For stellarator .node files rgn is column 6, for tokamaks column 4
59  });
60  }
61 
62  void convert_to_node_format(const View<int**, CLayout, HostType>& tmp_ele, View<Vertex*, HostType>& nodes){
63  Kokkos::parallel_for("transpose_x_and_rgn", Kokkos::RangePolicy<HostExSpace>(0,tmp_ele.extent(0)), KOKKOS_LAMBDA( const int i ){
64  // tmp_ele[i,0]; <-- first column is ignored
65  nodes(i).vertex[0] = tmp_ele(i,1);
66  nodes(i).vertex[1] = tmp_ele(i,2);
67  nodes(i).vertex[2] = tmp_ele(i,3);
68  });
69  }
70 
71  View<int*, HostType> pack_scalars(int max_nodes_on_surface, int n_xpts, int n_separatrices) const{
72  int n_scalars = 11;
73  View<int*, HostType> scalars_s("scalars", n_scalars);
74  scalars_s(0) = nnodes;
75  scalars_s(1) = ntriangles;
76  scalars_s(2) = nsurfaces;
77  scalars_s(3) = max_nodes_on_surface;
78  scalars_s(4) = num_non_aligned;
79  scalars_s(5) = n_xpts;
80  scalars_s(6) = n_separatrices;
81  scalars_s(7) = nsurfaces1;
82  scalars_s(8) = nsurfaces2;
83  scalars_s(9) = nsurfaces3a;
84  scalars_s(10) = nsurfaces3b;
85  return scalars_s;
86  }
87 
88  void unpack_scalars(const View<int*, HostType>& scalars, int& max_nodes_on_surface, int& n_xpts, int& n_separatrices){
89  nnodes = scalars(0);
90  ntriangles = scalars(1);
91  nsurfaces = scalars(2);
92  max_nodes_on_surface = scalars(3);
93  num_non_aligned = scalars(4);
94  n_xpts = scalars(5);
95  n_separatrices = scalars(6);
96  nsurfaces1 = scalars(7);
97  nsurfaces2 = scalars(8);
98  nsurfaces3a = scalars(9);
99  nsurfaces3b = scalars(10);
100  }
101 
102  // Default constructor, to use to deallocate the file arrays after use
104 
105 // Should remove this USE_MPI at some point
106 #ifdef USE_MPI
107 
108  // Copy data from a PlaneFiles object on a neighboring rank
109  PlaneFiles(const MyMPI& my_mpi, const PlaneFiles& plane_files){
110  /*** Broadcast grid contents to other MPI ranks ***/
111 
112  // First broadcast scalars so we know the sizes of the arrays
113  int max_nodes_on_surface = plane_files.surf_idx.extent(1);
114  int n_xpts = plane_files.i_x.size();
115  int n_separatrices = plane_files.i_surf_separatrices.size();
116  auto scalars_s = plane_files.pack_scalars(max_nodes_on_surface, n_xpts, n_separatrices);
117  View<int*, HostType> scalars("scalars", scalars_s.layout());
118  send_right(my_mpi, scalars_s, scalars);
119  unpack_scalars(scalars, max_nodes_on_surface, n_xpts, n_separatrices);
120 
122  // node and ele views
123  x = View<RZPair*, HostType>(NoInit("x"),nnodes);
124  rgn = View<int*, HostType>(NoInit("rgn"),nnodes);
125  nodes = View<Vertex*, HostType>("nodes",ntriangles);
126 
127  // flx.aif views
128  i_x = View<int*, HostType>("i_x", n_xpts);
129  i_surf_separatrices = View<int*, HostType>(NoInit("i_surf_separatrices"), n_separatrices);
130  non_aligned_vert = View<int*, HostType>(NoInit("non_aligned_vert"), num_non_aligned);
131  non_aligned_nsurf = View<int*, HostType>(NoInit("non_aligned_nsurf"), num_non_aligned);
132  non_aligned_surf_idx = View<int**, CLayout, HostType>(NoInit("non_aligned_surf_idx"), num_non_aligned, 4);
133  nnodes_on_surface = View<int*, HostType>(NoInit("nnodes_on_surface"), nsurfaces);
134  surf_idx = View<int**, CLayout, HostType>(NoInit("surf_idx"),nsurfaces, max_nodes_on_surface);
135 
136  // Send grid arrays
137  send_right(my_mpi, plane_files.rgn, rgn);
138  // Cast these as double and int respectively for MPI simplicity (or could create MPI_Types)
139  View<double*,CLayout, HostType, Kokkos::MemoryTraits<Kokkos::Unmanaged>> x_dbl((double*)(x.data()), x.size()*2);
140  View<double*,CLayout, HostType, Kokkos::MemoryTraits<Kokkos::Unmanaged>> src_x_dbl((double*)(plane_files.x.data()), plane_files.x.size()*2);
141  send_right(my_mpi, src_x_dbl, x_dbl);
142  View<int*,CLayout, HostType, Kokkos::MemoryTraits<Kokkos::Unmanaged>> nodes_dbl((int*)(nodes.data()), nodes.size()*3);
143  View<int*,CLayout, HostType, Kokkos::MemoryTraits<Kokkos::Unmanaged>> src_nodes_dbl((int*)(plane_files.nodes.data()), plane_files.nodes.size()*3);
144  send_right(my_mpi, src_nodes_dbl, nodes_dbl);
145 
146  // These might be size 0 on one plane but not the other
147  send_right(my_mpi, plane_files.i_x, i_x);
148  send_right(my_mpi, plane_files.i_surf_separatrices, i_surf_separatrices);
149  send_right(my_mpi, plane_files.non_aligned_vert, non_aligned_vert);
150  send_right(my_mpi, plane_files.non_aligned_nsurf, non_aligned_nsurf);
151  send_right(my_mpi, plane_files.non_aligned_surf_idx, non_aligned_surf_idx);
152  send_right(my_mpi, plane_files.nnodes_on_surface, nnodes_on_surface);
153  send_right(my_mpi, plane_files.surf_idx, surf_idx);
154  }
155 
156  // Read the files
157  PlaneFiles(const std::string& node_file, const std::string& element_file, const std::string& flx_aif_file, bool is_stellarator, const MPI_Comm& comm){
158  // MPI setup
159  int my_rank;
160 #ifdef USE_MPI
161  MPI_Comm_rank(comm, &my_rank);
162 #else
163  my_rank = 0;
164 #endif
165  int ROOT_RANK = 0;
166 
167  int max_nodes_on_surface;
168  int n_xpts;
169  int n_separatrices;
170  // Only one MPI rank reads the files, then broadcasts to the others
171  if(my_rank==ROOT_RANK){
172 
173  FileReader file_reader;
174 
176  file_reader.open(node_file);
177  file_reader.read(&nnodes, 1); // Get nnodes from file
178  file_reader.skip(3); // Skip 3 unused values in file
179 
180  // Use dummy array since reformatting is needed
181  View<double**, CLayout, HostType> tmp_nodes(NoInit("tmp_nodes"),nnodes,4); // 4 entries per node: index, r, z, region
182  if(is_stellarator){
183  // AM: I think (psi/psi_a)^(1/2) and poloidal angle are only used for debugging in Fortran XGC-S so do not store them.
184  Kokkos::resize(tmp_nodes,nnodes,6); // 6 entries per node: index, r, z, (psi/psi_a)^(1/2), poloidal angle, region
185  }
186  file_reader.read(tmp_nodes.data(), tmp_nodes.size());
187 
188  // Reformat to x and rgn
189  x = View<RZPair*, HostType>(NoInit("x"),nnodes);
190  rgn = View<int*, HostType>(NoInit("rgn"),nnodes);
191  convert_to_x_rgn_format(tmp_nodes, x, rgn, is_stellarator);
192 
194  file_reader.open(element_file);
195  file_reader.read(&ntriangles, 1); // Get ntriangles from file
196  file_reader.skip(2); // Skip 2 unused values in file
197 
198  // Use dummy array since reformatting is needed
199  View<int**, CLayout, HostType> tmp_ele(NoInit("tmp_ele"),ntriangles,4); // 4 entries per node: index, 3 vertices of the triangle
200  file_reader.read(tmp_ele.data(), tmp_ele.size());
201 
202  // Reformat since first column is ignored
203  nodes = View<Vertex*, HostType>("nodes",ntriangles);
204  convert_to_node_format(tmp_ele, nodes);
205 
207  if(nnodes==1) exit_XGC("\nGrid must have more than 1 nnode!\n"); // This is assumed in the flx format check
208  bool old_flx_aif_format = is_stellarator ? false : is_old_flx_format(flx_aif_file);
209  set_fortran_flx_aif_format(old_flx_aif_format ? 1 : 0);
210 
211  file_reader.open(flx_aif_file);
212 
213  if(!is_stellarator){
214  if(old_flx_aif_format){
215  file_reader.read(&nsurfaces, 1); // Get number of surfaces
216  n_xpts = 0;
217  n_separatrices = 0;
218  i_x = View<int*, HostType>("i_x",0); // No xpoints specified
219  i_surf_separatrices = View<int*, HostType>(NoInit("i_surf_separatrices"), 0); // No separatrices specified
220 
221  // These don't get read or used by the old format
222  nsurfaces1 = 0;
223  nsurfaces2 = 0;
224  nsurfaces3a = 0;
225  nsurfaces3b = 0;
226  }else{
227  file_reader.read(&n_xpts, 1); // Get number of x-points
228  i_x = View<int*, HostType>("i_x", n_xpts);
229 
230  // Current flx.aif files list 2 x-point indices, regardless of n_xpts
231  int i_xpt1, i_xpt2;
232  file_reader.read(&i_xpt1, 1); // Get index of 1st x-point
233  file_reader.read(&i_xpt2, 1); // Get index of 2nd x-point
234 
235  // Check for consistency
236  if(n_xpts==1 && i_xpt2 != -1) exit_XGC("Error in flx.aif file: n_xpts=1, but i_xpt2 is specified.\n");
237  if(n_xpts==2 && i_xpt2 <= 0) exit_XGC("Error in flx.aif file: n_xpts=2, but i_xpt2 is not specified.\n");
238 
239  if(n_xpts>=1) i_x(0) = i_xpt1 - 1; // switch to zero-index
240  if(n_xpts>=2) i_x(1) = i_xpt2 - 1;
241 
242  // Next are the nsurface variables
243  file_reader.read(&nsurfaces1, 1); // Get number of surfaces
244  file_reader.read(&nsurfaces2, 1); // Get number of surfaces
245  file_reader.read(&nsurfaces3a, 1); // Get number of surfaces
246  file_reader.read(&nsurfaces3b, 1); // Get number of surfaces
247 
248  // nsurfaces is actually the sum of these 4 inputs
249  nsurfaces = nsurfaces1 + nsurfaces2 + nsurfaces3a + nsurfaces3b;
250 
251  // Get surface index of separatrices
252  // Hard-coded to 2, but generalize the read just in case
253  n_separatrices = 2;
254  i_surf_separatrices = View<int*, HostType>(NoInit("i_surf_separatrices"), n_separatrices);
255  file_reader.read(&i_surf_separatrices(0), 1);
256  file_reader.read(&i_surf_separatrices(1), 1);
257  }
258  }else{ // Stellarator .grid file
259  n_xpts = 0;
260  n_separatrices = 0;
261  i_x = View<int*, HostType>("i_x",0); // No xpoints specified
262  i_surf_separatrices = View<int*, HostType>(NoInit("i_surf_separatrices"), 0); // No separatrices specified
263 
264  file_reader.read(&nsurfaces1, 1); // Get number of surfaces in the core
265  nsurfaces2 = 0; // No surfaces in the SOL
266  nsurfaces3a = 0; // No surfaces in the lower private region
267  nsurfaces3b = 0; // No surfaces in the upper private region
268 
269  nsurfaces = nsurfaces1 + nsurfaces2 + nsurfaces3a + nsurfaces3b;
270  }
271  // Read in nnodes on each surface
272  nnodes_on_surface = View<int*, HostType>(NoInit("nnodes_on_surface"), nsurfaces);
273  if(!is_stellarator){
274  file_reader.read(nnodes_on_surface.data(), nnodes_on_surface.size());
275  }else{ // Stellarator .grid file
276  psiN = View<double*, HostType>(NoInit("normalized_psi_on_surface"), nsurfaces);
277  int dummy_idx;
278  for(int i=0; i<nnodes_on_surface.size(); i++){
279  file_reader.read(&nnodes_on_surface(i), 1);
280  file_reader.read(&psiN(i), 1);
281  if(i<nnodes_on_surface.size() - 1) file_reader.read(&dummy_idx, 1);
282  }
283  }
284 
285  // Loop through surfaces to find max number of nodes in a surface
286  max_nodes_on_surface = 0;
287  for(int i=0; i<nnodes_on_surface.size(); i++) max_nodes_on_surface = std::max(nnodes_on_surface(i), max_nodes_on_surface);
288 
289  // Allocate surf_idx and initialize to 0
290  surf_idx = View<int**, CLayout, HostType>("surf_idx",nsurfaces, max_nodes_on_surface);
291 
292  // Read in values
293  if(!is_stellarator){
294  for(int i=0; i<nnodes_on_surface.size(); i++){
295  file_reader.read(&surf_idx(i,0), nnodes_on_surface(i));
296  }
297  }else{ // The stellarator .grid file has no enumeration of surface indices so we assume they are in order
298  int surf_idx_counter = 0;
299  for(int i=0; i<nnodes_on_surface.size(); i++){
300  for(int j=0; j<nnodes_on_surface(i); j++){
301  surf_idx_counter++;
302  surf_idx(i,j) = surf_idx_counter;
303  }
304  }
305  }
306 
307  int should_be_negative_one;
308  if(!is_stellarator){
309  // There is a "-1" in the file here to mark the end of the section
310  file_reader.read(&should_be_negative_one, 1);
311  if(should_be_negative_one != -1){
312  printf("\nExpected -1 at end of section in surface file, encountered %d instead\n", should_be_negative_one);
313  exit_XGC("Error reading surface file.\n");
314  }
315 
316  file_reader.read(&num_non_aligned, 1);
317  printf(" - reading num_non_aligned: %d",num_non_aligned);
318  }else{ // The stellarator .grid file has no non-aligned vertex data
319  num_non_aligned = 0;
320  }
321 
322  // Read the non-aligned vertex data
323  non_aligned_vert = View<int*, HostType>(NoInit("non_aligned_vert"), num_non_aligned);
324  non_aligned_nsurf = View<int*, HostType>(NoInit("non_aligned_nsurf"), num_non_aligned);
325  // Initialize to 0
326  non_aligned_surf_idx = View<int**, CLayout, HostType>("non_aligned_surf_idx", num_non_aligned, 4);
327  for(int i=0; i<num_non_aligned; i++){
328  file_reader.read(&non_aligned_vert(i), 1);
329  if(old_flx_aif_format){
330  file_reader.read(&non_aligned_nsurf(i), 1);
331  }else{
332  // New format does not contain non_aligned_nsurf
333  non_aligned_nsurf(i)=2;
334  }
335  file_reader.read(&non_aligned_surf_idx(i,0), non_aligned_nsurf(i));
336  }
337 
338  if(!is_stellarator){
339  // Another "-1" at the end of the file
340  file_reader.read(&should_be_negative_one, 1);
341  if(should_be_negative_one != -1){
342  printf("\nExpected -1 at end of the surface file, encountered %d instead\n", should_be_negative_one);
343  exit_XGC("Error reading surface file.\n");
344  }
345  }
346  }
347 
348  /*** Broadcast grid contents to other MPI ranks ***/
349 
350  // First broadcast scalars so we know the sizes of the arrays
351  auto scalars = pack_scalars(max_nodes_on_surface, n_xpts, n_separatrices);
352 
353  // Broadcast file size to other processes
354 #ifdef USE_MPI
355  MPI_Bcast(scalars.data(), scalars.size(), MPI_INT, ROOT_RANK, comm);
356 #endif
357 
358  if(my_rank!=ROOT_RANK){
359  // Set scalars on non-root ranks
360  unpack_scalars(scalars, max_nodes_on_surface, n_xpts, n_separatrices);
361 
363  // node and ele views
364  x = View<RZPair*, HostType>(NoInit("x"),nnodes);
365  rgn = View<int*, HostType>(NoInit("rgn"),nnodes);
366  nodes = View<Vertex*, HostType>("nodes",ntriangles);
367 
368  // flx.aif views
369  i_x = View<int*, HostType>("i_x", n_xpts);
370  i_surf_separatrices = View<int*, HostType>(NoInit("i_surf_separatrices"), n_separatrices);
371  non_aligned_vert = View<int*, HostType>(NoInit("non_aligned_vert"), num_non_aligned);
372  non_aligned_nsurf = View<int*, HostType>(NoInit("non_aligned_nsurf"), num_non_aligned);
373  non_aligned_surf_idx = View<int**, CLayout, HostType>(NoInit("non_aligned_surf_idx"), num_non_aligned, 4);
374  nnodes_on_surface = View<int*, HostType>(NoInit("nnodes_on_surface"), nsurfaces);
375  surf_idx = View<int**, CLayout, HostType>(NoInit("surf_idx"),nsurfaces, max_nodes_on_surface);
376  }
377 
378  // Broadcast grid arrays
379 #ifdef USE_MPI
380  MPI_Bcast( x.data(), x.size()*2, MPI_DOUBLE, ROOT_RANK, comm);
381  MPI_Bcast( rgn.data(), rgn.size(), MPI_INT, ROOT_RANK, comm);
382  MPI_Bcast( nodes.data(), nodes.size()*3, MPI_INT, ROOT_RANK, comm);
383 
384  if(i_x.size()>0) MPI_Bcast( i_x.data(), i_x.size(), MPI_INT, ROOT_RANK, comm);
385  if(i_surf_separatrices.size()>0) MPI_Bcast( i_surf_separatrices.data(), i_surf_separatrices.size(), MPI_INT, ROOT_RANK, comm);
386  if(non_aligned_vert.size()>0) MPI_Bcast( non_aligned_vert.data(), non_aligned_vert.size(), MPI_INT, ROOT_RANK, comm);
387  if(non_aligned_nsurf.size()>0) MPI_Bcast( non_aligned_nsurf.data(), non_aligned_nsurf.size(), MPI_INT, ROOT_RANK, comm);
388  if(non_aligned_surf_idx.size()>0) MPI_Bcast( non_aligned_surf_idx.data(), non_aligned_surf_idx.size(), MPI_INT, ROOT_RANK, comm);
389  MPI_Bcast( nnodes_on_surface.data(), nnodes_on_surface.size(), MPI_INT, ROOT_RANK, comm);
390  MPI_Bcast( surf_idx.data(), surf_idx.size(), MPI_INT, ROOT_RANK, comm);
391 #endif
392  }
393 #endif
394 
395  // Options: hex vs circular
399  };
400 
410  void construct_hex_grid(int nshells, double raxis, double zaxis, double rscale, double zscale) {
411  // First, build construction of what .node and .ele files would look like
412 
413  // Build hexagonal sample grid of concentric shells
414  // Minimum 1 shell, which produces 7 nodes and 6 triangles
415 
416  // Grid construction uses this enum for directions
417  enum { NE = 0, N, NW, SW, S, SE };
418 
419  // Node indices of each surface
420  int max_nodes_on_surface = 6*nshells;
421  surf_idx = View<int**, CLayout, HostType>(NoInit("surf_idx"),nsurfaces, max_nodes_on_surface);
422 
423  // Work from center outward
424  x(0).r = raxis;
425  x(0).z = zaxis;
426  rgn(0) = 0; // not wall
427  nnodes_on_surface(0) = 1;
428  surf_idx(0, 0) = 1; // Since surf_idx is 1-indexed
429  int inode = 1;
430  int itri = 0;
431  for(int ish = 1; ish<=nshells; ish++){
432  int nodes_in_shell = 6*ish;
433  int first_inner_shell_node = ish==1 ? 0 : inode - 6*(ish-1);
434  int first_node = inode;
435  int inner_shell_node = first_inner_shell_node;
436  // Starting position:
437  double r_track = raxis + rscale*ish;
438  double z_track = zaxis;
439  for (int inodesh = 0; inodesh < nodes_in_shell; inodesh++){
440  x(inode).r = r_track;
441  x(inode).z = z_track;
442  rgn(inode) = (ish==nshells? OnWall : Inside); // last shell is wall
443 
444  surf_idx(ish, inodesh) = inode+1; // Since surf_idx is 1-indexed
445 
446  if (inodesh>0){
447  nodes(itri).vertex[0] = inner_shell_node;
448  nodes(itri).vertex[1] = inode-1;
449  nodes(itri).vertex[2] = inode;
450  itri++;
451  if (inodesh%ish>0){ // corners only add one triangle
452  nodes(itri).vertex[0] = inner_shell_node;
453  nodes(itri).vertex[1] = (inodesh == nodes_in_shell - 1) ? first_inner_shell_node : (inner_shell_node+1);
454  nodes(itri).vertex[2] = inode;
455  inner_shell_node++;
456  itri++;
457  }
458  }
459 
460  // Determine position of next node
461  int hex_region = inodesh/ish;
462  double r_move, z_move;
463  if(hex_region == NE)
464  {r_move = -0.5; z_move = 1.0;}
465  else if (hex_region == N)
466  {r_move = -1.0; z_move = 0.0;}
467  else if (hex_region == NW)
468  {r_move = -0.5; z_move = -1.0;}
469  else if (hex_region == SW)
470  {r_move = 0.5; z_move = -1.0;}
471  else if (hex_region == S)
472  {r_move = 1.0; z_move = 0.0;}
473  else if (hex_region == SE)
474  {r_move = 0.5; z_move = 1.0;}
475  r_track += rscale*r_move;
476  z_track += zscale*z_move;
477  inode++;
478  }
479 
480  // Add last element to close ring
481  nodes(itri).vertex[0] = first_inner_shell_node;
482  nodes(itri).vertex[1] = inode-1;
483  nodes(itri).vertex[2] = first_node;
484  itri++;
485 
486  nnodes_on_surface(ish) = nodes_in_shell;
487  }
488 
489  // Shift node list to 1-indexing (input files are 1-indexed)
490  for(int i = 0; i<ntriangles; i++){
491  nodes(i).vertex[0]++;
492  nodes(i).vertex[1]++;
493  nodes(i).vertex[2]++;
494  }
495 
496  // Field geometry
497  int n_xpts = 0; // Assume no X-points for now
498  int n_separatrices = 0; // Assume no separatrices for now
499  i_x = View<int*, HostType>("i_x", n_xpts);
500  i_surf_separatrices = View<int*, HostType>(NoInit("i_surf_separatrices"), n_separatrices);
501 
502  // Assume no non-aligned surfaces for now
503  num_non_aligned = 0;
504  non_aligned_nsurf = View<int*, HostType>(NoInit("non_aligned_nsurf"), num_non_aligned);
505  non_aligned_surf_idx = View<int**, CLayout, HostType>(NoInit("non_aligned_surf_idx"), num_non_aligned, 4);
506  }
507 
514  return sqrt((a.r-b.r)*(a.r-b.r) + (a.z-b.z)*(a.z-b.z));
515  }
516 
526  void construct_circular_grid(int nshells, double raxis, double zaxis, double rscale, double zscale) {
527  // First, build construction of what .node and .ele files would look like
528 
529  // Build circular sample grid of concentric shells
530  // Minimum 1 shell, which produces 7 nodes and 6 triangles
531 
532  // Work from center outward
533  x(0).r = raxis;
534  x(0).z = zaxis;
535  rgn(0) = 0; // not wall
536 
537  // Node indices of each surface
538  int max_nodes_on_surface = 6*nshells;
539  surf_idx = View<int**, CLayout, HostType>(NoInit("surf_idx"),nsurfaces, max_nodes_on_surface);
540 
541  nnodes_on_surface(0) = 1;
542  surf_idx(0, 0) = 1; // Since surf_idx is 1-indexed
543  int inode = 1;
544  int itri = 0;
545  for(int ish = 1; ish<=nshells; ish++){
546  int nodes_in_shell = 6*ish;
547  int nodes_in_inner_shell = ish==1 ? 1 : 6*(ish-1);
548  int first_inner_shell_node = ish==1 ? 0 : inode - 6*(ish-1);
549  int first_node = inode;
550  int inner_shell_node = first_inner_shell_node;
551  // Starting position:
552  for (int inodesh = 0; inodesh <= nodes_in_shell; inodesh++){
553  // Create node
554  if(inodesh<nodes_in_shell){
555  double angle = TWOPI*double(inodesh)/double(nodes_in_shell);
556  x(inode).r = raxis + ish*rscale*cos(angle);
557  x(inode).z = zaxis + ish*zscale*sin(angle);
558  rgn(inode) = (ish==nshells? OnWall : Inside); // last shell is wall
559 
560  surf_idx(ish, inodesh) = inode+1; // Since surf_idx is 1-indexed
561  }
562 
563  // Create element(s)
564  if (inodesh>0){
565  int outer_node = (inodesh < nodes_in_shell) ? inode : first_node;
566  int prev_outer_node = inode-1;
567 
568  // First, check cross-distances
569  int next_inner_shell_node = inner_shell_node+1;
570  if(next_inner_shell_node-first_inner_shell_node==nodes_in_inner_shell) next_inner_shell_node = first_inner_shell_node;
571  double distance1 = node_to_node_dist(x(outer_node), x(inner_shell_node));
572  double distance2 = node_to_node_dist(x(prev_outer_node), x(next_inner_shell_node));
573 
574  // The next node is far away enough that we need to add two elements
575  if(distance1>distance2 && next_inner_shell_node!=inner_shell_node){
576  if(itri>=nodes.size()) {printf("\nMore elements than expected in circular grid generation. Dynamically resize?"); exit(1); }
577  nodes(itri).vertex[0] = inner_shell_node;
578  nodes(itri).vertex[1] = next_inner_shell_node;
579  nodes(itri).vertex[2] = prev_outer_node;
580  inner_shell_node = next_inner_shell_node;
581  itri++;
582  }
583 
584  if(itri>=nodes.size()) {printf("\nMore elements than expected in circular grid generation. Dynamically resize?"); exit(1); }
585  nodes(itri).vertex[0] = inner_shell_node;
586  nodes(itri).vertex[1] = prev_outer_node;
587  nodes(itri).vertex[2] = outer_node;
588  itri++;
589  }
590 
591  if(inodesh<nodes_in_shell){
592  inode++;
593  }
594  }
595 
596  nnodes_on_surface(ish) = nodes_in_shell;
597  }
598 
599  if(itri!=nodes.size()) {printf("\nFewer elements than expected in circular grid generation. Dynamically resize?"); exit(1); }
600 
601  // Shift node list to 1-indexing (input files are 1-indexed)
602  for(int i = 0; i<ntriangles; i++){
603  nodes(i).vertex[0]++;
604  nodes(i).vertex[1]++;
605  nodes(i).vertex[2]++;
606  }
607 
608  // Field geometry
609  int n_xpts = 0; // Assume no X-points for now
610  int n_separatrices = 0; // Assume no separatrices for now
611  i_x = View<int*, HostType>("i_x", n_xpts);
612  i_surf_separatrices = View<int*, HostType>(NoInit("i_surf_separatrices"), n_separatrices);
613 
615  nsurfaces2 = 0;
616  nsurfaces3a = 0;
617  nsurfaces3b = 0;
618 
619  // Assume no non-aligned surfaces for now
620  num_non_aligned = 0;
621  non_aligned_nsurf = View<int*, HostType>(NoInit("non_aligned_nsurf"), num_non_aligned);
622  non_aligned_surf_idx = View<int**, CLayout, HostType>(NoInit("non_aligned_surf_idx"), num_non_aligned, 4);
623  non_aligned_vert = View<int*, HostType>(NoInit("non_aligned_vert"), num_non_aligned);
624  }
625 
626  // Returns the number of shells needed to create a grid with at least nnodes nodes.
627  static int nshells_from_nnodes(int nnodes){
628  // nnodes = 1 + 3*(nshells+1)*nshells
629  // 0 = 3s^2 + 3s + 1 - n
630  // s = (-3 +/- sqrt(9 + 4*3*(n-1)))/6
631  // s = sqrt(9 + 12*(n-1))/6 - 0.5
632  if(nnodes<1) return 0;
633  else return std::ceil(sqrt(12*nnodes - 3)/6 - 0.5);
634  }
635 
646  PlaneFiles(int nshells, double raxis, double zaxis, double rscale, double zscale, GridCreateOpt grid_opt=Hexagonal)
647  : nnodes(1 + 3*(nshells+1)*nshells),
648  x(NoInit("x"),nnodes),
649  rgn(NoInit("rgn"),nnodes),
650  ntriangles(6*nshells*nshells),
651  nodes(NoInit("nodes"),ntriangles),
652  nsurfaces(nshells+1),
653  nnodes_on_surface(NoInit("nnodes_on_surface"),nsurfaces),
655  {
656  if(grid_opt==Hexagonal){
657  construct_hex_grid(nshells, raxis, zaxis, rscale, zscale);
658  }else{
659  construct_circular_grid(nshells, raxis, zaxis, rscale, zscale);
660  }
661  }
662 
671  : nnodes(1 + 3*(nshells+1)*nshells),
672  x(NoInit("x"),nnodes),
673  rgn(NoInit("rgn"),nnodes),
674  ntriangles(6*nshells*nshells),
675  nodes(NoInit("nodes"),ntriangles),
676  nsurfaces(nshells+1),
677  nnodes_on_surface(NoInit("nnodes_on_surface"),nsurfaces),
679  {
680  // Magnetic axis will be the center of the grid
681  double raxis = magnetic_field.equil.axis.r;
682  double zaxis = magnetic_field.equil.axis.z;
683 
684  // Since the axis is not necessarily centered in the field map, take the shorter
685  // distance to the edge to use as the dimensions of the grid
686  double width = 2*std::min(magnetic_field.bounds.max_r - magnetic_field.equil.axis.r,
687  magnetic_field.equil.axis.r - magnetic_field.bounds.min_r);
688 
689  double height = 2*std::min(magnetic_field.bounds.max_z - magnetic_field.equil.axis.z,
690  magnetic_field.equil.axis.z - magnetic_field.bounds.min_z);
691 
692  // Make sure the grid is within the psi bounds
693  width *= 0.99;
694  height *= 0.99;
695 
696  if(grid_opt==Hexagonal){
697  double triangle_width = width/(2*nshells);
698  double triangle_height = height/(2*nshells);
699  construct_hex_grid(nshells, raxis, zaxis, triangle_width, triangle_height);
700  }else{
701  // Circular grid; take the smaller of the width or height to make sure it fits
702  double scale = std::min(width,height)/(2*nshells);
703  construct_circular_grid(nshells, raxis, zaxis, scale, scale);
704  }
705  }
706 
707  int n_specified_xpts() const{
708  return i_x.size();
709  }
710 
711  View<Equil::XPoint*, HostType> get_xpts() const{
712  View<Equil::XPoint*, HostType> xpts(NoInit("xpts"), i_x.size());
713  for(int i = 0; i<i_x.size(); i++){
714  xpts(i).r = x(i_x(i)).r;
715  xpts(i).z = x(i_x(i)).z;
716  }
717  return xpts;
718  }
719 
720  // Return first vertex coordinates. Assumes this is the axis
721  RZPair get_axis() const{
722  return x(0);
723  }
724 
725  // Returns last vertex coordinates. Assumes this is on the LCFS for stellarators
727  return x(x.size()-1);
728  }
729 
730  // Returns normalized psi coordinates of surfaces temporarily needed for the stellarator version
731  View<double*, HostType> get_psiN() const{
732  return psiN;
733  }
734 };
735 
736 struct GridFiles{
740 
741  std::vector<PlaneFiles> plane_files_vec;
742 
743  // Allows array-like access pattern
744  const PlaneFiles& operator [](int i) const {return plane_files_vec[i];}
745 
746  // For test grids
747  GridFiles(int nshells, double raxis, double zaxis, double rscale, double zscale, GridCreateOpt grid_opt=Hexagonal){
748  plane_files_vec.push_back(PlaneFiles(nshells, raxis, zaxis, rscale, zscale, grid_opt));
749  }
750 
752  plane_files_vec.push_back(PlaneFiles(magnetic_field, nshells, grid_opt));
753  }
754 
755 #ifdef USE_MPI
756  /* Returns PlaneFiles object using the file-reader constructor. This may be updated so the input file
757  * can specify an analytic grid.
758  * */
759  GridFiles(NLReader::NamelistReader& nlr, bool is_stellarator, const MyMPI& my_mpi){
760  // Read file names and path from input file
761  nlr.use_namelist("sml_param");
762 
763  if(nlr.get<bool>("sml_create_analytic_grid", false)==true){ // Option to use a generated grid instead of reading a grid from file. The grid is composed of concentric circles of vertices. analytic_grid_param can be used to control basic grid properties.
764  // Keep analytic grid options simple: number of surfaces, and axis_r
765  int nshells = 30;
766  double axis_r = 1.5;
767  if(nlr.namelist_present("analytic_grid_param")){
768  nlr.use_namelist("analytic_grid_param");
769 
770  nshells = nlr.get<int>("nsurfaces", 30); // Number of surfaces in the analytic circular grid.
771  axis_r = nlr.get<double>("axis_r", 1.5);
772  }
773 
774  if(is_rank_zero()) printf("\nCreating a circular analytic grid with %d surfaces and %d nodes.\n", nshells, 1 + 3*(nshells+1)*nshells);
775 
776  double axis_z = 0.0;
777 
778  // Leave a margin so that all vertices are inside equilibrium
779  double scale = 0.999/nshells; // Height of each triangle
780  double rscale = scale; // rscale = zscale for circular
781  double zscale = scale;
782  plane_files_vec.push_back(PlaneFiles(nshells, axis_r, axis_z, rscale, zscale, GridCreateOpt::Circular));
783 
784  if(is_stellarator){
785  // Push back a copy of the grid since stellarator needs two
786  plane_files_vec.push_back(PlaneFiles(nshells, axis_r, axis_z, rscale, zscale, GridCreateOpt::Circular));
787  }
788  }else{
789  nlr.use_namelist("sml_param");
790  std::string node_file = nlr.get<std::string>("sml_node_file","example_file.node"); // The grid file
791  std::string ele_file = nlr.get<std::string>("sml_ele_file","example_file.ele"); // The element file
792  std::string surf_file = nlr.get<std::string>("sml_surf_file","example_file.flx.aif"); // The flux surface file
793  std::string input_file_dir = nlr.get<std::string>("sml_input_file_dir","./"); // Directory containing the grid files
794 
795  if(is_stellarator){
796  surf_file = nlr.get<std::string>("sml_surf_file","example_file.grid"); // The surface file for stellarators
797  // Add path to file names
798  // Remove the extensions ".node" and ".ele" since we have the add plane numbers
799  node_file = input_file_dir + "/equilibrium/" + node_file.substr(0, node_file.size() - 5);
800  ele_file = input_file_dir + "/equilibrium/" + ele_file.substr(0, ele_file.size() - 4);
801  // At the moment the surf_file (.grid) is identical for all stellarator planes
802  surf_file = input_file_dir + "/equilibrium/" + surf_file;
803 
804  int iplane = my_mpi.my_intpl_rank;
805  int ifile = iplane*2 + 1; // Midplane is .1, Edge is .2, next midplane is .3, etc.
806 
807  size_t n_char = STELLARATOR_PLANE_LIMIT;
808  if(std::to_string(ifile + 1).length() > n_char)
809  exit_XGC("Error, plane number " + std::to_string(ifile + 1) + " is larger than the maximum number of allowed planes.\n");
810 
811  std::string ifile_str = std::to_string(ifile);
812  ifile_str = std::string(n_char - std::min(n_char, ifile_str.length()), '0') + ifile_str;
813 
814  std::string midplane_node_file = node_file + "_" + ifile_str + ".node";
815  std::string midplane_ele_file = ele_file + "_" + ifile_str + ".ele";
816  std::string midplane_surf_file = surf_file;
817 
818  // Read in midplane grid files (one rank per plane reads)
819  plane_files_vec.push_back(PlaneFiles(midplane_node_file, midplane_ele_file, midplane_surf_file, is_stellarator, my_mpi.plane_comm));
820 
821  ifile += 1; // shift to edge plane
822  ifile_str = std::to_string(ifile);
823  ifile_str = std::string(n_char - std::min(n_char, ifile_str.length()), '0') + ifile_str;
824  node_file = node_file + "_" + ifile_str + ".node";
825  ele_file = ele_file + "_" + ifile_str + ".ele";
826 
827  // Read in edgeplane grid files (one rank per plane reads)
828  plane_files_vec.push_back(PlaneFiles(node_file, ele_file, surf_file, is_stellarator, my_mpi.plane_comm));
829 
830  // Send right edgeplane to next MPI plane
831  plane_files_vec.push_back(PlaneFiles(my_mpi, plane_files_vec[1]));
832  }else{
833  // Add path to file names
834  node_file = input_file_dir + "/" + node_file;
835  ele_file = input_file_dir + "/" + ele_file;
836  surf_file = input_file_dir + "/" + surf_file;
837 
838  // Read in grid files of single plane
839  plane_files_vec.push_back(PlaneFiles(node_file, ele_file, surf_file, is_stellarator, SML_COMM_WORLD));
840  }
841  }
842  }
843 #endif
844 
845  int n_planes() const{
846  return plane_files_vec.size();
847  }
848 
849  void clear(){
850  plane_files_vec.clear();
851  }
852 
853  // Returns the number of shells needed to create a grid with at least nnodes nodes.
854  static int nshells_from_nnodes(int nnodes){
855  return PlaneFiles::nshells_from_nnodes(nnodes);
856  }
857 
858 };
859 
860 #endif
View< double *, HostType > get_psiN() const
Definition: grid_files.hpp:731
RZPair get_reference_point() const
Definition: grid_files.hpp:726
View< int *, HostType > pack_scalars(int max_nodes_on_surface, int n_xpts, int n_separatrices) const
Definition: grid_files.hpp:71
Definition: grid_files.hpp:398
PlaneFiles(int nshells, double raxis, double zaxis, double rscale, double zscale, GridCreateOpt grid_opt=Hexagonal)
Definition: grid_files.hpp:646
FileRegion
Definition: grid_files.hpp:23
static constexpr GridCreateOpt Hexagonal
Definition: grid_files.hpp:738
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:49
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:627
static constexpr GridCreateOpt Circular
Definition: grid_files.hpp:739
int n_planes() const
Definition: grid_files.hpp:845
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:103
View< int *, HostType > non_aligned_vert
Definition: grid_files.hpp:46
RZPair get_axis() const
Definition: grid_files.hpp:721
Definition: NamelistReader.hpp:193
Definition: magnetic_field.hpp:12
View< int *, HostType > non_aligned_nsurf
Definition: grid_files.hpp:47
int ntriangles
Definition: grid_files.hpp:34
RZBounds bounds
Simulation boundary.
Definition: magnetic_field.hpp:41
Definition: grid_files.hpp:736
Equilibrium equil
The object containing information about the magnetic equilibrium.
Definition: magnetic_field.hpp:44
void read(T *var, int n)
Definition: file_reader.hpp:28
int num_non_aligned
Definition: grid_files.hpp:45
void construct_circular_grid(int nshells, double raxis, double zaxis, double rscale, double zscale)
Definition: grid_files.hpp:526
View< Equil::XPoint *, HostType > get_xpts() const
Definition: grid_files.hpp:711
int nsurfaces3a
Definition: grid_files.hpp:41
View< int *, HostType > nnodes_on_surface
Definition: grid_files.hpp:43
int n_specified_xpts() const
Definition: grid_files.hpp:707
Definition: grid_files.hpp:397
GridFiles(int nshells, double raxis, double zaxis, double rscale, double zscale, GridCreateOpt grid_opt=Hexagonal)
Definition: grid_files.hpp:747
double z
Definition: grid_structs.hpp:30
Definition: grid_structs.hpp:28
int nsurfaces2
Definition: grid_files.hpp:41
std::vector< PlaneFiles > plane_files_vec
Definition: grid_files.hpp:741
Definition: my_mpi.hpp:19
GridFiles(MagneticField< HostType > &magnetic_field, int nshells, GridCreateOpt grid_opt=Hexagonal)
Definition: grid_files.hpp:751
void use_namelist(const string &namelist)
Definition: NamelistReader.hpp:355
View< int *, HostType > i_x
Definition: grid_files.hpp:38
double max_z
Definition: rz_bounds.hpp:8
Definition: grid_files.hpp:25
#define STELLARATOR_PLANE_LIMIT
Definition: grid_files.hpp:6
double min_r
Definition: rz_bounds.hpp:5
Definition: col_grid.hpp:35
int nsurfaces1
Definition: grid_files.hpp:41
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:62
const PlaneFiles & operator[](int i) const
Definition: grid_files.hpp:744
View< double *, HostType > psiN
An array of normalized psi coordinates temporarily needed for the stellarator version.
Definition: grid_files.hpp:51
void unpack_scalars(const View< int *, HostType > &scalars, int &max_nodes_on_surface, int &n_xpts, int &n_separatrices)
Definition: grid_files.hpp:88
GridCreateOpt
Definition: grid_files.hpp:396
View< int **, CLayout, HostType > non_aligned_surf_idx
Definition: grid_files.hpp:48
bool namelist_present(const string &namelist)
Definition: NamelistReader.hpp:351
int nnodes
Definition: grid_files.hpp:29
void construct_hex_grid(int nshells, double raxis, double zaxis, double rscale, double zscale)
Definition: grid_files.hpp:410
int nsurfaces3b
Definition: grid_files.hpp:41
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:670
static int nshells_from_nnodes(int nnodes)
Definition: grid_files.hpp:854
void convert_to_x_rgn_format(const View< double **, CLayout, HostType > &tmp_nodes, View< RZPair *, HostType > &x, View< int *, HostType > &rgn, bool is_stellarator)
Definition: grid_files.hpp:53
Definition: grid_files.hpp:21
double node_to_node_dist(RZPair a, RZPair b)
Definition: grid_files.hpp:513
void clear()
Definition: grid_files.hpp:849
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
RZPair axis
coordinates of axis
Definition: equil.hpp:96
View< int *, HostType > rgn
Definition: grid_files.hpp:31
int nsurfaces
Definition: grid_files.hpp:40
void skip(int n)
Definition: file_reader.hpp:45
View< RZPair *, HostType > x
Definition: grid_files.hpp:30
Kokkos::ViewAllocateWithoutInitializing NoInit
Definition: space_settings.hpp:69
constexpr double TWOPI
Definition: constants.hpp:9
Definition: grid_files.hpp:24
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:42
View< Vertex *, HostType > nodes
Definition: grid_files.hpp:35
int sum_nnodes_on_surfaces
Definition: grid_files.hpp:44
void open(const std::string &input)
Definition: file_reader.hpp:13
double max_r
Definition: rz_bounds.hpp:6