v0.14.0
Public Types | Public Member Functions | Public Attributes | Friends | List of all members
ADOLCPlasticity::ClosestPointProjection Struct Referenceabstract

Closest point projection algorithm. More...

#include <users_modules/adolc-plasticity/src/ADOLCPlasticity.hpp>

Inheritance diagram for ADOLCPlasticity::ClosestPointProjection:
[legend]
Collaboration diagram for ADOLCPlasticity::ClosestPointProjection:
[legend]

Public Types

enum  TypesTags { W = 0, Y, H, LAST_TAPE }
 

Public Member Functions

VectorAdaptor getPlasticStrain ()
 
VectorAdaptor getTotalStrain ()
 
VectorAdaptor getInternalVariables ()
 
VectorAdaptor getActiveVariablesYH ()
 
VectorAdaptor getStress ()
 
VectorAdaptor getInternalFluxes ()
 
 ClosestPointProjection ()
 
MoFEMErrorCode recordW ()
 Record strain energy. More...
 
MoFEMErrorCode recordY ()
 Record yield function. More...
 
MoFEMErrorCode recordH ()
 Record flow potential. More...
 
MoFEMErrorCode playW ()
 
MoFEMErrorCode playW_NoHessian ()
 
MoFEMErrorCode playY ()
 
MoFEMErrorCode playY_NoGradient ()
 
MoFEMErrorCode playH ()
 
MoFEMErrorCode playH_NoHessian ()
 
MoFEMErrorCode createMatAVecR ()
 
MoFEMErrorCode evaluatePotentials ()
 
MoFEMErrorCode playPotentials ()
 
MoFEMErrorCode playPotentials_NoHessian ()
 
MoFEMErrorCode calculateR (Vec R)
 
MoFEMErrorCode calculateA ()
 
MoFEMErrorCode snesCreate ()
 Create nested snes. More...
 
MoFEMErrorCode solveClosetProjection ()
 Solve nonlinear system of equations to find stress, internal fluxes, and Lagrange plastic multiplier. More...
 
MoFEMErrorCode consistentTangent ()
 Calculate consistent tangent matrix. More...
 
MoFEMErrorCode recordTapes ()
 Record tapes. More...
 
virtual MoFEMErrorCode addMatBlockOps (MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string block_name, Sev sev)
 Get model parameters from blocks. More...
 
virtual MoFEMErrorCode setParams (short tag, bool &recalculate_elastic_tangent)
 Set parameters for ADOL-C of constitutive relations. More...
 
virtual MoFEMErrorCode freeHemholtzEnergy ()=0
 Set Hemholtz energy. More...
 
virtual MoFEMErrorCode yieldFunction ()=0
 Set yield function. More...
 
virtual MoFEMErrorCode flowPotential ()=0
 Set flow potential. More...
 

Public Attributes

int nbInternalVariables
 
boost::function< int(int, int, int)> integrationRule
 
VectorDouble internalVariables0
 
VectorDouble plasticStrain0
 
std::array< int, LAST_TAPEtapesTags
 
VectorAdaptor activeVariablesW
 
VectorAdaptor gradientW
 
double w
 
double y
 
double h
 
double deltaGamma
 
MatrixDouble Ep
 
MatrixDouble Cp
 
MatrixDouble Cep
 
ublas::symmetric_matrix< double, ublas::lower > C
 
ublas::symmetric_matrix< double, ublas::lower > D
 
MatrixDouble hessianW
 
VectorDouble gradientY
 
VectorDouble gradientH
 
MatrixDouble hessianH
 
MatrixDouble partialWStrainPlasticStrain
 
VectorAdaptor partialYSigma
 
VectorAdaptor partialYFlux
 
VectorAdaptor partialHSigma
 
VectorAdaptor partialHFlux
 
ublas::symmetric_matrix< double, ublas::lower > partial2HSigma
 
ublas::symmetric_matrix< double, ublas::lower > partial2HFlux
 
MatrixDouble partial2HSigmaFlux
 
SmartPetscObj< Mat > A
 
SmartPetscObj< Vec > R
 
SmartPetscObj< Vec > Chi
 
SmartPetscObj< Vec > dChi
 
ublas::matrix< double, ublas::column_major > dataA
 
SmartPetscObj< SNES > sNes
 
ublas::vector< adoublea_sTrain
 
ublas::vector< adoublea_plasticStrain
 
ublas::vector< adoublea_internalVariables
 
ublas::vector< adoublea_sTress
 
ublas::vector< adoublea_internalFluxes
 
adouble a_w
 
adouble a_y
 
adouble a_h
 

Friends

MoFEMErrorCode ADOLCPlasticityRes (SNES, Vec, Vec, void *ctx)
 Function executed by nested SENES at evaluationg residual. More...
 
MoFEMErrorCode ADOLCPlasticityJac (SNES, Vec, Mat, Mat, void *ctx)
 Function executed by nested SENES at evaluationg tangent matrix. More...
 

Detailed Description

Closest point projection algorithm.

Definition at line 139 of file ADOLCPlasticity.hpp.

Member Enumeration Documentation

◆ TypesTags

Enumerator
LAST_TAPE 
Examples
ADOLCPlasticity.hpp.

Definition at line 168 of file ADOLCPlasticity.hpp.

168 { W = 0, Y, H, LAST_TAPE }; //< W - energy, Y - yield, H - flow

Constructor & Destructor Documentation

◆ ClosestPointProjection()

ADOLCPlasticity::ClosestPointProjection::ClosestPointProjection ( )
Examples
ADOLCPlasticity.hpp, and ADOLCPlasticityMaterialModels.hpp.

Definition at line 15 of file ClosetPointProjection.cpp.

15  {
16 
17  if (!LogManager::checkIfChannelExist("ADOLCPlasticityWold")) {
18  auto core_log = logging::core::get();
19 
20  core_log->add_sink(LogManager::createSink(LogManager::getStrmWorld(),
21  "ADOLCPlasticityWold"));
22  core_log->add_sink(LogManager::createSink(LogManager::getStrmSync(),
23  "ADOLCPlasticitySync"));
24  core_log->add_sink(LogManager::createSink(LogManager::getStrmSelf(),
25  "ADOLCPlasticitySelf"));
26 
27  LogManager::setLog("ADOLCPlasticityWold");
28  LogManager::setLog("ADOLCPlasticitySync");
29  LogManager::setLog("ADOLCPlasticitySelf");
30 
31  MOFEM_LOG_TAG("ADOLCPlasticityWold", "ADOL-C Plasticity");
32  MOFEM_LOG_TAG("ADOLCPlasticitySync", "ADOL-C Plasticity");
33  MOFEM_LOG_TAG("ADOLCPlasticitySelf", "ADOL-C Plasticity");
34  }
35 
36  tapesTags = {0, 1, 2}; //< set tapes default tags
37 
38 }

Member Function Documentation

