pyopensn.solver.DiscreteOrdinatesCurvilinearProblem

class pyopensn.solver.DiscreteOrdinatesCurvilinearProblem

Base class for discrete ordinates problems in curvilinear geometry.

Wrapper of opensn::DiscreteOrdinatesCurvilinearProblem.

ComputeFissionProduction(self: pyopensn.solver.LBSProblem, scalar_flux_iterate: str) float

Computes the total fission production (nu*fission).

Parameters:

scalar_flux_iterate ({'old', 'new'}) –

Specifies which scalar flux vector to use in the calculation.
  • ’old’: Use the previous scalar flux iterate.

  • ’new’: Use the current scalar flux iterate.

Returns:

The total fission production.

Return type:

float

Raises:

ValueError – If scalar_flux_iterate is not ‘old’ or ‘new’.

ComputeFissionRate(self: pyopensn.solver.LBSProblem, scalar_flux_iterate: str) float

Computes the total fission rate.

Parameters:

scalar_flux_iterate ({'old', 'new'}) –

Specifies which scalar flux vector to use in the calculation.
  • ’old’: Use the previous scalar flux iterate.

  • ’new’: Use the current scalar flux iterate.

Returns:

The total fission rate.

Return type:

float

Raises:

ValueError – If scalar_flux_iterate is not ‘old’ or ‘new’.

ComputeLeakage(self: pyopensn.solver.DiscreteOrdinatesProblem, bnd_names: list) dict

Compute leakage for the problem.

Parameters:

bnd_names (List[str]) – A list of boundary names for which leakage should be computed.

Returns:

A dictionary mapping boundary names to group-wise leakage vectors. Each array contains the outgoing angular flux (per group) integrated over the corresponding boundary surface.

Return type:

Dict[str, numpy.ndarray]

Raises:
  • RuntimeError – If save_angular_flux option was not enabled during problem setup.

  • ValueError – If one or more boundary ids are not present on the current mesh.

CreateAndWriteSourceMoments(self: pyopensn.solver.LBSProblem, file_base: str) None

Write source moments from latest flux iterate to file.

Parameters:

file_base (str) – File basename.

GetAngularFieldFunctionList(self: pyopensn.solver.DiscreteOrdinatesProblem, groups: list, angles: list) list

Return field functions for selected angular flux components.

This requires Initialize() to have completed; it throws if the problem is not yet initialized. Note: You must enable angular flux storage (save_angular_flux=True) in the problem options, or use a transient/time-dependent solver that retains angular fluxes, otherwise the field functions will remain zero.

Example

`python solver.Initialize() solver.Execute() ang_ff = phys.GetAngularFieldFunctionList(groups=[0], angles=[0]) `

For transient/time-dependent runs, call this after each timestep. Two common patterns are:

1) Use TransientSolver.Execute() with a post-advance callback: ```python solver = TransientSolver(problem=phys) solver.Initialize()

def post_advance():

ang_ff = phys.GetAngularFieldFunctionList(groups=[0], angles=[0]) FieldFunctionGridBased.ExportMultipleToPVTU(ang_ff, “angular_flux_t”)

solver.SetPostAdvanceCallback(post_advance) solver.Execute() ```

2) Use a custom Python loop with TransientSolver.Advance(): ```python solver = TransientSolver(problem=phys) solver.Initialize() for _ in range(num_steps):

solver.Advance() ang_ff = phys.GetAngularFieldFunctionList(groups=[0], angles=[0]) FieldFunctionGridBased.ExportMultipleToPVTU(ang_ff, “angular_flux_t”)

```

Parameters:
  • groups (List[int]) – Global group indices to export.

  • angles (List[int]) – Angle indices within the groupset quadrature for each group.

Returns:

Field functions for the requested (group, angle) pairs.

Return type:

List[pyopensn.fieldfunc.FieldFunctionGridBased]

GetFieldFunctions(self: pyopensn.solver.Problem) list

Get the list of field functions.

Returns:

List of grid-based field functions representing solution data such as scalar fluxes.

Return type:

List[pyopensn.fieldfunc.FieldFunctionGridBased]

GetPhiNewLocal(self: pyopensn.solver.LBSProblem) memoryview

Get the current scalar flux iterate (local vector).

Returns:

Memory view of the local new scalar flux vector.

Return type:

memoryview

GetPhiOldLocal(self: pyopensn.solver.LBSProblem) memoryview

Get the previous scalar flux iterate (local vector).

Returns:

Memory view of the local old scalar flux vector.

Return type:

