NVIDIA Modulus Sym (Latest Release)
Sym (Latest Release)

Modulus Sym Equations

Advection diffusion equation Reference: https://en.wikipedia.org/wiki/Convection%E2%80%93diffusion_equation

class modulus.sym.eq.pdes.advection_diffusion.AdvectionDiffusion(T='T', D='D', Q=0, rho='rho', dim=3, time=False, mixed_form=False)[source]

Bases: PDE

Advection diffusion equation

Parameters
  • T (str) – The dependent variable.

  • D (float, Sympy Symbol/Expr, str) – 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.

  • Q (float, Sympy Symbol/Expr, str) – 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.

  • rho (float, Sympy Symbol/Expr, str) – 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.

  • dim (int) – Dimension of the diffusion equation (1, 2, or 3). Default is 3.

  • time (bool) – If time-dependent equations or not. Default is False.

  • mixed_form (bool) – If True, use the mixed formulation of the wave equation.

Examples

Copy
Copied!
            

>>> 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.sym.eq.pdes.basic.Curl(vector, curl_name=['u', 'v', 'w'])[source]

Bases: PDE

del cross vector operator

Parameters
  • vector (tuple of 3 Sympy Exprs, floats, ints or strings) – This will be the vector to take the curl of.

  • curl_name (tuple of 3 strings) – These will be the output names of the curl operations.

Examples

Copy
Copied!
            

>>> c = Curl((0,0,'phi'), ('u','v','w')) >>> c.pprint() u: phi__y v: -phi__x w: 0

class modulus.sym.eq.pdes.basic.GradNormal(T, dim=3, time=True)[source]

Bases: PDE

Implementation of the gradient boundary condition

Parameters
  • T (str) – The dependent variable.

  • dim (int) – Dimension of the equations (1, 2, or 3). Default is 3.

  • time (bool) – If time-dependent equations or not. Default is True.

Examples

Copy
Copied!
            

>>> gn = ns = GradNormal(T='T') >>> gn.pprint() normal_gradient_T: normal_x*T__x + normal_y*T__y + normal_z*T__z

class modulus.sym.eq.pdes.basic.NormalDotVec(vec=['u', 'v', 'w'])[source]

Bases: PDE

Normal dot velocity

Parameters

dim (int) – Dimension of the equations (1, 2, or 3). Default is 3.

Diffusion equation

class modulus.sym.eq.pdes.diffusion.Diffusion(T='T', D='D', Q=0, dim=3, time=True, mixed_form=False)[source]

Bases: PDE

Diffusion equation

Parameters
  • T (str) – The dependent variable.

  • D (float, Sympy Symbol/Expr, str) – 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.

  • Q (float, Sympy Symbol/Expr, str) – 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.

  • dim (int) – Dimension of the diffusion equation (1, 2, or 3). Default is 3.

  • time (bool) – If time-dependent equations or not. Default is True.

  • mixed_form (bool) – If True, use the mixed formulation of the diffusion equations.

Examples

Copy
Copied!
            

>>> 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.sym.eq.pdes.diffusion.DiffusionInterface(T_1, T_2, D_1, D_2, dim=3, time=True)[source]

Bases: PDE

Matches the boundary conditions at an interface

Parameters
  • T_1 (str) – Dependent variables to match the boundary conditions at the interface.

  • T_2 (str) – Dependent variables to match the boundary conditions at the interface.

  • D_1 (float) – Diffusivity at the interface.

  • D_2 (float) – Diffusivity at the interface.

  • dim (int) – Dimension of the equations (1, 2, or 3). Default is 3.

  • time (bool) – If time-dependent equations or not. Default is True.

Example

Copy
Copied!
            

>>> 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.sym.eq.pdes.electromagnetic.MaxwellFreqReal(ux='ux', uy='uy', uz='uz', k=1.0, mixed_form=False)[source]

Bases: PDE

Frequency domain Maxwell’s equation

Parameters
  • ux (str) – Ex

  • uy (str) – Ey

  • uz (str) – Ez

  • k (float, Sympy Symbol/Expr, str) – 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.

  • mixed_form (bool) – If True, use the mixed formulation of the diffusion equations.

class modulus.sym.eq.pdes.electromagnetic.PEC(ux='ux', uy='uy', uz='uz', dim=3)[source]

Bases: PDE

Perfect Electric Conduct BC for

Parameters
  • ux (str) – Ex

  • uy (str) – Ey

  • uz (str) – Ez

  • dim (int) – Dimension of the equations (2, or 3). Default is 3.

class modulus.sym.eq.pdes.electromagnetic.SommerfeldBC(ux='ux', uy='uy', uz='uz')[source]

Bases: PDE

Frequency domain ABC, Sommerfeld radiation condition Only for real part Equation: ‘n x _curl(E) = 0’

