-
Notifications
You must be signed in to change notification settings - Fork 298
/
ParallelHexahedronFEMForceField.inl
270 lines (223 loc) · 11.2 KB
/
ParallelHexahedronFEMForceField.inl
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
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
/******************************************************************************
* SOFA, Simulation Open-Framework Architecture *
* (c) 2006 INRIA, USTL, UJF, CNRS, MGH *
* *
* This program is free software; you can redistribute it and/or modify it *
* under the terms of the GNU Lesser General Public License as published by *
* the Free Software Foundation; either version 2.1 of the License, or (at *
* your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, but WITHOUT *
* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or *
* FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License *
* for more details. *
* *
* You should have received a copy of the GNU Lesser General Public License *
* along with this program. If not, see <http://www.gnu.org/licenses/>. *
*******************************************************************************
* Authors: The SOFA Team and external contributors (see Authors.txt) *
* *
* Contact information: contact@sofa-framework.org *
******************************************************************************/
#pragma once
#include <MultiThreading/component/solidmechanics/fem/elastic/ParallelHexahedronFEMForceField.h>
#include <sofa/simulation/TaskScheduler.h>
#include <sofa/simulation/MainTaskSchedulerFactory.h>
#include <sofa/simulation/ParallelForEach.h>
namespace multithreading::component::forcefield::solidmechanics::fem::elastic
{
template<class DataTypes>
void ParallelHexahedronFEMForceField<DataTypes>::init()
{
Inherit1::init();
initTaskScheduler();
const auto indexedElements = *this->getIndexedElements();
m_vertexIdInAdjacentHexahedra.resize(this->l_topology->getNbPoints());
m_around.clear();
for (sofa::Size i = 0; i < this->l_topology->getNbPoints(); ++i)
{
const auto& around = this->l_topology->getHexahedraAroundVertex(i);
m_around.push_back(around);
sofa::Size j {};
for (const auto hexaId : around)
{
const auto element = indexedElements[hexaId];
sofa::Size indexInElement {};
for (const auto v : element)
{
if (v == i)
{
m_vertexIdInAdjacentHexahedra[i][j] = indexInElement;
}
++indexInElement;
}
++j;
}
}
}
template<class DataTypes>
void ParallelHexahedronFEMForceField<DataTypes>::addForce(const sofa::core::MechanicalParams* mparams, DataVecDeriv& f,
const DataVecCoord& p, const DataVecDeriv& v)
{
using namespace sofa::component::solidmechanics::fem::elastic;
if (this->method != HexahedronFEMForceField<DataTypes>::LARGE)
{
static bool firstTime = true;
msg_warning_when(firstTime) << "Multithreading is only partially supported for other methods than 'large'";
firstTime = false;
HexahedronFEMForceField<DataTypes>::addForce(mparams, f, p, v);
return;
}
WDataRefVecDeriv _f = f;
RDataRefVecCoord _p = p;
_f.resize(_p.size());
if (this->needUpdateTopology)
{
this->reinit();
this->needUpdateTopology = false;
}
const auto* indexedElements = this->getIndexedElements();
this->m_potentialEnergy = 0;
const auto& elementStiffnesses = this->_elementStiffnesses.getValue();
updateStiffnessMatrices = this->f_updateStiffnessMatrix.getValue();
if (updateStiffnessMatrices)
{
static bool first = true;
msg_warning_when(first) << "Updating the stiffness matrix is not fully supported in a multithreaded context. "
"It may lead to a crash. "
"Set 'updateStiffnessMatrix' to false to remove this warning.";
first = false;
}
std::mutex mutex;
sofa::simulation::parallelForEachRange(*m_taskScheduler,
indexedElements->begin(), indexedElements->end(),
[this, &_p, &elementStiffnesses, &mutex, &_f](const auto& range)
{
auto elementId = std::distance(this->getIndexedElements()->begin(), range.start);
std::vector<sofa::type::Vec<8, Deriv>> fElements;
fElements.reserve(std::distance(range.start, range.end));
SReal potentialEnergy { 0_sreal };
for (auto it = range.start; it != range.end; ++it, ++elementId)
{
sofa::type::Vec<8, Deriv> forceInElement;
this->computeTaskForceLarge(_p, elementId, *it, elementStiffnesses, potentialEnergy, forceInElement);
fElements.emplace_back(forceInElement);
}
std::lock_guard guard(mutex);
this->m_potentialEnergy += potentialEnergy;
auto it = range.start;
for (const auto& forceInElement : fElements)
{
for (int w = 0; w < 8; ++w)
{
_f[(*it)[w]] += forceInElement[w];
}
++it;
}
});
this->m_potentialEnergy/=-2.0;
}
template<class DataTypes>
void ParallelHexahedronFEMForceField<DataTypes>::computeTaskForceLarge(RDataRefVecCoord &p,
sofa::Index elementId,
const Element& elem,
const VecElementStiffness& elementStiffnesses,
SReal& OutPotentialEnery,
sofa::type::Vec<8, Deriv>& OutF)
{
sofa::type::Vec<8,Coord> nodes;
for(int w=0; w<8; ++w)
nodes[w] = p[elem[w]];
Coord horizontal;
horizontal = (nodes[1]-nodes[0] + nodes[2]-nodes[3] + nodes[5]-nodes[4] + nodes[6]-nodes[7])*.25;
Coord vertical;
vertical = (nodes[3]-nodes[0] + nodes[2]-nodes[1] + nodes[7]-nodes[4] + nodes[6]-nodes[5])*.25;
// Transformation R_0_2; // Rotation matrix (deformed and displaced Tetrahedron/world)
this->computeRotationLarge(this->_rotations[elementId], horizontal, vertical);
// positions of the deformed and displaced Tetrahedron in its frame
sofa::type::Vec<8,Coord> deformed;
for(int w=0; w<8; ++w)
deformed[w] = this->_rotations[elementId] * nodes[w];
// displacement
sofa::type::Vec<24, Real> D;
for(int k=0 ; k<8 ; ++k )
{
const int index = k*3;
for(int j=0 ; j<3 ; ++j )
D[index+j] = this->_rotatedInitialElements[elementId][k][j] - deformed[k][j];
}
if(updateStiffnessMatrices)
{
//Not thread safe: a message warned the user earlier that updating the stiffness matrices is not fully supported
this->computeElementStiffness((*this->_elementStiffnesses.beginEdit())[elementId], this->_materialsStiffnesses[elementId], deformed, elementId, this->_sparseGrid ? this->_sparseGrid->getStiffnessCoef(elementId) : 1.0 );
}
sofa::type::Vec<24, Real> F; //forces
this->computeForce( F, D, elementStiffnesses[elementId] ); // compute force on element
for(int w=0; w<8; ++w)
OutF[w] += this->_rotations[elementId].multTranspose(Deriv(F[w * 3], F[w * 3 + 1], F[w * 3 + 2] ) );
OutPotentialEnery += dot(Deriv( F[0], F[1], F[2] ) ,-Deriv( D[0], D[1], D[2]));
OutPotentialEnery += dot(Deriv( F[3], F[4], F[5] ) ,-Deriv( D[3], D[4], D[5] ));
OutPotentialEnery += dot(Deriv( F[6], F[7], F[8] ) ,-Deriv( D[6], D[7], D[8] ));
OutPotentialEnery += dot(Deriv( F[9], F[10], F[11]),-Deriv( D[9], D[10], D[11] ));
OutPotentialEnery += dot(Deriv( F[12], F[13], F[14]),-Deriv( D[12], D[13], D[14] ));
OutPotentialEnery += dot(Deriv( F[15], F[16], F[17]),-Deriv( D[15], D[16], D[17] ));
OutPotentialEnery += dot(Deriv( F[18], F[19], F[20]),-Deriv( D[18], D[19], D[20] ));
OutPotentialEnery += dot(Deriv( F[21], F[22], F[23]),-Deriv( D[21], D[22], D[23] ));
}
template<class DataTypes>
void ParallelHexahedronFEMForceField<DataTypes>::addDForce (const sofa::core::MechanicalParams *mparams, DataVecDeriv& v, const DataVecDeriv& x)
{
WDataRefVecDeriv _df = v;
RDataRefVecCoord _dx = x;
const Real kFactor = (Real)sofa::core::mechanicalparams::kFactorIncludingRayleighDamping(mparams, this->rayleighStiffness.getValue());
if (_df.size() != _dx.size())
_df.resize(_dx.size());
const auto& elementStiffnesses = this->_elementStiffnesses.getValue();
const auto& indexedElements = *this->getIndexedElements();
m_elementsDf.resize(indexedElements.size());
sofa::simulation::parallelForEachRange(*m_taskScheduler,
indexedElements.begin(), indexedElements.end(),
[this, &_dx, &elementStiffnesses, kFactor, &indexedElements](const auto& range)
{
auto elementId = std::distance(indexedElements.begin(), range.start);
auto elementsDfIt = m_elementsDf.begin() + elementId;
auto elementStiffnessesIt = elementStiffnesses.begin() + elementId;
auto rotationIt = this->_rotations.begin() + elementId;
for (auto it = range.start; it != range.end; ++it, ++elementId)
{
sofa::type::Vec<24, Real> X(sofa::type::NOINIT); //displacement
sofa::type::Vec<24, Real> F(sofa::type::NOINIT); //force
const auto& r = *rotationIt++;
const auto& element = *it;
for (sofa::Size w = 0; w < 8; ++w)
{
const Coord x_2 = r * _dx[element[w]];
X[w * 3] = x_2[0];
X[w * 3 + 1] = x_2[1];
X[w * 3 + 2] = x_2[2];
}
// F = K * X
this->computeForce(F, X, *elementStiffnessesIt++);
sofa::type::Vec<8, Deriv>& df = *elementsDfIt++;
for (sofa::Size w = 0; w < 8; ++w)
{
df[w] = -r.multTranspose(Deriv(F[w * 3], F[w * 3 + 1], F[w * 3 + 2])) * kFactor;
}
}
});
sofa::simulation::parallelForEachRange(*m_taskScheduler,
static_cast<std::size_t>(0), _df.size(), [&_df, this](const auto& range)
{
for (auto vertexId = range.start; vertexId < range.end; ++vertexId)
{
const auto& around = m_around[vertexId];
sofa::Size hexaAroundId {};
for (const auto hexaId : around)
{
_df[vertexId] += m_elementsDf[hexaId][m_vertexIdInAdjacentHexahedra[vertexId][hexaAroundId]];
hexaAroundId++;
}
}
});
}
} //namespace sofa::component::forcefield