def test_2d_ho_integration():
    # -------------------- Test Higher-Order Integration ----------------------
    level_sets = [-y, x - 1, y - 1, -x]
    nr_ls = len(level_sets)

    # ---------------------------- Background Mesh ----------------------------
    geo = SplineGeometry()
    geo.AddRectangle((-0.2, -0.2), (1.2, 1.2),
                     bcs=("bottom", "right", "top", "left"))
    mesh = Mesh(geo.GenerateMesh(maxh=0.4))

    # ------------------------------- LEVELSET --------------------------------
    level_sets_p1 = tuple(GridFunction(H1(mesh, order=1)) for i in range(nr_ls))
    for i, lsetp1 in enumerate(level_sets_p1):
        InterpolateToP1(level_sets[i], lsetp1)

    square = DomainTypeArray([(NEG, NEG, NEG, NEG)])
    lset_dom_square = {"levelset": level_sets_p1, "domain_type": square}

    # ------------------------------- Integrate -------------------------------
    result = Integrate(levelset_domain=lset_dom_square, mesh=mesh, 
                       cf=x * (1 - x) * y * (1 - y), order=4)
    assert abs(result - 1/36) < 1e-12

    del mesh, level_sets, level_sets_p1, square 
Beispiel #2
0
def test_new_integrateX_via_straight_cutted_quad2D_polynomial(
        order, domain, quad_dominated, alpha, dim):
    square = SplineGeometry()
    square.AddRectangle([0, 0], [1, 1], bc=1)
    mesh = Mesh(square.GenerateMesh(maxh=100, quad_dominated=quad_dominated))

    levelset = 1 - 2 * x - 2 * y
    val_pos = pow(2, -alpha - 2) / (alpha * alpha + 3 * alpha + 2)
    referencevals = {POS: val_pos, NEG: 1. / (alpha + 1) - val_pos}

    lset_approx = GridFunction(H1(mesh, order=1))
    InterpolateToP1(levelset, lset_approx)

    f = pow(dim, alpha)

    integral = Integrate(levelset_domain={
        "levelset": lset_approx,
        "domain_type": domain
    },
                         cf=f,
                         mesh=mesh,
                         order=order)
    error = abs(integral - referencevals[domain])

    assert error < 5e-15 * (order + 1) * (order + 1)
Beispiel #3
0
def create_mesh(mesh_size):
    geo = SplineGeometry()
    geo.AddRectangle((0, 0), (2, 0.41), bcs=(
        "wall", "outlet", "wall", "inlet"))
    geo.AddCircle((0.2, 0.2), r=0.05, leftdomain=0, rightdomain=1, bc="cyl")
    mesh = Mesh(geo.GenerateMesh(maxh=mesh_size))
    mesh.Curve(3)
    return mesh
Beispiel #4
0
def area_of_a_sphere_ST_error(n_steps = 64, i=3, structured_mesh= True):
    if structured_mesh:
        length = 1
        mesh = MakeStructured2DMesh(quads=False,nx=2**(i),ny=2**(i),mapping= lambda x,y : (2*length*x-length,2*length*y-length))
    else:
        square = SplineGeometry()
        square.AddRectangle([-1,-1],[1,1])
        ngmesh = square.GenerateMesh(maxh=(1/2)**(i-1), quad_dominated=False)
        mesh = Mesh (ngmesh)

    coef_told = Parameter(0)
    coef_delta_t = Parameter(0)
    tref = ReferenceTimeVariable()
    t = coef_told + coef_delta_t*tref
    
    r0 = 0.9
    r = sqrt(x**2+y**2+t**2)
    
    # level set
    levelset= r - r0
    
    time_order = 1
    fes1 = H1(mesh, order=1)
    tfe = ScalarTimeFE(time_order)
    st_fes = SpaceTimeFESpace(fes1,tfe)
    
    tend = 1
    delta_t = tend/n_steps
    coef_delta_t.Set(delta_t)
    told = 0

    lset_p1 = GridFunction(st_fes)
    
    sum_vol = 0
    sum_int = 0
    for i in range(n_steps):
        SpaceTimeInterpolateToP1(levelset,tref,0.,delta_t,lset_p1) # call for the master spacetime_weihack -- 0 and tend are ununsed parameter
    
        val_vol = Integrate({ "levelset" : lset_p1, "domain_type" : NEG}, CoefficientFunction(1.0), mesh, time_order = time_order)
        val_int = Integrate({ "levelset" : lset_p1, "domain_type" : IF}, CoefficientFunction(1.0), mesh, time_order = time_order)
        #print(val_vol, val_int)
        sum_vol += val_vol*delta_t
        sum_int += val_int*delta_t
        
        told = told + delta_t
        coef_told.Set(told)

    print("SUM VOL: ", sum_vol)
    print("VOL: ", 2/3*pi*r0**3)
    vol_err = abs(sum_vol - 2/3*pi*r0**3)
    print("\t\tDIFF: ", vol_err)
    
    print("SUM INT: ", sum_int)
    print("AREA: ", 0.5*pi**2*r0**2)
    int_err = abs(sum_int - 0.5*pi**2*r0**2)
    print("\t\tDIFF: ",int_err)
    return (vol_err, int_err)
def test_spacetime_lsetcurving_maxdist(imax, order):

    ngsglobals.msg_level = 1

    maxdist_at_reflevels = []
    # polynomial order in time
    k_t = order
    # polynomial order in space
    k_s = k_t

    for i in range(imax):
        square = SplineGeometry()
        square.AddRectangle([-0.6, -0.6], [0.6, 0.6])
        ngmesh = square.GenerateMesh(maxh=0.5**(i + 1), quad_dominated=False)
        mesh = Mesh(ngmesh)

        tstart = 0
        tend = 0.5
        delta_t = (tend - tstart) / (2**(i))

        told = Parameter(tstart)
        t = told + delta_t * tref

        lset_adap_st = LevelSetMeshAdaptation_Spacetime(mesh,
                                                        order_space=k_s,
                                                        order_time=k_t,
                                                        threshold=0.5,
                                                        discontinuous_qn=True)

        r0 = 0.5

        rho = CoefficientFunction((1 / (pi)) * sin(2 * pi * t))
        r = sqrt(x**2 + (y - rho)**2)
        levelset = CoefficientFunction(r - r0)

        maxdists = []
        while tend - told.Get() > delta_t / 2:
            dfm = lset_adap_st.CalcDeformation(levelset)
            mesh.deformation = dfm
            maxdist = lset_adap_st.CalcMaxDistance(levelset)
            mesh.deformation = None
            maxdists.append(maxdist)
            told.Set(told.Get() + delta_t)

        print("i: ", i, "\tdist max: ", max(maxdists))
        maxdist_at_reflevels.append(max(maxdists))

    eocs = [
        log(maxdist_at_reflevels[i - 1] / maxdist_at_reflevels[i]) / log(2)
        for i in range(1, len(maxdist_at_reflevels))
    ]
    print("eocs: ", eocs)
    avg_eoc = sum(eocs[-2:]) / len(eocs[-2:])
    print("avg: ", avg_eoc)

    assert avg_eoc > order + 0.75