Parameters
  • ux (str) – Ex

  • uy (str) – Ey

  • uz (str) – Ez

Equations related to linear elasticity

class modulus.sym.eq.pdes.linear_elasticity.LinearElasticity(E=None, nu=None, lambda_=None, mu=None, rho=1, dim=3, time=False)[source]

Bases: PDE

Linear elasticity equations. Use either (E, nu) or (lambda_, mu) to define the material properties.

Parameters
  • E (float, Sympy Symbol/Expr, str) – The Young’s modulus

  • nu (float, Sympy Symbol/Expr, str) – The Poisson’s ratio

  • lambda (float, Sympy Symbol/Expr, str) – Lamé’s first parameter

  • mu (float, Sympy Symbol/Expr, str) – Lamé’s second parameter (shear modulus)

  • rho (float, Sympy Symbol/Expr, str) – Mass density.

  • dim (int) – Dimension of the linear elasticity (2 or 3). Default is 3.

Example

Copy
Copied!
            

>>> 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.sym.eq.pdes.linear_elasticity.LinearElasticityPlaneStress(E=None, nu=None, lambda_=None, mu=None, rho=1, time=False)[source]

Bases: PDE

Linear elasticity plane stress equations. Use either (E, nu) or (lambda_, mu) to define the material properties.

Parameters
  • E (float, Sympy Symbol/Expr, str) – The Young’s modulus

  • nu (float, Sympy Symbol/Expr, str) – The Poisson’s ratio

  • lambda (float, Sympy Symbol/Expr, str) – Lamé’s first parameter.

  • mu (float, Sympy Symbol/Expr, str) – Lamé’s second parameter (shear modulus)

  • rho (float, Sympy Symbol/Expr, str) – Mass density.

Example

Copy
Copied!
            

>>> 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

class modulus.sym.eq.pdes.navier_stokes.CompressibleIntegralContinuity(rho=1, vec=['u', 'v', 'w'])[source]

Bases: PDE

Compressible Integral Continuity

Parameters
  • rho (float, Sympy Symbol/Expr, str) – 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 compressibility. Default is 1.

  • dim (int) – Dimension of the equations (1, 2, or 3). Default is 3.

class modulus.sym.eq.pdes.navier_stokes.Curl(vector, curl_name=['u', 'v', 'w'])[source]

Bases: PDE

del cross vector operator

Parameters
  • vector (tuple of 3 Sympy Exprs, floats, ints or strings) – This will be the vector to take the curl of.

  • curl_name (tuple of 3 strings) – These will be the output names of the curl operations.

Examples

Copy
Copied!
            

>>> c = Curl((0,0,'phi'), ('u','v','w')) >>> c.pprint() u: phi__y v: -phi__x w: 0

class modulus.sym.eq.pdes.navier_stokes.FluxContinuity(T='T', D='D', rho=1, vec=['u', 'v', 'w'])[source]

Bases: PDE

Flux Continuity for arbitrary variable. Includes advective and diffusive flux

Parameters
  • T (str) – The dependent variable.

  • rho (float, Sympy Symbol/Expr, str) – 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 compressibility. Default is 1.

  • dim (int) – Dimension of the equations (1, 2, or 3). Default is 3.

class modulus.sym.eq.pdes.navier_stokes.GradNormal(T, dim=3, time=True)[source]

Bases: PDE

Implementation of the gradient boundary condition

Parameters
  • T (str) – The dependent variable.

  • dim (int) – Dimension of the equations (1, 2, or 3). Default is 3.

  • time (bool) – If time-dependent equations or not. Default is True.

Examples

Copy
Copied!
            

>>> gn = ns = GradNormal(T='T') >>> gn.pprint() normal_gradient_T: normal_x*T__x + normal_y*T__y + normal_z*T__z

class modulus.sym.eq.pdes.navier_stokes.NavierStokes(nu, rho=1, dim=3, time=True, mixed_form=False)[source]

Bases: PDE

Compressible Navier Stokes equations Reference: https://turbmodels.larc.nasa.gov/implementrans.html

Parameters
  • nu (float, Sympy Symbol/Expr, str) – 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.

  • rho (float, Sympy Symbol/Expr, str) – 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.

  • dim (int) – Dimension of the Navier Stokes (2 or 3). Default is 3.

  • time (bool) – If time-dependent equations or not. Default is True.

  • mixed_form (bool) – If True, use the mixed formulation of the Navier-Stokes equations.

Examples

Copy
Copied!
            

>>> 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 - 2*nu__x*u__x - nu__y*u__y - nu__y*v__x + p__x momentum_y: -nu*v__x__x - nu*v__y__y + u*v__x + v*v__y - nu__x*u__y - nu__x*v__x - 2*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.sym.eq.pdes.signed_distance_function.ScreenedPoissonDistance(distance='normal_distance', tau=0.1, dim=3)[source]