memoryview

GetPowerFieldFunction(self: pyopensn.solver.LBSProblem) pyopensn.fieldfunc.FieldFunctionGridBased

Returns the power generation field function, if enabled.

GetPsi(self: pyopensn.solver.DiscreteOrdinatesProblem) list

Return psi as a list of NumPy arrays (float64), using zero-copy views into the underlying data.

GetScalarFieldFunctionList(self: pyopensn.solver.LBSProblem, only_scalar_flux: bool = True) list

Return field functions grouped by energy group and, optionally, by moment.

Parameters:

only_scalar_flux (bool, default=True) –

If True, returns only the zeroth moment (scalar flux) field function for each group. The result is a flat list of field functions, one per group.

If False, returns all moment field functions for each group. The result is a nested list where each entry corresponds to a group and contains a list of field functions for all moments (e.g., scalar flux, higher-order moments).

Returns:

The structure of the returned list depends on the only_scalar_flux flag.

Return type:

Union[List[pyopensn.fieldfunc.FieldFunctionGridBased], List[List[pyopensn.fieldfunc.FieldFunctionGridBased]]]

Notes

The moment index varies more rapidly than the group index when only_scalar_flux is False.

GetTime(self: pyopensn.solver.LBSProblem) float

Get the current simulation time in seconds.

GetTimeStep(self: pyopensn.solver.LBSProblem) float

Get the current timestep size.

IsAdjoint(self: pyopensn.solver.LBSProblem) bool

Return True if the problem is in adjoint mode, otherwise False.

IsTimeDependent(self: pyopensn.solver.DiscreteOrdinatesProblem) bool

Return True if the problem is currently in time-dependent mode.

ReadAngularFluxes(self: opensn::DiscreteOrdinatesProblem, file_base: str) None

Read angular fluxes from file.

Parameters:

file_base (str) – File basename.

ReadFluxMoments(self: pyopensn.solver.LBSProblem, file_base: str, single_file_flag: bool) None

Read flux moment data.

Parameters:
  • file_base (str) – File basename.

  • single_file_flag (bool) – True if all flux moments are in a single file.

ReadFluxMomentsAndMakeSourceMoments(self: pyopensn.solver.LBSProblem, file_base: str, single_file_flag: bool) None

Read flux moments and compute corresponding source moments.

Parameters:
  • file_base (str) – File basename.

  • single_file_flag (bool) – True if all flux moments are in a single file.

ReadSourceMoments(self: pyopensn.solver.LBSProblem, file_base: str, single_file_flag: bool) None

Read source moments from file.

Parameters:
  • file_base (str) – File basename.

  • single_file_flag (bool) – True if all source moments are in a single file.

SetAdjoint(self: pyopensn.solver.LBSProblem, arg0: bool) None

Set forward/adjoint transport mode.

Parameters:

adjoint (bool) – True enables adjoint mode and False enables forward mode.

Notes

This is one of three supported mode-setting paths in Python:
  1. options={'adjoint': ...} in the problem constructor.

  2. SetOptions(adjoint=...).

  3. SetAdjoint(...) (this method).

In typical Python workflows, use constructor options or SetOptions.

When this changes mode after the problem has been initialized, OpenSn performs a full mode-transition reset:

  • Materials are reinitialized in the selected mode.

  • Point and volumetric sources are cleared.

  • Boundary conditions are cleared.

  • Scalar and angular flux vectors are zeroed.

The block-id to cross-section map is preserved across the transition. However, the mode change is applied to the mapped MultiGroupXS objects themselves. If those objects are shared with other problems, they observe the same mode toggle.

After a mode change, reapply the desired driving terms before solving, typically:

This routine is intentionally destructive with respect to source/boundary/flux state to avoid hidden coupling between forward and adjoint setups.

SetBoundaryOptions(self: pyopensn.solver.DiscreteOrdinatesProblem, **kwargs) None

Set or clear boundary conditions.

Parameters:
  • clear_boundary_conditions (bool, default=False) – If true, all current boundary conditions are deleted.

  • boundary_conditions (List[Dict]) –

    A list of boundary condition dictionaries. Each dictionary supports:
    • name: str (required)

      Boundary name that identifies the specific boundary.

    • type: {‘vacuum’, ‘isotropic’, ‘reflecting’, ‘arbitrary’} (required)

      Boundary type specification.

    • group_strength: List[float], optional

      Required when type='isotropic'. Isotropic strength per group.

    • function: AngularFluxFunction, optional

      Required when type='arbitrary'. Callable that returns incoming angular flux.

