1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466
|
// MFEM Example 37
//
// Compile with: make ex37
//
// Sample runs:
// ex37 -alpha 10
// ex37 -alpha 10 -pv
// ex37 -lambda 0.1 -mu 0.1
// ex37 -o 2 -alpha 5.0 -mi 50 -vf 0.4 -ntol 1e-5
// ex37 -r 6 -o 1 -alpha 25.0 -epsilon 0.02 -mi 50 -ntol 1e-5
//
// Description: This example code demonstrates the use of MFEM to solve a
// density-filtered [3] topology optimization problem. The
// objective is to minimize the compliance
//
// minimize ∫_Ω f⋅u dx over u ∈ [H¹(Ω)]² and ρ ∈ L¹(Ω)
//
// subject to
//
// -Div(r(ρ̃)Cε(u)) = f in Ω + BCs
// -ϵ²Δρ̃ + ρ̃ = ρ in Ω + Neumann BCs
// 0 ≤ ρ ≤ 1 in Ω
// ∫_Ω ρ dx = θ vol(Ω)
//
// Here, r(ρ̃) = ρ₀ + ρ̃³ (1-ρ₀) is the solid isotropic material
// penalization (SIMP) law, C is the elasticity tensor for an
// isotropic linearly elastic material, ϵ > 0 is the design
// length scale, and 0 < θ < 1 is the volume fraction.
//
// The problem is discretized and gradients are computing using
// finite elements [1]. The design is optimized using an entropic
// mirror descent algorithm introduced by Keith and Surowiec [2]
// that is tailored to the bound constraint 0 ≤ ρ ≤ 1.
//
// This example highlights the ability of MFEM to deliver high-
// order solutions to inverse design problems and showcases how
// to set up and solve PDE-constrained optimization problems
// using the so-called reduced space approach.
//
// [1] Andreassen, E., Clausen, A., Schevenels, M., Lazarov, B. S., & Sigmund, O.
// (2011). Efficient topology optimization in MATLAB using 88 lines of
// code. Structural and Multidisciplinary Optimization, 43(1), 1-16.
// [2] Keith, B. and Surowiec, T. (2023) Proximal Galerkin: A structure-
// preserving finite element method for pointwise bound constraints.
// arXiv:2307.12444 [math.NA]
// [3] Lazarov, B. S., & Sigmund, O. (2011). Filters in topology optimization
// based on Helmholtz‐type differential equations. International Journal
// for Numerical Methods in Engineering, 86(6), 765-781.
#include "mfem.hpp"
#include <iostream>
#include <fstream>
#include "ex37.hpp"
using namespace std;
using namespace mfem;
/**
* @brief Bregman projection of ρ = sigmoid(ψ) onto the subspace
* ∫_Ω ρ dx = θ vol(Ω) as follows:
*
* 1. Compute the root of the R → R function
* f(c) = ∫_Ω sigmoid(ψ + c) dx - θ vol(Ω)
* 2. Set ψ ← ψ + c.
*
* @param psi a GridFunction to be updated
* @param target_volume θ vol(Ω)
* @param tol Newton iteration tolerance
* @param max_its Newton maximum iteration number
* @return real_t Final volume, ∫_Ω sigmoid(ψ)
*/
real_t proj(GridFunction &psi, real_t target_volume, real_t tol=1e-12,
int max_its=10)
{
MappedGridFunctionCoefficient sigmoid_psi(&psi, sigmoid);
MappedGridFunctionCoefficient der_sigmoid_psi(&psi, der_sigmoid);
LinearForm int_sigmoid_psi(psi.FESpace());
int_sigmoid_psi.AddDomainIntegrator(new DomainLFIntegrator(sigmoid_psi));
LinearForm int_der_sigmoid_psi(psi.FESpace());
int_der_sigmoid_psi.AddDomainIntegrator(new DomainLFIntegrator(
der_sigmoid_psi));
bool done = false;
for (int k=0; k<max_its; k++) // Newton iteration
{
int_sigmoid_psi.Assemble(); // Recompute f(c) with updated ψ
const real_t f = int_sigmoid_psi.Sum() - target_volume;
int_der_sigmoid_psi.Assemble(); // Recompute df(c) with updated ψ
const real_t df = int_der_sigmoid_psi.Sum();
const real_t dc = -f/df;
psi += dc;
if (abs(dc) < tol) { done = true; break; }
}
if (!done)
{
mfem_warning("Projection reached maximum iteration without converging. "
"Result may not be accurate.");
}
int_sigmoid_psi.Assemble();
return int_sigmoid_psi.Sum();
}
/*
* ---------------------------------------------------------------
* ALGORITHM PREAMBLE
* ---------------------------------------------------------------
*
* The Lagrangian for this problem is
*
* L(u,ρ,ρ̃,w,w̃) = (f,u) - (r(ρ̃) C ε(u),ε(w)) + (f,w)
* - (ϵ² ∇ρ̃,∇w̃) - (ρ̃,w̃) + (ρ,w̃)
*
* where
*
* r(ρ̃) = ρ₀ + ρ̃³ (1 - ρ₀) (SIMP rule)
*
* ε(u) = (∇u + ∇uᵀ)/2 (symmetric gradient)
*
* C e = λtr(e)I + 2μe (isotropic material)
*
* NOTE: The Lame parameters can be computed from Young's modulus E
* and Poisson's ratio ν as follows:
*
* λ = E ν/((1+ν)(1-2ν)), μ = E/(2(1+ν))
*
* ---------------------------------------------------------------
*
* Discretization choices:
*
* u ∈ V ⊂ (H¹)ᵈ (order p)
* ψ ∈ L² (order p - 1), ρ = sigmoid(ψ)
* ρ̃ ∈ H¹ (order p)
* w ∈ V (order p)
* w̃ ∈ H¹ (order p)
*
* ---------------------------------------------------------------
* ALGORITHM
* ---------------------------------------------------------------
*
* Update ρ with projected mirror descent via the following algorithm.
*
* 1. Initialize ψ = inv_sigmoid(vol_fraction) so that ∫ sigmoid(ψ) = θ vol(Ω)
*
* While not converged:
*
* 2. Solve filter equation ∂_w̃ L = 0; i.e.,
*
* (ϵ² ∇ ρ̃, ∇ v ) + (ρ̃,v) = (ρ,v) ∀ v ∈ H¹.
*
* 3. Solve primal problem ∂_w L = 0; i.e.,
*
* (λ r(ρ̃) ∇⋅u, ∇⋅v) + (2 μ r(ρ̃) ε(u), ε(v)) = (f,v) ∀ v ∈ V.
*
* NB. The dual problem ∂_u L = 0 is the negative of the primal problem due to symmetry.
*
* 4. Solve for filtered gradient ∂_ρ̃ L = 0; i.e.,
*
* (ϵ² ∇ w̃ , ∇ v ) + (w̃ ,v) = (-r'(ρ̃) ( λ |∇⋅u|² + 2 μ |ε(u)|²),v) ∀ v ∈ H¹.
*
* 5. Project the gradient onto the discrete latent space; i.e., solve
*
* (G,v) = (w̃,v) ∀ v ∈ L².
*
* 6. Bregman proximal gradient update; i.e.,
*
* ψ ← ψ - αG + c,
*
* where α > 0 is a step size parameter and c ∈ R is a constant ensuring
*
* ∫_Ω sigmoid(ψ - αG + c) dx = θ vol(Ω).
*
* end
*/
int main(int argc, char *argv[])
{
// 1. Parse command-line options.
int ref_levels = 5;
int order = 2;
real_t alpha = 1.0;
real_t epsilon = 0.01;
real_t vol_fraction = 0.5;
int max_it = 1e3;
real_t itol = 1e-1;
real_t ntol = 1e-4;
real_t rho_min = 1e-6;
real_t lambda = 1.0;
real_t mu = 1.0;
bool glvis_visualization = true;
bool paraview_output = false;
OptionsParser args(argc, argv);
args.AddOption(&ref_levels, "-r", "--refine",
"Number of times to refine the mesh uniformly.");
args.AddOption(&order, "-o", "--order",
"Order (degree) of the finite elements.");
args.AddOption(&alpha, "-alpha", "--alpha-step-length",
"Step length for gradient descent.");
args.AddOption(&epsilon, "-epsilon", "--epsilon-thickness",
"Length scale for ρ.");
args.AddOption(&max_it, "-mi", "--max-it",
"Maximum number of gradient descent iterations.");
args.AddOption(&ntol, "-ntol", "--rel-tol",
"Normalized exit tolerance.");
args.AddOption(&itol, "-itol", "--abs-tol",
"Increment exit tolerance.");
args.AddOption(&vol_fraction, "-vf", "--volume-fraction",
"Volume fraction for the material density.");
args.AddOption(&lambda, "-lambda", "--lambda",
"Lamé constant λ.");
args.AddOption(&mu, "-mu", "--mu",
"Lamé constant μ.");
args.AddOption(&rho_min, "-rmin", "--psi-min",
"Minimum of density coefficient.");
args.AddOption(&glvis_visualization, "-vis", "--visualization", "-no-vis",
"--no-visualization",
"Enable or disable GLVis visualization.");
args.AddOption(¶view_output, "-pv", "--paraview", "-no-pv",
"--no-paraview",
"Enable or disable ParaView output.");
args.Parse();
if (!args.Good())
{
args.PrintUsage(mfem::out);
return 1;
}
args.PrintOptions(mfem::out);
Mesh mesh = Mesh::MakeCartesian2D(3, 1, mfem::Element::Type::QUADRILATERAL,
true, 3.0, 1.0);
int dim = mesh.Dimension();
// 2. Set BCs.
for (int i = 0; i<mesh.GetNBE(); i++)
{
Element * be = mesh.GetBdrElement(i);
Array<int> vertices;
be->GetVertices(vertices);
real_t * coords1 = mesh.GetVertex(vertices[0]);
real_t * coords2 = mesh.GetVertex(vertices[1]);
Vector center(2);
center(0) = 0.5*(coords1[0] + coords2[0]);
center(1) = 0.5*(coords1[1] + coords2[1]);
if (abs(center(0) - 0.0) < 1e-10)
{
// the left edge
be->SetAttribute(1);
}
else
{
// all other boundaries
be->SetAttribute(2);
}
}
mesh.SetAttributes();
// 3. Refine the mesh.
for (int lev = 0; lev < ref_levels; lev++)
{
mesh.UniformRefinement();
}
// 4. Define the necessary finite element spaces on the mesh.
H1_FECollection state_fec(order, dim); // space for u
H1_FECollection filter_fec(order, dim); // space for ρ̃
L2_FECollection control_fec(order-1, dim,
BasisType::GaussLobatto); // space for ψ
FiniteElementSpace state_fes(&mesh, &state_fec,dim);
FiniteElementSpace filter_fes(&mesh, &filter_fec);
FiniteElementSpace control_fes(&mesh, &control_fec);
int state_size = state_fes.GetTrueVSize();
int control_size = control_fes.GetTrueVSize();
int filter_size = filter_fes.GetTrueVSize();
mfem::out << "Number of state unknowns: " << state_size << std::endl;
mfem::out << "Number of filter unknowns: " << filter_size << std::endl;
mfem::out << "Number of control unknowns: " << control_size << std::endl;
// 5. Set the initial guess for ρ.
GridFunction u(&state_fes);
GridFunction psi(&control_fes);
GridFunction psi_old(&control_fes);
GridFunction rho_filter(&filter_fes);
u = 0.0;
rho_filter = vol_fraction;
psi = inv_sigmoid(vol_fraction);
psi_old = inv_sigmoid(vol_fraction);
// ρ = sigmoid(ψ)
MappedGridFunctionCoefficient rho(&psi, sigmoid);
// Interpolation of ρ = sigmoid(ψ) in control fes (for ParaView output)
GridFunction rho_gf(&control_fes);
// ρ - ρ_old = sigmoid(ψ) - sigmoid(ψ_old)
DiffMappedGridFunctionCoefficient succ_diff_rho(&psi, &psi_old, sigmoid);
// 6. Set-up the physics solver.
int maxat = mesh.bdr_attributes.Max();
Array<int> ess_bdr(maxat);
ess_bdr = 0;
ess_bdr[0] = 1;
ConstantCoefficient one(1.0);
ConstantCoefficient lambda_cf(lambda);
ConstantCoefficient mu_cf(mu);
LinearElasticitySolver * ElasticitySolver = new LinearElasticitySolver();
ElasticitySolver->SetMesh(&mesh);
ElasticitySolver->SetOrder(state_fec.GetOrder());
ElasticitySolver->SetupFEM();
Vector center(2); center(0) = 2.9; center(1) = 0.5;
Vector force(2); force(0) = 0.0; force(1) = -1.0;
real_t r = 0.05;
VolumeForceCoefficient vforce_cf(r,center,force);
ElasticitySolver->SetRHSCoefficient(&vforce_cf);
ElasticitySolver->SetEssentialBoundary(ess_bdr);
// 7. Set-up the filter solver.
ConstantCoefficient eps2_cf(epsilon*epsilon);
DiffusionSolver * FilterSolver = new DiffusionSolver();
FilterSolver->SetMesh(&mesh);
FilterSolver->SetOrder(filter_fec.GetOrder());
FilterSolver->SetDiffusionCoefficient(&eps2_cf);
FilterSolver->SetMassCoefficient(&one);
Array<int> ess_bdr_filter;
if (mesh.bdr_attributes.Size())
{
ess_bdr_filter.SetSize(mesh.bdr_attributes.Max());
ess_bdr_filter = 0;
}
FilterSolver->SetEssentialBoundary(ess_bdr_filter);
FilterSolver->SetupFEM();
BilinearForm mass(&control_fes);
mass.AddDomainIntegrator(new InverseIntegrator(new MassIntegrator(one)));
mass.Assemble();
SparseMatrix M;
Array<int> empty;
mass.FormSystemMatrix(empty,M);
// 8. Define the Lagrange multiplier and gradient functions.
GridFunction grad(&control_fes);
GridFunction w_filter(&filter_fes);
// 9. Define some tools for later.
ConstantCoefficient zero(0.0);
GridFunction onegf(&control_fes);
onegf = 1.0;
GridFunction zerogf(&control_fes);
zerogf = 0.0;
LinearForm vol_form(&control_fes);
vol_form.AddDomainIntegrator(new DomainLFIntegrator(one));
vol_form.Assemble();
real_t domain_volume = vol_form(onegf);
const real_t target_volume = domain_volume * vol_fraction;
// 10. Connect to GLVis. Prepare for VisIt output.
char vishost[] = "localhost";
int visport = 19916;
socketstream sout_r;
if (glvis_visualization)
{
sout_r.open(vishost, visport);
sout_r.precision(8);
}
mfem::ParaViewDataCollection paraview_dc("ex37", &mesh);
if (paraview_output)
{
rho_gf.ProjectCoefficient(rho);
paraview_dc.SetPrefixPath("ParaView");
paraview_dc.SetLevelsOfDetail(order);
paraview_dc.SetDataFormat(VTKFormat::BINARY);
paraview_dc.SetHighOrderOutput(true);
paraview_dc.SetCycle(0);
paraview_dc.SetTime(0.0);
paraview_dc.RegisterField("displacement",&u);
paraview_dc.RegisterField("density",&rho_gf);
paraview_dc.RegisterField("filtered_density",&rho_filter);
paraview_dc.Save();
}
// 11. Iterate:
for (int k = 1; k <= max_it; k++)
{
if (k > 1) { alpha *= ((real_t) k) / ((real_t) k-1); }
mfem::out << "\nStep = " << k << std::endl;
// Step 1 - Filter solve
// Solve (ϵ^2 ∇ ρ̃, ∇ v ) + (ρ̃,v) = (ρ,v)
FilterSolver->SetRHSCoefficient(&rho);
FilterSolver->Solve();
rho_filter = *FilterSolver->GetFEMSolution();
// Step 2 - State solve
// Solve (λ r(ρ̃) ∇⋅u, ∇⋅v) + (2 μ r(ρ̃) ε(u), ε(v)) = (f,v)
SIMPInterpolationCoefficient SIMP_cf(&rho_filter,rho_min, 1.0);
ProductCoefficient lambda_SIMP_cf(lambda_cf,SIMP_cf);
ProductCoefficient mu_SIMP_cf(mu_cf,SIMP_cf);
ElasticitySolver->SetLameCoefficients(&lambda_SIMP_cf,&mu_SIMP_cf);
ElasticitySolver->Solve();
u = *ElasticitySolver->GetFEMSolution();
// Step 3 - Adjoint filter solve
// Solve (ϵ² ∇ w̃, ∇ v) + (w̃ ,v) = (-r'(ρ̃) ( λ |∇⋅u|² + 2 μ |ε(u)|²),v)
StrainEnergyDensityCoefficient rhs_cf(&lambda_cf,&mu_cf,&u, &rho_filter,
rho_min);
FilterSolver->SetRHSCoefficient(&rhs_cf);
FilterSolver->Solve();
w_filter = *FilterSolver->GetFEMSolution();
// Step 4 - Compute gradient
// Solve G = M⁻¹w̃
GridFunctionCoefficient w_cf(&w_filter);
LinearForm w_rhs(&control_fes);
w_rhs.AddDomainIntegrator(new DomainLFIntegrator(w_cf));
w_rhs.Assemble();
M.Mult(w_rhs,grad);
// Step 5 - Update design variable ψ ← proj(ψ - αG)
psi.Add(-alpha, grad);
const real_t material_volume = proj(psi, target_volume);
// Compute ||ρ - ρ_old|| in control fes.
real_t norm_increment = zerogf.ComputeL1Error(succ_diff_rho);
real_t norm_reduced_gradient = norm_increment/alpha;
psi_old = psi;
real_t compliance = (*(ElasticitySolver->GetLinearForm()))(u);
mfem::out << "norm of the reduced gradient = " << norm_reduced_gradient <<
std::endl;
mfem::out << "norm of the increment = " << norm_increment << endl;
mfem::out << "compliance = " << compliance << std::endl;
mfem::out << "volume fraction = " << material_volume / domain_volume <<
std::endl;
if (glvis_visualization)
{
GridFunction r_gf(&filter_fes);
r_gf.ProjectCoefficient(SIMP_cf);
sout_r << "solution\n" << mesh << r_gf
<< "window_title 'Design density r(ρ̃)'" << flush;
}
if (paraview_output)
{
rho_gf.ProjectCoefficient(rho);
paraview_dc.SetCycle(k);
paraview_dc.SetTime((real_t)k);
paraview_dc.Save();
}
if (norm_reduced_gradient < ntol && norm_increment < itol)
{
break;
}
}
delete ElasticitySolver;
delete FilterSolver;
return 0;
}
|