Second-Order Linear Solver (solvcon.parcel.linear)

This package implements a basic second-order, three-dimensional CESE solver that uses the CFL-insensitive c\mbox{-}\tau scheme. The basic solver can be extended to solver for any linear system of first-order hyperbolic PDEs. An example is the velocity-stress equation solver in .velstress.

Numerical Implementation (._algorithm)

Let

  • u_m be the m-th solution variable and m = 1, \ldots,
M. M is the total number of variables.
  • u_{mx_{\mu}} be the component of the gradient of u_m along the x_{\mu} axis in a Cartesian coordinate system. \mu = 1,
2 in two-dimensional space and \mu = 1, 2, 3 in three-dimensional space.

Common Data Structure

sc_linear_algorithm_t

Basic Information of the Solver

int neq
double time
double time_increment

sc_linear_algorithm_t.neq is the number of equations or the number of variables in a mesh cell. sc_linear_algorithm_t.time is the current time of the solver. sc_linear_algorithm_t.time_increment is \Delta t.

Parameters to the c\mbox{-}\tau Scheme

int alpha
double sigma0
double taylor
double cnbfac
double sftfac
double taumin
double tauscale

Metric Arrays for the CESE Method

double *cecnd
double *cevol
double *sfmrc

Group Data

int ngroup
int gdlen
double *grpda

Scalar Parameters

int nsca
double *amsca

Vector Parameters

int nvec
double *amvec

Solution Arrays

double *sol
double *soln
double *solt
double *dsol
double *dsoln
double *stm
double *cfl
double *ocfl

Metric for CEs & SEs

void sc_linear_prepare_ce_3d(sc_mesh_t *msd, sc_linear_algorithm_t *alg)
void sc_linear_prepare_ce_2d(sc_mesh_t *msd, sc_linear_algorithm_t *alg)

Calculate the volume and centroid of conservation elements.

void sc_linear_prepare_sf_3d(sc_mesh_t *msd, sc_linear_algorithm_t *alg)
void sc_linear_prepare_sf_2d(sc_mesh_t *msd, sc_linear_algorithm_t *alg)

Calculate the centroid and normal of hyperplanes of conservation elements and solution elements.

Time Marching

void sc_linear_calc_soln_3d(sc_mesh_t *msd, sc_linear_algorithm_t *alg)
void sc_linear_calc_soln_2d(sc_mesh_t *msd, sc_linear_algorithm_t *alg)

Calculate the solutions of the next half time step ((u_m)_j^{n+1/2}) based on the informaiton at the current time step (n).

void sc_linear_calc_solt_3d(sc_mesh_t *msd, sc_linear_algorithm_t *alg)
void sc_linear_calc_solt_2d(sc_mesh_t *msd, sc_linear_algorithm_t *alg)

Calculate the changing rate of solutions ((u_mt)_j^n).

void sc_linear_calc_jaco_3d(sc_mesh_t *msd, sc_linear_algorithm_t *alg, int icl, double fcn[NEQ][NDIM], double jacos[NEQ][NEQ][NDIM])
void sc_linear_calc_jaco_2d(sc_mesh_t *msd, sc_linear_algorithm_t *alg, int icl, double fcn[NEQ][NDIM], double jacos[NEQ][NEQ][NDIM])

Calculate the Jacobian matrices.

void sc_linear_calc_dsoln_3d(sc_mesh_t *msd, sc_linear_algorithm_t *alg)
void sc_linear_calc_dsoln_2d(sc_mesh_t *msd, sc_linear_algorithm_t *alg)

Calculate the gradient of solutions of the next half time step ((u_{mx_{\mu}})_j^{n+1/2}) based on the information at the current time step (n).

void sc_linear_calc_cfl_3d(sc_mesh_t *msd, sc_linear_algorithm_t *alg)
void sc_linear_calc_cfl_2d(sc_mesh_t *msd, sc_linear_algorithm_t *alg)

Calculate the CFL number based on the eigenvalues of the linear Jacobian matrices.

Plane Wave Solution

void sc_linear_calc_planewave_3d(sc_mesh_t *msd, sc_linear_algorithm_t *alg, double *asol, double *adsol, double *amp, double *ctr, double *wvec, double afreq)
void sc_linear_calc_planewave_2d(sc_mesh_t *msd, sc_linear_algorithm_t *alg, double *asol, double *adsol, double *amp, double *ctr, double *wvec, double afreq)