Notes

Mode transitions via SetOptions(adjoint=...) (or the low-level LBSProblem.SetAdjoint()) clear all boundary conditions. Reapply boundaries with this method before solving in the new mode.

SetOptions(self: pyopensn.solver.DiscreteOrdinatesProblem, **kwargs) None

Set problem options from a large list of parameters.

Parameters:

**kwargs – Same option keys and semantics as LBSProblem.SetOptions().

Notes

This call is additive: only options explicitly supplied are updated.

As with LBSProblem.SetOptions(), if adjoint is explicitly supplied and differs from the current mode, the same mode-transition behavior as LBSProblem.SetAdjoint() is used.

SetPointSources(self: pyopensn.solver.LBSProblem, **kwargs) None

Set or clear point sources.

Parameters:
  • clear_point_sources (bool, default=False) – If true, all current the point sources of the problem are deleted.

  • point_sources (List[pyopensn.source.PointSource]) – List of new point sources to be added to the problem.

SetSaveAngularFlux(self: pyopensn.solver.LBSProblem, save: bool) None
SetSteadyStateMode(self: pyopensn.solver.DiscreteOrdinatesProblem) None

Set the problem to steady-state mode.

Notes

You only need to call this if the problem was constructed in time-dependent mode (time_dependent=True) and you want to switch to steady-state mode.

If called before initialization, this only sets the mode flag used at initialization time; no reset is performed because runtime state does not yet exist. Runtime state will be set when Initialize is called. This does not force save_angular_flux in the pre-initialize path.

If called after initialization, this updates problem internals (sweep chunk mode and source-function) while preserving user boundary conditions and fixed sources. If save_angular_flux was auto-enabled when entering time-dependent mode, this call restores the previous value.

SetTimeDependentMode(self: pyopensn.solver.DiscreteOrdinatesProblem) None

Set the problem to time-dependent mode.

Notes

You only need to call this if the problem was constructed in steady-state mode (the default) and you want to switch to time-dependent mode.

If called before initialization, this only sets the mode flag used at initialization time; no reset is performed because runtime state does not yet exist. It may also force save_angular_flux=True because transient mode requires angular-flux storage. Runtime state will be set when Initialize is called.

If called after initialization, this updates problem internals (sweep chunk mode and source-function) while preserving user boundary conditions and fixed sources. If save_angular_flux is currently false, this transition enables it while in time-dependent mode.

SetVolumetricSources(self: pyopensn.solver.LBSProblem, **kwargs) None

Set or clear volumetric sources.

Parameters:
  • clear_volumetric_sources (bool, default=False) – If true, all current the volumetric sources of the problem are deleted.

  • volumetric_sources (List[pyopensn.source.VolumetricSource]) – List of new volumetric sources to be added to the problem.

SetXSMap(self: pyopensn.solver.LBSProblem, **kwargs) None

Replace the block-id to cross-section map.

Parameters:

xs_map (List[Dict]) –

A list of block-id to cross-section mapping dictionaries. Each dictionary supports:
  • block_ids: List[int] (required)

    Mesh block ids to associate with the cross section.

  • xs: pyopensn.xs.MultiGroupXS (required)

    Cross section object.

Notes

Forward/adjoint mode toggles via SetOptions(adjoint=...) (or the low-level LBSProblem.SetAdjoint()) do not change this map. The MultiGroupXS objects themselves are mutable and shared by pointer. If the same MultiGroupXS handle is shared across multiple problems, toggling adjoint mode in one problem also changes the transport mode seen by the others.

WriteAngularFluxes(self: opensn::DiscreteOrdinatesProblem, file_base: str) None

Write angular flux data to file.

Parameters:

file_base (str) – File basename.

WriteFluxMoments(self: pyopensn.solver.LBSProblem, file_base: str) None

Write flux moments to file.

Parameters:

file_base (str) – File basename.

ZeroPhi(self: pyopensn.solver.LBSProblem) None

Zero the scalar-flux vectors (phi_old and phi_new) in-place.

ZeroPsi(self: pyopensn.solver.LBSProblem) None

Zero the angular-flux storage (if present for this problem type).

__init__(self: pyopensn.solver.DiscreteOrdinatesCurvilinearProblem, **kwargs) None

Construct a discrete ordinates problem for curvilinear geometry.

Warning

DiscreteOrdinatesCurvilinearProblem is experimental and should be used with caution!

