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 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431
|
/******************************************************************************
* SOFA, Simulation Open-Framework Architecture, version 1.0 beta 4 *
* (c) 2006-2009 MGH, INRIA, USTL, UJF, CNRS *
* *
* This library 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 library 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 library; if not, write to the Free Software Foundation, *
* Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. *
*******************************************************************************
* SOFA :: Modules *
* *
* Authors: The SOFA Team and external contributors (see Authors.txt) *
* *
* Contact information: contact@sofa-framework.org *
******************************************************************************/
#ifndef SOFA_GPU_CUDA_CUDASPRINGFORCEFIELD_INL
#define SOFA_GPU_CUDA_CUDASPRINGFORCEFIELD_INL
#include "CudaSpringForceField.h"
#include <sofa/component/forcefield/SpringForceField.inl>
#include <sofa/component/forcefield/StiffSpringForceField.inl>
#include <sofa/component/forcefield/MeshSpringForceField.inl>
#include <sofa/component/forcefield/TriangleBendingSprings.inl>
#include <sofa/component/forcefield/QuadBendingSprings.inl>
namespace sofa
{
namespace gpu
{
namespace cuda
{
extern "C"
{
void SpringForceFieldCuda3f_addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v);
void SpringForceFieldCuda3f_addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2);
void StiffSpringForceFieldCuda3f_addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v, void* dfdx);
void StiffSpringForceFieldCuda3f_addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2, void* dfdx);
void StiffSpringForceFieldCuda3f_addDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* dx, const void* x, const void* dfdx, double factor);
void StiffSpringForceFieldCuda3f_addExternalDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* dx1, const void* x1, const void* dx2, const void* x2, const void* dfdx, double factor);
void SpringForceFieldCuda3f1_addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v);
void SpringForceFieldCuda3f1_addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2);
void StiffSpringForceFieldCuda3f1_addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v, void* dfdx);
void StiffSpringForceFieldCuda3f1_addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2, void* dfdx);
void StiffSpringForceFieldCuda3f1_addDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* dx, const void* x, const void* dfdx, double factor);
void StiffSpringForceFieldCuda3f1_addExternalDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* dx1, const void* x1, const void* dx2, const void* x2, const void* dfdx, double factor);
#ifdef SOFA_GPU_CUDA_DOUBLE
void SpringForceFieldCuda3d_addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v);
void SpringForceFieldCuda3d_addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2);
void StiffSpringForceFieldCuda3d_addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v, void* dfdx);
void StiffSpringForceFieldCuda3d_addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2, void* dfdx);
void StiffSpringForceFieldCuda3d_addDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* dx, const void* x, const void* dfdx, double factor);
void StiffSpringForceFieldCuda3d_addExternalDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* dx1, const void* x1, const void* dx2, const void* x2, const void* dfdx, double factor);
void SpringForceFieldCuda3d1_addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v);
void SpringForceFieldCuda3d1_addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2);
void StiffSpringForceFieldCuda3d1_addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v, void* dfdx);
void StiffSpringForceFieldCuda3d1_addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2, void* dfdx);
void StiffSpringForceFieldCuda3d1_addDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* dx, const void* x, const void* dfdx, double factor);
void StiffSpringForceFieldCuda3d1_addExternalDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* dx1, const void* x1, const void* dx2, const void* x2, const void* dfdx, double factor);
#endif // SOFA_GPU_CUDA_DOUBLE
} // extern "C"
template<>
class CudaKernelsSpringForceField<CudaVec3fTypes>
{
public:
static void addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v)
{ SpringForceFieldCuda3f_addForce(nbVertex, nbSpringPerVertex, springs, f, x, v); }
static void addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2)
{ SpringForceFieldCuda3f_addExternalForce(nbVertex, nbSpringPerVertex, springs, f1, x1, v1, x2, v2); }
static void addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v, void* dfdx)
{ StiffSpringForceFieldCuda3f_addForce(nbVertex, nbSpringPerVertex, springs, f, x, v, dfdx); }
static void addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2, void* dfdx)
{ StiffSpringForceFieldCuda3f_addExternalForce(nbVertex, nbSpringPerVertex, springs, f1, x1, v1, x2, v2, dfdx); }
static void addDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* dx, const void* x, const void* dfdx, double factor)
{ StiffSpringForceFieldCuda3f_addDForce(nbVertex, nbSpringPerVertex, springs, f, dx, x, dfdx, factor); }
static void addExternalDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* dx1, const void* x1, const void* dx2, const void* x2, const void* dfdx, double factor)
{ StiffSpringForceFieldCuda3f_addExternalDForce(nbVertex, nbSpringPerVertex, springs, f1, dx1, x1, dx2, x2, dfdx, factor); }
};
template<>
class CudaKernelsSpringForceField<CudaVec3f1Types>
{
public:
static void addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v)
{ SpringForceFieldCuda3f1_addForce(nbVertex, nbSpringPerVertex, springs, f, x, v); }
static void addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2)
{ SpringForceFieldCuda3f1_addExternalForce(nbVertex, nbSpringPerVertex, springs, f1, x1, v1, x2, v2); }
static void addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v, void* dfdx)
{ StiffSpringForceFieldCuda3f1_addForce(nbVertex, nbSpringPerVertex, springs, f, x, v, dfdx); }
static void addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2, void* dfdx)
{ StiffSpringForceFieldCuda3f1_addExternalForce(nbVertex, nbSpringPerVertex, springs, f1, x1, v1, x2, v2, dfdx); }
static void addDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* dx, const void* x, const void* dfdx, double factor)
{ StiffSpringForceFieldCuda3f1_addDForce(nbVertex, nbSpringPerVertex, springs, f, dx, x, dfdx, factor); }
static void addExternalDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* dx1, const void* x1, const void* dx2, const void* x2, const void* dfdx, double factor)
{ StiffSpringForceFieldCuda3f1_addExternalDForce(nbVertex, nbSpringPerVertex, springs, f1, dx1, x1, dx2, x2, dfdx, factor); }
};
#ifdef SOFA_GPU_CUDA_DOUBLE
template<>
class CudaKernelsSpringForceField<CudaVec3dTypes>
{
public:
static void addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v)
{ SpringForceFieldCuda3d_addForce(nbVertex, nbSpringPerVertex, springs, f, x, v); }
static void addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2)
{ SpringForceFieldCuda3d_addExternalForce(nbVertex, nbSpringPerVertex, springs, f1, x1, v1, x2, v2); }
static void addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v, void* dfdx)
{ StiffSpringForceFieldCuda3d_addForce(nbVertex, nbSpringPerVertex, springs, f, x, v, dfdx); }
static void addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2, void* dfdx)
{ StiffSpringForceFieldCuda3d_addExternalForce(nbVertex, nbSpringPerVertex, springs, f1, x1, v1, x2, v2, dfdx); }
static void addDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* dx, const void* x, const void* dfdx, double factor)
{ StiffSpringForceFieldCuda3d_addDForce(nbVertex, nbSpringPerVertex, springs, f, dx, x, dfdx, factor); }
static void addExternalDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* dx1, const void* x1, const void* dx2, const void* x2, const void* dfdx, double factor)
{ StiffSpringForceFieldCuda3d_addExternalDForce(nbVertex, nbSpringPerVertex, springs, f1, dx1, x1, dx2, x2, dfdx, factor); }
};
template<>
class CudaKernelsSpringForceField<CudaVec3d1Types>
{
public:
static void addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v)
{ SpringForceFieldCuda3d1_addForce(nbVertex, nbSpringPerVertex, springs, f, x, v); }
static void addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2)
{ SpringForceFieldCuda3d1_addExternalForce(nbVertex, nbSpringPerVertex, springs, f1, x1, v1, x2, v2); }
static void addForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* x, const void* v, void* dfdx)
{ StiffSpringForceFieldCuda3d1_addForce(nbVertex, nbSpringPerVertex, springs, f, x, v, dfdx); }
static void addExternalForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* x1, const void* v1, const void* x2, const void* v2, void* dfdx)
{ StiffSpringForceFieldCuda3d1_addExternalForce(nbVertex, nbSpringPerVertex, springs, f1, x1, v1, x2, v2, dfdx); }
static void addDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f, const void* dx, const void* x, const void* dfdx, double factor)
{ StiffSpringForceFieldCuda3d1_addDForce(nbVertex, nbSpringPerVertex, springs, f, dx, x, dfdx, factor); }
static void addExternalDForce(unsigned int nbVertex, unsigned int nbSpringPerVertex, const void* springs, void* f1, const void* dx1, const void* x1, const void* dx2, const void* x2, const void* dfdx, double factor)
{ StiffSpringForceFieldCuda3d1_addExternalDForce(nbVertex, nbSpringPerVertex, springs, f1, dx1, x1, dx2, x2, dfdx, factor); }
};
#endif // SOFA_GPU_CUDA_DOUBLE
} // namespace cuda
} // namespace gpu
namespace component
{
namespace forcefield
{
using namespace gpu::cuda;
template<class TCoord, class TDeriv, class TReal>
void SpringForceFieldInternalData< gpu::cuda::CudaVectorTypes<TCoord,TDeriv,TReal> >::init(Main* m, bool stiff)
{
Data& data = m->data;
m->Inherit::init();
const sofa::helper::vector<Spring>& springs = m->springs.getValue();
if (!springs.empty())
{
bool external = (m->mstate1!=m->mstate2);
if (external)
{
std::map<int,int> nsprings1;
std::map<int,int> nsprings2;
for (unsigned int i=0;i<springs.size();i++)
{
nsprings1[springs[i].m1]++;
nsprings2[springs[i].m2]++;
}
int nmax1 = 0;
for (std::map<int,int>::const_iterator it = nsprings1.begin(); it != nsprings1.end(); ++it)
if (it->second > nmax1)
nmax1 = it->second;
data.springs1.init(nsprings1.begin()->first, nsprings1.rbegin()->first - nsprings1.begin()->first + 1, nmax1);
int nmax2 = 0;
for (std::map<int,int>::const_iterator it = nsprings2.begin(); it != nsprings2.end(); ++it)
if (it->second > nmax2)
nmax2 = it->second;
data.springs2.init(nsprings2.begin()->first, nsprings2.rbegin()->first - nsprings2.begin()->first + 1, nmax2);
nsprings1.clear();
nsprings2.clear();
for (unsigned int i=0;i<springs.size();i++)
{
int m1 = springs[i].m1 - data.springs1.vertex0;
int m2 = springs[i].m2 - data.springs2.vertex0;
data.springs1.set(m1, nsprings1[m1]++, m2,
(float)springs[i].initpos,
(float)springs[i].ks,
(float)springs[i].kd);
data.springs2.set(m2, nsprings2[m2]++, m1,
(float)springs[i].initpos,
(float)springs[i].ks,
(float)springs[i].kd);
}
}
else
{
std::map<int,int> nsprings;
for (unsigned int i=0;i<springs.size();i++)
{
nsprings[springs[i].m1]++;
nsprings[springs[i].m2]++;
}
int nmax = 0;
for (std::map<int,int>::const_iterator it = nsprings.begin(); it != nsprings.end(); ++it)
if (it->second > nmax)
nmax = it->second;
data.springs1.init(nsprings.begin()->first, nsprings.rbegin()->first - nsprings.begin()->first + 1, nmax);
std::cout << "CUDA SpringForceField: "<<springs.size()<<" springs, "<<data.springs1.nbVertex<<" attached points, max "<<data.springs1.nbSpringPerVertex<<" springs per point."<<std::endl;
nsprings.clear();
for (unsigned int i=0;i<springs.size();i++)
{
int m1 = springs[i].m1 - data.springs1.vertex0;
int m2 = springs[i].m2 - data.springs1.vertex0;
data.springs1.set(m1, nsprings[m1]++, m2,
(float)springs[i].initpos,
(float)springs[i].ks,
(float)springs[i].kd);
data.springs1.set(m2, nsprings[m2]++, m1,
(float)springs[i].initpos,
(float)springs[i].ks,
(float)springs[i].kd);
}
}
}
if (stiff)
{
data.springs1.dfdx.resize(data.springs1.springs.size());
data.springs2.dfdx.resize(data.springs2.springs.size());
}
}
// -- InteractionForceField interface
template<class TCoord, class TDeriv, class TReal>
void SpringForceFieldInternalData< gpu::cuda::CudaVectorTypes<TCoord,TDeriv,TReal> >::addForce(Main* m, bool stiff, VecDeriv& f1, VecDeriv& f2, const VecCoord& x1, const VecCoord& x2, const VecDeriv& v1, const VecDeriv& v2)
{
Data& data = m->data;
if (m->mstate1 == m->mstate2)
{
VecDeriv& f = f1;
const VecCoord& x = x1;
const VecDeriv& v = v1;
f.resize(x.size());
int d = data.springs1.vertex0;
if (data.springs1.nbSpringPerVertex > 0)
{
if (!stiff)
Kernels::addForce(data.springs1.nbVertex,
data.springs1.nbSpringPerVertex,
data.springs1.springs.deviceRead(),
( Deriv*)f.deviceWrite() + d,
(const Coord*)x.deviceRead() + d,
(const Deriv*)v.deviceRead() + d);
else
Kernels::addForce(data.springs1.nbVertex,
data.springs1.nbSpringPerVertex,
data.springs1.springs.deviceRead(),
( Deriv*)f.deviceWrite() + d,
(const Coord*)x.deviceRead() + d,
(const Deriv*)v.deviceRead() + d,
data.springs1.dfdx.deviceWrite());
}
}
else
{
f1.resize(x1.size());
f2.resize(x2.size());
int d1 = data.springs1.vertex0;
int d2 = data.springs2.vertex0;
if (data.springs1.nbSpringPerVertex > 0)
{
if (!stiff)
Kernels::addExternalForce(data.springs1.nbVertex,
data.springs1.nbSpringPerVertex,
data.springs1.springs.deviceRead(),
( Deriv*)f1.deviceWrite() + d1,
(const Coord*)x1.deviceRead() + d1,
(const Deriv*)v1.deviceRead() + d1,
(const Coord*)x2.deviceRead() + d2,
(const Deriv*)v2.deviceRead() + d2);
else
Kernels::addExternalForce(data.springs1.nbVertex,
data.springs1.nbSpringPerVertex,
data.springs1.springs.deviceRead(),
( Deriv*)f1.deviceWrite() + d1,
(const Coord*)x1.deviceRead() + d1,
(const Deriv*)v1.deviceRead() + d1,
(const Coord*)x2.deviceRead() + d2,
(const Deriv*)v2.deviceRead() + d2,
data.springs1.dfdx.deviceWrite());
}
if (data.springs2.nbSpringPerVertex > 0)
{
if (!stiff)
Kernels::addExternalForce(data.springs2.nbVertex,
data.springs2.nbSpringPerVertex,
data.springs2.springs.deviceRead(),
( Deriv*)f2.deviceWrite() + d2,
(const Coord*)x2.deviceRead() + d2,
(const Deriv*)v2.deviceRead() + d2,
(const Coord*)x1.deviceRead() + d1,
(const Deriv*)v1.deviceRead() + d1);
else
Kernels::addExternalForce(data.springs2.nbVertex,
data.springs2.nbSpringPerVertex,
data.springs2.springs.deviceRead(),
( Deriv*)f2.deviceWrite() + d2,
(const Coord*)x2.deviceRead() + d2,
(const Deriv*)v2.deviceRead() + d2,
(const Coord*)x1.deviceRead() + d1,
(const Deriv*)v1.deviceRead() + d1,
data.springs2.dfdx.deviceWrite());
}
}
}
template<class TCoord, class TDeriv, class TReal>
void SpringForceFieldInternalData< gpu::cuda::CudaVectorTypes<TCoord,TDeriv,TReal> >::addDForce(Main* m, bool stiff, VecDeriv& df1, VecDeriv& df2, const VecDeriv& dx1, const VecDeriv& dx2, double kFactor, double /*bFactor*/)
{
if (!stiff) return;
Data& data = m->data;
if (m->mstate1 == m->mstate2)
{
VecDeriv& df = df1;
const VecDeriv& dx = dx1;
const VecCoord& x = *m->mstate1->getX();
df.resize(x.size());
int d = data.springs1.vertex0;
if (data.springs1.nbSpringPerVertex > 0)
{
Kernels::addDForce(data.springs1.nbVertex,
data.springs1.nbSpringPerVertex,
data.springs1.springs.deviceRead(),
( Deriv*)df.deviceWrite() + d,
(const Deriv*)dx.deviceRead() + d,
(const Coord*)x.deviceRead() + d,
data.springs1.dfdx.deviceRead(),
kFactor);
}
}
else
{
const VecCoord& x1 = *m->mstate1->getX();
const VecCoord& x2 = *m->mstate2->getX();
df1.resize(x1.size());
df2.resize(x2.size());
int d1 = data.springs1.vertex0;
int d2 = data.springs2.vertex0;
if (data.springs1.nbSpringPerVertex > 0)
{
Kernels::addExternalDForce(data.springs1.nbVertex,
data.springs1.nbSpringPerVertex,
data.springs1.springs.deviceRead(),
( Deriv*)df1.deviceWrite() + d1,
(const Coord*)x1.deviceRead() + d1,
(const Deriv*)dx1.deviceRead() + d1,
(const Coord*)x2.deviceRead() + d2,
(const Deriv*)dx2.deviceRead() + d2,
data.springs1.dfdx.deviceRead(),
kFactor);
}
if (data.springs2.nbSpringPerVertex > 0)
{
Kernels::addExternalDForce(data.springs2.nbVertex,
data.springs2.nbSpringPerVertex,
data.springs2.springs.deviceRead(),
( Deriv*)df2.deviceWrite() + d2,
(const Deriv*)dx2.deviceRead() + d2,
(const Coord*)x2.deviceRead() + d2,
(const Deriv*)dx1.deviceRead() + d1,
(const Coord*)x1.deviceRead() + d1,
data.springs2.dfdx.deviceRead(),
kFactor);
}
}
}
// I know using macros is bad design but this is the only way not to repeat the code for all CUDA types
#define CudaSpringForceField_ImplMethods(T) \
template<> void SpringForceField< T >::init() \
{ data.init(this, false); } \
template<> void SpringForceField< T >::addForce(VecDeriv& f1, VecDeriv& f2, const VecCoord& x1, const VecCoord& x2, const VecDeriv& v1, const VecDeriv& v2) \
{ data.addForce(this, false, f1, f2, x1, x2, v1, v2); } \
template<> void StiffSpringForceField< T >::init() \
{ data.init(this, true); } \
template<> void StiffSpringForceField< T >::addForce(VecDeriv& f1, VecDeriv& f2, const VecCoord& x1, const VecCoord& x2, const VecDeriv& v1, const VecDeriv& v2) \
{ data.addForce(this, true, f1, f2, x1, x2, v1, v2); } \
template<> void StiffSpringForceField< T >::addDForce(VecDeriv& df1, VecDeriv& df2, const VecDeriv& dx1, const VecDeriv& dx2, double kFactor, double bFactor) \
{ data.addDForce(this, true, df1, df2, dx1, dx2, kFactor, bFactor); }
CudaSpringForceField_ImplMethods(gpu::cuda::CudaVec3fTypes);
CudaSpringForceField_ImplMethods(gpu::cuda::CudaVec3f1Types);
#ifdef SOFA_GPU_CUDA_DOUBLE
CudaSpringForceField_ImplMethods(gpu::cuda::CudaVec3dTypes);
CudaSpringForceField_ImplMethods(gpu::cuda::CudaVec3d1Types);
#endif // SOFA_GPU_CUDA_DOUBLE
#undef CudaSpringForceField_ImplMethods
} // namespace forcefield
} // namespace component
} // namespace sofa
#endif
|