◆ addMatBlockOps()

virtual MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::addMatBlockOps ( MoFEM::Interface m_field,
boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &  pip,
std::string  block_name,
Sev  sev 
)
inlinevirtual

Get model parameters from blocks.

Reimplemented in ADOLCPlasticity::ParaboloidalPlasticity, and ADOLCPlasticity::J2Plasticity< 3 >.

Examples
ADOLCPlasticity.hpp.

Definition at line 261 of file ADOLCPlasticity.hpp.

263  {
264  return 0;
265  }

◆ calculateA()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::calculateA ( )
Examples
ADOLCPlasticity.hpp.

Definition at line 382 of file ClosetPointProjection.cpp.

382  {
384  // matrix A
385  CHKERR MatZeroEntries(A);
386  // row 0 (Strain)
388  MatrixDouble a01 = prod(partial2HSigmaFlux, D);
389  for (int ii = 0; ii < 6; ii++) {
390  for (int jj = 0; jj < 6; jj++) {
391  dataA(ii, jj) = deltaGamma * a00(ii, jj);
392  }
393  for (unsigned int jj = 0; jj < nbInternalVariables; jj++) {
394  dataA(ii, 6 + jj) = deltaGamma * a01(ii, jj);
395  }
396  dataA(ii, ii) -= 1;
397  }
398  for (int ii = 0; ii < 6; ii++) {
399  dataA(ii, 6 + nbInternalVariables) = partialHSigma[ii];
400  }
401  // row 1 (Fluxes)
402  MatrixDouble a11 = prod(partial2HFlux, D);
403  MatrixDouble a10 =
405  for (unsigned int ii = 0; ii < nbInternalVariables; ii++) {
406  for (int jj = 0; jj < 6; jj++) {
407  dataA(6 + ii, jj) += deltaGamma * a10(ii, jj);
408  }
409  for (unsigned int jj = 0; jj < nbInternalVariables; jj++) {
410  dataA(6 + ii, 6 + jj) += deltaGamma * a11(ii, jj);
411  }
412  dataA(6 + ii, 6 + ii) -= 1;
413  }
414  for (unsigned int ii = 0; ii < nbInternalVariables; ii++) {
415  dataA(6 + ii, 6 + nbInternalVariables) = partialHFlux[ii];
416  }
417  // row 3 (Plastic multiplier)
418  VectorDouble partial_y_sigma_c =
420  VectorDouble partial_y_flux_d = prod(trans(D), partialYFlux);
421  for (unsigned int dd = 0; dd < partial_y_sigma_c.size(); dd++) {
422  dataA(6 + nbInternalVariables, dd) = partial_y_sigma_c[dd];
423  }
424  for (unsigned int dd = 0; dd < partial_y_flux_d.size(); dd++) {
425  dataA(6 + nbInternalVariables, 6 + dd) = partial_y_flux_d[dd];
426  }
428  CHKERR MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
429  CHKERR MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
431 }

◆ calculateR()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::calculateR ( Vec  R)
Examples
ADOLCPlasticity.hpp.

Definition at line 360 of file ClosetPointProjection.cpp.

360  {
362 
363  auto plastic_strain = getPlasticStrain();
364  auto internal_variables = getInternalVariables();
365 
366  // Residual
367  double *array;
368  CHKERR VecGetArray(R, &array);
369  for (int ii = 0; ii < 6; ii++) {
370  array[ii] = -plastic_strain[ii] + plasticStrain0[ii] +
372  }
373  for (unsigned int ii = 0; ii < nbInternalVariables; ii++) {
374  array[6 + ii] = -internal_variables[ii] + internalVariables0[ii] +
375  deltaGamma * partialHFlux[ii];
376  }
377  array[6 + nbInternalVariables] = y;
378  CHKERR VecRestoreArray(R, &array);
380 }

◆ consistentTangent()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::consistentTangent ( )

Calculate consistent tangent matrix.

Examples
ADOLCPlasticity.hpp.

Definition at line 509 of file ClosetPointProjection.cpp.

509  {
512  CHKERR calculateA();
513  Ep.resize(6, 6, false);
514  MatrixDouble ep_row; // inrement internal varaibles as result of incremen of
515  // total strain
516  ep_row = deltaGamma * prod(partial2HSigma, C);
517  MatrixDouble alpha_row;
518  alpha_row = deltaGamma * prod(trans(partial2HSigmaFlux), C);
519  VectorDouble y_row;
520  y_row = prod(C, partialYSigma);
521  const int n = 6 + nbInternalVariables;
522  // Iterate strains
523  for (int dd = 0; dd < 6; dd++) {
524  CHKERR VecZeroEntries(R);
525  double *array;
526  CHKERR VecGetArray(R, &array);
527  {
528  for (auto ii = 0; ii < 6; ii++) {
529  array[ii] = ep_row(ii, dd);
530  }
531  for (auto ii = 0; ii < nbInternalVariables; ii++) {
532  array[6 + ii] = alpha_row(ii, dd);
533  }
534  array[n] = y_row[dd];
535  }
536  CHKERR VecRestoreArray(R, &array);
537  CHKERR VecAssemblyBegin(R);
538  CHKERR VecAssemblyEnd(R);
539 
540  KSP ksp;
541  CHKERR SNESGetKSP(sNes, &ksp);
542  CHKERR KSPSolve(ksp, R, dChi);
543  CHKERR VecGetArray(dChi, &array);
544  {
545  for (auto ii = 0; ii < 6; ii++) {
546  Ep(ii, dd) = array[ii];
547  }
548  }
549  CHKERR VecRestoreArray(dChi, &array);
550  }
551  Cp.resize(6, 6, false);
552  noalias(Cp) = prod(C, Ep);
553  Cep.resize(6, 6, false);
554  noalias(Cep) = C + Cp;
556 }

◆ createMatAVecR()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::createMatAVecR ( )
Examples
ADOLCPlasticity.hpp.

Definition at line 314 of file ClosetPointProjection.cpp.

314  {
316  PetscInt n;
317  n = 1 + 6 + nbInternalVariables;
318  dataA.resize(n, n, false);
319  Mat A_tmp;
320  Vec R_tmp, Chi_tmp, dChi_tmp;
321  CHKERR MatCreateSeqDense(PETSC_COMM_SELF, n, n, &dataA(0, 0), &A_tmp);
322  CHKERR VecCreateSeq(PETSC_COMM_SELF, n, &R_tmp);
323  CHKERR VecCreateSeq(PETSC_COMM_SELF, n, &Chi_tmp);
324  CHKERR VecCreateSeq(PETSC_COMM_SELF, n, &dChi_tmp);
325 
326  A = SmartPetscObj<Mat>(A_tmp);
327  R = SmartPetscObj<Vec>(R_tmp);
328  Chi = SmartPetscObj<Vec>(Chi_tmp);
329  dChi = SmartPetscObj<Vec>(dChi_tmp);
330 
332 }