Beispiel #6
0
def mesh():
    geo = SplineGeometry()
    geo.AddRectangle( (-10, 0.0), (10, 1), leftdomain=1, bcs=["bottom","right","slave","left"])

    geo.AddCircle ( (-1, 5), r=1, leftdomain=2, bc="master")
    geo.SetMaterial(1, "brick")
    geo.SetMaterial(2, "ball")
    mesh = Mesh( geo.GenerateMesh(maxh=0.5))
    mesh.Curve(5)
    Draw(mesh)
    return mesh
def test_calc_linearized():
    square = SplineGeometry()
    square.AddRectangle([-1.5, -1.5], [1.5, 1.5], bc=1)
    mesh = Mesh(square.GenerateMesh(maxh=0.8))

    levelset = (sqrt(sqrt(x * x * x * x + y * y * y * y)) - 1.0)

    lsetp1 = GridFunction(H1(mesh))
    lsetp1.Set(levelset)
    lset_neg = {"levelset": lsetp1, "domain_type": NEG, "subdivlvl": 0}

    ci = CutInfo(mesh, lsetp1)
    hasneg = ci.GetElementsOfType(HASNEG)

    Vh = H1(mesh, order=2, dirichlet=[])
    active_dofs = GetDofsOfElements(Vh, hasneg)
    Vh = Compress(Vh, active_dofs)

    u, v = Vh.TnT()

    gfu = GridFunction(Vh)
    gfu.Set(sin(x))

    a1 = BilinearForm(Vh)
    a1 += SymbolicBFI(levelset_domain=lset_neg, form=u * u * v)
    a1.AssembleLinearization(gfu.vec)

    a2 = BilinearForm(Vh)
    a2 += SymbolicBFI(levelset_domain=lset_neg, form=2 * gfu * u * v)
    a2.Assemble()

    a3 = BilinearForm(Vh)
    a3 += SymbolicBFI(levelset_domain=lset_neg, form=gfu * u * v)

    w1 = gfu.vec.CreateVector()
    w2 = gfu.vec.CreateVector()

    w1.data = a1.mat * gfu.vec
    w2.data = a2.mat * gfu.vec

    w1.data -= w2
    diff = Norm(w1)
    print("diff : ", diff)
    assert diff < 1e-12

    a1.Apply(gfu.vec, w1)
    a3.Apply(gfu.vec, w2)

    w1.data -= w2
    diff = Norm(w1)
    print("diff : ", diff)
    assert diff < 1e-12
Beispiel #8
0
def test_spacetime_vtk(quad):
    tref = ReferenceTimeVariable()
    square = SplineGeometry()
    A = 1.25
    square.AddRectangle([-A, -A], [A, A])
    ngmesh = square.GenerateMesh(maxh=0.3, quad_dominated=False)
    mesh = Mesh(ngmesh)

    #### expression for the time variable:

    coef_told = Parameter(0)
    coef_delta_t = Parameter(0)
    tref = ReferenceTimeVariable()
    t = coef_told + coef_delta_t * tref

    #### the data:
    # radius of disk (the geometry)
    r0 = 0.5
    x0 = lambda t: r0 * cos(pi * t)
    y0 = lambda t: r0 * sin(pi * t)

    #levelset= x - t
    levelset = sqrt((x - x0(t))**2 + (y - y0(t))**2) - 0.4
    # spatial FESpace for solution
    fesp1 = H1(mesh, order=1)
    # polynomial order in time for level set approximation
    lset_order_time = 1
    # integration order in time
    time_order = 2
    # time finite element (nodal!)
    tfe = ScalarTimeFE(1)
    # space-time finite element space
    st_fes = SpaceTimeFESpace(fesp1, tfe)

    coef_delta_t.Set(0.25)

    lset_p1 = GridFunction(st_fes)
    SpaceTimeInterpolateToP1(levelset, tref, lset_p1)

    vtk = SpaceTimeVTKOutput(ma=mesh,
                             coefs=[levelset, lset_p1],
                             names=["levelset", "lset_p1"],
                             filename="spacetime_vtk_test",
                             subdivision_x=1,
                             subdivision_t=1)
    vtk.Do(t_start=coef_told.Get(), t_end=coef_told.Get() + coef_delta_t.Get())
    import os
    os.remove("spacetime_vtk_test_0.vtk")
