Skip to content

Commit

Permalink
update look up table
Browse files Browse the repository at this point in the history
  • Loading branch information
maxnezdyur committed May 16, 2024
1 parent 652634e commit a2f5738
Show file tree
Hide file tree
Showing 2 changed files with 38 additions and 34 deletions.
4 changes: 4 additions & 0 deletions framework/include/variables/MooseVariableBase.h
Original file line number Diff line number Diff line change
Expand Up @@ -188,6 +188,8 @@ class MooseVariableBase : public MooseObject,
*/
bool doDerivatives() const;

void setStride() const;

/// System this variable is part of
SystemBase & _sys;

Expand Down Expand Up @@ -244,6 +246,8 @@ class MooseVariableBase : public MooseObject,

/// Whether this variable lives on lower dimensional blocks
bool _is_lower_d;

mutable std::vector<dof_id_type> _stride;
};

inline void
Expand Down
68 changes: 34 additions & 34 deletions framework/src/variables/MooseVariableBase.C
Original file line number Diff line number Diff line change
Expand Up @@ -17,11 +17,14 @@
#include "InputParameterWarehouse.h"
#include "BlockRestrictable.h"

#include "libmesh/id_types.h"
#include "libmesh/int_range.h"
#include "libmesh/variable.h"
#include "libmesh/dof_map.h"
#include "libmesh/system.h"
#include "libmesh/fe_type.h"
#include "libmesh/string_to_enum.h"
#include <algorithm>

// Users should never actually create this object
registerMooseObject("MooseApp", MooseVariableBase);
Expand Down Expand Up @@ -156,7 +159,7 @@ MooseVariableBase::allDofIndices() const
if (it != _sys.subproblem()._var_dof_map.end())
return it->second;
else
mooseError("VariableAllDoFMap not prepared for ",
mooseError("VariableAllDOFMap not prepared for ",
name(),
" . Check nonlocal coupling requirement for the variable.");
}
Expand All @@ -171,51 +174,48 @@ std::vector<dof_id_type>
MooseVariableBase::componentDofIndices(const std::vector<dof_id_type> & dof_indices,
unsigned int component) const
{

// Initialize the new_dof_indices vector with the input dof_indices
std::vector<dof_id_type> new_dof_indices(dof_indices);

// Only perform mapping if the component is not 0
if (component != 0)
{
// Grab all DoFs for component 0; this will be used as a reference
std::vector<dof_id_type> comp_0_dofs;
_dof_map.local_variable_indices(comp_0_dofs, _mesh, _var_num);
// Set once. Need to have redo stride when mesh or dof map changes.
setStride();

// Determine the variable number for the specified component
auto comp_num = _var_num + component;

// Grab all DoFs for the specified component
std::vector<dof_id_type> comp_n_dofs;
_dof_map.local_variable_indices(comp_n_dofs, _mesh, comp_num);

// Create a map to store the index of each DoF in comp_0_dofs
std::unordered_map<dof_id_type, size_t> index_map;
for (size_t i = 0; i < comp_0_dofs.size(); ++i)
{
index_map[comp_0_dofs[i]] = i;
}
// We expect that the ordering of the dofs from one compoent to another
// is the same, relatively.
// Map the input dof_indices to the corresponding DoFs in comp_n_dofs
for (size_t i = 0; i < dof_indices.size(); ++i)
{
// Find the index of the current DoF in comp_0_dofs using the map
auto it = index_map.find(dof_indices[i]);
if (it != index_map.end())
{
// Use the found index to get the corresponding DoF in comp_n_dofs
new_dof_indices[i] = comp_n_dofs[it->second];
}
else
{
mooseError("DoF not found in comp_0_dofs");
}
}
for (auto & id : new_dof_indices)
id += component * _stride[id];
}

return new_dof_indices;
}

void
MooseVariableBase::setStride() const
{
if (_stride.empty())
{
std::vector<dof_id_type> comp_0_dofs;
_dof_map.local_variable_indices(comp_0_dofs, _mesh, _var_num);
// Determine the variable number for the specified component
auto comp_num = _var_num + 1;

// Grab all DOFs for the specified component
std::vector<dof_id_type> comp_1_dofs;
_dof_map.local_variable_indices(comp_1_dofs, _mesh, comp_num);

// This hangs when a processor doesn't have the variable.
_communicator.allgather(comp_0_dofs);
_communicator.allgather(comp_1_dofs);

_stride.resize(*std::max_element(comp_0_dofs.begin(), comp_0_dofs.end()) + 1,
std::numeric_limits<dof_id_type>::max());
for (const auto i : index_range(comp_0_dofs))
_stride[comp_0_dofs[i]] = comp_1_dofs[i] - comp_0_dofs[i];
}
}

void
MooseVariableBase::scalingFactor(const std::vector<Real> & factor)
{
Expand Down

0 comments on commit a2f5738

Please sign in to comment.