◆ evaluatePotentials()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::evaluatePotentials ( )
Examples
ADOLCPlasticity.hpp.

Definition at line 335 of file ClosetPointProjection.cpp.

335  {
337  CHKERR recordW();
338  CHKERR recordY();
339  CHKERR recordH();
341 }

◆ flowPotential()

virtual MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::flowPotential ( )
pure virtual

◆ freeHemholtzEnergy()

virtual MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::freeHemholtzEnergy ( )
pure virtual

◆ getActiveVariablesYH()

VectorAdaptor ADOLCPlasticity::ClosestPointProjection::getActiveVariablesYH ( )
inline
Examples
ADOLCPlasticity.hpp.

Definition at line 158 of file ADOLCPlasticity.hpp.

158  {
160  }

◆ getInternalFluxes()

VectorAdaptor ADOLCPlasticity::ClosestPointProjection::getInternalFluxes ( )
inline
Examples
ADOLCPlasticity.hpp.

Definition at line 164 of file ADOLCPlasticity.hpp.

164  {
166  }

◆ getInternalVariables()

VectorAdaptor ADOLCPlasticity::ClosestPointProjection::getInternalVariables ( )
inline
Examples
ADOLCPlasticity.hpp.

Definition at line 155 of file ADOLCPlasticity.hpp.

155  {
157  }

◆ getPlasticStrain()

VectorAdaptor ADOLCPlasticity::ClosestPointProjection::getPlasticStrain ( )
inline
Examples
ADOLCPlasticity.hpp.

Definition at line 149 of file ADOLCPlasticity.hpp.

149  {
150  return getVectorAdaptor(&(activeVariablesW[0]), 6);
151  }

◆ getStress()

VectorAdaptor ADOLCPlasticity::ClosestPointProjection::getStress ( )
inline
Examples
ADOLCPlasticity.hpp.

Definition at line 161 of file ADOLCPlasticity.hpp.

161  {
162  return getVectorAdaptor(&(gradientW[6]), 6);
163  }

◆ getTotalStrain()

VectorAdaptor ADOLCPlasticity::ClosestPointProjection::getTotalStrain ( )
inline
Examples
ADOLCPlasticity.hpp.

Definition at line 152 of file ADOLCPlasticity.hpp.

152  {
153  return getVectorAdaptor(&(activeVariablesW[6]), 6);
154  }

◆ playH()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::playH ( )
Examples
ADOLCPlasticity.hpp.

Definition at line 230 of file ClosetPointProjection.cpp.

230  {
232  auto active_variables_yh = getActiveVariablesYH();
233  int r;
234  int adloc_return_value = 0;
235  r = ::function(tapesTags[ClosestPointProjection::H], 1,
236  active_variables_yh.size(), &active_variables_yh[0], &h);
237  if (r < adloc_return_value) {
238  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
239  "ADOL-C function evaluation with error");
240  }
241  gradientH.resize(active_variables_yh.size());
242  r = gradient(tapesTags[ClosestPointProjection::H], active_variables_yh.size(),
243  &active_variables_yh[0], &gradientH[0]);
244  if (r < adloc_return_value) {
245  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
246  "ADOL-C function evaluation with error");
247  }
248  hessianH.resize(6 + nbInternalVariables, 6 + nbInternalVariables, false);
249  hessianH.clear();
250  vector<double *> hessian_h(active_variables_yh.size());
251  for (int dd = 0; dd < active_variables_yh.size(); dd++) {
252  hessian_h[dd] = &hessianH(dd, 0);
253  }
254  r = hessian(tapesTags[ClosestPointProjection::H], active_variables_yh.size(),
255  &active_variables_yh[0], &hessian_h[0]);
256  if (r < adloc_return_value) {
257  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
258  "ADOL-C function evaluation with error");
259  }
260  partialHSigma =
261  VectorAdaptor(6, ublas::shallow_array_adaptor<double>(6, &gradientH[0]));
264  ublas::shallow_array_adaptor<double>(nbInternalVariables, &gradientH[6]));
265  partialHFlux *= -1;
266  partial2HSigma.resize(6, 6, false);
267  for (int ii = 0; ii < 6; ii++) {
268  for (int jj = 0; jj <= ii; jj++) {
269  partial2HSigma(ii, jj) = hessianH(ii, jj);
270  }
271  }
273  for (unsigned int ii = 0; ii < nbInternalVariables; ii++) {
274  for (unsigned int jj = 0; jj <= ii; jj++) {
275  partial2HFlux(ii, jj) = -hessianH(6 + ii, 6 + jj);
276  }
277  }
278  partial2HSigmaFlux.resize(6, nbInternalVariables, false);
279  partial2HSigmaFlux.clear();
280  for (int ii = 0; ii < 6; ii++) {
281  for (unsigned int jj = 0; jj < nbInternalVariables; jj++) {
282  partial2HSigmaFlux(ii, jj) = -hessianH(6 + jj, ii);
283  }
284  }
286 }

◆ playH_NoHessian()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::playH_NoHessian ( )
Examples
ADOLCPlasticity.hpp.

Definition at line 287 of file ClosetPointProjection.cpp.

287  {
289  auto active_variables_yh = getActiveVariablesYH();
290  int r;
291  int adloc_return_value = 0;
292  r = ::function(tapesTags[ClosestPointProjection::H], 1,
293  active_variables_yh.size(), &active_variables_yh[0], &h);
294  if (r < adloc_return_value) {
295  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
296  "ADOL-C function evaluation with error");
297  }
298  gradientH.resize(active_variables_yh.size());
299  r = gradient(tapesTags[ClosestPointProjection::H], active_variables_yh.size(),
300  &active_variables_yh[0], &gradientH[0]);
301  if (r < adloc_return_value) {
302  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
303  "ADOL-C function evaluation with error");
304  }
305  partialHSigma =
306  VectorAdaptor(6, ublas::shallow_array_adaptor<double>(6, &gradientH[0]));
309  ublas::shallow_array_adaptor<double>(nbInternalVariables, &gradientH[6]));
310  partialHFlux *= -1;
312 }

◆ playPotentials()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::playPotentials ( )
Examples
ADOLCPlasticity.hpp.

Definition at line 343 of file ClosetPointProjection.cpp.

343  {
345  CHKERR playW();
346  CHKERR playY();
347  CHKERR playH();
349 }

◆ playPotentials_NoHessian()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::playPotentials_NoHessian ( )
Examples
ADOLCPlasticity.hpp.

Definition at line 351 of file ClosetPointProjection.cpp.

351  {
357 }

◆ playW()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::playW ( )
Examples
ADOLCPlasticity.hpp.

Definition at line 122 of file ClosetPointProjection.cpp.

