Chi-Tech
Simplified LDFESQ

Functions

void chi_lua::chiLocallyRefineSLDFESQAngularQuadrature (in handle, vec3 reference_direction, double cone_size, bool invert_logic)
 
void chi_lua::chiCreateSLDFESQAngularQuadrature (int initial_refinement_level)
 
void chi_lua::chiPrintToPythonSLDFESQAngularQuadrature (int handle, string file_name_prefix)
 

Detailed Description

This quadrature is a simplified implementation of the quadrature defined in the paper by Cheuk Lau and Marvin Adams: "Discrete Ordinates Quadratures Based on Linear and Quadratic Discontinuous Finite Elements over Spherical Quadrilaterals", Nuclear Science and Engineering, 185:1, pages 36-52, 2017.

LDFE Spherical Quadrilaterals at a refinement level of n=40

Synopsis

Specific notes for the simplified version

Initial sub-division of the inscribed cube face

The orthogonal grid in the xy-tilde reference frame is dependent on a a diagonal spacing of grid points. The paper by Lau and Adams did not specify exactly how to generate this spacing but only stated that their spacing minimizes the max/avg and max/min ratio of projected SQ to 1.1 and 1.3 respectively. We found that the application of a weighting function, $ \alpha (\cos \beta \frac{\pi}{2} x - \cos \beta \frac{\pi}{2})$ produces the optimal spacing. See the code for additional clarification.

Placement of the quadrature points

Instead of using a multi-variate secant method to place the quadrature points we found that placing the points on the 4 Gauss-Legendre points for a quadrilateral ( $[\pm \frac{1}{\sqrt{3}},\pm \frac{1}{\sqrt{3}} ]$)still produced 4th order convergence and added a lot of speed to the algorithm.

Integration of basis functions

The determinant of the Jacobian can be hard to derive. It can be done in the two ways: using a cross-product or in angle-space.

Local refinement

The data structures employed allows the SQs to easily be refined in certain regions.

Function Documentation

◆ chiCreateSLDFESQAngularQuadrature()

void chi_lua::chiCreateSLDFESQAngularQuadrature ( int  initial_refinement_level)

Creates a Simplified Linear Discontinuous Finite Element (SLDFE) quadrature based on Spherical Quadrilaterals (SQ). Hence SLDFE-SQ.

Parameters
initial_refinement_levelint Initial refinement level, $n$ to be used. The total number of angles will be $ 8{\times}12(n+1)^2 $.

_

Example:

Example with refinement level 2.

void chiCreateSLDFESQAngularQuadrature(int initial_refinement_level)

With direction points:

Author
Jan

Definition at line 2531 of file lua_functions.c.

◆ chiLocallyRefineSLDFESQAngularQuadrature()

void chi_lua::chiLocallyRefineSLDFESQAngularQuadrature ( in  handle,
vec3  reference_direction,
double  cone_size,
bool  invert_logic 
)

Applies a local refinement of angles.

Parameters
handleint. Handle to the reference quadrature.
reference_directionvec3 Reference vector. $ \vec{r} $
cone_sizedouble Cone size in radians. $ \theta $
invert_logicbool Optional[Default:false]. If supplied, interprets SQ-splitting as when $|\omega \cdot \vec{r}| < \sin(\theta) $. Otherwise, SQs will be split if $ \omega \cdot \vec{r} > \cos(\theta)$

_

Example:

Example with refinement level 2 and a triple directional refinement:

chiLocallyRefineSLDFESQAngularQuadrature(pquad,{1,0,0},45.0*math.pi/180,false)
chiLocallyRefineSLDFESQAngularQuadrature(pquad,{1,0,0},23.0*math.pi/180,false)
chiLocallyRefineSLDFESQAngularQuadrature(pquad,{1,0,0},12.0*math.pi/180,false)
void chiLocallyRefineSLDFESQAngularQuadrature(in handle, vec3 reference_direction, double cone_size, bool invert_logic)

Example with refinement level 2 and a triple planar refinement:

chiLocallyRefineSLDFESQAngularQuadrature(pquad,{1,0,0},22.50*math.pi/180,true)
chiLocallyRefineSLDFESQAngularQuadrature(pquad,{1,0,0},11.75*math.pi/180,true)
chiLocallyRefineSLDFESQAngularQuadrature(pquad,{1,0,0},5.000*math.pi/180,true)
Author
Jan

Definition at line 2510 of file lua_functions.c.

◆ chiPrintToPythonSLDFESQAngularQuadrature()

void chi_lua::chiPrintToPythonSLDFESQAngularQuadrature ( int  handle,
string  file_name_prefix 
)

Outputs the quadrature information to python format.

Parameters
handleint Handle to the reference quadrature.
file_name_prefixstring Prefix to be used in front of file.

_

Example:

Example of printing a quadrature: Example with refinement level 2 and a triple directional refinement:

chiLocallyRefineSLDFESQAngularQuadrature(pquad,{1,0,0},45.0*math.pi/180,false)
chiLocallyRefineSLDFESQAngularQuadrature(pquad,{1,0,0},23.0*math.pi/180,false)
chiLocallyRefineSLDFESQAngularQuadrature(pquad,{1,0,0},12.0*math.pi/180,false)
void chiPrintToPythonSLDFESQAngularQuadrature(int handle, string file_name_prefix)
Author
Jan

Definition at line 2554 of file lua_functions.c.