Browse Source

before modify arbitrary mesh

multiple_top
cflin 2 years ago
parent
commit
a071b40ed8
  1. 5
      CMakeLists.txt
  2. 15
      main.cpp
  3. 2
      src/Boundary.h
  4. 17
      src/IrregularMesh.h
  5. 151
      src/Top3d.cpp
  6. 148
      src/Top3d.h
  7. 11
      src/Util.h

5
CMakeLists.txt

@ -70,4 +70,7 @@ 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}_lib PUBLIC CMAKE_SOURCE_DIR="${CMAKE_SOURCE_DIR}")
target_compile_definitions(${PROJECT_NAME} PUBLIC WRITE_TENSOR) target_compile_definitions(${PROJECT_NAME} PUBLIC
WRITE_TENSOR
DEBUG
)

15
main.cpp

@ -17,20 +17,21 @@ void print_tensor(const Tensor3f &t3) {
int main() { int main() {
auto para = std::make_shared<top::CtrlPara>(); auto para = std::make_shared<top::CtrlPara>();
para->max_loop=200; para->max_loop=200;
para->volfrac=0.2; para->volfrac=0.15;
para->penal=3.0; para->penal=3.0;
para->r_min=1.5; para->r_min=1.5;
double E = 1.0; double E = 1.0;
double Poisson_ratio = 0.3; double Poisson_ratio = 0.3;
auto material = std::make_shared<SIM::Material>(E, Poisson_ratio); auto material = std::make_shared<SIM::Material>(E, Poisson_ratio);
auto mesh = std::make_shared<top::Mesh>(20, 120,20); auto mesh = std::make_shared<top::Mesh>(30, 120,40);
top::Top3d top3d(para, material, mesh); top::Top3d top3d(para, material, mesh);
top::Boundary bdr(mesh); top::Boundary bdr(mesh);
top3d.AddDBC(bdr.GetCornerPoint(0),{1,1,1}); // top3d.AddDBC(bdr.GetCornerPoint(0),{1,1,1});
top3d.AddDBC(bdr.GetCornerPoint(1),{1,1,1}); // top3d.AddDBC(bdr.GetCornerPoint(1),{1,1,1});
top3d.AddDBC(bdr.GetCornerPoint(2),{1,1,1}); // top3d.AddDBC(bdr.GetCornerPoint(2),{1,1,1});
top3d.AddDBC(bdr.GetCornerPoint(3),{1,1,1}); // top3d.AddDBC(bdr.GetCornerPoint(3),{1,1,1});
top3d.AddNBC(bdr.GetTopSurfaceCenter(),{0,0,-1}); top3d.AddDBC(bdr.GetLeftBoundary(),{1,1,1});
top3d.AddNBC(bdr.GetRightBottomMidPoint(),{0,0,-1});
top3d.TopOptMainLoop(); top3d.TopOptMainLoop();

2
src/Boundary.h

@ -50,7 +50,7 @@ namespace top {
return coord; return coord;
} }
Eigen::MatrixXi GetCornerPoint(int i){ Eigen::MatrixXi GetCornerPoint(int i){
bool is_top=i%4; bool is_top=i/4==1;
i%=4; i%=4;
Eigen::MatrixXi coord(1,3); Eigen::MatrixXi coord(1,3);
int lx=sp_mesh_->GetLx(),ly=sp_mesh_->GetLy(),lz=sp_mesh_->GetLz(); int lx=sp_mesh_->GetLx(),ly=sp_mesh_->GetLy(),lz=sp_mesh_->GetLz();

17
src/IrregularMesh.h

@ -17,12 +17,17 @@ namespace top {
public: public:
IrregularMesh(const ModelMesh &arbitrary_mesh) : Mesh() { IrregularMesh(const ModelMesh &arbitrary_mesh) : Mesh() {
// get num_pixels // get num_pixels
const double percentage_of_min_len = 0.10;// 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;
#else
min_point_box_ = arbitrary_mesh.points.colwise().minCoeff(); min_point_box_ = arbitrary_mesh.points.colwise().minCoeff();
Eigen::Vector3d box_max_point = arbitrary_mesh.points.colwise().maxCoeff(); Eigen::Vector3d box_max_point = arbitrary_mesh.points.colwise().maxCoeff();
Eigen::Vector3d box_len = box_max_point - min_point_box_; #endif
const double percentage_of_min_len = 0.20; len_box_ = box_max_point - min_point_box_;
len_pixel_ = box_len.minCoeff() * percentage_of_min_len; len_pixel_ = len_box_.minCoeff() * percentage_of_min_len;
Eigen::Vector3d d_num_pixels = box_len / len_pixel_; Eigen::Vector3d d_num_pixels = len_box_ / len_pixel_;
Eigen::Vector3i num_pixels(std::ceil(d_num_pixels(0)), std::ceil(d_num_pixels(1)), Eigen::Vector3i num_pixels(std::ceil(d_num_pixels(0)), std::ceil(d_num_pixels(1)),
std::ceil(d_num_pixels(2))); std::ceil(d_num_pixels(2)));
@ -61,6 +66,7 @@ namespace top {
1, 1, 1, 1, 1, 1,
0, 1, 1 0, 1, 1
).finished(); ).finished();
// get num_pixel_ && fill pixel id // get num_pixel_ && fill pixel id
int cnt_pixel = 0; int cnt_pixel = 0;
for (int k = 0; k < lz_; ++k) { for (int k = 0; k < lz_; ++k) {
@ -78,6 +84,7 @@ namespace top {
} }
} }
num_pixel_ = cnt_pixel; num_pixel_ = cnt_pixel;
// get_num_node_pixel && fill node_id // get_num_node_pixel && fill node_id
int cnt_node_pixel = 0; int cnt_node_pixel = 0;
for (int k = 0; k < lz_ + 1; ++k) { for (int k = 0; k < lz_ + 1; ++k) {
@ -90,6 +97,7 @@ namespace top {
} }
} }
num_node_pixel_ = cnt_node_pixel; num_node_pixel_ = cnt_node_pixel;
// fill mat_ele_id2dofs_ // fill mat_ele_id2dofs_
mat_ele_id2dofs_.resize(num_pixel_, NUM_NODES_EACH_ELE * 3); mat_ele_id2dofs_.resize(num_pixel_, NUM_NODES_EACH_ELE * 3);
for (int k = 0; k < lz_; ++k) { for (int k = 0; k < lz_; ++k) {
@ -129,6 +137,7 @@ namespace top {
double num_pixel_; double num_pixel_;
double num_node_pixel_; double num_node_pixel_;
Eigen::Vector3d min_point_box_; Eigen::Vector3d min_point_box_;
Eigen::Vector3d len_box_;
}; };

151
src/Top3d.cpp

@ -5,4 +5,155 @@
#include "Top3d.h" #include "Top3d.h"
namespace top { 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());
int loop = 0;
double change = 1.0;
double E0 = sp_material_->YM;
double Emin = sp_material_->YM * sp_para_->E_factor;
// precompute
Eigen::VectorXd dv(sp_mesh_->GetNumEles());
dv.setOnes();
dv = H_ * (dv.array() / Hs_.array()).matrix().eval();
spdlog::info("end precompute");
// clear output dir
clear_dir(CMAKE_SOURCE_DIR "/output");
// start iteration
while (change > sp_para_->tol_x && loop < sp_para_->max_loop) {
++loop;
Eigen::VectorXd sK = (sKe_ * (Emin + xPhys_col.array().pow(sp_para_->penal) *
(E0 - Emin)).matrix().transpose()).reshaped();
auto v_tri = Vec2Triplet(iK_, jK_, sK);
K_.setFromTriplets(v_tri.begin(), v_tri.end());
IntroduceFixedDofs(K_, F_);
#ifdef USE_SUITESPARSE
Eigen::CholmodSupernodalLLT<Eigen::SparseMatrix<double>> solver;
solver.compute(K_);
Eigen::VectorXd U = solver.solve(F_);
#else
// dense direct method
Eigen::MatrixXd denseK = Eigen::MatrixXd(K_);
Eigen::VectorXd U = denseK.fullPivLu().solve(Eigen::VectorXd(F_));
#endif
// compliance
Eigen::VectorXd ce(sp_mesh_->GetNumEles());
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);
ce(i) = Ue.transpose() * Ke_ * Ue;
}
double c = ce.transpose() * (Emin + xPhys_col.array().pow(sp_para_->penal) * (E0 - Emin)).matrix();
double v = xPhys_col.sum();
Eigen::VectorXd dc =
-sp_para_->penal * (E0 - Emin) * xPhys_col.array().pow(sp_para_->penal - 1.0) * ce.array();
// mma solver
size_t num_constrants = 1;
size_t num_variables = sp_mesh_->GetNumEles();
auto mma = std::make_shared<MMASolver>(num_variables, num_constrants);
Eigen::VectorXd variables_tmp = xPhys_col;
double f0val = c;
Eigen::VectorXd df0dx = dc;
double fval = v - sp_mesh_->GetNumEles() * sp_para_->volfrac;
Eigen::VectorXd dfdx = dv;
static Eigen::VectorXd low_bounds = Eigen::VectorXd::Zero(sp_mesh_->GetNumEles());
static Eigen::VectorXd up_bounds = Eigen::VectorXd::Ones(sp_mesh_->GetNumEles());
// spdlog::info("mma update");
mma->Update(variables_tmp.data(), df0dx.data(), &fval, dfdx.data(),
low_bounds.data(), up_bounds.data());
change = (variables_tmp - xPhys_col).cwiseAbs().maxCoeff();
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);
}
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_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);
}
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_to_vtk(CMAKE_SOURCE_DIR "/output/vtk/iter_" + std::to_string(loop) + ".vtk",
ten_xPhys_to_write);
}
void Top3d::Precompute() {
Eigen::MatrixXi mat_ele2dofs = sp_mesh_->GetEleId2DofsMap();
iK_ = Eigen::KroneckerProduct(mat_ele2dofs, Eigen::VectorXi::Ones(24)).transpose().reshaped();
jK_ = Eigen::KroneckerProduct(mat_ele2dofs, Eigen::RowVectorXi::Ones(24)).transpose().reshaped();
sp_material_->computeKe(0.5, 0.5, 0.5, sp_material_->D, Ke_);
sKe_ = Ke_.reshaped();
// precompute filter
Eigen::VectorXi iH = Eigen::VectorXi::Ones(
sp_mesh_->GetNumEles() * std::pow(2.0 * (std::ceil(sp_para_->r_min) - 1.0) + 1, 3));
Eigen::VectorXi jH = iH;
Eigen::VectorXd sH(iH.size());
sH.setZero();
int cnt = 0;
Hs_ = Eigen::VectorXd(sp_mesh_->GetNumEles());
int delta = std::ceil(sp_para_->r_min) - 1;
for (int k = 0; k < sp_mesh_->GetLz(); ++k) {
for (int j = 0; j < sp_mesh_->GetLy(); ++j) {
for (int i = 0; i < sp_mesh_->GetLx(); ++i) {
int ele_id0 = sp_mesh_->MapEleCoord2Id((Eigen::MatrixXi(1, 3) << i, j, k).finished())(0);
if(ele_id0==-1){
continue;
}
for (int k2 = std::max(k - delta, 0); k2 <= std::min(k + delta, sp_mesh_->GetLz() - 1); ++k2) {
for (int j2 = std::max(j - delta, 0);
j2 <= std::min(j + delta, sp_mesh_->GetLy() - 1); ++j2) {
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);
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());
Hs_(ele_id0)+=sH(cnt);
++cnt;
}
}
}
}
}
}
std::vector<Eigen::Triplet<double>> v_tri;
for (int i = 0; i < cnt; ++i) {
v_tri.push_back({iH(i), jH(i), sH(i)});
}
H_ = SpMat(sp_mesh_->GetNumEles(), sp_mesh_->GetNumEles());
H_.setFromTriplets(v_tri.begin(), v_tri.end());
}
} // top } // top

