Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fix the finite element dual assignment operators #1032

Merged
merged 6 commits into from
Dec 20, 2023
Merged
Show file tree
Hide file tree
Changes from 2 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
7 changes: 2 additions & 5 deletions src/serac/physics/solid_mechanics.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -1105,11 +1105,8 @@ class SolidMechanics<order, dim, Parameters<parameter_space...>, std::integer_se
// the material state buffers to be updated
residual_->update_qdata = true;

// this seems like the wrong way to be doing this assignment, but
// reactions_ = residual(displacement, ...);
// isn't currently supported
reactions_.Vector::operator=(
(*residual_)(displacement_, acceleration_, shape_displacement_, *parameters_[parameter_indices].state...));
reactions_ =
(*residual_)(displacement_, acceleration_, shape_displacement_, *parameters_[parameter_indices].state...);

residual_->update_qdata = false;
}
Expand Down
53 changes: 52 additions & 1 deletion src/serac/physics/state/finite_element_dual.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -28,7 +28,6 @@ namespace serac {
class FiniteElementDual : public FiniteElementVector {
public:
using FiniteElementVector::FiniteElementVector;
using FiniteElementVector::operator=;
using mfem::Vector::Print;

/**
Expand Down Expand Up @@ -57,6 +56,58 @@ class FiniteElementDual : public FiniteElementVector {
return *this;
}

/**
* @brief Move assignment
*
* @param rhs The right hand side input Dual
* @return The assigned FiniteElementDual
*/
FiniteElementDual& operator=(FiniteElementDual&& rhs)
{
std::unique_ptr<mfem::ParLinearForm> lf(std::move(rhs.linear_form_));

FiniteElementVector::operator=(rhs);

this->linear_form_ = std::move(lf);
return *this;
}

/**
* @brief Copy assignment with HypreParVector
*
* @param rhs The right hand side input HypreParVector
* @return The assigned FiniteElementDual
*/
FiniteElementDual& operator=(const mfem::HypreParVector& rhs)
{
FiniteElementVector::operator=(rhs);
return *this;
}

/**
* @brief Copy assignment with mfem::Vector
*
* @param rhs The right hand side input Dual
* @return The assigned FiniteElementDual
*/
FiniteElementDual& operator=(const mfem::Vector& rhs)
{
FiniteElementVector::operator=(rhs);
return *this;
}

/**
* @brief Copy assignment with double
*
* @param rhs The right hand side input double
* @return The assigned FiniteElementDual
*/
FiniteElementDual& operator=(double rhs)
{
FiniteElementVector::operator=(rhs);
return *this;
}

/**
* @brief Fill a user-provided linear form based on the underlying true vector
*
Expand Down
53 changes: 50 additions & 3 deletions src/serac/physics/state/finite_element_state.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -48,7 +48,6 @@ inline bool is_vector_valued(const GeneralCoefficient& coef)
class FiniteElementState : public FiniteElementVector {
public:
using FiniteElementVector::FiniteElementVector;
using FiniteElementVector::operator=;
using mfem::Vector::Print;

/**
Expand Down Expand Up @@ -77,13 +76,61 @@ class FiniteElementState : public FiniteElementVector {
return *this;
}

/**
* @brief Move assignment
*
* @param rhs The right hand side input State
* @return The assigned FiniteElementState
*/
FiniteElementState& operator=(FiniteElementState&& rhs)
{
FiniteElementVector::operator=(rhs);
return *this;
}

/**
* @brief Copy assignment with HypreParVector
*
* @param rhs The right hand side input HypreParVector
* @return The assigned FiniteElementState
*/
FiniteElementState& operator=(const mfem::HypreParVector& rhs)
{
FiniteElementVector::operator=(rhs);
return *this;
}

/**
* @brief Copy assignment with mfem::Vector
*
* @param rhs The right hand side input State
* @return The assigned FiniteElementState
*/
FiniteElementState& operator=(const mfem::Vector& rhs)
{
FiniteElementVector::operator=(rhs);
return *this;
}

/**
* @brief Copy assignment with double
*
* @param rhs The right hand side input double
* @return The assigned FiniteElementState
*/
FiniteElementState& operator=(double rhs)
{
FiniteElementVector::operator=(rhs);
return *this;
}

/**
* @brief Fill a user-provided grid function based on the underlying true vector
*
* This distributes true vector dofs to the finite element (local) dofs by multiplying the true dofs
* by the prolongation operator.
*
* @see <a href="https://mfem.org/pri-dual-vec/">MFEM documentation</a> for details
* @see <a href="https://mfem.org/pri-State-vec/">MFEM documentation</a> for details
jamiebramwell marked this conversation as resolved.
Show resolved Hide resolved
*
*/
void fillGridFunction(mfem::ParGridFunction& grid_function) const { grid_function.SetFromTrueDofs(*this); }
Expand All @@ -94,7 +141,7 @@ class FiniteElementState : public FiniteElementVector {
* This distributes the grid function dofs to the true vector dofs by multiplying by the
* restriction operator.
*
* @see <a href="https://mfem.org/pri-dual-vec/">MFEM documentation</a> for details
* @see <a href="https://mfem.org/pri-State-vec/">MFEM documentation</a> for details
jamiebramwell marked this conversation as resolved.
Show resolved Hide resolved
*
* @param grid_function The grid function used to initialize the underlying true vector.
*/
Expand Down
6 changes: 6 additions & 0 deletions src/serac/physics/state/finite_element_vector.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -50,6 +50,12 @@ FiniteElementVector& FiniteElementVector::operator=(const mfem::HypreParVector&
return *this;
}

FiniteElementVector& FiniteElementVector::operator=(const mfem::Vector& rhs)
{
Vector::operator=(rhs);
return *this;
}

FiniteElementVector& FiniteElementVector::operator=(const FiniteElementVector& rhs)
{
SLIC_ERROR_IF(Size() != rhs.Size(),
Expand Down
8 changes: 8 additions & 0 deletions src/serac/physics/state/finite_element_vector.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -142,6 +142,14 @@ class FiniteElementVector : public mfem::HypreParVector {
*/
FiniteElementVector& operator=(const mfem::HypreParVector& rhs);

/**
* @brief Copy assignment from a hypre par vector
*
* @param rhs The rhs input hypre par vector
* @return The copy assigned input vector
*/
FiniteElementVector& operator=(const mfem::Vector& rhs);

/**
* @brief Returns the MPI communicator for the state
* @return The underlying MPI communicator
Expand Down