122  {
124  int r;
125  int adloc_return_value = 0;
126  r = ::function(tapesTags[ClosestPointProjection::W], 1,
127  activeVariablesW.size(), &activeVariablesW[0], &w);
128  if (r < adloc_return_value) {
129  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
130  "ADOL-C function evaluation with error");
131  }
133  &activeVariablesW[0], &gradientW[0]);
134  if (r < adloc_return_value) {
135  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
136  "ADOL-C function evaluation with error");
137  }
138  hessianW.resize(activeVariablesW.size(), activeVariablesW.size(), false);
139  hessianW.clear();
140  vector<double *> hessian_w(activeVariablesW.size());
141  for (unsigned int dd = 0; dd < activeVariablesW.size(); dd++) {
142  hessian_w[dd] = &hessianW(dd, 0);
143  }
145  &activeVariablesW[0], &hessian_w[0]);
146  if (r < adloc_return_value) {
147  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
148  "ADOL-C function evaluation with error");
149  }
150  partialWStrainPlasticStrain.resize(6, 6, false);
151  for (int ii = 0; ii < 6; ii++) {
152  for (int jj = 0; jj < 6; jj++) {
153  partialWStrainPlasticStrain(ii, jj) = hessianW(6 + ii, jj);
154  }
155  }
156  C.resize(6, 6, false);
157  for (int ii = 0; ii < 6; ii++) {
158  for (int jj = 0; jj <= ii; jj++) {
159  C(ii, jj) = hessianW(6 + ii, 6 + jj);
160  }
161  }
162  D.resize(nbInternalVariables, nbInternalVariables, false);
163  for (unsigned int ii = 0; ii < nbInternalVariables; ii++) {
164  for (unsigned int jj = 0; jj <= ii; jj++) {
165  D(ii, jj) = hessianW(12 + ii, 12 + jj);
166  }
167  }
169 }

◆ playW_NoHessian()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::playW_NoHessian ( )
Examples
ADOLCPlasticity.hpp.

Definition at line 172 of file ClosetPointProjection.cpp.

172  {
174  int r;
175  int adloc_return_value = 0;
176  r = ::function(tapesTags[ClosestPointProjection::W], 1,
177  activeVariablesW.size(), &activeVariablesW[0], &w);
178  if (r < adloc_return_value) {
179  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
180  "ADOL-C function evaluation with error");
181  }
183  &activeVariablesW[0], &gradientW[0]);
184  if (r < adloc_return_value) {
185  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
186  "ADOL-C function evaluation with error");
187  }
189 }

◆ playY()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::playY ( )
Examples
ADOLCPlasticity.hpp.

Definition at line 191 of file ClosetPointProjection.cpp.

191  {
193  auto active_variables_yh = getActiveVariablesYH();
194  int r;
195  int adloc_return_value = 0;
196  r = ::function(tapesTags[ClosestPointProjection::Y], 1,
197  active_variables_yh.size(), &active_variables_yh[0], &y);
198  if (r < adloc_return_value) {
199  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
200  "ADOL-C function evaluation with error");
201  }
202  gradientY.resize(active_variables_yh.size());
203  r = gradient(tapesTags[ClosestPointProjection::Y], active_variables_yh.size(),
204  &active_variables_yh[0], &gradientY[0]);
205  if (r < adloc_return_value) {
206  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
207  "ADOL-C function evaluation with error");
208  }
209  partialYSigma =
210  VectorAdaptor(6, ublas::shallow_array_adaptor<double>(6, &gradientY[0]));
213  ublas::shallow_array_adaptor<double>(nbInternalVariables, &gradientY[6]));
215 }

◆ playY_NoGradient()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::playY_NoGradient ( )
Examples
ADOLCPlasticity.hpp.

Definition at line 217 of file ClosetPointProjection.cpp.

217  {
219  auto active_variables_yh = getActiveVariablesYH();
220  int r;
221  int adloc_return_value = 0;
222  r = ::function(tapesTags[ClosestPointProjection::Y], 1,
223  active_variables_yh.size(), &active_variables_yh[0], &y);
224  if (r < adloc_return_value) {
225  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
226  "ADOL-C function evaluation with error");
227  }
229 }

◆ recordH()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::recordH ( )

Record flow potential.

Examples
ADOLCPlasticity.hpp.

Definition at line 96 of file ClosetPointProjection.cpp.

96  {
98 
99  auto stress = getStress();
100  auto internal_fluxes = getInternalFluxes();
101 
103  {
104  // active variables
105  a_sTress.resize(6, false);
106  for (int ii = 0; ii < 6; ii++) {
107  a_sTress[ii] <<= stress[ii];
108  }
109  a_internalFluxes.resize(nbInternalVariables, false);
110  for (unsigned int ii = 0; ii < nbInternalVariables; ii++) {
111  a_internalFluxes[ii] <<= internal_fluxes[ii];
112  }
113  // evaluete functions
115  // return variables
116  a_h >>= h;
117  }
118  trace_off();
120 }

◆ recordTapes()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::recordTapes ( )

Record tapes.

Examples
ADOLCPlasticity.hpp.

Definition at line 610 of file ClosetPointProjection.cpp.

610  {
612  VectorDouble active(12 + nbInternalVariables);
613  VectorDouble gradient(12 + nbInternalVariables);
614  auto tmp_active = getVectorAdaptor(&(active[0]), 12 + nbInternalVariables);
615  auto tmp_gradient =
616  getVectorAdaptor(&(gradient[0]), 12 + nbInternalVariables);
617  activeVariablesW.swap(tmp_active);
618  gradientW.swap(tmp_gradient);
621 }

◆ recordW()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::recordW ( )

Record strain energy.

Examples
ADOLCPlasticity.hpp.

Definition at line 40 of file ClosetPointProjection.cpp.

40  {
42 
43  auto plastic_strain = getPlasticStrain();
44  auto total_strain = getTotalStrain();
45  auto internal_variables = getInternalVariables();
46 
48  {
49  // active variables
50  a_plasticStrain.resize(6, false);
51  for (int ii = 0; ii < 6; ii++) {
52  a_plasticStrain[ii] <<= plastic_strain[ii];
53  }
54  a_sTrain.resize(6, false);
55  for (int ii = 0; ii < 6; ii++) {
56  a_sTrain[ii] <<= total_strain[ii];
57  }
59  for (unsigned int ii = 0; ii < nbInternalVariables; ii++) {
60  a_internalVariables[ii] <<= internal_variables[ii];
61  }
62  // evaluete functions
64  a_w >>= w;
65  }
66  trace_off();
68 }

◆ recordY()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::recordY ( )

Record yield function.

Examples
ADOLCPlasticity.hpp.

Definition at line 70 of file ClosetPointProjection.cpp.