def test_multilevelsetcutinfo():
    # ---------------------------- Background Mesh ----------------------------
    geo = SplineGeometry()
    geo.AddRectangle((-1, -1), (1, 1), bc=1)
    mesh = Mesh(geo.GenerateMesh(maxh=0.2))

    # ------------------------------ Test Update ------------------------------
    ba = [BitArray(mesh.ne) for i in range(4)]
    for ba_ in ba:
        ba_.Clear()

    P1 = H1(mesh, order=1)
    lsets = tuple(GridFunction(P1) for i in range(2))

    InterpolateToP1(x + 0.5, lsets[0])
    InterpolateToP1(x - 0.5, lsets[1])

    mlci = MultiLevelsetCutInfo(mesh, lsets)
    ba[0] |= mlci.GetElementsOfType((POS, NEG))

    InterpolateToP1(y + 0.5, lsets[0])
    InterpolateToP1(y - 0.5, lsets[1])    

    ba[1] |= mlci.GetElementsOfType((POS, NEG))
    assert MarkersEqual(ba[0], ba[1])

    mlci.Update(lsets)
    ba[2] |= mlci.GetElementsOfType((POS, NEG))
    assert MarkersEqual(ba[0], ba[2]) == False

    # -------------------------- Test Safety Checks ---------------------------
    with pytest.raises(netgen.libngpy._meshing.NgException):
        mlci2 = MultiLevelsetCutInfo(mesh, (x, y))

    with pytest.raises(netgen.libngpy._meshing.NgException):
        a1 = mlci.GetElementsOfType((NEG, NEG, NEG))
    with pytest.raises(netgen.libngpy._meshing.NgException):
        a2 = mlci.GetElementsOfType([(NEG, NEG, NEG)])
    with pytest.raises(netgen.libngpy._meshing.NgException):
        a3 = mlci.GetElementsWithContribution((NEG, NEG, NEG))
    with pytest.raises(netgen.libngpy._meshing.NgException):
        a4 = mlci.GetElementsWithContribution([(NEG, NEG, NEG)])
Beispiel #10
0
def test_new_integrateX_via_straight_cutted_quad2D(order, domain, quad):
    square = SplineGeometry()
    square.AddRectangle([0, 0], [1, 1], bc=1)
    mesh = Mesh(square.GenerateMesh(maxh=100, quad_dominated=quad))

    levelset = 1 - 2 * x - 2 * y
    referencevals = {NEG: 7 / 8, POS: 1 / 8, IF: 1 / sqrt(2)}

    f = CoefficientFunction(1)

    integral = Integrate(levelset_domain={
        "levelset": levelset,
        "domain_type": domain
    },
                         cf=f,
                         mesh=mesh,
                         order=order)
    error = abs(integral - referencevals[domain])

    assert error < 5e-15 * (order + 1) * (order + 1)
Beispiel #11
0
def test_new_integrateX_via_circle_geom(quad_dominated, order, domain):
    square = SplineGeometry()
    square.AddRectangle([0, 0], [1, 1], bc=1)
    mesh = Mesh(square.GenerateMesh(maxh=100, quad_dominated=quad_dominated))
    r = 0.6

    levelset = sqrt(x * x + y * y) - r
    referencevals = {
        POS: 1 - pi * r * r / 4,
        NEG: pi * r * r / 4,
        IF: r * pi / 2
    }

    n_ref = 8
    errors = []

    for i in range(n_ref):
        f = CoefficientFunction(1)

        integral = Integrate(levelset_domain={
            "levelset": levelset,
            "domain_type": domain
        },
                             cf=f,
                             mesh=mesh,
                             order=order)
        print("Result of Integration Reflevel ", i, ", Key ", domain, " : ",
              integral)
        errors.append(abs(integral - referencevals[domain]))

        if i < n_ref - 1:
            mesh.Refine()

    eoc = [log(errors[i + 1] / errors[i]) / log(0.5) for i in range(n_ref - 1)]

    print("L2-errors:", errors)
    print("experimental order of convergence (L2):", eoc)

    mean_eoc_array = eoc[1:]
    mean_eoc = sum(mean_eoc_array) / len(mean_eoc_array)
    assert mean_eoc > 1.75
Beispiel #12
0
def test_apply():
    square = SplineGeometry()
    square.AddRectangle([-1.5, -1.5], [1.5, 1.5], bc=1)
    mesh = Mesh(square.GenerateMesh(maxh=0.8))

    levelset = (sqrt(sqrt(x * x * x * x + y * y * y * y)) - 1.0)

    lsetp1 = GridFunction(H1(mesh))
    lsetp1.Set(levelset)
    lset_neg = {"levelset": lsetp1, "domain_type": NEG, "subdivlvl": 0}

    ci = CutInfo(mesh, lsetp1)
    hasneg = ci.GetElementsOfType(HASNEG)

    Vh = H1(mesh, order=2, dirichlet=[])
    active_dofs = GetDofsOfElements(Vh, hasneg)
    Vh = Compress(Vh, active_dofs)

    u, v = Vh.TnT()

    a = BilinearForm(Vh)
    a += SymbolicBFI(levelset_domain=lset_neg, form=u * v)
    a += SymbolicBFI(levelset_domain=lset_neg, form=grad(u)[0] * v)

    gfu = GridFunction(Vh)
    gfu.Set(1)

    a.Assemble()

    Au1 = gfu.vec.CreateVector()
    Au2 = gfu.vec.CreateVector()

    Au1.data = a.mat * gfu.vec

    a.Apply(gfu.vec, Au2)

    Au1.data -= Au2
    diff = Norm(Au1)
    print("diff : ", diff)
    assert diff < 1e-12
# with Dirichlet boundary condition u = 0

from ngsolve import *
from ngsolve.krylovspace import CG
from netgen.geom2d import unit_square
from netgen.geom2d import SplineGeometry
import ngsolve
from xfem import *
from xfem.cutmg import *  #MultiGridCL, CutFemSmoother, LinearMGIterator
from xfem.lsetcurv import *

ngsglobals.msg_level = 1

# generate a triangular mesh of mesh-size 0.2
square = SplineGeometry()
square.AddRectangle([-1.5, -1.5], [1.5, 1.5], bc=1)
mesh = Mesh(square.GenerateMesh(maxh=0.5, quad_dominated=False))

mu = [1e-3, 1]
order = 2
kappac = "harmonic"
gamma_stab = 10

rhsscal = -4 * mu[0] * mu[1]
coef_f = [rhsscal, rhsscal]
distsq = x * x + y * y - 1
solution = [mu[1] * distsq, mu[0] * distsq]

levelset = (sqrt(x * x + y * y) - 1)
subdivlvl = 0
Beispiel #14
0
from ngsolve import *

# viscosity
nu = 0.001

# timestepping parameters
tau = 0.001
tend = 10

# mesh = Mesh("cylinder.vol")
from netgen.geom2d import SplineGeometry
geo = SplineGeometry()
geo.AddRectangle((0, 0), (2, 0.41), bcs=("wall", "outlet", "wall", "inlet"))
geo.AddCircle((0.2, 0.2), r=0.05, leftdomain=0, rightdomain=1, bc="cyl")
mesh = Mesh(geo.GenerateMesh(maxh=0.08))

