Advection diffusion equation Reference: https://en.wikipedia.org/wiki/Convection%E2%80%93diffusion_equation
- class modulus.PDES.advection_diffusion.AdvectionDiffusion(T='T', D='D', Q=0, rho='rho', dim=3, time=False)
Bases:
modulus.pdes.PDES
Advection diffusion equation
- Tstr
- Dfloat, Sympy Symbol/Expr, str
- Qfloat, Sympy Symbol/Expr, str
- rhofloat, Sympy Symbol/Expr, str
- dimint
- timebool
The dependent variable.
Diffusivity. If D is a str then it is converted to Sympy Function of form ‘D(x,y,z,t)’. If ‘D’ is a Sympy Symbol or Expression then this is substituted into the equation.
The source term. If Q is a str then it is converted to Sympy Function of form ‘Q(x,y,z,t)’. If ‘Q’ is a Sympy Symbol or Expression then this is substituted into the equation. Default is 0.
The density. If rho is a str then it is converted to Sympy Function of form ‘rho(x,y,z,t)’. If ‘rho’ is a Sympy Symbol or Expression then this is substituted into the equation to allow for compressible Navier Stokes.
Dimension of the diffusion equation (1, 2, or 3). Default is 3.
If time-dependent equations or not. Default is False.
>>> ad = AdvectionDiffusion(D=0.1, rho=1.) >>> ad.pprint() advection_diffusion: u*T__x + v*T__y + w*T__z - 0.1*T__x__x - 0.1*T__y__y - 0.1*T__z__z >>> ad = AdvectionDiffusion(D='D', rho=1, dim=2, time=True) >>> ad.pprint() advection_diffusion: -D*T__x__x - D*T__y__y + u*T__x + v*T__y - D__x*T__x - D__y*T__y + T__t
Basic equations
- class modulus.PDES.basic.Curl(vector, curl_name=['u', 'v', 'w'])
Bases:
modulus.pdes.PDES
del cross vector operator
- vectortuple of 3 Sympy Exprs, floats, ints or strings
- curl_nametuple of 3 strings
This will be the vector to take the curl of.
These will be the output names of the curl operations.
>>> c = Curl((0,0,'phi'), ('u','v','w')) >>> c.pprint() u: phi__y v: -phi__x w: 0
- class modulus.PDES.basic.GradNormal(T, dim=3, time=True)
Bases:
modulus.pdes.PDES
Implementation of the gradient boundary condition
- Tstr
- dimint
- timebool
The dependent variable.
Dimension of the equations (1, 2, or 3). Default is 3.
If time-dependent equations or not. Default is True.
>>> gn = ns = GradNormal(T='T') >>> gn.pprint() normal_gradient_T: normal_x*T__x + normal_y*T__y + normal_z*T__z
- class modulus.PDES.basic.NormalDotVec(vec=['u', 'v', 'w'])
Bases:
modulus.pdes.PDES
Normal dot velocity
- dimint
Dimension of the equations (1, 2, or 3). Default is 3.
Diffusion equation
- class modulus.PDES.diffusion.Diffusion(T='T', D='D', Q=0, dim=3, time=True)
Bases:
modulus.pdes.PDES
Diffusion equation
- Tstr
- Dfloat, Sympy Symbol/Expr, str
- Qfloat, Sympy Symbol/Expr, str
- dimint
- timebool
The dependent variable.
Diffusivity. If D is a str then it is converted to Sympy Function of form ‘D(x,y,z,t)’. If ‘D’ is a Sympy Symbol or Expression then this is substituted into the equation.
The source term. If Q is a str then it is converted to Sympy Function of form ‘Q(x,y,z,t)’. If ‘Q’ is a Sympy Symbol or Expression then this is substituted into the equation. Default is 0.
Dimension of the diffusion equation (1, 2, or 3). Default is 3.
If time-dependent equations or not. Default is True.
>>> diff = Diffusion(D=0.1, Q=1, dim=2) >>> diff.pprint() diffusion_T: T__t - 0.1*T__x__x - 0.1*T__y__y - 1 >>> diff = Diffusion(T='u', D='D', Q='Q', dim=3, time=False) >>> diff.pprint() diffusion_u: -D*u__x__x - D*u__y__y - D*u__z__z - Q - D__x*u__x - D__y*u__y - D__z*u__z
- class modulus.PDES.diffusion.DiffusionInterface(T_1, T_2, D_1, D_2, dim=3, time=True)
Bases:
modulus.pdes.PDES
Matches the boundary conditions at an interface
- T_1, T_2str
- D_1, D_2float
- dimint
- timebool
Dependent variables to match the boundary conditions at the interface.
Diffusivity at the interface.
Dimension of the equations (1, 2, or 3). Default is 3.
If time-dependent equations or not. Default is True.
>>> diff = DiffusionInterface('theta_s', 'theta_f', 0.1, 0.05, dim=2) >>> diff.pprint() diffusion_interface_dirichlet_theta_s_theta_f: -theta_f + theta_s diffusion_interface_neumann_theta_s_theta_f: -0.05*normal_x*theta_f__x + 0.1*normal_x*theta_s__x - 0.05*normal_y*theta_f__y + 0.1*normal_y*theta_s__y
Maxwell’s equation
- class modulus.PDES.electromagnetic.MaxwellFreqReal(ux='ux', uy='uy', uz='uz', k=1.0)
Bases:
modulus.pdes.PDES
Frequency domain Maxwell’s equation
- uxstr
- uystr
- uzstr
- kfloat, Sympy Symbol/Expr, str
Ex
Ey
Ez
Wave number. If k is a str then it is converted to Sympy Function of form ‘k(x,y,z,t)’. If ‘k’ is a Sympy Symbol or Expression then this is substituted into the equation.
- class modulus.PDES.electromagnetic.PEC(ux='ux', uy='uy', uz='uz', dim=3)
Bases:
modulus.pdes.PDES
Perfect Electric Conduct BC for
- uxstr
- uystr
- uzstr
- dimint
Ex
Ey
Ez
Dimension of the equations (2, or 3). Default is 3.
- class modulus.PDES.electromagnetic.SommerfeldBC(ux='ux', uy='uy', uz='uz')
Bases:
modulus.pdes.PDES
Frequency domain ABC, Sommerfeld radiation condition Only for real part Equation: ‘n x _curl(E) = 0’
- uxstr
- uystr
- uzstr
Ex
Ey
Ez
Energy equation Reference: https://www.comsol.com/multiphysics/heat-transfer-conservation-of-energy http://dl.icdst.org/pdfs/files1/2fe68e957cdf09a4862088ed279f00b0.pdf http://farside.ph.utexas.edu/teaching/336L/Fluidhtml/node14.html#e4.67
- class modulus.PDES.energy_equation.EnergyFluid(cp='cp', kappa='kappa', rho='rho', nu='nu', visc_heating=False, dim=3, time=False)
Bases:
modulus.pdes.PDES
Energy equation Supports compressible flow. For Ideal gases only (uses the assumption that beta*T = 1). No heat/energy source added.
- cpstr
- kappastr
- rhoSympy Symbol/Expr, str
- nuSympy Symbol/Expr, str
- visc_heatingbool
- dimint
- timebool
The specific heat.
The conductivity.
The density. If rho is a str then it is converted to Sympy Function of form ‘rho(x,y,z,t)’. If ‘rho’ is a Sympy Symbol or Expression then this is substituted into the equation.
The kinematic viscosity. If nu is a str then it is converted to Sympy Function of form ‘nu(x,y,z,t)’. If ‘nu’ is a Sympy Symbol or Expression then this is substituted into the equation.
If viscous heating is applied or not. Default is False.
Dimension of the energy equation (2 or 3). Default is 3.
If time-dependent equations or not. Default is False.
>>> ene = EnergyFluid(nu=0.1, rho='rho', cp=2.0, kappa=5, dim=2, time=False, visc_heating=False) >>> ene.pprint() temperauture_fluid: 2.0*(u(x, y)*Derivative(T(x, y), x) + v(x, y)*Derivative(T(x, y), y))*rho(x, y) - u(x, y)*Derivative(p(x, y), x) - v(x, y)*Derivative(p(x, y), y) - 5*Derivative(T(x, y), (x, 2)) - 5*Derivative(T(x, y), (y, 2))
Equations related to linear elasticity
- class modulus.PDES.linear_elasticity.LinearElasticity(E=None, nu=None, lambda_=None, mu=None, rho=1, dim=3, time=False)
Bases:
modulus.pdes.PDES
Linear elasticity equations. Use either (E, nu) or (lambda_, mu) to define the material properties.
- Efloat, Sympy Symbol/Expr, str
- nufloat, Sympy Symbol/Expr, str
- lambda_: float, Sympy Symbol/Expr, str
- mu: float, Sympy Symbol/Expr, str
- rho: float, Sympy Symbol/Expr, str
- dimint
The Young’s modulus
The Poisson’s ratio
Lamé’s first parameter
Lamé’s second parameter (shear modulus)
Mass density.
Dimension of the linear elasticity (2 or 3). Default is 3.
>>> elasticity_equations = LinearElasticity(E=10, nu=0.3, dim=2) >>> elasticity_equations.pprint() navier_x: -13.4615384615385*u__x__x - 3.84615384615385*u__y__y - 9.61538461538461*v__x__y navier_y: -3.84615384615385*v__x__x - 13.4615384615385*v__y__y - 9.61538461538461*u__x__y stress_disp_xx: -sigma_xx + 13.4615384615385*u__x + 5.76923076923077*v__y stress_disp_yy: -sigma_yy + 5.76923076923077*u__x + 13.4615384615385*v__y stress_disp_xy: -sigma_xy + 3.84615384615385*u__y + 3.84615384615385*v__x equilibrium_x: -sigma_xx__x - sigma_xy__y equilibrium_y: -sigma_xy__x - sigma_yy__y traction_x: normal_x*sigma_xx + normal_y*sigma_xy traction_y: normal_x*sigma_xy + normal_y*sigma_yy
- class modulus.PDES.linear_elasticity.LinearElasticityPlaneStress(E=None, nu=None, lambda_=None, mu=None, rho=1, time=False)
Bases:
modulus.pdes.PDES
Linear elasticity plane stress equations. Use either (E, nu) or (lambda_, mu) to define the material properties.
- Efloat, Sympy Symbol/Expr, str
- nufloat, Sympy Symbol/Expr, str
- lambda_: float, Sympy Symbol/Expr, str
- mu: float, Sympy Symbol/Expr, str
- rho: float, Sympy Symbol/Expr, str
The Young’s modulus
The Poisson’s ratio
Lamé’s first parameter.
Lamé’s second parameter (shear modulus).
Mass density.
>>> plane_stress_equations = LinearElasticityPlaneStress(E=10, nu=0.3) >>> plane_stress_equations.pprint() stress_disp_xx: -sigma_xx + 10.989010989011*u__x + 3.2967032967033*v__y stress_disp_yy: -sigma_yy + 3.2967032967033*u__x + 10.989010989011*v__y stress_disp_xy: -sigma_xy + 3.84615384615385*u__y + 3.84615384615385*v__x equilibrium_x: -sigma_xx__x - sigma_xy__y equilibrium_y: -sigma_xy__x - sigma_yy__y traction_x: normal_x*sigma_xx + normal_y*sigma_xy traction_y: normal_x*sigma_xy + normal_y*sigma_yy
Equations related to Navier Stokes Equations
Bases:
modulus.pdes.PDES
del cross vector operator
- vectortuple of 3 Sympy Exprs, floats, ints or strings
- curl_nametuple of 3 strings
This will be the vector to take the curl of.
These will be the output names of the curl operations.
>>> c = Curl((0,0,'phi'), ('u','v','w')) >>> c.pprint() u: phi__y v: -phi__x w: 0
Bases:
modulus.pdes.PDES
Implementation of the gradient boundary condition
- Tstr
- dimint
- timebool
The dependent variable.
Dimension of the equations (1, 2, or 3). Default is 3.
If time-dependent equations or not. Default is True.
>>> gn = ns = GradNormal(T='T') >>> gn.pprint() normal_gradient_T: normal_x*T__x + normal_y*T__y + normal_z*T__z
Bases:
modulus.pdes.PDES
Compressible Navier Stokes equations
- nufloat, Sympy Symbol/Expr, str
- rhofloat, Sympy Symbol/Expr, str
- dimint
- timebool
The kinematic viscosity. If nu is a str then it is converted to Sympy Function of form nu(x,y,z,t). If nu is a Sympy Symbol or Expression then this is substituted into the equation. This allows for variable viscosity.
The density of the fluid. If rho is a str then it is converted to Sympy Function of form ‘rho(x,y,z,t)’. If ‘rho’ is a Sympy Symbol or Expression then this is substituted into the equation to allow for compressible Navier Stokes. Default is 1.
Dimension of the Navier Stokes (2 or 3). Default is 3.
If time-dependent equations or not. Default is True.
>>> ns = NavierStokes(nu=0.01, rho=1, dim=2) >>> ns.pprint() continuity: u__x + v__y momentum_x: u*u__x + v*u__y + p__x + u__t - 0.01*u__x__x - 0.01*u__y__y momentum_y: u*v__x + v*v__y + p__y + v__t - 0.01*v__x__x - 0.01*v__y__y >>> ns = NavierStokes(nu='nu', rho=1, dim=2, time=False) >>> ns.pprint() continuity: u__x + v__y momentum_x: -nu*u__x__x - nu*u__y__y + u*u__x + v*u__y - nu__x*u__x - nu__y*u__y + p__x momentum_y: -nu*v__x__x - nu*v__y__y + u*v__x + v*v__y - nu__x*v__x - nu__y*v__y + p__y
Screened Poisson Distance Equation taken from, https://www.researchgate.net/publication/266149392_Dynamic_Distance-Based_Shape_Features_for_Gait_Recognition, Equation 6 in paper.
- class modulus.PDES.signed_distance_function.ScreenedPoissonDistance(distance='normal_distance', tau=0.1, dim=3)
Bases:
modulus.pdes.PDES
Screened Poisson Distance
- distancestr
- taufloat
- dimint
A user-defined variable for distance. Default is “normal_distance”.
A small, positive parameter. Default is 0.1.
Dimension of the Screened Poisson Distance (1, 2, or 3). Default is 3.
>>> s = ScreenedPoissonDistance(tau=0.1, dim=2) >>> s.pprint() screened_poisson_normal_distance: -normal_distance__x**2 + 0.316227766016838*normal_distance__x__x - normal_distance__y**2 + 0.316227766016838*normal_distance__y__y + 1
Zero Equation Turbulence model References: https://www.eureka.im/954.html https://knowledge.autodesk.com/support/cfd/learn-explore/caas/CloudHelp/cloudhelp/2019/ENU/SimCFD-Learning/files/GUID-BBA4E008-8346-465B-9FD3-D193CF108AF0-htm.html
- class modulus.PDES.turbulence_zero_eq.ZeroEquation(nu, max_distance, rho=1, dim=3, time=True)
Bases:
modulus.pdes.PDES
Zero Equation Turbulence model
- nufloat
- max_distancefloat
- rhofloat, Sympy Symbol/Expr, str
- dimint
- timebool
The kinematic viscosity of the fluid.
The maximum wall distance in the flow field.
The density. If rho is a str then it is converted to Sympy Function of form ‘rho(x,y,z,t)’. If ‘rho’ is a Sympy Symbol or Expression then this is substituted into the equation. Default is 1.
Dimension of the Zero Equation Turbulence model (2 or 3). Default is 3.
If time-dependent equations or not. Default is True.
>>> zeroEq = ZeroEquation(nu=0.1, max_distance=2.0, dim=2) >>> kEp.pprint() nu: sqrt((u__y + v__x)**2 + 2*u__x**2 + 2*v__y**2) *Min(0.18, 0.419*normal_distance)**2 + 0.1
Wave equation Reference: https://en.wikipedia.org/wiki/Wave_equation
- class modulus.PDES.wave_equation.HelmholtzEquation(u, k, dim=3)
Bases: modulus.pdes.PDES
- class modulus.PDES.wave_equation.WaveEquation(u='u', c='c', dim=3, time=True)
Bases:
modulus.pdes.PDES
Wave equation
- ustr
- cfloat, Sympy Symbol/Expr, str
- dimint
- timebool
The dependent variable.
Wave speed coefficient. If c is a str then it is converted to Sympy Function of form ‘c(x,y,z,t)’. If ‘c’ is a Sympy Symbol or Expression then this is substituted into the equation.
Dimension of the wave equation (1, 2, or 3). Default is 2.
If time-dependent equations or not. Default is True.
>>> we = WaveEquation(c=0.8, dim=3) >>> we.pprint() wave_equation: u__t__t - 0.64*u__x__x - 0.64*u__y__y - 0.64*u__z__z >>> we = WaveEquation(c='c', dim=2, time=False) >>> we.pprint() wave_equation: -c**2*u__x__x - c**2*u__y__y - 2*c*c__x*u__x - 2*c*c__y*u__y