70  {
72 
73  auto stress = getStress();
74  auto internal_fluxes = getInternalFluxes();
75 
77  {
78  // active variables
79  a_sTress.resize(6, false);
80  for (int ii = 0; ii < 6; ii++) {
81  a_sTress[ii] <<= stress[ii];
82  }
84  for (unsigned int ii = 0; ii < nbInternalVariables; ii++) {
85  a_internalFluxes[ii] <<= internal_fluxes[ii];
86  }
87  // evaluete functions
89  // return variables
90  a_y >>= y;
91  }
92  trace_off();
94 }

◆ setParams()

virtual MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::setParams ( short  tag,
bool recalculate_elastic_tangent 
)
inlinevirtual

Set parameters for ADOL-C of constitutive relations.

Parameters
tag[in] - tag of the tape
recalculate_elastic_tangent[out] - if setParam set it to true, tangent matrix for elastic domain should be recalculated

Reimplemented in ADOLCPlasticity::ParaboloidalPlasticity, and ADOLCPlasticity::J2Plasticity< 3 >.

Examples
ADOLCPlasticity.hpp.

Definition at line 274 of file ADOLCPlasticity.hpp.

275  {
276  return 0;
277  }

◆ snesCreate()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::snesCreate ( )

Create nested snes.

Examples
ADOLCPlasticity.hpp.

Definition at line 433 of file ClosetPointProjection.cpp.

433  {
435  sNes = createSNES(PETSC_COMM_SELF);
436  CHKERR SNESSetFunction(sNes, R, ADOLCPlasticityRes, (void *)this);
437  CHKERR SNESSetJacobian(sNes, A, A, ADOLCPlasticityJac, (void *)this);
438 
439  KSP ksp;
440  PC pc;
441  SNESLineSearch line_search;
442  CHKERR SNESGetKSP(sNes, &ksp);
443  CHKERR KSPGetPC(ksp, &pc);
444  CHKERR KSPSetType(ksp, KSPPREONLY);
445  CHKERR PCSetType(pc, PCLU);
446  CHKERR SNESSetTolerances(sNes, 1e-8, 1e-6, 1e-12, 20, 1000);
447  CHKERR SNESGetLineSearch(sNes, &line_search);
448  // CHKERR SNESLineSearchSetType(line_search, SNESLINESEARCHBASIC);
449  CHKERR SNESLineSearchSetType(line_search, SNESLINESEARCHBT);
450  // CHKERR SNESLineSearchSetType(line_search, SNESLINESEARCHL2);
451 
452  CHKERR SNESAppendOptionsPrefix(sNes, "cp_");
453  CHKERR SNESSetFromOptions(sNes);
454 
456 }

◆ solveClosetProjection()

MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::solveClosetProjection ( )

Solve nonlinear system of equations to find stress, internal fluxes, and Lagrange plastic multiplier.

Examples
ADOLCPlasticity.hpp.

Definition at line 459 of file ClosetPointProjection.cpp.

459  {
461 
462  auto plastic_strain = getPlasticStrain();
463  auto internal_variables = getInternalVariables();
464 
465  {
466  double *array;
467  CHKERR VecGetArray(Chi, &array);
468  for (int ii = 0; ii < 6; ii++) {
469  array[ii] = plastic_strain[ii];
470  }
471  int nb_internal_variables = nbInternalVariables;
472  for (unsigned int ii = 0; ii < nb_internal_variables; ii++) {
473  array[6 + ii] = internal_variables[ii];
474  }
475  deltaGamma = 0;
476  CHKERR VecRestoreArray(Chi, &array);
477  }
478 
479  CHKERR SNESSolve(sNes, PETSC_NULL, Chi);
480  SNESConvergedReason reason;
481  CHKERR SNESGetConvergedReason(sNes, &reason);
482  if (reason < 0) {
483  int its;
484  CHKERR SNESGetIterationNumber(sNes, &its);
485  MOFEM_LOG_C(
486  "ADOLCPlasticitySelf", Sev::warning,
487  "Plasticity Closet Point Projection - number of Newton iterations = %d",
488  its);
489  plastic_strain = plasticStrain0;
490  internal_variables = internalVariables0;
491  deltaGamma = 0;
492  } else {
493  double *array;
494  CHKERR VecGetArray(Chi, &array);
495  for (int ii = 0; ii < 6; ii++) {
496  plastic_strain[ii] = array[ii];
497  }
498  int nb_internal_variables = nbInternalVariables;
499  for (unsigned int ii = 0; ii < nb_internal_variables; ii++) {
500  internal_variables[ii] = array[6 + ii];
501  }
502  deltaGamma = array[6 + nb_internal_variables];
503  CHKERR VecRestoreArray(Chi, &array);
504  }
506 }

◆ yieldFunction()

virtual MoFEMErrorCode ADOLCPlasticity::ClosestPointProjection::yieldFunction ( )
pure virtual

Friends And Related Function Documentation

◆ ADOLCPlasticityJac

MoFEMErrorCode ADOLCPlasticityJac ( SNES  snes,
Vec  chi,
Mat  A,
Mat  ,
void *  ctx 
)
friend

Function executed by nested SENES at evaluationg tangent matrix.

Examples
ADOLCPlasticity.hpp.

Definition at line 584 of file ClosetPointProjection.cpp.

585  {
588  cp = (ClosestPointProjection *)ctx;
589 
590  auto plastic_strain = cp->getPlasticStrain();
591  auto internal_variables = cp->getInternalVariables();
592 
593  {
594  const double *array;
595  CHKERR VecGetArrayRead(chi, &array);
596  for (auto ii = 0; ii < 6; ii++) {
597  plastic_strain[ii] = array[ii];
598  }
599  for (auto ii = 0; ii < cp->nbInternalVariables; ii++) {
600  internal_variables[ii] = array[6 + ii];
601  }
602  cp->deltaGamma = array[6 + cp->nbInternalVariables];
603  CHKERR VecRestoreArrayRead(chi, &array);
604  }
605  CHKERR cp->playPotentials();
606  CHKERR cp->calculateA();
608 }

◆ ADOLCPlasticityRes

MoFEMErrorCode ADOLCPlasticityRes ( SNES  snes,
Vec  chi,
Vec  r,
void *  ctx 
)
friend

Function executed by nested SENES at evaluationg residual.

Examples
ADOLCPlasticity.hpp.

Definition at line 558 of file ClosetPointProjection.cpp.

558  {
561  cp = (ClosestPointProjection *)ctx;
562 
563  auto plastic_strain = cp->getPlasticStrain();
564  auto internal_variables = cp->getInternalVariables();
565 
566  {
567  const double *array;
568  CHKERR VecGetArrayRead(chi, &array);
569  for (auto ii = 0; ii < 6; ii++) {
570  plastic_strain[ii] = array[ii];
571  }
572  for (auto ii = 0; ii < cp->nbInternalVariables; ii++) {
573  internal_variables[ii] = array[6 + ii];
574  }
575  cp->deltaGamma = array[6 + cp->nbInternalVariables];
576  CHKERR VecRestoreArrayRead(chi, &array);
577  }
578  CHKERR cp->playPotentials_NoHessian();
579  CHKERR cp->calculateR(r);
580 
582 }

