KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
apply_write_result_scalar_process.hpp
Go to the documentation of this file.
1 // KRATOS___
2 // // ) )
3 // // ___ ___
4 // // ____ //___) ) // ) )
5 // // / / // // / /
6 // ((____/ / ((____ ((___/ / MECHANICS
7 //
8 // License: geo_mechanics_application/license.txt
9 //
10 // Main authors: Vahid Galavi
11 //
12 
13 #pragma once
14 
15 #include "includes/table.h"
16 #include "includes/kratos_flags.h"
18 #include "processes/process.h"
19 
21 
22 namespace Kratos
23 {
24 
26 {
27 
28 public:
29 
31 
34  using IndexType = std::size_t;
35 
37  Parameters rParameters
38  ) : Process(Flags()) , mrModelPart(model_part)
39  {
41 
42  //only include validation with c++11 since raw_literals do not exist in c++03
43  Parameters default_parameters( R"(
44  {
45  "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
46  "variable_name": "PLEASE_PRESCRIBE_VARIABLE_NAME",
47  "append_file" : false
48  } )" );
49 
50  // Some values need to be mandatorily prescribed since no meaningful default value exist. For this reason try accessing to them
51  // So that an error is thrown if they don't exist
52  rParameters["variable_name"];
53  rParameters["model_part_name"];
54  rParameters["append_file"];
55 
56  // Now validate agains defaults -- this also ensures no type mismatch
57  rParameters.ValidateAndAssignDefaults(default_parameters);
58 
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];
63 
64  KRATOS_CATCH("")
65  }
66 
67  ~ApplyWriteScalarProcess() override = default;
72 
74  void Execute() override
75  {
76  }
77 
80  void ExecuteInitialize() override
81  {
83 
84  const SizeType nNodes = mrModelPart.NumberOfNodes();
85 
86  if (nNodes > 0) {
87  const Variable<double> &var = KratosComponents< Variable<double> >::Get(mVariableName);
88  const double Time = mrModelPart.GetProcessInfo()[TIME]/mTimeUnitConverter;
89 
90  ModelPart::NodesContainerType::iterator it_begin = mrModelPart.NodesBegin();
91  mOutFile.resize(nNodes);
92 
93  for (IndexType i = 0; i<nNodes; ++i) {
94  ModelPart::NodesContainerType::iterator it = it_begin + i;
95 
96  const IndexType nodeId = it->Id();
97  std::string fileName = mModelPartName + "_" + std::to_string(nodeId) + "_" + mVariableName + ".res";
98 
99  if (mAppendFile) {
100  // append instead of overwrite
101  mOutFile[i].open(fileName, std::ios::app);
102  } else {
103  // open a new file and overwrite
104  mOutFile[i].open(fileName, std::ios::trunc); // overwrite
105  mOutFile[i] << "Time" << " " << mVariableName << "\n";
106  const double value = it->FastGetSolutionStepValue(var);
107  mOutFile[i] << Time << " " << value << "\n";
108  }
109  }
110  }
111 
112  KRATOS_CATCH("")
113  }
114 
119  {
120  KRATOS_TRY
121 
122  const IndexType nNodes = mrModelPart.NumberOfNodes();
123 
124  if (nNodes > 0) {
125  const Variable<double> &var = KratosComponents< Variable<double> >::Get(mVariableName);
126  const double Time = mrModelPart.GetProcessInfo()[TIME]/mTimeUnitConverter;
127  ModelPart::NodesContainerType::iterator it_begin = mrModelPart.NodesBegin();
128 
129  for (IndexType i = 0; i<nNodes; ++i) {
130  ModelPart::NodesContainerType::iterator it = it_begin + i;
131 
132  const double value = it->FastGetSolutionStepValue(var);
133  mOutFile[i] << Time << " " << value << "\n";
134  }
135  }
136 
137  KRATOS_CATCH("")
138  }
139 
143  void ExecuteFinalize() override
144  {
145  KRATOS_TRY
146 
147  for (auto& outFile : mOutFile) {
148  outFile.close();
149  }
150 
151  KRATOS_CATCH("")
152  }
153 
155  std::string Info() const override
156  {
157  return "ApplyWriteScalarProcess";
158  }
159 
160 private:
162  ModelPart& mrModelPart;
163  std::string mVariableName;
164  std::string mModelPartName;
165  bool mAppendFile;
166  double mTimeUnitConverter;
167  std::vector<std::ofstream> mOutFile;
168 
169 }; // Class ApplyWriteScalarProcess
170 
172 inline std::istream& operator >> (std::istream& rIStream,
173  ApplyWriteScalarProcess& rThis);
174 
176 inline std::ostream& operator << (std::ostream& rOStream,
177  const ApplyWriteScalarProcess& rThis)
178 {
179  rThis.PrintInfo(rOStream);
180  rOStream << std::endl;
181  rThis.PrintData(rOStream);
182 
183  return rOStream;
184 }
185 
186 }
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
Definition: flags.h:58
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
Definition: table.h:435
#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