Skip to content

Commit

Permalink
Undid last commit and implemented better way of storing domain materi…
Browse files Browse the repository at this point in the history
…al axes.
  • Loading branch information
SteveMaas1978 committed Apr 30, 2024
1 parent 16622f1 commit 2c081ff
Show file tree
Hide file tree
Showing 10 changed files with 15 additions and 49 deletions.
6 changes: 0 additions & 6 deletions FEBioLib/input.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -96,12 +96,6 @@ void FEBioModel::print_parameter(FEParam& p, int level)
feLog("%s : %lg,%lg,%lg,%lg,%lg,%lg,%lg,%lg,%lg\n", sz, m(0,0), m(0,1), m(0,2), m(1,0), m(1,1), m(1,2), m(2,0), m(2,1), m(2,2));
}
break;
case FE_PARAM_QUAT:
{
quatd v = p.value<quatd>();
feLog("%s : %lg,%lg,%lg,%lg\n", sz, v.x, v.y, v.z, v.w);
}
break;
case FE_PARAM_TENS3DRS:
{
tens3drs m = p.value<tens3drs>();
Expand Down
35 changes: 0 additions & 35 deletions FEBioXML/FileImport.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -671,41 +671,6 @@ bool FEFileSection::ReadParameter(XMLTag& tag, FEParameterList& pl, const char*
case FE_PARAM_TENS3DRS: value(tag, pp->value<tens3drs>()); break;
case FE_PARAM_STRING: value(tag, pp->cvalue()); break;
case FE_PARAM_STD_STRING: value(tag, pp->value<string>()); break;
case FE_PARAM_QUAT:
{
const char* sztype = tag.AttributeValue("type", true);
if (sztype == nullptr) sztype = "vector";
if (strcmp(sztype, "vector") == 0)
{
vec3d a(1,0,0), d(0,1,0);
if (tag.isleaf()) throw XMLReader::InvalidValue(tag);

++tag;
do {
if (tag == "a") value(tag, a);
else if (tag == "d") value(tag, d);
++tag;
} while (!tag.isend());

mat3d Q(a, d);
pp->value<quatd>() = quatd(Q);
}
else if (strcmp(sztype, "mat3") == 0)
{
mat3d Q;
value(tag, Q);
pp->value<quatd>() = quatd(Q);
}
else if (strcmp(sztype, "quat") == 0)
{
double v[4] = { 0,0,0,1 };
int nr = tag.value(v, 4);
if (nr != 4) throw XMLReader::InvalidValue(tag);
pp->value<quatd>() = quatd(v[0], v[1], v[2], v[3]);
}
else throw XMLReader::InvalidAttributeValue(tag, "type");
}
break;
case FE_PARAM_DATA_ARRAY:
{
// get the surface map
Expand Down
8 changes: 6 additions & 2 deletions FECore/FEDomain.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -33,7 +33,7 @@ SOFTWARE.*/
//-----------------------------------------------------------------------------
FEDomain::FEDomain(int nclass, FEModel* fem) : FEMeshPartition(nclass, fem)
{

m_matAxis = nullptr;
}

//-----------------------------------------------------------------------------
Expand All @@ -55,6 +55,10 @@ void FEDomain::SetMatID(int mid)
// determines how many integration points an element gets).
void FEDomain::CreateMaterialPointData()
{
// This function is called before Init is called, so we need to do the
// initialization of the mat_axis here.
if (m_matAxis) m_matAxis->Init();

FEMaterial* pmat = GetMaterial();
FEMesh* mesh = GetMesh();
if (pmat) ForEachElement([=](FEElement& el) {
Expand All @@ -66,7 +70,7 @@ void FEDomain::CreateMaterialPointData()
for (int k = 0; k < el.GaussPoints(); ++k)
{
FEMaterialPoint* mp = new FEMaterialPoint(pmat->CreateMaterialPointData());
mp->m_Q = m_matAxis;
mp->m_Q = (m_matAxis ? m_matAxis->operator()(*mp) : mat3d::identity());
mp->m_r0 = el.Evaluate(r, k);
mp->m_index = k;
el.SetMaterialPointData(mp, k);
Expand Down
3 changes: 2 additions & 1 deletion FECore/FEDomain.h
Original file line number Diff line number Diff line change
Expand Up @@ -28,6 +28,7 @@ SOFTWARE.*/

#pragma once
#include "FEMeshPartition.h"
#include "FEMat3dValuator.h"

// forward declaration of material class
class FEMaterial;
Expand Down Expand Up @@ -90,5 +91,5 @@ class FECORE_API FEDomain : public FEMeshPartition
void UnpackLM(FEElement& el, const FEDofList& dof, vector<int>& lm);

protected:
quatd m_matAxis; // initial material axis
FEMat3dValuator* m_matAxis; // initial material axis
};
1 change: 0 additions & 1 deletion FECore/FEParam.h
Original file line number Diff line number Diff line change
Expand Up @@ -64,7 +64,6 @@ enum FEParamType {
FE_PARAM_MAT3D_MAPPED,
FE_PARAM_MAT3DS_MAPPED,
FE_PARAM_MATERIALPOINT,
FE_PARAM_QUAT
};

//-----------------------------------------------------------------------------
Expand Down
1 change: 0 additions & 1 deletion FECore/FEParameterList.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -369,7 +369,6 @@ FEParam* FEParamContainer::AddParameter(vec2d& v, const char*
FEParam* FEParamContainer::AddParameter(vec3d& v, const char* sz) { return AddParameter(&v, FE_PARAM_VEC3D, 1, sz); }
FEParam* FEParamContainer::AddParameter(mat3d& v, const char* sz) { return AddParameter(&v, FE_PARAM_MAT3D, 1, sz); }
FEParam* FEParamContainer::AddParameter(mat3ds& v, const char* sz) { return AddParameter(&v, FE_PARAM_MAT3DS, 1, sz); }
FEParam* FEParamContainer::AddParameter(quatd& v, const char* sz) { return AddParameter(&v, FE_PARAM_QUAT, 1, sz); }
FEParam* FEParamContainer::AddParameter(FEParamDouble& v, const char* sz) { return AddParameter(&v, FE_PARAM_DOUBLE_MAPPED, 1, sz); }
FEParam* FEParamContainer::AddParameter(FEParamVec3& v, const char* sz) { return AddParameter(&v, FE_PARAM_VEC3D_MAPPED, 1, sz); }
FEParam* FEParamContainer::AddParameter(FEParamMat3d& v, const char* sz) { return AddParameter(&v, FE_PARAM_MAT3D_MAPPED, 1, sz); }
Expand Down
2 changes: 0 additions & 2 deletions FECore/FEParameterList.h
Original file line number Diff line number Diff line change
Expand Up @@ -30,7 +30,6 @@ SOFTWARE.*/
#include "vec2d.h"
#include "vec3d.h"
#include "mat3d.h"
#include "quatd.h"
#include <assert.h>
#include <list>
#include <memory>
Expand Down Expand Up @@ -167,7 +166,6 @@ class FECORE_API FEParamContainer
FEParam* AddParameter(vec3d& v, const char* sz);
FEParam* AddParameter(mat3d& v, const char* sz);
FEParam* AddParameter(mat3ds& v, const char* sz);
FEParam* AddParameter(quatd& v, const char* sz);
FEParam* AddParameter(FEParamDouble& v, const char* sz);
FEParam* AddParameter(FEParamVec3& v, const char* sz);
FEParam* AddParameter(FEParamMat3d& v, const char* sz);
Expand Down
4 changes: 4 additions & 0 deletions FECore/FEShellDomain.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -31,6 +31,10 @@ SOFTWARE.*/
#include "FEMesh.h"
#include "FEMaterial.h"

BEGIN_FECORE_CLASS(FEShellDomain, FEDomain)
ADD_PROPERTY(m_matAxis, "mat_axis");
END_FECORE_CLASS();

//-----------------------------------------------------------------------------
//! constructor
FEShellDomain::FEShellDomain(FEModel* fem) : FEDomain(FE_DOMAIN_SHELL, fem)
Expand Down
2 changes: 2 additions & 0 deletions FECore/FEShellDomain.h
Original file line number Diff line number Diff line change
Expand Up @@ -75,6 +75,8 @@ class FECORE_API FEShellDomain : public FEDomain

public:
void ForEachShellElement(std::function<void(FEShellElement& el)> f);

DECLARE_FECORE_CLASS();
};

//-----------------------------------------------------------------------------
Expand Down
2 changes: 1 addition & 1 deletion FECore/FESolidDomain.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -34,7 +34,7 @@ SOFTWARE.*/
#include "FEModel.h"

BEGIN_FECORE_CLASS(FESolidDomain, FEDomain)
ADD_PARAMETER(m_matAxis, "mat_axis");
ADD_PROPERTY(m_matAxis, "mat_axis");
END_FECORE_CLASS();

//-----------------------------------------------------------------------------
Expand Down

0 comments on commit 2c081ff

Please sign in to comment.