An open source method of characteristics neutron transport code.
|
This is an abstract base class which different Solver subclasses implement for different architectures or source iteration algorithms. More...
#include "src/Solver.h"
Public Member Functions | |
Solver (TrackGenerator *track_generator=NULL) | |
Constructor initializes an empty Solver class with array pointers set to NULL. More... | |
virtual | ~Solver () |
Destructor deletes arrays of boundary angular fluxes, scalar fluxes and sources for each FSR and energy group. More... | |
void | setGeometry (Geometry *geometry) |
Sets the Geometry for the Solver. More... | |
Geometry * | getGeometry () |
Returns a pointer to the Geometry. More... | |
TrackGenerator * | getTrackGenerator () |
Returns a pointer to the TrackGenerator. More... | |
FP_PRECISION | getFSRVolume (long fsr_id) |
Returns the calculated volume for a flat source region. More... | |
int | getNumPolarAngles () |
Returns the number of angles used for the polar quadrature. More... | |
int | getNumIterations () |
Returns the number of source iterations to converge the source. More... | |
double | getTotalTime () |
Returns the total time to converge the source (seconds). More... | |
double | getKeff () |
Returns the converged eigenvalue . More... | |
double | getConvergenceThreshold () |
Returns the threshold for source/flux convergence. More... | |
FP_PRECISION | getMaxOpticalLength () |
Get the maximum allowable optical length for a track segment. More... | |
bool | isUsingDoublePrecision () |
Returns whether the solver is using double floating point precision. More... | |
bool | isUsingExponentialInterpolation () |
Returns whether the Solver uses linear interpolation to compute exponentials. More... | |
bool | is3D () |
Returns whether the Solver is tackling a 3D problem. More... | |
void | initializeSolver (solverMode solver_mode) |
Initializes most components of Solver. Mostly needed from the Python side. | |
virtual void | initializeFixedSources () |
Assigns fixed sources assigned by Cell, Material to FSRs. | |
void | printFissionRates (std::string fname, int nx, int ny, int nz) |
Prints fission rates to a binary file. More... | |
virtual void | printInputParamsSummary () |
A function that prints a summary of the input parameters. | |
void | setKeffFromNeutronBalance () |
Sets computation method of k-eff from fission, absorption, and leakage rates rather than from fission rates. keff = fission/(absorption + leakage) | |
void | setResidualByReference (std::string fname) |
Sets residuals to be computed a error relative to a reference. More... | |
void | dumpFSRFluxes (std::string fname) |
Prints scalar fluxes to a binary file. More... | |
void | loadInitialFSRFluxes (std::string fname) |
Load the initial scalar flux distribution from a binary file. More... | |
void | loadFSRFluxes (std::string fname, bool assign_k_eff=false, double tolerance=0.01) |
Load scalar fluxes from a binary file. More... | |
double | getFlux (long fsr_id, int group) |
Returns the scalar flux for some FSR and energy group. More... | |
virtual void | getFluxes (FP_PRECISION *out_fluxes, int num_fluxes)=0 |
double | getFSRSource (long fsr_id, int group) |
Returns the source for some energy group for a flat source region. More... | |
void | setTrackGenerator (TrackGenerator *track_generator) |
Sets the Solver's TrackGenerator with characteristic Tracks. More... | |
void | setConvergenceThreshold (double threshold) |
Sets the threshold for source/flux convergence. More... | |
virtual void | setFluxes (FP_PRECISION *in_fluxes, int num_fluxes)=0 |
virtual void | setFixedSourceByFSR (long fsr_id, int group, double source) |
Assign a fixed source for a flat source region and energy group. More... | |
void | setFixedSourceByCell (Cell *cell, int group, double source) |
Assign a fixed source for a Cell and energy group. More... | |
void | setFixedSourceByMaterial (Material *material, int group, double source) |
Assign a fixed source for a Material and energy group. More... | |
void | setMaxOpticalLength (FP_PRECISION max_optical_length) |
Set the maximum allowable optical length for a track segment. More... | |
void | setExpPrecision (double precision) |
Set the precision, or maximum allowable approximation error, of the the exponential interpolation table. More... | |
void | useExponentialInterpolation () |
Informs the Solver to use linear interpolation to compute the exponential in the transport equation. | |
void | useExponentialIntrinsic () |
Informs the Solver to use the exponential intrinsic exp(...) function to compute the exponential in the transport equation. | |
void | setSolverMode (solverMode solver_mode) |
Choose between direct and adjoint mode. More... | |
void | setRestartStatus (bool is_restart) |
Informs the Solver that this is a 'restart' calculation and therefore k_eff, track angular and region scalar fluxes should not be reset. More... | |
void | allowNegativeFluxes (bool negative_fluxes_on) |
Informs the Solver that this calculation may involve negative fluxes for computing higher eigenmodes for example. More... | |
void | correctXS () |
Directs OpenMOC to correct unphysical cross-sections. More... | |
void | stabilizeTransport (double stabilization_factor, stabilizationType stabilization_type=DIAGONAL) |
Directs OpenMOC to use the diagonal stabilizing correction to the source iteration transport sweep. More... | |
void | setInitialSpectrumCalculation (double threshold) |
Instructs OpenMOC to perform an initial spectrum calculation. More... | |
void | setCheckXSLogLevel (logLevel log_level) |
Determines which log level to set cross-section warnings. More... | |
void | setChiSpectrumMaterial (Material *material) |
Sets the chi spectrum for use as an inital flux guess. More... | |
void | resetMaterials (solverMode mode) |
Returns the Material data to its original state. More... | |
void | fissionTransportSweep () |
This method performs one transport sweep using the fission source. More... | |
void | scatterTransportSweep () |
This method performs one transport sweep using the scatter source. More... | |
void | computeFlux (int max_iters=1000, bool only_fixed_source=true) |
Computes the scalar flux distribution by performing a series of transport sweeps. More... | |
void | computeSource (int max_iters=1000, double k_eff=1.0, residualType res_type=TOTAL_SOURCE) |
Computes the total source distribution by performing a series of transport sweep and source updates. More... | |
void | computeEigenvalue (int max_iters=1000, residualType res_type=FISSION_SOURCE) |
Computes keff by performing a series of transport sweep and source updates. More... | |
void | printBGQMemory () |
virtual void | computeFSRFissionRates (double *fission_rates, long num_FSRs, bool nu=false)=0 |
Computes the volume-weighted, energy integrated fission rate in each FSR and stores them in an array indexed by FSR ID. More... | |
float * | getBoundaryFlux (long track_id, bool fwd) |
Returns the boundary flux array at the requested indexes. More... | |
void | setVerboseIterationReport () |
Sets the solver to print extra information for each iteration. | |
void | printTimerReport () |
Prints a report of the timing statistics to the console. | |
FP_PRECISION * | getFluxesArray () |
A function that returns the array of scalar fluxes. More... | |
void | limitXS () |
Limits cross-sections so that there are no negative cross-sections. More... | |
void | setLimitingXSMaterials (std::vector< int > material_ids, int reset_iteration) |
Instructs MOC to limit negative cross-sections for early iterations. More... | |
void | checkLimitXS (int iteration) |
Checks to see if limited XS should be reset. More... | |
void | setOTFTransport () |
Activate On-The-Fly transport, to OTF ray-trace and propagate the track angular fluxes at the same time. | |
Protected Member Functions | |
virtual void | initializeFluxArrays ()=0 |
Initializes Track boundary angular flux and leakage and FSR scalar flux arrays. | |
virtual void | initializeSourceArrays ()=0 |
Allocates memory for FSR source arrays. | |
virtual void | initializeExpEvaluators () |
Initializes new ExpEvaluator object to compute exponentials. | |
void | initializeMaterials (solverMode mode) |
Initializes the Material's production matrices. More... | |
virtual void | initializeFSRs () |
Initializes the FSR volumes and Materials array. More... | |
void | countFissionableFSRs () |
Counts the number of fissionable flat source regions. More... | |
void | checkXS () |
All material cross-sections in the geometry are checked for consistency. More... | |
virtual void | initializeCmfd () |
Initializes a Cmfd object for acceleration prior to source iteration. More... | |
void | calculateInitialSpectrum (double threshold) |
Performs a spectrum calculation to update the scalar fluxes. More... | |
virtual void | zeroTrackFluxes ()=0 |
Zero each Track's boundary fluxes for each energy group and polar angle in the "forward" and "reverse" directions. | |
virtual void | flattenFSRFluxes (FP_PRECISION value)=0 |
Set the scalar flux for each FSR and energy group to some value. More... | |
virtual void | flattenFSRFluxesChiSpectrum ()=0 |
Set the scalar flux for each FSR to a chi spectrum. | |
virtual void | storeFSRFluxes ()=0 |
Stores the current scalar fluxes in the old scalar flux array. | |
virtual double | normalizeFluxes ()=0 |
Normalizes all FSR scalar fluxes and Track boundary angular fluxes to the total fission source (times ). | |
virtual void | computeStabilizingFlux ()=0 |
Computes the stabilizing flux for transport stabilization. | |
virtual void | stabilizeFlux ()=0 |
Adjusts the scalar flux for transport stabilization. | |
virtual void | computeFSRSources (int iteration)=0 |
Computes the total source (fission, scattering, fixed) for each FSR and energy group. | |
virtual void | computeFSRFissionSources ()=0 |
Computes the total fission source for each FSR and energy group. | |
virtual void | computeFSRScatterSources ()=0 |
Computes the total scattering source for each FSR and energy group. | |
virtual double | computeResidual (residualType res_type)=0 |
Computes the residual between successive flux/source iterations. More... | |
virtual void | computeKeff ()=0 |
Compute from total fission and absorption rates in each FSR and energy group. | |
virtual void | addSourceToScalarFlux ()=0 |
Add the source term contribution in the transport equation to the FSR scalar flux. | |
virtual void | transportSweep ()=0 |
This method performs one transport sweep of all azimuthal angles, Tracks, segments, polar angles and energy groups. | |
void | clearTimerSplits () |
Deletes the Timer's timing entries for each timed code section code in the source convergence loop. More... | |
This is an abstract base class which different Solver subclasses implement for different architectures or source iteration algorithms.
Solver::Solver | ( | TrackGenerator * | track_generator = NULL | ) |
Constructor initializes an empty Solver class with array pointers set to NULL.
track_generator | an optional pointer to a TrackGenerator object |
|
virtual |
Destructor deletes arrays of boundary angular fluxes, scalar fluxes and sources for each FSR and energy group.
Deallocates memory for all arrays allocated for the Solver, including fluxes, sources, quadrature weights, and exponential linear interpolation table.
Delete exponential evaluators
void Solver::allowNegativeFluxes | ( | bool | negative_fluxes_on | ) |
Informs the Solver that this calculation may involve negative fluxes for computing higher eigenmodes for example.
negative_fluxes_on | whether to allow negative fluxes |
|
protected |
Performs a spectrum calculation to update the scalar fluxes.
This function is meant to be used before transport sweeps in an eigenvalue calculation in order to gain a better initial guess on the flux shape. It is equivalent to performing a CMFD update with no current tallies (pure diffusion solve) across a coarse mesh with one mesh cell per domain.
threshold | The convergence threshold of the calculation |
void Solver::checkLimitXS | ( | int | iteration | ) |
Checks to see if limited XS should be reset.
iteration | The MOC iteration number |
|
protected |
All material cross-sections in the geometry are checked for consistency.
Each cross-section is checked to ensure that the total cross-section is greater than or equal to the scattering cross-section for each energy group and that all cross-sections are positive.
|
protected |
Deletes the Timer's timing entries for each timed code section code in the source convergence loop.
To stop and reset all timer splits
void Solver::computeEigenvalue | ( | int | max_iters = 1000 , |
residualType | res_type = FISSION_SOURCE |
||
) |
Computes keff by performing a series of transport sweep and source updates.
This is the main method exposed to the user through the Python interface to perform an eigenvalue calculation. The method makes an initial guess for the scalar and boundary fluxes and performs transport sweeps and source updates until convergence.
By default, this method will perform a maximum of 1000 transport sweeps with a 1E-5 threshold on the integrated FSR fission source. These values may be freely modified by the user at runtime.
The res_type parameter may be used to control the convergence criterion - SCALAR_FLUX, TOTAL_SOURCE and FISSION_SOURCE (default) are all supported options in OpenMOC at this time.
max_iters | the maximum number of source iterations to allow |
res_type | the type of residual used for the convergence criterion |
void Solver::computeFlux | ( | int | max_iters = 1000 , |
bool | only_fixed_source = true |
||
) |
Computes the scalar flux distribution by performing a series of transport sweeps.
This is the main method exposed to the user through the Python interface to compute the scalar flux distribution, e.g., for a fixed source calculation. This routine makes an initial guess for scalar and boundary fluxes and performs transport sweep until convergence.
By default, this method will perform a maximum of 1000 transport sweeps with a 1E-5 threshold on the average FSR scalar flux. These values may be freely modified by the user at runtime.
The only_fixed_source runtime parameter may be used to control the type of source distribution used in the calculation. By default, this paramter is true and only the fixed sources specified by the user will be considered. Alternatively, when the parameter is false, the source will be computed as the scattering and fission sources resulting from a previously computed flux distribution (e.g., an eigenvalue calculation) in addition to any user-defined fixed sources.
This method may be called by the user to compute the scalar flux for a fixed source distribution from Python as follows:
Alternatively, as described above, this method may be called by the user in Python to compute the flux from a superposition of fixed and / or eigenvalue sources as follows:
max_iters | the maximum number of source iterations to allow |
only_fixed_source | use only fixed sources (true by default) |
|
pure virtual |
Computes the volume-weighted, energy integrated fission rate in each FSR and stores them in an array indexed by FSR ID.
This is a helper method for SWIG to allow users to retrieve FSR fission rates as a NumPy array. An example of how this method can be called from Python is as follows:
fission_rates | an array to store the fission rates (implicitly passed in as a NumPy array from Python) |
num_FSRs | the number of FSRs passed in from Python |
nu | whether to return nu-fission rates instead of fission rates |
Implemented in CPUSolver.
|
protectedpure virtual |
void Solver::computeSource | ( | int | max_iters = 1000 , |
double | k_eff = 1.0 , |
||
residualType | res_type = TOTAL_SOURCE |
||
) |
Computes the total source distribution by performing a series of transport sweep and source updates.
This is the main method exposed to the user through the Python interface to compute the source distribution, e.g., for a fixed and/or external source calculation. This routine makes an initial guess for the scalar and boundary fluxes and performs transport sweeps and source updates until convergence.
By default, this method will perform a maximum of 1000 transport sweeps with a 1E-5 threshold on the integrated FSR total source. These values may be freely modified by the user at runtime.
The k_eff parameter may be used for fixed source calculations with fissionable material (e.g., start-up in a reactor from a fixed external source). In this case, the user must "guess" the critical eigenvalue to be be used to scale the fission source.
The res_type parameter may be used to control the convergence criterion - SCALAR_FLUX, TOTAL_SOURCE (default) and FISSION_SOURCE are all supported options in OpenMOC at this time.
This method may be called by the user from Python as follows:
max_iters | the maximum number of source iterations to allow |
k_eff | the sub/super-critical eigenvalue (default 1.0) |
res_type | the type of residual used for the convergence criterion |
void Solver::correctXS | ( | ) |
Directs OpenMOC to correct unphysical cross-sections.
If a material is found with greater total scattering cross-section than total cross-section, the total cross-section is set to the scattering cross-section.
|
protected |
Counts the number of fissionable flat source regions.
This routine is used by the Solver::computeEigenvalue(...) routine which uses the number of fissionable FSRs to normalize the residual on the fission source distribution.
void Solver::dumpFSRFluxes | ( | std::string | fname | ) |
Prints scalar fluxes to a binary file.
fname | the name of the file to dump the fluxes to |
void Solver::fissionTransportSweep | ( | ) |
This method performs one transport sweep using the fission source.
This is a helper routine used for Krylov subspace methods.
|
protectedpure virtual |
Set the scalar flux for each FSR and energy group to some value.
value | the value to assign to each FSR scalar flux |
Implemented in GPUSolver, CPUSolver, and CPULSSolver.
|
inline |
Returns the boundary flux array at the requested indexes.
track_id | The Track's Unique ID |
fwd | Whether the direction of the angular flux along the track is forward (True) or backward (False) |
double Solver::getConvergenceThreshold | ( | ) |
Returns the threshold for source/flux convergence.
double Solver::getFlux | ( | long | fsr_id, |
int | group | ||
) |
Returns the scalar flux for some FSR and energy group.
fsr_id | the ID for the FSR of interest |
group | the energy group of interest |
FP_PRECISION * Solver::getFluxesArray | ( | ) |
A function that returns the array of scalar fluxes.
double Solver::getFSRSource | ( | long | fsr_id, |
int | group | ||
) |
Returns the source for some energy group for a flat source region.
This is a helper routine used by the openmoc.process module.
fsr_id | the ID for the FSR of interest |
group | the energy group of interest |
FP_PRECISION Solver::getFSRVolume | ( | long | fsr_id | ) |
Returns the calculated volume for a flat source region.
fsr_id | the flat source region ID of interest |
Geometry * Solver::getGeometry | ( | ) |
double Solver::getKeff | ( | ) |
Returns the converged eigenvalue .
FP_PRECISION Solver::getMaxOpticalLength | ( | ) |
Get the maximum allowable optical length for a track segment.
int Solver::getNumIterations | ( | ) |
Returns the number of source iterations to converge the source.
int Solver::getNumPolarAngles | ( | ) |
Returns the number of angles used for the polar quadrature.
double Solver::getTotalTime | ( | ) |
Returns the total time to converge the source (seconds).
TrackGenerator * Solver::getTrackGenerator | ( | ) |
Returns a pointer to the TrackGenerator.
|
protectedvirtual |
Initializes a Cmfd object for acceleration prior to source iteration.
Instantiates a dummy Cmfd object if one was not assigned to the Solver by the user and initializes FSRs, materials, fluxes and the Mesh object. This method is for internal use only and should not be called directly by the user.
Reimplemented in CPULSSolver.
|
protectedvirtual |
Initializes the FSR volumes and Materials array.
This method assigns each FSR a unique, monotonically increasing ID, sets the Material for each FSR, and assigns a volume based on the cumulative length of all of the segments inside the FSR.
Reimplemented in GPUSolver, CPUSolver, and CPULSSolver.
|
protected |
Initializes the Material's production matrices.
In an adjoint calculation, this routine will transpose the scattering and fission matrices in each material.
mode | the solution type (FORWARD or ADJOINT) |
bool Solver::is3D | ( | ) |
Returns whether the Solver is tackling a 3D problem.
bool Solver::isUsingDoublePrecision | ( | ) |
Returns whether the solver is using double floating point precision.
bool Solver::isUsingExponentialInterpolation | ( | ) |
Returns whether the Solver uses linear interpolation to compute exponentials.
void Solver::limitXS | ( | ) |
Limits cross-sections so that there are no negative cross-sections.
A copy of the original cross-section is saved
void Solver::loadFSRFluxes | ( | std::string | fname, |
bool | assign_k_eff = false , |
||
double | tolerance = 0.01 |
||
) |
Load scalar fluxes from a binary file.
The matching source regions between the current calculation and those in the loaded file are determined by comparing centroids
fname | The file containing the scalar fluxes |
assign_k_eff | Whether to set k-eff to the one loaded in the binary file |
tolerance | The width of the region in which to search for the matching centroid |
void Solver::loadInitialFSRFluxes | ( | std::string | fname | ) |
Load the initial scalar flux distribution from a binary file.
fname | The file containing the scalar fluxes |
void Solver::printFissionRates | ( | std::string | fname, |
int | nx, | ||
int | ny, | ||
int | nz | ||
) |
Prints fission rates to a binary file.
fname | the name of the file to dump the fission rates to |
nx | number of mesh cells in the x-direction |
ny | number of mesh cells in the y-direction |
nz | number of mesh cells in the z-direction |
void Solver::resetMaterials | ( | solverMode | mode | ) |
Returns the Material data to its original state.
In an adjoint calculation, the scattering and fission matrices in each material are transposed during initialization. This routine returns both matrices to their original (FORWARD) state at the end of a calculation.
mode | the solution type (FORWARD or ADJOINT) |
void Solver::scatterTransportSweep | ( | ) |
This method performs one transport sweep using the scatter source.
This is a helper routine used for Krylov subspace methods.
void Solver::setCheckXSLogLevel | ( | logLevel | log_level | ) |
Determines which log level to set cross-section warnings.
The default log level is ERROR
log_level | The log level for outputing cross-section inconsistencies |
void Solver::setChiSpectrumMaterial | ( | Material * | material | ) |
Sets the chi spectrum for use as an inital flux guess.
material | The material used for obtaining the chi spectrum |
void Solver::setConvergenceThreshold | ( | double | threshold | ) |
Sets the threshold for source/flux convergence.
The default threshold for convergence is 1E-5.
threshold | the threshold for source/flux convergence |
void Solver::setExpPrecision | ( | double | precision | ) |
Set the precision, or maximum allowable approximation error, of the the exponential interpolation table.
By default, the precision is 1E-5 based on the analysis in Yamamoto's 2003 paper.
precision | the precision of the exponential interpolation table, |
void Solver::setFixedSourceByCell | ( | Cell * | cell, |
int | group, | ||
double | source | ||
) |
|
virtual |
Assign a fixed source for a flat source region and energy group.
This is a helper routine to perform error checking for the subclasses which store the source in the appropriate array.
fsr_id | the flat source region ID |
group | the energy group |
source | the volume-averaged source in this group |
void Solver::setFixedSourceByMaterial | ( | Material * | material, |
int | group, | ||
double | source | ||
) |
Assign a fixed source for a Material and energy group.
This routine will add the fixed source to all instances of the Material in the geometry (e.g., all FSRs with this Material).
material | a pointer to the Material of interest |
group | the energy group |
source | the volume-averaged source in this group |
void Solver::setGeometry | ( | Geometry * | geometry | ) |
void Solver::setInitialSpectrumCalculation | ( | double | threshold | ) |
Instructs OpenMOC to perform an initial spectrum calculation.
threshold | The convergence threshold of the spectrum calculation |
void Solver::setLimitingXSMaterials | ( | std::vector< int > | material_ids, |
int | reset_iteration | ||
) |
Instructs MOC to limit negative cross-sections for early iterations.
material_ids | The material IDs of the cross-sections to limit |
reset_iteration | The iteration to reset cross-sections to their defaults |
void Solver::setMaxOpticalLength | ( | FP_PRECISION | max_optical_length | ) |
Set the maximum allowable optical length for a track segment.
max_optical_length | The max optical length |
void Solver::setResidualByReference | ( | std::string | fname | ) |
Sets residuals to be computed a error relative to a reference.
fname | The file containing the flux solution of the reference |
void Solver::setRestartStatus | ( | bool | is_restart | ) |
Informs the Solver that this is a 'restart' calculation and therefore k_eff, track angular and region scalar fluxes should not be reset.
is_restart | whether solver run is a restart run |
void Solver::setSolverMode | ( | solverMode | solver_mode | ) |
Choose between direct and adjoint mode.
solver_mode | openmoc.FORWARD or .ADJOINT |
void Solver::setTrackGenerator | ( | TrackGenerator * | track_generator | ) |
Sets the Solver's TrackGenerator with characteristic Tracks.
The TrackGenerator must already have generated Tracks and have used ray tracing to segmentize them across the Geometry. This should be initated in Python prior to assigning the TrackGenerator to the Solver:
track_generator | a pointer to a TrackGenerator object |
void Solver::stabilizeTransport | ( | double | stabilization_factor, |
stabilizationType | stabilization_type = DIAGONAL |
||
) |
Directs OpenMOC to use the diagonal stabilizing correction to the source iteration transport sweep.
The source iteration process which MOC uses can be unstable if negative cross-sections arise from transport correction. This instability causes issues in convergence. The stabilizing correction fixes this by adding a diagonal matrix to both sides of the discretized transport equation which introduces no bias but transforms the iteration matrix into one that is stable.
Three stabilization options exist: DIAGONAL, YAMAMOTO, and GLOBAL.
DIAGONAL: The stabilization is only applied to fluxes where the associated in-scatter cross-section is negative. The added stabilizing flux is equal to the magnitude of the in-scatter cross-section divided by the total cross-section and scaled by the stabilization factor. YAMAMOTO: This is the same as DIAGONAL except that the largest stabilization is applied to all regions, not just those containing negative in-scatter cross-sections. GLOBAL: This method applies a global stabilization factor to all fluxes defined by the user. In addition, the stabilization factor in this option refers to a damping factor, not the magnitude of the stabilizing correction.
stabilization_factor | The factor applied to the stabilizing correction |
stabilization_type | The type of stabilization to use |
|
protected |
|
protected |
The angular leakages for each Track. This array stores the weighted outgoing angular fluxes for use in non-CMFD eigenvalue calculations.
|
protected |
Boolean for whether to perform a spectrum calculation for the initial flux guess
|
protected |
Boolean for whether to calculate residuals from reference flux
|
protected |
Material to be used for calculating the initial flux guess from chi
|
protected |
A pointer to a Coarse Mesh Finite Difference (CMFD) acceleration object
|
protected |
The tolerance for converging the source/flux
|
protected |
Boolean for whether to correct unphysical cross-sections
|
protected |
A matrix of ExpEvaluators to compute exponentials in the transport equation. The matrix is indexed by azimuthal index and polar index
|
protected |
A mapping of fixed sources keyed by the pair (Cell*, energy group)
|
protected |
A mapping of fixed sources keyed by the pair (FSR ID, energy group)
|
protected |
A mapping of fixed sources keyed by the pair (Material*, energy group)
|
protected |
Optional user-specified fixed sources in each FSR and energy group
|
protected |
Boolean to indicate whether there are any fixed sources
|
protected |
The number of flux varies stored per track in each direction
|
protected |
The FSR "volumes" (i.e., areas) indexed by FSR UID
|
protected |
Temporary scratch pads for intermediate storage of computing steps
|
protected |
File to load initial FSR fluxes from
|
protected |
Convergence threshold for the initial spectrum calculation
|
protected |
Indicator of whether the flux array has been defined by the user
|
protected |
The current iteration's approximation to k-effective
|
protected |
How to compute the k-effective when not using CMFD
|
protected |
Map of the materials with limited cross sections
|
protected |
Whether to limit cross sections or not
|
protected |
Ids of materials that will have their cross section limited
|
protected |
Boolean for whether to load initial FSR flux profile from file
|
protected |
Boolean for whether the solver allows negative fluxes
|
protected |
The number of azimuthal angles
|
protected |
The number of exponential evaluators in the azimuthal direction
|
protected |
The number of exponential evaluators in the polar direction
|
protected |
The number of fissionable flat source regions
|
protected |
The number of flat source regions
|
protected |
The number of energy groups
|
protected |
The number of source iterations needed to reach convergence
|
protected |
The number of Materials
|
protected |
The number of polar angles
|
protected |
The old scalar flux for each energy group in each FSR
|
protected |
Map of the original materials
|
protected |
A pointer to a polar quadrature
|
protected |
Ratios of source to total cross-section for each FSR and energy group
|
protected |
File to load reference FSR fluxes from
|
protected |
The reference scalar flux for each energy group in each FSR
|
protected |
Iteration at which to revert to original cross sections
|
protected |
The scalar flux for each energy group in each FSR
|
protected |
Determines the type of track segmentation to use for 3D MOC
|
protected |
Boolean for whether to solve in 3D (true) or 2D (false)
|
protected |
Whether the solver is using the direct or adjoint mode
|
protected |
A string indicating the type of source approximation
|
protected |
The factor applied to the source iteration stabilization
|
protected |
The type of source iteration stabilization
|
protected |
Boolean for whether to apply the stabilizing correction to the source iteration (transport sweep) process
|
protected |
The stabilizing flux for each energy group in each FSR
|
protected |
A timer to record timing data for a simulation
|
protected |
The total number of Tracks
|
protected |
A pointer to a TrackGenerator which contains Tracks
|
protected |
A pointer to the Tracks array, only for GPUSolver and VectorizedSolver
|
protected |
A pointer to an array with the number of Tracks per azimuthal angle
|
protected |
Boolean to indicate whether a user loaded his own fluxes
|
protected |
Boolean for whether to print verbose iteration reports
|
protected |
The log level for outputting cross-section inconsistencies