Parameters:
  • mesh (MeshContinuum) – The spatial mesh.

  • coord_system (int) – Coordinate system to use. Must be set to 2 (cylindrical coordinates).

  • num_groups (int) – The total number of energy groups.

  • groupsets (list of dict) –

    A list of input parameter blocks, each block provides the iterative properties for a groupset. Each dictionary supports:

    • groups_from_to: List[int] (required)

      Two-entry list with the first and last group id for the groupset, e.g. [0, 3].

    • angular_quadrature: pyopensn.aquad.AngularQuadrature, optional

      Handle to an angular quadrature.

    • angle_aggregation_type: {‘polar’, ‘single’, ‘azimuthal’}, default=’polar’

      Angle aggregation method to use during sweeping.

    • angle_aggregation_num_subsets: int, default=1

      Number of angle subsets used for aggregation.

    • inner_linear_method: {‘classic_richardson’, ‘petsc_richardson’, ‘petsc_gmres’, ‘petsc_bicgstab’}, default=’petsc_richardson’

      Iterative method used for inner linear solves.

    • l_abs_tol: float, default=1.0e-6

      Inner linear solver absolute residual tolerance.

    • l_max_its: int, default=200

      Inner linear solver maximum iterations.

    • gmres_restart_interval: int, default=30

      GMRES restart interval, if GMRES is used.

    • allow_cycles: bool, default=True

      Whether cyclic dependencies are allowed in sweeps.

    • apply_wgdsa: bool, default=False

      Enable within-group DSA for this groupset.

    • wgdsa_l_abs_tol: float, default=1.0e-4

      WGDSA linear absolute tolerance.

    • wgdsa_l_max_its: int, default=30

      WGDSA maximum iterations.

    • wgdsa_verbose: bool, default=False

      Verbose WGDSA output.

    • wgdsa_petsc_options: str, default=’’

      PETSc options string for the WGDSA solver.

    • apply_tgdsa: bool, default=False

      Enable two-grid DSA for this groupset.

    • tgdsa_l_abs_tol: float, default=1.0e-4

      TGDSA linear absolute tolerance.

    • tgdsa_l_max_its: int, default=30

      TGDSA maximum iterations.

    • tgdsa_verbose: bool, default=False

      Verbose TGDSA output.

    • tgdsa_petsc_options: str, default=’’

  • xs_map (list of dict) –

    A list of mappings from block ids to cross-section definitions. Each dictionary supports:
    • block_ids: List[int] (required)

      Mesh block IDs to associate with the cross section.

    • xs: pyopensn.xs.MultiGroupXS (required)

      Cross-section object to assign to the specified blocks.

  • boundary_conditions (List[Dict], default=[]) –

    A list containing tables for each boundary specification. Each dictionary supports:
    • name: str (required)

      Boundary name that identifies the specific boundary.

    • type: {‘vacuum’, ‘isotropic’, ‘reflecting’, ‘arbitrary’} (required)

      Boundary type specification.

    • group_strength: List[float], optional

      Required when type='isotropic'. Isotropic strength per group.

    • function: AngularFluxFunction, optional

      Required when type='arbitrary'. Callable that returns incoming angular flux.

  • point_sources (List[pyopensn.source.PointSource], default=[]) – A list of point sources.

  • volumetric_sources (List[pyopensn.source.VolumetricSource], default=[]) – A list of volumetric sources.

  • options (dict, optional) –

    A block of optional configuration parameters. Each dictionary supports the same keys as LBSProblem.SetOptions(), including:

    • max_mpi_message_size: int, default=32768

    • restart_writes_enabled: bool, default=False

    • write_delayed_psi_to_restart: bool, default=True

    • read_restart_path: str, default=’’

    • write_restart_path: str, default=’’

    • write_restart_time_interval: int, default=0

    • use_precursors: bool, default=False

    • use_source_moments: bool, default=False

    • save_angular_flux: bool, default=False

    • verbose_inner_iterations: bool, default=True

    • verbose_outer_iterations: bool, default=True

    • max_ags_iterations: int, default=100

    • ags_tolerance: float, default=1.0e-6

    • ags_convergence_check: {‘l2’, ‘pointwise’}, default=’l2’

    • verbose_ags_iterations: bool, default=True

    • power_field_function_on: bool, default=False

    • power_default_kappa: float, default=3.20435e-11

    • power_normalization: float, default=-1.0

    • field_function_prefix_option: {‘prefix’, ‘solver_name’}, default=’prefix’

    • field_function_prefix: str, default=’’

  • sweep_type (str, optional) – The sweep type to use. Must be one of AAH or CBC. Defaults to AAH.