- Thank you received: 0
Trialfunction does not support BND-forms, maybe a Trace() operator is misssing
6 years 1 month ago - 6 years 1 month ago #1198
by noname
Trialfunction does not support BND-forms, maybe a Trace() operator is misssing was created by noname
Hello everyone ,
I am trying to calculate an integral on the dirichlet boundary, however for some reason I keep getting an error saying "Trialfunction does not support BND-forms, maybe a Trace() operator is misssing"
My geometry is and function space is:
and I want to impose Dirichlet boundary conditions only on the "t" boundary weakly. For that I need to calculate following:
but this calculation gives me the above-mentioned error.
Anyone can see what is going on here?
Thanks in advance.
I am trying to calculate an integral on the dirichlet boundary, however for some reason I keep getting an error saying "Trialfunction does not support BND-forms, maybe a Trace() operator is misssing"
My geometry is and function space is:
Code:
from netgen.geom2d import SplineGeometry
geo = SplineGeometry()
geo.AddRectangle( (0, 0), (0.41, 2), bcs = ("b", "r", "t", "l"))
geo.AddCircle ( (0.2, 1.8), r=0.05, leftdomain=0, rightdomain=1, bc="cyl")
mesh = Mesh(geo.GenerateMesh(maxh=0.08))
mesh.Curve(3)
V = VectorH1(mesh,order=3, dirichlet="t|cyl|l|r")
Q = H1(mesh,order=2)
X = FESpace([V,Q])
u,p = X.TrialFunction()
v,q = X.TestFunction()
h = specialcf.mesh_size
n = specialcf.normal(2)
and I want to impose Dirichlet boundary conditions only on the "t" boundary weakly. For that I need to calculate following:
Code:
a += SymbolicBFI(-nu*InnerProduct(grad(u)*n,v)
-nu*InnerProduct(grad(v)*n,u)
+nu*alpha/h*InnerProduct(u,v),
definedon=mesh.Boundaries("t"))
but this calculation gives me the above-mentioned error.
Anyone can see what is going on here?
Thanks in advance.
Last edit: 6 years 1 month ago by noname.
6 years 1 month ago #1200
by schruste
Replied by schruste on topic Trialfunction does not support BND-forms, maybe a Trace() operator is misssing
Hi noname,
You want the gradient of volume elements (aligned to boundary elements) evaluated. The SymbolicBFI on the boundary however only has the information of the boundary elements. Use the option "skeleton=True" combined with "VOL_or_BND=BND" to obtain access to the gradients of the aligned volume elements.
Best,
Christoph
You want the gradient of volume elements (aligned to boundary elements) evaluated. The SymbolicBFI on the boundary however only has the information of the boundary elements. Use the option "skeleton=True" combined with "VOL_or_BND=BND" to obtain access to the gradients of the aligned volume elements.
Best,
Christoph
The following user(s) said Thank You: noname
Time to create page: 0.093 seconds