mesh.Curve(3)

V = H1(mesh, order=3, dirichlet="wall|cyl|inlet")
Q = H1(mesh, order=2)

X = FESpace([V, V, Q])

ux, uy, p = X.TrialFunction()
vx, vy, q = X.TestFunction()

div_u = grad(ux)[0] + grad(uy)[1]
div_v = grad(vx)[0] + grad(vy)[1]

stokes = nu * grad(ux) * grad(vx) + nu * grad(uy) * grad(
    vy) + div_u * q + div_v * p - 1e-10 * p * q
Beispiel #15
0
from netgen.meshing import MeshingParameters
from ngsolve import *
from xfem import *
# For LevelSetAdaptationMachinery
from xfem.lsetcurv import *

r2 = 3 / 4  # outer radius
r1 = 1 / 4  # inner radius
rc = (r1 + r2) / 2.0
rr = (r2 - r1) / 2.0
r = sqrt(x * x + y * y)
levelset = IfPos(r - rc, r - rc - rr, rc - r - rr)

# Geometry
square = SplineGeometry()
square.AddRectangle([-1, -1], [1, 1], bc=1)


def meshonlvl(lvl):
    ngmesh = square.GenerateMesh(maxh=0.1, quad_dominated=False)
    mesh = Mesh(ngmesh)
    for i in range(lvl):
        mesh.Refine()
    return mesh


#trig version:
ngmesh = square.GenerateMesh(maxh=0.1, quad_dominated=False)
mesh = Mesh(ngmesh)
p1prolong = P2CutProlongation(mesh)
Beispiel #16
0
aneg = 1.0 / mu[0]
apos = 1.0 / mu[1] + (1.0 / mu[0] - 1.0 / mu[1]) * exp(x * x + y * y - R * R)

# Discretisation parameters

# Stabilization parameters for ghost-penalty
gamma_stab_v = 0.05  # if set to zero: no GP-stabilization for velocity
gamma_stab_p = 0.05
# Stabilization parameter for Nitsche
lambda_nitsche = 0.5 * (mu[0] + mu[1]) * 20 * order * order

# ----------------------------------- MAIN ------------------------------------
# Generate the background mesh
square = SplineGeometry()
square.AddRectangle((-1, -1), (1, 1), bcs=[1, 2, 3, 4])
mesh = Mesh(square.GenerateMesh(maxh=maxh, quad_dominated=False))

d = mesh.dim

# Construct the exact level set function, radius and radius squared:
rsqr = x**2 + y**2
r = sqrt(rsqr)

levelset = r - R

# Construct the exact solution:
gammaf = 0.5  # surface tension = pressure jump
vel_exact = [
    a * CoefficientFunction((-y, x)) * exp(-1 * rsqr) for a in [aneg, apos]
]
Beispiel #17
0
from ngsolve import *
from time import sleep
from netgen.geom2d import unit_square
from netgen.geom2d import SplineGeometry
from netgen.meshing import MeshingParameters
from ngsolve.internal import *
from xfem import *
from numpy import pi
from xfem.lset_spacetime import *


square = SplineGeometry()
square.AddRectangle([0,0],[2,2],bc=1)
maxh = 0.3
mesh = Mesh (square.GenerateMesh(maxh=maxh, quad_dominated=False))


h = specialcf.mesh_size

# FE Spaces
k_t = 1
k_s = 1
fes1 = H1(mesh,order=k_s,dirichlet=[1,2,3,4])

lset_order_time = 2
tfe = ScalarTimeFE(k_t) 

st_fes = SpaceTimeFESpace(fes1,tfe)
#visoptions.deformation = 1

tend = 1
Beispiel #18
0
# Finite element space order
order = 2

# Diffusion coefficients for the sub-domains (NEG/POS):
alpha = [1.0, 2.0]
# Nitsche penalty parameter
lambda_nitsche = 20

formulation = "XFEM"  # or "CUTFEM"

# ----------------------------------- MAIN ------------------------------------
# We generate the background mesh of the domain and use a simplicity
# triangulation to obtain a mesh with quadrilaterals use
# 'quad_dominated=True'
square = SplineGeometry()
square.AddRectangle(ll, ur, bc=1)
mesh = Mesh(square.GenerateMesh(maxh=maxh, quad_dominated=False))

# Manufactured solution and corresponding r.h.s. data CoefficientFunctions:
r22 = x**2 + y**2
r44 = x**4 + y**4
r66 = x**6 + y**6
r41 = sqrt(sqrt(r44))
r4m3 = 1.0 / r41**3
solution = [1 + pi / 2 - sqrt(2.0) * cos(pi / 4 * r44), pi / 2 * r41]
coef_f = [
    -alpha[i] * (solution[i].Diff(x).Diff(x) + solution[i].Diff(y).Diff(y))
    for i in range(2)
]

# Level set function of the domain (phi = ||x||_4 - 1) and its interpolation:
Beispiel #19
0
# Mesh diameter
maxh = 0.1
# Finite element space order
order = 2

# Stabilization parameter for ghost-penalty
gamma_stab = 0.1
# Stabilization parameter for Nitsche
lambda_nitsche = 10 * order * order
# Stabilization parameter for interior penalty
lambda_dg = 10 * order * order

# ----------------------------------- MAIN ------------------------------------
# Geometry and Mesh
square = SplineGeometry()
square.AddRectangle((-1, -1), (1, 1), bc=1)

ngmesh = square.GenerateMesh(maxh=maxh, quad_dominated=quad_mesh)
mesh = Mesh(ngmesh)

# Manufactured exact solution for monitoring the error
r2 = 3 / 4  # outer radius
r1 = 1 / 4  # inner radius
rc = (r1 + r2) / 2.0
rr = (r2 - r1) / 2.0
r = sqrt(x**2 + y**2)
levelset = IfPos(r - rc, r - rc - rr, rc - r - rr)

