Files
SERiF/src/poly/utils/public/polyMFEMUtils.h

190 lines
6.5 KiB
C++

/* ***********************************************************************
//
// Copyright (C) 2025 -- The 4D-STAR Collaboration
// File Author: Emily Boudreaux
// Last Modified: March 19, 2025
//
// 4DSSE is free software; you can use it and/or modify
// it under the terms and restrictions the GNU General Library Public
// License version 3 (GPLv3) as published by the Free Software Foundation.
//
// 4DSSE is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
// See the GNU Library General Public License for more details.
//
// You should have received a copy of the GNU Library General Public License
// along with this software; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
//
// *********************************************************************** */
#ifndef POLYMFEMUTILS_H
#define POLYMFEMUTILS_H
#include "mfem.hpp"
#include <string>
#include <vector>
#include "config.h"
#include "probe.h"
#include <unordered_map>
/**
* @file polyMFEMUtils.h
* @brief A collection of utilities for working with MFEM and solving the lane-emden equation.
*/
/**
* @namespace polyMFEMUtils
* @brief A namespace for utilities for working with MFEM and solving the lane-emden equation.
*/
namespace polyMFEMUtils {
/**
* @brief A class for nonlinear power integrator.
*/
class NonlinearPowerIntegrator: public mfem::NonlinearFormIntegrator {
private:
Config& config = Config::getInstance();
mfem::Coefficient &coeff_;
double polytropicIndex;
public:
/**
* @brief Constructor for NonlinearPowerIntegrator.
*
* @param coeff The function coefficient.
* @param n The polytropic index.
*/
NonlinearPowerIntegrator(mfem::Coefficient &coeff, double n);
/**
* @brief Assembles the element vector.
*
* @param el The finite element.
* @param Trans The element transformation.
* @param elfun The element function.
* @param elvect The element vector to be assembled.
*/
virtual void AssembleElementVector(const mfem::FiniteElement &el, mfem::ElementTransformation &Trans, const mfem::Vector &elfun, mfem::Vector &elvect) override;
/**
* @brief Assembles the element gradient.
*
* @param el The finite element.
* @param Trans The element transformation.
* @param elfun The element function.
* @param elmat The element matrix to be assembled.
*/
virtual void AssembleElementGrad (const mfem::FiniteElement &el, mfem::ElementTransformation &Trans, const mfem::Vector &elfun, mfem::DenseMatrix &elmat) override;
};
/**
* @brief A wrapper class for bilinear integrator.
*/
class BilinearIntegratorWrapper : public mfem::NonlinearFormIntegrator
{
private:
mfem::BilinearFormIntegrator *integrator;
public:
/**
* @brief Constructor for BilinearIntegratorWrapper.
*
* @param integratorInput The bilinear form integrator input.
*/
BilinearIntegratorWrapper(mfem::BilinearFormIntegrator *integratorInput);
/**
* @brief Destructor for BilinearIntegratorWrapper.
*/
virtual ~BilinearIntegratorWrapper();
/**
* @brief Assembles the element vector.
*
* @param el The finite element.
* @param Trans The element transformation.
* @param elfun The element function.
* @param elvect The element vector to be assembled.
*/
virtual void AssembleElementVector(const mfem::FiniteElement &el, mfem::ElementTransformation &Trans, const mfem::Vector &elfun, mfem::Vector &elvect) override;
/**
* @brief Assembles the element gradient.
*
* @param el The finite element.
* @param Trans The element transformation.
* @param elfun The element function.
* @param elmat The element matrix to be assembled.
*/
virtual void AssembleElementGrad(const mfem::FiniteElement &el, mfem::ElementTransformation &Trans, const mfem::Vector &elfun, mfem::DenseMatrix &elmat) override;
};
/**
* @brief A class for composite nonlinear integrator.
*/
class CompositeNonlinearIntegrator: public mfem::NonlinearFormIntegrator {
private:
std::vector<mfem::NonlinearFormIntegrator*> integrators;
public:
/**
* @brief Constructor for CompositeNonlinearIntegrator.
*/
CompositeNonlinearIntegrator();
/**
* @brief Destructor for CompositeNonlinearIntegrator.
*/
virtual ~CompositeNonlinearIntegrator();
/**
* @brief Adds an integrator to the composite integrator.
*
* @param integrator The nonlinear form integrator to add.
*/
void add_integrator(mfem::NonlinearFormIntegrator *integrator);
/**
* @brief Assembles the element vector.
*
* @param el The finite element.
* @param Trans The element transformation.
* @param elfun The element function.
* @param elvect The element vector to be assembled.
*/
virtual void AssembleElementVector(const mfem::FiniteElement &el, mfem::ElementTransformation &Trans, const mfem::Vector &elfun, mfem::Vector &elvect) override;
/**
* @brief Assembles the element gradient.
*
* @param el The finite element.
* @param Trans The element transformation.
* @param elfun The element function.
* @param elmat The element matrix to be assembled.
*/
virtual void AssembleElementGrad(const mfem::FiniteElement &el, mfem::ElementTransformation &Trans, const mfem::Vector &elfun, mfem::DenseMatrix &elmat) override;
};
class ConstraintIntegrator: public mfem::BilinearFormIntegrator {
private:
Config& m_config = Config::getInstance();
Probe::LogManager& logManager = Probe::LogManager::getInstance();
quill::Logger* m_logger = logManager.getLogger("log");
const double m_gamma;
mfem::Array<int> m_originElementIDs;
mfem::Array<mfem::IntegrationPoint> m_originIntegrationPoints;
mfem::DenseMatrix m_originCoordinateMatrix;
mfem::Vector m_originCoordinateVector;
mfem::Mesh* m_mesh;
mfem::Array<int> m_originConnectedElementIds;
public:
ConstraintIntegrator(const double gamma, mfem::Mesh* mesh);
~ConstraintIntegrator() = default;
void AssembleElementMatrix(const mfem::FiniteElement &el, mfem::ElementTransformation &Trans, mfem::DenseMatrix &elmat) override;
};
} // namespace polyMFEMUtils
#endif // POLYMFEMUTILS_H