Calculate the plane-wave solutions.

Wrapper for Numerical Code

class solvcon.parcel.linear.LinearAlgorithm

This class wraps around the C functions for the second-order CESE method.

Numerical Controller (.solver)

class solvcon.parcel.linear.solver.LinearSolver(blk, **kw)

This class controls the underneath algorithm LinearAlgorithm.

Inheritance diagram of LinearSolver

class solvcon.parcel.linear.solver.LinearPeriodic(**kw)

General periodic boundary condition for sequential runs.

Inheritance diagram of LinearPeriodic

Simulation Controller (.case)

class solvcon.parcel.linear.case.LinearCase(**kw)

Basic case with linear CESE method.

Inheritance diagram of LinearCase

Helpers for Plane Wave (.planewave)

class solvcon.parcel.linear.planewave.PlaneWaveSolution(**kw)
class solvcon.parcel.linear.planewave.PlaneWaveAnchor(svr, planewaves=None, **kw)

Use PlaneWaveSolution to calculate plane-wave solutions for LinearSolver.

Inheritance diagram of PlaneWaveAnchor

planewaves = None

Sequence of PlaneWaveSolution objects.

class solvcon.parcel.linear.planewave.PlaneWaveHook(svr, planewaves=None, **kw)

Inheritance diagram of PlaneWaveHook

planewaves = None

Sequence of PlaneWaveSolution objects.

norm = None

A dict holding the calculated norm.

Helpers for I/O (.inout)

class solvcon.parcel.linear.inout.MeshInfoHook(cse, show_bclist=False, perffn=None, **kw)

Print mesh information.

Inheritance diagram of MeshInfoHook

show_bclist = None

Flag to show the list of boundary conditions. Default is False.

perffn = None

Performance file name.

class solvcon.parcel.linear.inout.ProgressHook(cse, linewidth=50, **kw)

Print simulation progess.

Inheritance diagram of ProgressHook

linewidth = None

The maximum width for progress mark.

class solvcon.parcel.linear.inout.FillAnchor(svr, mappers=None, **kw)

Fill the specified arrays of a LinearSolver with corresponding value.

Inheritance diagram of FillAnchor

mappers = None

A dict maps the names of attributes of the MeshAnchor.svr to the filling value.

class solvcon.parcel.linear.inout.CflAnchor(svr, rsteps=None, **kw)

Counting CFL numbers. Use MeshSolver.marchret to return results. Implements postmarch() method.

Inheritance diagram of CflAnchor

rsteps = None

Steps to run (int).

class solvcon.parcel.linear.inout.CflHook(cse, name='cfl', cflmin=0.0, cflmax=1.0, fullstop=True, rsteps=None, **kw)

Makes sure CFL number is bounded and print averaged CFL number over time. Reports CFL information per time step and on finishing. Implements (i) postmarch() and (ii) postloop() methods.

Inheritance diagram of CflHook

rsteps = None

Steps to run.

name = None

Name of the CFL tool.

cflmin = None

Miminum CFL value.

cflmax = None

Maximum CFL value.

fullstop = None

Flag to stop when CFL is out of bound. Default is True.

aCFL = None

Accumulated CFL.

mCFL = None

Mean CFL.

hnCFL = None

Hereditary minimum CFL.

hxCFL = None

Hereditary maximum CFL.

aadj = None

Number of adjusted CFL accumulated since last report.

haadj = None

Total number of adjusted CFL since simulation started.

class solvcon.parcel.linear.inout.MarchSaveAnchor(svr, anames=None, compressor=None, fpdtype=None, psteps=None, vtkfn_tmpl=None, **kw)

Save solution data into VTK XML format for a solver.

Inheritance diagram of MarchSaveAnchor

anames = None

The arrays in LinearSolver or MeshSolver.der to be saved.

compressor = None

Compressor for binary data. Can be either 'gz' or ''.

fpdtype = None

String for floating point data type (NumPy convention).

psteps = None

The interval in step to save data.

vtkfn_tmpl = None

The template string for the VTK file.

class solvcon.parcel.linear.inout.PMarchSave(cse, anames=None, compressor='gz', fpdtype=None, altdir='', altsym='', vtkfn_tmpl=None, **kw)

