Chi-Tech
test/modules/LinearBoltzmannSolvers/MGDiffusion_Steady/MIPDiffusion3D_1b_Ortho.lua
-- 3D Transport test with Vacuum and Incident-isotropic BC.
-- SDM: PWLD
-- Test: Max-value=5.28310e-01 and 8.04576e-04
num_procs = 4
if (reflecting == nil) then reflecting = true end
--############################################### Check num_procs
if (check_num_procs==nil and chi_number_of_processes ~= num_procs) then
chiLog(LOG_0ERROR,"Incorrect amount of processors. " ..
"Expected "..tostring(num_procs)..
". Pass check_num_procs=false to override if possible.")
os.exit(false)
end
--############################################### Setup mesh
nodes={}
N=10
L=5
xmin = -L/2
dx = L/N
for i=1,(N+1) do
k=i-1
nodes[i] = xmin + k*dx
end
znodes={}
for i=1,(N/2+1) do
k=i-1
znodes[i] = xmin + k*dx
end
if (reflecting) then
meshgen1 = chi_mesh.OrthogonalMeshGenerator.Create({ node_sets = {nodes,nodes,znodes} })
else
meshgen1 = chi_mesh.OrthogonalMeshGenerator.Create({ node_sets = {nodes,nodes,nodes} })
end
--############################################### Set Material IDs
vol0 = chi_mesh.RPPLogicalVolume.Create({infx=true, infy=true, infz=true})
--############################################### Add materials
materials = {}
materials[1] = chiPhysicsAddMaterial("Test Material");
chiPhysicsMaterialAddProperty(materials[1],TRANSPORT_XSECTIONS)
chiPhysicsMaterialAddProperty(materials[1],ISOTROPIC_MG_SOURCE)
num_groups = 21
chiPhysicsMaterialSetProperty(materials[1],TRANSPORT_XSECTIONS,
CHI_XSFILE,"../Transport_Steady/xs_graphite_pure.cxs")
src={}
for g=1,num_groups do
src[g] = 0.0
end
src[1] = 1.0
chiPhysicsMaterialSetProperty(materials[1],ISOTROPIC_MG_SOURCE,FROM_ARRAY,src)
--############################################### Setup Physics
pquad0 = chiCreateProductQuadrature(GAUSS_LEGENDRE_CHEBYSHEV,2, 2)
lbs_block =
{
num_groups = num_groups,
groupsets =
{
{
groups_from_to = {0, 20},
angular_quadrature_handle = pquad0,
angle_aggregation_type = "single",
angle_aggregation_num_subsets = 1,
groupset_num_subsets = 1,
inner_linear_method = "gmres",
l_abs_tol = 1.0e-6,
l_max_its = 300,
gmres_restart_interval = 100,
},
}
}
lbs_options =
{
scattering_order = 1,
}
if (reflecting) then
lbs_options.boundary_conditions = {{name = "zmax", type = "reflecting"}}
end
phys1 = lbs.DiffusionDFEMSolver.Create(lbs_block)
lbs.SetOptions(phys1, lbs_options)
--############################################### Initialize and Execute Solver
ss_solver = lbs.SteadyStateSolver.Create({lbs_solver_handle = phys1})
chiSolverExecute(ss_solver)
--############################################### Get field functions
fflist,count = chiLBSGetScalarFieldFunctionList(phys1)
--############################################### Slice plot
--slices = {}
--for k=1,count do
-- slices[k] = chiFFInterpolationCreate(SLICE)
-- chiFFInterpolationSetProperty(slices[k],SLICE_POINT,0.0,0.0,0.8001)
-- chiFFInterpolationSetProperty(slices[k],ADD_FIELDFUNCTION,fflist[k])
-- --chiFFInterpolationSetProperty(slices[k],SLICE_TANGENT,0.393,1.0-0.393,0)
-- --chiFFInterpolationSetProperty(slices[k],SLICE_NORMAL,-(1.0-0.393),-0.393,0.0)
-- --chiFFInterpolationSetProperty(slices[k],SLICE_BINORM,0.0,0.0,1.0)
--end
--############################################### Volume integrations
curffi = ffi1
chiFFInterpolationSetProperty(curffi,OPERATION,OP_MAX)
chiFFInterpolationSetProperty(curffi,LOGICAL_VOLUME,vol0)
chiFFInterpolationSetProperty(curffi,ADD_FIELDFUNCTION,fflist[1])
chiLog(LOG_0,string.format("Max-value1=%.5e", maxval))
curffi = ffi1
chiFFInterpolationSetProperty(curffi,OPERATION,OP_MAX)
chiFFInterpolationSetProperty(curffi,LOGICAL_VOLUME,vol0)
chiFFInterpolationSetProperty(curffi,ADD_FIELDFUNCTION,fflist[20])
chiLog(LOG_0,string.format("Max-value2=%.5e", maxval))
--############################################### Exports
if (master_export == nil) then
if (reflecting) then
chiExportMultiFieldFunctionToVTK(fflist,"ZPhi3DReflected")
else
end
end
--############################################### Plots
if (chi_location_id == 0 and master_export == nil) then
--os.execute("python ZPFFI00.py")
----os.execute("python ZPFFI11.py")
--local handle = io.popen("python ZPFFI00.py")
print("Execution completed")
end
-- DO
--[0] Max-value1=2.52092e+00
--[0] Max-value2=5.79100e-03
-- MGDiffusion
--[0] Max-value1=2.74873e-01
--[0] Max-value2=9.47508e-05
virtual void Execute()
Pair chiLBSGetScalarFieldFunctionList(int SolverIndex)
void chiFFInterpolationExportPython(int FFIHandle, char BaseName)
Handle chiFFInterpolationCreate(int FFITypeIndex)
void chiFFInterpolationExecute(int FFIHandle)
void chiFFInterpolationGetValue(int FFIHandle)
void chiFFInterpolationInitialize(int FFIHandle)
Handle chiFFInterpolationSetProperty(int FFIHandle, int PropertyIndex)
void chiLog(int LogType, char LogMsg)
MaterialHandle chiPhysicsAddMaterial(char Name)
void chiPhysicsMaterialAddProperty(int MaterialHandle, int PropertyIndex)
void chiPhysicsMaterialSetProperty(int MaterialHandle, int PropertyIndex, int OperationIndex, varying Information)
Returns chiCreateProductQuadrature(int QuadratureType, varying values)
void chiVolumeMesherSetMatIDToAll(int material_id)
void chiSolverExecute(int solver_handle)
void chiSolverInitialize(int solver_handle)
Definition: lua_functions.c:92
void chiExportMultiFieldFunctionToVTK(table listFFHandles, char BaseName)
void Set(VecDbl &x, const double &val)