45 "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
46 "variable_name": "PLEASE_PRESCRIBE_VARIABLE_NAME",
52 rParameters[
"variable_name"];
53 rParameters[
"model_part_name"];
54 rParameters[
"append_file"];
59 mModelPartName = rParameters[
"model_part_name"].
GetString();
60 mVariableName = rParameters[
"variable_name"].
GetString();
61 mAppendFile = rParameters[
"append_file"].
GetBool();
62 mTimeUnitConverter =
model_part.GetProcessInfo()[TIME_UNIT_CONVERTER];
88 const double Time = mrModelPart.
GetProcessInfo()[TIME]/mTimeUnitConverter;
90 ModelPart::NodesContainerType::iterator it_begin = mrModelPart.
NodesBegin();
91 mOutFile.resize(nNodes);
94 ModelPart::NodesContainerType::iterator it = it_begin +
i;
97 std::string fileName = mModelPartName +
"_" + std::to_string(nodeId) +
"_" + mVariableName +
".res";
101 mOutFile[
i].open(fileName, std::ios::app);
104 mOutFile[
i].open(fileName, std::ios::trunc);
105 mOutFile[
i] <<
"Time" <<
" " << mVariableName <<
"\n";
106 const double value = it->FastGetSolutionStepValue(var);
107 mOutFile[
i] << Time <<
" " << value <<
"\n";
126 const double Time = mrModelPart.
GetProcessInfo()[TIME]/mTimeUnitConverter;
127 ModelPart::NodesContainerType::iterator it_begin = mrModelPart.
NodesBegin();
130 ModelPart::NodesContainerType::iterator it = it_begin +
i;
132 const double value = it->FastGetSolutionStepValue(var);
133 mOutFile[
i] << Time <<
" " << value <<
"\n";
147 for (
auto& outFile : mOutFile) {
155 std::string
Info()
const override
157 return "ApplyWriteScalarProcess";
163 std::string mVariableName;
164 std::string mModelPartName;
166 double mTimeUnitConverter;
167 std::vector<std::ofstream> mOutFile;
180 rOStream << std::endl;
Definition: apply_write_result_scalar_process.hpp:26
void ExecuteInitialize() override
Definition: apply_write_result_scalar_process.hpp:80
void ExecuteFinalize() override
This function is designed for being called at the end of the computations.
Definition: apply_write_result_scalar_process.hpp:143
void Execute() override
Execute method is used to execute the ApplyWriteScalarProcess algorithms.
Definition: apply_write_result_scalar_process.hpp:74
ApplyWriteScalarProcess(ModelPart &model_part, Parameters rParameters)
Definition: apply_write_result_scalar_process.hpp:36
ApplyWriteScalarProcess(const ApplyWriteScalarProcess &)=delete
ApplyWriteScalarProcess(ApplyWriteScalarProcess &&)=delete
~ApplyWriteScalarProcess() override=default
ApplyWriteScalarProcess & operator=(ApplyWriteScalarProcess &&)=delete
void ExecuteFinalizeSolutionStep() override
This function will be executed at every time step AFTER performing the solve phase.
Definition: apply_write_result_scalar_process.hpp:118
KRATOS_CLASS_POINTER_DEFINITION(ApplyWriteScalarProcess)
ApplyWriteScalarProcess & operator=(const ApplyWriteScalarProcess &)=delete
std::string Info() const override
Turn back information as a string.
Definition: apply_write_result_scalar_process.hpp:155
std::size_t IndexType
Definition: flags.h:74
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
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
SizeType NumberOfNodes(IndexType ThisIndex=0) const
Definition: model_part.h:341
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
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
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: process.h:204
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: process.h:210
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
std::size_t SizeType
The definition of the size type.
Definition: mortar_classes.h:43
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
model_part
Definition: face_heat.py:14
integer i
Definition: TensorModule.f:17