xfem.xfem module
- class xfem.xfem.BitArrayCF
Bases:
CoefficientFunctionCoefficientFunction that evaluates a BitArray. On elements with an index i where the BitArray evaluates to true the CoefficientFunction will evaluate as 1, otherwise as 0.
Similar functionality (also for facets) can be obtained with IndicatorCF.
- class xfem.xfem.COMBINED_DOMAIN_TYPE
Bases:
pybind11_objectMembers:
NO
CDOM_NEG
CDOM_POS
UNCUT
CDOM_IF
HASNEG
HASPOS
ANY
- ANY = <COMBINED_DOMAIN_TYPE.ANY: 7>
- CDOM_IF = <COMBINED_DOMAIN_TYPE.CDOM_IF: 4>
- CDOM_NEG = <COMBINED_DOMAIN_TYPE.CDOM_NEG: 1>
- CDOM_POS = <COMBINED_DOMAIN_TYPE.CDOM_POS: 2>
- HASNEG = <COMBINED_DOMAIN_TYPE.HASNEG: 5>
- HASPOS = <COMBINED_DOMAIN_TYPE.HASPOS: 6>
- NO = <COMBINED_DOMAIN_TYPE.NO: 0>
- UNCUT = <COMBINED_DOMAIN_TYPE.UNCUT: 3>
- property name
- property value
- class xfem.xfem.CSpaceTimeFESpace
Bases:
FESpace- IsTimeNodeActive(self: xfem.xfem.CSpaceTimeFESpace, arg0: int) bool
Return bool whether node is active
- SetOverrideTime(self: xfem.xfem.CSpaceTimeFESpace, arg0: bool) None
Set flag to or not to override the time variable
- SetTime(self: xfem.xfem.CSpaceTimeFESpace, arg0: float) None
- Set the time variable
Also sets override time
- TimeFE_nodes(self: xfem.xfem.CSpaceTimeFESpace) list
Return nodes of time FE
- k_t(self: xfem.xfem.CSpaceTimeFESpace) int
Return order of the time FE
- property spaceFES
get space FESpace
- class xfem.xfem.CXFESpace
Bases:
FESpaceXFESpace-class [For documentation of the XFESpace-constructor see help(XFESpace)]:
Extended finite element space. Takes a basis FESpace and creates an enrichment space based on cut information. The cut information is provided by a CutInfo object or - if a level set function is only provided - a CutInfo object is created. The enrichment doubles the unknowns on all cut elements and assigns to them a sign (NEG/POS). One of the differential operators neg(…) or pos(…) evaluates like the basis function of the origin space, the other one as zero for every basis function. Away from cut elements no basis function is supported.
- BaseDofOfXDof(self: xfem.xfem.CXFESpace, arg0: int) int
To an unknown of the extended space, get the corresponding unknown of the base FESpace.
Parameters
- iint
degree of freedom
- GetCutInfo(self: xfem.xfem.CXFESpace) xfem.xfem.CutInfo
Get Information of cut geometry
- GetDomainNrs(self: xfem.xfem.CXFESpace, arg0: int) ngcore::Array<DOMAIN_TYPE, unsigned long>
Get Array of Domains (Array of NEG/POS) of degrees of freedom of the extended FESpace on one element.
Parameters
- elnrint
element number
- GetDomainOfDof(self: xfem.xfem.CXFESpace, arg0: int) xfem.xfem.DOMAIN_TYPE
Get Domain (NEG/POS) of a degree of freedom of the extended FESpace.
Parameters
- iint
degree of freedom
- xfem.xfem.CompoundBitArray(balist: list) pyngcore.pyngcore.BitArray
Takes a list of BitArrays and merges them to one larger BitArray. Can be useful for CompoundFESpaces.
- class xfem.xfem.CompoundProlongation
Bases:
Prolongationprolongation for compound spaces
- AddProlongation(self: xfem.xfem.CompoundProlongation, p1prol: ngsolve.comp.Prolongation) None
- Prolongate(self: xfem.xfem.CompoundProlongation, finelevel: int, vec: ngsolve.la.BaseVector) None
- Restrict(self: xfem.xfem.CompoundProlongation, finelevel: int, vec: ngsolve.la.BaseVector) None
- Update(self: xfem.xfem.CompoundProlongation, fespace: ngsolve.comp.FESpace) None
- xfem.xfem.CreateTimeRestrictedGF(gf: ngsolve.comp.GridFunction, reference_time: float = 0.0) ngsolve.comp.GridFunction
Create spatial-only Gridfunction corresponding to a fixed time.
- class xfem.xfem.CutDifferentialSymbol
Bases:
DifferentialSymbolCutDifferentialSymbol that allows to formulate linear, bilinear forms and integrals on level set domains in an intuitive form:
Example use case:
dCut = CutDifferentialSymbol(VOL) dx = dCut(lset,NEG) a = BilinearForm(…) a += u * v * dx
Note that the most important options are set in the second line when the basic CutDifferentialSymbol is further specified.
- order(self: xfem.xfem.CutDifferentialSymbol, order: int) xfem.xfem.CutDifferentialSymbol
- property vb
Volume of boundary?
- class xfem.xfem.CutInfo
Bases:
pybind11_objectA CutInfo stores and organizes cut informations in the mesh with respect to a level set function. Elements (BND / VOL) and facets can be either cut elements or in the positive (POS) or negative (NEG) part of the domain. A CutInfo provides information about the cut configuration in terms of BitArrays and Vectors of Ratios. (Internally also domain_types for different mesh nodes are stored.)
- GetCutRatios(self: xfem.xfem.CutInfo, VOL_or_BND: ngsolve.comp.VorB = <VorB.VOL: 0>) ngsolve.la.BaseVector
Returns Vector of the ratios between the measure of the NEG domain on a (boundary) element and the full (boundary) element
- GetElementsOfType(self: xfem.xfem.CutInfo, domain_type: object = <DOMAIN_TYPE.IF: 2>, VOL_or_BND: ngsolve.comp.VorB = <VorB.VOL: 0>) pyngcore.pyngcore.BitArray
Returns BitArray that is true for every element that has the corresponding combined domain type (NO/NEG/POS/UNCUT/IF/HASNEG/HASPOS/ANY)
- GetElementsWithThresholdContribution(self: xfem.xfem.CutInfo, domain_type: object = <DOMAIN_TYPE.NEG: 0>, threshold: float = 1.0, VOL_or_BND: ngsolve.comp.VorB = <VorB.VOL: 0>) pyngcore.pyngcore.BitArray
Returns BitArray marking the elements where the cut ratio is greater or equal to the given threshold.
Parameters
- domain_typeENUM
Check POS or NEG elements.
- thresholdfloat
Mark elements with cut ratio (volume of domain_type / volume background mesh) greater or equal to threshold.
- VOL_or_BNDngsolve.comp.VorB
input VOL, BND, ..
- GetFacetsOfType(self: xfem.xfem.CutInfo, domain_type: object = <DOMAIN_TYPE.IF: 2>) pyngcore.pyngcore.BitArray
Returns BitArray that is true for every facet that has the corresponding combined domain type (NO/NEG/POS/UNCUT/IF/HASNEG/HASPOS/ANY)
- Mesh(self: xfem.xfem.CutInfo) ngsolve.comp.Mesh
Returns mesh of CutInfo
- Update(self: xfem.xfem.CutInfo, levelset: ngsolve.fem.CoefficientFunction, subdivlvl: int = 0, time_order: int = -1, heapsize: int = 1000000) None
Updates a CutInfo based on a level set function.
Parameters
- levelsetngsolve.CoefficientFunction
level set function w.r.t. which the CutInfo is generated
- subdivlvlint
subdivision for numerical integration
- time_orderint
order in time that is used in the integration in time to check for cuts and the ratios. This is only relevant for space-time discretizations.
- class xfem.xfem.DOMAIN_TYPE
Bases:
pybind11_objectMembers:
POS
NEG
IF
- IF = <DOMAIN_TYPE.IF: 2>
- NEG = <DOMAIN_TYPE.NEG: 0>
- POS = <DOMAIN_TYPE.POS: 1>
- property name
- property value
- class xfem.xfem.ElementAggregation
Bases:
pybind11_object- ElementAggregation does the following:
It collects patches of elements that allow to stabilize bad cut elements by at least one good element (the root element).
)
- Update(self: xfem.xfem.ElementAggregation, root_elements: pyngcore.pyngcore.BitArray, bad_elements: pyngcore.pyngcore.BitArray, heapsize: int = 1000000) None
Updates a Element Aggregation based …
- property element_to_patch
vector mapping elements to (non-trivial) patches
- property els_in_nontrivial_patch
BitArray that is true for every element that is part of a (non-trivial) patch
- property els_in_trivial_patch
BitArray that is true for every element that is not part of a (non-trivial) patch
- property facet_to_patch
vector mapping facets to (non-trivial) patches
- property patch_interior_facets
BitArray that is true for every facet that is inside an aggregation cluster
- xfem.xfem.ExtensionEmbedding(elagg: xfem.xfem.ElementAggregation, fes: ngsolve.comp.FESpace, bf: ngsolve.comp.SumOfIntegrals, heapsize: int = 1000000) ngsolve.la.SparseMatrixd
Computes the embedding matrix for an extension minimizing a described energy on patched (followed by averaging if some dofs are shared by multiple patches)
Parameters:
- elaggElementAggregation
ElementAggregation instace defining the patches which are aggregated into a single element
- fesngsolve.FESpace
The finite element space which is aggregated.
- bfngsolve.SumOfIntegrals
The bilinear form describing the energy to be minized
- class xfem.xfem.FacetPatchDifferentialSymbol
Bases:
DifferentialSymbolFacetPatchDifferentialSymbol that allows to formulate integrals on facet patches. Example use case:
dFacetPatch = FacetPatchDifferentialSymbol(VOL) dw = dFacetPatch(definedonelements = …) a = BilinearForm(…) a += (u-u.Other()) * (v-v.Other()) * dw
- class xfem.xfem.GCC3FE
Bases:
ScalarTimeFE
- xfem.xfem.GetDofsOfElements(space: ngsolve.comp.FESpace, a: pyngcore.pyngcore.BitArray, heapsize: int = 1000000) pyngcore.pyngcore.BitArray
Given a BitArray marking some elements in a mesh extract all unknowns that are supported on these elements as a BitArray.
Parameters:
- spacengsolve.FESpace
finite element space from which the corresponding dofs should be extracted
- angsolve.BitArray
BitArray for marked elements
- heapsizeint
heapsize of local computations.
- xfem.xfem.GetDofsOfFacets(space: ngsolve.comp.FESpace, a: pyngcore.pyngcore.BitArray, heapsize: int = 1000000) pyngcore.pyngcore.BitArray
Given a BitArray marking some facets in a mesh extract all unknowns that are associated to these facets as a BitArray.
Parameters:
- spacengsolve.FESpace
finite element space from which the corresponding dofs should be extracted
- angsolve.BitArray
BitArray for marked Facets
- heapsizeint
heapsize of local computations.
- xfem.xfem.GetElementsWithNeighborFacets(mesh: ngsolve.comp.Mesh, a: pyngcore.pyngcore.BitArray, heapsize: int = 1000000) pyngcore.pyngcore.BitArray
Given a BitArray marking some facets extract a BitArray of elements that are neighboring these facets
Parameters:
- mesh
mesh
- angsolve.BitArray
BitArray for marked facets
- heapsizeint
heapsize of local computations.
- xfem.xfem.GetFacetsWithNeighborTypes(mesh: ngsolve.comp.Mesh, a: pyngcore.pyngcore.BitArray, bnd_val_a: bool = True, bnd_val_b: bool = True, use_and: bool = True, b: object = None, heapsize: int = 1000000) pyngcore.pyngcore.BitArray
Given a mesh and two BitArrays (if only one is provided these are set to be equal) facets will be marked (in terms of BitArrays) depending on the BitArray-values on the neighboring elements. The BitArrays are complemented with flags for potential boundary values for the BitArrays. The decision on every facet is now based on the values a and b (left and right) where a or b can also be obtained from the BitArray boundary values. The result is:
- result = (a(left) and b(right))
or (b(left) and a(right))
- or
- result = (a(left) or b(right))
or (b(left) or a(right))
Parameters:
- mesh
mesh
- angsolve.BitArray
first BitArray
- bngsolve.BitArray / None
second BitArray. If None, b=a
- bnd_val_aboolean
BitArray-replacement for a if a(left) or a(right) is not valid (at the boundary)
- bnd_val_aboolean
BitArray-replacement for b if b(left) or b(right) is not valid (at the boundary)
- use_andboolean
use ‘and’-relation to evaluate the result. Otherwise use ‘or’-relation
- heapsizeint
heapsize of local computations.
- class xfem.xfem.GlobalNgsxfemVariables
Bases:
pybind11_objectThe class GlobalNgsxfemVariables provides Python-access to several internal parameters and options used by different subprocedures of ngsxfem. For “mainstream” application cases, it should not be required to change parameters here. Most cases where this class is practically relevant will be debugging or special applications, like investigations in a regime of total error below ~1e-8.
Properties:
- eps_spacetime_lset_perturbationdouble
When handling cut topologies, it is sometimes cumbersome to include the case of a lset value of exactly 0. Hence, the value will be set to eps_spacetime_lset_perturbation in the routine for generating space-time quadrature rules in case its absolute value is smaller. Default: 1e-14
- eps_spacetime_cutrule_bisectiondouble
For high temporal orders, the space-time quadrature rule will apply a bisection method to find those time points with topology changes. This parameters controls how small 2 times the value must be in order to be counted as a root. Default: 1e-15
- eps_P1_perturbationdouble
Similar to eps_spacetime_lset_perturbation, but for the P1 interpolation routine. Default: 1e-14
- eps_spacetime_fes_nodedouble
When a Gridfunction is restricted, the given time point is compared to the nodes of the finite element, such that those node values can be extracted directly in a matching case. This parameters controlls how far a deviation will still be counted as coincidence. Default: 1e-9
- MultiplyAllEps(self: xfem.xfem.GlobalNgsxfemVariables, arg0: float) None
- Output(self: xfem.xfem.GlobalNgsxfemVariables) None
- SetDefaults(self: xfem.xfem.GlobalNgsxfemVariables) None
- SwitchSIMD(self: xfem.xfem.GlobalNgsxfemVariables, arg0: bool) None
- property do_naive_timeint
- property eps_P1_perturbation
- property eps_facetpatch_ips
- property eps_shifted_eval
- property eps_spacetime_cutrule_bisection
- property eps_spacetime_fes_node
- property eps_spacetime_lset_perturbation
- property fixed_point_maxiter_shifted_eval
- property max_dist_newton
- property naive_timeint_order
- property naive_timeint_subdivs
- property newton_maxiter
- property non_conv_warn_msg_lvl
- property simd_eval
- xfem.xfem.IntegrateX(levelset_domain: dict, mesh: ngsolve.comp.Mesh, cf: ngsolve.fem.CoefficientFunction = <ngsolve.fem.CoefficientFunction object at 0x7fba7da40d70>, deformation: object = None, ip_container: object = None, element_wise: bool = False, heapsize: int = 1000000) object
Integrate on a level set domains. The accuracy of the integration is ‘order’ w.r.t. a (multi-)linear approximation of the level set function. At first, this implies that the accuracy will, in general, only be second order. However, if the isoparametric approach is used (cf. lsetcurving functionality) this will be improved.
Parameters
- levelset_domaindictionary which provides levelsets, domain_types and integration specifica:
important keys are “levelset”, “domain_type”, “order”, the remainder are additional:
- “levelset”ngsolve.CoefficientFunction or a list thereof
CoefficientFunction that describes the geometry. In the best case lset is a GridFunction of an FESpace with scalar continuous piecewise (multi-) linear basis functions.
- “order”int
integration order.
- “domain_type”{NEG,POS,IF} (ENUM) or a list (of lists) thereof
Integration on the domain where either: * the level set function is negative (NEG) * the level set function is positive (POS) * the level set function is zero (IF )
- “subdivlvl”int
On simplex meshes a subtriangulation is created on which the level set function lset is interpolated piecewise linearly. Based on this approximation, the integration rule is constructed. Note: this argument only works on simplices without space-time and without multiple levelsets.
- “time_order”int
integration order in time for space-time integration
- “quad_dir_policy”int
policy for the selection of the order of integration directions
- mesh
Mesh to integrate on (on some part)
- cfngsolve.CoefficientFunction
the integrand
- deformationgridfunction (or None)
deformation of the mesh
- ip_containerlist (or None)
a list to store integration points (for debugging or visualization purposes)
- element_wisebool
result will return the integral w.r.t. each element individually.
- heapsizeint
heapsize for local computations.
- xfem.xfem.IntegrationPointExtrema(levelset_domain: dict, mesh: ngsolve.comp.Mesh, cf: ngsolve.fem.CoefficientFunction = <ngsolve.fem.CoefficientFunction object at 0x7fba7da41010>, heapsize: int = 1000000) tuple
Determine minimum and maximum on integration points on a level set domain. The sampling uses the same integration rule as in Integrate and is determined by ‘order’ w.r.t. a (multi-)linear approximation of the level set function. At first, this implies that the accuracy will, in general, only be second order. However, if the isoparametric approach is used (cf. lsetcurving functionality) this will be improved.
Parameters
- levelset_domaindictionary which provides levelsets, domain_types and integration specifica:
important keys are “levelset”, “domain_type”, “order”, the remainder are additional:
- “levelset”ngsolve.CoefficientFunction or a list thereof
CoefficientFunction that describes the geometry. In the best case lset is a GridFunction of an FESpace with scalar continuous piecewise (multi-) linear basis functions.
- “order”int
integration order.
- “domain_type”{NEG,POS,IF} (ENUM) or a list (of lists) thereof
Integration on the domain where either: * the level set function is negative (NEG) * the level set function is positive (POS) * the level set function is zero (IF )
- “subdivlvl”int
On simplex meshes a subtriangulation is created on which the level set function lset is interpolated piecewise linearly. Based on this approximation, the integration rule is constructed. Note: this argument only works on simplices without space-time and without multiple levelsets.
- “time_order”int
integration order in time for space-time integration
- “quad_dir_policy”int
policy for the selection of the order of integration directions
- mesh
Mesh to integrate on (on some part)
- cfngsolve.CoefficientFunction
the integrand
- heapsizeint
heapsize for local computations.
- xfem.xfem.InterpolateToP1(*args, **kwargs)
Overloaded function.
InterpolateToP1(gf_ho: ngsolve.comp.GridFunction = 0, gf_p1: ngsolve.comp.GridFunction = 0, eps_perturbation: float = 1e-14, heapsize: int = 1000000) -> None
Takes the vertex values of a GridFunction (also possible with a CoefficentFunction) and puts them into a piecewise (multi-) linear function.
Parameters
- gf_hongsolve.GridFunction
Function to interpolate
- gf_p1ngsolve.GridFunction
Function to interpolate to (should be P1)
- eps_perturbationfloat
If the absolute value if the function is smaller than eps_perturbation, it will be set to eps_perturbation. Thereby, exact and close-to zeros at vertices are avoided (Useful to reduce cut configurations for level set based methods).
- heapsizeint
heapsize of local computations.
InterpolateToP1(coef: ngsolve.fem.CoefficientFunction, gf: ngsolve.comp.GridFunction, eps_perturbation: float = 1e-14, heapsize: int = 1000000) -> None
Takes the vertex values of a CoefficentFunction) and puts them into a piecewise (multi-) linear function.
Parameters
- coefngsolve.CoefficientFunction
Function to interpolate
- gf_p1ngsolve.GridFunction
Function to interpolate to (should be P1)
- eps_perturbationfloat
If the absolute value if the function is smaller than eps_perturbation, it will be set to eps_perturbation. Thereby, exact and close-to zeros at vertices are avoided (Useful to reduce cut configurations for level set based methods).
- heapsizeint
heapsize of local computations.
- class xfem.xfem.MultiLevelsetCutInfo
Bases:
pybind11_objectA minimal version of a CutInfo that allows for several levelsets and a list of tuples of domain_types.
- GetElementsOfType(self: xfem.xfem.MultiLevelsetCutInfo, domain_type: object, VOL_or_BND: ngsolve.comp.VorB = <VorB.VOL: 0>, heapsize: int = 1000000) pyngcore.pyngcore.BitArray
Returns BitArray that is true for every element that has the corresponding domain type. This BitArray remains attached to the mlci class instance and is updated on mlci.Update(lsets).
Parameters
- domain_type{tuple(ENUM), list(tuple(ENUM)), DomainTypeArray}
Description of the domain.
- heapsizeint = 1000000
heapsize of local computations.
- GetElementsWithContribution(self: xfem.xfem.MultiLevelsetCutInfo, domain_type: object, VOL_or_BND: ngsolve.comp.VorB = <VorB.VOL: 0>, heapsize: int = 1000000) pyngcore.pyngcore.BitArray
Returns BitArray that is true for every element that has the a contribution to the corresponding level set domain. This BitArray remains attached to the mlci class instance and is updated on mlci.Update(lsets).
Parameters
- domain_type{tuple(ENUM), list(tuple(ENUM)), DomainTypeArray}
Description of the domain.
- heapsizeint = 1000000
heapsize of local computations.
- Mesh(self: xfem.xfem.MultiLevelsetCutInfo) ngsolve.comp.Mesh
Returns mesh of CutInfo.
- Update(self: xfem.xfem.MultiLevelsetCutInfo, levelsets: object, heapsize: int = 1000000) None
Updates the tuple of levelsets behind the MultiLevelsetCutInfo and recomputes any element marker arrays which have been created with this instance.
Parameters
- levelsetstuple(ngsolve.GridFunction)
tuple of GridFunctions w.r.t. which elements are marked.
- heapsizeint = 1000000
heapsize of local computations
- class xfem.xfem.P1Prolongation
Bases:
ProlongationProlongation for P1-type spaces (with possibly inactive dofs) — As is asks the fespace for dofs to vertices at several occasions the current implementation is not very fast and should be primarily used for prototype and testing…
- Update(self: xfem.xfem.P1Prolongation, space: ngsolve.comp.FESpace) None
- class xfem.xfem.P2CutProlongation
Bases:
ProlongationProlongation for P2 spaces (with possibly inactive dofs) — As is asks the fespace for dofs to vertices at several occasions the current implementation is not very fast and should be primarily used for prototype and testing…
- Update(self: xfem.xfem.P2CutProlongation, space: ngsolve.comp.FESpace) None
- class xfem.xfem.P2Prolongation
Bases:
ProlongationProlongation for P2 spaces (with possibly inactive dofs) — As is asks the fespace for dofs to vertices at several occasions the current implementation is not very fast and should be primarily used for prototype and testing…
- Update(self: xfem.xfem.P2Prolongation, space: ngsolve.comp.FESpace) None
- xfem.xfem.PatchwiseSolve(elagg: xfem.xfem.ElementAggregation, fes: ngsolve.comp.FESpace, bf: ngsolve.comp.SumOfIntegrals, lf: ngsolve.comp.SumOfIntegrals, heapsize: int = 1000000) ngsolve.la.BaseVector
Solve patch-wise problem based on the patches provided by the element aggregation input.
Parameters
- elagg: ElementAggregatetion
The instance defining the patches
- fesFESpace
The finite element space on which the local solve is performed.
- bfSumOfIntegrals
Integrators defining the matrix problem.
- lfSumOfIntegrals
Integrators defining the right-hand side.
- xfem.xfem.ProjectShift(lset_ho: ngsolve.comp.GridFunction = 0, lset_p1: ngsolve.comp.GridFunction = 0, deform: ngsolve.comp.GridFunction = 0, qn: ngsolve.fem.CoefficientFunction = 0, active_elements: object = None, blending: ngsolve.fem.CoefficientFunction = 0, lower: float = 0.0, upper: float = 0.0, threshold: float = 1.0, heapsize: int = 1000000) None
- class xfem.xfem.QUAD_DIRECTION_POLICY
Bases:
pybind11_objectMembers:
FIRST
OPTIMAL
FALLBACK
- FALLBACK = <QUAD_DIRECTION_POLICY.FALLBACK: 2>
- FIRST = <QUAD_DIRECTION_POLICY.FIRST: 0>
- OPTIMAL = <QUAD_DIRECTION_POLICY.OPTIMAL: 1>
- property name
- property value
- xfem.xfem.ReferenceTimeVariable() xfem.xfem.TimeVariableCoefficientFunction
This is the time variable. Call tref = ReferenceTimeVariable() to have a symbolic variable for the time like x,y,z for space. That can be used e.g. in lset functions for unfitted methods. Note that one would typically use tref in [0,1] as one time slab, leading to a call like t = told + delta_t * tref, when tref is our ReferenceTimeVariable. ngsxfem.__init__ defines tref.
- xfem.xfem.RefineAtLevelSet(gf: ngsolve.comp.GridFunction = 0, lower: float = 0.0, upper: float = 0.0, heapsize: int = 1000000) None
Mark mesh for refinement on all elements where the piecewise linear level set function lset_p1 has values in the interval [lower,upper] (default [0,0]).
Parameters
- gfngsolve.GridFunction
Scalar piecewise (multi-)linear Gridfunction
- lowerfloat
smallest level set value of interest
- upperfloat
largest level set value of interest
- heapsizeint
heapsize of local computations.
- class xfem.xfem.Restrict
Bases:
CompressWrapper Finite Element Spaces. The restricted fespace is a wrapper around a standard fespace which removes dofs from marked elements.
Parameters:
- fespacengsolve.comp.FESpace
finite element space
- active_elsBitArray or None
Only use dofs from these elements
- GetBaseSpace(self: xfem.xfem.Restrict) ngsolve.comp.FESpace
- property active_elements
active elements
- xfem.xfem.RestrictGFInTime(spacetime_gf: ngsolve.comp.GridFunction, reference_time: float = 0.0, space_gf: ngsolve.comp.GridFunction) None
Extract Gridfunction corresponding to a fixed time from a space-time GridFunction.
- class xfem.xfem.RestrictedBilinearFormComplex
Bases:
BilinearFormBilinearForm restricted on a set of elements and facets.
- property element_restriction
element restriction
- property facet_restriction
facet restriction
- class xfem.xfem.RestrictedBilinearFormDouble
Bases:
BilinearFormBilinearForm restricted on a set of elements and facets.
- property element_restriction
element restriction
- property facet_restriction
facet restriction
- xfem.xfem.SFESpace(arg0: ngsolve.comp.Mesh, arg1: ngsolve.fem.CoefficientFunction, arg2: int, arg3: dict) ngsolve.comp.FESpace
This is a special finite elemetn space which is a 1D polynomial along the zero level of the linearly approximated level set function lset and constantly extended in normal directions to this.
- class xfem.xfem.ScalarTimeFE
Bases:
FiniteElement
- xfem.xfem.SpaceTimeFESpace(spacefes: ngsolve.comp.FESpace, timefe: ngsolve.fem.FiniteElement, dirichlet: object = None, heapsize: int = 1000000, **kwargs) ngcomp::SpaceTimeFESpace
This function creates a SpaceTimeFiniteElementSpace based on a spacial FE space and a time Finite element Roughly, this is the tensor product between those two arguments. Further arguments specify several details.
Parameters
- spacefesngsolve.FESpace
This is the spacial finite element used for the space-time discretisation. Both scalar and vector valued spaces might be used. An example would be spacefes = H1(mesh, order=order) for given mesh and order.
- timefengsolve.FiniteElement
This is the time finite element for the space-time discretisation. That is essentially a simple finite element on the unit interval. There is a class ScalarTimeFE to create something fitting here. For example, one could call timefe = ScalarTimeFE(order) to create a time finite element of order order.
- dirichletlist or string
The boundary of the space domain which should have Dirichlet boundary values. Specification policy is the same as with the usual space finite element spaces.
- heapsizeint
Size of the local heap of this class. Increase this if you observe errors which look like a heap overflow.
dgjumps : bool
- xfem.xfem.SpaceTimeInterpolateToP1(spacetime_cf: ngsolve.fem.CoefficientFunction, time: ngsolve.fem.CoefficientFunction, spacetime_gf: ngsolve.comp.GridFunction) None
Interpolate nodal in time (possible high order) and nodal in space (P1).
- class xfem.xfem.SpaceTimeVTKOutput
Bases:
pybind11_object- Do(*args, **kwargs)
Overloaded function.
Do(self: xfem.xfem.SpaceTimeVTKOutput, vb: ngsolve.comp.VorB = <VorB.VOL: 0>, t_start: float = 0, t_end: float = 1) -> None
Do(self: xfem.xfem.SpaceTimeVTKOutput, vb: ngsolve.comp.VorB = <VorB.VOL: 0>, t_start: float = 0, t_end: float = 1, drawelems: pyngcore.pyngcore.BitArray) -> None
- xfem.xfem.SymbolicCutBFI(levelset_domain: dict, form: ngsolve.fem.CoefficientFunction, VOL_or_BND: ngsolve.comp.VorB = <VorB.VOL: 0>, element_boundary: bool = False, skeleton: bool = False, definedon: object = None, definedonelements: object = None, deformation: object = None) ngsolve.fem.BFI
see documentation of SymbolicBFI (which is a wrapper)
- xfem.xfem.SymbolicCutLFI(levelset_domain: dict, form: ngsolve.fem.CoefficientFunction, VOL_or_BND: ngsolve.comp.VorB = <VorB.VOL: 0>, element_boundary: bool = False, skeleton: bool = False, definedon: object = None, definedonelements: object = None, deformation: object = None) ngsolve.fem.LFI
see documentation of SymbolicLFI (which is a wrapper)
- xfem.xfem.SymbolicFacetPatchBFI(form: ngsolve.fem.CoefficientFunction, force_intorder: int = -1, time_order: int = -1, skeleton: bool = True, definedonelements: object = None, deformation: object = None, downscale: object = None) ngsolve.fem.BFI
Integrator on facet patches. Two versions are possible: * Either (skeleton=False) an integration on the element patch consisting of two neighboring elements is applied, * or (skeleton=True) the integration is applied on the facet.
Parameters
- formngsolve.CoefficientFunction
var form to integrate
- force_intorderint
(only active in the facet patch case (skeleton=False)) use this integration order in the integration
- skeletonboolean
decider on facet patch vs facet integration
- definedonelementsngsolve.BitArray/None
array which decides on which facets the integrator should be applied
- time_orderint
order in time that is used in the space-time integration. time_order=-1 means that no space-time rule will be applied. This is only relevant for space-time discretizations.
- downscaleNone | double
Downscale facet patch around facets.
- class xfem.xfem.TIME_DOMAIN_TYPE
Bases:
pybind11_objectMembers:
BOTTOM
TOP
INTERVAL
- BOTTOM = <TIME_DOMAIN_TYPE.BOTTOM: 0>
- INTERVAL = <TIME_DOMAIN_TYPE.INTERVAL: 2>
- TOP = <TIME_DOMAIN_TYPE.TOP: 1>
- property name
- property value
- class xfem.xfem.TimeVariableCoefficientFunction
Bases:
CoefficientFunction- FixTime(self: xfem.xfem.TimeVariableCoefficientFunction, arg0: float) None
- IsFixed(self: xfem.xfem.TimeVariableCoefficientFunction) bool
- UnfixTime(self: xfem.xfem.TimeVariableCoefficientFunction) None
- xfem.xfem.XFESpace(basefes: ngsolve.comp.FESpace, cutinfo: object = None, lset: object = None, flags: dict = {}, heapsize: int = 1000000) ngcomp::XFESpace
Constructor for XFESpace [For documentation of XFESpace-class see help(CXFESpace)]:
Extended finite element space. Takes a basis FESpace and creates an enrichment space based on cut information. The cut information is provided by a CutInfo object or - if a level set function is only provided - a CutInfo object is created. The enrichment doubles the unknowns on all cut elements and assigns to them a sign (NEG/POS). One of the differential operators neg(…) or pos(…) evaluates like the basis function of the origin space, the other one as zero for every basis function. Away from cut elements no basis function is supported.
Parameters
- basefesngsolve.FESpace
basic FESpace to be extended
- cutinfoxfem.CutInfo / None
Information on the cut configurations (cut elements, sign of vertices….)
- lsetngsolve.CoefficientFunction / None
level set function to construct own CutInfo (if no CutInfo is provided)
- flagsFlags
additional FESpace-flags
- heapsizeint
heapsize of local computations.
- xfem.xfem.XToNegPos(arg0: ngsolve.comp.GridFunction, arg1: ngsolve.comp.GridFunction) None
Takes a GridFunction of an extended FESpace, i.e. a compound space of V and VX = XFESpace(V) and interpretes it as a function in the CompoundFESpace of V and V. Updates the values of the vector of the corresponding second GridFunction.
- xfem.xfem.dn(*args, **kwargs)
Overloaded function.
dn(proxy: ngsolve.comp.ProxyFunction, order: int, comp: object = -1, hdiv: bool = False) -> ngsolve.comp.ProxyFunction
Normal derivative of higher order. This is evaluated via numerical differentiation which offers only limited accuracy (~ 1e-7).
Parameters
- proxyngsolve.ProxyFunction
test / trialfunction to the the normal derivative of
- orderint
order of derivative (in normal direction)
- compint
component of proxy if test / trialfunction is a component of a compound or vector test / trialfunction
- hdivboolean
assumes scalar FEs if false, otherwise assumes hdiv
dn(gf: ngsolve.comp.GridFunction, order: int) -> ngsolve.fem.CoefficientFunction
Normal derivative of higher order for a GridFunction. This is evaluated via numerical differentiation which offers only limited accuracy (~ 1e-7).
Parameters
- gfngsolve.GridFunction
(scalar) GridFunction to the the normal derivative of
- orderint
order of derivative (in normal direction)
- xfem.xfem.fix_tref_coef(arg0: ngsolve.fem.CoefficientFunction, arg1: object) ngsolve.fem.CoefficientFunction
fix_t fixes the evaluated time to a fixed value.
Parameters
- self: ngsolve.CoefficientFunction
CoefficientFunction in which the time should be fixed
- time: Parameter or double
Value the time should become (if Parameter, the value can be adjusted later on)
- xfem.xfem.fix_tref_gf(gf: ngsolve.comp.GridFunction, time: float = 0.0) ngsolve.fem.CoefficientFunction
fix_t fixes the time (ReferenceTimeVariable) of a given expression. This is the variant for a gridfunction.
Parameters
- self: ngsolve.GridFunction
Gridfunction in which the time should be fixed
- time: double
Value the time should become
- xfem.xfem.fix_tref_proxy(proxy: ngsolve.comp.ProxyFunction, time: float, comp: object = -1, use_FixAnyTime: bool = False) ngsolve.comp.ProxyFunction
- xfem.xfem.shifted_eval(gf: ngsolve.comp.GridFunction, back: object = None, forth: object = None) ngsolve.fem.CoefficientFunction
Returns a CoefficientFunction that evaluates Gridfunction gf at a shifted location, s.t. the original function to gf, gf: x -> f(x) is changed to cf: x -> f(s(x)) where z = s(x) is the shifted location that is computed ( pointwise ) from:
Psi_back(z) = Psi_forth(x),
< = > z = Inv(Psi_back)( Psi_forth(x) ) < = > s = Inv(Psi_back) o Psi_forth(x)
To compute z = s(x) a fixed point iteration is used.
ATTENTION:
If s(x) leaves the the element that the integration point x is defined on, it will NOT change the element but result in an integration point that lies outside of the physical element.
Parameters
- backngsolve.GridFunction
transformation describing Psi_back as I + d_back where d_back is the deformation (can be None).
- forthngsolve.GridFunction
transformation describing Psi_forth as I + d_forth where d_forth is the deformation (can be None).
ASSUMPTIONS:
2D or 3D mesh