Skip to content

Commit

Permalink
[MechanicalLoad] Update TrianglePressureForceField to remove plan sel…
Browse files Browse the repository at this point in the history
…ection and add callback on Data pressure and triangleList (sofa-framework#4743)

* [MechanicalLoad] Add callback on TrianglePressureForceField pressure Data

* [MechanicalLoad] Update TrianglePressureFF to remove plan selection and add callback on main Data

* Fix tests

* Update simpleBoundaryConditions.scn
  • Loading branch information
epernod authored Jun 11, 2024
1 parent 58d0ab0 commit ca70b0f
Show file tree
Hide file tree
Showing 5 changed files with 55 additions and 105 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -54,17 +54,16 @@ class TrianglePressureForceField : public core::behavior::ForceField<DataTypes>

Data<sofa::type::vector<Index> > triangleList; ///< Indices of triangles separated with commas where a pressure is applied

/// the normal used to define the edge subjected to the pressure force.
Data<Deriv> normal;

Data<Real> dmin; ///< Minimum distance from the origin along the normal direction
Data<Real> dmax; ///< Maximum distance from the origin along the normal direction
Data<bool> p_showForces; ///< draw triangles which have a given pressure
Data<bool> p_useConstantForce; ///< applied force is computed as the pressure vector times the area at rest

/// Link to be set to the topology container in the component graph.
SingleLink<TrianglePressureForceField<DataTypes>, sofa::core::topology::BaseMeshTopology, BaseLink::FLAG_STOREPATH | BaseLink::FLAG_STRONGLINK> l_topology;


core::objectmodel::lifecycle::DeprecatedData normal{ this, "v24.06", "v24.12", "normal", "Plan selection using normal, dmin, dmax has been removed. Triangles should be selected using an Engine.Select and passed using Data triangleList" };
core::objectmodel::lifecycle::DeprecatedData dmin{ this, "v24.06", "v24.12", "dmin", "Plan selection using normal, dmin, dmax has been removed. Triangles should be selected using an Engine.Select and passed using Data triangleList" };
core::objectmodel::lifecycle::DeprecatedData dmax{ this, "v24.06", "v24.12", "dmax", "Plan selection using normal, dmin, dmax has been removed. Triangles should be selected using an Engine.Select and passed using Data triangleList" };

protected:

class TrianglePressureInformation
Expand Down Expand Up @@ -125,16 +124,11 @@ class TrianglePressureForceField : public core::behavior::ForceField<DataTypes>
SReal getPotentialEnergy(const core::MechanicalParams* /*mparams*/, const DataVecCoord& /* x */) const override;
void draw(const core::visual::VisualParams* vparams) override;

void setDminAndDmax(const SReal _dmin, const SReal _dmax){dmin.setValue((Real)_dmin); dmax.setValue((Real)_dmax);}
void setNormal(const Coord n) { normal.setValue(n);}
void setPressure(Deriv _pressure) { this->pressure = _pressure; updateTriangleInformation(); }

protected :
void selectTrianglesAlongPlane();
void selectTrianglesFromString();
void updateTriangleInformation();
void initTriangleInformation();
bool isPointInPlane(Coord p);
};

#if !defined(SOFA_COMPONENT_FORCEFIELD_TRIANGLEPRESSUREFORCEFIELD_CPP)
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -40,19 +40,32 @@ template <class DataTypes> TrianglePressureForceField<DataTypes>::TrianglePress
: pressure(initData(&pressure, "pressure", "Pressure force per unit area"))
, cauchyStress(initData(&cauchyStress, MatSym3(),"cauchyStress", "Cauchy Stress applied on the normal of each triangle"))
, triangleList(initData(&triangleList,"triangleList", "Indices of triangles separated with commas where a pressure is applied"))
, normal(initData(&normal,"normal", "Normal direction for the plane selection of triangles"))
, dmin(initData(&dmin,(Real)0.0, "dmin", "Minimum distance from the origin along the normal direction"))
, dmax(initData(&dmax,(Real)0.0, "dmax", "Maximum distance from the origin along the normal direction"))
, p_showForces(initData(&p_showForces, (bool)false, "showForces", "draw triangles which have a given pressure"))
, p_useConstantForce(initData(&p_useConstantForce, (bool)true, "useConstantForce", "applied force is computed as the pressure vector times the area at rest"))
, l_topology(initLink("topology", "link to the topology container"))
, trianglePressureMap(initData(&trianglePressureMap, "trianglePressureMap", "Map between triangle indices and their pressure"))
, m_topology(nullptr)
{}
{
this->addUpdateCallback("pressure_change", { &pressure }, [this](const core::DataTracker& t)
{
SOFA_UNUSED(t);
updateTriangleInformation();
return sofa::core::objectmodel::ComponentState::Valid;
}, {});


this->addUpdateCallback("triangles_change", { &triangleList }, [this](const core::DataTracker& t)
{
SOFA_UNUSED(t);
initTriangleInformation();
return sofa::core::objectmodel::ComponentState::Valid;
}, {});

}


