Browse Source

add irregular model

multiple_top
cflin 2 years ago
parent
commit
7a8a7d0311
  1. 9
      CMakeLists.txt
  2. 16
      main.cpp
  3. 15
      src/Boundary.h
  4. 38
      src/IrregularMesh.h
  5. 32
      src/Mesh.h
  6. 34
      src/Top3d.cpp
  7. 1
      src/Top3d.h
  8. 7
      src/Util.h

9
CMakeLists.txt

@ -54,7 +54,7 @@ target_link_libraries(${PROJECT_NAME}_lib PUBLIC TBB::tbb)
if (USE_SUITESPARSE)
# SuiteSparse
find_package(SuiteSparse REQUIRED)
target_compile_definitions(${PROJECT_NAME} PUBLIC USE_SUITESPARSE)
target_compile_definitions(${PROJECT_NAME}_lib PUBLIC USE_SUITESPARSE)
target_link_libraries(${PROJECT_NAME}_lib PUBLIC ${SUITESPARSE_LIBRARIES})
target_include_directories(${PROJECT_NAME}_lib PUBLIC ${SUITESPARSE_INCLUDE_DIRS})
endif ()
@ -69,8 +69,5 @@ add_subdirectory(3rd/mma)
target_link_libraries(${PROJECT_NAME}_lib PUBLIC mma::mma mma::gcmma)
target_compile_definitions(${PROJECT_NAME}_lib PUBLIC CMAKE_SOURCE_DIR="${CMAKE_SOURCE_DIR}")
target_compile_definitions(${PROJECT_NAME} PUBLIC
WRITE_TENSOR
DEBUG
)
target_compile_definitions(${PROJECT_NAME}_lib PUBLIC WRITE_TENSOR)
target_compile_definitions(${PROJECT_NAME}_lib PUBLIC DEBUG)

16
main.cpp

@ -4,6 +4,8 @@
#include <Eigen/Eigen>
#include "Top3d.h"
#include "Boundary.h"
#include "IrregularMesh.h"
using Eigen::Tensor;
using Tensor3f = Tensor<float, 3>;
@ -23,15 +25,23 @@ int main() {
double E = 1.0;
double Poisson_ratio = 0.3;
auto material = std::make_shared<SIM::Material>(E, Poisson_ratio);
auto mesh = std::make_shared<top::Mesh>(30, 120,40);
#ifdef DEBUG
auto mesh = std::make_shared<top::IrregularMesh>(top::ModelMesh());
#else
auto mesh = std::make_shared<top::Mesh>(15, 60,20);
#endif
top::Top3d top3d(para, material, mesh);
top::Boundary bdr(mesh);
// top3d.AddDBC(bdr.GetCornerPoint(0),{1,1,1});
// top3d.AddDBC(bdr.GetCornerPoint(1),{1,1,1});
// top3d.AddDBC(bdr.GetCornerPoint(2),{1,1,1});
// top3d.AddDBC(bdr.GetCornerPoint(3),{1,1,1});
top3d.AddDBC(bdr.GetLeftBoundary(),{1,1,1});
top3d.AddNBC(bdr.GetRightBottomMidPoint(),{0,0,-1});
Eigen::VectorXi x_seq = Eigen::VectorXi::LinSpaced(21, 0, 20);
Eigen::VectorXi y_seq=x_seq;
Eigen::VectorXi z_seq=Eigen::Vector<int,2>(0,1);
Eigen::VectorXi z_seq_neu=Eigen::Vector<int,2>(19,20);
top3d.AddDBC(bdr.GetChosenCoords(x_seq,y_seq,z_seq), {1, 1, 1});
top3d.AddNBC(bdr.GetChosenCoords(x_seq,y_seq,z_seq_neu), {0, 0, -1});
top3d.TopOptMainLoop();

15
src/Boundary.h

@ -97,10 +97,23 @@ namespace top {
z_sequence.size(), 1);
chosen_coords.col(2) = z_sequence.transpose().replicate(x_sequence.size() * y_sequence.size(),
1).reshaped();
return chosen_coords;
return ReduceInvalidCoords(chosen_coords);
}
private:
Eigen::MatrixXi ReduceInvalidCoords(const Eigen::MatrixXi & coords){
Eigen::VectorXi node_ids=sp_mesh_->MapNodeCoord2Id(coords);
Eigen::MatrixXi masked_coord=coords;
int j=0;
for(int i=0;i<node_ids.rows();++i){
if(node_ids[i]!=-1){
masked_coord.row(j++)=coords.row(i);
}
}
assert(j>0);
return masked_coord.topRows(j);
}
std::shared_ptr<Mesh> sp_mesh_;
};

