fix(poly): bug fixing in block form
currently derivitive constraint is not working
This commit is contained in:
@@ -34,6 +34,8 @@
|
||||
#include "resourceManagerTypes.h"
|
||||
#include "operator.h"
|
||||
|
||||
#include "debug.h"
|
||||
|
||||
#include "quill/LogMacros.h"
|
||||
|
||||
|
||||
@@ -149,14 +151,14 @@ void PolySolver::assembleBlockSystem() {
|
||||
// A full derivation of the weak form can be found in the 4DSSE documentation
|
||||
|
||||
// --- Assemble the MixedBilinear and Bilinear forms (M, D, and Q) ---
|
||||
auto Mform = std::make_unique<mfem::MixedBilinearForm>(m_feTheta.get(), m_fePhi.get());
|
||||
auto Qform = std::make_unique<mfem::MixedBilinearForm>(m_fePhi.get(), m_feTheta.get());
|
||||
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(negOneCoeff));
|
||||
Qform->AddDomainIntegrator(new mfem::MixedVectorGradientIntegrator(negOneVCoeff));
|
||||
Dform->AddDomainIntegrator(new mfem::VectorFEMassIntegrator(oneCoeff));
|
||||
Mform->AddDomainIntegrator(new mfem::MixedVectorWeakDivergenceIntegrator());
|
||||
Qform->AddDomainIntegrator(new mfem::MixedVectorGradientIntegrator());
|
||||
Dform->AddDomainIntegrator(new mfem::VectorFEMassIntegrator());
|
||||
|
||||
Mform->Assemble();
|
||||
Mform->Finalize();
|
||||
@@ -187,47 +189,24 @@ void PolySolver::solve(){
|
||||
// --- Set the initial guess for the solution ---
|
||||
setInitialGuess();
|
||||
|
||||
// --- Set the essential true dofs for the operator ---
|
||||
mfem::Array<int> theta_ess_tdof_list, phi_ess_tdof_list;
|
||||
std::tie(theta_ess_tdof_list, phi_ess_tdof_list) = getEssentialTrueDof();
|
||||
m_polytropOperator->SetEssentialTrueDofs(theta_ess_tdof_list, phi_ess_tdof_list);
|
||||
setupOperator();
|
||||
|
||||
// --- Load configuration parameters ---
|
||||
double newtonRelTol = m_config.get<double>("Poly:Solver:Newton:RelTol", 1e-7);
|
||||
double newtonAbsTol = m_config.get<double>("Poly:Solver:Newton:AbsTol", 1e-7);
|
||||
int newtonMaxIter = m_config.get<int>("Poly:Solver:Newton:MaxIter", 200);
|
||||
int newtonPrintLevel = m_config.get<int>("Poly:Solver:Newton:PrintLevel", 1);
|
||||
// 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);
|
||||
state_vector.GetBlock(0) = *m_theta;
|
||||
state_vector.GetBlock(1) = *m_phi;
|
||||
|
||||
double gmresRelTol = m_config.get<double>("Poly:Solver:GMRES:RelTol", 1e-10);
|
||||
double gmresAbsTol = m_config.get<double>("Poly:Solver:GMRES:AbsTol", 1e-12);
|
||||
int gmresMaxIter = m_config.get<int>("Poly:Solver:GMRES:MaxIter", 2000);
|
||||
int gmresPrintLevel = m_config.get<int>("Poly:Solver:GMRES:PrintLevel", 0);
|
||||
mfem::Vector zero_rhs(block_offsets.Last());
|
||||
zero_rhs = 0.0;
|
||||
|
||||
LOG_DEBUG(m_logger, "Newton Solver (relTol: {:0.2E}, absTol: {:0.2E}, maxIter: {}, printLevel: {})", newtonRelTol, newtonAbsTol, newtonMaxIter, newtonPrintLevel);
|
||||
LOG_DEBUG(m_logger, "GMRES Solver (relTol: {:0.2E}, absTol: {:0.2E}, maxIter: {}, printLevel: {})", gmresRelTol, gmresAbsTol, gmresMaxIter, gmresPrintLevel);
|
||||
|
||||
// --- Set up the Newton solver ---
|
||||
mfem::NewtonSolver newtonSolver;
|
||||
newtonSolver.SetRelTol(newtonRelTol);
|
||||
newtonSolver.SetAbsTol(newtonAbsTol);
|
||||
newtonSolver.SetMaxIter(newtonMaxIter);
|
||||
newtonSolver.SetPrintLevel(newtonPrintLevel);
|
||||
newtonSolver.SetOperator(*m_polytropOperator);
|
||||
mfem::GMRESSolver gmresSolver;
|
||||
gmresSolver.SetRelTol(gmresRelTol);
|
||||
gmresSolver.SetAbsTol(gmresAbsTol);
|
||||
gmresSolver.SetMaxIter(gmresMaxIter);
|
||||
gmresSolver.SetPrintLevel(gmresPrintLevel);
|
||||
newtonSolver.SetSolver(gmresSolver);
|
||||
// newtonSolver.SetAdaptiveLinRtol();
|
||||
mfem::NewtonSolver newtonSolver = setupNewtonSolver();
|
||||
|
||||
mfem::Vector B(m_feTheta->GetTrueVSize());
|
||||
B = 0.0;
|
||||
|
||||
newtonSolver.Mult(B, *m_theta);
|
||||
newtonSolver.Mult(zero_rhs, state_vector);
|
||||
|
||||
// --- Save and view the solution ---
|
||||
saveAndViewSolution();
|
||||
saveAndViewSolution(state_vector);
|
||||
|
||||
|
||||
}
|
||||
@@ -295,22 +274,77 @@ void PolySolver::setInitialGuess() {
|
||||
|
||||
}
|
||||
|
||||
void PolySolver::saveAndViewSolution() {
|
||||
void PolySolver::saveAndViewSolution(const mfem::BlockVector& state_vector) {
|
||||
mfem::BlockVector x_block(const_cast<mfem::BlockVector&>(state_vector), m_polytropOperator->GetBlockOffsets());
|
||||
mfem::Vector& x_theta = x_block.GetBlock(0);
|
||||
mfem::Vector& x_phi = x_block.GetBlock(1);
|
||||
|
||||
bool doView = m_config.get<bool>("Poly:Output:View", false);
|
||||
if (doView) {
|
||||
Probe::glVisView(*m_theta, *m_mesh, "solution");
|
||||
Probe::glVisView(x_theta, *m_feTheta, "θ Solution");
|
||||
Probe::glVisView(x_phi, *m_fePhi, "ɸ Solution");
|
||||
}
|
||||
|
||||
// --- Extract the Solution ---
|
||||
bool write11DSolution = m_config.get<bool>("Poly:Output:1D:Save", true);
|
||||
if (write11DSolution) {
|
||||
std::string solutionPath = m_config.get<std::string>("Poly:Output:1D:Path", "polytropeSolution_1D.csv");
|
||||
std::string derivSolPath = "d" + solutionPath;
|
||||
|
||||
double rayCoLatitude = m_config.get<double>("Poly:Output:1D:RayCoLatitude", 0.0);
|
||||
double rayLongitude = m_config.get<double>("Poly:Output:1D:RayLongitude", 0.0);
|
||||
int raySamples = m_config.get<int>("Poly:Output:1D:RaySamples", 100);
|
||||
|
||||
std::vector rayDirection = {rayCoLatitude, rayLongitude};
|
||||
|
||||
Probe::getRaySolution(*m_theta, *m_feTheta, rayDirection, raySamples, solutionPath);
|
||||
Probe::getRaySolution(x_theta, *m_feTheta, rayDirection, raySamples, solutionPath);
|
||||
// Probe::getRaySolution(x_phi, *m_fePhi, rayDirection, raySamples, derivSolPath);
|
||||
}
|
||||
}
|
||||
|
||||
void PolySolver::setupOperator() {
|
||||
mfem::Array<int> theta_ess_tdof_list, phi_ess_tdof_list;
|
||||
std::tie(theta_ess_tdof_list, phi_ess_tdof_list) = getEssentialTrueDof();
|
||||
m_polytropOperator->SetEssentialTrueDofs(theta_ess_tdof_list, phi_ess_tdof_list);
|
||||
|
||||
// -- Finalize the operator --
|
||||
m_polytropOperator->finalize();
|
||||
|
||||
if (!m_polytropOperator->isFinalized()) {
|
||||
LOG_ERROR(m_logger, "PolytropeOperator is not finalized. Cannot solve.");
|
||||
throw std::runtime_error("PolytropeOperator is not finalized. Cannot solve.");
|
||||
}
|
||||
}
|
||||
|
||||
mfem::NewtonSolver PolySolver::setupNewtonSolver(){
|
||||
// --- Load configuration parameters ---
|
||||
double newtonRelTol = m_config.get<double>("Poly:Solver:Newton:RelTol", 1e-7);
|
||||
double newtonAbsTol = m_config.get<double>("Poly:Solver:Newton:AbsTol", 1e-7);
|
||||
int newtonMaxIter = m_config.get<int>("Poly:Solver:Newton:MaxIter", 200);
|
||||
int newtonPrintLevel = m_config.get<int>("Poly:Solver:Newton:PrintLevel", 1);
|
||||
|
||||
double gmresRelTol = m_config.get<double>("Poly:Solver:GMRES:RelTol", 1e-10);
|
||||
double gmresAbsTol = m_config.get<double>("Poly:Solver:GMRES:AbsTol", 1e-12);
|
||||
int gmresMaxIter = m_config.get<int>("Poly:Solver:GMRES:MaxIter", 2000);
|
||||
int gmresPrintLevel = m_config.get<int>("Poly:Solver:GMRES:PrintLevel", 0);
|
||||
|
||||
LOG_DEBUG(m_logger, "Newton Solver (relTol: {:0.2E}, absTol: {:0.2E}, maxIter: {}, printLevel: {})", newtonRelTol, newtonAbsTol, newtonMaxIter, newtonPrintLevel);
|
||||
LOG_DEBUG(m_logger, "GMRES Solver (relTol: {:0.2E}, absTol: {:0.2E}, maxIter: {}, printLevel: {})", gmresRelTol, gmresAbsTol, gmresMaxIter, gmresPrintLevel);
|
||||
|
||||
// --- Set up the Newton solver ---
|
||||
mfem::NewtonSolver newtonSolver;
|
||||
newtonSolver.SetRelTol(newtonRelTol);
|
||||
newtonSolver.SetAbsTol(newtonAbsTol);
|
||||
newtonSolver.SetMaxIter(newtonMaxIter);
|
||||
newtonSolver.SetPrintLevel(newtonPrintLevel);
|
||||
newtonSolver.SetOperator(*m_polytropOperator);
|
||||
mfem::GMRESSolver gmresSolver;
|
||||
gmresSolver.SetRelTol(gmresRelTol);
|
||||
gmresSolver.SetAbsTol(gmresAbsTol);
|
||||
gmresSolver.SetMaxIter(gmresMaxIter);
|
||||
gmresSolver.SetPrintLevel(gmresPrintLevel);
|
||||
newtonSolver.SetSolver(gmresSolver);
|
||||
// newtonSolver.SetAdaptiveLinRtol();
|
||||
|
||||
return newtonSolver;
|
||||
}
|
||||
@@ -1,6 +1,7 @@
|
||||
#ifndef POLYSOLVER_H
|
||||
#define POLYSOLVER_H
|
||||
|
||||
#include "linalg/solvers.hpp"
|
||||
#include "mfem.hpp"
|
||||
#include <memory>
|
||||
#include <utility>
|
||||
@@ -53,7 +54,9 @@ private: // Private methods
|
||||
std::pair<mfem::Array<int>, mfem::Array<int>> getEssentialTrueDof();
|
||||
mfem::Array<int> findCenterElement();
|
||||
void setInitialGuess();
|
||||
void saveAndViewSolution();
|
||||
void saveAndViewSolution(const mfem::BlockVector& state_vector);
|
||||
mfem::NewtonSolver setupNewtonSolver();
|
||||
void setupOperator();
|
||||
|
||||
|
||||
};
|
||||
|
||||
@@ -60,7 +60,8 @@ namespace polyMFEMUtils {
|
||||
double u_safe = std::max(u_val, 0.0);
|
||||
double u_nl = std::pow(u_safe, m_polytropicIndex);
|
||||
|
||||
double coeff_val = m_coeff.Eval(Trans, ip);
|
||||
// double coeff_val = m_coeff.Eval(Trans, ip);
|
||||
double coeff_val = 1.0;
|
||||
double x2_u_nl = coeff_val * u_nl;
|
||||
|
||||
for (int i = 0; i < dof; i++){
|
||||
@@ -94,7 +95,8 @@ namespace polyMFEMUtils {
|
||||
for (int j = 0; j < dof; j++) {
|
||||
u_val += elfun(j) * shape(j);
|
||||
}
|
||||
double coeff_val = m_coeff.Eval(Trans, ip);
|
||||
// double coeff_val = m_coeff.Eval(Trans, ip);
|
||||
double coeff_val = 1.0;
|
||||
|
||||
|
||||
// Calculate the Jacobian
|
||||
|
||||
@@ -3,6 +3,8 @@
|
||||
#include "linalg/vector.hpp"
|
||||
#include <memory>
|
||||
|
||||
#include "debug.h"
|
||||
|
||||
PolytropeOperator::PolytropeOperator(
|
||||
|
||||
std::unique_ptr<mfem::MixedBilinearForm> M,
|
||||
@@ -19,8 +21,12 @@ PolytropeOperator::PolytropeOperator(
|
||||
m_Q = std::move(Q);
|
||||
m_D = std::move(D);
|
||||
m_f = std::move(f);
|
||||
}
|
||||
|
||||
|
||||
void PolytropeOperator::finalize() {
|
||||
if (m_isFinalized) {
|
||||
return;
|
||||
}
|
||||
m_Mmat = std::make_unique<mfem::SparseMatrix>(m_M->SpMat());
|
||||
m_Qmat = std::make_unique<mfem::SparseMatrix>(m_Q->SpMat());
|
||||
m_Dmat = std::make_unique<mfem::SparseMatrix>(m_D->SpMat());
|
||||
@@ -28,14 +34,14 @@ PolytropeOperator::PolytropeOperator(
|
||||
m_negM_op = std::make_unique<mfem::ScaledOperator>(m_Mmat.get(), -1.0);
|
||||
m_negQ_op = std::make_unique<mfem::ScaledOperator>(m_Qmat.get(), -1.0);
|
||||
|
||||
MFEM_ASSERT(m_Mmat.get() != nullptr, "Matrix m_Mmat is null in PolytropeOperator constructor");
|
||||
MFEM_ASSERT(m_Qmat.get() != nullptr, "Matrix m_Qmat is null in PolytropeOperator constructor");
|
||||
MFEM_ASSERT(m_Dmat.get() != nullptr, "Matrix m_Dmat is null in PolytropeOperator constructor");
|
||||
MFEM_ASSERT(m_f.get() != nullptr, "NonlinearForm m_f is null in PolytropeOperator constructor");
|
||||
m_isFinalized = true;
|
||||
}
|
||||
|
||||
|
||||
void PolytropeOperator::Mult(const mfem::Vector &x, mfem::Vector &y) const {
|
||||
if (!m_isFinalized) {
|
||||
MFEM_ABORT("PolytropeOperator::Mult called before finalize");
|
||||
}
|
||||
// -- Create BlockVector views for input x and output y
|
||||
mfem::BlockVector x_block(const_cast<mfem::Vector&>(x), m_blockOffsets);
|
||||
mfem::BlockVector y_block(y, m_blockOffsets);
|
||||
@@ -89,6 +95,9 @@ void PolytropeOperator::Mult(const mfem::Vector &x, mfem::Vector &y) const {
|
||||
|
||||
}
|
||||
mfem::Operator& PolytropeOperator::GetGradient(const mfem::Vector &x) const {
|
||||
if (!m_isFinalized) {
|
||||
MFEM_ABORT("PolytropeOperator::GetGradient called before finalize");
|
||||
}
|
||||
// -- Get the gradient of f --
|
||||
mfem::BlockVector x_block(const_cast<mfem::Vector&>(x), m_blockOffsets);
|
||||
const mfem::Vector& x_theta = x_block.GetBlock(0);
|
||||
@@ -111,6 +120,7 @@ mfem::Operator& PolytropeOperator::GetGradient(const mfem::Vector &x) const {
|
||||
|
||||
void PolytropeOperator::SetEssentialTrueDofs(const mfem::Array<int> &theta_ess_tofs,
|
||||
const mfem::Array<int> &phi_ess_tofs) {
|
||||
m_isFinalized = false;
|
||||
m_theta_ess_tofs = theta_ess_tofs;
|
||||
m_phi_ess_tofs = phi_ess_tofs;
|
||||
|
||||
|
||||
@@ -20,13 +20,19 @@ public:
|
||||
void SetEssentialTrueDofs(const mfem::Array<int> &theta_ess_tofs,
|
||||
const mfem::Array<int> &phi_ess_tofs);
|
||||
|
||||
bool isFinalized() const { return m_isFinalized; }
|
||||
|
||||
void finalize();
|
||||
|
||||
const mfem::Array<int>& GetBlockOffsets() const { return m_blockOffsets; }
|
||||
|
||||
private:
|
||||
std::unique_ptr<mfem::MixedBilinearForm> m_M;
|
||||
std::unique_ptr<mfem::MixedBilinearForm> m_Q;
|
||||
std::unique_ptr<mfem::BilinearForm> m_D;
|
||||
std::unique_ptr<mfem::NonlinearForm> m_f;
|
||||
|
||||
const mfem::Array<int> &m_blockOffsets;
|
||||
const mfem::Array<int> m_blockOffsets;
|
||||
|
||||
mfem::Array<int> m_theta_ess_tofs;
|
||||
mfem::Array<int> m_phi_ess_tofs;
|
||||
@@ -39,6 +45,8 @@ private:
|
||||
std::unique_ptr<mfem::ScaledOperator> m_negM_op;
|
||||
std::unique_ptr<mfem::ScaledOperator> m_negQ_op;
|
||||
mutable std::unique_ptr<mfem::BlockOperator> m_jacobian;
|
||||
|
||||
bool m_isFinalized = false;
|
||||
};
|
||||
|
||||
|
||||
|
||||
Reference in New Issue
Block a user