Member Data Documentation

◆ A

SmartPetscObj<Mat> ADOLCPlasticity::ClosestPointProjection::A
Examples
ADOLCPlasticity.hpp.

Definition at line 317 of file ADOLCPlasticity.hpp.

◆ a_h

adouble ADOLCPlasticity::ClosestPointProjection::a_h

◆ a_internalFluxes

ublas::vector<adouble> ADOLCPlasticity::ClosestPointProjection::a_internalFluxes

◆ a_internalVariables

ublas::vector<adouble> ADOLCPlasticity::ClosestPointProjection::a_internalVariables

◆ a_plasticStrain

ublas::vector<adouble> ADOLCPlasticity::ClosestPointProjection::a_plasticStrain

◆ a_sTrain

ublas::vector<adouble> ADOLCPlasticity::ClosestPointProjection::a_sTrain

◆ a_sTress

ublas::vector<adouble> ADOLCPlasticity::ClosestPointProjection::a_sTress

◆ a_w

adouble ADOLCPlasticity::ClosestPointProjection::a_w

◆ a_y

adouble ADOLCPlasticity::ClosestPointProjection::a_y

◆ activeVariablesW

VectorAdaptor ADOLCPlasticity::ClosestPointProjection::activeVariablesW

◆ C

ublas::symmetric_matrix<double, ublas::lower> ADOLCPlasticity::ClosestPointProjection::C
Examples
ADOLCPlasticity.hpp.

Definition at line 182 of file ADOLCPlasticity.hpp.

◆ Cep

MatrixDouble ADOLCPlasticity::ClosestPointProjection::Cep
Examples
ADOLCPlasticity.hpp.

Definition at line 181 of file ADOLCPlasticity.hpp.

◆ Chi

SmartPetscObj<Vec> ADOLCPlasticity::ClosestPointProjection::Chi
Examples
ADOLCPlasticity.hpp.

Definition at line 319 of file ADOLCPlasticity.hpp.

◆ Cp

MatrixDouble ADOLCPlasticity::ClosestPointProjection::Cp
Examples
ADOLCPlasticity.hpp.

Definition at line 181 of file ADOLCPlasticity.hpp.

◆ D

ublas::symmetric_matrix<double, ublas::lower> ADOLCPlasticity::ClosestPointProjection::D
Examples
ADOLCPlasticity.hpp.

Definition at line 182 of file ADOLCPlasticity.hpp.

◆ dataA

ublas::matrix<double, ublas::column_major> ADOLCPlasticity::ClosestPointProjection::dataA
Examples
ADOLCPlasticity.hpp.

Definition at line 321 of file ADOLCPlasticity.hpp.

◆ dChi

SmartPetscObj<Vec> ADOLCPlasticity::ClosestPointProjection::dChi
Examples
ADOLCPlasticity.hpp.

Definition at line 320 of file ADOLCPlasticity.hpp.

◆ deltaGamma

double ADOLCPlasticity::ClosestPointProjection::deltaGamma
Examples
ADOLCPlasticity.hpp.

Definition at line 180 of file ADOLCPlasticity.hpp.

◆ Ep

MatrixDouble ADOLCPlasticity::ClosestPointProjection::Ep
Examples
ADOLCPlasticity.hpp.

Definition at line 181 of file ADOLCPlasticity.hpp.

◆ gradientH

VectorDouble ADOLCPlasticity::ClosestPointProjection::gradientH
Examples
ADOLCPlasticity.hpp.

Definition at line 207 of file ADOLCPlasticity.hpp.

◆ gradientW

VectorAdaptor ADOLCPlasticity::ClosestPointProjection::gradientW
Examples
ADOLCPlasticity.hpp.

Definition at line 174 of file ADOLCPlasticity.hpp.

◆ gradientY

VectorDouble ADOLCPlasticity::ClosestPointProjection::gradientY
Examples
ADOLCPlasticity.hpp.

Definition at line 203 of file ADOLCPlasticity.hpp.

◆ h

double ADOLCPlasticity::ClosestPointProjection::h
Examples
ADOLCPlasticity.hpp.

Definition at line 178 of file ADOLCPlasticity.hpp.

◆ hessianH

MatrixDouble ADOLCPlasticity::ClosestPointProjection::hessianH
Examples
ADOLCPlasticity.hpp.

Definition at line 208 of file ADOLCPlasticity.hpp.

◆ hessianW

MatrixDouble ADOLCPlasticity::ClosestPointProjection::hessianW
Examples
ADOLCPlasticity.hpp.

Definition at line 199 of file ADOLCPlasticity.hpp.

◆ integrationRule

boost::function<int(int, int, int)> ADOLCPlasticity::ClosestPointProjection::integrationRule
Initial value:
= [](int, int, int p) {
return 2 * (p - 1);
}
Examples
ADOLCPlasticity.hpp.

Definition at line 142 of file ADOLCPlasticity.hpp.

◆ internalVariables0

VectorDouble ADOLCPlasticity::ClosestPointProjection::internalVariables0
Examples
ADOLCPlasticity.hpp.

Definition at line 146 of file ADOLCPlasticity.hpp.

◆ nbInternalVariables

int ADOLCPlasticity::ClosestPointProjection::nbInternalVariables

◆ partial2HFlux

ublas::symmetric_matrix<double, ublas::lower> ADOLCPlasticity::ClosestPointProjection::partial2HFlux

Second partial derivative of flow potential with respect to internal fluxes

Examples
ADOLCPlasticity.hpp.

Definition at line 312 of file ADOLCPlasticity.hpp.

◆ partial2HSigma

ublas::symmetric_matrix<double, ublas::lower> ADOLCPlasticity::ClosestPointProjection::partial2HSigma

Second partial derivative of flow potential with respect to stresses and internal

Examples
ADOLCPlasticity.hpp.

Definition at line 309 of file ADOLCPlasticity.hpp.

◆ partial2HSigmaFlux

MatrixDouble ADOLCPlasticity::ClosestPointProjection::partial2HSigmaFlux

Mixed decond partial derivative of flow potential with respect to stresses and internal fluxes

Examples
ADOLCPlasticity.hpp.

Definition at line 315 of file ADOLCPlasticity.hpp.

◆ partialHFlux

VectorAdaptor ADOLCPlasticity::ClosestPointProjection::partialHFlux
Examples
ADOLCPlasticity.hpp.

Definition at line 304 of file ADOLCPlasticity.hpp.

◆ partialHSigma

VectorAdaptor ADOLCPlasticity::ClosestPointProjection::partialHSigma
Examples
ADOLCPlasticity.hpp.

Definition at line 302 of file ADOLCPlasticity.hpp.

◆ partialWStrainPlasticStrain

