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

[Core] Add distribute variable method to variable utils #7146

Merged
merged 15 commits into from
Jun 28, 2020
Merged
Show file tree
Hide file tree
Changes from all 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
8 changes: 8 additions & 0 deletions kratos/python/add_variable_utils_to_python.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -325,6 +325,14 @@ void AddVariableUtilsToPython(pybind11::module &m)
.def("ClearNonHistoricalData", &VariableUtils::ClearNonHistoricalData<ModelPart::NodesContainerType>)
.def("ClearNonHistoricalData", &VariableUtils::ClearNonHistoricalData<ModelPart::ConditionsContainerType>)
.def("ClearNonHistoricalData", &VariableUtils::ClearNonHistoricalData<ModelPart::ElementsContainerType>)
.def("WeightedAccumulateConditionVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes<double, ModelPart::ConditionsContainerType, int>)
.def("WeightedAccumulateConditionVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes<array_1d<double, 3>, ModelPart::ConditionsContainerType, int>)
.def("WeightedAccumulateElementVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes<double, ModelPart::ElementsContainerType, int>)
.def("WeightedAccumulateElementVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes<array_1d<double, 3>, ModelPart::ElementsContainerType, int>)
.def("WeightedAccumulateConditionVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes<double, ModelPart::ConditionsContainerType, double>)
.def("WeightedAccumulateConditionVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes<array_1d<double, 3>, ModelPart::ConditionsContainerType, double>)
.def("WeightedAccumulateElementVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes<double, ModelPart::ElementsContainerType, double>)
.def("WeightedAccumulateElementVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes<array_1d<double, 3>, ModelPart::ElementsContainerType, double>)
.def("SetFlag", &VariableUtils::SetFlag<ModelPart::NodesContainerType>)
.def("SetFlag", &VariableUtils::SetFlag<ModelPart::ConditionsContainerType>)
.def("SetFlag", &VariableUtils::SetFlag<ModelPart::ElementsContainerType>)
Expand Down
83 changes: 83 additions & 0 deletions kratos/tests/test_variable_utils.py
Original file line number Diff line number Diff line change
Expand Up @@ -761,6 +761,89 @@ def test_UpdateCurrentPosition(self):
self.assertAlmostEqual(node.Y, node.Y0 + 4.0 * float(node.Id))
self.assertAlmostEqual(node.Z, node.Z0 + 5.0 * float(node.Id))

def test_distribute_condition_variable(self):
current_model = KratosMultiphysics.Model()

##set the model part
model_part = current_model.CreateModelPart("Main")
model_part.AddNodalSolutionStepVariable(KratosMultiphysics.VISCOSITY)
model_part.AddNodalSolutionStepVariable(KratosMultiphysics.DISPLACEMENT)
model_part_io = KratosMultiphysics.ModelPartIO(GetFilePath("auxiliar_files_for_python_unittest/mdpa_files/test_model_part_io_read"))
model_part_io.ReadModelPart(model_part)

for node in model_part.Nodes:
node.SetValue(KratosMultiphysics.AUX_MESH_VAR, node.Id)

for condition in model_part.Conditions:
condition.SetValue(KratosMultiphysics.DISTANCE, condition.Id)
vector = KratosMultiphysics.Vector(3)
vector[0] = condition.Id * 3
vector[1] = condition.Id * 3 + 1
vector[2] = condition.Id * 3 + 2
condition.SetValue(KratosMultiphysics.VELOCITY, vector)

variable_utils = KratosMultiphysics.VariableUtils()
variable_utils.WeightedAccumulateConditionVariableOnNodes(model_part, KratosMultiphysics.DISTANCE, KratosMultiphysics.AUX_MESH_VAR, False)
variable_utils.WeightedAccumulateConditionVariableOnNodes(model_part, KratosMultiphysics.VELOCITY, KratosMultiphysics.AUX_MESH_VAR, False)

distance_vector = [
1.0, 3602.0, 10803.0, 3643056.0, 3789835.0, 1897352.0
]
velocity_vector = [
3.0, 4.0, 5.0, 10806.0, 10810.0, 10814.0, 32409.0, 32415.0, 32421.0, 10929168.0, 10931112.0, 10933056.0, 11369505.0, 11371451.0, 11373397.0, 5692056.0, 5693030.0, 5694004.0
]

local_index = 0
for node in model_part.Nodes:
self.assertEqual(node.GetValue(KratosMultiphysics.DISTANCE), distance_vector[local_index])
self.assertEqual(node.GetValue(KratosMultiphysics.VELOCITY)[0], velocity_vector[local_index * 3])
self.assertEqual(node.GetValue(KratosMultiphysics.VELOCITY)[1], velocity_vector[local_index * 3 + 1])
self.assertEqual(node.GetValue(KratosMultiphysics.VELOCITY)[2], velocity_vector[local_index * 3 + 2])

local_index += 1

def test_distribute_condition_variable_inverse(self):
current_model = KratosMultiphysics.Model()

##set the model part
model_part = current_model.CreateModelPart("Main")
model_part.AddNodalSolutionStepVariable(KratosMultiphysics.VISCOSITY)
model_part.AddNodalSolutionStepVariable(KratosMultiphysics.DISPLACEMENT)
model_part_io = KratosMultiphysics.ModelPartIO(GetFilePath("auxiliar_files_for_python_unittest/mdpa_files/test_model_part_io_read"))
model_part_io.ReadModelPart(model_part)

for node in model_part.Nodes:
node.SetValue(KratosMultiphysics.AUX_MESH_VAR, 1.0 / node.Id)

for condition in model_part.Conditions:
condition.SetValue(KratosMultiphysics.DISTANCE, condition.Id)
vector = KratosMultiphysics.Vector(3)
vector[0] = condition.Id * 3
vector[1] = condition.Id * 3 + 1
vector[2] = condition.Id * 3 + 2
condition.SetValue(KratosMultiphysics.VELOCITY, vector)

variable_utils = KratosMultiphysics.VariableUtils()
variable_utils.WeightedAccumulateConditionVariableOnNodes(model_part, KratosMultiphysics.DISTANCE, KratosMultiphysics.AUX_MESH_VAR, True)
variable_utils.WeightedAccumulateConditionVariableOnNodes(model_part, KratosMultiphysics.VELOCITY, KratosMultiphysics.AUX_MESH_VAR, True)

distance_vector = [
1.0, 3602.0, 10803.0, 3643056.0, 3789835.0, 1897352.0
]
velocity_vector = [
3.0, 4.0, 5.0, 10806.0, 10810.0, 10814.0, 32409.0, 32415.0, 32421.0, 10929168.0, 10931112.0, 10933056.0, 11369505.0, 11371451.0, 11373397.0, 5692056.0, 5693030.0, 5694004.0
]

local_index = 0
for node in model_part.Nodes:
self.assertAlmostEqual(node.GetValue(KratosMultiphysics.DISTANCE), distance_vector[local_index])
self.assertAlmostEqual(node.GetValue(KratosMultiphysics.VELOCITY)[0], velocity_vector[local_index * 3])
self.assertAlmostEqual(node.GetValue(KratosMultiphysics.VELOCITY)[1], velocity_vector[local_index * 3 + 1])
self.assertAlmostEqual(node.GetValue(KratosMultiphysics.VELOCITY)[2], velocity_vector[local_index * 3 + 2])

local_index += 1


if __name__ == '__main__':
KratosMultiphysics.Logger.GetDefaultOutput().SetSeverity(KratosMultiphysics.Logger.Severity.WARNING)
KratosUnittest.main()
82 changes: 82 additions & 0 deletions kratos/utilities/variable_utils.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -14,6 +14,7 @@
//

/* System includes */
#include <functional>

/* External includes */

Expand Down Expand Up @@ -486,4 +487,85 @@ void VariableUtils::AuxiliaryAtomicAdd(
rSumValue[2] += rPrivateValue[2];
}

template <>
ModelPart::ElementsContainerType& VariableUtils::GetContainer<ModelPart::ElementsContainerType>(ModelPart& rModelPart)
{
return rModelPart.Elements();
}

template <>
ModelPart::ConditionsContainerType& VariableUtils::GetContainer<ModelPart::ConditionsContainerType>(ModelPart& rModelPart)
{
return rModelPart.Conditions();
}

template <class TDataType, class TContainerType, class TWeightDataType>
void VariableUtils::WeightedAccumulateVariableOnNodes(
ModelPart& rModelPart,
rubenzorrilla marked this conversation as resolved.
Show resolved Hide resolved
const Variable<TDataType>& rVariable,
const Variable<TWeightDataType>& rWeightVariable,
const bool IsInverseWeightProvided)
{
KRATOS_TRY

SetNonHistoricalVariableToZero(rVariable, rModelPart.Nodes());

auto& r_entities = GetContainer<TContainerType>(rModelPart);
const int n_entities = r_entities.size();

const std::function<double(const Node<3>&)>& r_weight_method =
(IsInverseWeightProvided) ?
static_cast<std::function<double(const Node<3>&)>>([rWeightVariable](const Node<3>& rNode) -> double {return 1.0 / rNode.GetValue(rWeightVariable);}) :
static_cast<std::function<double(const Node<3>&)>>([rWeightVariable](const Node<3>& rNode) -> double {return rNode.GetValue(rWeightVariable);});
Comment on lines +516 to +519
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

  1. I think using std function has some overhead, use with caution in performance critical parts
  2. There is no check for rNode.GetValue(rWeightVariable). If it is 0 then it will crash. This happens if this var was not allocate

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

  1. I added the check in full debug only. Should I change it to release as well? I didnt do it considering the performance impact. (
    KRATOS_DEBUG_ERROR_IF(!r_node.Has(rWeightVariable))
    << "Non-historical nodal " << rWeightVariable.Name() << " at "
    << r_node << " is not initialized in " << rModelPart.Name()
    << ". Please initialize it first.";
    )
  2. As for the zeros, I can add a check in full debug again. Would it be sufficient?


#pragma omp parallel for
for (int i_entity = 0; i_entity < n_entities; ++i_entity)
{
auto it_entity = r_entities.begin() + i_entity;
auto& r_geometry = it_entity->GetGeometry();

const auto& r_value = it_entity->GetValue(rVariable);
for (int i_node = 0; i_node < static_cast<int>(r_geometry.PointsNumber()); ++i_node)
{
auto& r_node = r_geometry[i_node];

KRATOS_DEBUG_ERROR_IF(!r_node.Has(rWeightVariable))
<< "Non-historical nodal " << rWeightVariable.Name() << " at "
<< r_node << " is not initialized in " << rModelPart.Name()
<< ". Please initialize it first.";

const double weight = r_weight_method(r_node);

r_node.SetLock();
r_node.GetValue(rVariable) += r_value * weight;
r_node.UnSetLock();
}
}

rModelPart.GetCommunicator().AssembleNonHistoricalData(rVariable);

KRATOS_CATCH("");
}

// template instantiations
template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes<double, ModelPart::ConditionsContainerType, int>(
ModelPart&, const Variable<double>&, const Variable<int>&, const bool);
template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes<array_1d<double, 3>, ModelPart::ConditionsContainerType, int>(
ModelPart&, const Variable<array_1d<double, 3>>&, const Variable<int>&, const bool);

template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes<double, ModelPart::ElementsContainerType, int>(
ModelPart&, const Variable<double>&, const Variable<int>&, const bool);
template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes<array_1d<double, 3>, ModelPart::ElementsContainerType, int>(
ModelPart&, const Variable<array_1d<double, 3>>&, const Variable<int>&, const bool);

template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes<double, ModelPart::ConditionsContainerType, double>(
ModelPart&, const Variable<double>&, const Variable<double>&, const bool);
template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes<array_1d<double, 3>, ModelPart::ConditionsContainerType, double>(
ModelPart&, const Variable<array_1d<double, 3>>&, const Variable<double>&, const bool);

template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes<double, ModelPart::ElementsContainerType, double>(
ModelPart&, const Variable<double>&, const Variable<double>&, const bool);
template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes<array_1d<double, 3>, ModelPart::ElementsContainerType, double>(
ModelPart&, const Variable<array_1d<double, 3>>&, const Variable<double>&, const bool);

} /* namespace Kratos.*/
28 changes: 28 additions & 0 deletions kratos/utilities/variable_utils.h
Original file line number Diff line number Diff line change
Expand Up @@ -504,6 +504,31 @@ class KRATOS_API(KRATOS_CORE) VariableUtils
KRATOS_CATCH("")
}

/**
* @brief Distributes variable values in TContainerType container to nodes
*
* This method distributes variables values stored in TContainerType data value container in rModelPart
* to nodes. Constant weighting is used for each node based on rWeightVariable value. The result
* is stored in nodal non-historical data value container under the same rVariable. If IsInverseWeightProvided
* is true, then the weights provided by rWeightVariable is inverted to get nodal weight. Otherwise, the value
* given by rWeightVariable is used as weight.
*
*
* @tparam TDataType Data type
* @tparam TContainerType ContainerType of model part
* @tparam TWeightDataType Data type of weight variable (this should be either int or double)
* @param rModelPart Model part
* @param rVariable Variable to be distributed
* @param rWeightVariable Variable which holds weight to distribute entity values to nodes
* @param IsInverseWeightProvided Whether the weight is provided as inverse or not.
*/
template <class TDataType, class TContainerType, class TWeightDataType>
void WeightedAccumulateVariableOnNodes(
ModelPart& rModelPart,
const Variable<TDataType>& rVariable,
const Variable<TWeightDataType>& rWeightVariable,
rubenzorrilla marked this conversation as resolved.
Show resolved Hide resolved
const bool IsInverseWeightProvided = false);

/**
* @brief Sets a flag according to a given status over a given container
* @param rFlag flag to be set
Expand Down Expand Up @@ -1244,6 +1269,9 @@ class KRATOS_API(KRATOS_CORE) VariableUtils
KRATOS_CATCH("")
}

template <class TContainerType>
TContainerType& GetContainer(ModelPart& rModelPart);

///@}
///@name Private Acces
///@{
Expand Down