exact = (20 * (r2 - sqrt(x**2 + y**2)) * (sqrt(x**2 + y**2) - r1)).Compile()
coeff_f = -(exact.Diff(x).Diff(x) + exact.Diff(y).Diff(y)).Compile()
Beispiel #20
0
def test_spacetime_spaceP1_timeDGP1():
    ngsglobals.msg_level = 1

    square = SplineGeometry()
    square.AddRectangle([-1, -1], [1, 1])
    ngmesh = square.GenerateMesh(maxh=0.08, quad_dominated=False)
    mesh = Mesh(ngmesh)

    coef_told = Parameter(0)
    coef_delta_t = Parameter(0)
    tref = ReferenceTimeVariable()
    t = coef_told + coef_delta_t * tref

    r0 = 0.5

    # position shift of the geometry in time
    rho = CoefficientFunction((1 / (pi)) * sin(2 * pi * t))
    rhoL = lambda t: CoefficientFunction((1 / (pi)) * sin(2 * pi * t))
    #convection velocity:
    d_rho = CoefficientFunction(2 * cos(2 * pi * t))
    w = CoefficientFunction((0, d_rho))

    # level set
    r = sqrt(x**2 + (y - rho)**2)
    levelset = r - r0

    # diffusion coefficient
    alpha = 1

    # solution and r.h.s.
    Q = pi / r0
    u_exact = cos(Q * r) * sin(pi * t)
    u_exactL = lambda t: cos(Q * sqrt(x**2 + (y - rhoL(t))**2)) * sin(pi * t)
    coeff_f = (Q / r * sin(Q * r) + (Q**2) * cos(Q * r)) * sin(
        pi * t) + pi * cos(Q * r) * cos(pi * t)
    u_init = u_exact

    # polynomial order in time
    k_t = 1
    # polynomial order in space
    k_s = 1
    # spatial FESpace for solution
    fes1 = H1(mesh, order=k_s)
    # polynomial order in time for level set approximation
    lset_order_time = 1
    # integration order in time
    time_order = 2
    # time finite element (nodal!)
    tfe = ScalarTimeFE(k_t)
    # space-time finite element space
    st_fes = SpaceTimeFESpace(fes1, tfe, flags={"dgjumps": True})

    #Unfitted heat equation example
    tend = 1
    delta_t = tend / 32
    coef_delta_t.Set(delta_t)
    tnew = 0
    told = 0

    lset_p1 = GridFunction(st_fes)

    SpaceTimeInterpolateToP1(levelset, tref, lset_p1)

    lset_top = CreateTimeRestrictedGF(lset_p1, 1.0)
    lset_bottom = CreateTimeRestrictedGF(lset_p1, 0.0)

    gfu = GridFunction(st_fes)

    u_last = CreateTimeRestrictedGF(gfu, 0)
    u_last.Set(u_exactL(0.))

    u, v = st_fes.TnT()

    h = specialcf.mesh_size

    lset_neg = {"levelset": lset_p1, "domain_type": NEG, "subdivlvl": 0}
    lset_neg_bottom = {
        "levelset": lset_bottom,
        "domain_type": NEG,
        "subdivlvl": 0
    }
    lset_neg_top = {"levelset": lset_top, "domain_type": NEG, "subdivlvl": 0}

    def SpaceTimeNegBFI(form):
        return SymbolicBFI(levelset_domain=lset_neg,
                           form=form,
                           time_order=time_order)

    ci = CutInfo(mesh, time_order=time_order)

    hasneg_integrators_a = []
    hasneg_integrators_f = []
    patch_integrators_a = []
    hasneg_integrators_a.append(
        SpaceTimeNegBFI(form=delta_t * alpha * grad(u) * grad(v)))
    hasneg_integrators_a.append(
        SymbolicBFI(levelset_domain=lset_neg_top,
                    form=fix_tref(u, 1) * fix_tref(v, 1)))
    hasneg_integrators_a.append(SpaceTimeNegBFI(form=-u * dt(v)))
    hasneg_integrators_a.append(
        SpaceTimeNegBFI(form=-delta_t * u * InnerProduct(w, grad(v))))
    patch_integrators_a.append(
        SymbolicFacetPatchBFI(form=delta_t * 1.05 * h**(-2) * (u - u.Other()) *
                              (v - v.Other()),
                              skeleton=False,
                              time_order=time_order))
    hasneg_integrators_f.append(
        SymbolicLFI(levelset_domain=lset_neg,
                    form=delta_t * coeff_f * v,
                    time_order=time_order))
    hasneg_integrators_f.append(
        SymbolicLFI(levelset_domain=lset_neg_bottom,
                    form=u_last * fix_tref(v, 0)))

    a = BilinearForm(st_fes, check_unused=False, symmetric=False)
    for integrator in hasneg_integrators_a + patch_integrators_a:
        a += integrator

    f = LinearForm(st_fes)

    for integrator in hasneg_integrators_f:
        f += integrator

    while tend - told > delta_t / 2:
        SpaceTimeInterpolateToP1(levelset, tref, lset_p1)
        RestrictGFInTime(spacetime_gf=lset_p1,
                         reference_time=0.0,
                         space_gf=lset_bottom)
        RestrictGFInTime(spacetime_gf=lset_p1,
                         reference_time=1.0,
                         space_gf=lset_top)

        # update markers in (space-time) mesh
        ci.Update(lset_p1, time_order=time_order)

        # re-compute the facets for stabilization:
        ba_facets = GetFacetsWithNeighborTypes(mesh,
                                               a=ci.GetElementsOfType(HASNEG),
                                               b=ci.GetElementsOfType(IF))
        # re-evaluate the "active dofs" in the space time slab
        active_dofs = GetDofsOfElements(st_fes, ci.GetElementsOfType(HASNEG))

        # re-set definedonelements-markers according to new markings:
        for integrator in hasneg_integrators_a + hasneg_integrators_f:
            integrator.SetDefinedOnElements(ci.GetElementsOfType(HASNEG))
        for integrator in patch_integrators_a:
            integrator.SetDefinedOnElements(ba_facets)

        # assemble linear system
        a.Assemble()
        f.Assemble()

        # solve linear system
        inv = a.mat.Inverse(active_dofs, inverse="umfpack")
        gfu.vec.data = inv * f.vec

        # evaluate upper trace of solution for
        #  * for error evaluation
        #  * upwind-coupling to next time slab
        RestrictGFInTime(spacetime_gf=gfu, reference_time=1.0, space_gf=u_last)

        # update time variable (float and ParameterCL)
        told = told + delta_t
        coef_told.Set(told)

        # compute error at end of time slab
        l2error = sqrt(
            Integrate(lset_neg_top, (u_exactL(told) - u_last)**2, mesh))
        # print time and error
        print("t = {0:10}, l2error = {1:20}".format(told, l2error), end="\n")
        assert (l2error < 0.085)
