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.
dgeosettlement.h
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: Anne van de Graaf
11 //
12 
13 #pragma once
14 
15 #include <filesystem>
16 #include <functional>
17 #include <string>
18 
19 #include "includes/kernel.h"
21 
27 
28 namespace Kratos
29 {
30 
31 class InputUtility;
32 class ProcessInfoParser;
33 class Process;
34 class TimeLoopExecutorInterface;
35 class TimeIncrementor;
36 class StrategyWrapper;
37 
38 class KRATOS_API(GEO_MECHANICS_APPLICATION) KratosGeoSettlement
39 {
40 public:
41  KratosGeoSettlement(std::unique_ptr<InputUtility> pInputUtility,
42  std::unique_ptr<ProcessInfoParser> pProcessInfoParser,
43  std::unique_ptr<TimeLoopExecutorInterface> pTimeLoopExecutorInterface);
44 
46 
47  int RunStage(const std::filesystem::path& rWorkingDirectory,
48  const std::filesystem::path& rProjectParametersFile,
49  const std::function<void(const char*)>& rLogCallback,
50  const std::function<void(double)>& rReportProgress,
51  const std::function<void(const char*)>& rReportTextualProgress,
52  const std::function<bool()>& rShouldCancel);
53 
54  const InputUtility* GetInterfaceInputUtility() const;
55 
56 private:
57  ModelPart& AddNewModelPart(const std::string& rModelPartName);
58  void PrepareModelPart(const Parameters& rSolverSettings);
59 
60  ModelPart& GetMainModelPart();
61  ModelPart& GetComputationalModelPart();
62 
63  static void AddNodalSolutionStepVariablesTo(ModelPart& rModelPart);
64  static void AddDegreesOfFreedomTo(ModelPart& rModelPart);
65  void InitializeProcessFactory();
66  std::vector<std::shared_ptr<Process>> GetProcesses(const Parameters& project_parameters) const;
67  static std::unique_ptr<TimeIncrementor> MakeTimeIncrementor(const Parameters& rProjectParameters);
68  std::shared_ptr<StrategyWrapper> MakeStrategyWrapper(const Parameters& rProjectParameters,
69  const std::filesystem::path& rWorkingDirectory);
70  LoggerOutput::Pointer CreateLoggingOutput(std::stringstream& rKratosLogBuffer) const;
71  void FlushLoggingOutput(const std::function<void(const char*)>& rLogCallback,
72  LoggerOutput::Pointer pLoggerOutput,
73  const std::stringstream& rKratosLogBuffer) const;
74 
75  template <typename TVariableType>
76  void RestoreValuesOfNodalVariable(const TVariableType& rVariable, Node::IndexType SourceIndex, Node::IndexType DestinationIndex)
77  {
78  if (!GetComputationalModelPart().HasNodalSolutionStepVariable(rVariable)) return;
79 
80  VariableUtils{}.SetHistoricalVariableToZero(rVariable, GetComputationalModelPart().Nodes());
81 
82  block_for_each(GetComputationalModelPart().Nodes(),
83  [&rVariable, SourceIndex, DestinationIndex](auto& node) {
84  node.GetSolutionStepValue(rVariable, DestinationIndex) =
85  node.GetSolutionStepValue(rVariable, SourceIndex);
86  });
87  }
88 
89  template <typename ProcessType>
90  std::function<ProcessFactory::ProductType(const Parameters&)> MakeCreatorFor()
91  {
92  return [&model = mModel](const Parameters& rProcessSettings) {
93  auto& model_part = model.GetModelPart(rProcessSettings["model_part_name"].GetString());
94  return std::make_unique<ProcessType>(model_part, rProcessSettings);
95  };
96  }
97 
98  Kernel mKernel;
99  Model mModel;
100  std::string mModelPartName;
101  KratosGeoMechanicsApplication::Pointer mpGeoApp;
102  KratosLinearSolversApplication::Pointer mpLinearSolversApp;
103  KratosStructuralMechanicsApplication::Pointer mpStructuralMechanicsApp;
104  std::unique_ptr<ProcessFactory> mProcessFactory = std::make_unique<ProcessFactory>();
105  std::unique_ptr<InputUtility> mpInputUtility;
106  std::unique_ptr<ProcessInfoParser> mpProcessInfoParser;
107  std::unique_ptr<TimeLoopExecutorInterface> mpTimeLoopExecutor;
108  const std::string mComputationalSubModelPartName{"settlement_computational_model_part"};
109 };
110 
111 } // namespace Kratos
Definition: input_utility.h:21
Definition: dgeosettlement.h:39
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
std::size_t IndexType
The index type.
Definition: node.h:86
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
std::unique_ptr< Process > ProductType
Definition: process_factory.hpp:31
This class implements a set of auxiliar, already parallelized, methods to perform some common tasks r...
Definition: variable_utils.h:63
void SetHistoricalVariableToZero(const Variable< TType > &rVariable, NodesContainerType &rNodes)
Sets the nodal value of any variable to zero.
Definition: variable_utils.h:757
string path
Definition: DEM_run_all_benchmarks_analysis.py:10
bool HasNodalSolutionStepVariable(ModelPart &rModelPart, Variable< TDataType > const &rThisVariable)
Definition: add_model_part_to_python.cpp:39
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
void block_for_each(TIterator itBegin, TIterator itEnd, TFunction &&rFunction)
Execute a functor on all items of a range in parallel.
Definition: parallel_utilities.h:299
project_parameters
Definition: concentric_element_size_coarsener.py:171
model_part
Definition: face_heat.py:14
model
Definition: fluid_chimera_analysis.py:37
dictionary Model
TODO replace this "model" for real one once available in kratos core.
Definition: script.py:94
Definition: mesh_converter.cpp:38