feat(poly): lagrangian constrained weak form of 3D lane-Emden

added a basic implimentation of the 3D lane emden equation using a lagrangian multiplier to constrain the value at the center of a spherical domain
This commit is contained in:
2025-02-20 15:28:00 -05:00
parent deab5be0c1
commit 1fd1e624f2
4 changed files with 89 additions and 85 deletions

View File

@@ -132,11 +132,7 @@ namespace polyMFEMUtils {
CompositeNonlinearIntegrator::CompositeNonlinearIntegrator() { }
CompositeNonlinearIntegrator::~CompositeNonlinearIntegrator() {
for (size_t i = 0; i < integrators.size(); i++) {
delete integrators[i];
}
}
CompositeNonlinearIntegrator::~CompositeNonlinearIntegrator() { }
void CompositeNonlinearIntegrator::add_integrator(mfem::NonlinearFormIntegrator *integrator) {
integrators.push_back(integrator);
@@ -186,7 +182,7 @@ namespace polyMFEMUtils {
mfem::Vector shape(dof);
const int intOrder = 2 * el.GetOrder() + 3;
mfem::IntegrationRule ir(el.GetGeomType(), intOrder);
const mfem::IntegrationRule &ir = mfem::IntRules.Get(el.GetGeomType(), intOrder);
for (int i = 0; i < ir.GetNPoints(); i++) {
const mfem::IntegrationPoint &ip = ir.IntPoint(i);
@@ -210,7 +206,7 @@ namespace polyMFEMUtils {
: stdDev(stdDev_),
norm_coeff(1.0/(std::pow(std::sqrt(2*std::numbers::pi)*stdDev_,3))) {}
GaussianCoefficient::Eval(mfem::ElementTransformation &T, const mfem::IntegrationPoint &ip) {
double GaussianCoefficient::Eval(mfem::ElementTransformation &T, const mfem::IntegrationPoint &ip) {
mfem::Vector r;
T.Transform(ip, r);
double r2 = r * r;
@@ -221,11 +217,11 @@ namespace polyMFEMUtils {
AugmentedOperator::AugmentedOperator(mfem::NonlinearForm &nfl_, mfem::LinearForm &C_, int lambdaDofOffset_)
:
mfem::Operator( // Call the base class constructor
nfl_->FESpace()->GetTrueVSize()+1, // Sets the height attribute (+1 for the lambda component)
nfl_->FESpace()->GetTrueVSize()+1 // Sets the width attribute (+1 for the lambda component)
nfl_.FESpace()->GetTrueVSize()+1, // Sets the height attribute (+1 for the lambda component)
nfl_.FESpace()->GetTrueVSize()+1 // Sets the width attribute (+1 for the lambda component)
),
nfl(&nfl_),
C(&C_),
nfl(nfl_),
C(C_),
lambdaDofOffset(lambdaDofOffset_),
lastJacobian(nullptr) {}
@@ -236,51 +232,59 @@ namespace polyMFEMUtils {
// Compute the residual of the nonlinear form (F(u) - λC)
mfem::Vector F(lambdaDofOffset);
nfl->Mult(u, F); // This now includes the -λ∫vη(r) dΩ term
nfl.Mult(u, F); // This now includes the -λ∫vη(r) dΩ term
// Compute the transpose of C for the off diagonal terms of the augmented operator
y.SetSize(height);
y = 0.0;
y.SetBlock(0, F);
y[lambdaDofOffset] = (*C)(u); // Cᵀ×u. This is equivalent to ∫ η(r)θ dΩ
for (int i = 0; i < lambdaDofOffset; i++) {
y[i] = F[i];
}
mfem::GridFunction u_gf(C.FESpace());
u_gf.SetData(u);
y[lambdaDofOffset] = C.operator()(u_gf);
// add -lambda * C to the residual
mfem::Vector lambda_C(lambdaDofOffset);
C->GetVector(lambda_C) // Get the coefficient vector for C. I.e. ∫ vη(r) dΩ
lambda_C *= -lambda; // Multiply by -λ (this is now the term −λ ∫ vη(r)dΩ)
mfem::GridFunction constraint_gf(C.FESpace());
constraint_gf = 0.0;
y.AddBlockVector(0, lambda_C); // Add the term to the residual
mfem::Vector constraint_tdofs;
C.Assemble();
lambda_C = C.GetData();
lambda_C *= -lambda; // Multiply by -λ (this is now the term −λ ∫ vη(r)dΩ)
for (int i = 0; i < lambdaDofOffset; i++) {
y[i] += lambda_C[i];
}
}
mfem::Operator &AugmentedOperator::GetGradient(const mfem::Vector &x) const {
// Select the lambda component of the input vector and seperate it from the θ component
mfem::Vector u(x.GetData(), lambdaDofOffset);
double lambda = x[lambdaDofOffset];
// --- Create the augmented Jacobian matrix ---
mfem::Operator *Jnfl = nfl->GetGradient(u); // Get the gradient of the nonlinear form
mfem::SparseMatrix *J_aug = new mfem::SparseMatrix(height, width);
// Fill in the blocks of the augmented Jacobian matrix
// Top-Left: Jacobian of the nonlinear form
mfem::SparseMatrix *Jnfl_sparse = dynamic_cast<mfem::SparseMatrix*>(Jnfl);
// Copy the original Jacobian into the augmented Jacobian
if (Jnfl_sparse) { // Checks if the dynamic cast was successful
for (int i = 0; i < Jnfl_sparse->Height(); i++) {
for (int j = 0; j < Jnfl_sparse->Width(); j++) {
J_aug->Set(i, j, Jnfl_sparse->Get(i, j));
}
}
} else {
mfem::SparseMatrix *Jnfl_sparse = dynamic_cast<mfem::SparseMatrix*>(&nfl.GetGradient(u));
if (!Jnfl_sparse) {
MFEM_ABORT("GetGradient did not return a SparseMatrix");
}
mfem::SparseMatrix *J_aug = new mfem::SparseMatrix(height, width);
// Copy the original Jacobian into the augmented Jacobian
for (int i = 0; i < Jnfl_sparse->Height(); i++) {
for (int j = 0; j < Jnfl_sparse->Width(); j++) {
J_aug->Set(i, j, Jnfl_sparse->Elem(i, j));
}
}
// Bottom-left C (the constraint matrix)
mfem::Vector CVec(lambdaDofOffset);
C->GetVector(CVec);
mfem::GridFunction tempCGrid(C.FESpace());
C.Assemble();
CVec = C.GetData();
for (int i = 0; i < CVec.Size(); i++) {
J_aug->Set(lambdaDofOffset, i, CVec[i]);
}
@@ -297,7 +301,7 @@ namespace polyMFEMUtils {
return *lastJacobian;
}
AugmentedOperator::~AugAugmentedOperator() {
AugmentedOperator::~AugmentedOperator() {
delete lastJacobian;
}
} // namespace polyMFEMUtils