Chi-Tech
|
Modules | |
lbs.SetOptions | |
Functions | |
void | chi_lua::chiLBSInitializeMaterials (int SolverIndex) |
void | chi_lua::chiLBSAddPointSource (int SolverIndex, double Location_x, double Location_y, double Location_z, table Strength) |
void | chi_lua::chiLBSClearPointSources (int SolverIndex) |
void | chi_lua::chiLBSInitializePointSources (int SolverIndex) |
void | chi_lua::chiLBSSetProperty (int SolverIndex, int PropertyIndex) |
Pair | chi_lua::chiLBSGetScalarFieldFunctionList (int SolverIndex) |
double | chi_lua::chiLBSComputeFissionRate (int SolverIndex, string OldNewOption) |
The | chi_lua::chiLBSComputeLeakage (int SolverIndex, int GroupSetHandle, int BoundaryID) |
void | chi_lua::chiLBSComputeBalance (int SolverIndex) |
Lua functions
void chi_lua::chiLBSAddPointSource | ( | int | SolverIndex, |
double | Location_x, | ||
double | Location_y, | ||
double | Location_z, | ||
table | Strength | ||
) |
Adds a point source to an LBS solver.
SolverIndex | int Handle to the solver. |
Location_x | double X-location. |
Location_y | double Y-location. |
Location_z | double Z-location. |
Strength | table Source strength as a multigroup vector. |
test/framework/tutorials/tutorial_93_raytracing.lua
test/modules/LinearBoltzmannSolvers/Transport_Adjoint/Adjoint2D_3c_response.lua
test/modules/LinearBoltzmannSolvers/Transport_Adjoint/Adjoint2D_2a_forward.lua
test/modules/LinearBoltzmannSolvers/Transport_Adjoint/Adjoint2D_2c_response.lua
test/modules/LinearBoltzmannSolvers/Transport_Adjoint/Adjoint2D_3a_forward.lua
Definition at line 3315 of file lua_functions.c.
void chi_lua::chiLBSClearPointSources | ( | int | SolverIndex | ) |
Clears all the point sources from the solver. This is mostly useful for adjoint response calculations.
SolverIndex | int Handle to the solver. |
Definition at line 3322 of file lua_functions.c.
void chi_lua::chiLBSComputeBalance | ( | int | SolverIndex | ) |
Computes balance tables and prints it to the console.
SolverIndex | int Handle to the solver for which the list is to be obtained. |
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport2D_1Poly_Balance.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport2D_1Poly_BalanceMG.lua
Definition at line 4266 of file lua_functions.c.
double chi_lua::chiLBSComputeFissionRate | ( | int | SolverIndex, |
string | OldNewOption | ||
) |
Computes and returns the fission rate.
SolverIndex | int Handle to the solver maintaining the information. |
OldNewOption | string "NEW" or "OLD". For steady state solvers, the "OLD" option would give the fission rate for the previous iterate. [Default="NEW"] |
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport1D_1.lua
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport1D_3.lua
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport2D_3.lua
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport2D_2.lua
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport1D_2.lua
Definition at line 4101 of file lua_functions.c.
The chi_lua::chiLBSComputeLeakage | ( | int | SolverIndex, |
int | GroupSetHandle, | ||
int | BoundaryID | ||
) |
Computes the leakage for the specified groupset and boundary id.
SolverIndex | int Handle to the solver. |
GroupSetHandle | int Handle to the groupset. |
BoundaryID | int Id of the boundary for which leakage is to be computed. |
Definition at line 4254 of file lua_functions.c.
Pair chi_lua::chiLBSGetScalarFieldFunctionList | ( | int | SolverIndex | ) |
Obtains a list of field functions, related only to scalar flux, from the transport solver.
SolverIndex | int Handle to the solver for which the list is to be obtained. |
test/framework/tutorials/tutorial_93_raytracing.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady_Cyl/Transport2DCyl_1Monoenergetic.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady_Cyl/Transport2DCyl_2Multigroup.lua
test/modules/LinearBoltzmannSolvers/MGDiffusion_KEigen/KEigenvalueMIP2D_1a_QBlock.lua
test/modules/LinearBoltzmannSolvers/MGDiffusion_KEigen/KEigenvalueMIP2D_1b_QBlock.lua
test/modules/LinearBoltzmannSolvers/Transport_SteadyCBC/Transport2D_1Poly.lua
test/modules/LinearBoltzmannSolvers/Transport_SteadyCBC/Transport1D_1.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport2D_3Poly_quad_mod.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport2D_5PolyA_AniHeteroBndry.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport2D_1Poly_Balance.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport3D_5Cycles2.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport2D_4b_DSA_ortho.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport2D_1Poly_BalanceMG.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport2D_1Poly.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport3D_1a_Extruder.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport2D_4a_DSA_ortho.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport3D_6ASplitMesh.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport3D_3a_DSA_ortho.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport3D_1Poly_qmom_part1.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport3D_2Unstructured.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport3D_6BSplitMesh.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport1D_1.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport1D_3a_DSA_ortho.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport3D_1Poly_qmom_part2.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport1D_4_DSA_ortho_inf.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport3D_4Cycles1.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport3D_1b_Ortho.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport2D_2Unstructured.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady/Transport3D_1Poly_parmetis.lua
test/modules/LinearBoltzmannSolvers/Transport_Keigen/KEigenvalueTransport2D_1c_QBlock_CBC.lua
test/modules/LinearBoltzmannSolvers/Transport_Keigen/KEigenvalueTransport2D_1c_QBlock.lua
test/modules/LinearBoltzmannSolvers/Transport_Keigen/KEigenvalueTransport2D_1a_QBlock_CBC.lua
test/modules/LinearBoltzmannSolvers/Transport_Keigen/KEigenvalueTransport2D_1b_QBlock.lua
test/modules/LinearBoltzmannSolvers/Transport_Keigen/KEigenvalueTransport2D_1b_QBlock_CBC.lua
test/modules/LinearBoltzmannSolvers/Transport_Keigen/KEigenvalueTransport2D_1a_QBlock.lua
test/modules/LinearBoltzmannSolvers/Transport_Keigen/c5g7/c5g7.lua
test/modules/LinearBoltzmannSolvers/MGDiffusion_Steady/MIPDiffusion3D_3a_DSA_ortho.lua
test/modules/LinearBoltzmannSolvers/MGDiffusion_Steady/MIPDiffusion3D_1b_Ortho.lua
elements in it (indexed from 1). \
Definition at line 3603 of file lua_functions.c.
void chi_lua::chiLBSInitializeMaterials | ( | int | SolverIndex | ) |
Initializes or reinitializes the materials. This normally happens automatically during solver initialization but if the user wants to swap/change XSs during the run then this will allow the material structures to now deal with the new/changed materials.
SolverIndex | int Handle to the solver maintaining the information. |
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport1D_3.lua
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport2D_3.lua
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport2D_2.lua
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport1D_2.lua
Definition at line 3298 of file lua_functions.c.
void chi_lua::chiLBSInitializePointSources | ( | int | SolverIndex | ) |
Initializes the point sources. This is mostly useful for adjoint response calculations.
SolverIndex | int Handle to the solver. |
Definition at line 3329 of file lua_functions.c.
void chi_lua::chiLBSSetProperty | ( | int | SolverIndex, |
int | PropertyIndex | ||
) |
Set LBS property.
SolverIndex | int Handle to the solver for which the set is to be created. |
PropertyIndex | int Code for a specific property. |
DISCRETIZATION_METHOD
Discretization method.
BOUNDARY_CONDITION
Boundary condition type. See BoundaryIdentify.
SCATTERING_ORDER
Defines the level of harmonic expansion for the scattering source.Default 1. Expects to be followed by an integer.
SWEEP_EAGER_LIMIT
The eager limit to be used in message size during sweep initialization. This expects to be followed by a size in bytes (Max 64,0000).Default 32,000. See note below.
READ_RESTART_DATA
Indicates the reading of restart data from restart file. The value can be followed by two optional strings. The first is the folder name which can be relative or absolute, and the second is the file base name. These are defaulted to "YRestart" and "restart" respectively.
SAVE_ANGULAR_FLUX
Sets the flag for saving the angular flux. Expects to be followed by true/false. [Default=false]
USE_SOURCE_MOMENTS
Flag for using a vector of source moments instead the regular material/boundary source. Default false. This expects to be followed by a boolean.
VERBOSE_INNER_ITERATIONS
Flag for printing inner iteration information. This is primarily used for printing information related to group-set-level iterative methods. Default true. Expects to be followed by a boolean.
VERBOSE_OUTER_ITERATIONS
Flag for printing outer iteration information. This is primarily used for printing information aggregated over group sets such as k-eigenvalue iterations. Default true. Expects to be followed by a boolean.
USE_PRECURSORS
Flag for using delayed neutron precursors. Default false. This expects to be followed by a boolean.
WRITE_RESTART_DATA
Indicates the writing of restart data to restart files. The value can be followed by two optional strings and a number optional strings. The first string is the folder name which can be relative or absolute, and the second string is the file base name. The number is the time interval (in minutes) for a restart write to be triggered (apart from GMRES restarts and the conclusion of groupset completions) .These are defaulted to "YRestart", "restart" and 30 minutes respectively.
PWLD = Piecewise Linear Finite Element.
This value follows the argument BOUNDARY_CONDITION and identifies which boundary is under consideration. Right now only boundaries aligned with cartesian axes are considered. Followed by LBSBoundaryType.
XMAX = Right boundary
XMIN = Left boundary
YMAX = Front boundary
YMIN = Back boundary
ZMAX = Top boundary
ZMIN = Bottom boundary
Specifies the type of boundary. Depending on the type this argument needs to be followed by one or more values. Note: By default all boundaries are type VACUUM.
LBSBoundaryTypes.VACUUM
Specifies a vaccuum boundary condition. It is not followed by any value.
LBSBoundaryTypes.INCIDENT_ISOTROPIC
Incident isotropic flux. This argument needs to be followed by a lua table index 1 to G where G is the amount of energy groups. Note internally this is mapped as 0 to G-1.
LBSBoundaryTypes.REFLECTING
Reflecting boundary condition. Beware, when opposing reflecting boundary conditions are used this enduces a cyclic dependency which will increase the iteration convergence behavior.
LBSBoundaryTypes.INCIDENT_ANISTROPIC_HETEROGENEOUS
Expects to be followed by the name of a lua function. The lua function will get called with the following parameters:
and must return a 1D array of data-values ordered first by angle index, then by group index, e.g., n0g0, n0g1, n0g2, n1g0, n1g1, n1g2, etc.
Example lua function:
The eager limit is the message size limit before which non-blocking MPI send calls will execute without waiting for a matching receive call. The limit is platform dependent but in general 64 kb. Some systems have 32 kb as a limit and therefore we use that as a default limit in ChiTech. There is a fine interplay between message size and the shear amount of messages that will be sent. In general smaller messages tend to be more efficient, however, when there are too many small messages being sent around the communication system on the given platform will start to suffer. One can gain a small amount of parallel efficiency by lowering this limit, however, there is a point where the parallel efficiency will actually get worse so use with caution.
test/framework/tutorials/tutorial_93_raytracing.lua
test/modules/LinearBoltzmannSolvers/Transport_Steady_Cyl/Transport2DCyl_1Monoenergetic.lua
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport1D_1.lua
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport1D_3.lua
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport2D_3.lua
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport2D_2.lua
test/modules/LinearBoltzmannSolvers/Transport_Transient/TransientTransport1D_2.lua
Definition at line 3554 of file lua_functions.c.