proteus.mprans.MCorr module

class proteus.mprans.MCorr.Coefficients(applyCorrection=True, epsFactHeaviside=0.0, epsFactDirac=1.0, epsFactDiffusion=2.0, LSModel_index=3, V_model=2, me_model=5, VOFModel_index=4, checkMass=True, sd=True, nd=None, applyCorrectionToDOF=True, useMetrics=0.0, useConstantH=False, useQuadraticRegularization=False, edgeBasedStabilizationMethods=False, nullSpace='NoNullSpace', useExact=False, initialize=True)[source]

Bases: proteus.TransportCoefficients.TC_base

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.

levelSetConservationCoefficientsEvaluate()[source]
levelSetConservationCoefficientsEvaluate_sd()[source]
initialize()[source]
initializeMesh(mesh)[source]

Give the TC object access to the mesh for any mesh-dependent information.

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)

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

preStep(t, firstStep=False)[source]

Give the TC object an opportunity to modify itself before the time step.

postStep(t, firstStep=False)[source]

Give the TC object an opportunity to modify itself after the time step.

postAdaptStep()[source]
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.MCorr.LevelModel(uDict, phiDict, testSpaceDict, matType, dofBoundaryConditionsDict, dofBoundaryConditionsSetterDict, coefficients, elementQuadrature, elementBoundaryQuadrature, fluxBoundaryConditionsDict=None, advectiveFluxBoundaryConditionsSetterDict=None, diffusiveFluxBoundaryConditionsSetterDictDict=None, stressTraceBoundaryConditionsSetterDict=None, stabilization=None, shockCapturing=None, conservativeFluxDict=None, numericalFluxType=None, TimeIntegrationClass=None, massLumping=False, reactionLumping=False, options=None, name='defaultName', reuse_trial_and_test_quadrature=True, sd=True, movingDomain=False, bdyNullSpace=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

nCalls = 0[source]
FCTStep()[source]
calculateCoefficients()[source]
calculateElementResidual()[source]

Calculate all the element residuals

getResidual(u, r)[source]

Calculate the element residuals and add in to the global residual

getMassMatrix()[source]
getJacobian(jacobian)[source]
elementSolve(u, r)[source]
elementConstantSolve(u, r)[source]
globalConstantRJ(u, r, U)[source]
globalConstantSolve(u, r)[source]
calculateElementQuadrature()[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()[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.

calculateExteriorElementBoundaryQuadrature()[source]

Calculate the physical location and weights of the quadrature rules and the shape information at the quadrature points on global element boundaries.

This function should be called only when the mesh changes.

estimate_mt()[source]
calculateAuxiliaryQuantitiesAfterStep()[source]
calculateMass(q_phi)[source]
setMassQuadratureEdgeBasedStabilizationMethods()[source]
setMassQuadrature()[source]
calculateSolutionAtQuadrature()[source]
updateAfterMeshMotion()[source]
class proteus.mprans.MCorr.DummyNewton(linearSolver, F, J=None, du=None, par_du=None, rtol_r=0.0001, atol_r=1e-16, rtol_du=0.0001, atol_du=1e-16, maxIts=100, norm=<function l2Norm>, convergenceTest='r', computeRates=True, printInfo=True, fullNewton=True, directSolver=False, EWtol=True, maxLSits=100)[source]

Bases: proteus.NonlinearSolvers.NonlinearSolver

info()[source]
solve(u, r=None, b=None, par_u=None, par_r=None)[source]
class proteus.mprans.MCorr.ElementNewton(linearSolver, F, J=None, du=None, par_du=None, rtol_r=0.0001, atol_r=1e-16, rtol_du=0.0001, atol_du=1e-16, maxIts=100, norm=<function l2Norm>, convergenceTest='r', computeRates=True, printInfo=True, fullNewton=True, directSolver=False, EWtol=True, maxLSits=100)[source]

Bases: proteus.NonlinearSolvers.NonlinearSolver

info()[source]
solve(u, r=None, b=None, par_u=None, par_r=None)[source]
class proteus.mprans.MCorr.ElementConstantNewton(linearSolver, F, J=None, du=None, par_du=None, rtol_r=0.0001, atol_r=1e-16, rtol_du=0.0001, atol_du=1e-16, maxIts=100, norm=<function l2Norm>, convergenceTest='r', computeRates=True, printInfo=True, fullNewton=True, directSolver=False, EWtol=True, maxLSits=100)[source]

Bases: proteus.NonlinearSolvers.NonlinearSolver

info()[source]
solve(u, r=None, b=None, par_u=None, par_r=None)[source]
class proteus.mprans.MCorr.GlobalConstantNewton(linearSolver, F, J=None, du=None, par_du=None, rtol_r=0.0001, atol_r=1e-16, rtol_du=0.0001, atol_du=1e-16, maxIts=100, norm=<function l2Norm>, convergenceTest='r', computeRates=True, printInfo=True, fullNewton=True, directSolver=False, EWtol=True, maxLSits=100)[source]

Bases: proteus.NonlinearSolvers.NonlinearSolver

info()[source]
solve(u, r=None, b=None, par_u=None, par_r=None)[source]
proteus.mprans.MCorr.conservationNorm(x)[source]
class proteus.mprans.MCorr.Newton_controller(model, nOptions)[source]

Bases: proteus.StepControl.Newton_controller

initializeTimeHistory()[source]