template <class DataTypes> void TrianglePressureForceField<DataTypes>::init()
{

this->core::behavior::ForceField<DataTypes>::init();

if (l_topology.empty())
Expand All @@ -71,15 +84,6 @@ template <class DataTypes> void TrianglePressureForceField<DataTypes>::init()
return;
}

if (dmin.getValue()!=dmax.getValue())
{
selectTrianglesAlongPlane();
}
if (triangleList.getValue().size()>0)
{
selectTrianglesFromString();
}

trianglePressureMap.createTopologyHandler(m_topology);

initTriangleInformation();
Expand Down Expand Up @@ -138,101 +142,52 @@ void TrianglePressureForceField<DataTypes>::addDForce(const core::MechanicalPara

template<class DataTypes>
void TrianglePressureForceField<DataTypes>::initTriangleInformation()
{
const sofa::type::vector<Index>& my_map = trianglePressureMap.getMap2Elements();
{
if (triangleList.getValue().empty())
return;

// Get list of input triangle indices
type::vector<Index> _triangleList = triangleList.getValue();

// Get write access to TopologySubset Data storing pressure information per triangle
auto my_subset = sofa::helper::getWriteOnlyAccessor(trianglePressureMap);

// Set the list of triangles indices as map of this TopologySubset Data
trianglePressureMap.setMap2Elements(_triangleList);

// Fill pressure data
const VecCoord& x0 = this->mstate->read(core::ConstVecCoordId::restPosition())->getValue();
const Deriv& my_pressure = pressure.getValue();

for (unsigned int i=0; i<my_map.size(); ++i)
for (unsigned int i = 0; i < _triangleList.size(); ++i)
{
const auto& t = this->m_topology->getTriangle(my_map[i]);
const auto& t = this->m_topology->getTriangle(_triangleList[i]);

const auto& n0 = DataTypes::getCPos(x0[t[0]]);
const auto& n1 = DataTypes::getCPos(x0[t[1]]);
const auto& n2 = DataTypes::getCPos(x0[t[2]]);

my_subset[i].area = sofa::geometry::Triangle::area(n0, n1, n2);
my_subset[i].force=pressure.getValue()*my_subset[i].area;
TrianglePressureInformation tInfo;
tInfo.area = sofa::geometry::Triangle::area(n0, n1, n2);
tInfo.force = my_pressure * tInfo.area;
my_subset.push_back(tInfo);
}
}

template<class DataTypes>
bool TrianglePressureForceField<DataTypes>::isPointInPlane(Coord p)
{
Real d=dot(p,normal.getValue());
if ((d>dmin.getValue())&& (d<dmax.getValue()))
return true;
else
return false;
}

template<class DataTypes>
void TrianglePressureForceField<DataTypes>::updateTriangleInformation()
{
sofa::type::vector<TrianglePressureInformation>& my_subset = *(trianglePressureMap).beginEdit();

const Deriv& my_pressure = pressure.getValue();
for (unsigned int i = 0; i < my_subset.size(); ++i)
my_subset[i].force = (my_pressure * my_subset[i].area);

for (unsigned int i=0; i<my_subset.size(); ++i)
my_subset[i].force=(pressure.getValue()*my_subset[i].area);

trianglePressureMap.endEdit();
}


template <class DataTypes>
void TrianglePressureForceField<DataTypes>::selectTrianglesAlongPlane()
{
const VecCoord& x = this->mstate->read(core::ConstVecCoordId::restPosition())->getValue();
std::vector<bool> vArray;

vArray.resize(x.size());

for( unsigned int i=0; i<x.size(); ++i)
{
vArray[i]=isPointInPlane(x[i]);
}

sofa::type::vector<TrianglePressureInformation>& my_subset = *(trianglePressureMap).beginEdit();
type::vector<Index> inputTriangles;

for (size_t n=0; n<m_topology->getNbTriangles(); ++n)
{
if ((vArray[m_topology->getTriangle(n)[0]]) && (vArray[m_topology->getTriangle(n)[1]])&& (vArray[m_topology->getTriangle(n)[2]]) )
{
// insert a dummy element : computation of pressure done later
TrianglePressureInformation t;
t.area = 0;
my_subset.push_back(t);
inputTriangles.push_back(n);
}
}
trianglePressureMap.endEdit();
trianglePressureMap.setMap2Elements(inputTriangles);

return;
}


template <class DataTypes>
void TrianglePressureForceField<DataTypes>::selectTrianglesFromString()
{
sofa::type::vector<TrianglePressureInformation>& my_subset = *(trianglePressureMap).beginEdit();
type::vector<Index> _triangleList = triangleList.getValue();

trianglePressureMap.setMap2Elements(_triangleList);

for (unsigned int i = 0; i < _triangleList.size(); ++i)
{
TrianglePressureInformation t;
t.area = 0;
my_subset.push_back(t);
}

trianglePressureMap.endEdit();

return;
}


template<class DataTypes>
void TrianglePressureForceField<DataTypes>::draw(const core::visual::VisualParams* vparams)
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -65,16 +65,15 @@ struct TrianglePressureForceField_test : public ForceField_test<_TrianglePressur
DataTypes::set( v[2], 0,0,0);

//Force
f.resize(3);
f.resize(3);
Vec3 f0(0,0,0.1);
DataTypes::set( f[0], f0[0], f0[1], f0[2]);
DataTypes::set( f[1], f0[0], f0[1], f0[2]);
DataTypes::set( f[2], f0[0], f0[1], f0[2]);

// Set the properties of the force field
Inherited::force->normal.setValue(Deriv(0,0,1));
Inherited::force->dmin.setValue(-0.01);
Inherited::force->dmax.setValue(0.01);
sofa::type::vector<Index> indices = {0};
Inherited::force->triangleList.setValue(indices);
Inherited::force->pressure=Coord(0,0,0.6);
}

Expand Down
Original file line number Diff line number Diff line change
@@ -1,3 +1,4 @@
<?xml version="1.0"?>
<!-- TrianglePressureForceField example scene -->
<Node name="root" dt="0.05" showBoundingTree="0" gravity="0 0 0">
<RequiredPlugin name="Sofa.Component.Collision.Detection.Algorithm"/> <!-- Needed to use components [BVHNarrowPhase BruteForceBroadPhase CollisionPipeline] -->
Expand Down Expand Up @@ -32,13 +33,13 @@
<MechanicalObject src="@loader" name="Volume" />
<include href="Objects/TetrahedronSetTopology.xml" src="@loader" />
<DiagonalMass massDensity="0.5" />
<FixedPlaneProjectiveConstraint direction="0 0 1" dmin="-0.1" dmax="0.1" />
<FixedProjectiveConstraint indices="0" />
<FixedProjectiveConstraint indices="0-50" />
<TetrahedronFEMForceField name="FEM" youngModulus="60" poissonRatio="0.3" computeGlobalMatrix="false" method="large" />
<Node name="T">
<include href="Objects/TriangleSetTopology.xml" src="@../Container"/>
<Tetra2TriangleTopologicalMapping input="@../Container" output="@Container" />
<TrianglePressureForceField name="PFF" normal="0 0 1" dmin="0.9" dmax="1.1" pressure="1 0 0" />
<BoxROI name="Roi" position="@../Volume.position" triangles="@Container.triangles" box="-0.5 -0.5 0.9 0.5 0.5 1.1"/>
<TrianglePressureForceField name="PFF" triangleList="@Roi.triangleIndices" pressure="0 10 10" showForces="1"/>
<TriangleCollisionModel />
<Node name="Visu">
<OglModel name="Visual" color="yellow" />
Expand Down
3 changes: 2 additions & 1 deletion examples/Demos/simpleBoundaryConditions.scn
Original file line number Diff line number Diff line change
Expand Up @@ -54,7 +54,8 @@
<Tetra2TriangleTopologicalMapping input="@TetraContainer" output="@TriangleContainer" />

<!--Apply a pressure to the triangles in the x=1 plane -->
<TrianglePressureForceField normal="1 0 0" dmin="0.99" dmax="1.01" pressure="100 0 0" showForces="1"/>
<BoxROI name="pressureROI" position="@../../DOF.position" triangles="@TriangleContainer.triangles" box="0.99 -0.01 -0.01 1.01 1.01 1.01" drawBoxes="1"/>
<TrianglePressureForceField name="PFF" triangleList="@pressureROI.triangleIndices" pressure="100 0 0" showForces="1"/>
</Node>
</Node>
</Node>
Expand Down

0 comments on commit ca70b0f

Please sign in to comment.