Post processing#

General#

class struphy.post_processing.post_processing_tools.ParamsIn(env: EnvironmentOptions | None = None, base_units: BaseUnits | None = None, time_opts: Time | None = None, domain=None, equil=None, grid: TensorProductGrid | None = None, derham_opts=None, model: StruphyModel | None = None)[source]#

Bases: object

Holds the input parameters of a Struphy simulation as attributes.

struphy.post_processing.post_processing_tools.get_params_of_run(path: str) ParamsIn[source]#

Retrieve parameters of finished Struphy run.

Parameters:

path (str) – Absolute path of simulation output folder.

struphy.post_processing.post_processing_tools.create_femfields(path: str, params_in: ParamsIn, *, step: int = 1)[source]#

Creates instances of SplineFunction from distributed Struphy data.

Parameters:
  • path (str) – Absolute path of simulation output folder.

  • params_in (ParamsIn) – Simulation parameters.

  • step (int) – Whether to create FEM fields at every time step (step=1, default), every second time step (step=2), etc.

Returns:

  • fields (dict) – Nested dictionary holding SplineFunction: fields[t][name] contains the Field with the name “name” in the hdf5 file at time t.

  • t_grid (xp.ndarray) – Time grid.

struphy.post_processing.post_processing_tools.eval_femfields(params_in: ParamsIn, fields: dict, *, celldivide: list = [1, 1, 1], physical: bool = False)[source]#

Evaluate FEM fields obtained from struphy.post_processing.post_processing_tools.create_femfields().

Parameters:
  • params_in (ParamsIn) – Simulation parameters.

  • fields (dict) – Obtained from struphy.diagnostics.post_processing.create_femfields.

  • celldivide (list of ints) – Grid refinement in each eta direction.

  • physical (bool) – Wether to do post-processing into push-forwarded physical (xyz) components of fields.

Returns:

  • point_data (dict) – Nested dictionary holding values of FemFields on the grid as list of 3d xp.arrays: point_data[name][t] contains the values of the field with name “name” in fields[t].keys() at time t.

    If physical is True, physical components of fields are saved. Otherwise, logical components (differential n-forms) are saved.

  • grids_log (3-list) – 1d logical grids in each eta-direction with Nel[i]*cell_divide[i] + 1 entries in each direction.

  • grids_phy (3-list) – Mapped (physical) grids obtained by domain(*grids_log).

struphy.post_processing.post_processing_tools.create_vtk(path: str, t_grid: ndarray, grids_phy: list, point_data: dict, *, physical: bool = False)[source]#

Creates structured virtual toolkit files (.vts) for Paraview from evaluated field data.

Parameters:
  • path (str) – Absolute path of where to store the .vts files. Will then be in path/vtk/step_<step>.vts.

  • t_grid (xp.ndarray) – Time grid.

  • grids_phy (3-list) – Mapped (physical) grids obtained from struphy.diagnostics.post_processing.eval_femfields.

  • point_data (dict) – Field data obtained from struphy.diagnostics.post_processing.eval_femfields.

  • physical (bool) – Wether to create vtk for push-forwarded physical (xyz) components of fields.

struphy.post_processing.post_processing_tools.post_process_markers(path_in: str, path_out: str, species: str, domain: Domain, kind: str = 'Particles6D', step: int = 1)[source]#

Computes the Cartesian (x, y, z) coordinates of saved markers during a simulation and writes them to a .npy files and to .txt files. Also saves the weights.

  • .npy files:

    • Particles6D:

      index

      0 |
      1 | 2 | 3 |
      4 | 5 | 6 |
      7 |

      value

      ID

      position (xyz)

      velocities

      weight

    • Particles5D:

      index

      0 |
      1 | 2 | | 3 |

      4

      5

      6 |

      7

      value

      ID

      guiding_center

      v_parallel

      v_perp

      weight

      magn. moment

    • Particles3D:

      index

      0 |
      1 | 2 | 3 |
      4 |

      value

      ID

      position (xyz)

      weight

  • .txt files :

    index

    0 |
    1 | 2 | 3 |
    4 |

    value

    ID

    position (xyz)

    weight

.txt files can be imported to e.g. Paraview, see 08 - Kinetic data for details.

Parameters:
  • path_in (str) – Absolute path of simulation output folder.

  • path_out (str) – Absolute path of where to store the .txt files. Will be in path_out/orbits.

  • species (str) – Name of the species for which the post processing should be performed.

  • domain (Domain) – Domain object.

  • kind (str) – Name of the kinetic kind (Particles6D, Particles5D or Particles3D).

  • step (int, optional) – Whether to do post-processing at every time step (step=1, default), every second time step (step=2), etc.

struphy.post_processing.post_processing_tools.post_process_f(path_in, params_in: ParamsIn, path_out, species, step=1, compute_bckgr=False)[source]#

Computes and saves distribution functions of saved binning data during a simulation.

