diff --git a/kratos/python/add_variable_utils_to_python.cpp b/kratos/python/add_variable_utils_to_python.cpp index d29fe3efec1b..478cc6865dea 100644 --- a/kratos/python/add_variable_utils_to_python.cpp +++ b/kratos/python/add_variable_utils_to_python.cpp @@ -325,6 +325,14 @@ void AddVariableUtilsToPython(pybind11::module &m) .def("ClearNonHistoricalData", &VariableUtils::ClearNonHistoricalData) .def("ClearNonHistoricalData", &VariableUtils::ClearNonHistoricalData) .def("ClearNonHistoricalData", &VariableUtils::ClearNonHistoricalData) + .def("WeightedAccumulateConditionVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes) + .def("WeightedAccumulateConditionVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes, ModelPart::ConditionsContainerType, int>) + .def("WeightedAccumulateElementVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes) + .def("WeightedAccumulateElementVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes, ModelPart::ElementsContainerType, int>) + .def("WeightedAccumulateConditionVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes) + .def("WeightedAccumulateConditionVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes, ModelPart::ConditionsContainerType, double>) + .def("WeightedAccumulateElementVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes) + .def("WeightedAccumulateElementVariableOnNodes", &VariableUtils::WeightedAccumulateVariableOnNodes, ModelPart::ElementsContainerType, double>) .def("SetFlag", &VariableUtils::SetFlag) .def("SetFlag", &VariableUtils::SetFlag) .def("SetFlag", &VariableUtils::SetFlag) diff --git a/kratos/tests/test_variable_utils.py b/kratos/tests/test_variable_utils.py index 43cc9d0fbea1..43ab61bdb0f9 100644 --- a/kratos/tests/test_variable_utils.py +++ b/kratos/tests/test_variable_utils.py @@ -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() diff --git a/kratos/utilities/variable_utils.cpp b/kratos/utilities/variable_utils.cpp index 7e5959a90a9b..15a62ecd554c 100644 --- a/kratos/utilities/variable_utils.cpp +++ b/kratos/utilities/variable_utils.cpp @@ -14,6 +14,7 @@ // /* System includes */ +#include /* External includes */ @@ -486,4 +487,85 @@ void VariableUtils::AuxiliaryAtomicAdd( rSumValue[2] += rPrivateValue[2]; } +template <> +ModelPart::ElementsContainerType& VariableUtils::GetContainer(ModelPart& rModelPart) +{ + return rModelPart.Elements(); +} + +template <> +ModelPart::ConditionsContainerType& VariableUtils::GetContainer(ModelPart& rModelPart) +{ + return rModelPart.Conditions(); +} + +template +void VariableUtils::WeightedAccumulateVariableOnNodes( + ModelPart& rModelPart, + const Variable& rVariable, + const Variable& rWeightVariable, + const bool IsInverseWeightProvided) +{ + KRATOS_TRY + + SetNonHistoricalVariableToZero(rVariable, rModelPart.Nodes()); + + auto& r_entities = GetContainer(rModelPart); + const int n_entities = r_entities.size(); + + const std::function&)>& r_weight_method = + (IsInverseWeightProvided) ? + static_cast&)>>([rWeightVariable](const Node<3>& rNode) -> double {return 1.0 / rNode.GetValue(rWeightVariable);}) : + static_cast&)>>([rWeightVariable](const Node<3>& rNode) -> double {return rNode.GetValue(rWeightVariable);}); + +#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(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( + ModelPart&, const Variable&, const Variable&, const bool); +template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes, ModelPart::ConditionsContainerType, int>( + ModelPart&, const Variable>&, const Variable&, const bool); + +template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes( + ModelPart&, const Variable&, const Variable&, const bool); +template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes, ModelPart::ElementsContainerType, int>( + ModelPart&, const Variable>&, const Variable&, const bool); + +template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes( + ModelPart&, const Variable&, const Variable&, const bool); +template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes, ModelPart::ConditionsContainerType, double>( + ModelPart&, const Variable>&, const Variable&, const bool); + +template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes( + ModelPart&, const Variable&, const Variable&, const bool); +template KRATOS_API(KRATOS_CORE) void VariableUtils::WeightedAccumulateVariableOnNodes, ModelPart::ElementsContainerType, double>( + ModelPart&, const Variable>&, const Variable&, const bool); + } /* namespace Kratos.*/ diff --git a/kratos/utilities/variable_utils.h b/kratos/utilities/variable_utils.h index e29aad8351db..592a1aa45992 100644 --- a/kratos/utilities/variable_utils.h +++ b/kratos/utilities/variable_utils.h @@ -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 + void WeightedAccumulateVariableOnNodes( + ModelPart& rModelPart, + const Variable& rVariable, + const Variable& rWeightVariable, + const bool IsInverseWeightProvided = false); + /** * @brief Sets a flag according to a given status over a given container * @param rFlag flag to be set @@ -1244,6 +1269,9 @@ class KRATOS_API(KRATOS_CORE) VariableUtils KRATOS_CATCH("") } + template + TContainerType& GetContainer(ModelPart& rModelPart); + ///@} ///@name Private Acces ///@{