Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Algorithm.hpp
Go to the documentation of this file.
1 //==============================================================================
2 //
3 // Seismo Virtual Laboratory
4 // Module for Serial and Parallel Analysis of seismic
5 // wave propagation and soil-structure interaction simulation
6 // Copyright (C) 2018-2021, The California Institute of Technology
7 // All Rights Reserved.
8 //
9 // Commercial use of this program without express permission of the California
10 // Institute of Technology, is strictly prohibited. See file "COPYRIGHT" in
11 // main directory for information on usage and redistribution, and for a
12 // DISCLAIMER OF ALL WARRANTIES.
13 //
14 //==============================================================================
15 //
16 // Written by:
17 // Danilo S. Kusanovic (dkusanov@caltech.edu)
18 // Elnaz E. Seylabi (elnaze@unr.edu)
19 //
20 // Supervised by:
21 // Domniki M. Asimaki (domniki@caltech.edu)
22 //
23 // References :
24 // [1]
25 //
26 // Description:
29 //------------------------------------------------------------------------------
30 
31 #ifndef _ALGORITHM_HPP_
32 #define _ALGORITHM_HPP_
33 
34 #include <map>
35 #include <mpi.h>
36 #include <memory>
37 #include <Eigen/Dense>
38 #include <Eigen/SparseCore>
39 
40 #include "Mesh.hpp"
41 #include "LoadCombo.hpp"
42 
43 class Integrator;
44 
52 class Algorithm{
53 
54  public:
60  Algorithm(const std::shared_ptr<Mesh> &mesh, unsigned int flag=1, double NormFactor=1.0);
61 
63  virtual ~Algorithm() = 0;
64 
69  virtual bool ComputeNewIncrement(std::shared_ptr<Mesh> &mesh, unsigned int i) = 0;
70 
73  virtual const Eigen::VectorXd& GetDisplacementIncrement() = 0;
74 
77  virtual void SetIntegrator(std::shared_ptr<Integrator> &integrator) = 0;
78 
81  virtual void SetLoadFactor(double factor) = 0;
82 
87  void ReverseStatesIncrements(std::shared_ptr<Mesh> &mesh);
88 
94  void UpdateStatesIncrements(std::shared_ptr<Mesh> &mesh, const Eigen::VectorXd &dU);
95 
100  void ReducedParallelResidual(const Eigen::VectorXd &Feff, double &Residual);
101 
115  double ComputeConvergence(const Eigen::VectorXd &Force, const Eigen::VectorXd &Delta, const Eigen::VectorXd &delta, unsigned int k);
116 
117  protected:
120  Eigen::SparseMatrix<double> Total2FreeMatrix;
121 
122  private:
124  unsigned int flag;
125 
127  double NormFactor;
128 };
129 
130 #endif
Virtual class for defining how the solution between two steps is going to be carried out...
Definition: Algorithm.hpp:52
void ReducedParallelResidual(const Eigen::VectorXd &Feff, double &Residual)
Construct the residual vector force from each processor.
unsigned int flag
Convergence test to be performed.
Definition: Algorithm.hpp:124
void ReverseStatesIncrements(std::shared_ptr< Mesh > &mesh)
Reverse the incremental state variables to its pervious converged one in the mesh.
This file contains the "Load Combination" class declarations, which defines how the loads are going t...
virtual ~Algorithm()=0
Destroys this Algorithm object.
Class for defining how the static/dynamic analysis between two steps is going to be performed...
Definition: Integrator.hpp:53
Eigen::SparseMatrix< double > Total2FreeMatrix
Operator that enforced restrain/constraint.
Definition: Algorithm.hpp:120
void UpdateStatesIncrements(std::shared_ptr< Mesh > &mesh, const Eigen::VectorXd &dU)
Update the incremental state variables in the mesh.
This file contains the "Mesh object" declarations, which stores nodes, materials, elements...
Algorithm(const std::shared_ptr< Mesh > &mesh, unsigned int flag=1, double NormFactor=1.0)
Creates a Algorithm object.
virtual const Eigen::VectorXd & GetDisplacementIncrement()=0
Gets the displacement increment.
virtual void SetLoadFactor(double factor)=0
Set the load factor.
double ComputeConvergence(const Eigen::VectorXd &Force, const Eigen::VectorXd &Delta, const Eigen::VectorXd &delta, unsigned int k)
Computes convergence tests for this algorithm.
virtual void SetIntegrator(std::shared_ptr< Integrator > &integrator)=0
Sets the integrator for this algorithm.
virtual bool ComputeNewIncrement(std::shared_ptr< Mesh > &mesh, unsigned int i)=0
Computes a new incremental solution.
double NormFactor
Relative Unbalanced Factor.
Definition: Algorithm.hpp:127