Save the geometry and variables in a case when time marching in parallel VTK XML format.

Inheritance diagram of PMarchSave

anames = None

The arrays in LinearSolver or MeshSolver.der to be saved. Format is (name, inder, ndim), (name, inder, ndim) ... For ndim > 0 the array is a spatial vector, for ndim == 0 a simple scalar, and ndim < 0 a list of scalar.

compressor = None

Compressor for binary data. Can be either 'gz' or ''.

fpdtype = None

String for floating point data type (NumPy convention).

altdir = None

The alternate directory to save the VTK files.

altsym = None

The symbolic link in basedir pointing to the alternate directory to save the VTK files.

vtkfn_tmpl = None

The template string for the VTK file.

pextmpl = None

Velocity-Stress Equation Solver (.velstress)

See [1] for the basic formulation.

Solver Logic (.velstress.logic)

class solvcon.parcel.linear.velstress.logic.VslinPWSolution(**kw)

Plane-wave solutions for the velocity-stress equations.

class solvcon.parcel.linear.velstress.logic.VslinSolver(blk, mtrldict=None, **kw)

Basic elastic solver.

mtrldict = None

A dict that maps names to Material object.

mtrllist = None

A list of all Material objects.

class solvcon.parcel.linear.velstress.logic.VslinCase(**kw)

Case for anisotropic elastic solids.

Material Definition (.velstress.material)

solvcon.parcel.linear.velstress.material.mltregy = {'GaAs': <class 'solvcon.parcel.linear.velstress.material.GaAs'>, 'Acmite': <class 'solvcon.parcel.linear.velstress.material.Acmite'>, 'Trigonal': <class 'solvcon.parcel.linear.velstress.material.Trigonal'>, 'Monoclinic': <class 'solvcon.parcel.linear.velstress.material.Monoclinic'>, 'Cubic': <class 'solvcon.parcel.linear.velstress.material.Cubic'>, 'BariumTitanate': <class 'solvcon.parcel.linear.velstress.material.BariumTitanate'>, 'Orthorhombic': <class 'solvcon.parcel.linear.velstress.material.Orthorhombic'>, 'Material': <class 'solvcon.parcel.linear.velstress.material.Material'>, 'ZnO': <class 'solvcon.parcel.linear.velstress.material.ZnO'>, 'Albite': <class 'solvcon.parcel.linear.velstress.material.Albite'>, 'RickerSample': <class 'solvcon.parcel.linear.velstress.material.RickerSample'>, 'Beryl': <class 'solvcon.parcel.linear.velstress.material.Beryl'>, 'AlphaUranium': <class 'solvcon.parcel.linear.velstress.material.AlphaUranium'>, 'Isotropic': <class 'solvcon.parcel.linear.velstress.material.Isotropic'>, 'Zinc': <class 'solvcon.parcel.linear.velstress.material.Zinc'>, 'CdS': <class 'solvcon.parcel.linear.velstress.material.CdS'>, 'Tetragonal': <class 'solvcon.parcel.linear.velstress.material.Tetragonal'>, 'Triclinic': <class 'solvcon.parcel.linear.velstress.material.Triclinic'>, 'Hexagonal': <class 'solvcon.parcel.linear.velstress.material.Hexagonal'>, 'RickerSampleLight': <class 'solvcon.parcel.linear.velstress.material.RickerSampleLight'>, 'AlphaQuartz': <class 'solvcon.parcel.linear.velstress.material.AlphaQuartz'>}

Registry class for the name of types.

class solvcon.parcel.linear.velstress.material.Material(rho=None, al=None, be=None, ga=None, **kw)

Material properties. The constitutive relation needs not be symmetric.

rho = None

Density.

al = None

Alpha angle.

be = None

Beta angle.

ga = None

Gamma angle.

origstiff = None

Stiffness matrix in the crystal coordinate.

stiff = None

Stiffness matrix in the transformed global coordinate.

Bibliography

[1]Yung-Yu Chen, Lixiang Yang, and Sheng-Tao John Yu, “Hyperbolicity of Velocity-Stress Equations for Waves in Anisotropic Elastic Solids”, Journal of Elasticity, Volume 106, Number 2, Feb. 2012, Page 149-164. doi: 10.1007/s10659-011-9315-8 <http://dx.doi.org/10.1007/s10659-011-9315-8>.