xfem.lset_spacetime module
- class xfem.lset_spacetime.LevelSetMeshAdaptation_Spacetime(**kwargs)
Bases:
objectClass to compute a proper mesh deformation to improve a piecewise (multi-) linear (in space!) level set approximation to obtain a higher order accurate approximation.
Computes a continuous function that describes a shift between points that are on the (P1-in-space ) approximated level set function and its higher order accurate approximation. The transformation is only applied on elements where a level value inside a certain interval (lower,upper) exists.
The result is a space-time finite element deformation (deform). For each fixed time the behavior is as for an LevelSetMeshAdaptation object, i.e.
(1-b)*phi_h(x,t) + b*phi_lin( x,t ) = phi_h(Psi(x,t))
with Psi(x,t) = x + d(x,t) qn(x,t) =: D(x,t)
for all x on ‘cut’ elements
with
- phi_hself.lset_ho
the higher order space-time level set function
- phi_linself.lset_p1
the P1-in-space space-time level set function
- PsiId + self.deform
the resulting deformation
- qnself.qn
normal direction field
- b: self.smooth_blend
an optional space-time CoefficientFunction to ensure a smooth transition between curved and uncurved elements. This class offers the options of 1) a finite element blending, in which case b is assumed as b=0 and shall not be specified in the constructor, and 2) a smooth blending, in which the sufficiently regular function b ranging from 0 at cut locations to 1 outside with zero deformation is taken into consideration. Check out https://arxiv.org/abs/2311.02348 for a detailed mathematical description, the constructor of this class for specification, and the class LevelSet_SmoothBlending for generating function b in a pre-defined way.
- This class holds its own members for the higher order and lower order
(P1-in-space) space-time approximation of the level set function and only depends on the input of a mesh and a CoefficientFunction (the levelset function) (and options).
A LevelSetMeshAdaptation_Spacetime can also be used as a context manager. In this case, the mesh deformation is applied on the mesh inside the context.
- BFI(domain_type, time_type, form, time_order=None, definedonelements=None)
Convenience function to construct Symbolic(Cut)LFI based on the domain_type and the known space-time level set function.
- Args:
domain_type (DOMAIN_TYPE: NEG/POS/IF): domain type to integrate on time_type (BOTTOM/TOP/INTERVAL): type of integral on the space-time domain form (CoefficientFunction): integrand time_order (int, optional): integration order in time. Defaults to None. defineonelement (BitArray, optional): elements to define the bilinear form on. Defaults to None.
- CalcDeformation(**kwargs)
- CalcMaxDistance(**kwargs)
- Integrate(**kwargs)
- Integrator(SymbolicFI, domain_type, time_type, form, time_order=None, definedonelements=None)
Convenience function to construct space-time Symbolic(Cut)LFI/Symbolic(Cut)BFI
- Args:
SymbolicFI (SymbolicLFI/SymbolicBFI): type of integrator (LFI vs. BFI) domain_type (DOMAIN_TYPE: NEG/POS/IF): domain type to integrate on time_type (BOTTOM/TOP/INTERVAL): type of integral on the space-time domain form (CoefficientFunction): integrand time_order (int, optional): integration order in time. Defaults to None. defineonelement (BitArray, optional): elements to define the linear form on. Defaults to None.
- LFI(domain_type, time_type, form, time_order=None, definedonelements=None)
Convenience function to construct Symbolic(Cut)LFI based on the domain_type and the known space-time level set function.
- Args:
domain_type (DOMAIN_TYPE: NEG/POS/IF): domain type to integrate on time_type (BOTTOM/TOP/INTERVAL): type of integral on the space-time domain form (CoefficientFunction): integrand time_order (int, optional): integration order in time. Defaults to None. defineonelement (BitArray, optional): elements to define the linear form on. Defaults to None.
- ProjectGFs(**kwargs)
- ProjectOnUpdate(gf, update_domain=None)
When the LevelsetMeshAdaptation class generates a new deformation (due to a new level set function) all GridFunction that have been stored through ProjectOnUpdate will be projected from the previous to the new mesh by an essentially local projection (Oswald projection of the shifted evaluation)
- gfuGridFunction (or list of GridFunctions)
GridFunction(s) to store for later deformation updates.
- interpol_ho(**kwargs)
- interpol_p1(**kwargs)
- levelset_domain(domain_type=<DOMAIN_TYPE.IF: 2>, time_type=<TIME_DOMAIN_TYPE.INTERVAL: 2>)
- Return a levelset_domain dictionary.
- Args:
domain_type (DOMAIN_TYPE: NEG/POS/IF, optional): domain type for level set doamin. Defaults to IF. time_type (BOTTOM/TOP/INTERVAL): type of integral on the space-time domain
- Returns:
dict: levelset domain
- order_deform = 2
- order_lset = 2
- order_qn = 2