|
|
@ -7,9 +7,10 @@ |
|
|
|
#include <cassert> |
|
|
|
#include "Eigen/src/Core/Matrix.h" |
|
|
|
#include "Util.h" |
|
|
|
|
|
|
|
namespace da::sha { |
|
|
|
namespace top { |
|
|
|
Tensor3d Top3d::TopOptMainLoop() { |
|
|
|
namespace top { |
|
|
|
Tensor3d Top3d::TopOptMainLoop() { |
|
|
|
Eigen::VectorXd xPhys_col(sp_mesh_->GetNumEles()); |
|
|
|
Eigen::VectorXi chosen_ele_id(sp_mesh_->GetChosenEleIdx()); |
|
|
|
bool flg_chosen = chosen_ele_id.size() != 0; |
|
|
@ -47,6 +48,8 @@ Tensor3d Top3d::TopOptMainLoop() { |
|
|
|
// start iteration
|
|
|
|
while (change > sp_para_->tol_x && loop < sp_para_->max_loop) { |
|
|
|
++loop; |
|
|
|
// filter
|
|
|
|
xPhys_col = H_ * (xPhys_col.array() / Hs_.array()).matrix().eval(); |
|
|
|
Eigen::VectorXd sK = |
|
|
|
(sKe_ * (Emin + xPhys_col.array().pow(sp_para_->penal) * (E0 - Emin)).matrix().transpose()) |
|
|
|
.reshaped(); |
|
|
@ -71,8 +74,10 @@ Tensor3d Top3d::TopOptMainLoop() { |
|
|
|
ce.transpose() * (Emin + xPhys_col.array().pow(sp_para_->penal) * (E0 - Emin)).matrix(); |
|
|
|
double v = flg_chosen ? xPhys_col(chosen_ele_id).sum() : xPhys_col.sum(); |
|
|
|
|
|
|
|
Eigen::VectorXd dc = |
|
|
|
// sensitivity
|
|
|
|
Eigen::VectorXd dc_drho = |
|
|
|
-sp_para_->penal * (E0 - Emin) * xPhys_col.array().pow(sp_para_->penal - 1.0) * ce.array(); |
|
|
|
Eigen::VectorXd dc_dx = drho_dx_ * dc_drho; |
|
|
|
|
|
|
|
// mma solver
|
|
|
|
size_t num_constrants = 1; |
|
|
@ -81,8 +86,8 @@ Tensor3d Top3d::TopOptMainLoop() { |
|
|
|
Eigen::VectorXd variables_tmp = flg_chosen ? xPhys_col(chosen_ele_id) : xPhys_col; |
|
|
|
double f0val = c; |
|
|
|
Eigen::VectorXd df0dx = flg_chosen |
|
|
|
? dc(chosen_ele_id).eval() / dc(chosen_ele_id).cwiseAbs().maxCoeff() |
|
|
|
: dc / dc.cwiseAbs().maxCoeff(); |
|
|
|
? dc_dx(chosen_ele_id).eval() / dc_dx(chosen_ele_id).cwiseAbs().maxCoeff() |
|
|
|
: dc_dx / dc_dx.cwiseAbs().maxCoeff(); |
|
|
|
double fval = v - num_variables * sp_para_->volfrac; |
|
|
|
Eigen::VectorXd dfdx = flg_chosen ? dv(chosen_ele_id) : dv; |
|
|
|
static Eigen::VectorXd low_bounds = Eigen::VectorXd::Zero(num_variables); |
|
|
@ -120,7 +125,7 @@ Tensor3d Top3d::TopOptMainLoop() { |
|
|
|
assert(xPhys_col.size()); |
|
|
|
Eigen::VectorXi continue_idx = |
|
|
|
Eigen::VectorXi::LinSpaced(xPhys_col.size(), 0, xPhys_col.size() - 1); |
|
|
|
Eigen::VectorXi unchosen_idx =flg_chosen? SetDifference(continue_idx, chosen_ele_id): Eigen::VectorXi(); |
|
|
|
Eigen::VectorXi unchosen_idx = flg_chosen ? SetDifference(continue_idx, chosen_ele_id) : Eigen::VectorXi(); |
|
|
|
{ |
|
|
|
xPhys_col(unchosen_idx).setZero(); |
|
|
|
ele_to_write(pixel_idx) = xPhys_col; |
|
|
@ -146,15 +151,15 @@ Tensor3d Top3d::TopOptMainLoop() { |
|
|
|
} |
|
|
|
|
|
|
|
return rho_field_zero_filled_; |
|
|
|
} |
|
|
|
} |
|
|
|
|
|
|
|
std::vector<Tensor3d> Top3d::GetTensorOfStress(const Eigen::VectorXd &which_col_of_stress) { |
|
|
|
std::vector<Tensor3d> Top3d::GetTensorOfStress(const Eigen::VectorXd &which_col_of_stress) { |
|
|
|
Eigen::VectorXd ele_to_write = |
|
|
|
Eigen::VectorXd::Zero(sp_mesh_->GetLx() * sp_mesh_->GetLy() * sp_mesh_->GetLz()); |
|
|
|
Eigen::VectorXi pixel_idx = sp_mesh_->GetPixelIdx(); |
|
|
|
// stress
|
|
|
|
Eigen::MatrixXd mat_stress(sp_mesh_->GetNumEles(), 6); |
|
|
|
Eigen::MatrixXd B=sp_fea_->computeBe({0, 0, 0}); |
|
|
|
Eigen::MatrixXd B = sp_fea_->computeBe({0, 0, 0}); |
|
|
|
for (int i = 0; i < sp_mesh_->GetNumEles(); ++i) { |
|
|
|
Eigen::VectorXi dofs_in_ele_i = sp_mesh_->MapEleId2Dofs(i); |
|
|
|
Eigen::VectorXd Ue = U_(dofs_in_ele_i); |
|
|
@ -167,9 +172,9 @@ std::vector<Tensor3d> Top3d::GetTensorOfStress(const Eigen::VectorXd &which_col_ |
|
|
|
vt.push_back(GetTensorFromCol(ele_to_write)); |
|
|
|
} |
|
|
|
return vt; |
|
|
|
} |
|
|
|
} |
|
|
|
|
|
|
|
Tensor3d Top3d::GetTensorFromCol(const Eigen::VectorXd &proprty_col) { |
|
|
|
Tensor3d Top3d::GetTensorFromCol(const Eigen::VectorXd &proprty_col) { |
|
|
|
Tensor3d ten_prop_to_write(sp_mesh_->GetLx() * sp_mesh_->GetLy() * sp_mesh_->GetLz(), 1, 1); |
|
|
|
assert(proprty_col.size() == ten_prop_to_write.size()); |
|
|
|
for (int i = 0; i < proprty_col.size(); ++i) { |
|
|
@ -178,9 +183,9 @@ Tensor3d Top3d::GetTensorFromCol(const Eigen::VectorXd &proprty_col) { |
|
|
|
ten_prop_to_write = ten_prop_to_write.reshape( |
|
|
|
Eigen::array<Eigen::DenseIndex, 3>{sp_mesh_->GetLx(), sp_mesh_->GetLy(), sp_mesh_->GetLz()}); |
|
|
|
return ten_prop_to_write; |
|
|
|
} |
|
|
|
} |
|
|
|
|
|
|
|
void Top3d::Precompute() { |
|
|
|
void Top3d::Precompute() { |
|
|
|
Eigen::MatrixXi mat_ele2dofs = sp_mesh_->GetEleId2DofsMap(); |
|
|
|
int dofs_each_ele = sp_mesh_->Get_NUM_NODES_EACH_ELE() * |
|
|
|
sp_mesh_->Get_DOFS_EACH_NODE(); // 24 for mathe; 8 for heat
|
|
|
@ -190,7 +195,7 @@ void Top3d::Precompute() { |
|
|
|
jK_ = Eigen::KroneckerProduct(mat_ele2dofs, Eigen::RowVectorXi::Ones(dofs_each_ele)) |
|
|
|
.transpose() |
|
|
|
.reshaped(); |
|
|
|
Ke_=sp_fea_->computeKe(1.0); |
|
|
|
Ke_ = sp_fea_->computeKe(1.0); |
|
|
|
sKe_ = Ke_.reshaped(); |
|
|
|
|
|
|
|
// precompute filter
|
|
|
@ -217,14 +222,16 @@ void Top3d::Precompute() { |
|
|
|
for (int i2 = std::max(i - delta, 0); i2 <= std::min(i + delta, sp_mesh_->GetLx() - 1); |
|
|
|
++i2) { |
|
|
|
int ele_id1 = |
|
|
|
sp_mesh_->MapEleCoord2Id((Eigen::MatrixXi(1, 3) << i2, j2, k2).finished())(0); |
|
|
|
sp_mesh_->MapEleCoord2Id((Eigen::MatrixXi(1, 3) << i2, j2, k2).finished())( |
|
|
|
0); |
|
|
|
if (ele_id1 == -1) { |
|
|
|
continue; |
|
|
|
} |
|
|
|
iH(cnt) = ele_id0; |
|
|
|
jH(cnt) = ele_id1; |
|
|
|
sH(cnt) = |
|
|
|
std::max(0.0, sp_para_->r_min - Eigen::Vector3d(i - i2, j - j2, k - k2).norm()); |
|
|
|
std::max(0.0, |
|
|
|
sp_para_->r_min - Eigen::Vector3d(i - i2, j - j2, k - k2).norm()); |
|
|
|
Hs_(ele_id0) += sH(cnt); |
|
|
|
++cnt; |
|
|
|
} |
|
|
@ -236,6 +243,14 @@ void Top3d::Precompute() { |
|
|
|
std::vector<Eigen::Triplet<double>> v_tri = Vec2Triplet(iH, jH, sH); |
|
|
|
H_ = SpMat(sp_mesh_->GetNumEles(), sp_mesh_->GetNumEles()); |
|
|
|
H_.setFromTriplets(v_tri.begin(), v_tri.end()); |
|
|
|
} |
|
|
|
} // namespace top
|
|
|
|
|
|
|
|
Eigen::VectorXi i_Hs = Eigen::VectorXi::LinSpaced(Hs_.size(), 0, |
|
|
|
Hs_.size()); |
|
|
|
Eigen::SparseMatrix<double> sp_inv_Hs(i_Hs.size(), i_Hs.size()); |
|
|
|
auto v_inv_Hs_tri = Vec2Triplet(i_Hs, i_Hs, Eigen::VectorXd(1.0 / Hs_.array())); |
|
|
|
sp_inv_Hs.setFromTriplets(v_inv_Hs_tri.begin(), v_inv_Hs_tri.end()); |
|
|
|
drho_dx_ = sp_inv_Hs * H_; |
|
|
|
|
|
|
|
} |
|
|
|
} // namespace top
|
|
|
|
} // namespace da::sha
|