Beispiel #21
0
from ngsolve import *
from netgen.geom2d import SplineGeometry

geo = SplineGeometry()
geo.AddRectangle((-0.1, -0.25), (0.1, 0.25),
                 leftdomain=0,
                 rightdomain=1,
                 bc="scatterer")
geo.AddCircle((0, 0), r=1, leftdomain=1, rightdomain=2)
geo.AddCircle((0, 0), r=1.4, leftdomain=2, rightdomain=0)

ngmesh = geo.GenerateMesh(maxh=0.04)
# ngmesh.Save("scattering.vol")
mesh = Mesh(ngmesh)
# mesh = Mesh ("scattering.vol")

mesh.SetPML(pml.Radial(origin=(0, 0), rad=1, alpha=0.1j), definedon=2)

kx = 50
ky = 20
k = sqrt(kx * kx + ky * ky)

uin = exp(1J * kx * x + 1J * ky * y)

fes = H1(mesh, complex=True, order=5, dirichlet="scatterer")
u = fes.TrialFunction()
v = fes.TestFunction()

uscat = GridFunction(fes)
uscat.Set(uin, definedon=mesh.Boundaries("scatterer"))
Beispiel #22
0
def Make2DProblem(maxh=2):
    from netgen.geom2d import SplineGeometry
    square = SplineGeometry()
    square.AddRectangle([-1, -1], [1, 1], bc=1)
    mesh = Mesh(square.GenerateMesh(maxh=maxh, quad_dominated=False))
    return mesh
            if IF_INT:
                ts = Symbol('ts')
                C = Curve([ts, -(a * ts + d) / (b + c * ts)], (ts, x0, x1))
                print(C)
                referencevals[IF] += line_integrate(f_py, C, [xs, ys])

    for i in range(0, len(part) - 1):
        add_integration_on_interval(part[i], part[i + 1])
    return referencevals


if __name__ == "__main__":
    from netgen.geom2d import SplineGeometry
    square = SplineGeometry()
    square.AddRectangle([0, 0], [1, 1], bc=1)
    #mesh = Mesh (square.GenerateMesh(maxh=100, quad_dominated=False))
    mesh = Mesh(square.GenerateMesh(maxh=100, quad_dominated=True))

    lsetvals_list = [[-0.18687, 0.324987, 0.765764, 0.48983],
                     [0.765764, 0.324987, -0.18687, -0.48983],
                     [1, 2 / 3, -1, -2 / 3]]
    #lsetvals_list = [[1.,-1.,-3.,-1.]]
    #lsetvals_list = [[1,-1,-4,-2]]
    lsetvals_list.append([3, -1, 1, -1.023123])

    f = lambda x, y: 1 + 0 * x + 0 * y
    f_ngs = f(x, y)
    V = H1(mesh, order=1)
    lset_approx = GridFunction(V)
Beispiel #24
0
# Finite element space polynomial order
order = 2
# diffusion coefficients for the sub-domains (NEG/POS):
alpha = [1.0, 2.0]
# Nitsche penalty parameter
lambda_nitsche = 20
# Write VTK output for visualisation
do_vtk = True


# ----------------------------------- MAIN ------------------------------------
# We generate the background mesh of the domain and use a simplicity
# triangulation to obtain a mesh with quadrilaterals use
# 'quad_dominated=True'
geo = SplineGeometry()
geo.AddRectangle(ll, ur, bc=1)

if rank == 0:
    ngmesh = geo.GenerateMesh(maxh=maxh)
    if np > 1:
        ngmesh.Distribute(comm)
else:
    ngmesh = netgen.meshing.Mesh.Receive(comm)
    ngmesh.SetGeometry(geo)

mesh = Mesh(ngmesh)

# Level set function of the domain and higher-order machinery
levelset = (sqrt(sqrt(x**4 + y**4)) - 1)

