def test_polynomial_ET_Segm(domain, alpha): order = alpha m = ngm.Mesh() m.dim = 1 nel = 1 pnums = [] for i in range(0, nel+1): pnums.append (m.Add (ngm.MeshPoint (ngm.Pnt(i/nel, 0, 0)))) for i in range(0,nel): m.Add (ngm.Element1D ([pnums[i],pnums[i+1]], index=1)) m.Add (ngm.Element0D (pnums[0], index=1)) m.Add (ngm.Element0D (pnums[nel], index=2)) mesh = Mesh(m) x_ast = 0.78522 levelset = x_ast-x referencevals = {POS:pow(x_ast, alpha+1)/(alpha+1), NEG:(1-pow(x_ast, alpha+1))/(alpha+1), IF: pow(x_ast, alpha)} V = H1(mesh, order=1) lset_approx = GridFunction(V) #InterpolateToP1(levelset,lset_approx) lset_approx.Set(levelset) f = x**alpha integral = Integrate(levelset_domain = { "levelset" : lset_approx, "domain_type" : domain}, cf=f, mesh=mesh, order = order) print("Result of Integration Key ",domain," : ", integral) error = abs(integral - referencevals[domain]) print("Error: ", error) assert error < 5e-15*(order+1)*(order+1)
def test_convection1d_dg(): m = meshing.Mesh() m.dim = 1 nel = 20 pnums = [] for i in range(0, nel + 1): pnums.append(m.Add(meshing.MeshPoint(Pnt(i / nel, 0, 0)))) for i in range(0, nel): m.Add(meshing.Element1D([pnums[i], pnums[i + 1]], index=1)) m.Add(meshing.Element0D(pnums[0], index=1)) m.Add(meshing.Element0D(pnums[nel], index=2)) m.AddPointIdentification(pnums[0], pnums[nel], identnr=1, type=meshing.IdentificationType.PERIODIC) mesh = Mesh(m) fes = L2(mesh, order=4) u = fes.TrialFunction() v = fes.TestFunction() b = CoefficientFunction(1) bn = b * specialcf.normal(1) a = BilinearForm(fes) a += SymbolicBFI(-u * b * grad(v)) a += SymbolicBFI(bn * IfPos(bn, u, u.Other()) * v, element_boundary=True) u = GridFunction(fes) pos = 0.5 u0 = exp(-100 * (x - pos) * (x - pos)) u.Set(u0) w = u.vec.CreateVector() t = 0 tau = 1e-4 tend = 1 with TaskManager(): while t < tend - tau / 2: a.Apply(u.vec, w) fes.SolveM(rho=CoefficientFunction(1), vec=w) u.vec.data -= tau * w t += tau l2error = sqrt(Integrate((u - u0) * (u - u0), mesh)) print(l2error) assert l2error < 1e-2
def test_1dlaplace(): m = meshing.Mesh() m.dim = 1 nel = 20 pnums = [] for i in range(0, nel + 1): pnums.append(m.Add(meshing.MeshPoint(meshing.Pnt(i / nel, 0, 0)))) for i in range(0, nel): m.Add(meshing.Element1D([pnums[i], pnums[i + 1]], index=1)) m.Add(meshing.Element0D(pnums[0], index=1)) m.Add(meshing.Element0D(pnums[nel], index=2)) mesh = Mesh(m) order = 10 fes = H1(mesh, order=order, dirichlet=[1, 2]) u = fes.TrialFunction() v = fes.TestFunction() n = specialcf.normal(mesh.dim) h = specialcf.mesh_size ugf = GridFunction(fes) uex = GridFunction(fes) uex.Set(sin(pi * x)) a = BilinearForm(fes) a += SymbolicBFI(grad(u) * grad(v)) l = LinearForm(fes) l += SymbolicLFI(pi * pi * uex * v) l.Assemble() a.Assemble() invmat = a.mat.Inverse(fes.FreeDofs()) ugf.vec.data = invmat * l.vec error = sqrt(Integrate((ugf - uex) * (ugf - uex), mesh)) assert error <= 1e-14
def MakeMesh1D(): a = 0 b = 1 nel = 10 m = msh.Mesh() m.dim = 1 pnums = [] for i in range(nel + 1): pnums.append(m.Add(msh.MeshPoint(msh.Pnt(a + (b - a) * i / nel, 0, 0)))) for i in range(2): m.Add(msh.Element1D([pnums[i], pnums[i + 1]], index=1)) for i in range(2, 5): m.Add(msh.Element1D([pnums[i], pnums[i + 1]], index=3)) for i in range(5, 7): m.Add(msh.Element1D([pnums[i], pnums[i + 1]], index=1)) for i in range(7, nel): m.Add(msh.Element1D([pnums[i], pnums[i + 1]], index=2)) m.SetMaterial(1, 'base') m.SetMaterial(2, 'chip') m.SetMaterial(3, 'top') # add points m.Add(msh.Element0D(pnums[0], index=1)) m.Add(msh.Element0D(pnums[2], index=2)) m.Add(msh.Element0D(pnums[5], index=2)) m.Add(msh.Element0D(pnums[7], index=2)) m.Add(msh.Element0D(pnums[nel], index=2)) # set boundary condition names m.SetBCName(0, 'bottom') m.SetBCName(1, 'other') m.Save('temp.vol') return m
x_min = -50 x_max = 50 length = x_max - x_min pnums = [] xs = [] for i in range(0, n_elems + 1): pnt_x = x_min + length * i / n_elems xs.append(pnt_x) pnums.append(m.Add(ngm.MeshPoint(ngm.Pnt(pnt_x, 0, 0)))) for i in range(0, n_elems): m.Add(ngm.Element1D([pnums[i], pnums[i + 1]], index=1)) m.SetMaterial(1, 'material') m.Add(ngm.Element0D(pnums[0], index=1)) m.Add(ngm.Element0D(pnums[n_elems], index=2)) ## NGSolve mesh = ngs.Mesh(m) fes = ngs.H1(mesh, order=1, dirichlet=[1, 2], complex=True) u = fes.TrialFunction() v = fes.TestFunction() ## Potentials ### Potential barrier barrier_w = 2 barrier_h = 1 potential = ngs.CoefficientFunction( ngs.IfPos(x, barrier_h, 0) - ngs.IfPos(x - barrier_w, barrier_h, 0))