38
src/IrregularMesh.h

@ -17,10 +17,10 @@ namespace top {
public:
IrregularMesh(const ModelMesh &arbitrary_mesh) : Mesh() {
// get num_pixels
const double percentage_of_min_len = 0.10;// TODO: fixme
const double percentage_of_min_len = 0.05;// TODO: fixme
#ifdef DEBUG
min_point_box_=Eigen::Vector3d(1,1,1)-5.0;
Eigen::Vector3d box_max_point=Eigen::Vector3d(1,1,1)+5.0;
min_point_box_=Eigen::Vector3d(1,1,1).array()-5.0;
Eigen::Vector3d box_max_point=Eigen::Vector3d(1,1,1).array()+5.0;
#else
min_point_box_ = arbitrary_mesh.points.colwise().minCoeff();
Eigen::Vector3d box_max_point = arbitrary_mesh.points.colwise().maxCoeff();
@ -37,10 +37,10 @@ namespace top {
num_node_ = (lx_ + 1) * (ly_ + 1) * (lz_ + 1);
num_ele_ = lx_ * ly_ * lz_;
ten_ele_coord2ele_id = Tensor3i(lx_, ly_, lz_);
ten_node_coord2node_id.setConstant(-1);
ten_node_coord2node_id = Tensor3i(lx_ + 1, ly_ + 1, lz_ + 1);
ten_node_coord2node_id.setConstant(-1);
ten_ele_coord2ele_id_ = Tensor3i(lx_, ly_, lz_);
ten_ele_coord2ele_id_.setConstant(-1);
ten_node_coord2node_id_ = Tensor3i(lx_ + 1, ly_ + 1, lz_ + 1);
ten_node_coord2node_id_.setConstant(-1);
auto EvaluateSDF = [](const Eigen::Vector3d point) {
Eigen::Vector3d center(1, 1, 1);
@ -73,10 +73,10 @@ namespace top {
for (int j = 0; j < ly_; ++j) {
for (int i = 0; i < lx_; ++i) {
if (Is_in_model({i, j, k})) {
ten_ele_coord2ele_id(i, j, k) = cnt_pixel++;
ten_ele_coord2ele_id_(i, j, k) = cnt_pixel++;
for (int di = 0; di < delta_coord.rows(); ++di) {
Eigen::Vector3i cur_delta_coord = delta_coord.row(di);
ten_node_coord2node_id(i + cur_delta_coord(0), j + cur_delta_coord(1),
ten_node_coord2node_id_(i + cur_delta_coord(0), j + cur_delta_coord(1),
k + cur_delta_coord(2)) = 1;
}
}
@ -90,8 +90,8 @@ namespace top {
for (int k = 0; k < lz_ + 1; ++k) {
for (int j = 0; j < ly_ + 1; ++j) {
for (int i = 0; i < lx_ + 1; ++i) {
if (ten_node_coord2node_id(i, j, k) == 1) {
ten_node_coord2node_id(i, j, k) = cnt_node_pixel++;
if (ten_node_coord2node_id_(i, j, k) == 1) {
ten_node_coord2node_id_(i, j, k) = cnt_node_pixel++;
}
}
}
@ -103,13 +103,13 @@ namespace top {
for (int k = 0; k < lz_; ++k) {
for (int j = 0; j < ly_; ++j) {
for (int i = 0; i < lx_; ++i) {
int cur_ele_id = ten_ele_coord2ele_id(i, j, k);
int cur_ele_id = ten_ele_coord2ele_id_(i, j, k);
if (cur_ele_id == -1) {
continue;
}
Eigen::MatrixXi world_node_coords = delta_coord.rowwise() + Eigen::RowVector3i(i, j, k);
assert(world_node_coords.rows() == 8 && world_node_coords.cols() == 3);
Eigen::Vector<int, 8> node_ids = ten_node_coord2node_id(world_node_coords);
Eigen::Vector<int, 8> node_ids = ten_node_coord2node_id_(world_node_coords);
for (int nodi = 0; nodi < 8; ++nodi) {
mat_ele_id2dofs_(cur_ele_id, {3 * nodi, 3 * nodi + 1, 3 * nodi + 2}) = Eigen::Vector3i(
3 * node_ids(nodi), 3 * node_ids(nodi) + 1, 3 * node_ids(nodi) + 2);
@ -117,7 +117,17 @@ namespace top {
}
}
}
// fill valid_ele_idx_
valid_ele_idx_.resize(num_pixel_);
Tensor3i ten_tmp_col=ten_ele_coord2ele_id_.reshape(Eigen::array<int,3>{num_ele_,1,1});
cnt_pixel=0;
for(int i=0;i<ten_tmp_col.size();++i){
if(ten_tmp_col(i,0,0)!=-1){
valid_ele_idx_(cnt_pixel)=i;
++cnt_pixel;
}
}
assert(cnt_pixel==num_pixel_);
}
int GetNumDofs() const override {

32
src/Mesh.h

@ -14,17 +14,18 @@ namespace top {
Mesh(int len_x, int len_y, int len_z) : lx_(len_x), ly_(len_y), lz_(len_z),
num_node_((lx_ + 1) * (ly_ + 1) * (lz_ + 1)),
num_ele_(lx_ * ly_ * lz_) {
ten_node_coord2node_id = Tensor3i(num_node_, 1, 1);
valid_ele_idx_=Eigen::VectorXi::LinSpaced(num_ele_,0,num_ele_-1);
ten_node_coord2node_id_ = Tensor3i(num_node_, 1, 1);
for (int i = 0; i < num_node_; ++i) {
ten_node_coord2node_id(i, 0, 0) = i;
ten_node_coord2node_id_(i, 0, 0) = i;
}
ten_node_coord2node_id=ten_node_coord2node_id.reshape(Eigen::array<int, 3>{lx_ + 1, ly_ + 1, lz_ + 1});
ten_node_coord2node_id_=ten_node_coord2node_id_.reshape(Eigen::array<int, 3>{lx_ + 1, ly_ + 1, lz_ + 1});
ten_ele_coord2ele_id = Tensor3i(num_ele_, 1, 1);
ten_ele_coord2ele_id_ = Tensor3i(num_ele_, 1, 1);
for (int i = 0; i < num_ele_; ++i) {
ten_ele_coord2ele_id(i, 0, 0) = i;
ten_ele_coord2ele_id_(i, 0, 0) = i;
}
ten_ele_coord2ele_id=ten_ele_coord2ele_id.reshape(Eigen::array<int, 3>{lx_, ly_, lz_});
ten_ele_coord2ele_id_=ten_ele_coord2ele_id_.reshape(Eigen::array<int, 3>{lx_, ly_, lz_});
mat_ele_id2dofs_.resize(num_ele_, NUM_NODES_EACH_ELE*3);
static const Eigen::MatrixXi delta_coord = (Eigen::MatrixXi(8, 3) <<
@ -40,10 +41,10 @@ namespace top {
for (int k = 0; k < lz_; ++k) {
for (int j = 0; j < ly_; ++j) {
for (int i = 0; i < lx_; ++i) {
int cur_ele_id = ten_ele_coord2ele_id(i, j, k);
int cur_ele_id = ten_ele_coord2ele_id_(i, j, k);
Eigen::MatrixXi world_node_coords=delta_coord.rowwise()+Eigen::RowVector3i(i,j,k);
assert(world_node_coords.rows()==8 && world_node_coords.cols()==3);
Eigen::Vector<int,8> node_ids=ten_node_coord2node_id(world_node_coords);
Eigen::Vector<int,8> node_ids=ten_node_coord2node_id_(world_node_coords);
for(int nodi=0;nodi<8;++nodi){
mat_ele_id2dofs_(cur_ele_id, {3*nodi,3*nodi+1,3*nodi+2}) = Eigen::Vector3i(3*node_ids(nodi),3*node_ids(nodi)+1,3*node_ids(nodi)+2);
}
@ -64,6 +65,9 @@ namespace top {
int GetLz()const{
return lz_;
}
Eigen::VectorXi GetPixelIdx()const{
return valid_ele_idx_;
}
virtual int GetNumDofs() const{
return num_node_*3;
@ -77,12 +81,15 @@ namespace top {
Eigen::MatrixXi GetEleId2DofsMap()const{
return mat_ele_id2dofs_;
}
Tensor3i GetEleCoord2IdTensor()const{
return ten_ele_coord2ele_id_;
}
Eigen::VectorXi MapNodeCoord2Id(const Eigen::MatrixXi & node_coords)const{
return ten_node_coord2node_id(node_coords);
return ten_node_coord2node_id_(node_coords);
}
Eigen::VectorXi MapEleCoord2Id(const Eigen::MatrixXi& ele_coords)const{
return ten_ele_coord2ele_id(ele_coords);
return ten_ele_coord2ele_id_(ele_coords);
}
Eigen::MatrixXi MapEleId2Dofs(const Eigen::VectorXi &ele_ids)const{
return mat_ele_id2dofs_(ele_ids, all);
@ -98,8 +105,9 @@ namespace top {
int lz_;
int num_ele_;
int num_node_;
Tensor3i ten_node_coord2node_id;
Tensor3i ten_ele_coord2ele_id;
Eigen::VectorXi valid_ele_idx_;
Tensor3i ten_node_coord2node_id_;
Tensor3i ten_ele_coord2ele_id_;
Eigen::MatrixXi mat_ele_id2dofs_;// num_ele_x8
};

34
src/Top3d.cpp

@ -6,10 +6,9 @@
namespace top {
void Top3d::TopOptMainLoop() {
Tensor3d x = Tensor3d(sp_mesh_->GetLx(), sp_mesh_->GetLy(), sp_mesh_->GetLz());
x.setConstant(sp_para_->volfrac);
Tensor3d xPhys = x;
Eigen::VectorXd xPhys_col = Eigen::Map<Eigen::VectorXd>(xPhys.data(), xPhys.size());
Eigen::VectorXd xPhys_col(sp_mesh_->GetNumEles());
xPhys_col.setConstant(sp_para_->volfrac);
int loop = 0;
double change = 1.0;
double E0 = sp_material_->YM;
@ -19,6 +18,9 @@ namespace top {
Eigen::VectorXd dv(sp_mesh_->GetNumEles());
dv.setOnes();
dv = H_ * (dv.array() / Hs_.array()).matrix().eval();
Eigen::VectorXd ele_to_write = Eigen::VectorXd::Zero(sp_mesh_->GetLx() * sp_mesh_->GetLy() * sp_mesh_->GetLz());
Eigen::VectorXi pixel_idx = sp_mesh_->GetPixelIdx();
spdlog::info("end precompute");
// clear output dir
clear_dir(CMAKE_SOURCE_DIR "/output");
@ -74,22 +76,24 @@ namespace top {
xPhys_col = variables_tmp;
spdlog::critical("Iter: {:3d}, Comp: {:.3f}, Vol: {:.2f}, Change: {:.3f}", loop, c, v, change);
#ifdef WRITE_TENSOR
Tensor3d ten_xPhys_to_write(sp_mesh_->GetNumEles(), 1, 1);
for (int i = 0; i < xPhys_col.size(); ++i) {
ten_xPhys_to_write(i, 0, 0) = xPhys_col(i);
// result
ele_to_write(pixel_idx)=xPhys_col;
Tensor3d ten_xPhys_to_write(sp_mesh_->GetLx()*sp_mesh_->GetLy()*sp_mesh_->GetLz(),1,1);
for(int i=0;i<ele_to_write.size();++i){
ten_xPhys_to_write(i,0,0)=ele_to_write(i);
}
ten_xPhys_to_write = ten_xPhys_to_write.reshape(
Eigen::array<int, 3>{sp_mesh_->GetLx(), sp_mesh_->GetLy(), sp_mesh_->GetLz()});
write_tensor3d(CMAKE_SOURCE_DIR "/output/txt/iter_" + std::to_string(loop) + ".txt",
ten_xPhys_to_write);
write_tensor3d(CMAKE_SOURCE_DIR "/output/txt/iter_" + std::to_string(loop) + ".txt", ten_xPhys_to_write);
write_tensor3d_to_vtk(CMAKE_SOURCE_DIR "/output/vtk/iter_" + std::to_string(loop) + ".vtk",
ten_xPhys_to_write);
#endif
}
// result
Tensor3d ten_xPhys_to_write(sp_mesh_->GetNumEles(), 1, 1);
for (int i = 0; i < xPhys_col.size(); ++i) {
ten_xPhys_to_write(i, 0, 0) = xPhys_col(i);
ele_to_write(pixel_idx) = xPhys_col;
Tensor3d ten_xPhys_to_write(sp_mesh_->GetLx() * sp_mesh_->GetLy() * sp_mesh_->GetLz(), 1, 1);
for (int i = 0; i < ele_to_write.size(); ++i) {
ten_xPhys_to_write(i, 0, 0) = ele_to_write(i);
}
ten_xPhys_to_write = ten_xPhys_to_write.reshape(
Eigen::array<int, 3>{sp_mesh_->GetLx(), sp_mesh_->GetLy(), sp_mesh_->GetLz()});
@ -148,11 +152,7 @@ namespace top {
}
}
std::vector<Eigen::Triplet<double>> v_tri;
for (int i = 0; i < cnt; ++i) {
v_tri.push_back({iH(i), jH(i), sH(i)});
}
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());
}

1
src/Top3d.h

@ -51,7 +51,6 @@ namespace top {
if (directions[dir])
dofs_to_fixed_.insert(node_id_to_fix[i] * 3 + dir);
}
auto ttt = dofs_to_fixed_;
}
void AddNBC(const Eigen::MatrixXi &NBC_coords, const Eigen::Vector3d &forces) {

7
src/Util.h

@ -4,6 +4,7 @@
#ifndef TOP3D_UTIL_H
#define TOP3D_UTIL_H
#include <unsupported/Eigen/CXX11/Tensor>
#include <Eigen/Eigen>
#include "TensorWrapper.h"
@ -20,6 +21,7 @@
#include <spdlog/spdlog.h>
#include <iostream>
#include <filesystem>
namespace top {
using Tensor3d = TensorWrapper<double, 3>;
using Tensor3i = TensorWrapper<int, 3>;
@ -27,7 +29,8 @@ namespace top {
using SpMat = Eigen::SparseMatrix<double>;
template<typename Scalar>
inline std::vector<Eigen::Triplet<Scalar>> Vec2Triplet(const Eigen::VectorXi& I,const Eigen::VectorXi& J,const Eigen::Matrix<Scalar,-1,1>& V){
inline std::vector<Eigen::Triplet<Scalar>>
Vec2Triplet(const Eigen::VectorXi &I, const Eigen::VectorXi &J, const Eigen::Matrix<Scalar, -1, 1> &V) {
std::vector<Eigen::Triplet<Scalar>> v_tri;
for (int i = 0; i < I.size(); ++i) {
v_tri.push_back({I(i), J(i), V(i)});
@ -101,7 +104,7 @@ namespace top {
}
void clear_dir(const std::filesystem::path& path){
inline void clear_dir(const std::filesystem::path &path) {
for (const auto &entry: std::filesystem::recursive_directory_iterator(path)) {
if (entry.is_regular_file()) {
std::filesystem::remove(entry.path());

Loading…
Cancel
Save