Parameters:
  • path_in (str) – Absolute path of simulation output folder.

  • params_in (ParamsIn) – Simulation parameters.

  • path_out (str) – Absolute path of where to store the .txt files. Will be in path_out/orbits.

  • species (str) – Name of the species for which the post processing should be performed.

  • step (int, optional) – Whether to do post-processing at every time step (step=1, default), every second time step (step=2), etc.

  • compute_bckgr (bool) – Whether to compute the kinetic background values and add them to the binning data. This is used if non-standard weights are binned.

struphy.post_processing.post_processing_tools.post_process_n_sph(path_in, params_in: ParamsIn, path_out, species, step=1)[source]#

Computes and saves the density n of saved sph data during a simulation.

Parameters:
  • path_in (str) – Absolute path of simulation output folder.

  • params_in (ParamsIn) – Simulation parameters.

  • path_out (str) – Absolute path of where to store the .txt files. Will be in path_out/orbits.

  • species (str) – Name of the species for which the post processing should be performed.

  • step (int, optional) – Whether to do post-processing at every time step (step=1, default), every second time step (step=2), etc.

Particle orbits#

struphy.post_processing.orbits.orbits_tools.post_process_orbit_guiding_center(path_in, path_kinetics_species, species)[source]#

Computes the Cartesian guiding center from saved full-orbit marker orbits (Particles6D) and writes them to a .npy files and to .txt files.

  • .npy files :

    index

    0 |
    1 | 2 | | 3 |

    4

    5

    6

    value

    ID

    guiding_center

    v_parallel

    v_perp

    magn. moment

  • .txt files :

    index

    0 |
    1 | 2 | 3 |

    value

    ID

    guiding_center

.txt file can be imported to e.g. Paraview, see Tutorial 08 - Kinetic data for details..

Parameters:
  • path_in (str) – Absolute path of simulation output folder.

  • path_kinetics_species (str) – Absolute path of where to store the .txt files. Will be saved in path_kinetics_species/guiding_center.

  • species (str) – Name of the species for which the post processing should be performed.

struphy.post_processing.orbits.orbits_tools.post_process_orbit_classification(path_kinetics_species, species)[source]#

Classify guiding center orbits as “passing”, “trapped” or “lost”.

Classification data (0 for “passing”, 1 for “trapped” and -1 for “lost”) is added at the last column(7) of .npy files in a directory “kinetic_data/<name_of_species>/guiding_center/”.

.npy files :

index

0 |
1 | 2 | | 3 |

4

5

6

7

value

ID

guiding_center

v_parallel

v_perp

magn. moment

classification

Parameters:
  • path_kinetics_species (str) – Absolute path of where to store the .txt files. Will be saved in path_kinetics_species/guiding_center.

  • species (str) – Name of the species for which the post processing should be performed.

  • kind (str) – Name of the kinetic kind (Particles6D, Particles5D or Particles3D).

struphy.post_processing.orbits.orbits_kernels.calculate_guiding_center_from_6d(markers: float[:, :], B_cart: float[:, :])[source]#

Calculate guiding center positions \(\overline{\mathbf X}\) from 6d Cartesian phase-space coordinates \((\mathbf x, \mathbf v)\):

\[\overline{\mathbf X} = \mathbf x - \boldsymbol{\rho}_\textrm{L}(\mathbf x, \mathbf v) \,.\]

where \(\boldsymbol{\rho}_\textrm{L}\) is a normalized Larmor-radius vector defined as

\[\boldsymbol{\rho}_\textrm{L}(\mathbf x, \mathbf v) = \frac{1}{B_0(\mathbf x)} \mathbf{b}_0(\mathbf x) \times \mathbf v\,,\]

where \(\mathbf{b}_0\) is the unit magnetic field vector.

Profiling#

struphy.post_processing.cprofile_analyser.get_cprofile_data(path, print_callers=None)[source]#

Prepare Cprofile data and save to “profile_dict.sav”.

Parameters:
  • path (str) – Path to file “profile_tmp” (usually in output folder).

  • print_callers (str) – Part of function name for which to show calling functions.

struphy.post_processing.cprofile_analyser.compare_cprofile_data(path, list_of_funcs=None)[source]#

Print Cprofile data from “profile_dict.sav” to screen (see get_cprofile_data).

Parameters:
  • path (str) – Path to file “profile_dict.sav” (usually in output folder).

  • list_of_funcs (list) – Strings to watch for in “function name” of Cprofile data, allows to look at data of specific functions. If “None”, the 50 functions with the longest cumtime are listed.

struphy.post_processing.cprofile_analyser.replace_keys(d)[source]#

Replace keys from cprofile data with corresponding class names.

Parameters:
  • d (dict) – Dictionary with keys from cprofile_analyser.get_cprofile_data().

  • list_of_funcs (list[str]) – Names for keyword search.