10 #if !defined(KRATOS_ASSIGN_SCALAR_VARIABLES_TO_PFEM_ENTITIES_PROCESS_H_INCLUDED)
11 #define KRATOS_ASSIGN_SCALAR_VARIABLES_TO_PFEM_ENTITIES_PROCESS_H_INCLUDED
44 enum class AssignmentType { DIRECT, ADDITION, SUBTRACTION, MULTIPLICATION, DIVISION };
59 "model_part_name":"MODEL_PART_NAME",
60 "variable_name": "VARIABLE_NAME",
61 "entity_type": "NODES",
63 "compound_assignment": "direct"
70 mvariable_name = rParameters[
"variable_name"].
GetString();
72 if( rParameters[
"entity_type"].GetString() ==
"NODES" ){
73 mEntity = EntityType::NODES;
75 else if( rParameters[
"entity_type"].GetString() ==
"CONDITIONS" ){
76 mEntity = EntityType::CONDITIONS;
78 else if( rParameters[
"entity_type"].GetString() ==
"ELEMENTS" ){
79 mEntity = EntityType::ELEMENTS;
82 KRATOS_ERROR <<
" Entity type "<< rParameters[
"entity_type"].
GetString() <<
" is not supported " << std::endl;
85 if( mEntity == EntityType::NODES ){
91 KRATOS_ERROR <<
"trying to set a variable that is not in the model_part - variable name is " << mvariable_name << std::endl;
94 mdouble_value = rParameters[
"value"].
GetDouble();
100 KRATOS_ERROR <<
"trying to set a variable that is not in the model_part - variable name is " << mvariable_name << std::endl;
103 mint_value = rParameters[
"value"].
GetInt();
110 KRATOS_ERROR <<
"trying to set a variable that is not in the model_part - variable name is " << mvariable_name << std::endl;
113 mbool_value = rParameters[
"value"].
GetBool();
117 else if( mEntity == EntityType::CONDITIONS || mEntity == EntityType::ELEMENTS ){
121 mdouble_value = rParameters[
"value"].
GetDouble();
125 mint_value = rParameters[
"value"].
GetInt();
129 mbool_value = rParameters[
"value"].
GetBool();
132 KRATOS_ERROR <<
"trying to set a variable that is not in the model_part - variable name is " << mvariable_name << std::endl;
137 this->SetAssignmentType(rParameters[
"compound_assignment"].GetString(), mAssignment);
164 void Execute()
override;
205 void ExecuteFinalize()
override;
223 std::string
Info()
const override
225 return "AssignScalarVariableToPfemEntitiesProcess";
231 rOStream <<
"AssignScalarVariableToPfemEntitiesProcess";
283 if( method ==
"direct" ){
284 rAssignment = AssignmentType::DIRECT;
286 else if( method ==
"addition" ){
287 rAssignment = AssignmentType::ADDITION;
289 else if( method ==
"subtraction" ){
290 rAssignment = AssignmentType::SUBTRACTION;
292 else if( method ==
"multiplication" ){
293 rAssignment = AssignmentType::MULTIPLICATION;
295 else if( method ==
"division" ){
296 rAssignment = AssignmentType::DIVISION;
299 KRATOS_ERROR <<
" Assignment type "<< method <<
" is not supported " << std::endl;
306 template<
class TVarType,
class TDataType >
312 template<
class TVarType,
class TDataType >
318 template<
class TVarType,
class TDataType >
324 template<
class TVarType,
class TDataType >
330 template<
class TVarType,
class TDataType >
361 template<
class TEntityType,
class TVarType,
class TDataType >
364 rEntity.SetValue(rVariable,value);
367 template<
class TEntityType,
class TVarType,
class TDataType >
368 void AddAssignValue(TEntityType& rEntity,
const TVarType& rVariable,
const TDataType& value)
370 TDataType AddedValue = rEntity.GetValue(rVariable)+value;
371 rEntity.SetValue(rVariable,AddedValue);
374 template<
class TEntityType,
class TVarType,
class TDataType >
377 TDataType SubtractedValue = rEntity.GetValue(rVariable)-value;
378 rEntity.SetValue(rVariable,SubtractedValue);
381 template<
class TEntityType,
class TVarType,
class TDataType >
384 TDataType MultipliedValue = rEntity.GetValue(rVariable)*value;
385 rEntity.SetValue(rVariable,MultipliedValue);
388 template<
class TEntityType,
class TVarType,
class TDataType >
391 TDataType DividedValue = rEntity.GetValue(rVariable);
394 rEntity.SetValue(rVariable,DividedValue);
397 template<
class TEntityType >
400 Vector MultipliedValue = rEntity.GetValue(rVariable);
401 for(
unsigned int i=0;
i<MultipliedValue.size(); ++
i)
403 MultipliedValue[
i]*=value[
i];
405 rEntity.SetValue(rVariable,MultipliedValue);
408 template<
class TEntityType >
411 Vector DividedValue = rEntity.GetValue(rVariable);
412 for(
unsigned int i=0;
i<DividedValue.size(); ++
i)
415 DividedValue[
i]/=value[
i];
417 rEntity.SetValue(rVariable,DividedValue);
420 template<
class TEntityType >
423 Vector MultipliedValue = rEntity.GetValue(rVariable);
424 for(
unsigned int i=0;
i<3; ++
i)
426 MultipliedValue[
i]*=value[
i];
428 rEntity.SetValue(rVariable,MultipliedValue);
431 template<
class TEntityType >
434 Vector DividedValue = rEntity.GetValue(rVariable);
435 for(
unsigned int i=0;
i<3; ++
i)
438 DividedValue[
i]/=value[
i];
440 rEntity.SetValue(rVariable,DividedValue);
468 template<
class TMethodPo
interType >
471 TMethodPointerType AssignmentMethod =
nullptr;
472 switch( mAssignment )
474 case AssignmentType::DIRECT:
477 case AssignmentType::ADDITION:
480 case AssignmentType::SUBTRACTION:
483 case AssignmentType::MULTIPLICATION:
486 case AssignmentType::DIVISION:
490 KRATOS_ERROR <<
"Unexpected value for Assignment method " << std::endl;
492 return AssignmentMethod;
516 double mdouble_value;
522 template<
class TVarType,
class TDataType >
523 void AssignValueToNodes(
const TVarType& rVariable,
const TDataType value);
525 template<
class TVarType,
class TDataType >
526 void AssignValueToConditions(
const TVarType& rVariable,
const TDataType value);
528 template<
class TVarType,
class TDataType >
529 void AssignValueToElements(
const TVarType& rVariable,
const TDataType value);
574 rOStream << std::endl;
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:34
AssignScalarVariableToPfemEntitiesProcess(ModelPart &rModelPart, Parameters rParameters)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:52
TMethodPointerType GetAssignmentMethod()
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:469
void DivideAssignValue(TEntityType &rEntity, const Variable< array_1d< double, 3 >> &rVariable, const array_1d< double, 3 > &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:432
void SubtractAssignValue(ModelPart::NodeType &rNode, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:319
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:186
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:229
std::string Info() const override
Turn back information as a string.
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:223
void ExecuteInitialize() override
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:169
AssignScalarVariableToPfemEntitiesProcess(AssignScalarVariableToPfemEntitiesProcess const &rOther)
Copy constructor.
ModelPart & mrModelPart
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:250
void MultiplyAssignValue(TEntityType &rEntity, const TVarType &rVariable, const TDataType value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:382
EntityType mEntity
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:253
void ExecuteBeforeSolutionLoop() override
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:175
void DirectAssignValue(TEntityType &rEntity, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:362
void SetAssignmentType(std::string method, AssignmentType &rAssignment)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:272
void MultiplyAssignValue(TEntityType &rEntity, const Variable< array_1d< double, 3 >> &rVariable, const array_1d< double, 3 > &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:421
void DivideAssignValue(TEntityType &rEntity, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:389
void MultiplyAssignValue(TEntityType &rEntity, const Variable< Vector > &rVariable, const Vector &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:398
~AssignScalarVariableToPfemEntitiesProcess() override
Destructor.
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:144
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:152
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:181
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:192
void AddAssignValue(TEntityType &rEntity, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:368
KRATOS_CLASS_POINTER_DEFINITION(AssignScalarVariableToPfemEntitiesProcess)
Pointer definition of AssignScalarVariableToPfemEntitiesProcess.
EntityType
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:42
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:235
AssignmentType
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:44
std::string mvariable_name
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:251
AssignScalarVariableToPfemEntitiesProcess(ModelPart &rModelPart)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:50
void AddAssignValue(ModelPart::NodeType &rNode, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:313
void SubtractAssignValue(TEntityType &rEntity, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:375
void DivideAssignValue(TEntityType &rEntity, const Variable< Vector > &rVariable, const Vector &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:409
void DivideAssignValue(ModelPart::NodeType &rNode, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:331
void MultiplyAssignValue(ModelPart::NodeType &rNode, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:325
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:198
void DirectAssignValue(ModelPart::NodeType &rNode, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:307
AssignmentType mAssignment
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:254
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
VariablesList & GetNodalSolutionStepVariablesList()
Definition: model_part.h:549
This class defines the node.
Definition: node.h:65
TVariableType::Type & FastGetSolutionStepValue(const TVariableType &rThisVariable)
Definition: node.h:435
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
double GetDouble() const
This method returns the double contained in the current Parameter.
Definition: kratos_parameters.cpp:657
int GetInt() const
This method returns the integer contained in the current Parameter.
Definition: kratos_parameters.cpp:666
void ValidateAndAssignDefaults(const Parameters &rDefaultParameters)
This function is designed to verify that the parameters under testing match the form prescribed by th...
Definition: kratos_parameters.cpp:1306
std::string GetString() const
This method returns the string contained in the current Parameter.
Definition: kratos_parameters.cpp:684
bool GetBool() const
This method returns the boolean contained in the current Parameter.
Definition: kratos_parameters.cpp:675
The base class for all processes in Kratos.
Definition: process.h:49
bool Has(const VariableData &rThisVariable) const
Definition: variables_list.h:372
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR
Definition: exception.h:161
bool Has(const std::string &ModelerName)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:24
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
integer i
Definition: TensorModule.f:17