proteus.mprans.RANS2P module
Optimized Two-Phase Reynolds Averaged Navier-Stokes
- class proteus.mprans.RANS2P.SubgridError(coefficients, nd, lag=False, nStepsToDelay=1, hFactor=1.0, noPressureStabilization=False)[source]
Bases:
proteus.SubgridError.SGE_base
Create a SubgridError object for two-phase incompressible flow
The VMS subgrid error approximation
- Parameters
coefficients (proteus.TransportCoefficients.TC_base) – The coefficients object
nd (int) – Number of space dimensions
lag (bool) – Use prior time step to calculate
nStepsToDelay (int) – Lag only after nSteps
hFactor (float) – scaling factor based on order
noPressureStabilization (bool) – turn off pressure stab
- initializeElementQuadrature(mesh, t, cq)[source]
Allocated or set additional arrays for values at element quadrature
- Parameters
mesh (proteus.MeshTools.Mesh) – The mesh for the domain
t (float) – The current time
cq (dict) – The dictionary of element quadrature arrrays
- class proteus.mprans.RANS2P.NumericalFlux(vt, getPointwiseBoundaryConditions, getAdvectiveFluxBoundaryConditions, getDiffusiveFluxBoundaryConditions, getPeriodicBoundaryConditions=None)[source]
Bases:
proteus.NumericalFlux.NavierStokes_Advection_DiagonalUpwind_Diffusion_SIPG_exterior
- class proteus.mprans.RANS2P.ShockCapturing(coefficients, nd, shockCapturingFactor=0.25, lag=False, nStepsToDelay=1)[source]
- class proteus.mprans.RANS2P.Coefficients(epsFact=1.5, sigma=72.8, rho_0=998.2, nu_0=1.004e-06, rho_1=1.205, nu_1=1.5e-05, g=[0.0, 0.0, - 9.8], nd=3, ME_model=0, CLSVOF_model=None, LS_model=None, VF_model=None, KN_model=None, Closure_0_model=None, Closure_1_model=None, epsFact_density=None, stokes=False, sd=True, movingDomain=False, useVF=0.0, useRBLES=0.0, useMetrics=0.0, useConstant_he=False, dragAlpha=0.0, dragBeta=0.0, setParamsFunc=None, dragAlphaTypes=None, dragBetaTypes=None, porosityTypes=None, killNonlinearDrag=False, epsFact_solid=None, epsFact_porous=None, eb_adjoint_sigma=1.0, eb_penalty_constant=100.0, forceStrongDirichlet=False, turbulenceClosureModel=0, smagorinskyConstant=0.1, barycenters=None, NONCONSERVATIVE_FORM=1.0, MOMENTUM_SGE=1.0, PRESSURE_SGE=1.0, VELOCITY_SGE=1.0, PRESSURE_PROJECTION_STABILIZATION=0.0, phaseFunction=None, LAG_LES=1.0, use_ball_as_particle=1, ball_center=None, ball_mass=None, ball_stiffness=None, ball_force_range=None, particle_cfl=0.001, particle_box=[1.0, 1.0, 1.0], ball_radius=None, ball_velocity=None, ball_angular_velocity=None, ball_center_acceleration=None, ball_angular_acceleration=None, ball_density=None, particle_velocities=None, particle_centroids=None, particle_sdfList=None, particle_velocityList=None, nParticles=0, particle_epsFact=3.0, particle_alpha=1000.0, particle_beta=1000.0, particle_penalty_constant=100.0, ghost_penalty_constant=0.1, particle_nitsche=1.0, nullSpace='NoNullSpace', useExact=False, analyticalSolution=None, initialize=True, force_x=None, force_y=None, force_z=None, normalize_pressure=False, useInternalParticleSolver=False)[source]
Bases:
proteus.TransportCoefficients.TC_base
The coefficients for two incompresslble fluids governed by the Navier-Stokes equations and separated by a sharp interface represented by a level set function
Notes
The PRESSURE_PROJECTION_STABILIZATION flag allows the user to use the Bochev-Dohrmann-Gunzburger stabilization introduced in Stabilization of Low-Order Mixed Finite Elements for the Stokes Equation. This option should be turned off for most problems, but in some instances it may produced better preconditioning results than the full SGE approach.
Set the number of components (equations) of the PDE and initialize the dicitionaries describing the form of the coefficients. Strings naming each component (used for viewing and archiving) and a structure defining the sparsity pattern of diffusion tensors may also be provided.
- attachModels(modelList)[source]
Give the TC object access to other models in a loosely coupled split operator formulation (e.g. a transport equation for concentration might get velocity from a flow equation)
- initializeMesh(mesh)[source]
Give the TC object access to the mesh for any mesh-dependent information.
- initializeElementQuadrature(t, cq)[source]
Give the TC object access to the element quadrature storage
- initializeElementBoundaryQuadrature(t, cebq, cebq_global)[source]
Give the TC object access to the element boundary quadrature storage
- initializeGlobalExteriorElementBoundaryQuadrature(t, cebqe)[source]
Give the TC object access to the exterior element boundary quadrature storage
- evaluate(t, c)[source]
Evaluate the coefficients at a given time, t, using the coefficient storage passed in as the dictionary c.
- class proteus.mprans.RANS2P.LevelModel(uDict, phiDict, testSpaceDict, matType, dofBoundaryConditionsDict, dofBoundaryConditionsSetterDict, coefficients, elementQuadrature, elementBoundaryQuadrature, fluxBoundaryConditionsDict=None, advectiveFluxBoundaryConditionsSetterDict=None, diffusiveFluxBoundaryConditionsSetterDictDict=None, stressTraceBoundaryConditionsSetterDictDict=None, stabilization=None, shockCapturing=None, conservativeFluxDict=None, numericalFluxType=None, TimeIntegrationClass=None, massLumping=False, reactionLumping=False, options=None, name='RANS2P', reuse_trial_and_test_quadrature=False, sd=True, movingDomain=False)[source]
Bases:
proteus.Transport.OneLevelTransport
Allocate storage and initialize some variables.
- Parameters
uDict (dict) – Dictionary of
proteus.FemTools.FiniteElementFunction
objects.phiDict (dict) – Dictionary of
proteus.FemTools.FiniteElementFunction
objects.testSpaceDict (dict) – Dictionary of FiniteElementSpace objects
dofBoundaryConditionsDict (dict) – Dictionary of DOFBoundaryConditions objects for the Dirichlet conditions.
coefficients (
proteus.TransportCoefficients.TC_base
) – Problem’s Transport Coefficients class.elementQuadratureDict (dict) – Dictionary of dictionaries of quadrature rules for each element integral in each component equation.
elementBoundaryQuadratureDict (dict) – Dictionary of dictionaries of quadrature rules for each element boundary integral in each component equation
stabilization (bool) –
shockCapturing (bool) –
numericalFlux (bool) –
bdyNullSpace (bool) – Indicates whether the boundary conditions create a global null space.
Notes
The constructor sets the input arguments, calculates dimensions, and allocates storage. The meanings of variable suffixes are
global – per physical domain
element – per element
elementBoundary – per element boundary
The prefix n means ‘number of’.
Storage is divided into quantities required at different sets of points or geometric entities. Each type of storage has a dictionary for all the quantities of that type. The names and dimensions of the storage dictionaries are
e – at element
q – at element quadrature, unique to elements
ebq – at element boundary quadrature, unique to elements
ebq_global – at element boundary quadrature, unique to element boundary
ebqe – at element boundary quadrature, unique to global, exterior element boundary
phi_ip – at the generalized interpolation points required to build a nonlinear phi
- getResidual(u, r)[source]
Calculate the element residuals and add in to the global residual
- Parameters
u (
numpy.ndarray
) –r (
numpy.ndarray
) – Stores the calculated residual vector.
- calculateElementQuadrature(domainMoved=False)[source]
Calculate the physical location and weights of the quadrature rules and the shape information at the quadrature points.
This function should be called only when the mesh changes.
- calculateElementBoundaryQuadrature(domainMoved=False)[source]
Calculate the physical location and weights of the quadrature rules and the shape information at the quadrature points on element boundaries.
This function should be called only when the mesh changes.