MatrixDouble ADOLCPlasticity::ClosestPointProjection::partialWStrainPlasticStrain
Examples
ADOLCPlasticity.hpp.

Definition at line 296 of file ADOLCPlasticity.hpp.

◆ partialYFlux

VectorAdaptor ADOLCPlasticity::ClosestPointProjection::partialYFlux
Examples
ADOLCPlasticity.hpp.

Definition at line 300 of file ADOLCPlasticity.hpp.

◆ partialYSigma

VectorAdaptor ADOLCPlasticity::ClosestPointProjection::partialYSigma
Examples
ADOLCPlasticity.hpp.

Definition at line 298 of file ADOLCPlasticity.hpp.

◆ plasticStrain0

VectorDouble ADOLCPlasticity::ClosestPointProjection::plasticStrain0
Examples
ADOLCPlasticity.hpp.

Definition at line 147 of file ADOLCPlasticity.hpp.

◆ R

SmartPetscObj<Vec> ADOLCPlasticity::ClosestPointProjection::R
Examples
ADOLCPlasticity.hpp.

Definition at line 318 of file ADOLCPlasticity.hpp.

◆ sNes

SmartPetscObj<SNES> ADOLCPlasticity::ClosestPointProjection::sNes
Examples
ADOLCPlasticity.hpp.

Definition at line 322 of file ADOLCPlasticity.hpp.

◆ tapesTags

std::array<int, LAST_TAPE> ADOLCPlasticity::ClosestPointProjection::tapesTags
Examples
ADOLCPlasticity.hpp.

Definition at line 169 of file ADOLCPlasticity.hpp.

◆ w

double ADOLCPlasticity::ClosestPointProjection::w
Examples
ADOLCPlasticity.hpp.

Definition at line 176 of file ADOLCPlasticity.hpp.

◆ y

double ADOLCPlasticity::ClosestPointProjection::y
Examples
ADOLCPlasticity.hpp.

Definition at line 177 of file ADOLCPlasticity.hpp.


