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.
solution_scheme.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosSolidMechanicsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: March 2018 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_SOLUTION_SCHEME_H_INCLUDED)
11 #define KRATOS_SOLUTION_SCHEME_H_INCLUDED
12 
13 // System includes
14 //#include <set>
15 
16 // External includes
17 
18 // Project includes
19 #include "includes/checks.h"
20 #include "includes/model_part.h"
23 
24 namespace Kratos
25 {
26 
29 
33 
37 
41 
45 
50 template<class TSparseSpace,
51  class TDenseSpace //= DenseSpace<double>
52  >
53 class SolutionScheme : public Flags
54 {
55  public:
56 
60  typedef typename SolutionSchemeType::Pointer SolutionSchemePointerType;
63 
68 
72 
77  typedef typename IntegrationVectorType::Pointer IntegrationVectorPointerType;
78  typedef std::vector<IntegrationVectorPointerType> IntegrationMethodsVectorType;
79 
82  typedef typename IntegrationScalarType::Pointer IntegrationScalarPointerType;
83  typedef std::vector<IntegrationScalarPointerType> IntegrationMethodsScalarType;
84 
85  typedef typename SolverProcess::Pointer ProcessPointerType;
86  typedef std::vector<ProcessPointerType> ProcessPointerVectorType;
87 
90 
94 
97 
99  SolutionScheme(Flags& rOptions) : Flags(), mOptions(rOptions) {SetDefaultFlags();}
100 
102  SolutionScheme(IntegrationMethodsVectorType& rTimeVectorIntegrationMethods, Flags& rOptions) : Flags(), mOptions(rOptions), mTimeVectorIntegrationMethods(rTimeVectorIntegrationMethods) {SetDefaultFlags();}
103 
105  SolutionScheme(IntegrationMethodsVectorType& rTimeVectorIntegrationMethods) : Flags(), mTimeVectorIntegrationMethods(rTimeVectorIntegrationMethods) {SetDefaultFlags();}
106 
108  SolutionScheme(IntegrationMethodsScalarType& rTimeScalarIntegrationMethods, Flags& rOptions) : Flags(), mOptions(rOptions), mTimeScalarIntegrationMethods(rTimeScalarIntegrationMethods) {SetDefaultFlags();}
109 
111  SolutionScheme(IntegrationMethodsScalarType& rTimeScalarIntegrationMethods) : Flags(), mTimeScalarIntegrationMethods(rTimeScalarIntegrationMethods) {SetDefaultFlags();}
112 
114  SolutionScheme(IntegrationMethodsVectorType& rTimeVectorIntegrationMethods,
115  IntegrationMethodsScalarType& rTimeScalarIntegrationMethods,
116  Flags& rOptions)
117  : Flags(), mOptions(rOptions), mTimeVectorIntegrationMethods(rTimeVectorIntegrationMethods), mTimeScalarIntegrationMethods(rTimeScalarIntegrationMethods) {SetDefaultFlags();}
118 
120  SolutionScheme(IntegrationMethodsVectorType& rTimeVectorIntegrationMethods,
121  IntegrationMethodsScalarType& rTimeScalarIntegrationMethods)
122  : Flags(), mTimeVectorIntegrationMethods(rTimeVectorIntegrationMethods), mTimeScalarIntegrationMethods(rTimeScalarIntegrationMethods) {SetDefaultFlags();}
123 
126  , mProcesses(rOther.mProcesses)
127  {
128  std::copy(std::begin(rOther.mTimeVectorIntegrationMethods), std::end(rOther.mTimeVectorIntegrationMethods), std::back_inserter(mTimeVectorIntegrationMethods));
129  }
130 
133  {
134  return SolutionSchemePointerType( new SolutionScheme(*this) );
135  }
136 
138  ~SolutionScheme() override {}
139 
143 
147 
153  {
154  KRATOS_TRY
155 
156  if( this->mOptions.IsNotDefined(LocalFlagType::MOVE_MESH) )
157  mOptions.Set(LocalFlagType::MOVE_MESH,true); //default : lagrangian mesh update
158 
159  if( this->mOptions.IsNotDefined(LocalFlagType::UPDATE_VARIABLES) )
160  mOptions.Set(LocalFlagType::UPDATE_VARIABLES,true); //default : derivatives update
161 
162  if( this->mOptions.IsNotDefined(LocalFlagType::INCREMENTAL_SOLUTION) )
163  mOptions.Set(LocalFlagType::INCREMENTAL_SOLUTION,true); //default : dof is the variable increment
164 
165  KRATOS_CATCH("")
166  }
167 
172  virtual void Initialize(ModelPart& rModelPart)
173  {
174  KRATOS_TRY
175 
176  for(typename IntegrationMethodsVectorType::iterator it=mTimeVectorIntegrationMethods.begin();
177  it!=mTimeVectorIntegrationMethods.end(); ++it)
178  (*it)->SetParameters(rModelPart.GetProcessInfo());
179 
180  for(typename IntegrationMethodsScalarType::iterator it=mTimeScalarIntegrationMethods.begin();
181  it!=mTimeScalarIntegrationMethods.end(); ++it)
182  (*it)->SetParameters(rModelPart.GetProcessInfo());
183 
184  this->InitializeElements(rModelPart);
185 
186  this->InitializeConditions(rModelPart);
187 
188  this->Set(LocalFlagType::INITIALIZED, true);
189 
190  KRATOS_CATCH("")
191  }
192 
197  virtual void InitializeSolutionStep(ModelPart& rModelPart)
198  {
199  KRATOS_TRY
200 
201  ProcessInfo& rCurrentProcessInfo = rModelPart.GetProcessInfo();
202 
203  for(typename ProcessPointerVectorType::iterator it=mProcesses.begin(); it!=mProcesses.end(); ++it)
204  (*it)->ExecuteInitializeSolutionStep();
205 
206 #pragma omp parallel for
207  for(int i=0; i<static_cast<int>(rModelPart.Elements().size()); i++)
208  {
209  auto itElem = rModelPart.ElementsBegin() + i;
210  itElem->InitializeSolutionStep(rCurrentProcessInfo);
211  }
212 
213 #pragma omp parallel for
214  for(int i=0; i<static_cast<int>(rModelPart.Conditions().size()); i++)
215  {
216  auto itCond = rModelPart.ConditionsBegin() + i;
217  itCond->InitializeSolutionStep(rCurrentProcessInfo);
218  }
219 
220  KRATOS_CATCH("")
221  }
222 
223 
229  virtual void FinalizeSolutionStep(ModelPart& rModelPart)
230  {
231  KRATOS_TRY
232 
233  ProcessInfo& rCurrentProcessInfo = rModelPart.GetProcessInfo();
234 
235  for(typename ProcessPointerVectorType::iterator it=mProcesses.begin(); it!=mProcesses.end(); ++it)
236  (*it)->ExecuteFinalizeSolutionStep();
237 
238 
239 #pragma omp parallel for
240  for(int i=0; i<static_cast<int>(rModelPart.Elements().size()); i++)
241  {
242  auto itElem = rModelPart.ElementsBegin() + i;
243  itElem->FinalizeSolutionStep(rCurrentProcessInfo);
244  }
245 
246 
247 #pragma omp parallel for
248  for(int i=0; i<static_cast<int>(rModelPart.Conditions().size()); i++)
249  {
250  auto itCond = rModelPart.ConditionsBegin() + i;
251  itCond->FinalizeSolutionStep(rCurrentProcessInfo);
252  }
253 
254  KRATOS_CATCH("")
255  }
256 
261  virtual void InitializeNonLinearIteration(ModelPart& rModelPart)
262  {
263  KRATOS_TRY
264 
265  ProcessInfo& rCurrentProcessInfo = rModelPart.GetProcessInfo();
266 
267  for(typename ProcessPointerVectorType::iterator it=mProcesses.begin(); it!=mProcesses.end(); ++it)
268  (*it)->ExecuteInitialize(); //corresponds to ExecuteInitializeNonLinearIteration()
269 
270 #pragma omp parallel for
271  for(int i=0; i<static_cast<int>(rModelPart.Elements().size()); i++)
272  {
273  auto itElem = rModelPart.ElementsBegin() + i;
274  itElem->InitializeNonLinearIteration(rCurrentProcessInfo);
275  }
276 
277 
278 #pragma omp parallel for
279  for(int i=0; i<static_cast<int>(rModelPart.Conditions().size()); i++)
280  {
281  auto itCond = rModelPart.ConditionsBegin() + i;
282  itCond->InitializeNonLinearIteration(rCurrentProcessInfo);
283  }
284 
285 
286  KRATOS_CATCH("")
287  }
288 
289 
294  virtual void FinalizeNonLinearIteration(ModelPart& rModelPart)
295  {
296  KRATOS_TRY
297 
298  for(typename ProcessPointerVectorType::iterator it=mProcesses.begin(); it!=mProcesses.end(); ++it)
299  (*it)->ExecuteFinalize(); //corresponds to ExecuteFinalizeNonLinearIteration()
300 
301  ProcessInfo& rCurrentProcessInfo = rModelPart.GetProcessInfo();
302 
303 #pragma omp parallel for
304  for(int i=0; i<static_cast<int>(rModelPart.Elements().size()); i++)
305  {
306  auto itElem = rModelPart.ElementsBegin() + i;
307  itElem->FinalizeNonLinearIteration(rCurrentProcessInfo);
308  }
309 
310 
311 #pragma omp parallel for
312  for(int i=0; i<static_cast<int>(rModelPart.Conditions().size()); i++)
313  {
314  auto itCond = rModelPart.ConditionsBegin() + i;
315  itCond->FinalizeNonLinearIteration(rCurrentProcessInfo);
316  }
317 
318  KRATOS_CATCH("")
319  }
320 
321 
326  virtual void Predict(ModelPart& rModelPart,
327  DofsArrayType& rDofSet,
328  SystemVectorType& rDx)
329  {
330  KRATOS_TRY
331 
332  KRATOS_CATCH("")
333  }
334 
339  virtual void Update(ModelPart& rModelPart,
340  DofsArrayType& rDofSet,
341  SystemVectorType& rDx)
342  {
343  KRATOS_TRY
344 
345 
346  KRATOS_CATCH("")
347  }
348 
349 
354  static inline void SetSolution(ModelPart& rModelPart,
355  DofsArrayType& rDofSet,
356  SystemVectorType& rDx)
357  {
358  KRATOS_TRY
359 
360  const unsigned int NumThreads = ParallelUtilities::GetNumThreads();
361 
362  // Update of displacement (by DOF)
363  OpenMPUtils::PartitionVector DofPartition;
364  OpenMPUtils::DivideInPartitions(rDofSet.size(), NumThreads, DofPartition);
365 
366  const int ndof = static_cast<int>(rDofSet.size());
367  typename DofsArrayType::iterator DofBegin = rDofSet.begin();
368 
369  #pragma omp parallel for firstprivate(DofBegin)
370  for(int i = 0; i < ndof; i++)
371  {
372  typename DofsArrayType::iterator itDof = DofBegin + i;
373 
374  if (itDof->IsFree() )
375  {
376  itDof->GetSolutionStepValue() = TSparseSpace::GetValue(rDx,itDof->EquationId());
377  }
378  }
379 
380  KRATOS_CATCH("")
381  }
382 
383 
388  static inline void AddSolution(ModelPart& rModelPart,
389  DofsArrayType& rDofSet,
390  SystemVectorType& rDx)
391  {
392  KRATOS_TRY
393 
394  const unsigned int NumThreads = ParallelUtilities::GetNumThreads();
395 
396  // Update of displacement (by DOF)
397  OpenMPUtils::PartitionVector DofPartition;
398  OpenMPUtils::DivideInPartitions(rDofSet.size(), NumThreads, DofPartition);
399 
400  const int ndof = static_cast<int>(rDofSet.size());
401  typename DofsArrayType::iterator DofBegin = rDofSet.begin();
402 
403  #pragma omp parallel for firstprivate(DofBegin)
404  for(int i = 0; i < ndof; i++)
405  {
406  typename DofsArrayType::iterator itDof = DofBegin + i;
407 
408  if (itDof->IsFree() )
409  {
410  itDof->GetSolutionStepValue() += TSparseSpace::GetValue(rDx,itDof->EquationId());
411  }
412  }
413 
414  KRATOS_CATCH("")
415  }
416 
421  virtual void UpdateDofs(ModelPart& rModelPart,
422  DofsArrayType& rDofSet,
423  SystemVectorType& rDx)
424  {
425  KRATOS_TRY
426 
427  if( mOptions.Is(LocalFlagType::INCREMENTAL_SOLUTION) )
428  AddSolution(rModelPart,rDofSet,rDx); //dof = incremental variable
429  else
430  SetSolution(rModelPart,rDofSet,rDx); //dof = total variable
431 
432  KRATOS_CATCH("")
433  }
434 
435 
436 
441  virtual void UpdateVariables(ModelPart& rModelPart)
442  {
443  KRATOS_TRY
444 
445  if( this->mOptions.Is(LocalFlagType::UPDATE_VARIABLES) ){
446 
447  // Updating time derivatives (nodally for efficiency)
448  const unsigned int NumThreads = ParallelUtilities::GetNumThreads();
449  OpenMPUtils::PartitionVector NodePartition;
450  OpenMPUtils::DivideInPartitions(rModelPart.Nodes().size(), NumThreads, NodePartition);
451 
452  const int nnodes = static_cast<int>(rModelPart.Nodes().size());
453  NodesContainerType::iterator NodeBegin = rModelPart.Nodes().begin();
454 
455  #pragma omp parallel for firstprivate(NodeBegin)
456  for(int i = 0; i < nnodes; i++)
457  {
458  NodesContainerType::iterator itNode = NodeBegin + i;
459 
460  this->IntegrationMethodUpdate(*itNode);
461  }
462  }
463 
464  KRATOS_CATCH("")
465  }
466 
471  virtual void PredictVariables(ModelPart& rModelPart)
472  {
473  KRATOS_TRY
474 
475  // Updating time derivatives (nodally for efficiency)
476  const unsigned int NumThreads = ParallelUtilities::GetNumThreads();
477  OpenMPUtils::PartitionVector NodePartition;
478  OpenMPUtils::DivideInPartitions(rModelPart.Nodes().size(), NumThreads, NodePartition);
479 
480  const int nnodes = static_cast<int>(rModelPart.Nodes().size());
481  NodesContainerType::iterator NodeBegin = rModelPart.Nodes().begin();
482 
483  #pragma omp parallel for firstprivate(NodeBegin)
484  for(int i = 0; i < nnodes; i++)
485  {
486  NodesContainerType::iterator itNode = NodeBegin + i;
487 
488  this->IntegrationMethodPredict(*itNode);
489  }
490 
491  KRATOS_CATCH("")
492  }
493 
494 
499  virtual void MoveMesh(ModelPart& rModelPart)
500  {
501  KRATOS_TRY
502 
503  if( this->mOptions.Is(LocalFlagType::MOVE_MESH) ){
504 
505  if (rModelPart.NodesBegin()->SolutionStepsDataHas(DISPLACEMENT_X) == false)
506  {
507  KRATOS_ERROR << "It is impossible to move the mesh since the DISPLACEMENT variable is not in the Model Part. Add DISPLACEMENT to the list of variables" << std::endl;
508  }
509 
510  bool DisplacementIntegration = false;
511  for(typename IntegrationMethodsVectorType::iterator it=mTimeVectorIntegrationMethods.begin();
512  it!=mTimeVectorIntegrationMethods.end(); ++it)
513  {
514  if( "DISPLACEMENT" == (*it)->GetVariableName() ){
515  DisplacementIntegration = true;
516  break;
517  }
518  }
519 
520  if(DisplacementIntegration == true){
521 
522  // Update mesh positions : node coordinates
523  const int nnodes = rModelPart.NumberOfNodes();
524  ModelPart::NodesContainerType::iterator it_begin = rModelPart.NodesBegin();
525 
526 #pragma omp parallel for
527  for(int i = 0; i<nnodes; i++)
528  {
529  ModelPart::NodesContainerType::iterator it_node = it_begin + i;
530 
531  noalias(it_node->Coordinates()) = it_node->GetInitialPosition().Coordinates();
532  noalias(it_node->Coordinates()) += it_node->FastGetSolutionStepValue(DISPLACEMENT);
533  }
534 
535  }
536 
537  }
538 
539  KRATOS_CATCH("")
540  }
541 
545  virtual void Clear(Element::Pointer rCurrentElement)
546  {
547  }
548 
552  virtual void Clear(Condition::Pointer rCurrentCondition)
553  {
554  }
555 
559  virtual void Clear() {}
560 
565  virtual int Check(ModelPart& rModelPart)
566  {
567  KRATOS_TRY
568 
569  for(ModelPart::ElementsContainerType::iterator it=rModelPart.ElementsBegin();
570  it!=rModelPart.ElementsEnd(); ++it)
571  {
572  it->Check(rModelPart.GetProcessInfo());
573  }
574 
575  for(ModelPart::ConditionsContainerType::iterator it=rModelPart.ConditionsBegin();
576  it!=rModelPart.ConditionsEnd(); ++it)
577  {
578  it->Check(rModelPart.GetProcessInfo());
579  }
580 
581  for(typename IntegrationMethodsVectorType::iterator it=mTimeVectorIntegrationMethods.begin();
582  it!=mTimeVectorIntegrationMethods.end(); ++it)
583  (*it)->Check(rModelPart.GetProcessInfo());
584 
585  for(typename IntegrationMethodsScalarType::iterator it=mTimeScalarIntegrationMethods.begin();
586  it!=mTimeScalarIntegrationMethods.end(); ++it)
587  (*it)->Check(rModelPart.GetProcessInfo());
588 
589  KRATOS_CATCH("")
590 
591  return 0;
592  }
593 
609  virtual void CalculateSystemContributions(Element::Pointer pCurrentElement,
610  LocalSystemMatrixType& rLHS_Contribution,
611  LocalSystemVectorType& rRHS_Contribution,
612  Element::EquationIdVectorType& rEquationId,
613  ProcessInfo& rCurrentProcessInfo)
614  {
615  pCurrentElement->CalculateLocalSystem(rLHS_Contribution, rRHS_Contribution, rCurrentProcessInfo);
616  pCurrentElement->EquationIdVector(rEquationId, rCurrentProcessInfo);
617  }
618 
626  virtual void Calculate_RHS_Contribution(Element::Pointer pCurrentElement,
627  LocalSystemVectorType& rRHS_Contribution,
628  Element::EquationIdVectorType& rEquationId,
629  ProcessInfo& rCurrentProcessInfo)
630  {
631  pCurrentElement->CalculateRightHandSide(rRHS_Contribution, rCurrentProcessInfo);
632  pCurrentElement->EquationIdVector(rEquationId, rCurrentProcessInfo);
633  }
641  virtual void Calculate_LHS_Contribution(Element::Pointer pCurrentElement,
642  LocalSystemMatrixType& rLHS_Contribution,
643  Element::EquationIdVectorType& rEquationId,
644  ProcessInfo& rCurrentProcessInfo)
645  {
646  std::cout<< " it is C_LHS "<<std::endl;
647  pCurrentElement->CalculateLeftHandSide(rLHS_Contribution, rCurrentProcessInfo);
648  pCurrentElement->EquationIdVector(rEquationId, rCurrentProcessInfo);
649  }
650 
651  virtual void EquationId(Element::Pointer pCurrentElement,
652  Element::EquationIdVectorType& rEquationId,
653  ProcessInfo& rCurrentProcessInfo)
654  {
655  (pCurrentElement)->EquationIdVector(rEquationId, rCurrentProcessInfo);
656  }
657 
666  virtual void Condition_CalculateSystemContributions(Condition::Pointer pCurrentCondition,
667  LocalSystemMatrixType& rLHS_Contribution,
668  LocalSystemVectorType& rRHS_Contribution,
669  Element::EquationIdVectorType& rEquationId,
670  ProcessInfo& rCurrentProcessInfo)
671  {
672  pCurrentCondition->CalculateLocalSystem(rLHS_Contribution, rRHS_Contribution, rCurrentProcessInfo);
673  pCurrentCondition->EquationIdVector(rEquationId, rCurrentProcessInfo);
674  }
683  virtual void Condition_Calculate_RHS_Contribution(Condition::Pointer pCurrentCondition,
684  LocalSystemVectorType& rRHS_Contribution,
685  Element::EquationIdVectorType& rEquationId,
686  ProcessInfo& rCurrentProcessInfo)
687  {
688  pCurrentCondition->CalculateRightHandSide(rRHS_Contribution, rCurrentProcessInfo);
689  pCurrentCondition->EquationIdVector(rEquationId, rCurrentProcessInfo);
690  }
691 
699  virtual void Condition_Calculate_LHS_Contribution(Condition::Pointer pCurrentCondition,
700  LocalSystemMatrixType& rLHS_Contribution,
701  Element::EquationIdVectorType& rEquationId,
702  ProcessInfo& rCurrentProcessInfo)
703  {
704  pCurrentCondition->CalculateLeftHandSide(rLHS_Contribution, rCurrentProcessInfo);
705  pCurrentCondition->EquationIdVector(rEquationId, rCurrentProcessInfo);
706  }
707 
708  virtual void Condition_EquationId(Condition::Pointer pCurrentCondition,
709  Element::EquationIdVectorType& rEquationId,
710  ProcessInfo& rCurrentProcessInfo)
711  {
712  (pCurrentCondition)->EquationIdVector(rEquationId, rCurrentProcessInfo);
713  }
714 
721  virtual void GetElementalDofList(Element::Pointer pCurrentElement,
722  Element::DofsVectorType& rElementalDofList,
723  ProcessInfo& rCurrentProcessInfo)
724  {
725  pCurrentElement->GetDofList(rElementalDofList, rCurrentProcessInfo);
726  }
727 
734  virtual void GetConditionDofList(Condition::Pointer pCurrentCondition,
735  Element::DofsVectorType& rConditionDofList,
736  ProcessInfo& rCurrentProcessInfo)
737  {
738  pCurrentCondition->GetDofList(rConditionDofList, rCurrentProcessInfo);
739  }
740 
744 
745 
749  void SetOptions(Flags& rOptions)
750  {
751  mOptions = rOptions;
752  }
753 
759  {
760  return mOptions;
761  }
762 
767  {
768  mProcesses.push_back(pProcess); //NOTE: order set = order of execution
769  }
770 
775  {
776  mProcesses = rProcessVector;
777  }
778 
782 
786 
788  protected:
791 
795 
796  // Flags to set options
798 
799  // Time integration methods
802 
803  // Processes called after move mesh is called
805 
809 
813 
818  virtual void InitializeElements(ModelPart& rModelPart)
819  {
820  KRATOS_TRY
821 
822  if( this->IsNot(LocalFlagType::ELEMENTS_INITIALIZED) ){
823 
824 #pragma omp parallel for
825  for(int i=0; i<static_cast<int>(rModelPart.Elements().size()); i++)
826  {
827  auto itElem = rModelPart.ElementsBegin() + i;
828  itElem->Initialize(rModelPart.GetProcessInfo());
829  }
830 
831  this->Set(LocalFlagType::ELEMENTS_INITIALIZED, true);
832 
833  }
834 
835  KRATOS_CATCH("")
836  }
837 
842  virtual void InitializeConditions(ModelPart& rModelPart)
843  {
844  KRATOS_TRY
845 
846  if( this->IsNot(LocalFlagType::ELEMENTS_INITIALIZED) )
847  KRATOS_ERROR << "Before initilizing Conditions, initialize Elements FIRST" << std::endl;
848 
849  if( this->IsNot(LocalFlagType::CONDITIONS_INITIALIZED) ){
850 
851 #pragma omp parallel for
852  for(int i=0; i<static_cast<int>(rModelPart.Conditions().size()); i++)
853  {
854  auto itCond = rModelPart.ConditionsBegin() + i;
855  itCond->Initialize(rModelPart.GetProcessInfo());
856  }
857 
858  this->Set(LocalFlagType::CONDITIONS_INITIALIZED, true);
859  }
860 
861  KRATOS_CATCH("")
862  }
863 
868  virtual void InitializeNonLinearIteration(Condition::Pointer rCurrentCondition,
869  ProcessInfo& rCurrentProcessInfo)
870  {
871  KRATOS_TRY
872 
873  rCurrentCondition->InitializeNonLinearIteration(rCurrentProcessInfo);
874 
875  KRATOS_CATCH("")
876  }
877 
882  virtual void InitializeNonLinearIteration(Element::Pointer rCurrentElement,
883  ProcessInfo& rCurrentProcessInfo)
884  {
885  KRATOS_TRY
886 
887  rCurrentElement->InitializeNonLinearIteration(rCurrentProcessInfo);
888 
889  KRATOS_CATCH("")
890  }
891 
892 
893  virtual void IntegrationMethodUpdate(NodeType& rNode)
894  {
895  for(typename IntegrationMethodsVectorType::iterator it=mTimeVectorIntegrationMethods.begin();
896  it!=mTimeVectorIntegrationMethods.end(); ++it)
897  (*it)->Update(rNode);
898  for(typename IntegrationMethodsScalarType::iterator it=mTimeScalarIntegrationMethods.begin();
899  it!=mTimeScalarIntegrationMethods.end(); ++it)
900  (*it)->Update(rNode);
901  }
902 
903  virtual void IntegrationMethodPredict(NodeType& rNode)
904  {
905  for(typename IntegrationMethodsVectorType::iterator it=mTimeVectorIntegrationMethods.begin();
906  it!=mTimeVectorIntegrationMethods.end(); ++it)
907  (*it)->Predict(rNode);
908  for(typename IntegrationMethodsScalarType::iterator it=mTimeScalarIntegrationMethods.begin();
909  it!=mTimeScalarIntegrationMethods.end(); ++it)
910  (*it)->Predict(rNode);
911  }
912 
913 
917 
921 
925 
927  private:
930 
934 
938 
942 
946 
950 
954 
958 
960 }; // Class SolutionScheme
962 
965 
966 
970 
971 
973 
975 
976 } // namespace Kratos.
977 
978 #endif // KRATOS_SOLUTION_SCHEME_H_INCLUDED defined
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
Definition: flags.h:58
void Set(const Flags ThisFlag)
Definition: flags.cpp:33
bool Is(Flags const &rOther) const
Definition: flags.h:274
bool IsNot(Flags const &rOther) const
Definition: flags.h:291
bool IsNotDefined(Flags const &rOther) const
Definition: flags.h:296
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
ElementIterator ElementsBegin(IndexType ThisIndex=0)
Definition: model_part.h:1169
MeshType::ConditionsContainerType ConditionsContainerType
Condintions container. A vector set of Conditions with their Id's as key.
Definition: model_part.h:183
ConditionIterator ConditionsBegin(IndexType ThisIndex=0)
Definition: model_part.h:1361
MeshType::ElementsContainerType ElementsContainerType
Element container. A vector set of Elements with their Id's as key.
Definition: model_part.h:168
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
ConditionsContainerType & Conditions(IndexType ThisIndex=0)
Definition: model_part.h:1381
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
ElementsContainerType & Elements(IndexType ThisIndex=0)
Definition: model_part.h:1189
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
ElementIterator ElementsEnd(IndexType ThisIndex=0)
Definition: model_part.h:1179
SizeType NumberOfNodes(IndexType ThisIndex=0) const
Definition: model_part.h:341
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
ConditionIterator ConditionsEnd(IndexType ThisIndex=0)
Definition: model_part.h:1371
This class defines the node.
Definition: node.h:65
static void DivideInPartitions(const int NumTerms, const int NumThreads, PartitionVector &Partitions)
Divide an array of length NumTerms between NumThreads threads.
Definition: openmp_utils.h:158
std::vector< int > PartitionVector
Vector type for the output of DivideInPartitions method.
Definition: openmp_utils.h:53
static int GetNumThreads()
Returns the current number of threads.
Definition: parallel_utilities.cpp:34
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
boost::indirect_iterator< typename TContainerType::iterator > iterator
Definition: pointer_vector_set.h:95
size_type size() const
Returns the number of elements in the container.
Definition: pointer_vector_set.h:502
iterator begin()
Returns an iterator pointing to the beginning of the container.
Definition: pointer_vector_set.h:278
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Solution scheme base class.
Definition: solution_scheme.hpp:54
virtual void InitializeConditions(ModelPart &rModelPart)
Initialize the conditions.
Definition: solution_scheme.hpp:842
virtual void Condition_Calculate_RHS_Contribution(Condition::Pointer pCurrentCondition, LocalSystemVectorType &rRHS_Contribution, Element::EquationIdVectorType &rEquationId, ProcessInfo &rCurrentProcessInfo)
Definition: solution_scheme.hpp:683
SolutionScheme(IntegrationMethodsScalarType &rTimeScalarIntegrationMethods)
Constructor.
Definition: solution_scheme.hpp:111
SolutionScheme(IntegrationMethodsVectorType &rTimeVectorIntegrationMethods, Flags &rOptions)
Constructor.
Definition: solution_scheme.hpp:102
virtual void UpdateDofs(ModelPart &rModelPart, DofsArrayType &rDofSet, SystemVectorType &rDx)
Performing the update of the solution Dofs.
Definition: solution_scheme.hpp:421
IntegrationVectorType::Pointer IntegrationVectorPointerType
Definition: solution_scheme.hpp:77
virtual void InitializeNonLinearIteration(Condition::Pointer rCurrentCondition, ProcessInfo &rCurrentProcessInfo)
Initialize the conditions.
Definition: solution_scheme.hpp:868
std::vector< IntegrationScalarPointerType > IntegrationMethodsScalarType
Definition: solution_scheme.hpp:83
virtual int Check(ModelPart &rModelPart)
This function is designed to be called once to perform all the checks needed.
Definition: solution_scheme.hpp:565
Flags & GetOptions()
Get strategy options.
Definition: solution_scheme.hpp:758
virtual void Clear()
Liberates internal storage.
Definition: solution_scheme.hpp:559
KRATOS_CLASS_POINTER_DEFINITION(SolutionScheme)
Pointer definition of SolutionScheme.
ProcessPointerVectorType mProcesses
Definition: solution_scheme.hpp:804
SolutionScheme(SolutionScheme &rOther)
Copy contructor.
Definition: solution_scheme.hpp:125
virtual void Clear(Condition::Pointer rCurrentCondition)
Liberates internal storage for a condition.
Definition: solution_scheme.hpp:552
array_1d< double, 3 > VectorType
Definition: solution_scheme.hpp:74
virtual void FinalizeNonLinearIteration(ModelPart &rModelPart)
Performs all the required operations that should be done (for each iteration) after solving a solutio...
Definition: solution_scheme.hpp:294
virtual void Predict(ModelPart &rModelPart, DofsArrayType &rDofSet, SystemVectorType &rDx)
Performing the prediction of the solution.
Definition: solution_scheme.hpp:326
virtual void IntegrationMethodPredict(NodeType &rNode)
Definition: solution_scheme.hpp:903
SolutionSchemeType::Pointer SolutionSchemePointerType
Definition: solution_scheme.hpp:60
virtual void MoveMesh(ModelPart &rModelPart)
This function is designed to move the mesh.
Definition: solution_scheme.hpp:499
void SetProcess(ProcessPointerType pProcess)
Set process to execute after move_mesh.
Definition: solution_scheme.hpp:766
void SetOptions(Flags &rOptions)
Sets strategy options.
Definition: solution_scheme.hpp:749
virtual void Calculate_LHS_Contribution(Element::Pointer pCurrentElement, LocalSystemMatrixType &rLHS_Contribution, Element::EquationIdVectorType &rEquationId, ProcessInfo &rCurrentProcessInfo)
Definition: solution_scheme.hpp:641
SolutionScheme(IntegrationMethodsVectorType &rTimeVectorIntegrationMethods)
Constructor.
Definition: solution_scheme.hpp:105
TDenseSpace::MatrixType LocalSystemMatrixType
Definition: solution_scheme.hpp:66
virtual void Initialize(ModelPart &rModelPart)
Performs all the required operations that should be done (for each step) before solving the solution ...
Definition: solution_scheme.hpp:172
void SetProcessVector(ProcessPointerVectorType &rProcessVector)
Set list of processes to execute after move_mesh.
Definition: solution_scheme.hpp:774
ModelPart::ElementsContainerType ElementsContainerType
Definition: solution_scheme.hpp:70
virtual void PredictVariables(ModelPart &rModelPart)
Performing the prediction of the solution variables.
Definition: solution_scheme.hpp:471
virtual void InitializeNonLinearIteration(Element::Pointer rCurrentElement, ProcessInfo &rCurrentProcessInfo)
Initialize the elements.
Definition: solution_scheme.hpp:882
SolutionScheme(IntegrationMethodsScalarType &rTimeScalarIntegrationMethods, Flags &rOptions)
Constructor.
Definition: solution_scheme.hpp:108
TDenseSpace::VectorType LocalSystemVectorType
Definition: solution_scheme.hpp:67
std::vector< ProcessPointerType > ProcessPointerVectorType
Definition: solution_scheme.hpp:86
TSparseSpace::MatrixType SystemMatrixType
Definition: solution_scheme.hpp:64
Variable< VectorType > VariableVectorType
Definition: solution_scheme.hpp:75
SolutionScheme(Flags &rOptions)
Constructor.
Definition: solution_scheme.hpp:99
void SetDefaultFlags()
SetDefaultSchemeFlags.
Definition: solution_scheme.hpp:152
ModelPart::NodeType NodeType
Definition: solution_scheme.hpp:73
virtual void Condition_Calculate_LHS_Contribution(Condition::Pointer pCurrentCondition, LocalSystemMatrixType &rLHS_Contribution, Element::EquationIdVectorType &rEquationId, ProcessInfo &rCurrentProcessInfo)
Definition: solution_scheme.hpp:699
ModelPart::NodesContainerType NodesContainerType
Definition: solution_scheme.hpp:69
SolverProcess::Pointer ProcessPointerType
Definition: solution_scheme.hpp:85
static void SetSolution(ModelPart &rModelPart, DofsArrayType &rDofSet, SystemVectorType &rDx)
Performing the update of the solution Dofs (total solution)
Definition: solution_scheme.hpp:354
virtual void GetElementalDofList(Element::Pointer pCurrentElement, Element::DofsVectorType &rElementalDofList, ProcessInfo &rCurrentProcessInfo)
Definition: solution_scheme.hpp:721
virtual void InitializeElements(ModelPart &rModelPart)
Initialize the elements.
Definition: solution_scheme.hpp:818
ModelPart::DofsArrayType DofsArrayType
Definition: solution_scheme.hpp:62
static void AddSolution(ModelPart &rModelPart, DofsArrayType &rDofSet, SystemVectorType &rDx)
Performing the update of the solution Dofs (incremental solution)
Definition: solution_scheme.hpp:388
TimeIntegrationMethod< VariableVectorType, VectorType > IntegrationVectorType
Definition: solution_scheme.hpp:76
SolutionScheme< TSparseSpace, TDenseSpace > SolutionSchemeType
Definition: solution_scheme.hpp:59
virtual void GetConditionDofList(Condition::Pointer pCurrentCondition, Element::DofsVectorType &rConditionDofList, ProcessInfo &rCurrentProcessInfo)
Definition: solution_scheme.hpp:734
virtual SolutionSchemePointerType Clone()
Clone.
Definition: solution_scheme.hpp:132
virtual void InitializeSolutionStep(ModelPart &rModelPart)
Performs all the required operations that should be done (for each step) before solving the solution ...
Definition: solution_scheme.hpp:197
Flags mOptions
Definition: solution_scheme.hpp:797
virtual void Condition_EquationId(Condition::Pointer pCurrentCondition, Element::EquationIdVectorType &rEquationId, ProcessInfo &rCurrentProcessInfo)
Definition: solution_scheme.hpp:708
IntegrationMethodsVectorType mTimeVectorIntegrationMethods
Definition: solution_scheme.hpp:800
virtual void Condition_CalculateSystemContributions(Condition::Pointer pCurrentCondition, LocalSystemMatrixType &rLHS_Contribution, LocalSystemVectorType &rRHS_Contribution, Element::EquationIdVectorType &rEquationId, ProcessInfo &rCurrentProcessInfo)
Definition: solution_scheme.hpp:666
IntegrationScalarType::Pointer IntegrationScalarPointerType
Definition: solution_scheme.hpp:82
SolutionScheme(IntegrationMethodsVectorType &rTimeVectorIntegrationMethods, IntegrationMethodsScalarType &rTimeScalarIntegrationMethods, Flags &rOptions)
Constructor.
Definition: solution_scheme.hpp:114
IntegrationMethodsScalarType mTimeScalarIntegrationMethods
Definition: solution_scheme.hpp:801
virtual void FinalizeSolutionStep(ModelPart &rModelPart)
Performs all the required operations that should be done (for each step) after solving the solution s...
Definition: solution_scheme.hpp:229
virtual void EquationId(Element::Pointer pCurrentElement, Element::EquationIdVectorType &rEquationId, ProcessInfo &rCurrentProcessInfo)
Definition: solution_scheme.hpp:651
virtual void Update(ModelPart &rModelPart, DofsArrayType &rDofSet, SystemVectorType &rDx)
Performing the update of the solution.
Definition: solution_scheme.hpp:339
TSparseSpace::VectorType SystemVectorType
Definition: solution_scheme.hpp:65
virtual void UpdateVariables(ModelPart &rModelPart)
Performing the update of the solution variables.
Definition: solution_scheme.hpp:441
ModelPart::ConditionsContainerType ConditionsContainerType
Definition: solution_scheme.hpp:71
virtual void IntegrationMethodUpdate(NodeType &rNode)
Definition: solution_scheme.hpp:893
virtual void Calculate_RHS_Contribution(Element::Pointer pCurrentElement, LocalSystemVectorType &rRHS_Contribution, Element::EquationIdVectorType &rEquationId, ProcessInfo &rCurrentProcessInfo)
Definition: solution_scheme.hpp:626
SolutionScheme()
Default Constructor.
Definition: solution_scheme.hpp:96
~SolutionScheme() override
Destructor.
Definition: solution_scheme.hpp:138
SolverLocalFlags LocalFlagType
Definition: solution_scheme.hpp:61
TimeIntegrationMethod< VariableScalarType, double > IntegrationScalarType
Definition: solution_scheme.hpp:81
virtual void CalculateSystemContributions(Element::Pointer pCurrentElement, LocalSystemMatrixType &rLHS_Contribution, LocalSystemVectorType &rRHS_Contribution, Element::EquationIdVectorType &rEquationId, ProcessInfo &rCurrentProcessInfo)
Definition: solution_scheme.hpp:609
virtual void InitializeNonLinearIteration(ModelPart &rModelPart)
Performs all the required operations that should be done (for each iteration) before solving a soluti...
Definition: solution_scheme.hpp:261
Variable< double > VariableScalarType
Definition: solution_scheme.hpp:80
SolutionScheme(IntegrationMethodsVectorType &rTimeVectorIntegrationMethods, IntegrationMethodsScalarType &rTimeScalarIntegrationMethods)
Constructor.
Definition: solution_scheme.hpp:120
std::vector< IntegrationVectorPointerType > IntegrationMethodsVectorType
Definition: solution_scheme.hpp:78
virtual void Clear(Element::Pointer rCurrentElement)
Liberates internal storage for an element.
Definition: solution_scheme.hpp:545
Solver local flags class definition.
Definition: solution_local_flags.hpp:48
Short class definition.
Definition: time_integration_method.hpp:55
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR
Definition: exception.h:161
end
Definition: DEM_benchmarks.py:180
Vector VectorType
Definition: geometrical_transformation_utilities.h:56
Matrix MatrixType
Definition: geometrical_transformation_utilities.h:55
Parameters GetValue(Parameters &rParameters, const std::string &rEntry)
Definition: add_kratos_parameters_to_python.cpp:53
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
int nnodes
Definition: sensitivityMatrix.py:24
integer i
Definition: TensorModule.f:17