148
src/Top3d.h

@ -4,7 +4,7 @@
#ifndef TOP3D_TOP3D_H #ifndef TOP3D_TOP3D_H
#define TOP3D_TOP3D_H #define TOP3D_TOP3D_H
#include <filesystem>
#include "Util.h" #include "Util.h"
#include "Material.hpp" #include "Material.hpp"
#include "Mesh.h" #include "Mesh.h"
@ -39,7 +39,7 @@ namespace top {
), sp_material_(sp_material), sp_mesh_(sp_mesh) { ), sp_material_(sp_material), sp_mesh_(sp_mesh) {
F_ = SpMat(sp_mesh_->GetNumDofs(), 1); F_ = SpMat(sp_mesh_->GetNumDofs(), 1);
K_ = SpMat(sp_mesh_->GetNumDofs(), sp_mesh_->GetNumDofs()); K_ = SpMat(sp_mesh_->GetNumDofs(), sp_mesh_->GetNumDofs());
spdlog::info("sizeof K: {}",K_.rows()); spdlog::info("DOF: {}", K_.rows());
spdlog::info("start to precompute..."); spdlog::info("start to precompute...");
Precompute(); Precompute();
} }
@ -51,7 +51,7 @@ namespace top {
if (directions[dir]) if (directions[dir])
dofs_to_fixed_.insert(node_id_to_fix[i] * 3 + dir); dofs_to_fixed_.insert(node_id_to_fix[i] * 3 + dir);
} }
auto ttt=dofs_to_fixed_; auto ttt = dofs_to_fixed_;
} }
void AddNBC(const Eigen::MatrixXi &NBC_coords, const Eigen::Vector3d &forces) { void AddNBC(const Eigen::MatrixXi &NBC_coords, const Eigen::Vector3d &forces) {
@ -59,95 +59,11 @@ namespace top {
for (int i = 0; i < node_id_to_load.size(); ++i) { for (int i = 0; i < node_id_to_load.size(); ++i) {
for (int dir = 0; dir < 3; ++dir) for (int dir = 0; dir < 3; ++dir)
if (forces[dir]) if (forces[dir])
F_.coeffRef(GetDof(node_id_to_load[i], dir),0) += forces(dir); F_.coeffRef(GetDof(node_id_to_load[i], dir), 0) += forces(dir);
}
}
void 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());
int loop = 0;
double change = 1.0;
double E0 = sp_material_->YM;
double Emin = sp_material_->YM * sp_para_->E_factor;
// precompute
Eigen::VectorXd dv(sp_mesh_->GetNumEles());
dv.setOnes();
dv = H_ * (dv.array() / Hs_.array()).matrix().eval();
// start iteration
while (change > sp_para_->tol_x && loop < sp_para_->max_loop) {
++loop;
Eigen::VectorXd sK = (sKe_ * (Emin + xPhys_col.array().pow(sp_para_->penal) *
(E0 - Emin)).matrix().transpose()).reshaped();
auto v_tri = Vec2Triplet(iK_, jK_, sK);
K_.setFromTriplets(v_tri.begin(), v_tri.end());
IntroduceFixedDofs(K_, F_);
#ifdef USE_SUITESPARSE
Eigen::CholmodSupernodalLLT<Eigen::SparseMatrix<double>> solver;
solver.compute(K_);
Eigen::VectorXd U = solver.solve(F_);
#else
// dense direct method
Eigen::MatrixXd denseK = Eigen::MatrixXd(K_);
Eigen::VectorXd U = denseK.fullPivLu().solve(Eigen::VectorXd(F_));
#endif
// compliance
Eigen::VectorXd ce(sp_mesh_->GetNumEles());
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);
ce(i) = Ue.transpose() * Ke_ * Ue;
}
double c = ce.transpose() * (Emin + xPhys_col.array().pow(sp_para_->penal) * (E0 - Emin)).matrix();
double v = xPhys_col.sum();
Eigen::VectorXd dc =
-sp_para_->penal * (E0 - Emin) * xPhys_col.array().pow(sp_para_->penal - 1.0) * ce.array();
// mma solver
size_t num_constrants = 1;
size_t num_variables = sp_mesh_->GetNumEles();
auto mma = std::make_shared<MMASolver>(num_variables, num_constrants);
Eigen::VectorXd variables_tmp = xPhys_col;
double f0val = c;
Eigen::VectorXd df0dx = dc;
double fval = v-sp_mesh_->GetNumEles()*sp_para_->volfrac;
Eigen::VectorXd dfdx = dv;
static Eigen::VectorXd low_bounds = Eigen::VectorXd::Zero(sp_mesh_->GetNumEles());
static Eigen::VectorXd up_bounds = Eigen::VectorXd::Ones(sp_mesh_->GetNumEles());
// spdlog::info("mma update");
mma->Update(variables_tmp.data(), df0dx.data(), &fval, dfdx.data(),
low_bounds.data(), up_bounds.data());
change = (variables_tmp - xPhys_col).cwiseAbs().maxCoeff();
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);
} }
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_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);
}
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_to_vtk(CMAKE_SOURCE_DIR "/output/vtk/iter_" + std::to_string(loop)+".vtk", ten_xPhys_to_write);
} void TopOptMainLoop();
private: private:
@ -158,59 +74,7 @@ namespace top {
} }
} }
void Precompute() { void Precompute();
Eigen::MatrixXi mat_ele2dofs = sp_mesh_->GetEleId2DofsMap();
iK_ = Eigen::KroneckerProduct(mat_ele2dofs, Eigen::VectorXi::Ones(24)).transpose().reshaped();
jK_ = Eigen::KroneckerProduct(mat_ele2dofs, Eigen::RowVectorXi::Ones(24)).transpose().reshaped();
sp_material_->computeKe(0.5, 0.5, 0.5, sp_material_->D, Ke_);
sKe_ = Ke_.reshaped();
// precompute filter
Eigen::VectorXi iH = Eigen::VectorXi::Ones(
sp_mesh_->GetNumEles() * std::pow(2.0 * (std::ceil(sp_para_->r_min) - 1.0) + 1, 3));
Eigen::VectorXi jH = iH;
Eigen::VectorXd sH(iH.size());
sH.setZero();
int cnt = 0;
int delta = std::ceil(sp_para_->r_min) - 1;
for (int k = 0; k < sp_mesh_->GetLz(); ++k) {
for (int j = 0; j < sp_mesh_->GetLy(); ++j) {
for (int i = 0; i < sp_mesh_->GetLx(); ++i) {
int ele_id0 = sp_mesh_->MapEleCoord2Id((Eigen::MatrixXi(1, 3) << i, j, k).finished())(0);
for (int k2 = std::max(k - delta, 0); k2 <= std::min(k + delta, sp_mesh_->GetLz() - 1); ++k2) {
for (int j2 = std::max(j - delta, 0);
j2 <= std::min(j + delta, sp_mesh_->GetLy() - 1); ++j2) {
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);
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());
++cnt;
}
}
}
}
}
}
std::vector<Eigen::Triplet<double>> v_tri;
for (int i = 0; i < cnt; ++i) {
v_tri.push_back({iH(i), jH(i), sH(i)});
}
H_ = SpMat(sp_mesh_->GetNumEles(), sp_mesh_->GetNumEles());
H_.setFromTriplets(v_tri.begin(), v_tri.end());
Hs_ = Eigen::VectorXd(sp_mesh_->GetNumEles());
for (int i = 0; i < Hs_.size(); ++i) {
Hs_(i) = H_.row(i).sum();
}
}
int GetDof(int node_id, int dir) { int GetDof(int node_id, int dir) {
return node_id * 3 + dir; return node_id * 3 + dir;

11
src/Util.h

@ -19,6 +19,7 @@
#include <Eigen/Dense> #include <Eigen/Dense>
#include <spdlog/spdlog.h> #include <spdlog/spdlog.h>
#include <iostream> #include <iostream>
#include <filesystem>
namespace top { namespace top {
using Tensor3d=TensorWrapper<double,3>; using Tensor3d=TensorWrapper<double,3>;
using Tensor3i=TensorWrapper<int,3>; using Tensor3i=TensorWrapper<int,3>;
@ -100,6 +101,16 @@ namespace top {
} }
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());
}
}
std::string the_path=path;
spdlog::info("clear dir : {}",the_path);
}
} // top } // top

Loading…
Cancel
Save