feat(poly): constraint integrator

The NewtonSolver has been subclassed to try to auto enforce the zero boundary central condition by modifying the residual vector and the gradient matrix. This is a work in progress

BREAKING CHANGE:
This commit is contained in:
2025-03-05 12:55:53 -05:00
parent cd6da7065b
commit 59162a1a54
4 changed files with 435 additions and 202 deletions

View File

@@ -5,6 +5,8 @@
#include <numbers>
#include <csignal>
#include <fstream>
#include <array>
#include <vector>
#include "polyMFEMUtils.h"
#include "probe.h"
@@ -371,4 +373,323 @@ namespace polyMFEMUtils {
return gaussianIntegral(one_gf);
}
ZeroSlopeNewtonSolver::ZeroSlopeNewtonSolver(double alpha_, std::vector<double> zeroSlopeCoordinate_)
: alpha(alpha_), zeroSlopeCoordinate(zeroSlopeCoordinate_) {}
ZeroSlopeNewtonSolver::~ZeroSlopeNewtonSolver() {}
void ZeroSlopeNewtonSolver::SetOperator(const mfem::Operator &op) {
LOG_INFO(logger, "Setting operator for zero slope constraint...");
mfem::NewtonSolver::SetOperator(op); // Call the base class method
LOG_INFO(logger, "Setting operator for zero slope constraint...done");
LOG_INFO(logger, "Building location of zero slope constraint...");
mfem::NonlinearForm *nlf = dynamic_cast<mfem::NonlinearForm*>(const_cast<mfem::Operator*>(&op));
if (!nlf) {
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::SetOperator: input operator is not a NonlinearForm");
MFEM_ABORT("ZeroSlopeNewtonSolver::SetOperator: input operator is not a NonlinearForm");
}
mfem::FiniteElementSpace *fes = nlf->FESpace();
if (!fes) {
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::SetOperator: input operator does not have a finite element space");
MFEM_ABORT("ZeroSlopeNewtonSolver::SetOperator: input operator does not have a finite element space");
}
u_gf = std::make_unique<mfem::GridFunction>(fes);
mfem::Mesh *mesh = fes->GetMesh();
if (!mesh) {
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::SetOperator: input operator does not have a mesh");
MFEM_ABORT("ZeroSlopeNewtonSolver::SetOperator: input operator does not have a mesh");
}
if (mesh->SpaceDimension() != static_cast<int>(zeroSlopeCoordinate.size())) {
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::SetOperator: input operator mesh dimension does not match the zero slope coordinate dimension");
MFEM_ABORT("ZeroSlopeNewtonSolver::SetOperator: input operator mesh dimension does not match the zero slope coordinate dimension");
}
mfem::DenseMatrix zeroSlopeCoordinateMatrix(mesh->SpaceDimension(), 1);
for (int dimID = 0; dimID < mesh->SpaceDimension(); dimID++) {
zeroSlopeCoordinateMatrix(dimID, 0) = zeroSlopeCoordinate[dimID];
}
mfem::Array<int> elementsIDs;
mfem::Array<mfem::IntegrationPoint> ips;
mesh->FindPoints(zeroSlopeCoordinateMatrix, elementsIDs, ips);
zeroSlopeElemID = elementsIDs[0];
zeroSlopeIP = ips[0];
LOG_INFO(logger, "Getting element dofs for zero slope constraint...");
fes->GetElementDofs(zeroSlopeElemID, zeroSlopeDofs);
LOG_INFO(logger, "Getting element dofs for zero slope constraint...done");
LOG_INFO(logger, "Building location of zero slope constraint...done");
}
// void ZeroSlopeNewtonSolver::ProcessNewState(const mfem::Vector &x) const {
// LOG_INFO(logger, "Processing new state for zero slope constraint...");
// if (zeroSlopeElemID < 0) {
// LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ProcessNewState: zero slope element ID is not set");
// MFEM_ABORT("ZeroSlopeNewtonSolver::ProcessNewState: zero slope element ID is not set");
// }
// mfem::NonlinearForm *nlf = dynamic_cast<mfem::NonlinearForm*>(const_cast<mfem::Operator*>(oper));
// if (!nlf) {
// LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ProcessNewState: input operator is not a NonlinearForm");
// MFEM_ABORT("ZeroSlopeNewtonSolver::ProcessNewState: input operator is not a NonlinearForm");
// }
// mfem::FiniteElementSpace *fes = nlf->FESpace();
// if (!fes) {
// LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ProcessNewState: input operator does not have a finite element space");
// MFEM_ABORT("ZeroSlopeNewtonSolver::ProcessNewState: input operator does not have a finite element space");
// }
// mfem::Mesh *mesh = fes->GetMesh();
// if (!mesh) {
// LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ProcessNewState: input operator does not have a mesh");
// MFEM_ABORT("ZeroSlopeNewtonSolver::ProcessNewState: input operator does not have a mesh");
// }
// mfem::ElementTransformation *T = mesh->GetElementTransformation(zeroSlopeElemID);
// if (!T) {
// LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ProcessNewState: element transformation is not found");
// MFEM_ABORT("ZeroSlopeNewtonSolver::ProcessNewState: element transformation is not found");
// }
// mfem::Vector grad_u(3);
// mfem::GridFunction u_gf(fes);
// DEPRECATION_WARNING_OFF
// u_gf.SetData(x.GetData());
// DEPRECATION_WARNING_ON
// T->SetIntPoint(&zeroSlopeIP);
// u_gf.GetGradient(*T, grad_u);
// int dof;
// LOG_DEBUG(logger, "Adjusting the residual to enforce the zero slope constraint by {:0.4E}...", -alpha*grad_u[0]);
// double rNorm = r.Norml2();
// LOG_INFO(logger, "||r_B|| = {:0.4E}", rNorm);
// for (int i = 0; i < zeroSlopeDofs.Size(); i++) {
// dof = zeroSlopeDofs[i];
// r[dof] -= alpha * grad_u[0];
// r[dof] -= alpha * grad_u[1];
// r[dof] -= alpha * grad_u[2];
// }
// rNorm = r.Norml2();
// LOG_INFO(logger, "||r_A|| = {:0.4E}", rNorm);
// // This still is not working; however, I think I am close. I also need to modify the jacobain.
// }
void ZeroSlopeNewtonSolver::Mult(const mfem::Vector &b, mfem::Vector &x) const {
using namespace mfem;
using namespace std;
MFEM_VERIFY(oper != NULL, "the Operator is not set (use SetOperator).");
MFEM_VERIFY(prec != NULL, "the Solver is not set (use SetSolver).");
int it;
real_t norm0, norm, norm_goal;
const bool have_b = (b.Size() == Height());
if (!iterative_mode)
{
x = 0.0;
}
ProcessNewState(x);
oper->Mult(x, r);
if (have_b)
{
r -= b;
}
// ComputeConstrainedResidual(x, r);
norm0 = norm = initial_norm = Norm(r);
if (print_options.first_and_last && !print_options.iterations)
{
mfem::out << "Zero slope newton iteration " << setw(2) << 0
<< " : ||r|| = " << norm << "...\n";
}
norm_goal = std::max(rel_tol*norm, abs_tol);
prec->iterative_mode = false;
// x_{i+1} = x_i - [DF(x_i)]^{-1} [F(x_i)-b]
for (it = 0; true; it++)
{
MFEM_VERIFY(IsFinite(norm), "norm = " << norm);
if (print_options.iterations)
{
mfem::out << "Zero slope newton iteration " << setw(2) << it
<< " : ||r|| = " << norm;
if (it > 0)
{
mfem::out << ", ||r||/||r_0|| = " << norm/norm0;
}
mfem::out << '\n';
}
Monitor(it, norm, r, x);
if (norm <= norm_goal)
{
converged = true;
break;
}
if (it >= max_iter)
{
converged = false;
break;
}
grad = dynamic_cast<mfem::SparseMatrix*>(&oper->GetGradient(x));
if (!grad)
{
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::Mult: Operator does not return a SparseMatrix");
MFEM_ABORT("ZeroSlopeNewtonSolver::Mult: Operator does not return a SparseMatrix");
}
ComputeConstrainedGradient(x);
prec->SetOperator(*grad);
if (lin_rtol_type)
{
AdaptiveLinRtolPreSolve(x, it, norm);
}
prec->Mult(r, c); // c = [DF(x_i)]^{-1} [F(x_i)-b]
if (lin_rtol_type)
{
AdaptiveLinRtolPostSolve(c, r, it, norm);
}
const real_t c_scale = ComputeScalingFactor(x, b);
if (c_scale == 0.0)
{
converged = false;
break;
}
add(x, -c_scale, c, x);
ProcessNewState(x);
oper->Mult(x, r);
if (have_b)
{
r -= b;
}
// ComputeConstrainedResidual(x, r);
norm = Norm(r);
}
final_iter = it;
final_norm = norm;
if (print_options.summary || (!converged && print_options.warnings) ||
print_options.first_and_last)
{
mfem::out << "Newton: Number of iterations: " << final_iter << '\n'
<< " ||r|| = " << final_norm
<< ", ||r||/||r_0|| = " << final_norm/norm0 << '\n';
}
if (!converged && (print_options.summary || print_options.warnings))
{
mfem::out << "Newton: No convergence!\n";
}
}
void ZeroSlopeNewtonSolver::ComputeConstrainedResidual(const mfem::Vector &x, mfem::Vector &residual) const {
mfem::NonlinearForm *nlf = dynamic_cast<mfem::NonlinearForm*>(const_cast<mfem::Operator*>(oper));
if (!nlf) {
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ProcessNewState: input operator is not a NonlinearForm");
MFEM_ABORT("ZeroSlopeNewtonSolver::ProcessNewState: input operator is not a NonlinearForm");
}
mfem::FiniteElementSpace *fes = nlf->FESpace();
if (!fes) {
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ProcessNewState: input operator does not have a finite element space");
MFEM_ABORT("ZeroSlopeNewtonSolver::ProcessNewState: input operator does not have a finite element space");
}
mfem::Mesh *mesh = fes->GetMesh();
if (!mesh) {
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ProcessNewState: input operator does not have a mesh");
MFEM_ABORT("ZeroSlopeNewtonSolver::ProcessNewState: input operator does not have a mesh");
}
mfem::ElementTransformation *T = mesh->GetElementTransformation(zeroSlopeElemID);
if (!T) {
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ProcessNewState: element transformation is not found");
MFEM_ABORT("ZeroSlopeNewtonSolver::ProcessNewState: element transformation is not found");
}
DEPRECATION_WARNING_OFF
u_gf->SetData(x.GetData());
DEPRECATION_WARNING_ON
T->SetIntPoint(&zeroSlopeIP);
mfem::Vector grad_u(3); // TODO make this a unique pointer so it can be dimensionally adaptive
u_gf->GetGradient(*T, grad_u);
for (int i = 0; i < zeroSlopeDofs.Size(); i++) {
int dof = zeroSlopeDofs[i];
residual[dof] -= alpha * grad_u[0];
residual[dof] -= alpha * grad_u[1];
residual[dof] -= alpha * grad_u[2];
}
}
void ZeroSlopeNewtonSolver::ComputeConstrainedGradient(const mfem::Vector &x) const {
mfem::NonlinearForm *nlf = dynamic_cast<mfem::NonlinearForm*>(const_cast<mfem::Operator*>(oper));
if (!nlf) {
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ProcessNewState: input operator is not a NonlinearForm");
MFEM_ABORT("ZeroSlopeNewtonSolver::ProcessNewState: input operator is not a NonlinearForm");
}
mfem::FiniteElementSpace *fes = nlf->FESpace();
if (!fes) {
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ProcessNewState: input operator does not have a finite element space");
MFEM_ABORT("ZeroSlopeNewtonSolver::ProcessNewState: input operator does not have a finite element space");
}
mfem::Mesh *mesh = fes->GetMesh();
if (!mesh) {
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ProcessNewState: input operator does not have a mesh");
MFEM_ABORT("ZeroSlopeNewtonSolver::ProcessNewState: input operator does not have a mesh");
}
mfem::ElementTransformation *T = mesh->GetElementTransformation(zeroSlopeElemID);
if (!T) {
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ProcessNewState: element transformation is not found");
MFEM_ABORT("ZeroSlopeNewtonSolver::ProcessNewState: element transformation is not found");
}
const mfem::FiniteElement* fe = fes->GetFE(zeroSlopeElemID); // Get FE *once*.
mfem::DenseMatrix dshape; // For shape function derivatives.
dshape.SetSize(fe->GetDof(), mesh->Dimension());
T->SetIntPoint(&zeroSlopeIP);
fe->CalcDShape(zeroSlopeIP, dshape);
if (!grad) {
LOG_ERROR(logger, "ZeroSlopeNewtonSolver::ComputeConstrainedGradient: Grad is not set");
MFEM_ABORT("ZeroSlopeNewtonSolver::ComputeConstrainedGradient: Grad is not set");
}
// --- Modify Jacobian ---
LOG_INFO(logger, "Adjusting the Jacobian to enforce the zero slope constraint...");
for (int i = 0; i < zeroSlopeDofs.Size(); i++) {
for (int j = 0; j < zeroSlopeDofs.Size(); j++) {
grad->Add(zeroSlopeDofs[i], zeroSlopeDofs[j], alpha * dshape(j, 0));
grad->Add(zeroSlopeDofs[i], zeroSlopeDofs[j], alpha * dshape(j, 1));
grad->Add(zeroSlopeDofs[i], zeroSlopeDofs[j], alpha * dshape(j, 2));
}
}
LOG_INFO(logger, "Adjusting the Jacobian to enforce the zero slope constraint...done");
}
} // namespace polyMFEMUtils