v0.14.0
Functions | Variables
cloasest_point_projection.cpp File Reference
#include <MoFEM.hpp>
#include <boost/program_options.hpp>
#include <boost/numeric/ublas/vector_proxy.hpp>
#include <boost/numeric/ublas/matrix.hpp>
#include <boost/numeric/ublas/matrix_proxy.hpp>
#include <boost/numeric/ublas/vector.hpp>
#include <boost/numeric/ublas/symmetric.hpp>
#include <adolc/adolc.h>
#include <SmallStrainPlasticity.hpp>
#include <SmallStrainPlasticityMaterialModels.hpp>

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Variables

static char help []
 

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 44 of file cloasest_point_projection.cpp.

44  {
45 
46  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
47 
48  try {
49 
50  moab::Core mb_instance;
51  moab::Interface &moab = mb_instance;
52  MoFEM::Core core(moab);
53  MoFEM::Interface &m_field = core;
54 
56 
57  cp.tAgs.resize(3);
58  cp.tAgs[0] = 0;
59  cp.tAgs[1] = 1;
60  cp.tAgs[2] = 2;
61 
62  const double young_modulus = 1;
63  const double poisson_ratio = 0.0;
66  cp.sIgma_y = 1;
67  cp.H = 0.1;
68  cp.K = 0;
69  cp.phi = 1;
70 
71  cp.sTrain.resize(6, false);
72  cp.sTrain.clear();
73  cp.plasticStrain.resize(6, false);
74  cp.plasticStrain.clear();
75  cp.internalVariables.resize(7, false);
76  cp.internalVariables.clear();
77 
78  cp.createMatAVecR();
79  cp.snesCreate();
80  CHKERR SNESSetFromOptions(cp.sNes);
81 
82  const int idx = 3;
83  double alpha = 0.1;
84  int N = 40;
85  for (int ii = 0; ii < N; ii++) {
86  cout << "Step " << ii << endl;
87 
88  cp.gG = ii;
89 
90  if (ii <= 14) {
91  cp.sTrain[idx] += +alpha;
92  cp.sTrain[0] += +alpha;
93  } else if (ii <= 100) {
94  cp.sTrain[idx] += -alpha;
95  cp.sTrain[0] += -alpha;
96  } else {
97  cp.sTrain[idx] += +alpha;
98  }
99 
102  CHKERR VecSetValue(cp.Chi, 6 + cp.internalVariables.size(), 0,
103  INSERT_VALUES);
104  CHKERR VecAssemblyBegin(cp.Chi);
105  CHKERR VecAssemblyEnd(cp.Chi);
106  CHKERR SmallStrainPlasticityfRes(cp.sNes, cp.Chi, cp.R, &cp);
107 
108  double tol = 1e-8;
109  bool nonlinear = false;
110  if (cp.y > tol) {
111  nonlinear = true;
113  }
114 
115  cout << "plot " << cp.sTrain[idx] << " " << cp.sTress[idx] << " "
116  << cp.internalVariables[0] << " " << cp.plasticStrain[idx] << " "
117  << cp.deltaGamma << endl;
118 
119  if (nonlinear) {
121  cerr << endl;
122  cerr << "Cep " << cp.Cep << endl;
123  MatrixDouble e = cp.Cep;
124  VectorDouble strain0 = cp.sTrain;
125  const double eps = 1e-6;
126  for (int ii = 0; ii < 6; ii++) {
127  cp.sTrain = strain0;
128  cp.sTrain[ii] -= eps;
130  for (int dd = 0; dd < 6; dd++) {
131  e(dd, ii) += cp.sTress[dd] / (2 * eps);
132  }
133  cp.sTrain = strain0;
134  cp.sTrain[ii] += eps;
136  for (int dd = 0; dd < 6; dd++) {
137  e(dd, ii) -= cp.sTress[dd] / (2 * eps);
138  }
139  }
140  cp.sTrain = strain0;
141  cerr << "e " << e << endl;
142  }
143 
144  // std::string wait;
145  // std::cin >> wait;
146  }
147 
148  cp.destroyMatAVecR();
149  cp.snesDestroy();
150  }
151  CATCH_ERRORS;
152 
154 
155  return 0;
156 }

Variable Documentation

◆ help

char help[]
static
Initial value:
= "...\n"
"\n"

Definition at line 41 of file cloasest_point_projection.cpp.