The documentation for this struct was generated from the following files:
ADOLCPlasticity::ClosestPointProjection::getActiveVariablesYH
VectorAdaptor getActiveVariablesYH()
Definition: ADOLCPlasticity.hpp:158
ADOLCPlasticity::ClosestPointProjection::playW
MoFEMErrorCode playW()
Definition: ClosetPointProjection.cpp:122
ADOLCPlasticity::ClosestPointProjection::playW_NoHessian
MoFEMErrorCode playW_NoHessian()
Definition: ClosetPointProjection.cpp:172
ADOLCPlasticity::ClosestPointProjection::getInternalVariables
VectorAdaptor getInternalVariables()
Definition: ADOLCPlasticity.hpp:155
ADOLCPlasticity::ClosestPointProjection::flowPotential
virtual MoFEMErrorCode flowPotential()=0
Set flow potential.
ADOLCPlasticity::ClosestPointProjection::partial2HSigmaFlux
MatrixDouble partial2HSigmaFlux
Definition: ADOLCPlasticity.hpp:315
ADOLCPlasticity::ClosestPointProjection::plasticStrain0
VectorDouble plasticStrain0
Definition: ADOLCPlasticity.hpp:147
ADOLCPlasticity::ClosestPointProjection::a_h
adouble a_h
Definition: ADOLCPlasticity.hpp:330
ADOLCPlasticity::ClosestPointProjection::playH_NoHessian
MoFEMErrorCode playH_NoHessian()
Definition: ClosetPointProjection.cpp:287
ADOLCPlasticity::ClosestPointProjection::Ep
MatrixDouble Ep
Definition: ADOLCPlasticity.hpp:181
ADOLCPlasticity::ClosestPointProjection::deltaGamma
double deltaGamma
Definition: ADOLCPlasticity.hpp:180
ADOLCPlasticity::ClosestPointProjection::ClosestPointProjection
ClosestPointProjection()
Definition: ClosetPointProjection.cpp:15
MoFEM::createSNES
auto createSNES(MPI_Comm comm)
Definition: PetscSmartObj.hpp:255
ADOLCPlasticity::ClosestPointProjection::a_plasticStrain
ublas::vector< adouble > a_plasticStrain
Definition: ADOLCPlasticity.hpp:325
ADOLCPlasticity::ClosestPointProjection::evaluatePotentials
MoFEMErrorCode evaluatePotentials()
Definition: ClosetPointProjection.cpp:335
ADOLCPlasticity::ClosestPointProjection::A
SmartPetscObj< Mat > A
Definition: ADOLCPlasticity.hpp:317
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
ADOLCPlasticity::ClosestPointProjection::a_sTrain
ublas::vector< adouble > a_sTrain
Definition: ADOLCPlasticity.hpp:324
ADOLCPlasticity::ClosestPointProjection::partialYFlux
VectorAdaptor partialYFlux
Definition: ADOLCPlasticity.hpp:300
ADOLCPlasticity::ClosestPointProjection::partialWStrainPlasticStrain
MatrixDouble partialWStrainPlasticStrain
Definition: ADOLCPlasticity.hpp:296
ADOLCPlasticity::ClosestPointProjection::Cep
MatrixDouble Cep
Definition: ADOLCPlasticity.hpp:181
sdf.r
int r
Definition: sdf.py:8
ADOLCPlasticity::ClosestPointProjection::dChi
SmartPetscObj< Vec > dChi
Definition: ADOLCPlasticity.hpp:320
ADOLCPlasticity::ClosestPointProjection::hessianW
MatrixDouble hessianW
Definition: ADOLCPlasticity.hpp:199
ADOLCPlasticity::ClosestPointProjection::playY
MoFEMErrorCode playY()
Definition: ClosetPointProjection.cpp:191
ADOLCPlasticity::ClosestPointProjection::recordY
MoFEMErrorCode recordY()
Record yield function.
Definition: ClosetPointProjection.cpp:70
ADOLCPlasticity::ClosestPointProjection::gradientW
VectorAdaptor gradientW
Definition: ADOLCPlasticity.hpp:174
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
ADOLCPlasticity::ClosestPointProjection::H
@ H
Definition: ADOLCPlasticity.hpp:168
ADOLCPlasticity::ClosestPointProjection::getPlasticStrain
VectorAdaptor getPlasticStrain()
Definition: ADOLCPlasticity.hpp:149
ADOLCPlasticity::ClosestPointProjection::playPotentials
MoFEMErrorCode playPotentials()
Definition: ClosetPointProjection.cpp:343
MoFEM::getVectorAdaptor
auto getVectorAdaptor(T1 ptr, const size_t n)
Get Vector adaptor.
Definition: Templates.hpp:31
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
ADOLCPlasticity::ClosestPointProjection::getStress
VectorAdaptor getStress()
Definition: ADOLCPlasticity.hpp:161
ADOLCPlasticity::ClosestPointProjection::getInternalFluxes
VectorAdaptor getInternalFluxes()
Definition: ADOLCPlasticity.hpp:164
ADOLCPlasticity::ClosestPointProjection::h
double h
Definition: ADOLCPlasticity.hpp:178
ADOLCPlasticity::ClosestPointProjection::partial2HFlux
ublas::symmetric_matrix< double, ublas::lower > partial2HFlux
Definition: ADOLCPlasticity.hpp:312
ADOLCPlasticity::ClosestPointProjection::partialYSigma
VectorAdaptor partialYSigma
Definition: ADOLCPlasticity.hpp:298
ADOLCPlasticity::ClosestPointProjection::Chi
SmartPetscObj< Vec > Chi
Definition: ADOLCPlasticity.hpp:319
MOFEM_OPERATION_UNSUCCESSFUL
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition: definitions.h:34
ADOLCPlasticity::ClosestPointProjection::R
SmartPetscObj< Vec > R
Definition: ADOLCPlasticity.hpp:318
ADOLCPlasticity::ClosestPointProjection::a_w
adouble a_w
Definition: ADOLCPlasticity.hpp:328
ADOLCPlasticity::ClosestPointProjection::dataA
ublas::matrix< double, ublas::column_major > dataA
Definition: ADOLCPlasticity.hpp:321
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
ADOLCPlasticity::ClosestPointProjection::partialHSigma
VectorAdaptor partialHSigma
Definition: ADOLCPlasticity.hpp:302
MoFEM::Types::VectorAdaptor
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
ADOLCPlasticity::ClosestPointProjection::a_internalFluxes
ublas::vector< adouble > a_internalFluxes
Definition: ADOLCPlasticity.hpp:327
ADOLCPlasticity::ClosestPointProjection::nbInternalVariables
int nbInternalVariables
Definition: ADOLCPlasticity.hpp:141
ADOLCPlasticity::ClosestPointProjection::activeVariablesW
VectorAdaptor activeVariablesW
Definition: ADOLCPlasticity.hpp:173
ADOLCPlasticity::ClosestPointProjection::tapesTags
std::array< int, LAST_TAPE > tapesTags
Definition: ADOLCPlasticity.hpp:169
ADOLCPlasticity::ClosestPointProjection::partialHFlux
VectorAdaptor partialHFlux
Definition: ADOLCPlasticity.hpp:304
convert.n
n
Definition: convert.py:82
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
ADOLCPlasticity::ClosestPointProjection::ADOLCPlasticityRes
friend MoFEMErrorCode ADOLCPlasticityRes(SNES, Vec, Vec, void *ctx)
Function executed by nested SENES at evaluationg residual.
Definition: ClosetPointProjection.cpp:558
ADOLCPlasticity::ClosestPointProjection::w
double w
Definition: ADOLCPlasticity.hpp:176
ADOLCPlasticity::ClosestPointProjection::sNes
SmartPetscObj< SNES > sNes
Definition: ADOLCPlasticity.hpp:322
ADOLCPlasticity::ClosestPointProjection::getTotalStrain
VectorAdaptor getTotalStrain()
Definition: ADOLCPlasticity.hpp:152
ADOLCPlasticity::ClosestPointProjection::a_internalVariables
ublas::vector< adouble > a_internalVariables
Definition: ADOLCPlasticity.hpp:326
ADOLCPlasticity::ClosestPointProjection::recordW
MoFEMErrorCode recordW()
Record strain energy.
Definition: ClosetPointProjection.cpp:40
ADOLCPlasticity::ClosestPointProjection::playH
MoFEMErrorCode playH()
Definition: ClosetPointProjection.cpp:230
ADOLCPlasticity::ClosestPointProjection::a_y
adouble a_y
Definition: ADOLCPlasticity.hpp:329
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
ADOLCPlasticity::ClosestPointProjection::calculateA
MoFEMErrorCode calculateA()
Definition: ClosetPointProjection.cpp:382
ADOLCPlasticity::ClosestPointProjection::LAST_TAPE
@ LAST_TAPE
Definition: ADOLCPlasticity.hpp:168
ADOLCPlasticity::ClosestPointProjection::a_sTress
ublas::vector< adouble > a_sTress
Definition: ADOLCPlasticity.hpp:327
ADOLCPlasticity::ClosestPointProjection::gradientY
VectorDouble gradientY
Definition: ADOLCPlasticity.hpp:203
ADOLCPlasticity::ClosestPointProjection::D
ublas::symmetric_matrix< double, ublas::lower > D
Definition: ADOLCPlasticity.hpp:182
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
ADOLCPlasticity::ClosestPointProjection::gradientH
VectorDouble gradientH
Definition: ADOLCPlasticity.hpp:207
ADOLCPlasticity::ClosestPointProjection::C
ublas::symmetric_matrix< double, ublas::lower > C
Definition: ADOLCPlasticity.hpp:182
ADOLCPlasticity::ClosestPointProjection::ADOLCPlasticityJac
friend MoFEMErrorCode ADOLCPlasticityJac(SNES, Vec, Mat, Mat, void *ctx)
Function executed by nested SENES at evaluationg tangent matrix.
Definition: ClosetPointProjection.cpp:584
ADOLCPlasticity::ClosestPointProjection::partial2HSigma
ublas::symmetric_matrix< double, ublas::lower > partial2HSigma
Definition: ADOLCPlasticity.hpp:309
MoFEM::SmartPetscObj< Mat >
ADOLCPlasticity::ClosestPointProjection::freeHemholtzEnergy
virtual MoFEMErrorCode freeHemholtzEnergy()=0
Set Hemholtz energy.
ADOLCPlasticity::ClosestPointProjection::playY_NoGradient
MoFEMErrorCode playY_NoGradient()
Definition: ClosetPointProjection.cpp:217
convert.int
int
Definition: convert.py:64
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
ADOLCPlasticity::ClosestPointProjection::internalVariables0
VectorDouble internalVariables0
Definition: ADOLCPlasticity.hpp:146
ADOLCPlasticity::ClosestPointProjection::Cp
MatrixDouble Cp
Definition: ADOLCPlasticity.hpp:181
ADOLCPlasticity::ClosestPointProjection::Y
@ Y
Definition: ADOLCPlasticity.hpp:168
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
ADOLCPlasticity::ClosestPointProjection::yieldFunction
virtual MoFEMErrorCode yieldFunction()=0
Set yield function.
ADOLCPlasticity::ClosestPointProjection::hessianH
MatrixDouble hessianH
Definition: ADOLCPlasticity.hpp:208
ADOLCPlasticity::ClosestPointProjection::y
double y
Definition: ADOLCPlasticity.hpp:177
ADOLCPlasticity::ClosestPointProjection::recordH
MoFEMErrorCode recordH()
Record flow potential.
Definition: ClosetPointProjection.cpp:96
ADOLCPlasticity::ClosestPointProjection::W
@ W
Definition: ADOLCPlasticity.hpp:168