#ifndef __CONTACTOPS_HPP__
#define __CONTACTOPS_HPP__
struct CommonData :
public boost::enable_shared_from_this<CommonData> {
static SmartPetscObj<Vec>
constexpr int ghosts[] = {0, 1, 2, 3, 4};
}
const double *t_ptr;
"get array");
t_ptr[4]};
"restore array");
} else {
}
}
return boost::shared_ptr<MatrixDouble>(shared_from_this(),
}
return boost::shared_ptr<MatrixDouble>(shared_from_this(), &
contactDisp);
}
return boost::shared_ptr<MatrixDouble>(shared_from_this(),
}
return boost::shared_ptr<VectorDouble>(shared_from_this(), &
sdfVals);
}
return boost::shared_ptr<MatrixDouble>(shared_from_this(), &
gradsSdf);
}
return boost::shared_ptr<MatrixDouble>(shared_from_this(), &
hessSdf);
}
return boost::shared_ptr<VectorDouble>(shared_from_this(), &
constraintVals);
}
};
#ifdef PYTHON_SDF
struct SDFPython {
SDFPython() = default;
virtual ~SDFPython() = default;
try {
auto main_module = bp::import("__main__");
mainNamespace = main_module.attr("__dict__");
bp::exec_file(py_file.c_str(), mainNamespace, mainNamespace);
sdfFun = mainNamespace["sdf"];
sdfGradFun = mainNamespace["grad_sdf"];
sdfHessFun = mainNamespace["hess_sdf"];
} catch (bp::error_already_set const &) {
PyErr_Print();
}
};
template <typename T>
inline std::vector<T>
py_list_to_std_vector(const boost::python::object &iterable) {
return std::vector<T>(boost::python::stl_input_iterator<T>(iterable),
boost::python::stl_input_iterator<T>());
}
double delta_t,
double t, np::ndarray x, np::ndarray y, np::ndarray z,
np::ndarray tx, np::ndarray ty, np::ndarray tz, int block_id,
) {
try {
sdf = bp::extract<np::ndarray>(
sdfFun(delta_t,
t, x, y, z, tx, ty, tz, block_id));
} catch (bp::error_already_set const &) {
PyErr_Print();
}
}
double delta_t,
double t, np::ndarray x, np::ndarray y, np::ndarray z,
np::ndarray tx, np::ndarray ty, np::ndarray tz, int block_id,
) {
try {
sdfGradFun(delta_t,
t, x, y, z, tx, ty, tz, block_id));
} catch (bp::error_already_set const &) {
PyErr_Print();
}
}
double delta_t,
double t, np::ndarray x, np::ndarray y, np::ndarray z,
np::ndarray tx, np::ndarray ty, np::ndarray tz, int block_id,
) {
try {
sdfHessFun(delta_t,
t, x, y, z, tx, ty, tz, block_id));
} catch (bp::error_already_set const &) {
PyErr_Print();
}
}
private:
bp::object mainNamespace;
bp::object sdfFun;
bp::object sdfGradFun;
bp::object sdfHessFun;
};
static boost::weak_ptr<SDFPython> sdfPythonWeakPtr;
inline np::ndarray convert_to_numpy(
VectorDouble &data,
int nb_gauss_pts,
int id) {
auto dtype = np::dtype::get_builtin<double>();
auto size = bp::make_tuple(nb_gauss_pts);
auto stride = bp::make_tuple(3 * sizeof(double));
return (np::from_data(&data[id], dtype, size, stride, bp::object()));
};
#endif
double delta_t,
double t,
int nb_gauss_pts,
MatrixDouble &spatial_coords,
double delta_t,
double t,
int nb_gauss_pts,
MatrixDouble &spatial_coords,
double delta_t,
double t,
int nb_gauss_pts,
MatrixDouble &spatial_coords,
int nb_gauss_pts,
int block_id) {
#ifdef PYTHON_SDF
if (auto sdf_ptr = sdfPythonWeakPtr.lock()) {
bp::list python_coords;
bp::list python_normals;
for (int idx = 0; idx < 3; ++idx) {
python_coords.append(
convert_to_numpy(v_spatial_coords, nb_gauss_pts, idx));
python_normals.append(
convert_to_numpy(v_normal_at_pts, nb_gauss_pts, idx));
}
np::ndarray np_sdf = np::empty(bp::make_tuple(nb_gauss_pts),
np::dtype::get_builtin<double>());
bp::extract<np::ndarray>(python_coords[0]),
bp::extract<np::ndarray>(python_coords[1]),
bp::extract<np::ndarray>(python_coords[2]),
bp::extract<np::ndarray>(python_normals[0]),
bp::extract<np::ndarray>(python_normals[1]),
bp::extract<np::ndarray>(python_normals[2]),
block_id, np_sdf),
"Failed python call");
double *sdf_val_ptr = reinterpret_cast<double *>(np_sdf.get_data());
v_sdf.resize(nb_gauss_pts, false);
for (size_t gg = 0; gg < nb_gauss_pts; ++gg)
v_sdf[gg] = *(sdf_val_ptr + gg);
return v_sdf;
}
#endif
v_sdf.resize(nb_gauss_pts, false);
auto t_coords = getFTensor1FromPtr<3>(&m_spatial_coords(0, 0));
for (size_t gg = 0; gg < nb_gauss_pts; ++gg) {
v_sdf[gg] = -t_coords(2) - 0.1;
++t_coords;
}
return v_sdf;
}
#ifdef PYTHON_SDF
if (auto sdf_ptr = sdfPythonWeakPtr.lock()) {
bp::list python_coords;
bp::list python_normals;
for (int idx = 0; idx < 3; ++idx) {
python_coords.append(
convert_to_numpy(v_spatial_coords, nb_gauss_pts, idx));
python_normals.append(
convert_to_numpy(v_normal_at_pts, nb_gauss_pts, idx));
}
np::ndarray np_grad_sdf = np::empty(bp::make_tuple(nb_gauss_pts, 3),
np::dtype::get_builtin<double>());
delta_t,
t, bp::extract<np::ndarray>(python_coords[0]),
bp::extract<np::ndarray>(python_coords[1]),
bp::extract<np::ndarray>(python_coords[2]),
bp::extract<np::ndarray>(python_normals[0]),
bp::extract<np::ndarray>(python_normals[1]),
bp::extract<np::ndarray>(python_normals[2]), block_id,
np_grad_sdf),
"Failed python call");
double *grad_ptr = reinterpret_cast<double *>(np_grad_sdf.get_data());
m_grad_sdf.resize(3, nb_gauss_pts, false);
for (size_t gg = 0; gg < nb_gauss_pts; ++gg) {
for (int idx = 0; idx < 3; ++idx)
m_grad_sdf(idx, gg) = *(grad_ptr + (3 * gg + idx));
}
return m_grad_sdf;
}
#endif
m_grad_sdf.resize(3, nb_gauss_pts, false);
auto t_grad_sdf = getFTensor1FromMat<3>(m_grad_sdf);
for (size_t gg = 0; gg < nb_gauss_pts; ++gg) {
t_grad_sdf(
i) = t_grad_sdf_set(
i);
++t_grad_sdf;
}
return m_grad_sdf;
}
#ifdef PYTHON_SDF
if (auto sdf_ptr = sdfPythonWeakPtr.lock()) {
bp::list python_coords;
bp::list python_normals;
for (int idx = 0; idx < 3; ++idx) {
python_coords.append(
convert_to_numpy(v_spatial_coords, nb_gauss_pts, idx));
python_normals.append(
convert_to_numpy(v_normal_at_pts, nb_gauss_pts, idx));
};
np::ndarray np_hess_sdf = np::empty(bp::make_tuple(nb_gauss_pts, 6),
np::dtype::get_builtin<double>());
delta_t,
t, bp::extract<np::ndarray>(python_coords[0]),
bp::extract<np::ndarray>(python_coords[1]),
bp::extract<np::ndarray>(python_coords[2]),
bp::extract<np::ndarray>(python_normals[0]),
bp::extract<np::ndarray>(python_normals[1]),
bp::extract<np::ndarray>(python_normals[2]), block_id,
np_hess_sdf),
"Failed python call");
double *hess_ptr = reinterpret_cast<double *>(np_hess_sdf.get_data());
m_hess_sdf.resize(6, nb_gauss_pts, false);
for (size_t gg = 0; gg < nb_gauss_pts; ++gg) {
for (int idx = 0; idx < 6; ++idx)
m_hess_sdf(idx, gg) =
*(hess_ptr + (6 * gg + idx));
}
return m_hess_sdf;
}
#endif
m_hess_sdf.resize(6, nb_gauss_pts, false);
auto t_hess_sdf = getFTensor2SymmetricFromMat<3>(m_hess_sdf);
for (size_t gg = 0; gg < nb_gauss_pts; ++gg) {
t_hess_sdf(
i,
j) = t_hess_sdf_set(
i,
j);
++t_hess_sdf;
}
return m_hess_sdf;
}
template <int DIM, IntegrationType I, typename BoundaryEleOp>
struct OpAssembleTotalContactTractionImpl;
template <int DIM, IntegrationType I, typename BoundaryEleOp>
struct OpAssembleTotalContactAreaImpl;
template <int DIM, IntegrationType I, typename BoundaryEleOp>
struct OpEvaluateSDFImpl;
template <int DIM, IntegrationType I, typename AssemblyBoundaryEleOp>
struct OpConstrainBoundaryRhsImpl;
template <int DIM, IntegrationType I, typename AssemblyBoundaryEleOp>
struct OpConstrainBoundaryLhs_dUImpl;
template <int DIM, IntegrationType I, typename AssemblyBoundaryEleOp>
struct OpConstrainBoundaryLhs_dTractionImpl;
template <typename T1, typename T2, int DIM1, int DIM2>
size_t nb_gauss_pts) {
m_spatial_coords.clear();
auto t_spatial_coords = getFTensor1FromPtr<3>(&m_spatial_coords(0, 0));
for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
t_spatial_coords(
i) = t_coords(
i) + t_disp(
i);
++t_spatial_coords;
++t_coords;
++t_disp;
}
return m_spatial_coords;
}
template <typename T1, int DIM1>
size_t nb_gauss_pts) {
m_normals_at_pts.clear();
auto t_set_normal = getFTensor1FromMat<3>(m_normals_at_pts);
for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
t_set_normal(
i) = t_normal_at_pts(
i) / t_normal_at_pts.l2();
++t_set_normal;
++t_normal_at_pts;
}
return m_normals_at_pts;
}
template <int DIM, typename BoundaryEleOp>
OpAssembleTotalContactTractionImpl(
boost::shared_ptr<CommonData> common_data_ptr,
double scale = 1,
private:
boost::shared_ptr<CommonData> commonDataPtr;
const double scaleTraction;
bool isAxisymmetric;
};
template <int DIM, typename BoundaryEleOp>
OpAssembleTotalContactAreaImpl(
boost::shared_ptr<CommonData> common_data_ptr,
boost::shared_ptr<Range> contact_range_ptr = nullptr);
private:
boost::shared_ptr<CommonData> commonDataPtr;
bool isAxisymmetric;
boost::shared_ptr<Range> contactRange;
};
template <int DIM, typename BoundaryEleOp>
OpEvaluateSDFImpl(boost::shared_ptr<CommonData> common_data_ptr);
private:
boost::shared_ptr<CommonData> commonDataPtr;
};
template <int DIM, typename AssemblyBoundaryEleOp>
OpConstrainBoundaryRhsImpl(
const std::string
field_name,
boost::shared_ptr<CommonData> common_data_ptr,
private:
boost::shared_ptr<CommonData> commonDataPtr;
bool isAxisymmetric;
};
template <int DIM, typename AssemblyBoundaryEleOp>
OpConstrainBoundaryLhs_dUImpl(const std::string row_field_name,
const std::string col_field_name,
boost::shared_ptr<CommonData> common_data_ptr,
boost::shared_ptr<CommonData> commonDataPtr;
bool isAxisymmetric;
};
template <int DIM, typename AssemblyBoundaryEleOp>
OpConstrainBoundaryLhs_dTractionImpl(
const std::string row_field_name, const std::string col_field_name,
boost::shared_ptr<CommonData> common_data_ptr,
private:
boost::shared_ptr<CommonData> commonDataPtr;
bool isAxisymmetric;
};
template <typename BoundaryEleOp> struct ContactIntegrators {
template <int DIM, IntegrationType I>
OpAssembleTotalContactTractionImpl<DIM, I, BoundaryEleOp>;
template <int DIM, IntegrationType I>
OpAssembleTotalContactAreaImpl<DIM, I, BoundaryEleOp>;
template <int DIM, IntegrationType I>
template <AssemblyType A> struct Assembly {
template <int DIM, IntegrationType I>
OpConstrainBoundaryRhsImpl<DIM, I, AssemblyBoundaryEleOp>;
template <int DIM, IntegrationType I>
OpConstrainBoundaryLhs_dUImpl<DIM, I, AssemblyBoundaryEleOp>;
template <int DIM, IntegrationType I>
OpConstrainBoundaryLhs_dTractionImpl<DIM, I, AssemblyBoundaryEleOp>;
};
};
inline double sign(
double x) {
constexpr
auto eps = std::numeric_limits<float>::epsilon();
return 0;
return 1;
else
return -1;
};
inline double w(
const double sdf,
const double tn) {
}
return (1 - s) / 2;
}
template <int DIM, typename BoundaryEleOp>
OpAssembleTotalContactTractionImpl<DIM, GAUSS, BoundaryEleOp>::
OpAssembleTotalContactTractionImpl(
boost::shared_ptr<CommonData> common_data_ptr,
double scale,
commonDataPtr(common_data_ptr), scaleTraction(
scale),
template <int DIM, typename BoundaryEleOp>
auto t_w = BoundaryEleOp::getFTensor0IntegrationWeight();
auto t_traction = getFTensor1FromMat<DIM>(commonDataPtr->contactTraction);
auto t_coords = BoundaryEleOp::getFTensor1CoordsAtGaussPts();
const auto nb_gauss_pts = BoundaryEleOp::getGaussPts().size2();
for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
double jacobian = 1.;
if (isAxisymmetric) {
jacobian = 2. * M_PI * t_coords(0);
}
const double alpha = t_w * jacobian * BoundaryEleOp::getMeasure();
t_sum_t(
i) += alpha * t_traction(
i);
++t_w;
++t_traction;
++t_coords;
}
t_sum_t(
i) *= scaleTraction;
constexpr
int ind[] = {0, 1, 2};
ADD_VALUES);
}
template <int DIM, typename BoundaryEleOp>
OpAssembleTotalContactAreaImpl<DIM, GAUSS, BoundaryEleOp>::
OpAssembleTotalContactAreaImpl(
boost::shared_ptr<Range> contact_range_ptr)
contactRange(contact_range_ptr) {}
template <int DIM, typename BoundaryEleOp>
auto fe_type = BoundaryEleOp::getFEType();
const auto fe_ent = BoundaryEleOp::getFEEntityHandle();
if (contactRange->find(fe_ent) != contactRange->end()) {
auto t_w = BoundaryEleOp::getFTensor0IntegrationWeight();
auto t_traction = getFTensor1FromMat<DIM>(commonDataPtr->contactTraction);
auto t_coords = BoundaryEleOp::getFTensor1CoordsAtGaussPts();
auto t_grad = getFTensor2FromMat<DIM, DIM>(commonDataPtr->contactDispGrad);
auto t_normal_at_pts = BoundaryEleOp::getFTensor1NormalsAtGaussPts();
const auto nb_gauss_pts = BoundaryEleOp::getGaussPts().size2();
BoundaryEleOp::getFTensor1CoordsAtGaussPts(),
getFTensor1FromMat<DIM>(commonDataPtr->contactDisp), nb_gauss_pts);
BoundaryEleOp::getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
auto t_normal = getFTensor1FromMat<3>(m_normals_at_pts);
auto ts_time = BoundaryEleOp::getTStime();
auto ts_time_step = BoundaryEleOp::getTStimeStep();
int block_id = 0;
auto v_sdf =
surfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
m_spatial_coords, m_normals_at_pts, block_id);
auto m_grad_sdf = gradSurfaceDistanceFunction(
ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
block_id);
auto t_grad_sdf = getFTensor1FromMat<3>(m_grad_sdf);
for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
double jacobian = 1.;
if (isAxisymmetric) {
jacobian = 2. * M_PI * t_coords(0);
}
auto tn = -t_traction(
i) * t_grad_sdf(
i);
double alpha = t_w * jacobian;
t_normal_current(
i) = det * (invF(
j,
i) * t_normal_at_pts(
j));
alpha *= sqrt(t_normal_current(
i) * t_normal_current(
i));
if (fe_type == MBTRI) {
alpha /= 2;
}
t_sum_a(0) += alpha;
}
t_sum_a(1) += alpha;
++t_w;
++t_traction;
++t_coords;
++t_sdf;
++t_grad_sdf;
++t_grad;
++t_normal_at_pts;
}
constexpr
int ind[] = {3, 4};
ADD_VALUES);
}
}
template <int DIM, typename BoundaryEleOp>
boost::shared_ptr<CommonData> common_data_ptr)
commonDataPtr(common_data_ptr) {}
template <int DIM, typename BoundaryEleOp>
const auto nb_gauss_pts = BoundaryEleOp::getGaussPts().size2();
auto &sdf_vec = commonDataPtr->sdfVals;
auto &grad_mat = commonDataPtr->gradsSdf;
auto &hess_mat = commonDataPtr->hessSdf;
auto &constraint_vec = commonDataPtr->constraintVals;
auto &contactTraction_mat = commonDataPtr->contactTraction;
sdf_vec.resize(nb_gauss_pts, false);
grad_mat.resize(DIM, nb_gauss_pts, false);
hess_mat.resize((DIM * (DIM + 1)) / 2, nb_gauss_pts, false);
constraint_vec.resize(nb_gauss_pts, false);
auto t_traction = getFTensor1FromMat<DIM>(contactTraction_mat);
auto t_grad_sdf = getFTensor1FromMat<DIM>(grad_mat);
auto t_hess_sdf = getFTensor2SymmetricFromMat<DIM>(hess_mat);
auto t_disp = getFTensor1FromMat<DIM>(commonDataPtr->contactDisp);
auto t_coords = BoundaryEleOp::getFTensor1CoordsAtGaussPts();
auto t_normal_at_pts = BoundaryEleOp::getFTensor1NormalsAtGaussPts();
auto ts_time = BoundaryEleOp::getTStime();
auto ts_time_step = BoundaryEleOp::getTStimeStep();
BoundaryEleOp::getFTensor1CoordsAtGaussPts(),
getFTensor1FromMat<DIM>(commonDataPtr->contactDisp), nb_gauss_pts);
BoundaryEleOp::getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
int block_id = 0;
auto v_sdf =
surfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
m_spatial_coords, m_normals_at_pts, block_id);
auto m_grad_sdf =
gradSurfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
m_spatial_coords, m_normals_at_pts, block_id);
auto m_hess_sdf =
hessSurfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
m_spatial_coords, m_normals_at_pts, block_id);
auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
auto t_hess_sdf_v = getFTensor2SymmetricFromMat<3>(m_hess_sdf);
auto next = [&]() {
++t_sdf;
++t_sdf_v;
++t_grad_sdf;
++t_grad_sdf_v;
++t_hess_sdf;
++t_hess_sdf_v;
++t_disp;
++t_traction;
++t_constraint;
};
for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
auto tn = -t_traction(
i) * t_grad_sdf_v(
i);
t_sdf = t_sdf_v;
t_grad_sdf(
i) = t_grad_sdf_v(
i);
t_hess_sdf(
i,
j) = t_hess_sdf_v(
i,
j);
next();
}
}
template <int DIM, typename AssemblyBoundaryEleOp>
OpConstrainBoundaryRhsImpl<DIM, GAUSS, AssemblyBoundaryEleOp>::
OpConstrainBoundaryRhsImpl(
const std::string
field_name,
boost::shared_ptr<CommonData> common_data_ptr,
AssemblyBoundaryEleOp::OPROW),
template <int DIM, typename AssemblyBoundaryEleOp>
const size_t nb_gauss_pts = AssemblyBoundaryEleOp::getGaussPts().size2();
auto &nf = AssemblyBoundaryEleOp::locF;
auto t_normal_at_pts = AssemblyBoundaryEleOp::getFTensor1NormalsAtGaussPts();
auto t_w = AssemblyBoundaryEleOp::getFTensor0IntegrationWeight();
auto t_disp = getFTensor1FromMat<DIM>(commonDataPtr->contactDisp);
auto t_traction = getFTensor1FromMat<DIM>(commonDataPtr->contactTraction);
auto t_coords = AssemblyBoundaryEleOp::getFTensor1CoordsAtGaussPts();
size_t nb_base_functions = data.getN().size2() / 3;
auto t_base = data.getFTensor1N<3>();
BoundaryEleOp::getFTensor1CoordsAtGaussPts(),
getFTensor1FromMat<DIM>(commonDataPtr->contactDisp), nb_gauss_pts);
BoundaryEleOp::getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
auto t_normal = getFTensor1FromMat<3>(m_normals_at_pts);
auto ts_time = AssemblyBoundaryEleOp::getTStime();
auto ts_time_step = AssemblyBoundaryEleOp::getTStimeStep();
int block_id = 0;
auto v_sdf =
surfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
m_spatial_coords, m_normals_at_pts, block_id);
auto m_grad_sdf =
gradSurfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
m_spatial_coords, m_normals_at_pts, block_id);
auto t_grad_sdf = getFTensor1FromMat<3>(m_grad_sdf);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
auto t_nf = getFTensor1FromPtr<DIM>(&nf[0]);
double jacobian = 1.;
if (isAxisymmetric) {
jacobian = 2. * M_PI * t_coords(0);
}
const double alpha = t_w * jacobian * AssemblyBoundaryEleOp::getMeasure();
auto tn = -t_traction(
i) * t_grad_sdf(
i);
t_cP(
i,
j) = (
c * t_grad_sdf(
i)) * t_grad_sdf(
j);
+
c * (t_sdf * t_grad_sdf(
i));
size_t bb = 0;
for (; bb != AssemblyBoundaryEleOp::nbRows / DIM; ++bb) {
const double beta = alpha * (t_base(
i) * t_normal(
i));
t_nf(
i) -= beta * t_rhs(
i);
++t_nf;
++t_base;
}
for (; bb < nb_base_functions; ++bb)
++t_base;
++t_disp;
++t_traction;
++t_coords;
++t_w;
++t_normal;
++t_sdf;
++t_grad_sdf;
}
}
template <int DIM, typename AssemblyBoundaryEleOp>
OpConstrainBoundaryLhs_dUImpl<DIM, GAUSS, AssemblyBoundaryEleOp>::
OpConstrainBoundaryLhs_dUImpl(const std::string row_field_name,
const std::string col_field_name,
boost::shared_ptr<CommonData> common_data_ptr,
AssemblyBoundaryEleOp::OPROWCOL),
AssemblyBoundaryEleOp::sYmm = false;
}
template <int DIM, typename AssemblyBoundaryEleOp>
const size_t nb_gauss_pts = AssemblyBoundaryEleOp::getGaussPts().size2();
auto &locMat = AssemblyBoundaryEleOp::locMat;
auto t_normal_at_pts = AssemblyBoundaryEleOp::getFTensor1NormalsAtGaussPts();
auto t_traction = getFTensor1FromMat<DIM>(commonDataPtr->contactTraction);
auto t_coords = AssemblyBoundaryEleOp::getFTensor1CoordsAtGaussPts();
auto t_w = AssemblyBoundaryEleOp::getFTensor0IntegrationWeight();
auto t_row_base = row_data.getFTensor1N<3>();
size_t nb_face_functions = row_data.getN().size2() / 3;
BoundaryEleOp::getFTensor1CoordsAtGaussPts(),
getFTensor1FromMat<DIM>(commonDataPtr->contactDisp), nb_gauss_pts);
BoundaryEleOp::getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
auto t_normal = getFTensor1FromMat<3>(m_normals_at_pts);
auto ts_time = AssemblyBoundaryEleOp::getTStime();
auto ts_time_step = AssemblyBoundaryEleOp::getTStimeStep();
int block_id = 0;
auto v_sdf =
surfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
m_spatial_coords, m_normals_at_pts, block_id);
auto m_grad_sdf =
gradSurfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
m_spatial_coords, m_normals_at_pts, block_id);
auto m_hess_sdf =
hessSurfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
m_spatial_coords, m_normals_at_pts, block_id);
auto t_grad_sdf = getFTensor1FromMat<3>(m_grad_sdf);
auto t_hess_sdf = getFTensor2SymmetricFromMat<3>(m_hess_sdf);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
double jacobian = 1.;
if (isAxisymmetric) {
jacobian = 2. * M_PI * t_coords(0);
}
const double alpha = t_w * jacobian * AssemblyBoundaryEleOp::getMeasure();
auto tn = -t_traction(
i) * t_grad_sdf(
i);
t_cP(
i,
j) = (
c * t_grad_sdf(
i)) * t_grad_sdf(
j);
(t_hess_sdf(
i,
j) * (t_grad_sdf(
k) * t_traction(
k)) +
t_grad_sdf(
i) * t_hess_sdf(
k,
j) * t_traction(
k)) +
c * t_sdf * t_hess_sdf(
i,
j);
}
size_t rr = 0;
for (; rr != AssemblyBoundaryEleOp::nbRows / DIM; ++rr) {
auto t_mat = getFTensor2FromArray<DIM, DIM, DIM>(locMat, DIM * rr);
const double row_base = t_row_base(
i) * t_normal(
i);
auto t_col_base = col_data.getFTensor0N(gg, 0);
for (size_t cc = 0; cc != AssemblyBoundaryEleOp::nbCols / DIM; ++cc) {
const double beta = alpha * row_base * t_col_base;
t_mat(
i,
j) -= beta * t_res_dU(
i,
j);
++t_col_base;
++t_mat;
}
++t_row_base;
}
for (; rr < nb_face_functions; ++rr)
++t_row_base;
++t_traction;
++t_coords;
++t_w;
++t_normal;
++t_sdf;
++t_grad_sdf;
++t_hess_sdf;
}
}
template <int DIM, typename AssemblyBoundaryEleOp>
OpConstrainBoundaryLhs_dTractionImpl<DIM, GAUSS, AssemblyBoundaryEleOp>::
OpConstrainBoundaryLhs_dTractionImpl(
const std::string row_field_name, const std::string col_field_name,
AssemblyBoundaryEleOp::OPROWCOL),
AssemblyBoundaryEleOp::sYmm = false;
}
template <int DIM, typename AssemblyBoundaryEleOp>
OpConstrainBoundaryLhs_dTractionImpl<DIM, GAUSS, AssemblyBoundaryEleOp>::
const size_t nb_gauss_pts = AssemblyBoundaryEleOp::getGaussPts().size2();
auto &locMat = AssemblyBoundaryEleOp::locMat;
auto t_normal_at_pts = AssemblyBoundaryEleOp::getFTensor1NormalsAtGaussPts();
auto t_traction = getFTensor1FromMat<DIM>(commonDataPtr->contactTraction);
auto t_coords = AssemblyBoundaryEleOp::getFTensor1CoordsAtGaussPts();
auto t_w = AssemblyBoundaryEleOp::getFTensor0IntegrationWeight();
auto t_row_base = row_data.getFTensor1N<3>();
size_t nb_face_functions = row_data.getN().size2() / 3;
BoundaryEleOp::getFTensor1CoordsAtGaussPts(),
getFTensor1FromMat<DIM>(commonDataPtr->contactDisp), nb_gauss_pts);
BoundaryEleOp::getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
auto t_normal = getFTensor1FromMat<3>(m_normals_at_pts);
auto ts_time = AssemblyBoundaryEleOp::getTStime();
auto ts_time_step = AssemblyBoundaryEleOp::getTStimeStep();
int block_id = 0;
auto v_sdf =
surfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
m_spatial_coords, m_normals_at_pts, block_id);
auto m_grad_sdf =
gradSurfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
m_spatial_coords, m_normals_at_pts, block_id);
auto t_grad_sdf = getFTensor1FromMat<3>(m_grad_sdf);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
double jacobian = 1.;
if (isAxisymmetric) {
jacobian = 2. * M_PI * t_coords(0);
}
const double alpha = t_w * jacobian * AssemblyBoundaryEleOp::getMeasure();
auto tn = -t_traction(
i) * t_grad_sdf(
i);
t_cP(
i,
j) = (
c * t_grad_sdf(
i)) * t_grad_sdf(
j);
size_t rr = 0;
for (; rr != AssemblyBoundaryEleOp::nbRows / DIM; ++rr) {
auto t_mat = getFTensor2FromArray<DIM, DIM, DIM>(locMat, DIM * rr);
const double row_base = t_row_base(
i) * t_normal(
i);
auto t_col_base = col_data.getFTensor1N<3>(gg, 0);
for (size_t cc = 0; cc != AssemblyBoundaryEleOp::nbCols / DIM; ++cc) {
const double col_base = t_col_base(
i) * t_normal(
i);
const double beta = alpha * row_base * col_base;
t_mat(
i,
j) -= beta * t_res_dt(
i,
j);
++t_col_base;
++t_mat;
}
++t_row_base;
}
for (; rr < nb_face_functions; ++rr)
++t_row_base;
++t_traction;
++t_coords;
++t_w;
++t_normal;
++t_sdf;
++t_grad_sdf;
}
}
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
A>::template LinearForm<I>;
using OpMixDivURhs = typename B::template OpMixDivTimesU<3, DIM, DIM>;
using OpMixDivUCylRhs =
typename B::template OpMixDivTimesU<3, DIM, DIM, CYLINDRICAL>;
using OpMixLambdaGradURhs = typename B::template OpMixTensorTimesGradU<DIM>;
using OpMixUTimesDivLambdaRhs =
typename B::template OpMixVecTimesDivLambda<SPACE_DIM>;
using OpMixUTimesLambdaRhs =
auto common_data_ptr = boost::make_shared<ContactOps::CommonData>();
auto mat_grad_ptr = boost::make_shared<MatrixDouble>();
auto div_stress_ptr = boost::make_shared<MatrixDouble>();
auto contact_stress_ptr = boost::make_shared<MatrixDouble>();
else
return 1.;
};
pip.push_back(new OpCalculateVectorFieldValues<DIM>(
u, common_data_ptr->contactDispPtr()));
pip.push_back(
new OpCalculateHVecTensorField<DIM, DIM>(sigma, contact_stress_ptr));
pip.push_back(
new OpCalculateHVecTensorDivergence<DIM, DIM>(sigma, div_stress_ptr));
} else {
pip.push_back(new OpCalculateHVecTensorDivergence<DIM, DIM, CYLINDRICAL>(
sigma, div_stress_ptr));
}
pip.push_back(
new OpMixDivURhs(sigma, common_data_ptr->contactDispPtr(), jacobian));
} else {
pip.push_back(new OpMixDivUCylRhs(sigma, common_data_ptr->contactDispPtr(),
jacobian));
}
pip.push_back(new OpMixLambdaGradURhs(sigma, mat_grad_ptr, jacobian));
pip.push_back(new OpMixUTimesDivLambdaRhs(u, div_stress_ptr, jacobian));
pip.push_back(new OpMixUTimesLambdaRhs(u, contact_stress_ptr, jacobian));
}
template <
typename OpMixLhs>
struct OpMixLhsSide :
public OpMixLhs {
using OpMixLhs::OpMixLhs;
EntityType col_type,
auto side_fe_entity = OpMixLhs::getSidePtrFE()->getFEEntityHandle();
auto side_fe_data = OpMixLhs::getSideEntity(row_side, row_type);
if (side_fe_entity == side_fe_data) {
CHKERR OpMixLhs::doWork(row_side, col_side, row_type, col_type, row_data,
col_data);
}
}
};
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEle>
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string fe_domain_name, std::string sigma, std::string u,
std::string geom, ForcesAndSourcesCore::RuleHookFun rule,
auto op_loop_side = new OpLoopSide<DomainEle>(
m_field, fe_domain_name, DIM, Sev::noisy,
boost::make_shared<ForcesAndSourcesCore::UserDataOperator::AdjCache>());
pip.push_back(op_loop_side);
CHKERR AddHOOps<DIM, DIM, DIM>::add(op_loop_side->getOpPtrVector(),
{H1, HDIV}, geom);
using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
using OpMixDivULhs = typename B::template OpMixDivTimesVec<DIM>;
using OpMixDivUCylLhs =
typename B::template OpMixDivTimesVec<DIM, CYLINDRICAL>;
using OpLambdaGraULhs = typename B::template OpMixTensorTimesGrad<DIM>;
using OpMixDivULhsSide = OpMixLhsSide<OpMixDivULhs>;
using OpMixDivUCylLhsSide = OpMixLhsSide<OpMixDivUCylLhs>;
using OpLambdaGraULhsSide = OpMixLhsSide<OpLambdaGraULhs>;
auto unity = []() { return 1; };
else
return 1.;
};
op_loop_side->getOpPtrVector().push_back(
new OpMixDivULhsSide(sigma, u, unity, jacobian, true));
} else {
op_loop_side->getOpPtrVector().push_back(
new OpMixDivUCylLhsSide(sigma, u, unity, jacobian, true));
}
op_loop_side->getOpPtrVector().push_back(
new OpLambdaGraULhsSide(sigma, u, unity, jacobian, true));
op_loop_side->getSideFEPtr()->getRuleHook = rule;
}
template <int DIM, AssemblyType A, IntegrationType I, typename BoundaryEleOp>
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
using C = ContactIntegrators<BoundaryEleOp>;
auto common_data_ptr = boost::make_shared<ContactOps::CommonData>();
pip.push_back(new OpCalculateVectorFieldValues<DIM>(
u, common_data_ptr->contactDispPtr()));
pip.push_back(new OpCalculateHVecTensorTrace<DIM, BoundaryEleOp>(
sigma, common_data_ptr->contactTractionPtr()));
pip.push_back(
new typename C::template Assembly<A>::template OpConstrainBoundaryLhs_dU<
pip.push_back(new typename C::template Assembly<A>::
template OpConstrainBoundaryLhs_dTraction<DIM, GAUSS>(
}
template <int DIM, AssemblyType A, IntegrationType I, typename BoundaryEleOp>
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
using C = ContactIntegrators<BoundaryEleOp>;
auto common_data_ptr = boost::make_shared<ContactOps::CommonData>();
pip.push_back(new OpCalculateVectorFieldValues<DIM>(
u, common_data_ptr->contactDispPtr()));
pip.push_back(new OpCalculateHVecTensorTrace<DIM, BoundaryEleOp>(
sigma, common_data_ptr->contactTractionPtr()));
pip.push_back(
new typename C::template Assembly<A>::template OpConstrainBoundaryRhs<
}
template <int DIM, IntegrationType I, typename BoundaryEleOp>
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
using C = ContactIntegrators<BoundaryEleOp>;
auto common_data_ptr = boost::make_shared<ContactOps::CommonData>();
pip.push_back(new OpCalculateHVecTensorTrace<DIM, BoundaryEleOp>(
sigma, common_data_ptr->contactTractionPtr()));
pip.push_back(new typename C::template OpAssembleTotalContactTraction<DIM, I>(
}
template <int DIM, IntegrationType I, typename BoundaryEleOp>
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
OpLoopSide<SideEle> *op_loop_side, std::string sigma, std::string u,
boost::shared_ptr<Range> contact_range_ptr = nullptr) {
using C = ContactIntegrators<BoundaryEleOp>;
auto common_data_ptr = boost::make_shared<ContactOps::CommonData>();
op_loop_side->getOpPtrVector().push_back(
"U", common_data_ptr->contactDispGradPtr()));
if (contact_range_ptr) {
pip.push_back(new OpCalculateVectorFieldValues<DIM>(
u, common_data_ptr->contactDispPtr()));
pip.push_back(new OpCalculateHVecTensorTrace<DIM, BoundaryEleOp>(
sigma, common_data_ptr->contactTractionPtr()));
pip.push_back(op_loop_side);
pip.push_back(new typename C::template OpAssembleTotalContactArea<DIM, I>(
}
}
};
#endif // __CONTACTOPS_HPP__