SmallStrainPlasticity::ClosestPointProjection::plasticStrain0
VectorDouble plasticStrain0
Definition: SmallStrainPlasticity.hpp:488
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
SmallStrainPlasticity::ClosestPointProjection::internalVariables0
VectorDouble internalVariables0
Definition: SmallStrainPlasticity.hpp:487
young_modulus
double young_modulus
Young modulus.
Definition: plastic.cpp:172
SmallStrainJ2Plasticity::sIgma_y
double sIgma_y
Definition: SmallStrainPlasticityMaterialModels.hpp:45
SmallStrainPlasticity::ClosestPointProjection::destroyMatAVecR
virtual PetscErrorCode destroyMatAVecR()
Definition: SmallStrainPlasticity.cpp:936
SmallStrainPlasticity::ClosestPointProjection::plasticStrain
VectorDouble plasticStrain
Definition: SmallStrainPlasticity.hpp:489
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
SmallStrainPlasticity::ClosestPointProjection::R
Vec R
Definition: SmallStrainPlasticity.hpp:542
SmallStrainPlasticity::ClosestPointProjection::sTrain
VectorDouble sTrain
Definition: SmallStrainPlasticity.hpp:485
LAMBDA
#define LAMBDA(E, NU)
Definition: fem_tools.h:22
SmallStrainPlasticity::ClosestPointProjection::solveColasetProjection
PetscErrorCode solveColasetProjection()
Definition: SmallStrainPlasticity.cpp:1050
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
SmallStrainPlasticity::ClosestPointProjection::tAgs
vector< int > tAgs
Definition: SmallStrainPlasticity.hpp:506
SmallStrainPlasticity::ClosestPointProjection::snesCreate
PetscErrorCode snesCreate()
Definition: SmallStrainPlasticity.cpp:1027
SmallStrainPlasticity::ClosestPointProjection::gG
int gG
Definition: SmallStrainPlasticity.hpp:492
SmallStrainPlasticity::ClosestPointProjection::createMatAVecR
virtual PetscErrorCode createMatAVecR()
Definition: SmallStrainPlasticity.cpp:924
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
SmallStrainPlasticityfRes
PetscErrorCode SmallStrainPlasticityfRes(SNES snes, Vec chi, Vec r, void *ctx)
Definition: SmallStrainPlasticity.cpp:1154
SmallStrainPlasticity::ClosestPointProjection::sNes
SNES sNes
Definition: SmallStrainPlasticity.hpp:556
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
SmallStrainPlasticity::ClosestPointProjection::snesDestroy
PetscErrorCode snesDestroy()
Definition: SmallStrainPlasticity.cpp:1044
SmallStrainJ2Plasticity::phi
double phi
Definition: SmallStrainPlasticityMaterialModels.hpp:44
SmallStrainPlasticity::ClosestPointProjection::Chi
Vec Chi
Definition: SmallStrainPlasticity.hpp:542
SmallStrainPlasticity::ClosestPointProjection::y
double y
Definition: SmallStrainPlasticity.hpp:516
poisson_ratio
double poisson_ratio
Poisson ratio.
Definition: plastic.cpp:173
SmallStrainJ2Plasticity
J2 plasticity (Kinematic Isotropic (Linear) Hardening)
Definition: SmallStrainPlasticityMaterialModels.hpp:32
SmallStrainJ2Plasticity::mu
double mu
Definition: SmallStrainPlasticityMaterialModels.hpp:40
SmallStrainPlasticity::ClosestPointProjection::sTress
VectorAdaptor sTress
Definition: SmallStrainPlasticity.hpp:494
SmallStrainJ2Plasticity::lambda
double lambda
Definition: SmallStrainPlasticityMaterialModels.hpp:41
SmallStrainPlasticity::ClosestPointProjection::internalVariables
VectorDouble internalVariables
Definition: SmallStrainPlasticity.hpp:486
N
const int N
Definition: speed_test.cpp:3
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
FTensor::dd
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
SmallStrainJ2Plasticity::H
double H
Definition: SmallStrainPlasticityMaterialModels.hpp:42
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
SmallStrainJ2Plasticity::K
double K
Definition: SmallStrainPlasticityMaterialModels.hpp:43
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
SmallStrainPlasticity::ClosestPointProjection::consistentTangent
PetscErrorCode consistentTangent()
Definition: SmallStrainPlasticity.cpp:1095
SmallStrainPlasticity::ClosestPointProjection::deltaGamma
double deltaGamma
Definition: SmallStrainPlasticity.hpp:490
help
static char help[]
Definition: cloasest_point_projection.cpp:41
tol
double tol
Definition: mesh_smoothing.cpp:27
SmallStrainPlasticity::ClosestPointProjection::Cep
MatrixDouble Cep
Definition: SmallStrainPlasticity.hpp:566
MU
#define MU(E, NU)
Definition: fem_tools.h:23