2025-02-12 16:44:10 -05:00
|
|
|
#include "mfem.hpp"
|
|
|
|
|
#include <string>
|
|
|
|
|
#include <iostream>
|
|
|
|
|
#include <cmath>
|
2025-02-19 14:35:15 -05:00
|
|
|
#include <numbers>
|
2025-02-12 16:44:10 -05:00
|
|
|
|
|
|
|
|
#include "polyMFEMUtils.h"
|
|
|
|
|
|
2025-02-20 16:05:02 -05:00
|
|
|
#include "warning_control.h"
|
|
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
namespace polyMFEMUtils {
|
|
|
|
|
NonlinearPowerIntegrator::NonlinearPowerIntegrator(
|
|
|
|
|
mfem::Coefficient &coeff,
|
|
|
|
|
double n) : coeff_(coeff), polytropicIndex(n) {
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
void NonlinearPowerIntegrator::AssembleElementVector(
|
|
|
|
|
const mfem::FiniteElement &el,
|
|
|
|
|
mfem::ElementTransformation &Trans,
|
|
|
|
|
const mfem::Vector &elfun,
|
|
|
|
|
mfem::Vector &elvect) {
|
|
|
|
|
|
|
|
|
|
const mfem::IntegrationRule *ir = &mfem::IntRules.Get(el.GetGeomType(), 2 * el.GetOrder() + 3);
|
|
|
|
|
int dof = el.GetDof();
|
|
|
|
|
elvect.SetSize(dof);
|
|
|
|
|
elvect = 0.0;
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
mfem::Vector shape(dof);
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
for (int iqp = 0; iqp < ir->GetNPoints(); iqp++) {
|
|
|
|
|
mfem::IntegrationPoint ip = ir->IntPoint(iqp);
|
|
|
|
|
Trans.SetIntPoint(&ip);
|
|
|
|
|
double weight = ip.weight * Trans.Weight();
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
el.CalcShape(ip, shape);
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
double u_val = 0.0;
|
|
|
|
|
for (int j = 0; j < dof; j++) {
|
|
|
|
|
u_val += elfun(j) * shape(j);
|
|
|
|
|
}
|
|
|
|
|
double u_safe = std::max(u_val, 0.0);
|
|
|
|
|
double u_nl = std::pow(u_safe, polytropicIndex);
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
double coeff_val = coeff_.Eval(Trans, ip);
|
|
|
|
|
double x2_u_nl = coeff_val * u_nl;
|
|
|
|
|
|
|
|
|
|
for (int i = 0; i < dof; i++){
|
|
|
|
|
elvect(i) += shape(i) * x2_u_nl * weight;
|
|
|
|
|
}
|
2025-02-12 16:44:10 -05:00
|
|
|
}
|
2025-02-19 14:35:15 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
void NonlinearPowerIntegrator::AssembleElementGrad (
|
|
|
|
|
const mfem::FiniteElement &el,
|
|
|
|
|
mfem::ElementTransformation &Trans,
|
|
|
|
|
const mfem::Vector &elfun,
|
|
|
|
|
mfem::DenseMatrix &elmat) {
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
const mfem::IntegrationRule *ir = &mfem::IntRules.Get(el.GetGeomType(), 2 * el.GetOrder() + 3);
|
|
|
|
|
int dof = el.GetDof();
|
|
|
|
|
elmat.SetSize(dof);
|
|
|
|
|
elmat = 0.0;
|
|
|
|
|
mfem::Vector shape(dof);
|
|
|
|
|
|
|
|
|
|
for (int iqp = 0; iqp < ir->GetNPoints(); iqp++) {
|
|
|
|
|
const mfem::IntegrationPoint &ip = ir->IntPoint(iqp);
|
|
|
|
|
Trans.SetIntPoint(&ip);
|
|
|
|
|
double weight = ip.weight * Trans.Weight();
|
|
|
|
|
|
|
|
|
|
el.CalcShape(ip, shape);
|
|
|
|
|
|
|
|
|
|
double u_val = 0.0;
|
|
|
|
|
|
|
|
|
|
for (int j = 0; j < dof; j++) {
|
|
|
|
|
u_val += elfun(j) * shape(j);
|
|
|
|
|
}
|
|
|
|
|
double coeff_val = coeff_.Eval(Trans, ip);
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
// Calculate the Jacobian
|
|
|
|
|
double u_safe = std::max(u_val, 0.0);
|
|
|
|
|
double d_u_nl = coeff_val * polytropicIndex * std::pow(u_safe, polytropicIndex - 1);
|
|
|
|
|
double x2_d_u_nl = d_u_nl;
|
|
|
|
|
|
|
|
|
|
for (int i = 0; i < dof; i++) {
|
|
|
|
|
for (int j = 0; j < dof; j++) {
|
|
|
|
|
elmat(i, j) += shape(i) * x2_d_u_nl * shape(j) * weight;
|
|
|
|
|
}
|
|
|
|
|
}
|
2025-02-12 16:44:10 -05:00
|
|
|
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
BilinearIntegratorWrapper::BilinearIntegratorWrapper(
|
|
|
|
|
mfem::BilinearFormIntegrator *integratorInput
|
|
|
|
|
) : integrator(integratorInput) { }
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
BilinearIntegratorWrapper::~BilinearIntegratorWrapper() {
|
|
|
|
|
delete integrator;
|
|
|
|
|
}
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
void BilinearIntegratorWrapper::AssembleElementVector(
|
|
|
|
|
const mfem::FiniteElement &el,
|
|
|
|
|
mfem::ElementTransformation &Trans,
|
|
|
|
|
const mfem::Vector &elfun,
|
|
|
|
|
mfem::Vector &elvect) {
|
|
|
|
|
int dof = el.GetDof();
|
|
|
|
|
mfem::DenseMatrix elMat(dof);
|
|
|
|
|
integrator->AssembleElementMatrix(el, Trans, elMat);
|
|
|
|
|
elvect.SetSize(dof);
|
|
|
|
|
elvect = 0.0;
|
|
|
|
|
for (int i = 0; i < dof; i++)
|
|
|
|
|
{
|
|
|
|
|
double sum = 0.0;
|
|
|
|
|
for (int j = 0; j < dof; j++)
|
|
|
|
|
{
|
|
|
|
|
sum += elMat(i, j) * elfun(j);
|
|
|
|
|
}
|
|
|
|
|
elvect(i) = sum;
|
|
|
|
|
}
|
|
|
|
|
}
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
void BilinearIntegratorWrapper::AssembleElementGrad(const mfem::FiniteElement &el,
|
|
|
|
|
mfem::ElementTransformation &Trans,
|
|
|
|
|
const mfem::Vector &elfun,
|
|
|
|
|
mfem::DenseMatrix &elmat) {
|
|
|
|
|
int dof = el.GetDof();
|
|
|
|
|
elmat.SetSize(dof, dof);
|
|
|
|
|
elmat = 0.0;
|
|
|
|
|
integrator->AssembleElementMatrix(el, Trans, elmat);
|
|
|
|
|
}
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
CompositeNonlinearIntegrator::CompositeNonlinearIntegrator() { }
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
|
2025-02-20 15:28:00 -05:00
|
|
|
CompositeNonlinearIntegrator::~CompositeNonlinearIntegrator() { }
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
void CompositeNonlinearIntegrator::add_integrator(mfem::NonlinearFormIntegrator *integrator) {
|
|
|
|
|
integrators.push_back(integrator);
|
|
|
|
|
}
|
2025-02-12 16:44:10 -05:00
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
void CompositeNonlinearIntegrator::AssembleElementVector(
|
|
|
|
|
const mfem::FiniteElement &el,
|
|
|
|
|
mfem::ElementTransformation &Trans,
|
|
|
|
|
const mfem::Vector &elfun,
|
|
|
|
|
mfem::Vector &elvect) {
|
|
|
|
|
int dof = el.GetDof();
|
|
|
|
|
elvect.SetSize(dof);
|
|
|
|
|
elvect = 0.0;
|
|
|
|
|
mfem::Vector temp(dof);
|
|
|
|
|
|
|
|
|
|
for (size_t i = 0; i < integrators.size(); i++) {
|
|
|
|
|
temp= 0.0;
|
|
|
|
|
integrators[i]->AssembleElementVector(el, Trans, elfun, temp);
|
|
|
|
|
elvect.Add(1.0, temp);
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
void CompositeNonlinearIntegrator::AssembleElementGrad(
|
|
|
|
|
const mfem::FiniteElement &el,
|
|
|
|
|
mfem::ElementTransformation &Trans,
|
|
|
|
|
const mfem::Vector &elfun,
|
|
|
|
|
mfem::DenseMatrix &elmat) {
|
|
|
|
|
int dof = el.GetDof();
|
|
|
|
|
elmat.SetSize(dof, dof);
|
|
|
|
|
elmat = 0.0;
|
|
|
|
|
mfem::DenseMatrix temp(dof);
|
|
|
|
|
temp.SetSize(dof, dof);
|
|
|
|
|
for (size_t i = 0; i < integrators.size(); i++) {
|
|
|
|
|
temp = 0.0;
|
|
|
|
|
integrators[i] -> AssembleElementGrad(el, Trans, elfun, temp);
|
|
|
|
|
elmat.Add(1.0, temp);
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
ConstraintIntegrator::ConstraintIntegrator(mfem::Coefficient &eta_) : eta(eta_) {}
|
|
|
|
|
|
|
|
|
|
void ConstraintIntegrator::AssembleElementVector(const mfem::FiniteElement &el, mfem::ElementTransformation &Trans, const mfem::Vector &elfun, mfem::Vector &elvect) {
|
|
|
|
|
int dof = el.GetDof();
|
|
|
|
|
elvect.SetSize(dof);
|
|
|
|
|
elvect = 0.0;
|
|
|
|
|
|
|
|
|
|
mfem::Vector shape(dof);
|
|
|
|
|
const int intOrder = 2 * el.GetOrder() + 3;
|
|
|
|
|
|
2025-02-20 15:28:00 -05:00
|
|
|
const mfem::IntegrationRule &ir = mfem::IntRules.Get(el.GetGeomType(), intOrder);
|
2025-02-19 14:35:15 -05:00
|
|
|
|
|
|
|
|
for (int i = 0; i < ir.GetNPoints(); i++) {
|
|
|
|
|
const mfem::IntegrationPoint &ip = ir.IntPoint(i);
|
|
|
|
|
Trans.SetIntPoint(&ip);
|
|
|
|
|
el.CalcShape(ip, shape);
|
|
|
|
|
|
|
|
|
|
double eta_val = eta.Eval(Trans, ip);
|
|
|
|
|
|
|
|
|
|
double weight = ip.weight * Trans.Weight();
|
|
|
|
|
elvect.Add(eta_val * weight, shape);
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
void ConstraintIntegrator::AssembleElementGrad(const mfem::FiniteElement &el, mfem::ElementTransformation &Trans, const mfem::Vector &elfun, mfem::DenseMatrix &elmat) {
|
|
|
|
|
int dof = el.GetDof();
|
|
|
|
|
elmat.SetSize(dof);
|
|
|
|
|
elmat = 0.0;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
GaussianCoefficient::GaussianCoefficient(double stdDev_)
|
|
|
|
|
: stdDev(stdDev_),
|
|
|
|
|
norm_coeff(1.0/(std::pow(std::sqrt(2*std::numbers::pi)*stdDev_,3))) {}
|
|
|
|
|
|
2025-02-20 15:28:00 -05:00
|
|
|
double GaussianCoefficient::Eval(mfem::ElementTransformation &T, const mfem::IntegrationPoint &ip) {
|
2025-02-19 14:35:15 -05:00
|
|
|
mfem::Vector r;
|
|
|
|
|
T.Transform(ip, r);
|
|
|
|
|
double r2 = r * r;
|
|
|
|
|
return norm_coeff * std::exp(-r2/(2*std::pow(stdDev, 2)));
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
AugmentedOperator::AugmentedOperator(mfem::NonlinearForm &nfl_, mfem::LinearForm &C_, int lambdaDofOffset_)
|
|
|
|
|
:
|
|
|
|
|
mfem::Operator( // Call the base class constructor
|
2025-02-20 15:28:00 -05:00
|
|
|
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)
|
2025-02-19 14:35:15 -05:00
|
|
|
),
|
2025-02-20 15:28:00 -05:00
|
|
|
nfl(nfl_),
|
|
|
|
|
C(C_),
|
2025-02-19 14:35:15 -05:00
|
|
|
lambdaDofOffset(lambdaDofOffset_),
|
|
|
|
|
lastJacobian(nullptr) {}
|
|
|
|
|
|
|
|
|
|
void AugmentedOperator::Mult(const mfem::Vector &x, mfem::Vector &y) 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];
|
|
|
|
|
|
|
|
|
|
// Compute the residual of the nonlinear form (F(u) - λC)
|
|
|
|
|
mfem::Vector F(lambdaDofOffset);
|
2025-02-20 15:28:00 -05:00
|
|
|
nfl.Mult(u, F); // This now includes the -λ∫vη(r) dΩ term
|
2025-02-19 14:35:15 -05:00
|
|
|
|
|
|
|
|
// Compute the transpose of C for the off diagonal terms of the augmented operator
|
|
|
|
|
y.SetSize(height);
|
|
|
|
|
y = 0.0;
|
|
|
|
|
|
2025-02-20 15:28:00 -05:00
|
|
|
for (int i = 0; i < lambdaDofOffset; i++) {
|
|
|
|
|
y[i] = F[i];
|
|
|
|
|
}
|
|
|
|
|
mfem::GridFunction u_gf(C.FESpace());
|
2025-02-20 16:05:02 -05:00
|
|
|
|
|
|
|
|
DEPRECATION_WARNING_OFF
|
|
|
|
|
u_gf.SetData(u);
|
|
|
|
|
DEPRECATION_WARNING_ON
|
|
|
|
|
|
2025-02-20 15:28:00 -05:00
|
|
|
y[lambdaDofOffset] = C.operator()(u_gf);
|
2025-02-19 14:35:15 -05:00
|
|
|
|
|
|
|
|
// add -lambda * C to the residual
|
|
|
|
|
mfem::Vector lambda_C(lambdaDofOffset);
|
2025-02-20 15:28:00 -05:00
|
|
|
mfem::GridFunction constraint_gf(C.FESpace());
|
|
|
|
|
constraint_gf = 0.0;
|
2025-02-19 14:35:15 -05:00
|
|
|
|
2025-02-20 15:28:00 -05:00
|
|
|
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];
|
|
|
|
|
}
|
2025-02-19 14:35:15 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
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);
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
// Fill in the blocks of the augmented Jacobian matrix
|
|
|
|
|
// Top-Left: Jacobian of the nonlinear form
|
2025-02-20 15:28:00 -05:00
|
|
|
mfem::SparseMatrix *Jnfl_sparse = dynamic_cast<mfem::SparseMatrix*>(&nfl.GetGradient(u));
|
|
|
|
|
if (!Jnfl_sparse) {
|
|
|
|
|
MFEM_ABORT("GetGradient did not return a SparseMatrix");
|
|
|
|
|
}
|
2025-02-19 14:35:15 -05:00
|
|
|
|
2025-02-20 15:28:00 -05:00
|
|
|
mfem::SparseMatrix *J_aug = new mfem::SparseMatrix(height, width);
|
2025-02-19 14:35:15 -05:00
|
|
|
// Copy the original Jacobian into the augmented Jacobian
|
2025-02-20 15:28:00 -05:00
|
|
|
for (int i = 0; i < Jnfl_sparse->Height(); i++) {
|
2025-02-20 16:05:02 -05:00
|
|
|
const int *J_cols = Jnfl_sparse->GetRowColumns(i);
|
|
|
|
|
const double *J_vals = Jnfl_sparse->GetRowEntries(i);
|
|
|
|
|
|
|
|
|
|
for (int jj = 0; jj < Jnfl_sparse->RowSize(i); jj++) {
|
|
|
|
|
int j = J_cols[jj];
|
|
|
|
|
double val = J_vals[jj];
|
|
|
|
|
J_aug->Set(i, j, val);
|
2025-02-12 16:44:10 -05:00
|
|
|
}
|
2025-02-19 14:35:15 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// Bottom-left C (the constraint matrix)
|
|
|
|
|
mfem::Vector CVec(lambdaDofOffset);
|
2025-02-20 15:28:00 -05:00
|
|
|
mfem::GridFunction tempCGrid(C.FESpace());
|
|
|
|
|
C.Assemble();
|
|
|
|
|
CVec = C.GetData();
|
2025-02-19 14:35:15 -05:00
|
|
|
for (int i = 0; i < CVec.Size(); i++) {
|
|
|
|
|
J_aug->Set(lambdaDofOffset, i, CVec[i]);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// Top-right -Cᵀ (the negative transpose of the constraint matrix)
|
|
|
|
|
for (int i = 0; i < CVec.Size(); i++) {
|
|
|
|
|
J_aug->Set(i, lambdaDofOffset, -CVec[i]);
|
2025-02-12 16:44:10 -05:00
|
|
|
}
|
|
|
|
|
|
2025-02-19 14:35:15 -05:00
|
|
|
J_aug->Finalize();
|
|
|
|
|
|
|
|
|
|
delete lastJacobian;
|
|
|
|
|
lastJacobian = J_aug;
|
|
|
|
|
return *lastJacobian;
|
|
|
|
|
}
|
|
|
|
|
|
2025-02-20 15:28:00 -05:00
|
|
|
AugmentedOperator::~AugmentedOperator() {
|
2025-02-19 14:35:15 -05:00
|
|
|
delete lastJacobian;
|
2025-02-12 16:44:10 -05:00
|
|
|
}
|
2025-02-19 14:35:15 -05:00
|
|
|
} // namespace polyMFEMUtils
|