// // Created by cflin on 6/9/23. // #include #include #include "Boundary.h" #include "Mesh/HeatMesh.h" #include "Util.h" #include "ThermoelasticTop3d.h" #include "FEA/MechanicalLinearFEA.h" #include "FEA/ThermalLinearFEA.h" #include "TensorWrapper.h" int main() { using namespace da::sha; using top::fs_path; using std::string; top::fs_path output_dir(OUTPUT_DIR); top::fs_path config_file(CONFIG_FILE); top::fs_path assets_dir(ASSETS_DIR); spdlog::info("Algo read from '{}'", config_file.string()); spdlog::info("Algo output to '{}'", output_dir.string()); spdlog::info("assets dir: '{}'", assets_dir.string()); // read json std::ifstream f(config_file.c_str()); if (!f) { spdlog::critical("f open fail!"); exit(-7); } nlohmann::json j_config = nlohmann::json::parse(f); // read title std::string condition = j_config["TopologyOptimizationExample"]; spdlog::critical("TopologyOptimizationExample: {}", condition); // set topology parameters auto para = std::make_shared(); para->max_loop = j_config["topology"]["max_loop"]; para->volfrac = j_config["topology"]["volfrac"]; para->r_min = j_config["topology"]["r_min"]; para->penal = j_config["topology"]["penal"]; para->T_ref = j_config["topology"]["T_ref"]; para->T_limit = j_config["topology"]["T_limit"]; para->R_E = j_config["topology"]["R_E"]; para->R_lambda = j_config["topology"]["R_lambda"]; para->R_beta = j_config["topology"]["R_beta"]; // set material parameters double E = j_config["material"]["E"]; double Poisson_ratio = j_config["material"]["poisson_ratio"]; double thermal_conductivity = j_config["material"]["thermal_conductivity"]; double thermal_expansion_coefficient = j_config["material"]["thermal_expansion_coefficient"]; auto material = std::make_shared(E, Poisson_ratio, thermal_conductivity, thermal_expansion_coefficient); // set fea auto sp_mech_fea = std::make_shared(material);// for mechanical auto sp_thermal_fea = std::make_shared(material);// for thermal // set mesh(regular) int len_x = j_config["model"]["regular_model"]["lx"]; int len_y = j_config["model"]["regular_model"]["ly"]; int len_z = j_config["model"]["regular_model"]["lz"]; std::shared_ptr sp_mech_mesh; std::shared_ptr sp_thermal_mesh; sp_mech_mesh = std::make_shared(len_x, len_y, len_z); sp_thermal_mesh = std::make_shared(len_x, len_y, len_z); // initialize Top3d auto sp_mech_top3d = std::make_shared(para, sp_mech_fea, sp_mech_mesh); auto sp_thermal_top3d = std::make_shared(para, sp_thermal_fea, sp_thermal_mesh); // auxiliary class(Boundary) help to get boundary coordinates // see the comments in Boundary.h for more information top::Boundary mech_bdr(sp_mech_mesh); top::Boundary thermal_bdr(sp_thermal_mesh); { // add Dirichlet boundary, see the comments in Top3d::AddDBC for more information assert(j_config.count("mechanical_boundary_condition")); assert(j_config["mechanical_boundary_condition"].count("DBC")); assert(j_config["mechanical_boundary_condition"].count("NBC")); int DBCNum = j_config["mechanical_boundary_condition"]["DBC"].size(); for (int _i = 0; _i < DBCNum; ++_i) { const auto &DBCI = j_config["mechanical_boundary_condition"]["DBC"][_i]; Eigen::Vector3d minBBox(DBCI["min"][0], DBCI["min"][1], DBCI["min"][2]); Eigen::Vector3d maxBBox(DBCI["max"][0], DBCI["max"][1], DBCI["max"][2]); Eigen::Vector3i dir(DBCI["dir"][0], DBCI["dir"][1], DBCI["dir"][2]); top::Dir t_dir(minBBox, maxBBox, dir); sp_mech_top3d->AddDBC(mech_bdr.GetChosenCoordsByRelativeAlignedBox(t_dir.box), t_dir.dir); } // add Neumann boundary, see the comments in Top3d::AddNBC for more information int NBCNum = j_config["mechanical_boundary_condition"]["NBC"].size(); for (int _i = 0; _i < NBCNum; ++_i) { const auto &NBCI = j_config["mechanical_boundary_condition"]["NBC"][_i]; Eigen::Vector3d minBBox(NBCI["min"][0], NBCI["min"][1], NBCI["min"][2]); Eigen::Vector3d maxBBox(NBCI["max"][0], NBCI["max"][1], NBCI["max"][2]); Eigen::Vector3d val(NBCI["val"][0], NBCI["val"][1], NBCI["val"][2]); top::Neu t_neu(minBBox, maxBBox, val); Eigen::MatrixXi coords = mech_bdr.GetChosenCoordsByRelativeAlignedBox(t_neu.box); sp_mech_top3d->AddNBC(coords, t_neu.val / coords.rows()); } } { // add Dirichlet boundary, see the comments in Top3d::AddDBC for more information assert(j_config.count("thermal_boundary_condition")); assert(j_config["thermal_boundary_condition"].count("DBC")); assert(j_config["thermal_boundary_condition"].count("NBC")); int DBCNum = j_config["thermal_boundary_condition"]["DBC"].size(); for (int _i = 0; _i < DBCNum; ++_i) { const auto &DBCI = j_config["thermal_boundary_condition"]["DBC"][_i]; Eigen::Vector3d minBBox(DBCI["min"][0], DBCI["min"][1], DBCI["min"][2]); Eigen::Vector3d maxBBox(DBCI["max"][0], DBCI["max"][1], DBCI["max"][2]); top::Dir t_dir(minBBox, maxBBox, Eigen::Vector3i(1, 0, 0)); sp_thermal_top3d->AddDBC(mech_bdr.GetChosenCoordsByRelativeAlignedBox(t_dir.box), t_dir.dir, DBCI["temperature"]); } // add Neumann boundary, see the comments in Top3d::AddNBC for more information int NBCNum = j_config["thermal_boundary_condition"]["NBC"].size(); for (int _i = 0; _i < NBCNum; ++_i) { const auto &NBCI = j_config["thermal_boundary_condition"]["NBC"][_i]; Eigen::Vector3d minBBox(NBCI["min"][0], NBCI["min"][1], NBCI["min"][2]); Eigen::Vector3d maxBBox(NBCI["max"][0], NBCI["max"][1], NBCI["max"][2]); Eigen::Vector3d val(NBCI["heat_flux"], 0, 0); top::Neu t_neu(minBBox, maxBBox, val); Eigen::MatrixXi coords = mech_bdr.GetChosenCoordsByRelativeAlignedBox(t_neu.box); sp_thermal_top3d->AddNBC(coords, t_neu.val / coords.rows()); } } // init thermoelastic top3d top::ThermoelasticTop3d ther_top3d(sp_mech_top3d, sp_thermal_top3d); // process topology optimization top::Tensor3d ten_rho = ther_top3d.TopOptMainLoop(); // extract rho (txt or vtk) fs_path rho_txt_path = output_dir / "txt" / (condition + "_rho.txt"); write_tensor3d(rho_txt_path, ten_rho, sp_mech_mesh->GetOrigin(), sp_mech_mesh->GetOrigin() + sp_mech_mesh->GetLenBox()); spdlog::info("write density txt to: {}", rho_txt_path.c_str()); fs_path rho_vtk_path = output_dir / "vtk" / (condition + "_rho.vtk"); WriteTensorToVtk(rho_vtk_path, ten_rho, sp_mech_mesh); spdlog::info("write density vtk to: {}", rho_vtk_path.c_str()); // extract temperature(vtk) fs_path U_vtk_path = output_dir / "vtk" / (condition + "_U.vtk"); WriteUToVtk(U_vtk_path, ther_top3d.GetTemperature(), sp_mech_mesh); spdlog::info("write displacement norm vtk to: {}", U_vtk_path.c_str()); // extract stress field(txt or vtk) auto vonStress = ther_top3d.GetVonStress(); fs_path stressVon_txt_pah = output_dir / "txt" / (condition + "_stressVon.txt"); write_tensor3d(stressVon_txt_pah, vonStress, sp_mech_mesh->GetOrigin(), sp_mech_mesh->GetOrigin() + sp_mech_mesh->GetLenBox()); spdlog::info("write stressVon txt to: {}", stressVon_txt_pah.c_str()); fs_path stressVon_vtk_pah = output_dir / "vtk" / (condition + "_stressVon.vtk"); WriteTensorToVtk(stressVon_vtk_pah, vonStress, sp_mech_mesh); spdlog::info("write stressVon vtk to: {}", stressVon_vtk_pah.c_str()); auto v_tenosr = ther_top3d.GetTensorOfStress(Eigen::Vector3d{0, 1, 2}); fs_path stressX_txt_pah = output_dir / "txt" / (condition + "_stressX.txt"); write_tensor3d(stressX_txt_pah, v_tenosr[0], sp_mech_mesh->GetOrigin(), sp_mech_mesh->GetOrigin() + sp_mech_mesh->GetLenBox()); spdlog::info("write stressX txt to: {}", stressX_txt_pah.c_str()); fs_path stressX_vtk_path = output_dir / "vtk" / (condition + "_stressX.vtk"); WriteTensorToVtk(stressX_vtk_path, v_tenosr[0], sp_mech_mesh); spdlog::info("write stressX vtk to: {}", stressX_vtk_path.c_str()); fs_path stressY_txt_pah = output_dir / "txt" / (condition + "_stressY.txt"); write_tensor3d(stressY_txt_pah, v_tenosr[1], sp_mech_mesh->GetOrigin(), sp_mech_mesh->GetOrigin() + sp_mech_mesh->GetLenBox()); spdlog::info("write stressY txt to: {}", stressY_txt_pah.c_str()); fs_path stressY_vtk_path = output_dir / "vtk" / (condition + "_stressY.vtk"); WriteTensorToVtk(stressY_vtk_path, v_tenosr[1], sp_mech_mesh); spdlog::info("write stressY vtk to: {}", stressY_vtk_path.c_str()); fs_path stressZ_txt_pah = output_dir / "txt" / (condition + "_stressZ.txt"); write_tensor3d(stressZ_txt_pah, v_tenosr[2], sp_mech_mesh->GetOrigin(), sp_mech_mesh->GetOrigin() + sp_mech_mesh->GetLenBox()); spdlog::info("write stressZ txt to: {}", stressZ_txt_pah.c_str()); fs_path stressZ_vtk_path = output_dir / "vtk" / (condition + "_stressZ.vtk"); WriteTensorToVtk(stressZ_vtk_path, v_tenosr[2], sp_mech_mesh); spdlog::info("write stressZ vtk to: {}", stressZ_vtk_path.c_str()); }