lsetadap = LevelSetMeshAdaptation(mesh, order=order, threshold=1000,
Beispiel #25
0
def test_nxfem(integration_quadtrig_newold,order):
    quad_dominated, newintegration = integration_quadtrig_newold
    square = SplineGeometry()
    square.AddRectangle([-1.5,-1.5],[1.5,1.5],bc=1)
    mesh = Mesh (square.GenerateMesh(maxh=0.2, quad_dominated=quad_dominated))
    
    r44 = (x*x*x*x+y*y*y*y)
    r41 = sqrt(sqrt(x*x*x*x+y*y*y*y))
    r4m3 = (1.0/(r41*r41*r41))
    r66 = (x*x*x*x*x*x+y*y*y*y*y*y)
    r63 = sqrt(r66)
    r22 = (x*x+y*y)
    r21 = sqrt(r22)
    
    levelset = (sqrt(sqrt(x*x*x*x+y*y*y*y)) - 1.0)

    solution = [1.0+pi/2.0-sqrt(2.0)*cos(pi/4.0*r44),pi/2.0*r41]
    alpha = [1.0,2.0]
    coef_f = [ (-1.0*sqrt(2.0)*pi*(pi*cos(pi/4*(r44))*(r66)+3*sin(pi/4*(r44))*(r22))),
              (-2.0*pi*3/2*(r4m3)*(-(r66)/(r44)+(r22))) ]
    
    order = order

    lsetmeshadap = LevelSetMeshAdaptation(mesh, order=order, threshold=1000, discontinuous_qn=True)
    deformation = lsetmeshadap.CalcDeformation(levelset)
    if (newintegration):
        lsetp1 = lsetmeshadap.lset_p1
        gradlsetp1 = grad(lsetp1)
    else:
        lsetp1 = 1.0*lsetmeshadap.lset_p1
        gradlsetp1 = 1.0*grad(lsetmeshadap.lset_p1)

    
    # extended FESpace 
    
    Vh = H1(mesh, order=order, dirichlet=[1,2,3,4])
    Vhx = XFESpace(Vh, lsetp1)
    VhG = FESpace([Vh,Vhx])
    print("unknowns in extended FESpace:", VhG.ndof)
    
    # coefficients / parameters: 
    
    n = 1.0/gradlsetp1.Norm() * gradlsetp1
    h = specialcf.mesh_size
    
    kappa = [CutRatioGF(Vhx.GetCutInfo()),1.0-CutRatioGF(Vhx.GetCutInfo())]
    
    stab = 10*(alpha[1]+alpha[0])*(order+1)*order/h
    
    # expressions of test and trial functions:
    
    u_std, u_x = VhG.TrialFunction()
    v_std, v_x = VhG.TestFunction()
    
    u = [u_std + op(u_x) for op in [neg,pos]]
    v = [v_std + op(v_x) for op in [neg,pos]]
    
    gradu = [grad(u_std) + op(u_x) for op in [neg_grad,pos_grad]]
    gradv = [grad(v_std) + op(v_x) for op in [neg_grad,pos_grad]]
    
    average_flux_u = sum([- kappa[i] * alpha[i] * gradu[i] * n for i in [0,1]])
    average_flux_v = sum([- kappa[i] * alpha[i] * gradv[i] * n for i in [0,1]])
    
    # integration domains (and integration parameter "subdivlvl" and "force_intorder")
    
    lset_neg = { "levelset" : lsetp1, "domain_type" : NEG, "subdivlvl" : 0}
    lset_pos = { "levelset" : lsetp1, "domain_type" : POS, "subdivlvl" : 0}
    lset_if  = { "levelset" : lsetp1, "domain_type" : IF , "subdivlvl" : 0}
    
    # bilinear forms:
    
    a = BilinearForm(VhG, symmetric = True, flags = { })
    a += SymbolicBFI(levelset_domain = lset_neg, form = alpha[0] * gradu[0] * gradv[0])
    a += SymbolicBFI(levelset_domain = lset_pos, form = alpha[1] * gradu[1] * gradv[1])
    a += SymbolicBFI(levelset_domain = lset_if , form =     average_flux_u * (v[0]-v[1]))
    a += SymbolicBFI(levelset_domain = lset_if , form =     average_flux_v * (u[0]-u[1]))
    a += SymbolicBFI(levelset_domain = lset_if , form = stab * (u[0]-u[1]) * (v[0]-v[1]))
    
    f = LinearForm(VhG)
    f += SymbolicLFI(levelset_domain = lset_neg, form = coef_f[0] * v[0])
    f += SymbolicLFI(levelset_domain = lset_pos, form = coef_f[1] * v[1])
    
    gfu = GridFunction(VhG)
    
    gfu.components[0].Set(solution[1], BND)
    
    mesh.SetDeformation(deformation)
    
    a.Assemble();
    f.Assemble();
    rhs = gfu.vec.CreateVector()
    rhs.data = f.vec - a.mat * gfu.vec
    update = gfu.vec.CreateVector()
    update.data = a.mat.Inverse(VhG.FreeDofs()) * rhs;
    gfu.vec.data += update
    
    sol_coef = IfPos(lsetp1,solution[1],solution[0])
    u_coef = gfu.components[0] + IfPos(lsetp1, pos(gfu.components[1]), neg(gfu.components[1]))
    u = [gfu.components[0] + op(gfu.components[1]) for op in [neg,pos]]
    
    err_sqr_coefs = [ (u[i] - solution[i])*(u[i] - solution[i]) for i in [0,1] ]
    
    l2error = sqrt(   Integrate(levelset_domain = lset_neg, cf=err_sqr_coefs[0], mesh=mesh, order=order, heapsize=1000000)
                    + Integrate(levelset_domain = lset_pos, cf=err_sqr_coefs[1], mesh=mesh, order=order, heapsize=1000000) )

    print("L2 error : ",l2error)

    mesh.UnsetDeformation()
    
    if (order == 1):
        assert l2error < 0.06
    if (order == 2):
        assert l2error < 0.004
    if (order == 3):
        assert l2error < 0.0004
Beispiel #26
0
    stokes interface problems. PAMM, 2016 
    [3]: A.Hansbo, P.Hansbo, An unfitted finite element method, based on Nitsche's method, for
    elliptic interface problems, Comp. Meth. Appl. Mech. Eng., 2002

"""

from netgen.geom2d import SplineGeometry
from netgen.meshing import MeshingParameters
from ngsolve import *
from xfem import *
# For LevelSetAdaptationMachinery
from xfem.lsetcurv import *
from math import pi
# 2D: circle configuration
square = SplineGeometry()
square.AddRectangle([-1, -1], [1, 1], bcs=[1, 2, 3, 4])
mesh = Mesh(square.GenerateMesh(maxh=4, quad_dominated=False))
for i in range(3):
    mesh.Refine()

order = 2  #=k in P(k)P(k-1)-Taylor-Hood pair

mu1 = 1.0
mu2 = 10.0
mu = [mu1, mu2]

R = 2.0 / 3.0
aneg = 1.0 / mu1
apos = 1.0 / mu2 + (1.0 / mu1 - 1.0 / mu2) * exp(x * x + y * y - R * R)
gammaf = 0.5  # surface tension = pressure jump
#q = gammaf - pi*R*R/4.0*gammaf
Beispiel #27
0
def test_spacetime_model_spacetime(pitfal1, pitfal2, pitfal3):
    square = SplineGeometry()
    square.AddRectangle([0, 0], [1, 1], bc=1)
    ngmesh = square.GenerateMesh(maxh=0.05, quad_dominated=False)
    mesh = Mesh(ngmesh)

    fes1 = V = H1(mesh, order=1, dirichlet=[1, 2, 3, 4])
    k_t = 1
    tfe = ScalarTimeFE(k_t)

    st_fes = SpaceTimeFESpace(fes1, tfe)
    st_fes_ic = SpaceTimeFESpace(fes1, tfe)

    tend = 1.0
    delta_t = 1 / 32

    told = Parameter(0)
    tref = ReferenceTimeVariable()
    t = told + delta_t * tref

    u_exact = lambda t: CoefficientFunction(
        sin(pi * t) * sin(pi * x) * sin(pi * x) * sin(pi * y) * sin(pi * y))
    coeff_f = CoefficientFunction(
        pi * cos(pi * t) * sin(pi * x) * sin(pi * x) * sin(pi * y) *
        sin(pi * y) - 2 * pi * pi * sin(pi * t) *
        (cos(pi * x) * cos(pi * x) * sin(pi * y) * sin(pi * y) -
         sin(pi * x) * sin(pi * x) * sin(pi * y) * sin(pi * y) +
         cos(pi * y) * cos(pi * y) * sin(pi * x) * sin(pi * x) -
         sin(pi * x) * sin(pi * x) * sin(pi * y) * sin(pi * y)))

    u0 = GridFunction(st_fes)
    u0_ic = GridFunction(fes1)
    u = st_fes.TrialFunction()
    v = st_fes.TestFunction()

    # dummy lset domain to call symboliccutbfi instead of usual symbolicbfi...
    levelset = (sqrt(x * x + y * y) - 1000.5)
    lsetp1 = GridFunction(H1(mesh, order=1))
    InterpolateToP1(levelset, lsetp1)
    lset_neg = {"levelset": lsetp1, "domain_type": NEG, "subdivlvl": 0}

    a = BilinearForm(st_fes, symmetric=False)
    a += SymbolicBFI(levelset_domain=lset_neg,
                     form=delta_t * grad(u) * grad(v),
                     time_order=2)
    a += SymbolicBFI(form=fix_tref(u, 0) * fix_tref(v, 0))
    a += SymbolicBFI(levelset_domain=lset_neg, form=dt(u) * v, time_order=2)
    a.Assemble()

    t_old = 0
    u0_ic.Set(u_exact(0))
    if pitfal1:
        u0_ic.Set(u_exact(t))

    while tend - t_old > delta_t / 2:
        f = LinearForm(st_fes)
        f += SymbolicLFI(levelset_domain=lset_neg,
                         form=delta_t * coeff_f * v,
                         time_order=2)
        f += SymbolicLFI(form=u0_ic * fix_tref(v, 0))
        if pitfal2:
            f += SymbolicLFI(form=u0_ic * v)
        f.Assemble()

        u0.vec.data = a.mat.Inverse(st_fes.FreeDofs(), "umfpack") * f.vec

        # exploiting the nodal property of the time fe:
        #u0_ic.vec[:] = u0.vec[0:fes1.ndof]
        u0_ic.vec[:].data = u0.vec[fes1.ndof:2 * fes1.ndof]

        t_old = t_old + delta_t
        told.Set(t_old)

        l2error = sqrt(Integrate((u_exact(t_old) - u0_ic)**2, mesh))
        if pitfal3:
            l2error = sqrt(Integrate((u_exact(t) - u0_ic)**2, mesh))

        print("t = {0}, l2error = {1}".format(t_old, l2error))
        assert l2error < 5e-3
    assert l2error < 2e-4
from netgen import gui
from ngsolve import *
from netgen.geom2d import SplineGeometry

geo = SplineGeometry()
geo.AddRectangle((-1, -1), (1, 1), bcs=("bottom", "right", "top", "left"))
mesh = Mesh(geo.GenerateMesh(maxh=0.25))
Draw(mesh)
fes = H1(mesh, order=3, dirichlet="bottom|top")

u, v = fes.TnT()

time = 0.0
dt = 0.01

gfu = GridFunction(fes)
gfuold = GridFunction(fes)

a = BilinearForm(fes, symmetric=False)
a += (u * v + dt * 3 * u**2 * grad(u) * grad(v) +
      dt * u**3 * grad(u) * grad(v) - dt * 1 * v - gfuold * v) * dx

from math import pi
gfu = GridFunction(fes)
gfu.Set(sin(2 * pi * x))
Draw(gfu, mesh, "u")
SetVisualization(deformation=True)
t = 0


def SolveNonlinearMinProblem(a, gfu, tol=1e-13, maxits=25):
# unfitted Heat equation with Neumann b.c.
from ngsolve import *
from time import sleep
from netgen.geom2d import unit_square
from netgen.geom2d import SplineGeometry
from netgen.meshing import MeshingParameters

from ngsolve.internal import *
from xfem import *
from numpy import pi

from xfem.lset_spacetime import *

square = SplineGeometry()
square.AddRectangle([-1, -0.75], [1, 1.5], bc=1)
ngmesh = square.GenerateMesh(maxh=10, quad_dominated=True)
mesh = Mesh(ngmesh)

for i in range(6):
    mesh.Refine()

fes1 = H1(mesh, order=1, dirichlet=[])
fes_lset_slice = H1(mesh, order=1, dirichlet=[])
k_t = 1
k_s = 1
lset_order_time = 2
tfe = ScalarTimeFE(k_t)

st_fes = SpaceTimeFESpace(fes1, tfe)
visoptions.deformation = 1
Beispiel #30
0
from ngsolve import *
from netgen.geom2d import SplineGeometry
from netgen.meshing import MeshingParameters


geo = SplineGeometry()
geo.AddRectangle( (0, 0), (1, 0.1), bcs = ("bottom", "right", "top", "left"))

mp = MeshingParameters(maxh=0.05)
mp.RestrictH(0,0,0, h=0.01)
mp.RestrictH(0,0.1,0, h=0.01)
mesh = Mesh( geo.GenerateMesh(mp=mp))



from ngs_templates.Elasticity import *

loadfactor = Parameter(0)

model = Elasticity(mesh=mesh, materiallaw=NeoHookeMaterial(200,0.2), \
                       # dirichlet="left",
                       # volumeforce=loadfactor*CoefficientFunction((0,1)), \
                       boundarydisplacement = { "left" : (0,0) },
                       boundaryforce = { "right" : loadfactor*(0,1) },
                       # boundarydisplacement = { "left" : (0,0), "right" : loadfactor*(0,1) },
                       nonlinear=True, order=4)

                       
Draw (model.displacement)
Draw (model.stress, mesh, "stress")
SetVisualization (deformation=True)