-
Notifications
You must be signed in to change notification settings - Fork 1
/
fenics_manufactured_solution.py
146 lines (122 loc) · 7.19 KB
/
fenics_manufactured_solution.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
#!/usr/bin/python3
"""
Manufactured solution by explicitly setting the solution displacement functions to:
| 1e-3 * z *exp(x) |
u (x, y, z) = | 1e-3 * z *exp(y) |
| 1e-3 * z *exp(z) |
with (x, y, z) material coordinates (undeformed, at rest). The external forces applied to the domain becomes:
ext = integrate( Div(P), volume_domain) + integrate( P . N , neummann_domain)
where P is the first Piola-Kirchhoff stress tensor, N is the surface normal at rest. The geometry of the volumetric
domain is a cylinder having a radius of 1 and a length of 3.
To run with docker :
docker run --rm -v $PWD:/w -w /w jnbrunet/caribou_validation python manufactured_solution.py
"""
import Sofa, SofaRuntime, SofaCaribou
import meshio, numpy as np
from manufactured_solution import assemble, integrate, compute_solution, ConstantForceField
order = "quadratic"
element = "tetrahedron"
material = "SaintVenantKirchhoff"
if element == "tetrahedron":
if order == "linear":
# mesh = meshio.read('meshes/refined_meshes/beam/p1/beam_p1_9.vtu')
mesh = meshio.read('meshes/refined_meshes/cylinder/p1/cylinder_p1_5.0.vtu')
surface_template = "Triangle"
volume_template = "Tetrahedron"
fenics_volume_indices = mesh.cells[1].data
else:
# mesh = meshio.read('meshes/refined_meshes/beam/p2/beam_p2_9.vtu')
mesh = meshio.read('meshes/refined_meshes/cylinder/p2/cylinder_p2_4.0.vtu')
surface_template = "Triangle6"
volume_template = "Tetrahedron10"
fenics_volume_indices = mesh.cells[1].data[:, [0, 1, 2, 3, 9, 8, 5, 7, 6, 4]]
else:
if order == "linear":
mesh = meshio.read('meshes/refined_meshes/cylinder/q1/cylinder_q1_5.0.vtu')
surface_template = "Quad"
volume_template = "Hexahedron"
fenics_volume_indices = mesh.cells[1].data[:, [4, 5, 0, 1, 7, 6, 3, 2]]
else:
mesh = meshio.read('meshes/refined_meshes/cylinder/p2/cylinder_p2_5.0.vtu')
surface_template = "Triangle6"
volume_template = "Tetrahedron10"
fenics_volume_indices = mesh.cells[1].data[:, [0, 1, 2, 3, 9, 8, 5, 7, 6, 4]]
mu = 1.0
l = 1.25
rad = 7.5
length = 80
poisson_ratio = 1. / (2 * ((mu / l) + 1))
young_modulus = 2 * mu * (1 + poisson_ratio)
P, f, u_s = compute_solution(mu, l, rad, length)
def create_scene(root):
root.addObject('RequiredPlugin',
pluginName='SofaBaseMechanics SofaBoundaryCondition SofaImplicitOdeSolver SofaSparseSolver')
root.addObject('RequiredPlugin', name='SofaCaribou')
root.addObject('VisualStyle', displayFlags='showBehaviorModels showForceFields')
sofa_node = root.addChild("sofa_node")
sofa_node.addObject('StaticODESolver', newton_iterations=10, residual_tolerance_threshold=1e-10, printLog=True)
sofa_node.addObject('LDLTSolver', backend='Pardiso')
sofa_node.addObject('MechanicalObject', name='mo', position=mesh.points.tolist())
sofa_node.addObject('CaribouTopology', name='volume', template=volume_template, indices=mesh.cells[1].data.tolist())
sofa_node.addObject('CaribouTopology', name='dirichlet_boundary', template=surface_template,
indices=mesh.cells[0].data[
np.ma.masked_equal(mesh.cell_data['gmsh:physical'][0], 1).mask].tolist())
sofa_node.addObject('CaribouTopology', name='neumann_boundary', template=surface_template,
indices=mesh.cells[0].data[
np.ma.masked_inside(mesh.cell_data['gmsh:physical'][0], 2, 3).mask].tolist())
sofa_node.addObject(material + "Material", young_modulus=young_modulus, poisson_ratio=poisson_ratio)
sofa_node.addObject('HyperelasticForcefield', topology='@volume', printLog=True)
sofa_node.addObject('FixedConstraint',
indices=np.unique(root.sofa_node.dirichlet_boundary.indices.array()).tolist())
sofa_node.addObject(ConstantForceField(name='external_forces_sofa'))
fenics_node = root.addChild("fenics_node")
fenics_node.addObject('StaticODESolver', newton_iterations=10, residual_tolerance_threshold=1e-10, printLog=True)
fenics_node.addObject('LDLTSolver', backend='Pardiso')
fenics_node.addObject('MechanicalObject', name='mo', position=mesh.points.tolist())
fenics_node.addObject('CaribouTopology', name='volume', template=volume_template,
indices=fenics_volume_indices.tolist())
fenics_node.addObject('CaribouTopology', name='dirichlet_boundary', template=surface_template,
indices=mesh.cells[0].data[
np.ma.masked_equal(mesh.cell_data['gmsh:physical'][0], 1).mask].tolist())
fenics_node.addObject('CaribouTopology', name='neumann_boundary', template=surface_template,
indices=mesh.cells[0].data[
np.ma.masked_inside(mesh.cell_data['gmsh:physical'][0], 2, 3).mask].tolist())
fenics_node.addObject('FEniCS_Material', template=volume_template, young_modulus=young_modulus,
poisson_ratio=poisson_ratio, material_name=material)
fenics_node.addObject('HyperelasticForcefield_FEniCS', topology='@volume', printLog=True)
fenics_node.addObject('FixedConstraint',
indices=np.unique(root.fenics_node.dirichlet_boundary.indices.array()).tolist())
fenics_node.addObject(ConstantForceField(name='external_forces_fenics'))
if __name__ == '__main__':
root = Sofa.Core.Node()
create_scene(root)
Sofa.Simulation.init(root)
print('Assembling the external force vector...', end='', flush=True)
external_forces = \
assemble(root.sofa_node.volume.domain(), lambda x, y, z, _: f(x, y, z)) + \
assemble(root.sofa_node.neumann_boundary.domain(), lambda x, y, z, t: np.dot(P(x, y, z), t.normal()))
print(' Done.', flush=True)
exact_error = np.sqrt(
integrate(root.sofa_node.volume.domain(), lambda x, y, z, _: np.dot(u_s(x, y, z), u_s(x, y, z))))
print(f"Exact error is {exact_error}")
for load in [1e-3, 1e-2, 1e-1, 0.15, 0.5, 1.0]:
# for load in [1e-3]:
root.sofa_node.external_forces_sofa.forces = (external_forces * load)
root.fenics_node.external_forces_fenics.forces = (external_forces * load)
Sofa.Simulation.animate(root, 1)
u_h_sofa = (root.sofa_node.mo.position.array() - root.sofa_node.mo.rest_position.array())
error_L2_sofa = np.sqrt(integrate(
root.sofa_node.volume.domain(),
lambda x, y, z, u_g, _: np.dot((u_s(x, y, z) - u_g), (u_s(x, y, z) - u_g)),
u_h_sofa
))
print(f"SOFA relative L2 error at {int(load * 100)}% load: {error_L2_sofa / exact_error}")
u_h_fenics = (root.fenics_node.mo.position.array() - root.fenics_node.mo.rest_position.array())
error_L2_fenics = np.sqrt(integrate(
root.sofa_node.volume.domain(),
lambda x, y, z, u_g, _: np.dot((u_s(x, y, z) - u_g), (u_s(x, y, z) - u_g)),
u_h_fenics
))
print(f"FEniCS relative L2 error at {int(load * 100)}% load: {error_L2_fenics / exact_error}")
print(
f"Difference in relative L2 error at {int(load * 100)}% load: {np.abs((error_L2_fenics / exact_error) - (error_L2_sofa / exact_error))}")