MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex37p.cpp
Go to the documentation of this file.
1// MFEM Example 37 - Parallel Version
2//
3// Compile with: make ex37p
4//
5// Sample runs:
6// mpirun -np 4 ex37p -alpha 10 -pv
7// mpirun -np 4 ex37p -lambda 0.1 -mu 0.1
8// mpirun -np 4 ex37p -o 2 -alpha 5.0 -mi 50 -vf 0.4 -ntol 1e-5
9// mpirun -np 4 ex37p -r 6 -o 2 -alpha 10.0 -epsilon 0.02 -mi 50 -ntol 1e-5
10//
11// Description: This example code demonstrates the use of MFEM to solve a
12// density-filtered [3] topology optimization problem. The
13// objective is to minimize the compliance
14//
15// minimize ∫_Ω f⋅u dx over u ∈ [H¹(Ω)]² and ρ ∈ L¹(Ω)
16//
17// subject to
18//
19// -Div(r(ρ̃)Cε(u)) = f in Ω + BCs
20// -ϵ²Δρ̃ + ρ̃ = ρ in Ω + Neumann BCs
21// 0 ≤ ρ ≤ 1 in Ω
22// ∫_Ω ρ dx = θ vol(Ω)
23//
24// Here, r(ρ̃) = ρ₀ + ρ̃³ (1-ρ₀) is the solid isotropic material
25// penalization (SIMP) law, C is the elasticity tensor for an
26// isotropic linearly elastic material, ϵ > 0 is the design
27// length scale, and 0 < θ < 1 is the volume fraction.
28//
29// The problem is discretized and gradients are computing using
30// finite elements [1]. The design is optimized using an entropic
31// mirror descent algorithm introduced by Keith and Surowiec [2]
32// that is tailored to the bound constraint 0 ≤ ρ ≤ 1.
33//
34// This example highlights the ability of MFEM to deliver high-
35// order solutions to inverse design problems and showcases how
36// to set up and solve PDE-constrained optimization problems
37// using the so-called reduced space approach.
38//
39// [1] Andreassen, E., Clausen, A., Schevenels, M., Lazarov, B. S., & Sigmund, O.
40// (2011). Efficient topology optimization in MATLAB using 88 lines of
41// code. Structural and Multidisciplinary Optimization, 43(1), 1-16.
42// [2] Keith, B. and Surowiec, T. (2023) Proximal Galerkin: A structure-
43// preserving finite element method for pointwise bound constraints.
44// arXiv:2307.12444 [math.NA]
45// [3] Lazarov, B. S., & Sigmund, O. (2011). Filters in topology optimization
46// based on Helmholtz‐type differential equations. International Journal
47// for Numerical Methods in Engineering, 86(6), 765-781.
48
49#include "mfem.hpp"
50#include <iostream>
51#include <fstream>
52#include "ex37.hpp"
53
54using namespace std;
55using namespace mfem;
56
57/**
58 * @brief Bregman projection of ρ = sigmoid(ψ) onto the subspace
59 * ∫_Ω ρ dx = θ vol(Ω) as follows:
60 *
61 * 1. Compute the root of the R → R function
62 * f(c) = ∫_Ω sigmoid(ψ + c) dx - θ vol(Ω)
63 * 2. Set ψ ← ψ + c.
64 *
65 * @param psi a GridFunction to be updated
66 * @param target_volume θ vol(Ω)
67 * @param tol Newton iteration tolerance
68 * @param max_its Newton maximum iteration number
69 * @return real_t Final volume, ∫_Ω sigmoid(ψ)
70 */
71real_t proj(ParGridFunction &psi, real_t target_volume, real_t tol=1e-12,
72 int max_its=10)
73{
74 MappedGridFunctionCoefficient sigmoid_psi(&psi, sigmoid);
75 MappedGridFunctionCoefficient der_sigmoid_psi(&psi, der_sigmoid);
76
77 ParLinearForm int_sigmoid_psi(psi.ParFESpace());
78 int_sigmoid_psi.AddDomainIntegrator(new DomainLFIntegrator(sigmoid_psi));
79 ParLinearForm int_der_sigmoid_psi(psi.ParFESpace());
80 int_der_sigmoid_psi.AddDomainIntegrator(new DomainLFIntegrator(
81 der_sigmoid_psi));
82 bool done = false;
83 for (int k=0; k<max_its; k++) // Newton iteration
84 {
85 int_sigmoid_psi.Assemble(); // Recompute f(c) with updated ψ
86 real_t f = int_sigmoid_psi.Sum();
87 MPI_Allreduce(MPI_IN_PLACE, &f, 1, MPITypeMap<real_t>::mpi_type,
88 MPI_SUM, MPI_COMM_WORLD);
89 f -= target_volume;
90
91 int_der_sigmoid_psi.Assemble(); // Recompute df(c) with updated ψ
92 real_t df = int_der_sigmoid_psi.Sum();
93 MPI_Allreduce(MPI_IN_PLACE, &df, 1, MPITypeMap<real_t>::mpi_type,
94 MPI_SUM, MPI_COMM_WORLD);
95
96 const real_t dc = -f/df;
97 psi += dc;
98 if (abs(dc) < tol) { done = true; break; }
99 }
100 if (!done)
101 {
102 mfem_warning("Projection reached maximum iteration without converging. "
103 "Result may not be accurate.");
104 }
105 int_sigmoid_psi.Assemble();
106 real_t material_volume = int_sigmoid_psi.Sum();
107 MPI_Allreduce(MPI_IN_PLACE, &material_volume, 1,
108 MPITypeMap<real_t>::mpi_type, MPI_SUM, MPI_COMM_WORLD);
109 return material_volume;
110}
111
112/*
113 * ---------------------------------------------------------------
114 * ALGORITHM PREAMBLE
115 * ---------------------------------------------------------------
116 *
117 * The Lagrangian for this problem is
118 *
119 * L(u,ρ,ρ̃,w,w̃) = (f,u) - (r(ρ̃) C ε(u),ε(w)) + (f,w)
120 * - (ϵ² ∇ρ̃,∇w̃) - (ρ̃,w̃) + (ρ,w̃)
121 *
122 * where
123 *
124 * r(ρ̃) = ρ₀ + ρ̃³ (1 - ρ₀) (SIMP rule)
125 *
126 * ε(u) = (∇u + ∇uᵀ)/2 (symmetric gradient)
127 *
128 * C e = λtr(e)I + 2μe (isotropic material)
129 *
130 * NOTE: The Lame parameters can be computed from Young's modulus E
131 * and Poisson's ratio ν as follows:
132 *
133 * λ = E ν/((1+ν)(1-2ν)), μ = E/(2(1+ν))
134 *
135 * ---------------------------------------------------------------
136 *
137 * Discretization choices:
138 *
139 * u ∈ V ⊂ (H¹)ᵈ (order p)
140 * ψ ∈ L² (order p - 1), ρ = sigmoid(ψ)
141 * ρ̃ ∈ H¹ (order p)
142 * w ∈ V (order p)
143 * w̃ ∈ H¹ (order p)
144 *
145 * ---------------------------------------------------------------
146 * ALGORITHM
147 * ---------------------------------------------------------------
148 *
149 * Update ρ with projected mirror descent via the following algorithm.
150 *
151 * 1. Initialize ψ = inv_sigmoid(vol_fraction) so that ∫ sigmoid(ψ) = θ vol(Ω)
152 *
153 * While not converged:
154 *
155 * 2. Solve filter equation ∂_w̃ L = 0; i.e.,
156 *
157 * (ϵ² ∇ ρ̃, ∇ v ) + (ρ̃,v) = (ρ,v) ∀ v ∈ H¹.
158 *
159 * 3. Solve primal problem ∂_w L = 0; i.e.,
160 *
161 * (λ r(ρ̃) ∇⋅u, ∇⋅v) + (2 μ r(ρ̃) ε(u), ε(v)) = (f,v) ∀ v ∈ V.
162 *
163 * NB. The dual problem ∂_u L = 0 is the negative of the primal problem due to symmetry.
164 *
165 * 4. Solve for filtered gradient ∂_ρ̃ L = 0; i.e.,
166 *
167 * (ϵ² ∇ w̃ , ∇ v ) + (w̃ ,v) = (-r'(ρ̃) ( λ |∇⋅u|² + 2 μ |ε(u)|²),v) ∀ v ∈ H¹.
168 *
169 * 5. Project the gradient onto the discrete latent space; i.e., solve
170 *
171 * (G,v) = (w̃,v) ∀ v ∈ L².
172 *
173 * 6. Bregman proximal gradient update; i.e.,
174 *
175 * ψ ← ψ - αG + c,
176 *
177 * where α > 0 is a step size parameter and c ∈ R is a constant ensuring
178 *
179 * ∫_Ω sigmoid(ψ - αG + c) dx = θ vol(Ω).
180 *
181 * end
182 */
183
184int main(int argc, char *argv[])
185{
186 // 0. Initialize MPI and HYPRE.
187 Mpi::Init();
188 int num_procs = Mpi::WorldSize();
189 int myid = Mpi::WorldRank();
190 Hypre::Init();
191
192 // 1. Parse command-line options.
193 int ref_levels = 5;
194 int order = 2;
195 real_t alpha = 1.0;
196 real_t epsilon = 0.01;
197 real_t vol_fraction = 0.5;
198 int max_it = 1e3;
199 real_t itol = 1e-1;
200 real_t ntol = 1e-4;
201 real_t rho_min = 1e-6;
202 real_t lambda = 1.0;
203 real_t mu = 1.0;
204 bool glvis_visualization = true;
205 bool paraview_output = false;
206
207 OptionsParser args(argc, argv);
208 args.AddOption(&ref_levels, "-r", "--refine",
209 "Number of times to refine the mesh uniformly.");
210 args.AddOption(&order, "-o", "--order",
211 "Order (degree) of the finite elements.");
212 args.AddOption(&alpha, "-alpha", "--alpha-step-length",
213 "Step length for gradient descent.");
214 args.AddOption(&epsilon, "-epsilon", "--epsilon-thickness",
215 "Length scale for ρ.");
216 args.AddOption(&max_it, "-mi", "--max-it",
217 "Maximum number of gradient descent iterations.");
218 args.AddOption(&ntol, "-ntol", "--rel-tol",
219 "Normalized exit tolerance.");
220 args.AddOption(&itol, "-itol", "--abs-tol",
221 "Increment exit tolerance.");
222 args.AddOption(&vol_fraction, "-vf", "--volume-fraction",
223 "Volume fraction for the material density.");
224 args.AddOption(&lambda, "-lambda", "--lambda",
225 "Lamé constant λ.");
226 args.AddOption(&mu, "-mu", "--mu",
227 "Lamé constant μ.");
228 args.AddOption(&rho_min, "-rmin", "--psi-min",
229 "Minimum of density coefficient.");
230 args.AddOption(&glvis_visualization, "-vis", "--visualization", "-no-vis",
231 "--no-visualization",
232 "Enable or disable GLVis visualization.");
233 args.AddOption(&paraview_output, "-pv", "--paraview", "-no-pv",
234 "--no-paraview",
235 "Enable or disable ParaView output.");
236 args.Parse();
237 if (!args.Good())
238 {
239 if (myid == 0)
240 {
241 args.PrintUsage(cout);
242 }
243 MPI_Finalize();
244 return 1;
245 }
246 if (myid == 0)
247 {
248 mfem::out << num_procs << " number of process created.\n";
249 args.PrintOptions(cout);
250 }
251
253 true, 3.0, 1.0);
254 int dim = mesh.Dimension();
255
256 // 2. Set BCs.
257 for (int i = 0; i<mesh.GetNBE(); i++)
258 {
259 Element * be = mesh.GetBdrElement(i);
260 Array<int> vertices;
261 be->GetVertices(vertices);
262
263 real_t * coords1 = mesh.GetVertex(vertices[0]);
264 real_t * coords2 = mesh.GetVertex(vertices[1]);
265
266 Vector center(2);
267 center(0) = 0.5*(coords1[0] + coords2[0]);
268 center(1) = 0.5*(coords1[1] + coords2[1]);
269
270 if (abs(center(0) - 0.0) < 1e-10)
271 {
272 // the left edge
273 be->SetAttribute(1);
274 }
275 else
276 {
277 // all other boundaries
278 be->SetAttribute(2);
279 }
280 }
281 mesh.SetAttributes();
282
283 // 3. Refine the mesh.
284 for (int lev = 0; lev < ref_levels; lev++)
285 {
286 mesh.UniformRefinement();
287 }
288
289 ParMesh pmesh(MPI_COMM_WORLD, mesh);
290 mesh.Clear();
291
292 // 4. Define the necessary finite element spaces on the mesh.
293 H1_FECollection state_fec(order, dim); // space for u
294 H1_FECollection filter_fec(order, dim); // space for ρ̃
295 L2_FECollection control_fec(order-1, dim,
296 BasisType::GaussLobatto); // space for ψ
297 ParFiniteElementSpace state_fes(&pmesh, &state_fec,dim);
298 ParFiniteElementSpace filter_fes(&pmesh, &filter_fec);
299 ParFiniteElementSpace control_fes(&pmesh, &control_fec);
300
301 HYPRE_BigInt state_size = state_fes.GlobalTrueVSize();
302 HYPRE_BigInt control_size = control_fes.GlobalTrueVSize();
303 HYPRE_BigInt filter_size = filter_fes.GlobalTrueVSize();
304 if (myid==0)
305 {
306 cout << "Number of state unknowns: " << state_size << endl;
307 cout << "Number of filter unknowns: " << filter_size << endl;
308 cout << "Number of control unknowns: " << control_size << endl;
309 }
310
311 // 5. Set the initial guess for ρ.
312 ParGridFunction u(&state_fes);
313 ParGridFunction psi(&control_fes);
314 ParGridFunction psi_old(&control_fes);
315 ParGridFunction rho_filter(&filter_fes);
316 u = 0.0;
317 rho_filter = vol_fraction;
318 psi = inv_sigmoid(vol_fraction);
319 psi_old = inv_sigmoid(vol_fraction);
320
321 // ρ = sigmoid(ψ)
323 // Interpolation of ρ = sigmoid(ψ) in control fes (for ParaView output)
324 ParGridFunction rho_gf(&control_fes);
325 // ρ - ρ_old = sigmoid(ψ) - sigmoid(ψ_old)
326 DiffMappedGridFunctionCoefficient succ_diff_rho(&psi, &psi_old, sigmoid);
327
328 // 6. Set-up the physics solver.
329 int maxat = pmesh.bdr_attributes.Max();
330 Array<int> ess_bdr(maxat);
331 ess_bdr = 0;
332 ess_bdr[0] = 1;
333 ConstantCoefficient one(1.0);
334 ConstantCoefficient lambda_cf(lambda);
335 ConstantCoefficient mu_cf(mu);
336 LinearElasticitySolver * ElasticitySolver = new LinearElasticitySolver();
337 ElasticitySolver->SetMesh(&pmesh);
338 ElasticitySolver->SetOrder(state_fec.GetOrder());
339 ElasticitySolver->SetupFEM();
340 Vector center(2); center(0) = 2.9; center(1) = 0.5;
341 Vector force(2); force(0) = 0.0; force(1) = -1.0;
342 real_t r = 0.05;
343 VolumeForceCoefficient vforce_cf(r,center,force);
344 ElasticitySolver->SetRHSCoefficient(&vforce_cf);
345 ElasticitySolver->SetEssentialBoundary(ess_bdr);
346
347 // 7. Set-up the filter solver.
349 DiffusionSolver * FilterSolver = new DiffusionSolver();
350 FilterSolver->SetMesh(&pmesh);
351 FilterSolver->SetOrder(filter_fec.GetOrder());
352 FilterSolver->SetDiffusionCoefficient(&eps2_cf);
353 FilterSolver->SetMassCoefficient(&one);
354 Array<int> ess_bdr_filter;
355 if (pmesh.bdr_attributes.Size())
356 {
357 ess_bdr_filter.SetSize(pmesh.bdr_attributes.Max());
358 ess_bdr_filter = 0;
359 }
360 FilterSolver->SetEssentialBoundary(ess_bdr_filter);
361 FilterSolver->SetupFEM();
362
363 ParBilinearForm mass(&control_fes);
365 mass.Assemble();
367 Array<int> empty;
368 mass.FormSystemMatrix(empty,M);
369
370 // 8. Define the Lagrange multiplier and gradient functions.
371 ParGridFunction grad(&control_fes);
372 ParGridFunction w_filter(&filter_fes);
373
374 // 9. Define some tools for later.
375 ConstantCoefficient zero(0.0);
376 ParGridFunction onegf(&control_fes);
377 onegf = 1.0;
378 ParGridFunction zerogf(&control_fes);
379 zerogf = 0.0;
380 ParLinearForm vol_form(&control_fes);
381 vol_form.AddDomainIntegrator(new DomainLFIntegrator(one));
382 vol_form.Assemble();
383 real_t domain_volume = vol_form(onegf);
384 const real_t target_volume = domain_volume * vol_fraction;
385
386 // 10. Connect to GLVis. Prepare for VisIt output.
387 char vishost[] = "localhost";
388 int visport = 19916;
389 socketstream sout_r;
390 if (glvis_visualization)
391 {
392 sout_r.open(vishost, visport);
393 sout_r.precision(8);
394 }
395
396 mfem::ParaViewDataCollection paraview_dc("ex37p", &pmesh);
397 if (paraview_output)
398 {
399 rho_gf.ProjectCoefficient(rho);
400 paraview_dc.SetPrefixPath("ParaView");
401 paraview_dc.SetLevelsOfDetail(order);
402 paraview_dc.SetDataFormat(VTKFormat::BINARY);
403 paraview_dc.SetHighOrderOutput(true);
404 paraview_dc.SetCycle(0);
405 paraview_dc.SetTime(0.0);
406 paraview_dc.RegisterField("displacement",&u);
407 paraview_dc.RegisterField("density",&rho_gf);
408 paraview_dc.RegisterField("filtered_density",&rho_filter);
409 paraview_dc.Save();
410 }
411
412 // 11. Iterate:
413 for (int k = 1; k <= max_it; k++)
414 {
415 if (k > 1) { alpha *= ((real_t) k) / ((real_t) k-1); }
416
417 if (myid == 0)
418 {
419 cout << "\nStep = " << k << endl;
420 }
421
422 // Step 1 - Filter solve
423 // Solve (ϵ^2 ∇ ρ̃, ∇ v ) + (ρ̃,v) = (ρ,v)
424 FilterSolver->SetRHSCoefficient(&rho);
425 FilterSolver->Solve();
426 rho_filter = *FilterSolver->GetFEMSolution();
427
428 // Step 2 - State solve
429 // Solve (λ r(ρ̃) ∇⋅u, ∇⋅v) + (2 μ r(ρ̃) ε(u), ε(v)) = (f,v)
430 SIMPInterpolationCoefficient SIMP_cf(&rho_filter,rho_min, 1.0);
431 ProductCoefficient lambda_SIMP_cf(lambda_cf,SIMP_cf);
432 ProductCoefficient mu_SIMP_cf(mu_cf,SIMP_cf);
433 ElasticitySolver->SetLameCoefficients(&lambda_SIMP_cf,&mu_SIMP_cf);
434 ElasticitySolver->Solve();
435 u = *ElasticitySolver->GetFEMSolution();
436
437 // Step 3 - Adjoint filter solve
438 // Solve (ϵ² ∇ w̃, ∇ v) + (w̃ ,v) = (-r'(ρ̃) ( λ |∇⋅u|² + 2 μ |ε(u)|²),v)
439 StrainEnergyDensityCoefficient rhs_cf(&lambda_cf,&mu_cf,&u, &rho_filter,
440 rho_min);
441 FilterSolver->SetRHSCoefficient(&rhs_cf);
442 FilterSolver->Solve();
443 w_filter = *FilterSolver->GetFEMSolution();
444
445 // Step 4 - Compute gradient
446 // Solve G = M⁻¹w̃
447 GridFunctionCoefficient w_cf(&w_filter);
448 ParLinearForm w_rhs(&control_fes);
450 w_rhs.Assemble();
451 M.Mult(w_rhs,grad);
452
453 // Step 5 - Update design variable ψ ← proj(ψ - αG)
454 psi.Add(-alpha, grad);
455 const real_t material_volume = proj(psi, target_volume);
456
457 // Compute ||ρ - ρ_old|| in control fes.
458 real_t norm_increment = zerogf.ComputeL1Error(succ_diff_rho);
459 real_t norm_reduced_gradient = norm_increment/alpha;
460 psi_old = psi;
461
462 real_t compliance = (*(ElasticitySolver->GetLinearForm()))(u);
463 MPI_Allreduce(MPI_IN_PLACE, &compliance, 1, MPITypeMap<real_t>::mpi_type,
464 MPI_SUM, MPI_COMM_WORLD);
465 if (myid == 0)
466 {
467 mfem::out << "norm of the reduced gradient = " << norm_reduced_gradient << endl;
468 mfem::out << "norm of the increment = " << norm_increment << endl;
469 mfem::out << "compliance = " << compliance << endl;
470 mfem::out << "volume fraction = " << material_volume / domain_volume << endl;
471 }
472
473 if (glvis_visualization)
474 {
475 ParGridFunction r_gf(&filter_fes);
476 r_gf.ProjectCoefficient(SIMP_cf);
477 sout_r << "parallel " << num_procs << " " << myid << "\n";
478 sout_r << "solution\n" << pmesh << r_gf
479 << "window_title 'Design density r(ρ̃)'" << flush;
480 }
481
482 if (paraview_output)
483 {
484 rho_gf.ProjectCoefficient(rho);
485 paraview_dc.SetCycle(k);
486 paraview_dc.SetTime((real_t)k);
487 paraview_dc.Save();
488 }
489
490 if (norm_reduced_gradient < ntol && norm_increment < itol)
491 {
492 break;
493 }
494 }
495
496 delete ElasticitySolver;
497 delete FilterSolver;
498
499 return 0;
500}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
int Size() const
Return the logical size of the array.
Definition array.hpp:144
@ GaussLobatto
Closed type.
Definition fe_base.hpp:32
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
A coefficient that is constant across space and time.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetTime(real_t t)
Set physical time (for time-dependent simulations)
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
Returns f(u(x)) - f(v(x)) where u, v are scalar GridFunctions and f:R → R.
Definition ex37.hpp:66
Class for solving Poisson's equation:
Definition ex37.hpp:216
void SetDiffusionCoefficient(Coefficient *diffcf_)
Definition ex37.hpp:258
void SetRHSCoefficient(Coefficient *rhscf_)
Definition ex37.hpp:260
void SetOrder(int order_)
Definition ex37.hpp:257
void SetMesh(Mesh *mesh_)
Definition ex37.hpp:248
void SetMassCoefficient(Coefficient *masscf_)
Definition ex37.hpp:259
GridFunction * GetFEMSolution()
Definition ex37.hpp:534
void SetEssentialBoundary(const Array< int > &ess_bdr_)
Definition ex37.hpp:261
Class for domain integration .
Definition lininteg.hpp:109
Abstract data type element.
Definition element.hpp:29
virtual void GetVertices(Array< int > &v) const =0
Get the indices defining the vertices.
void SetAttribute(const int attr)
Set element's attribute.
Definition element.hpp:58
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Definition fe_coll.hpp:225
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, real_t alpha=1.0, real_t beta=0.0) const
Computes y = alpha * A * x + beta * y.
Definition hypre.cpp:1815
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.hpp:74
Integrator that inverts the matrix assembled by another integrator.
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
Class for solving linear elasticity:
Definition ex37.hpp:304
GridFunction * GetFEMSolution()
Definition ex37.hpp:716
void SetEssentialBoundary(const Array< int > &ess_bdr_)
Definition ex37.hpp:344
void SetRHSCoefficient(VectorCoefficient *rhs_cf_)
Definition ex37.hpp:343
void SetMesh(Mesh *mesh_)
Definition ex37.hpp:332
LinearForm * GetLinearForm()
Definition ex37.hpp:353
void SetOrder(int order_)
Definition ex37.hpp:341
void SetLameCoefficients(Coefficient *lambda_cf_, Coefficient *mu_cf_)
Definition ex37.hpp:342
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Returns f(u(x)) where u is a scalar GridFunction and f:R → R.
Definition ex37.hpp:41
Mesh data type.
Definition mesh.hpp:56
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:282
void Clear()
Clear the contents of the Mesh.
Definition mesh.hpp:730
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
const Element * GetBdrElement(int i) const
Return pointer to the i'th boundary element object.
Definition mesh.hpp:1298
int GetNBE() const
Returns number of boundary elements.
Definition mesh.hpp:1229
static Mesh MakeCartesian2D(int nx, int ny, Element::Type type, bool generate_edges=false, real_t sx=1.0, real_t sy=1.0, bool sfc_ordering=true)
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny quadrilaterals if type = QUADRILATER...
Definition mesh.cpp:4243
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
virtual void SetAttributes()
Determine the sets of unique attribute values in domain and boundary elements.
Definition mesh.cpp:1604
const real_t * GetVertex(int i) const
Return pointer to vertex i's coordinates.
Definition mesh.hpp:1265
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
Class for parallel bilinear form.
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
Abstract parallel finite element space.
Definition pfespace.hpp:29
HYPRE_BigInt GlobalTrueVSize() const
Definition pfespace.hpp:285
Class for parallel grid function.
Definition pgridfunc.hpp:33
real_t ComputeL1Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const override
ParFiniteElementSpace * ParFESpace() const
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Class for parallel linear form.
void Assemble()
Assembles the ParLinearForm i.e. sums over all domain/bdr integrators.
Class for parallel meshes.
Definition pmesh.hpp:34
Helper class for ParaView visualization data.
void SetHighOrderOutput(bool high_order_output_)
void SetLevelsOfDetail(int levels_of_detail_)
virtual void Save() override
void SetDataFormat(VTKFormat fmt)
Scalar coefficient defined as the product of two scalar coefficients or a scalar and a scalar coeffic...
Solid isotropic material penalization (SIMP) coefficient.
Definition ex37.hpp:98
Strain energy density coefficient.
Definition ex37.hpp:122
Vector data type.
Definition vector.hpp:80
real_t Sum() const
Return the sum of the vector entries.
Definition vector.cpp:1283
Vector & Add(const real_t a, const Vector &Va)
(*this) += a * Va
Definition vector.cpp:247
Volumetric force for linear elasticity.
Definition ex37.hpp:168
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
const real_t alpha
Definition ex15.cpp:369
int dim
Definition ex24.cpp:53
real_t mu
Definition ex25.cpp:140
real_t epsilon
Definition ex25.cpp:141
real_t proj(ParGridFunction &psi, real_t target_volume, real_t tol=1e-12, int max_its=10)
Bregman projection of ρ = sigmoid(ψ) onto the subspace ∫_Ω ρ dx = θ vol(Ω) as follows:
Definition ex37p.cpp:71
int main()
HYPRE_Int HYPRE_BigInt
const int visport
real_t sigmoid(real_t x)
Sigmoid function.
Definition ex37.hpp:20
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
void mfem_warning(const char *msg)
Definition error.cpp:187
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
real_t inv_sigmoid(real_t x)
Inverse sigmoid function.
Definition ex37.hpp:12
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
real_t der_sigmoid(real_t x)
Derivative of sigmoid function.
Definition ex37.hpp:33
const char vishost[]
Helper struct to convert a C++ type to an MPI type.