15 x = std::min(std::max(tol,x),
real_t(1.0)-tol);
16 return std::log(x/(1.0-x));
24 return 1.0/(1.0+std::exp(-x));
28 return std::exp(x)/(1.0+std::exp(x));
36 return tmp - std::pow(tmp,2);
91 return value1 - value2;
139 MFEM_ASSERT(rho_min_ >= 0.0,
"rho_min must be >= 0");
140 MFEM_ASSERT(rho_min_ < 1.0, "rho_min must be > 1
");
141 MFEM_ASSERT(u, "displacement field is not set
");
142 MFEM_ASSERT(rho_filter, "density field is not set
");
145 virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
147 real_t L = lambda->Eval(T, ip);
148 real_t M = mu->Eval(T, ip);
149 u->GetVectorGradient(T, grad);
150 real_t div_u = grad.Trace();
151 real_t density = L*div_u*div_u;
152 int dim = T.GetSpaceDim();
153 for (int i=0; i<dim; i++)
155 for (int j=0; j<dim; j++)
157 density += M*grad(i,j)*(grad(i,j)+grad(j,i));
160 real_t val = rho_filter->GetValue(T,ip);
162 return -exponent * pow(val, exponent-1.0) * (1-rho_min) * density;
167class VolumeForceCoefficient : public VectorCoefficient
174 VolumeForceCoefficient(real_t r_,Vector & center_, Vector & force_) :
175 VectorCoefficient(center_.Size()), r(r_), center(center_), force(force_) { }
177 using VectorCoefficient::Eval;
179 virtual void Eval(Vector &V, ElementTransformation &T,
180 const IntegrationPoint &ip)
182 Vector xx; xx.SetSize(T.GetDimension());
184 for (int i=0; i<xx.Size(); i++)
186 xx[i]=xx[i]-center[i];
189 real_t cr=xx.Norml2();
190 V.SetSize(T.GetDimension());
201 void Set(real_t r_,Vector & center_, Vector & force_)
218 Mesh * mesh = nullptr;
220 // diffusion coefficient
221 Coefficient * diffcf = nullptr;
223 Coefficient * masscf = nullptr;
224 Coefficient * rhscf = nullptr;
225 Coefficient * essbdr_cf = nullptr;
226 Coefficient * neumann_cf = nullptr;
227 VectorCoefficient * gradient_cf = nullptr;
231 FiniteElementCollection * fec = nullptr;
232 FiniteElementSpace * fes = nullptr;
234 Array<int> neumann_bdr;
235 GridFunction * u = nullptr;
236 LinearForm * b = nullptr;
239 ParMesh * pmesh = nullptr;
240 ParFiniteElementSpace * pfes = nullptr;
244 DiffusionSolver() { }
245 DiffusionSolver(Mesh * mesh_, int order_, Coefficient * diffcf_,
248 void SetMesh(Mesh * mesh_)
253 pmesh = dynamic_cast<ParMesh *>(mesh);
254 if (pmesh) { parallel = true; }
257 void SetOrder(int order_) { order = order_ ; }
258 void SetDiffusionCoefficient(Coefficient * diffcf_) { diffcf = diffcf_; }
259 void SetMassCoefficient(Coefficient * masscf_) { masscf = masscf_; }
260 void SetRHSCoefficient(Coefficient * rhscf_) { rhscf = rhscf_; }
261 void SetEssentialBoundary(const Array<int> & ess_bdr_) { ess_bdr = ess_bdr_;};
262 void SetNeumannBoundary(const Array<int> & neumann_bdr_) { neumann_bdr = neumann_bdr_;};
263 void SetNeumannData(Coefficient * neumann_cf_) {neumann_cf = neumann_cf_;}
264 void SetEssBdrData(Coefficient * essbdr_cf_) {essbdr_cf = essbdr_cf_;}
265 void SetGradientData(VectorCoefficient * gradient_cf_) {gradient_cf = gradient_cf_;}
271 GridFunction * GetFEMSolution();
272 LinearForm * GetLinearForm() {return b;}
274 ParGridFunction * GetParFEMSolution();
275 ParLinearForm * GetParLinearForm()
279 return dynamic_cast<ParLinearForm *>(b);
283 MFEM_ABORT("Wrong code path. Call GetLinearForm
");
303class LinearElasticitySolver
306 Mesh * mesh = nullptr;
308 Coefficient * lambda_cf = nullptr;
309 Coefficient * mu_cf = nullptr;
310 VectorCoefficient * essbdr_cf = nullptr;
311 VectorCoefficient * rhs_cf = nullptr;
315 FiniteElementCollection * fec = nullptr;
316 FiniteElementSpace * fes = nullptr;
318 Array<int> neumann_bdr;
319 GridFunction * u = nullptr;
320 LinearForm * b = nullptr;
323 ParMesh * pmesh = nullptr;
324 ParFiniteElementSpace * pfes = nullptr;
328 LinearElasticitySolver() { }
329 LinearElasticitySolver(Mesh * mesh_, int order_,
330 Coefficient * lambda_cf_, Coefficient * mu_cf_);
332 void SetMesh(Mesh * mesh_)
337 pmesh = dynamic_cast<ParMesh *>(mesh);
338 if (pmesh) { parallel = true; }
341 void SetOrder(int order_) { order = order_ ; }
342 void SetLameCoefficients(Coefficient * lambda_cf_, Coefficient * mu_cf_) { lambda_cf = lambda_cf_; mu_cf = mu_cf_; }
343 void SetRHSCoefficient(VectorCoefficient * rhs_cf_) { rhs_cf = rhs_cf_; }
344 void SetEssentialBoundary(const Array<int> & ess_bdr_) { ess_bdr = ess_bdr_;};
345 void SetNeumannBoundary(const Array<int> & neumann_bdr_) { neumann_bdr = neumann_bdr_;};
346 void SetEssBdrData(VectorCoefficient * essbdr_cf_) {essbdr_cf = essbdr_cf_;}
352 GridFunction * GetFEMSolution();
353 LinearForm * GetLinearForm() {return b;}
355 ParGridFunction * GetParFEMSolution();
356 ParLinearForm * GetParLinearForm()
360 return dynamic_cast<ParLinearForm *>(b);
364 MFEM_ABORT("Wrong code path. Call GetLinearForm
");
370 ~LinearElasticitySolver();
377DiffusionSolver::DiffusionSolver(Mesh * mesh_, int order_,
378 Coefficient * diffcf_, Coefficient * rhscf_)
379 : mesh(mesh_), order(order_), diffcf(diffcf_), rhscf(rhscf_)
383 pmesh = dynamic_cast<ParMesh *>(mesh);
384 if (pmesh) { parallel = true; }
390void DiffusionSolver::SetupFEM()
392 dim = mesh->Dimension();
393 fec = new H1_FECollection(order, dim);
398 pfes = new ParFiniteElementSpace(pmesh, fec);
399 u = new ParGridFunction(pfes);
400 b = new ParLinearForm(pfes);
404 fes = new FiniteElementSpace(mesh, fec);
405 u = new GridFunction(fes);
406 b = new LinearForm(fes);
409 fes = new FiniteElementSpace(mesh, fec);
410 u = new GridFunction(fes);
411 b = new LinearForm(fes);
417 if (mesh->bdr_attributes.Size())
419 ess_bdr.SetSize(mesh->bdr_attributes.Max());
425void DiffusionSolver::Solve()
429 Array<int> ess_tdof_list;
434 pfes->GetEssentialTrueDofs(ess_bdr,ess_tdof_list);
438 fes->GetEssentialTrueDofs(ess_bdr,ess_tdof_list);
441 fes->GetEssentialTrueDofs(ess_bdr,ess_tdof_list);
450 b = new ParLinearForm(pfes);
454 b = new LinearForm(fes);
457 b = new LinearForm(fes);
462 b->AddDomainIntegrator(new DomainLFIntegrator(*rhscf));
466 MFEM_VERIFY(neumann_bdr.Size(), "neumann_bdr attributes not provided
");
467 b->AddBoundaryIntegrator(new BoundaryLFIntegrator(*neumann_cf),neumann_bdr);
469 else if (gradient_cf)
471 MFEM_VERIFY(neumann_bdr.Size(), "neumann_bdr attributes not provided
");
472 b->AddBoundaryIntegrator(new BoundaryNormalLFIntegrator(*gradient_cf),
478 BilinearForm * a = nullptr;
483 a = new ParBilinearForm(pfes);
487 a = new BilinearForm(fes);
490 a = new BilinearForm(fes);
492 a->AddDomainIntegrator(new DiffusionIntegrator(*diffcf));
495 a->AddDomainIntegrator(new MassIntegrator(*masscf));
500 u->ProjectBdrCoefficient(*essbdr_cf,ess_bdr);
502 a->FormLinearSystem(ess_tdof_list, *u, *b, A, X, B);
504 CGSolver * cg = nullptr;
505 Solver * M = nullptr;
509 M = new HypreBoomerAMG;
510 dynamic_cast<HypreBoomerAMG*>(M)->SetPrintLevel(0);
511 cg = new CGSolver(pmesh->GetComm());
515 M = new GSSmoother((SparseMatrix&)(*A));
519 M = new GSSmoother((SparseMatrix&)(*A));
522 cg->SetRelTol(1e-12);
523 cg->SetMaxIter(10000);
524 cg->SetPrintLevel(0);
525 cg->SetPreconditioner(*M);
530 a->RecoverFEMSolution(X, *b, *u);
534GridFunction * DiffusionSolver::GetFEMSolution()
540ParGridFunction * DiffusionSolver::GetParFEMSolution()
544 return dynamic_cast<ParGridFunction*>(u);
548 MFEM_ABORT("Wrong code path. Call GetFEMSolution
");
554DiffusionSolver::~DiffusionSolver()
556 delete u; u = nullptr;
557 delete fes; fes = nullptr;
559 delete pfes; pfes=nullptr;
561 delete fec; fec = nullptr;
568LinearElasticitySolver::LinearElasticitySolver(Mesh * mesh_, int order_,
569 Coefficient * lambda_cf_, Coefficient * mu_cf_)
570 : mesh(mesh_), order(order_), lambda_cf(lambda_cf_), mu_cf(mu_cf_)
573 pmesh = dynamic_cast<ParMesh *>(mesh);
574 if (pmesh) { parallel = true; }
579void LinearElasticitySolver::SetupFEM()
581 dim = mesh->Dimension();
582 fec = new H1_FECollection(order, dim,BasisType::Positive);
587 pfes = new ParFiniteElementSpace(pmesh, fec, dim);
588 u = new ParGridFunction(pfes);
589 b = new ParLinearForm(pfes);
593 fes = new FiniteElementSpace(mesh, fec,dim);
594 u = new GridFunction(fes);
595 b = new LinearForm(fes);
598 fes = new FiniteElementSpace(mesh, fec, dim);
599 u = new GridFunction(fes);
600 b = new LinearForm(fes);
606 if (mesh->bdr_attributes.Size())
608 ess_bdr.SetSize(mesh->bdr_attributes.Max());
614void LinearElasticitySolver::Solve()
616 GridFunction * x = nullptr;
619 Array<int> ess_tdof_list;
624 x = new ParGridFunction(pfes);
625 pfes->GetEssentialTrueDofs(ess_bdr,ess_tdof_list);
629 x = new GridFunction(fes);
630 fes->GetEssentialTrueDofs(ess_bdr,ess_tdof_list);
633 x = new GridFunction(fes);
634 fes->GetEssentialTrueDofs(ess_bdr,ess_tdof_list);
643 b = new ParLinearForm(pfes);
647 b = new LinearForm(fes);
650 b = new LinearForm(fes);
655 b->AddDomainIntegrator(new VectorDomainLFIntegrator(*rhs_cf));
662 BilinearForm * a = nullptr;
667 a = new ParBilinearForm(pfes);
671 a = new BilinearForm(fes);
674 a = new BilinearForm(fes);
676 a->AddDomainIntegrator(new ElasticityIntegrator(*lambda_cf, *mu_cf));
680 u->ProjectBdrCoefficient(*essbdr_cf,ess_bdr);
682 a->FormLinearSystem(ess_tdof_list, *x, *b, A, X, B);
684 CGSolver * cg = nullptr;
685 Solver * M = nullptr;
689 M = new HypreBoomerAMG;
690 dynamic_cast<HypreBoomerAMG*>(M)->SetPrintLevel(0);
691 cg = new CGSolver(pmesh->GetComm());
695 M = new GSSmoother((SparseMatrix&)(*A));
699 M = new GSSmoother((SparseMatrix&)(*A));
702 cg->SetRelTol(1e-10);
703 cg->SetMaxIter(10000);
704 cg->SetPrintLevel(0);
705 cg->SetPreconditioner(*M);
710 a->RecoverFEMSolution(X, *b, *x);
716GridFunction * LinearElasticitySolver::GetFEMSolution()
722ParGridFunction * LinearElasticitySolver::GetParFEMSolution()
726 return dynamic_cast<ParGridFunction*>(u);
730 MFEM_ABORT("Wrong code path. Call GetFEMSolution
");
736LinearElasticitySolver::~LinearElasticitySolver()
738 delete u; u = nullptr;
739 delete fes; fes = nullptr;
741 delete pfes; pfes=nullptr;
743 delete fec; fec = nullptr;
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Data type dense matrix using column-major storage.
Returns f(u(x)) - f(v(x)) where u, v are scalar GridFunctions and f:R → R.
GridFunctionCoefficient OtherGridF_cf
void SetFunction(std::function< real_t(const real_t)> fun_)
DiffMappedGridFunctionCoefficient(const GridFunction *gf, const GridFunction *other_gf, std::function< real_t(const real_t)> fun_, int comp=1)
const GridFunction * OtherGridF
DiffMappedGridFunctionCoefficient()
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
std::function< real_t(const real_t)> fun
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
Class for grid function - Vector with associated FE space.
virtual real_t GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Class for integration point with weight.
Returns f(u(x)) where u is a scalar GridFunction and f:R → R.
MappedGridFunctionCoefficient()
void SetFunction(std::function< real_t(const real_t)> fun_)
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
std::function< real_t(const real_t)> fun
MappedGridFunctionCoefficient(const GridFunction *gf, std::function< real_t(const real_t)> fun_, int comp=1)
Solid isotropic material penalization (SIMP) coefficient.
SIMPInterpolationCoefficient(GridFunction *rho_filter_, real_t min_val_=1e-6, real_t max_val_=1.0, real_t exponent_=3)
GridFunction * rho_filter
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
Strain energy density coefficient.
StrainEnergyDensityCoefficient(Coefficient *lambda_, Coefficient *mu_, GridFunction *u_, GridFunction *rho_filter_, real_t rho_min_=1e-6, real_t exponent_=3.0)
GridFunction * rho_filter
real_t sigmoid(real_t x)
Sigmoid function.
real_t inv_sigmoid(real_t x)
Inverse sigmoid function.
real_t der_sigmoid(real_t x)
Derivative of sigmoid function.