fix(poly): fixed -M bug in form
MFEM MixedVectorWeakDivergenceIntegrator is actually already -M in our derivation, I have negated this so that Mform -> M directly
This commit is contained in:
@@ -106,74 +106,65 @@ PolySolver::PolySolver(const double n, const double order) {
|
||||
PolySolver::~PolySolver() = default;
|
||||
|
||||
void PolySolver::assembleBlockSystem() {
|
||||
mfem::Array<int> blockOffsets = computeBlockOffsets();
|
||||
|
||||
// Start by defining the block structure of the system
|
||||
// Block 0: Theta (scalar space, uses m_feTheta)
|
||||
// Block 1: Phi (vector space, uses m_fePhi)
|
||||
mfem::Array<mfem::FiniteElementSpace*> feSpaces;
|
||||
feSpaces.Append(m_feTheta.get());
|
||||
feSpaces.Append(m_fePhi.get());
|
||||
const std::unique_ptr<formBundle> forms = buildIndividualForms(blockOffsets);
|
||||
|
||||
// Create the block offsets. These define the start of each block in the combined vector.
|
||||
// Block offsets will be [0, thetaDofs, thetaDofs + phiDofs]
|
||||
// --- Build the BlockOperator ---
|
||||
m_polytropOperator = std::make_unique<PolytropeOperator>(
|
||||
std::move(forms->M),
|
||||
std::move(forms->Q),
|
||||
std::move(forms->D),
|
||||
std::move(forms->f),
|
||||
blockOffsets);
|
||||
}
|
||||
|
||||
mfem::Array<int> PolySolver::computeBlockOffsets() const {
|
||||
mfem::Array<int> blockOffsets;
|
||||
blockOffsets.SetSize(3);
|
||||
blockOffsets[0] = 0;
|
||||
blockOffsets[1] = feSpaces[0]->GetVSize();
|
||||
blockOffsets[2] = feSpaces[1]->GetVSize();
|
||||
blockOffsets[1] = m_feTheta->GetVSize(); // Get actual number of dofs *before* applying BCs
|
||||
blockOffsets[2] = m_fePhi->GetVSize();
|
||||
blockOffsets.PartialSum(); // Cumulative sum to get the offsets
|
||||
return blockOffsets;
|
||||
}
|
||||
|
||||
// Add integrators to block form one by one
|
||||
// We add integrators corresponding to each term in the weak form
|
||||
// The block form of the residual matrix
|
||||
// ⎡ 0 -M ⎤ ⎡ θ ⎤ + ⎡f(θ)⎤ = ⎡ 0 ⎤ = R(X)
|
||||
// ⎣ -Q D ⎦ ⎣ Φ ⎦ ⎣ 0 ⎦ ⎣ 0 ⎦
|
||||
// This then simplifies to
|
||||
// ⎡f(θ) - MΘ ⎤ = ⎡ 0 ⎤ = R
|
||||
// ⎣ -Qɸ Dθ ⎦ ⎣ 0 ⎦
|
||||
// Here M, Q, and D are
|
||||
// M = ∫∇ψᶿ·Nᵠ dV (MixedVectorWeakDivergenceIntegrator)
|
||||
// D = ∫ψᵠ·Nᵠ dV (VectorFEMassIntegrator)
|
||||
// Q = ∫ψᵠ·∇Nᶿ dV (MixedVectorGradientIntegrator)
|
||||
// f(θ) = ∫ψᶿ·θⁿ dV (NonlinearPowerIntegrator)
|
||||
// Here ψᶿ and ψᵠ are the test functions for the theta and phi spaces, respectively
|
||||
// Nᵠ and Nᶿ are the basis functions for the theta and phi spaces, respectively
|
||||
// A full derivation of the weak form can be found in the 4DSSE documentation
|
||||
|
||||
std::unique_ptr<formBundle> PolySolver::buildIndividualForms(const mfem::Array<int> &blockOffsets) {
|
||||
// --- Assemble the MixedBilinear and Bilinear forms (M, D, and Q) ---
|
||||
auto Mform = std::make_unique<mfem::MixedBilinearForm>(m_fePhi.get(), m_feTheta.get());
|
||||
auto Qform = std::make_unique<mfem::MixedBilinearForm>(m_feTheta.get(), m_fePhi.get());
|
||||
auto Dform = std::make_unique<mfem::BilinearForm>(m_fePhi.get());
|
||||
|
||||
// TODO: Check the sign on all of the integrators
|
||||
Mform->AddDomainIntegrator(new mfem::MixedVectorWeakDivergenceIntegrator());
|
||||
Qform->AddDomainIntegrator(new mfem::MixedVectorGradientIntegrator());
|
||||
Dform->AddDomainIntegrator(new mfem::VectorFEMassIntegrator());
|
||||
auto forms = std::make_unique<formBundle>(
|
||||
std::make_unique<mfem::MixedBilinearForm>(m_fePhi.get(), m_feTheta.get()),
|
||||
std::make_unique<mfem::MixedBilinearForm>(m_feTheta.get(), m_fePhi.get()),
|
||||
std::make_unique<mfem::BilinearForm>(m_fePhi.get()),
|
||||
std::make_unique<mfem::NonlinearForm>(m_feTheta.get())
|
||||
);
|
||||
|
||||
Mform->Assemble();
|
||||
Mform->Finalize();
|
||||
// --- -M negation -> M ---
|
||||
mfem::Vector negOneVec(m_fePhi->GetVDim());
|
||||
negOneVec = -1.0;
|
||||
|
||||
Qform->Assemble();
|
||||
Qform->Finalize();
|
||||
m_negationCoeff = std::make_unique<mfem::VectorConstantCoefficient>(negOneVec);
|
||||
|
||||
Dform->Assemble();
|
||||
Dform->Finalize();
|
||||
// --- Add the integrators to the forms ---
|
||||
forms->M->AddDomainIntegrator(new mfem::MixedVectorWeakDivergenceIntegrator(*m_negationCoeff));
|
||||
forms->Q->AddDomainIntegrator(new mfem::MixedVectorGradientIntegrator());
|
||||
forms->D->AddDomainIntegrator(new mfem::VectorFEMassIntegrator());
|
||||
|
||||
// --- Assemble the NonlinearForm (f) ---
|
||||
// Note that the nonlinear form is built here but its essential true dofs (boundary conditions) are
|
||||
// not set until later (when solve is called). They are applied through the operator rather than directly.
|
||||
auto fform = std::make_unique<mfem::NonlinearForm>(m_feTheta.get());
|
||||
fform->AddDomainIntegrator(new polyMFEMUtils::NonlinearPowerIntegrator(m_polytropicIndex));
|
||||
// --- Assemble and Finalize the forms ---
|
||||
assembleAndFinalizeForm(forms->M);
|
||||
assembleAndFinalizeForm(forms->Q);
|
||||
assembleAndFinalizeForm(forms->D);
|
||||
|
||||
// -- Build the BlockOperator --
|
||||
m_polytropOperator = std::make_unique<PolytropeOperator>(
|
||||
std::move(Mform),
|
||||
std::move(Qform),
|
||||
std::move(Dform),
|
||||
std::move(fform),
|
||||
blockOffsets
|
||||
);
|
||||
forms->f->AddDomainIntegrator(new polyMFEMUtils::NonlinearPowerIntegrator(m_polytropicIndex));
|
||||
|
||||
return std::move(forms);
|
||||
}
|
||||
|
||||
void PolySolver::assembleAndFinalizeForm(auto &f) {
|
||||
// This constructs / ensures the matrix representation for each form
|
||||
// Assemble => Computes the local element matrices across the domain. Adds these to the global matrix . Adds these to the global matrix.
|
||||
// Finalize => Builds the sparsity pattern and allows the SparseMatrix representation to be extracted.
|
||||
f->Assemble();
|
||||
f->Finalize();
|
||||
}
|
||||
|
||||
void PolySolver::solve() const {
|
||||
@@ -193,6 +184,9 @@ void PolySolver::solve() const {
|
||||
throw std::runtime_error("PolytropeOperator is not finalized. Cannot solve.");
|
||||
}
|
||||
|
||||
GetDofCoordinates(*m_feTheta, "dof2posTheta.csv");
|
||||
GetDofCoordinates(*m_fePhi, "dof2posPhi.csv");
|
||||
|
||||
// It's safer to get the offsets directly from the operator after finalization
|
||||
const mfem::Array<int>& block_offsets = m_polytropOperator->GetBlockOffsets(); // Assuming a getter exists or accessing member if public/friend
|
||||
mfem::BlockVector state_vector(block_offsets);
|
||||
@@ -213,6 +207,8 @@ void PolySolver::solve() const {
|
||||
// with updating the preconditioner at every newton step as the
|
||||
// changes to the jacobian are automatically propagated through the
|
||||
// solving chain. This is at least true with MFEM 4.8-rc0
|
||||
std::string windowTitle = "testWindow";
|
||||
std::string keyset = "";
|
||||
sb.newton.Mult(zero_rhs, state_vector);
|
||||
|
||||
// --- Save and view the solution ---
|
||||
@@ -279,7 +275,8 @@ void PolySolver::setInitialGuess() const {
|
||||
const double radius = Probe::getMeshRadius(*m_mesh);
|
||||
const double u = 1/radius;
|
||||
|
||||
return -std::pow((u*r), 2)+1.0; // The series expansion is a better guess; however, this is cheaper and ensures that the value at the surface is very close to zero in a way that the series expansion does not
|
||||
return (-1.0/radius) * r + 1;
|
||||
// return -std::pow((u*r), 2)+1.0; // The series expansion is a better guess; however, this is cheaper and ensures that the value at the surface is very close to zero in a way that the series expansion does not
|
||||
}
|
||||
);
|
||||
m_theta->ProjectCoefficient(thetaInitGuess);
|
||||
@@ -338,6 +335,43 @@ void PolySolver::LoadSolverUserParams(double &newtonRelTol, double &newtonAbsTol
|
||||
LOG_DEBUG(m_logger, "GMRES Solver (relTol: {:0.2E}, absTol: {:0.2E}, maxIter: {}, printLevel: {})", gmresRelTol, gmresAbsTol, gmresMaxIter, gmresPrintLevel);
|
||||
}
|
||||
|
||||
void PolySolver::GetDofCoordinates(mfem::FiniteElementSpace &fes, const std::string& filename) const {
|
||||
mfem::Mesh *mesh = fes.GetMesh();
|
||||
double r = Probe::getMeshRadius(*mesh);
|
||||
std::ofstream outputFile(filename, std::ios::out | std::ios::trunc);
|
||||
outputFile << "dof,R,r,x,y,z" << '\n';
|
||||
|
||||
const int nElements = mesh->GetNE();
|
||||
|
||||
mfem::Vector coords;
|
||||
mfem::IntegrationPoint ipZero;
|
||||
double p[3] = {0.0, 0.0, 0.0};
|
||||
int actual_idx;
|
||||
ipZero.Set3(p);
|
||||
for (int i = 0; i < nElements; i++) {
|
||||
mfem::Array<int> elemDofs;
|
||||
fes.GetElementDofs(i, elemDofs);
|
||||
mfem::ElementTransformation* T = mesh->GetElementTransformation(i);
|
||||
mfem::Vector physCoord(3);
|
||||
T->Transform(ipZero, physCoord);
|
||||
for (int dofID = 0; dofID < elemDofs.Size(); dofID++) {
|
||||
if (elemDofs[dofID] < 0) {
|
||||
actual_idx = -elemDofs[dofID] - 1;
|
||||
} else {
|
||||
actual_idx = elemDofs[dofID];
|
||||
}
|
||||
outputFile << actual_idx;
|
||||
if (dofID != elemDofs.Size() - 1) {
|
||||
outputFile << "|";
|
||||
} else {
|
||||
outputFile << ",";
|
||||
}
|
||||
}
|
||||
outputFile << r << "," << physCoord.Norml2() << "," << physCoord[0] << "," << physCoord[1] << "," << physCoord[2] << '\n';
|
||||
}
|
||||
outputFile.close();
|
||||
}
|
||||
|
||||
solverBundle PolySolver::setupNewtonSolver() const {
|
||||
// --- Load configuration parameters ---
|
||||
double newtonRelTol, newtonAbsTol, gmresRelTol, gmresAbsTol;
|
||||
|
||||
Reference in New Issue
Block a user