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.
relaxed_residualbased_newton_rapshon_strategy.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ \.
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Joaquin Gonzalez-Usua
11 //
12 
13 #if !defined(KRATOS_RELAXED_RESIDUALBASED_NEWTON_RAPHSON_STRATEGY)
14 #define KRATOS_RELAXED_RESIDUALBASED_NEWTON_RAPHSON_STRATEGY
15 
16 // System includes
17 #include <iostream>
18 
19 // External includes
20 
21 // Project includes
22 #include "includes/define.h"
28 
29 //default builder and solver
31 
32 namespace Kratos
33 {
34 
37 
41 
43 
46 
50 
54 
62 template <class TSparseSpace,
63  class TDenseSpace, // = DenseSpace<double>,
64  class TLinearSolver //= LinearSolver<TSparseSpace,TDenseSpace>
65  >
67  : public ResidualBasedNewtonRaphsonStrategy<TSparseSpace, TDenseSpace, TLinearSolver>
68 {
69  public:
73 
74  // Counted pointer of ClassName
76 
78 
80 
82 
84 
86 
87  typedef typename BaseType::TDataType TDataType;
88 
89  typedef TSparseSpace SparseSpaceType;
90 
92 
93  //typedef typename BaseType::DofSetType DofSetType;
94 
96 
98 
100 
102 
104 
106 
108 
112 
113  // /**
114  // * @brief Default constructor
115  // */
117  {
118  }
119 
120  // /**
121  // * @brief Default constructor. (with parameters)
122  // * @param rModelPart The model part of the problem
123  // * @param ThisParameters The configuration parameters
124  // */
127  {
128  }
129 
135  explicit RelaxedResidualBasedNewtonRaphsonStrategy(ModelPart& rModelPart, Parameters ThisParameters)
136  : BaseType(rModelPart, ThisParameters),
140  {
141  }
142 
155  ModelPart& rModelPart,
156  typename TSchemeType::Pointer pScheme,
157  typename TLinearSolver::Pointer pNewLinearSolver,
158  typename TConvergenceCriteriaType::Pointer pNewConvergenceCriteria,
159  int MaxIterations = 30,
160  bool CalculateReactions = false,
161  bool ReformDofSetAtEachStep = false,
162  bool MoveMeshFlag = false)
163  : BaseType(rModelPart, pScheme, pNewLinearSolver, pNewConvergenceCriteria, MaxIterations, CalculateReactions, ReformDofSetAtEachStep, MoveMeshFlag),
164  mpScheme(pScheme),
165  mpConvergenceCriteria(pNewConvergenceCriteria),
166  mReformDofSetAtEachStep(ReformDofSetAtEachStep),
167  mCalculateReactionsFlag(CalculateReactions),
169  mMaxIterationNumber(MaxIterations),
172  {
173  }
174 
187  ModelPart& rModelPart,
188  typename TSchemeType::Pointer pScheme,
189  typename TConvergenceCriteriaType::Pointer pNewConvergenceCriteria,
190  typename TBuilderAndSolverType::Pointer pNewBuilderAndSolver,
191  int MaxIterations = 30,
192  bool CalculateReactions = false,
193  bool ReformDofSetAtEachStep = false,
194  bool MoveMeshFlag = false)
195  : BaseType(rModelPart, pScheme, pNewConvergenceCriteria, pNewBuilderAndSolver, MaxIterations, CalculateReactions, ReformDofSetAtEachStep, MoveMeshFlag),
196  mpScheme(pScheme),
197  mpBuilderAndSolver(pNewBuilderAndSolver),
198  mpConvergenceCriteria(pNewConvergenceCriteria),
199  mReformDofSetAtEachStep(ReformDofSetAtEachStep),
200  mCalculateReactionsFlag(CalculateReactions),
202  mMaxIterationNumber(MaxIterations),
205  {
206  KRATOS_TRY
207  // Getting builder and solver
208  auto p_builder_and_solver = GetBuilderAndSolver();
209 
210  // Tells to the builder and solver if the reactions have to be Calculated or not
211  p_builder_and_solver->SetCalculateReactionsFlag(mCalculateReactionsFlag);
212 
213  // Tells to the Builder And Solver if the system matrix and vectors need to
214  //be reshaped at each step or not
215  p_builder_and_solver->SetReshapeMatrixFlag(mReformDofSetAtEachStep);
216 
217  // Set EchoLevel to the default value (only time is displayed)
218  SetEchoLevel(1);
219 
220  // By default the matrices are rebuilt at each iteration
221  this->SetRebuildLevel(2);
222 
226  KRATOS_CATCH("")
227  }
228 
229 
239  ModelPart& rModelPart,
240  typename TSchemeType::Pointer pScheme,
241  typename TLinearSolver::Pointer pNewLinearSolver,
242  typename TConvergenceCriteriaType::Pointer pNewConvergenceCriteria,
243  Parameters Settings)
244  : BaseType(rModelPart, pScheme, pNewLinearSolver, pNewConvergenceCriteria, Settings),
245  mpScheme(pScheme),
246  mpConvergenceCriteria(pNewConvergenceCriteria),
250  {
251  }
252 
263  ModelPart& rModelPart,
264  typename TSchemeType::Pointer pScheme,
265  typename TConvergenceCriteriaType::Pointer pNewConvergenceCriteria,
266  typename TBuilderAndSolverType::Pointer pNewBuilderAndSolver,
267  Parameters Settings)
268  : BaseType(rModelPart, pScheme, pNewConvergenceCriteria, pNewBuilderAndSolver, Settings),
269  mpScheme(pScheme),
270  mpBuilderAndSolver(pNewBuilderAndSolver),
271  mpConvergenceCriteria(pNewConvergenceCriteria),
275  {
276  }
277 
283  {
284  // If the linear solver has not been deallocated, clean it before
285  // deallocating mpA. This prevents a memory error with the the ML
286  // solver (which holds a reference to it).
287  // NOTE: The linear solver is hold by the B&S
288  auto p_builder_and_solver = this->GetBuilderAndSolver();
289  if (p_builder_and_solver != nullptr) {
290  p_builder_and_solver->Clear();
291  }
292 
293  // Deallocating system vectors to avoid errors in MPI. Clear calls
294  // TrilinosSpace::Clear for the vectors, which preserves the Map of
295  // current vectors, performing MPI calls in the process. Due to the
296  // way Python garbage collection works, this may happen after
297  // MPI_Finalize has already been called and is an error. Resetting
298  // the pointers here prevents Clear from operating with the
299  // (now deallocated) vectors.
300  mpA.reset();
301  mpDx.reset();
302  mpb.reset();
303 
304  Clear();
305  }
310  void SetScheme(typename TSchemeType::Pointer pScheme)
311  {
312  BaseType::SetScheme(pScheme);
313  };
314 
315  // /**
316  // * @brief Get method for the time scheme
317  // * @return mpScheme: The pointer to the time scheme considered
318  // */
319  typename TSchemeType::Pointer GetScheme()
320  {
321  return mpScheme;
322  };
323 
324  void SetBuilderAndSolver(typename TBuilderAndSolverType::Pointer pNewBuilderAndSolver)
325  {
326  mpBuilderAndSolver = pNewBuilderAndSolver;
327  };
328 
333  typename TBuilderAndSolverType::Pointer GetBuilderAndSolver()
334  {
335  return mpBuilderAndSolver;
336  };
337 
338 //*********************************************************************************
340  void SetInitializePerformedFlag(bool InitializePerformedFlag = true)
341  {
342  mInitializeWasPerformed = InitializePerformedFlag;
343  }
344 
350  {
352  }
353 
358  void SetCalculateReactionsFlag(bool CalculateReactionsFlag)
359  {
360  mCalculateReactionsFlag = CalculateReactionsFlag;
361  }
362 
368  {
370  }
371 
376  void SetUseOldStiffnessInFirstIterationFlag(bool UseOldStiffnessInFirstIterationFlag)
377  {
378  mUseOldStiffnessInFirstIteration = UseOldStiffnessInFirstIterationFlag;
379  }
380 
381 
382 
384  {
386  }
387 
393  {
395  GetBuilderAndSolver()->SetReshapeMatrixFlag(mReformDofSetAtEachStep);
396  }
397 
403  {
405  }
406 
411  void SetMaxIterationNumber(unsigned int MaxIterationNumber)
412  {
413  mMaxIterationNumber = MaxIterationNumber;
414  }
415 
420  unsigned int GetMaxIterationNumber()
421  {
422  return mMaxIterationNumber;
423  }
424 
425  void SetEchoLevel(int Level) override
426  {
427  BaseType::SetEchoLevel(Level);
428  }
434  typename SolvingStrategyType::Pointer Create(
435  ModelPart& rModelPart,
436  Parameters ThisParameters
437  ) const override
438  {
439  return Kratos::make_shared<ClassType>(rModelPart, ThisParameters);
440  }
444  void Clear() override
445  {
446  BaseType::Clear();
447  }
452  void InitializeSolutionStep() override
453  {
454  KRATOS_TRY;
455 
457  // Pointers needed in the solution
458  typename TSchemeType::Pointer p_scheme = GetScheme();
459  typename TBuilderAndSolverType::Pointer p_builder_and_solver = GetBuilderAndSolver();
460  ModelPart& r_model_part = BaseType::GetModelPart();
461 
462  //set up the system, operation performed just once unless it is required
463  //to reform the dof set at each iteration
464  BuiltinTimer system_construction_time;
465  if (p_builder_and_solver->GetDofSetIsInitializedFlag() == false ||
466  mReformDofSetAtEachStep == true)
467  {
468  //setting up the list of the DOFs to be solved
469  BuiltinTimer setup_dofs_time;
470  p_builder_and_solver->SetUpDofSet(p_scheme, r_model_part);
471  KRATOS_INFO_IF("ResidualBasedNewtonRaphsonStrategy", BaseType::GetEchoLevel() > 0) << "Setup Dofs Time: "
472  << setup_dofs_time.ElapsedSeconds() << std::endl;
473 
474  //shaping correctly the system
475  BuiltinTimer setup_system_time;
476  p_builder_and_solver->SetUpSystem(r_model_part);
477  KRATOS_INFO_IF("ResidualBasedNewtonRaphsonStrategy", BaseType::GetEchoLevel() > 0) << "Setup System Time: "
478  << setup_system_time.ElapsedSeconds() << std::endl;
479 
480  //setting up the Vectors involved to the correct size
481  BuiltinTimer system_matrix_resize_time;
482  p_builder_and_solver->ResizeAndInitializeVectors(p_scheme, mpA, mpDx, mpb,
483  r_model_part);
484  KRATOS_INFO_IF("ResidualBasedNewtonRaphsonStrategy", BaseType::GetEchoLevel() > 0) << "System Matrix Resize Time: "
485  << system_matrix_resize_time.ElapsedSeconds() << std::endl;
486  }
487 
488  KRATOS_INFO_IF("ResidualBasedNewtonRaphsonStrategy", BaseType::GetEchoLevel() > 0) << "System Construction Time: "
489  << system_construction_time.ElapsedSeconds() << std::endl;
490 
491  TSystemMatrixType& rA = *mpA;
492  TSystemVectorType& rDx = *mpDx;
493  TSystemVectorType& rb = *mpb;
494 
495  // Initial operations ... things that are constant over the Solution Step
496  p_builder_and_solver->InitializeSolutionStep(r_model_part, rA, rDx, rb);
497 
498  // Initial operations ... things that are constant over the Solution Step
499  p_scheme->InitializeSolutionStep(r_model_part, rA, rDx, rb);
500 
501  // Initialisation of the convergence criteria
502  if (mpConvergenceCriteria->GetActualizeRHSflag() == true)
503  {
504  TSparseSpace::SetToZero(rb);
505  p_builder_and_solver->BuildRHS(p_scheme, r_model_part, rb);
506  }
507 
508  mpConvergenceCriteria->InitializeSolutionStep(r_model_part, p_builder_and_solver->GetDofSet(), rA, rDx, rb);
509 
510  if (mpConvergenceCriteria->GetActualizeRHSflag() == true)
511  TSparseSpace::SetToZero(rb);
512 
514  }
515 
516  KRATOS_CATCH("");
517  }
518 
523  void FinalizeSolutionStep() override
524  {
525  KRATOS_TRY;
526 
527  ModelPart& r_model_part = BaseType::GetModelPart();
528 
529  typename TSchemeType::Pointer p_scheme = GetScheme();
530  typename TBuilderAndSolverType::Pointer p_builder_and_solver = GetBuilderAndSolver();
531 
532  TSystemMatrixType& rA = *mpA;
533  TSystemVectorType& rDx = *mpDx;
534  TSystemVectorType& rb = *mpb;
535 
536  //Finalisation of the solution step,
537  //operations to be done after achieving convergence, for example the
538  //Final Residual Vector (mb) has to be saved in there
539  //to avoid error accumulation
540 
541  p_scheme->FinalizeSolutionStep(r_model_part, rA, rDx, rb);
542  p_builder_and_solver->FinalizeSolutionStep(r_model_part, rA, rDx, rb);
543  mpConvergenceCriteria->FinalizeSolutionStep(r_model_part, p_builder_and_solver->GetDofSet(), rA, rDx, rb);
544 
545  //Cleaning memory after the solution
546  p_scheme->Clean();
547 
548  //reset flags for next step
550 
551  if (mReformDofSetAtEachStep == true) //deallocate the systemvectors
552  {
553  this->Clear();
554  }
555 
556  KRATOS_CATCH("");
557  }
558 
562  bool SolveSolutionStep() override
563  {
564  // Pointers needed in the solution
565  ModelPart& r_model_part = BaseType::GetModelPart();
566  typename TSchemeType::Pointer p_scheme = GetScheme();
567  typename TBuilderAndSolverType::Pointer p_builder_and_solver = GetBuilderAndSolver();
568  auto& r_dof_set = p_builder_and_solver->GetDofSet();
569 
570  TSystemMatrixType& rA = *mpA;
571  TSystemVectorType& rDx = *mpDx;
572  TSystemVectorType& rb = *mpb;
573  double alpha = r_model_part.GetProcessInfo()[RELAXATION_ALPHA];
574  //initializing the parameters of the Newton-Raphson cycle
575  unsigned int iteration_number = 1;
576  r_model_part.GetProcessInfo()[NL_ITERATION_NUMBER] = iteration_number;
577  bool residual_is_updated = false;
578  p_scheme->InitializeNonLinIteration(r_model_part, rA, rDx, rb);
579  mpConvergenceCriteria->InitializeNonLinearIteration(r_model_part, r_dof_set, rA, rDx, rb);
580  bool is_converged = mpConvergenceCriteria->PreCriteria(r_model_part, r_dof_set, rA, rDx, rb);
581 
582  // Function to perform the building and the solving phase.
584  TSparseSpace::SetToZero(rA);
585  TSparseSpace::SetToZero(rDx);
586  TSparseSpace::SetToZero(rb);
587 
589  p_builder_and_solver->BuildAndSolveLinearizedOnPreviousIteration(p_scheme, r_model_part, rA, rDx, rb,BaseType::MoveMeshFlag());
590  } else {
591  p_builder_and_solver->BuildAndSolve(p_scheme, r_model_part, rA, rDx, rb);
592  }
593  } else {
594  TSparseSpace::SetToZero(rDx); // Dx = 0.00;
595  TSparseSpace::SetToZero(rb);
596  p_builder_and_solver->BuildRHSAndSolve(p_scheme, r_model_part, rA, rDx, rb);
597  }
598  // Debugging info
599  EchoInfo(iteration_number);
600  TSparseSpace::InplaceMult(rDx, alpha);
601  // Updating the results stored in the database
602  UpdateDatabase(rA, rDx, rb, BaseType::MoveMeshFlag());
603  p_scheme->FinalizeNonLinIteration(r_model_part, rA, rDx, rb);
604  mpConvergenceCriteria->FinalizeNonLinearIteration(r_model_part, r_dof_set, rA, rDx, rb);
605 
606  if (is_converged) {
607  if (mpConvergenceCriteria->GetActualizeRHSflag()) {
608  TSparseSpace::SetToZero(rb);
609 
610  p_builder_and_solver->BuildRHS(p_scheme, r_model_part, rb);
611  }
612 
613  is_converged = mpConvergenceCriteria->PostCriteria(r_model_part, r_dof_set, rA, rDx, rb);
614  }
615 
616  //Iteration Cycle... performed only for NonLinearProblems
617  while (is_converged == false &&
618  iteration_number++ < mMaxIterationNumber)
619  {
620  //setting the number of iteration
621  r_model_part.GetProcessInfo()[NL_ITERATION_NUMBER] = iteration_number;
622 
623  p_scheme->InitializeNonLinIteration(r_model_part, rA, rDx, rb);
624  mpConvergenceCriteria->InitializeNonLinearIteration(r_model_part, r_dof_set, rA, rDx, rb);
625 
626  is_converged = mpConvergenceCriteria->PreCriteria(r_model_part, r_dof_set, rA, rDx, rb);
627 
628  //call the linear system solver to find the correction mDx for the
629  //it is not called if there is no system to solve
630  if (SparseSpaceType::Size(rDx) != 0)
631  {
633  {
635  {
636  //A = 0.00;
637  TSparseSpace::SetToZero(rA);
638  TSparseSpace::SetToZero(rDx);
639  TSparseSpace::SetToZero(rb);
640 
641  p_builder_and_solver->BuildAndSolve(p_scheme, r_model_part, rA, rDx, rb);
642  TSparseSpace::InplaceMult(rDx, alpha);
643  }
644  else
645  {
646  TSparseSpace::SetToZero(rDx);
647  TSparseSpace::SetToZero(rb);
648 
649  p_builder_and_solver->BuildRHSAndSolve(p_scheme, r_model_part, rA, rDx, rb);
650  }
651  }
652  else
653  {
654  TSparseSpace::SetToZero(rDx);
655  TSparseSpace::SetToZero(rb);
656 
657  p_builder_and_solver->BuildRHSAndSolve(p_scheme, r_model_part, rA, rDx, rb);
658  }
659  }
660  else
661  {
662  KRATOS_WARNING("NO DOFS") << "ATTENTION: no free DOFs!! " << std::endl;
663  }
664 
665  // Debugging info
666  EchoInfo(iteration_number);
667 
668  // Updating the results stored in the database
669  UpdateDatabase(rA, rDx, rb, BaseType::MoveMeshFlag());
670 
671  p_scheme->FinalizeNonLinIteration(r_model_part, rA, rDx, rb);
672  mpConvergenceCriteria->FinalizeNonLinearIteration(r_model_part, r_dof_set, rA, rDx, rb);
673 
674  residual_is_updated = false;
675 
676  if (is_converged == true)
677  {
678  if (mpConvergenceCriteria->GetActualizeRHSflag() == true)
679  {
680  TSparseSpace::SetToZero(rb);
681 
682  p_builder_and_solver->BuildRHS(p_scheme, r_model_part, rb);
683  residual_is_updated = true;
684  }
685 
686  is_converged = mpConvergenceCriteria->PostCriteria(r_model_part, r_dof_set, rA, rDx, rb);
687  }
688  }
689 
690  //plots a warning if the maximum number of iterations is exceeded
691  if (iteration_number >= mMaxIterationNumber) {
693  } else {
694  KRATOS_INFO_IF("RelaxedResidualBasedNewtonRaphsonStrategy", this->GetEchoLevel() > 0)
695  << "Convergence achieved after " << iteration_number << " / "
696  << mMaxIterationNumber << " iterations" << std::endl;
697  }
698 
699  //recalculate residual if needed
700  //(note that some convergence criteria need it to be recalculated)
701  if (residual_is_updated == false)
702  {
703  // NOTE:
704  // The following part will be commented because it is time consuming
705  // and there is no obvious reason to be here. If someone need this
706  // part please notify the community via mailing list before uncommenting it.
707  // Pooyan.
708 
709  // TSparseSpace::SetToZero(mb);
710  // p_builder_and_solver->BuildRHS(p_scheme, r_model_part, mb);
711  }
712 
713  //calculate reactions if required
714  if (mCalculateReactionsFlag == true)
715  p_builder_and_solver->CalculateReactions(p_scheme, r_model_part, rA, rDx, rb);
716 
717  return is_converged;
718  }
719 
724  static std::string Name()
725  {
726  return "relaxed_newton_raphson_strategy";
727  }
728 
731 
733 
737 
740 
742 
747  {
748  TSystemMatrixType &mA = *mpA;
749 
750  return mA;
751  }
752 
758  {
760 
761  return mb;
762  }
763 
769  {
770  TSystemVectorType& mDx = *mpDx;
771 
772  return mDx;
773  }
774 
780  {
782  }
783 
789  {
791  }
792 
796 
800 
802  std::string Info() const override
803  {
804  return "RelaxedResidualBasedNewtonRaphsonStrategy";
805  }
806 
808  void PrintInfo(std::ostream& rOStream) const override
809  {
810  rOStream << Info();
811  }
812 
814  void PrintData(std::ostream& rOStream) const override
815  {
816  rOStream << Info();
817  }
818 
822 
824 
825  private:
828 
832 
836 
840 
844 
848 
852 
854 
855  protected:
858 
862 
863  typename TSchemeType::Pointer mpScheme = nullptr;
864  typename TBuilderAndSolverType::Pointer mpBuilderAndSolver = nullptr;
865  typename TConvergenceCriteriaType::Pointer mpConvergenceCriteria = nullptr;
866 
870 
879 
885 
891 
893 
894  unsigned int mMaxIterationNumber;
895 
897 
898  bool mKeepSystemConstantDuringIterations; // Flag to allow keeping system matrix constant during iterations
899 
900  virtual void UpdateDatabase(
901  TSystemMatrixType& rA,
902  TSystemVectorType& rDx,
903  TSystemVectorType& rb,
904  const bool MoveMesh) override
905  {
907  }
908 
909 
914  virtual void EchoInfo(const unsigned int IterationNumber) override
915  {
916  BaseType::EchoInfo(IterationNumber);
917  }
918 
919  virtual void MaxIterationsExceeded() override
920  {
921  KRATOS_INFO_IF("RelaxedResidualBasedNewtonRaphsonStrategy", this->GetEchoLevel() > 0)
922  << "ATTENTION: max iterations ( " << mMaxIterationNumber
923  << " ) exceeded!" << std::endl;
924  }
925 
930  void AssignSettings(const Parameters ThisParameters) override
931  {
932  BaseType::AssignSettings(ThisParameters);
933  }
934 
935  void WriteDofInfo(std::string FileName, const TSystemVectorType& rDX)
936  {
938  }
939  /*
940  * @brief Flag telling if it is needed to reform the DofSet at each
941  solution step or if it is possible to form it just once
942  * @details Default = false
943  - true : Reforme at each time step
944  - false : Form just one (more efficient)
945  */
946 
950 
951 
955 
959 
963 
967 
968  /*
969  * Copy constructor.
970  */
971 
973 
975 
976 }; /* Class RelaxedResidualBasedNewtonRaphsonStrategy */
977 
979 
982 
984 
985 } /* namespace Kratos. */
986 
987 #endif /* KRATOS_RELAXED_RESIDUALBASED_NEWTON_RAPHSON_STRATEGY defined */
Definition: builtin_timer.h:26
double ElapsedSeconds() const
Definition: builtin_timer.h:40
This is the base class to define the different convergence criterion considered.
Definition: convergence_criteria.h:58
Implicit solving strategy base class This is the base class from which we will derive all the implici...
Definition: implicit_solving_strategy.h:61
BaseType::TSystemVectorType TSystemVectorType
Definition: implicit_solving_strategy.h:72
void SetRebuildLevel(int Level) override
This sets the build level.
Definition: implicit_solving_strategy.h:207
int mRebuildLevel
Definition: implicit_solving_strategy.h:263
bool mStiffnessMatrixIsBuilt
The current rebuild level.
Definition: implicit_solving_strategy.h:264
BaseType::TSystemVectorPointerType TSystemVectorPointerType
Definition: implicit_solving_strategy.h:76
BaseType::TSystemMatrixPointerType TSystemMatrixPointerType
Definition: implicit_solving_strategy.h:74
BaseType::TSystemMatrixType TSystemMatrixType
Definition: implicit_solving_strategy.h:70
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
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
This is the base Newton Raphson strategy.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:68
RelaxedResidualBasedNewtonRaphsonStrategy(ModelPart &rModelPart, typename TSchemeType::Pointer pScheme, typename TLinearSolver::Pointer pNewLinearSolver, typename TConvergenceCriteriaType::Pointer pNewConvergenceCriteria, int MaxIterations=30, bool CalculateReactions=false, bool ReformDofSetAtEachStep=false, bool MoveMeshFlag=false)
Definition: relaxed_residualbased_newton_rapshon_strategy.h:154
bool GetUseOldStiffnessInFirstIterationFlag()
Definition: relaxed_residualbased_newton_rapshon_strategy.h:383
void SetUseOldStiffnessInFirstIterationFlag(bool UseOldStiffnessInFirstIterationFlag)
This method sets the flag mFullUpdateFlag.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:376
RelaxedResidualBasedNewtonRaphsonStrategy(const RelaxedResidualBasedNewtonRaphsonStrategy &Other)
Definition: relaxed_residualbased_newton_rapshon_strategy.h:972
bool mKeepSystemConstantDuringIterations
Flag to set as initialized the strategy.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:898
BaseType::TSystemMatrixPointerType TSystemMatrixPointerType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:105
void FinalizeSolutionStep() override
Performs all the required operations that should be done (for each step) after solving the solution s...
Definition: relaxed_residualbased_newton_rapshon_strategy.h:523
void SetInitializePerformedFlag(bool InitializePerformedFlag=true)
Definition: relaxed_residualbased_newton_rapshon_strategy.h:340
bool mSolutionStepIsInitialized
Definition: relaxed_residualbased_newton_rapshon_strategy.h:892
TSystemVectorPointerType mpDx
The pointer to the convergence criteria employed.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:867
BaseType::LocalSystemVectorType LocalSystemVectorType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:101
TSystemVectorType & GetSystemVector() override
This method returns the RHS vector.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:757
SolvingStrategyType::Pointer Create(ModelPart &rModelPart, Parameters ThisParameters) const override
Create method.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:434
TConvergenceCriteriaType::Pointer mpConvergenceCriteria
The pointer to the builder and solver employed.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:865
static std::string Name()
Returns the name of the class as used in the settings (snake_case format)
Definition: relaxed_residualbased_newton_rapshon_strategy.h:724
TSchemeType::Pointer mpScheme
Definition: relaxed_residualbased_newton_rapshon_strategy.h:863
RelaxedResidualBasedNewtonRaphsonStrategy(ModelPart &rModelPart, Parameters ThisParameters)
Default constructor. (with parameters)
Definition: relaxed_residualbased_newton_rapshon_strategy.h:135
void InitializeSolutionStep() override
Performs all the required operations that should be done (for each step) before solving the solution ...
Definition: relaxed_residualbased_newton_rapshon_strategy.h:452
BaseType::TSystemVectorType TSystemVectorType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:99
bool mInitializeWasPerformed
The maximum number of iterations, 30 by default.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:896
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:814
RelaxedResidualBasedNewtonRaphsonStrategy(ModelPart &rModelPart)
Definition: relaxed_residualbased_newton_rapshon_strategy.h:125
BaseType::TSystemMatrixType TSystemMatrixType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:97
TBuilderAndSolverType::Pointer GetBuilderAndSolver()
Get method for the builder and solver.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:333
bool mReformDofSetAtEachStep
The LHS matrix of the system of equations.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:878
KRATOS_CLASS_POINTER_DEFINITION(RelaxedResidualBasedNewtonRaphsonStrategy)
void AssignSettings(const Parameters ThisParameters) override
This method assigns settings to member variables.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:930
ResidualBasedNewtonRaphsonStrategy< TSparseSpace, TDenseSpace, TLinearSolver > BaseType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:81
ImplicitSolvingStrategy< TSparseSpace, TDenseSpace, TLinearSolver > BaseStrategyType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:79
BaseType::TBuilderAndSolverType TBuilderAndSolverType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:85
TSystemVectorType & GetSolutionVector() override
This method returns the solution vector.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:768
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:808
void SetScheme(typename TSchemeType::Pointer pScheme)
Set method for the time scheme.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:310
bool GetReformDofSetAtEachStepFlag()
This method returns the flag mReformDofSetAtEachStep.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:402
unsigned int GetMaxIterationNumber()
This method gets the flag mMaxIterationNumber.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:420
TSparseSpace SparseSpaceType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:89
RelaxedResidualBasedNewtonRaphsonStrategy< TSparseSpace, TDenseSpace, TLinearSolver > ClassType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:83
bool GetCalculateReactionsFlag()
This method returns the flag mCalculateReactionsFlag.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:367
TSystemMatrixType & GetSystemMatrix() override
This method returns the LHS matrix.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:746
void SetCalculateReactionsFlag(bool CalculateReactionsFlag)
This method sets the flag mCalculateReactionsFlag.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:358
BaseType::TSchemeType TSchemeType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:91
TSystemMatrixPointerType mpA
The RHS vector of the system of equations.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:869
void SetEchoLevel(int Level) override
It sets the level of echo for the solving strategy.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:425
BaseType::LocalSystemMatrixType LocalSystemMatrixType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:103
~RelaxedResidualBasedNewtonRaphsonStrategy() override
Destructor.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:282
RelaxedResidualBasedNewtonRaphsonStrategy()
Definition: relaxed_residualbased_newton_rapshon_strategy.h:116
virtual void UpdateDatabase(TSystemMatrixType &rA, TSystemVectorType &rDx, TSystemVectorType &rb, const bool MoveMesh) override
Here the database is updated.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:900
TBuilderAndSolverType::Pointer mpBuilderAndSolver
The pointer to the time scheme employed.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:864
void SetKeepSystemConstantDuringIterations(bool Value)
Set method for the flag mKeepSystemConstantDuringIterations.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:779
bool mCalculateReactionsFlag
Flag telling if it is needed or not to compute the reactions.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:884
ConvergenceCriteria< TSparseSpace, TDenseSpace > TConvergenceCriteriaType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:72
bool SolveSolutionStep() override
Solves the current step. This function returns true if a solution has been found, false otherwise.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:562
RelaxedResidualBasedNewtonRaphsonStrategy(ModelPart &rModelPart, typename TSchemeType::Pointer pScheme, typename TConvergenceCriteriaType::Pointer pNewConvergenceCriteria, typename TBuilderAndSolverType::Pointer pNewBuilderAndSolver, int MaxIterations=30, bool CalculateReactions=false, bool ReformDofSetAtEachStep=false, bool MoveMeshFlag=false)
Constructor specifying the builder and solver.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:186
TSchemeType::Pointer GetScheme()
Definition: relaxed_residualbased_newton_rapshon_strategy.h:319
void SetReformDofSetAtEachStepFlag(bool Flag)
This method sets the flag mReformDofSetAtEachStep.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:392
TSystemVectorPointerType mpb
The increment in the solution.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:868
bool GetKeepSystemConstantDuringIterations()
Get method for the flag mKeepSystemConstantDuringIterations.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:788
BaseType::TSystemVectorPointerType TSystemVectorPointerType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:107
std::string Info() const override
Turn back information as a string.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:802
bool mUseOldStiffnessInFirstIteration
Flag telling if a full update of the database will be performed at the first iteration.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:890
RelaxedResidualBasedNewtonRaphsonStrategy(ModelPart &rModelPart, typename TSchemeType::Pointer pScheme, typename TConvergenceCriteriaType::Pointer pNewConvergenceCriteria, typename TBuilderAndSolverType::Pointer pNewBuilderAndSolver, Parameters Settings)
Constructor specifying the builder and solver and using Parameters.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:262
SolvingStrategy< TSparseSpace, TDenseSpace > SolvingStrategyType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:77
void Clear() override
Clears the internal storage.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:444
unsigned int mMaxIterationNumber
Flag to set as initialized the solution step.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:894
virtual void MaxIterationsExceeded() override
This method prints information after reach the max number of iterations.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:919
BaseType::DofsArrayType DofsArrayType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:95
bool GetInitializePerformedFlag()
This method gets the flag mInitializeWasPerformed.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:349
void SetMaxIterationNumber(unsigned int MaxIterationNumber)
This method sets the flag mMaxIterationNumber.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:411
void SetBuilderAndSolver(typename TBuilderAndSolverType::Pointer pNewBuilderAndSolver)
Definition: relaxed_residualbased_newton_rapshon_strategy.h:324
virtual void EchoInfo(const unsigned int IterationNumber) override
This method returns the components of the system of equations depending of the echo level.
Definition: relaxed_residualbased_newton_rapshon_strategy.h:914
BaseType::TDataType TDataType
Definition: relaxed_residualbased_newton_rapshon_strategy.h:87
void WriteDofInfo(std::string FileName, const TSystemVectorType &rDX)
Definition: relaxed_residualbased_newton_rapshon_strategy.h:935
RelaxedResidualBasedNewtonRaphsonStrategy(ModelPart &rModelPart, typename TSchemeType::Pointer pScheme, typename TLinearSolver::Pointer pNewLinearSolver, typename TConvergenceCriteriaType::Pointer pNewConvergenceCriteria, Parameters Settings)
Definition: relaxed_residualbased_newton_rapshon_strategy.h:238
This is the base Newton Raphson strategy.
Definition: residualbased_newton_raphson_strategy.h:66
virtual void UpdateDatabase(TSystemMatrixType &rA, TSystemVectorType &rDx, TSystemVectorType &rb, const bool MoveMesh)
Here the database is updated.
Definition: residualbased_newton_raphson_strategy.h:1278
Parameters GetDefaultParameters() const override
This method provides the defaults parameters to avoid conflicts between the different constructors.
Definition: residualbased_newton_raphson_strategy.h:1069
BaseType::TSchemeType TSchemeType
Definition: residualbased_newton_raphson_strategy.h:87
void AssignSettings(const Parameters ThisParameters) override
This method assigns settings to member variables.
Definition: residualbased_newton_raphson_strategy.h:1351
BaseType::TBuilderAndSolverType TBuilderAndSolverType
Definition: residualbased_newton_raphson_strategy.h:81
void SetEchoLevel(int Level) override
It sets the level of echo for the solving strategy.
Definition: residualbased_newton_raphson_strategy.h:599
void Clear() override
Clears the internal storage.
Definition: residualbased_newton_raphson_strategy.h:707
virtual void EchoInfo(const unsigned int IterationNumber)
This method returns the components of the system of equations depending of the echo level.
Definition: residualbased_newton_raphson_strategy.h:1298
void WriteDofInfo(std::string FileName, const TSystemVectorType &rDX)
Definition: residualbased_newton_raphson_strategy.h:1375
void SetScheme(typename TSchemeType::Pointer pScheme)
Set method for the time scheme.
Definition: residualbased_newton_raphson_strategy.h:466
Solving strategy base class This is the base class from which we will derive all the strategies (impl...
Definition: solving_strategy.h:64
TSparseSpace::VectorPointerType TSystemVectorPointerType
Definition: solving_strategy.h:77
TDenseSpace::MatrixType LocalSystemMatrixType
Definition: solving_strategy.h:79
TSparseSpace::DataType TDataType
Definition: solving_strategy.h:69
ModelPart & GetModelPart()
Operations to get the pointer to the model.
Definition: solving_strategy.h:350
int GetEchoLevel()
This returns the level of echo for the solving strategy.
Definition: solving_strategy.h:271
TSparseSpace::MatrixType TSystemMatrixType
Definition: solving_strategy.h:71
TSparseSpace::MatrixPointerType TSystemMatrixPointerType
Definition: solving_strategy.h:75
bool MoveMeshFlag()
This function returns the flag that says if the mesh is moved.
Definition: solving_strategy.h:290
TSparseSpace::VectorType TSystemVectorType
Definition: solving_strategy.h:73
virtual void MoveMesh()
This function is designed to move the mesh.
Definition: solving_strategy.h:330
TDenseSpace::VectorType LocalSystemVectorType
Definition: solving_strategy.h:81
static IndexType Size(VectorType const &rV)
return size of vector rV
Definition: ublas_space.h:190
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_INFO_IF(label, conditional)
Definition: logger.h:251
#define KRATOS_WARNING(label)
Definition: logger.h:265
string FileName
Export to vtk.
Definition: GenerateWind.py:175
TSpaceType::VectorPointerType CreateEmptyVectorPointer(TSpaceType &dummy)
Definition: add_strategies_to_python.cpp:187
TSpaceType::MatrixPointerType CreateEmptyMatrixPointer(TSpaceType &dummy)
Definition: add_strategies_to_python.cpp:181
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
alpha
Definition: generate_convection_diffusion_explicit_element.py:113
double precision function mb(a)
Definition: TensorModule.f:849