Bases: PDE

Screened Poisson Distance

Parameters
  • distance (str) – A user-defined variable for distance. Default is “normal_distance”.

  • tau (float) – A small, positive parameter. Default is 0.1.

  • dim (int) – Dimension of the Screened Poisson Distance (1, 2, or 3). Default is 3.

Example

Copy
Copied!
            

>>> 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.sym.eq.pdes.turbulence_zero_eq.ZeroEquation(nu, max_distance, rho=1, dim=3, time=True)[source]

Bases: PDE

Zero Equation Turbulence model

Parameters
  • nu (float) – The kinematic viscosity of the fluid.

  • max_distance (float) – The maximum wall distance in the flow field.

  • rho (float, Sympy Symbol/Expr, str) – 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.

  • dim (int) – Dimension of the Zero Equation Turbulence model (2 or 3). Default is 3.

  • time (bool) – If time-dependent equations or not. Default is True.

Example

Copy
Copied!
            

>>> 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.sym.eq.pdes.wave_equation.HelmholtzEquation(u, k, dim=3, mixed_form=False)[source]

Bases: PDE

class modulus.sym.eq.pdes.wave_equation.WaveEquation(u='u', c='c', dim=3, time=True, mixed_form=False)[source]

Bases: PDE

Wave equation

Parameters
  • u (str) – The dependent variable.

  • c (float, Sympy Symbol/Expr, str) – 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.

  • dim (int) – Dimension of the wave equation (1, 2, or 3). Default is 2.

  • time (bool) – If time-dependent equations or not. Default is True.

  • mixed_form (bool) – If True, use the mixed formulation of the wave equation.

Examples

Copy
Copied!
            

>>> 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

class modulus.sym.eq.derivatives.Derivative(bwd_derivative_dict: Dict[Key, List[Key]])[source]

Bases: Module

Module to compute derivatives using backward automatic differentiation

forward(input_var: Dict[str, Tensor]) → Dict[str, Tensor][source]

Define the computation performed at every call.

Should be overridden by all subclasses.

Note

Although the recipe for forward pass needs to be defined within this function, one should call the Module instance afterwards instead of this since the former takes care of running the registered hooks while the latter silently ignores them.

class modulus.sym.eq.derivatives.MeshlessFiniteDerivative(model: Module, derivatives: List[Key], dx: Union[float, Callable], order: int = 2, max_batch_size: Optional[int] = None, double_cast: bool = True, input_keys: Optional[List[Key]] = None)[source]

Bases: Module

Module to compute derivatives using meshless finite difference

Parameters
  • model (torch.nn.Module) – Forward torch module for calculating stencil values

  • derivatives (List[Key]) – List of derivative keys to calculate

  • dx (Union[float, Callable]) – Spatial discretization of all axis, can be function with parameter count which is the number of forward passes for dynamically adjusting dx

  • order (int, optional) – Order of derivative, by default 2

  • max_batch_size (Union[int, None], optional) – Max batch size of stencil calucations, by default uses batch size of inputs

  • double_cast (bool, optional) – Cast fields to double precision to calculate derivatives, by default True

  • jit (bool, optional) – Use torch script for finite deriv calcs, by default True

forward(inputs: Dict[str, Tensor]) → Dict[str, Tensor][source]

Define the computation performed at every call.

Should be overridden by all subclasses.

Note

Although the recipe for forward pass needs to be defined within this function, one should call the Module instance afterwards instead of this since the former takes care of running the registered hooks while the latter silently ignores them.

classmethod make_node(node_model: Union[Node, Module], derivatives: List[Key], dx: Union[float, Callable], order: int = 2, max_batch_size: Optional[int] = None, name: Optional[str] = None, double_cast: bool = True, input_keys: Optional[Union[List[Key], List[str]]] = None)[source]

Makes a meshless finite derivative node.

Parameters
  • node_model (Union[Node, torch.nn.Module]) – Node or torch.nn.Module for computing FD stencil values. Part of the inputs to this model should consist of the independent variables and output the functional value

  • derivatives (List[Key]) – List of derivatives to be computed

  • dx (Union[float, Callable]) – Spatial discretization for finite diff calcs, can be function

  • order (int, optional) – Order of accuracy of finite diff calcs, by default 2

  • max_batch_size (Union[int, None], optional) – Maximum batch size to used with the stenicl foward passes, by default None

  • name (str, optional) – Name of node, by default None

  • double_cast (bool, optional) – Cast tensors to double precision for derivatives, by default True

  • input_keys (Union[List[Key], List[str], None], optional) – List of input keys to be used for input of forward model. Should be used if node_model is not a Node, by default None

Previous Modulus Sym Geometry
Next Modulus Sym Models
© Copyright 2023, NVIDIA Modulus Team. Last updated on Jul 25, 2024.