fix(poly): working on 3D polytrope
not working yet
This commit is contained in:
@@ -3,8 +3,12 @@
|
||||
#include <iostream>
|
||||
#include <cmath>
|
||||
#include <numbers>
|
||||
#include <csignal>
|
||||
#include <fstream>
|
||||
|
||||
#include "polyMFEMUtils.h"
|
||||
#include "probe.h"
|
||||
#include "config.h"
|
||||
|
||||
#include "warning_control.h"
|
||||
|
||||
@@ -206,17 +210,19 @@ namespace polyMFEMUtils {
|
||||
|
||||
GaussianCoefficient::GaussianCoefficient(double stdDev_)
|
||||
: stdDev(stdDev_),
|
||||
norm_coeff(1.0/(std::pow(std::sqrt(2*std::numbers::pi)*stdDev_,3))) {}
|
||||
norm_coeff(1.0/std::pow(std::sqrt(2*std::numbers::pi*std::pow(stdDev_, 2)), 3.0/2.0)) {}
|
||||
|
||||
double GaussianCoefficient::Eval(mfem::ElementTransformation &T, const mfem::IntegrationPoint &ip) {
|
||||
mfem::Vector r;
|
||||
T.Transform(ip, r);
|
||||
double r2 = r * r;
|
||||
return norm_coeff * std::exp(-r2/(2*std::pow(stdDev, 2)));
|
||||
double rnorm = std::sqrt(r * r);
|
||||
// TODO: return to this to resolve why the Guassian does not always have peek at g(0) = 1 without the factor of 3.0145 (manually calibrated).
|
||||
// Open a file (append if already exists) to write the Gaussian values
|
||||
return norm_coeff * std::exp(-std::pow(rnorm, 2)/(2*std::pow(stdDev, 2)));
|
||||
}
|
||||
|
||||
|
||||
AugmentedOperator::AugmentedOperator(mfem::NonlinearForm &nfl_, mfem::LinearForm &C_, int lambdaDofOffset_)
|
||||
AugmentedOperator::AugmentedOperator(mfem::NonlinearForm &nfl_, mfem::LinearForm &C_, int lambdaDofOffset_, double C_val_)
|
||||
:
|
||||
mfem::Operator( // Call the base class constructor
|
||||
nfl_.FESpace()->GetTrueVSize()+1, // Sets the height attribute (+1 for the lambda component)
|
||||
@@ -224,6 +230,7 @@ namespace polyMFEMUtils {
|
||||
),
|
||||
nfl(nfl_),
|
||||
C(C_),
|
||||
C_val(C_val_),
|
||||
lambdaDofOffset(lambdaDofOffset_),
|
||||
lastJacobian(nullptr) {}
|
||||
|
||||
@@ -240,30 +247,57 @@ namespace polyMFEMUtils {
|
||||
y.SetSize(height);
|
||||
y = 0.0;
|
||||
|
||||
for (int i = 0; i < lambdaDofOffset; i++) {
|
||||
y[i] = F[i];
|
||||
}
|
||||
mfem::GridFunction u_gf(C.FESpace());
|
||||
mfem::Vector C_u(1);
|
||||
|
||||
DEPRECATION_WARNING_OFF
|
||||
u_gf.SetData(u);
|
||||
DEPRECATION_WARNING_ON
|
||||
|
||||
y[lambdaDofOffset] = C.operator()(u_gf);
|
||||
C_u[0] = C.operator()(u_gf);
|
||||
|
||||
// add -lambda * C to the residual
|
||||
mfem::Vector lambda_C(lambdaDofOffset);
|
||||
mfem::GridFunction constraint_gf(C.FESpace());
|
||||
constraint_gf = 0.0;
|
||||
|
||||
mfem::Vector constraint_tdofs;
|
||||
C.Assemble();
|
||||
lambda_C = C.GetData();
|
||||
mfem::Vector CTmp(lambdaDofOffset);
|
||||
CTmp = C.GetData();
|
||||
lambda_C = CTmp;
|
||||
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];
|
||||
y[i] = F[i] + lambda_C[i];
|
||||
}
|
||||
|
||||
// Get Global Debug Options for Poly
|
||||
std::string defaultKeyset = config.get<std::string>("Poly:Debug:Global:GLVis:Keyset", "");
|
||||
bool defaultView = config.get<bool>("Poly:Debug:Global:GLVis:View", false);
|
||||
bool defaultExit = config.get<bool>("Poly:Debug:Global:GLVis:Exit", false);
|
||||
|
||||
|
||||
if (config.get<bool>("Poly:Debug:GLVis:C_gf_View:View", defaultView)) {
|
||||
Probe::glVisView(CTmp, *C.FESpace(), "CTmp", config.get<std::string>("Poly:Debug:C_gf_View:Keyset", defaultKeyset));
|
||||
if (config.get<bool>("Poly:Debug:GLVis:C_gf_View:Exit", defaultExit)) {
|
||||
std::raise(SIGINT);
|
||||
}
|
||||
}
|
||||
|
||||
if (config.get<bool>("Poly:Debug:GLVis:F_gf_View:View", defaultView)) {
|
||||
Probe::glVisView(lambda_C, *nfl.FESpace(), "lambda_C", config.get<std::string>("Poly:Debug:F_gf_View:Keyset", defaultKeyset));
|
||||
if (config.get<bool>("Poly:Debug:GLVis:F_gf_View:Exit", defaultExit)) {
|
||||
std::raise(SIGINT);
|
||||
}
|
||||
}
|
||||
|
||||
if (config.get<bool>("Poly:Debug:GLVis:M_gf_View:View", defaultView)) {
|
||||
Probe::glVisView(y, *nfl.FESpace(), "y", config.get<std::string>("Poly:Debug:M_gf_View:Keyset", defaultKeyset));
|
||||
if (config.get<bool>("Poly:Debug:GLVis:M_gf_View:Exit", defaultExit)) {
|
||||
std::raise(SIGINT);
|
||||
}
|
||||
}
|
||||
|
||||
// Compute the constraint residual (C(u))
|
||||
y[lambdaDofOffset] = C_u[0] - C_val; // Enforce the constraint C(u) = C_val
|
||||
}
|
||||
|
||||
mfem::Operator &AugmentedOperator::GetGradient(const mfem::Vector &x) const {
|
||||
@@ -305,6 +339,8 @@ namespace polyMFEMUtils {
|
||||
J_aug->Set(i, lambdaDofOffset, -CVec[i]);
|
||||
}
|
||||
|
||||
J_aug->Set(lambdaDofOffset, lambdaDofOffset, 0.0); // The bottom right corner is zero
|
||||
|
||||
J_aug->Finalize();
|
||||
|
||||
delete lastJacobian;
|
||||
@@ -315,4 +351,24 @@ namespace polyMFEMUtils {
|
||||
AugmentedOperator::~AugmentedOperator() {
|
||||
delete lastJacobian;
|
||||
}
|
||||
|
||||
double calculateGaussianIntegral(mfem::Mesh &mesh, polyMFEMUtils::GaussianCoefficient &gaussianCoeff) {
|
||||
// Use a discontinuous L2 finite element space (order 0) for integrating the Gaussian.
|
||||
// We use L2 because the Gaussian is not necessarily continuous across element boundaries
|
||||
// if the Gaussian is narrow relative to the element size.
|
||||
mfem::L2_FECollection feCollection(0, mesh.SpaceDimension());
|
||||
mfem::FiniteElementSpace feSpace(&mesh, &feCollection);
|
||||
|
||||
mfem::LinearForm gaussianIntegral(&feSpace);
|
||||
gaussianIntegral.AddDomainIntegrator(new mfem::DomainLFIntegrator(gaussianCoeff));
|
||||
gaussianIntegral.Assemble();
|
||||
|
||||
// Create a GridFunction with a constant value of 1.0 on the L2 space.
|
||||
mfem::GridFunction one_gf(&feSpace);
|
||||
one_gf = 1.0;
|
||||
|
||||
// Evaluate the linear form on the constant GridFunction. This gives the integral.
|
||||
return gaussianIntegral(one_gf);
|
||||
|
||||
}
|
||||
} // namespace polyMFEMUtils
|
||||
Reference in New Issue
Block a user