v0.14.0
Loading...
Searching...
No Matches
closest_point_projection.cpp
Go to the documentation of this file.
1/** \file small_strain_plasticity.cpp
2 * \ingroup small_strain_plasticity
3 * \brief Small strain plasticity example
4 *
5 */
6
7/* This file is part of MoFEM.
8 * MoFEM is free software: you can redistribute it and/or modify it under
9 * the terms of the GNU Lesser General Public License as published by the
10 * Free Software Foundation, either version 3 of the License, or (at your
11 * option) any later version.
12 *
13 * MoFEM is distributed in the hope that it will be useful, but WITHOUT
14 * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
15 * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
16 * License for more details.
17 *
18 * You should have received a copy of the GNU Lesser General Public
19 * License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
20
21#include <MoFEM.hpp>
22using namespace MoFEM;
23
24#include <boost/program_options.hpp>
25using namespace std;
26namespace po = boost::program_options;
27
28#include <ADOLCPlasticity.hpp>
30using namespace ADOLCPlasticity;
31
32using namespace boost::numeric;
33
34static char help[] = "...\n"
35 "\n";
36
37int main(int argc, char *argv[]) {
38
39 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
40
41 try {
42
43 moab::Core mb_instance;
44 moab::Interface &moab = mb_instance;
45 MoFEM::Core core(moab);
46 MoFEM::Interface &m_field = core;
47
49
50 VectorDouble active_variables(12 + cp_ptr->nbInternalVariables);
51 VectorDouble grad_w(12 + cp_ptr->nbInternalVariables);
52
53 auto tmp_active = getVectorAdaptor(&(active_variables[0]),
54 12 + cp_ptr->nbInternalVariables);
55 cp_ptr->activeVariablesW.swap(tmp_active);
56 auto tmp_gradient =
57 getVectorAdaptor(&(grad_w[0]), 12 + cp_ptr->nbInternalVariables);
58 cp_ptr->gradientW.swap(tmp_gradient);
59
60 auto strain = cp_ptr->getTotalStrain();
61 auto internal_variables = cp_ptr->getInternalVariables();
62 auto stress = cp_ptr->getStress();
63 auto plastic_strain = cp_ptr->getPlasticStrain();
64
65 cp_ptr->createMatAVecR();
66 cp_ptr->snesCreate();
67
68 const int idx = 3;
69 double alpha = 0.1;
70 int N = 40;
71 for (int ii = 0; ii < N; ii++) {
72 cout << "Step " << ii << endl;
73
74 if (ii <= 14) {
75 strain[idx] += +alpha;
76 strain[0] += +alpha;
77 } else if (ii <= 100) {
78 strain[idx] += -alpha;
79 strain[0] += -alpha;
80 } else {
81 strain[idx] += +alpha;
82 }
83
84 cp_ptr->plasticStrain0 = plastic_strain;
85 cp_ptr->internalVariables0 = internal_variables;
86 CHKERR VecSetValue(cp_ptr->Chi, 6 + cp_ptr->nbInternalVariables, 0,
87 INSERT_VALUES);
88 CHKERR VecAssemblyBegin(cp_ptr->Chi);
89 CHKERR VecAssemblyEnd(cp_ptr->Chi);
90 CHKERR ADOLCPlasticityRes(cp_ptr->sNes, cp_ptr->Chi, cp_ptr->R,
91 cp_ptr.get());
92
93 double tol = 1e-8;
94 bool nonlinear = false;
95 if (cp_ptr->y > tol) {
96 nonlinear = true;
97 CHKERR cp_ptr->solveClosestProjection();
98 }
99
100 cout << "plot " << strain[idx] << " " << stress[idx] << " "
101 << internal_variables[0] << " " << plastic_strain[idx] << " "
102 << cp_ptr->deltaGamma << endl;
103
104 if (nonlinear) {
105 CHKERR cp_ptr->consistentTangent();
106 cerr << endl;
107 cerr << "Cep " << cp_ptr->Cep << endl;
108 MatrixDouble e = cp_ptr->Cep;
109 VectorDouble strain0 = strain;
110 const double eps = 1e-6;
111 for (int ii = 0; ii < 6; ii++) {
112 strain = strain0;
113 strain[ii] -= eps;
114 CHKERR cp_ptr->solveClosestProjection();
115 for (int dd = 0; dd < 6; dd++) {
116 e(dd, ii) += stress[dd] / (2 * eps);
117 }
118 strain = strain0;
119 strain[ii] += eps;
120 CHKERR cp_ptr->solveClosestProjection();
121 for (int dd = 0; dd < 6; dd++) {
122 e(dd, ii) -= stress[dd] / (2 * eps);
123 }
124 }
125 strain = strain0;
126 cerr << "e " << e << endl;
127 }
128
129 // std::string wait;
130 // std::cin >> wait;
131 }
132 }
134
136
137 return 0;
138}
Matetial models for plasticity.
int main()
static const double eps
static char help[]
#define CATCH_ERRORS
Catch errors.
#define CHKERR
Inline error check.
double tol
auto createMaterial(std::array< int, ClosestPointProjection::LAST_TAPE > tapes_tags={0, 1, 2})
MoFEMErrorCode ADOLCPlasticityRes(SNES snes, Vec chi, Vec r, void *ctx)
Internal SNES function used at integration points to calulate stress.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
auto getVectorAdaptor(T1 ptr, const size_t n)
Get Vector adaptor.
Definition Templates.hpp:31
const int N
Definition speed_test.cpp:3
Core (interface) class.
Definition Core.hpp:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:112
Deprecated interface functions.