You signed in with another tab or window. Reload to refresh your session.You signed out in another tab or window. Reload to refresh your session.You switched accounts on another tab or window. Reload to refresh your session.Dismiss alert
This is not a classical code-bug, but it could be a bug from a user perspective.
Summary (for v7.19.1)
Adopt the README code-block to the MINI-element formulation ✅
Force-displacement curve is realistic ✅
Stresses at quadrature points are straightforward ✅
Stress-projection is error-prone for MINI-elements ⚠️
Possible Solutions
Add a new argument to RegionTetraMINI(mesh, project=True) which raises the order of the quadrature scheme from 2 to 5!? This is a manual solution.
Raise the quadrature schemes of the template regions to always produce correct results for projection? This affects the performance (more computationally expensive assembly)!
Fall-back to tools.extrapolate(mean=True) for the projection of cell-means? --> This seems to be very promising compared to the results obtained with hexahedrons, see Stress Projection for TetraMINI is complicated #690 (comment)! mean=True should enforce extrapolate() independent of the shape of the input array. The project method should enforce regions with quadrature schemes which are good for projection and if the shape of the input array is incompatible then an error is raised. This is implemented in PR Fall-back to tools.extrapolate(mean=True) for project(mean=True) #691
Detailed Description
The first part is easy...
importfelupeasfemfrompypardisoimportspsolvemesh=fem.Cube(n=(16, 6, 6)).triangulate(3).add_midpoints_volumes(cell_type="tetra")
region=fem.RegionTetraMINI(mesh, fem.TetrahedronQuadrature(order=2))
field=fem.FieldsMixed(region)
boundaries, loadcase=fem.dof.uniaxial(field, clamped=True)
nh=fem.NeoHooke(mu=1)
solid=fem.SolidBody(fem.NearlyIncompressible(nh, bulk=500), field)
move=fem.math.linsteps([0, 2], num=3)
step=fem.Step(items=[solid], ramp={boundaries["move"]: move}, boundaries=boundaries)
job=fem.CharacteristicCurve(steps=[step], boundary=boundaries["move"])
job.evaluate(solver=spsolve, parallel=True)
fig, ax=job.plot(
xlabel="Displacement $u$ in mm $\longrightarrow$",
ylabel="Normal Force $F$ in N $\longrightarrow$",
)
...this is more or less the README code-block ✅. Obtaining the Cauchy-stress at quadrature points is also straightforward ✅. But the projection produces weird stresses ⚠️ at the midvolume-points.
importnumpyasnp# stress at quadrature pointscauchy=solid.evaluate.cauchy_stress(field)
# project the stress to mesh-pointsstress=fem.project(cauchy, region)
stress[-mesh.ncells:] =np.nan# view cell-means of stress xxview=solid.view(cell_data=dict(Stress=np.mean(cauchy, axis=-2).T))
view.plot("Stress", show_undeformed=False).show()
# view projected point-data of stress xxview=solid.view(point_data=dict(Stress=stress))
view.plot("Stress", show_undeformed=False).show()
So these points have to be ignored 🗑️. But the result is still very different to the cell-data (cell-means of stresses).
# ignore projected stresses at midvolume-pointsstress[-mesh.ncells:] =np.nan# set displacements of mid-volume points to zerofield[0].values[-mesh.ncells:] =0# view projected point-data of stress xxview=solid.view(point_data=dict(Stress=stress))
view.plot("Stress", show_undeformed=False).show()
For RegionTetraMINI a quadrature scheme with order=5 has to be used to project values to mesh-points (order=2 is too low and the next implemented order is fifth-order). Repeat the above code-blocks with:
Another option would be to project the cell-means to the mesh-points (also use a fifth-order quadrature scheme!). This is probably not a good idea.
# project the stress to mesh-pointsstress=fem.project(cauchy, region, mean=True)
stress[-mesh.ncells:] =np.nan# view projected point-data of stress xxview=solid.view(point_data=dict(Stress=stress))
view.plot("Stress", show_undeformed=False).show()
A similar task could also be done directly within the PyVista plotter. This seems to be the easiest operation. Same results may be obtained by extrapolate(mean=True).
view=solid.view(cell_data=dict(Stress=np.mean(cauchy, axis=-2).T))
view.mesh=view.mesh.cell_data_to_point_data()
view.plot("Stress", show_undeformed=False).show()
# alternative with equal resultstress=fem.tools.extrapolate(cauchy, region, mean=True) # this is the old project-methodview=solid.view(point_data=dict(Stress=stress))
view.plot("Stress", show_undeformed=False).show()
Results obtained by a finer mesh:
(extrapolate cell-mean to point-data)
I have absolutely no idea how to simplify or make the stress projection for MINI-elements less error-prone. In the long term, we could separate global DOF from private DOF and then only the global DOF are considered for projection / plotting. Triangle and Tetra - Regions should be initiated with a higher-order quadrature scheme if it is used for interpolation or projection (and not only gradient evaluations). Something like Region(project=True) instead of manually setting the quadrature scheme.
The text was updated successfully, but these errors were encountered:
This is not a classical code-bug, but it could be a bug from a user perspective.
Summary (for v7.19.1)
Possible Solutions
Add a new argument toRegionTetraMINI(mesh, project=True)
which raises the order of the quadrature scheme from 2 to 5!? This is a manual solution.Raise the quadrature schemes of the template regions to always produce correct results for projection? This affects the performance (more computationally expensive assembly)!tools.extrapolate(mean=True)
for the projection of cell-means?-->
This seems to be very promising compared to the results obtained with hexahedrons, see Stress Projection forTetraMINI
is complicated #690 (comment)!mean=True
should enforceextrapolate()
independent of the shape of the input array. The project method should enforce regions with quadrature schemes which are good for projection and if the shape of the input array is incompatible then an error is raised. This is implemented in PR Fall-back totools.extrapolate(mean=True)
forproject(mean=True)
#691Detailed Description
The first part is easy...
...this is more or less the README code-block ✅. Obtaining the Cauchy-stress at quadrature points is also straightforward ✅. But the projection produces weird stresses⚠️ at the midvolume-points.
So these points have to be ignored 🗑️. But the result is still very different to the cell-data (cell-means of stresses).
For
RegionTetraMINI
a quadrature scheme withorder=5
has to be used to project values to mesh-points (order=2
is too low and the next implemented order is fifth-order). Repeat the above code-blocks with:Another option would be to project the cell-means to the mesh-points (also use a fifth-order quadrature scheme!).This is probably not a good idea.A similar task could also be done directly within the PyVista plotter. This seems to be the easiest operation. Same results may be obtained by
extrapolate(mean=True)
.Results obtained by a finer mesh:
(extrapolate cell-mean to point-data)
I have absolutely no idea how to simplify or make the stress projection for MINI-elements less error-prone. In the long term, we could separate global DOF from private DOF and then only the global DOF are considered for projection / plotting. Triangle and Tetra - Regions should be initiated with a higher-order quadrature scheme if it is used for interpolation or projection (and not only gradient evaluations). Something like
Region(project=True)
instead of manually setting the quadrature scheme.The text was updated successfully, but these errors were encountered: