v0.14.0
Public Member Functions | Public Attributes | Friends | List of all members
SmallStrainPlasticity::ClosestPointProjection Struct Reference

Closest point projection algorithm. More...

#include <users_modules/small_strain_plasticity/src/SmallStrainPlasticity.hpp>

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

Public Member Functions

 ClosestPointProjection ()
 
virtual PetscErrorCode setActiveVariablesW ()
 
virtual PetscErrorCode setActiveVariablesYH ()
 
virtual PetscErrorCode rEcordW ()
 
virtual PetscErrorCode rEcordY ()
 
virtual PetscErrorCode rEcordH ()
 
virtual PetscErrorCode pLayW ()
 
virtual PetscErrorCode pLayW_NoHessian ()
 
virtual PetscErrorCode pLayY ()
 
virtual PetscErrorCode pLayY_NoGradient ()
 
virtual PetscErrorCode pLayH ()
 
virtual PetscErrorCode createMatAVecR ()
 
virtual PetscErrorCode destroyMatAVecR ()
 
virtual PetscErrorCode evaluatePotentials ()
 
virtual PetscErrorCode cAlculateR (Vec R)
 
virtual PetscErrorCode cAlculateA ()
 
PetscErrorCode snesCreate ()
 
PetscErrorCode snesDestroy ()
 
PetscErrorCode solveColasetProjection ()
 
PetscErrorCode consistentTangent ()
 
virtual PetscErrorCode freeHemholtzEnergy ()
 
virtual PetscErrorCode yieldFunction ()
 
virtual PetscErrorCode flowPotential ()
 

Public Attributes

VectorDouble sTrain
 
VectorDouble internalVariables
 
VectorDouble internalVariables0
 
VectorDouble plasticStrain0
 
VectorDouble plasticStrain
 
double deltaGamma
 
double tOl
 
int gG
 
VectorAdaptor sTress
 
VectorAdaptor internalFluxes
 
ublas::symmetric_matrix< double, ublas::lower > C
 
ublas::symmetric_matrix< double, ublas::lower > D
 
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
 
vector< int > tAgs
 
VectorDouble activeVariablesW
 
VectorDouble activeVariablesYH
 
double w
 
double y
 
double h
 
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
 
VectorDouble gradientW
 
MatrixDouble hessianW
 
VectorDouble gradientY
 
VectorDouble gradientH
 
MatrixDouble hessianH
 
Mat A
 
ublas::matrix< double, ublas::column_major > dataA
 
Vec R
 
Vec Chi
 
Vec dChi
 
SNES sNes
 
KSP kSp
 
PC pC
 
SNESLineSearch lineSearch
 
MatrixDouble Ep
 
MatrixDouble Lp
 
MatrixDouble Cp
 
MatrixDouble Cep
 
VectorDouble Yp
 

Friends

PetscErrorCode SmallStrainPlasticityfRes (SNES, Vec, Vec, void *ctx)
 
PetscErrorCode SmallStrainPlasticityfJac (SNES, Vec, Mat, Mat, void *ctx)
 

Detailed Description

Closest point projection algorithm.

Examples
SmallStrainPlasticityMaterialModels.hpp.

Definition at line 482 of file SmallStrainPlasticity.hpp.

Constructor & Destructor Documentation

◆ ClosestPointProjection()

SmallStrainPlasticity::ClosestPointProjection::ClosestPointProjection ( )
inline

Member Function Documentation

◆ cAlculateA()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::cAlculateA ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 978 of file SmallStrainPlasticity.cpp.

978  {
979  PetscFunctionBegin;
980  // matrix A
981  ierr = MatZeroEntries(A); CHKERRQ(ierr);
982  // row 0
984  MatrixDouble a01 = prod(partial2HSigmaFlux,D);
985  for(int ii = 0;ii<6;ii++) {
986  for(int jj = 0;jj<6;jj++) {
987  dataA(ii,jj) = deltaGamma*a00(ii,jj);
988  }
989  for(unsigned int jj = 0;jj<internalVariables.size();jj++) {
990  dataA(ii,6+jj) = deltaGamma*a01(ii,jj);
991  }
992  dataA(ii,ii) -= 1;
993  }
994  for(int ii = 0;ii<6;ii++) {
995  dataA(ii,6+internalVariables.size()) = partialHSigma[ii];
996  }
997  // row 1
998  MatrixDouble a11 = prod(partial2HFlux,D);
1000  for(unsigned int ii = 0;ii<internalVariables.size();ii++) {
1001  for(int jj = 0;jj<6;jj++) {
1002  dataA(6+ii,jj) += deltaGamma*a10(ii,jj);
1003  }
1004  for(unsigned int jj = 0;jj<internalVariables.size();jj++) {
1005  dataA(6+ii,6+jj) += deltaGamma*a11(ii,jj);
1006  }
1007  dataA(6+ii,6+ii) -= 1;
1008  }
1009  for(unsigned int ii = 0;ii<internalVariables.size();ii++) {
1010  dataA(6+ii,6+internalVariables.size()) = partialHFlux[ii];
1011  }
1012  // row 3
1013  VectorDouble partial_y_sigma_c = prod(trans(partialWStrainPlasticStrain),partialYSigma);
1014  VectorDouble partial_y_flux_d = prod(trans(D),partialYFlux);
1015  for(unsigned int dd = 0;dd<partial_y_sigma_c.size();dd++) {
1016  dataA(6+internalVariables.size(),dd) = partial_y_sigma_c[dd];
1017  }
1018  for(unsigned int dd = 0;dd<partial_y_flux_d.size();dd++) {
1019  dataA(6+internalVariables.size(),6+dd) = partial_y_flux_d[dd];
1020  }
1021  dataA(6+internalVariables.size(),6+internalVariables.size()) = 0;
1022  ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1023  ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1024  PetscFunctionReturn(0);
1025 }

◆ cAlculateR()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::cAlculateR ( Vec  R)
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 962 of file SmallStrainPlasticity.cpp.

962  {
963  PetscFunctionBegin;
964  // Residual
965  double *array;
966  ierr = VecGetArray(R,&array); CHKERRQ(ierr);
967  for(int ii = 0;ii<6;ii++) {
968  array[ii] = -plasticStrain[ii] + plasticStrain0[ii] + deltaGamma*partialHSigma[ii];
969  }
970  for(unsigned int ii = 0;ii<internalVariables.size();ii++) {
971  array[6+ii] = -internalVariables[ii] + internalVariables0[ii] + deltaGamma*partialHFlux[ii];
972  }
973  array[6+internalVariables.size()] = y;
974  ierr = VecRestoreArray(R,&array); CHKERRQ(ierr);
975  PetscFunctionReturn(0);
976 }

◆ consistentTangent()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::consistentTangent ( )
Examples
SmallStrainPlasticity.hpp.

Definition at line 1095 of file SmallStrainPlasticity.cpp.

1095  {
1096  PetscFunctionBegin;
1097  ierr = evaluatePotentials(); CHKERRQ(ierr);
1098  ierr = cAlculateA(); CHKERRQ(ierr);
1099  Ep.resize(6,6,false);
1100  Lp.resize(internalVariables.size(),6,false);
1101  Yp.resize(6,false);
1102  MatrixDouble ep_row;
1103  ep_row = deltaGamma*prod(partial2HSigma,C);
1104  MatrixDouble alpha_row;
1105  alpha_row = deltaGamma*prod(trans(partial2HSigmaFlux),C);
1106  VectorDouble y_row;
1107  y_row = prod(trans(C),partialYSigma);
1108  const int n = 6+internalVariables.size();
1109  for(int dd = 0;dd<6;dd++) {
1110  ierr = VecZeroEntries(R); CHKERRQ(ierr);
1111  double *array;
1112  ierr = VecGetArray(R,&array); CHKERRQ(ierr);
1113  {
1114  for(int ii = 0;ii<6;ii++) {
1115  array[ii] = ep_row(ii,dd);
1116  }
1117  for(unsigned int ii = 0;ii<internalVariables.size();ii++) {
1118  array[6+ii] = alpha_row(ii,dd);
1119  }
1120  array[n] = y_row[dd];
1121  }
1122  ierr = VecRestoreArray(R,&array); CHKERRQ(ierr);
1123  ierr = VecAssemblyBegin(R); CHKERRQ(ierr);
1124  ierr = VecAssemblyEnd(R); CHKERRQ(ierr);
1125  ierr = KSPSolve(kSp,R,dChi); CHKERRQ(ierr);
1126  ierr = VecGetArray(dChi,&array); CHKERRQ(ierr);
1127  {
1128  for(int ii = 0;ii<6;ii++) {
1129  Ep(ii,dd) = array[ii];
1130  }
1131  for(unsigned int ii = 0;ii<internalVariables.size();ii++) {
1132  Lp(ii,dd) = array[6+ii];
1133  }
1134  Yp[dd] = array[n];
1135  }
1136  ierr = VecRestoreArray(dChi,&array); CHKERRQ(ierr);
1137  }
1138  Cp.resize(6,6,false);
1139  noalias(Cp) = prod(C,Ep);
1140  Cep.resize(6,6,false);
1141  noalias(Cep) = C+Cp;
1142  // cout << endl;
1143  // cout << "A " << dataA << endl;
1144  // cout << "C " << C << endl;
1145  // cout << "Ep " << Ep << endl;
1146  // // cout << "Lp " << Lp << endl;
1147  // // cout << "Yp " << Yp << endl;
1148  // cout << "Cp " << Cp << endl;
1149  // cout << "Cep " << Cep << endl;
1150  PetscFunctionReturn(0);
1151 }

◆ createMatAVecR()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::createMatAVecR ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 924 of file SmallStrainPlasticity.cpp.

924  {
925  PetscFunctionBegin;
926  PetscInt n;
927  n = 1+6+internalVariables.size();
928  dataA.resize(n,n,false);
929  ierr = MatCreateSeqDense(PETSC_COMM_SELF,n,n,&dataA(0,0),&A); CHKERRQ(ierr);
930  ierr = VecCreateSeq(PETSC_COMM_SELF,n,&R); CHKERRQ(ierr);
931  ierr = VecCreateSeq(PETSC_COMM_SELF,n,&Chi); CHKERRQ(ierr);
932  ierr = VecCreateSeq(PETSC_COMM_SELF,n,&dChi); CHKERRQ(ierr);
933  PetscFunctionReturn(0);
934 }

◆ destroyMatAVecR()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::destroyMatAVecR ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 936 of file SmallStrainPlasticity.cpp.

936  {
937  PetscFunctionBegin;
938  ierr = MatDestroy(&A); CHKERRQ(ierr);
939  ierr = VecDestroy(&R); CHKERRQ(ierr);
940  ierr = VecDestroy(&Chi); CHKERRQ(ierr);
941  ierr = VecDestroy(&dChi); CHKERRQ(ierr);
942  PetscFunctionReturn(0);
943 }

◆ evaluatePotentials()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::evaluatePotentials ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 945 of file SmallStrainPlasticity.cpp.

945  {
946  PetscFunctionBegin;
947  if(gG == 0) {
948  ierr = rEcordW(); CHKERRQ(ierr);
949  }
950  ierr = setActiveVariablesW(); CHKERRQ(ierr);
951  ierr = pLayW(); CHKERRQ(ierr);
952  if(gG == 0) {
953  ierr = rEcordY(); CHKERRQ(ierr);
954  ierr = rEcordH(); CHKERRQ(ierr);
955  }
956  ierr = setActiveVariablesYH(); CHKERRQ(ierr);
957  ierr = pLayY(); CHKERRQ(ierr);
958  ierr = pLayH(); CHKERRQ(ierr);
959  PetscFunctionReturn(0);
960 }

◆ flowPotential()

virtual PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::flowPotential ( )
inlinevirtual

Reimplemented in SmallStrainParaboloidalPlasticity, and SmallStrainJ2Plasticity.

Examples
SmallStrainPlasticity.hpp.

Definition at line 591 of file SmallStrainPlasticity.hpp.

591  {
592  PetscFunctionBegin;
593  SETERRQ(
594  PETSC_COMM_SELF,
596  "Not implemented (overload class and implement this function)"
597  );
598  PetscFunctionReturn(0);
599  }

◆ freeHemholtzEnergy()

virtual PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::freeHemholtzEnergy ( )
inlinevirtual

Reimplemented in SmallStrainParaboloidalPlasticity, and SmallStrainJ2Plasticity.

Examples
SmallStrainPlasticity.hpp.

Definition at line 571 of file SmallStrainPlasticity.hpp.

571  {
572  PetscFunctionBegin;
573  SETERRQ(
574  PETSC_COMM_SELF,
576  "Not implemented (overload class and implement this function)"
577  );
578  PetscFunctionReturn(0);
579  }

◆ pLayH()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::pLayH ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 870 of file SmallStrainPlasticity.cpp.

870  {
871  PetscFunctionBegin;
872  int r;
873  int adloc_return_value = 0;
874  r = ::function(tAgs[2],1,activeVariablesYH.size(),&activeVariablesYH[0],&h);
875  if(r<adloc_return_value) {
876  SETERRQ(PETSC_COMM_SELF,MOFEM_OPERATION_UNSUCCESSFUL,"ADOL-C function evaluation with error");
877  }
878  gradientH.resize(activeVariablesYH.size());
879  r = gradient(tAgs[2],activeVariablesYH.size(),&activeVariablesYH[0],&gradientH[0]);
880  if(r<adloc_return_value) {
881  SETERRQ(PETSC_COMM_SELF,MOFEM_OPERATION_UNSUCCESSFUL,"ADOL-C function evaluation with error");
882  }
883  hessianH.resize(6+internalVariables.size(),6+internalVariables.size(),false);
884  hessianH.clear();
885  vector<double*> hessian_h(activeVariablesYH.size());
886  for(int dd = 0;dd<activeVariablesYH.size();dd++) {
887  hessian_h[dd] = &hessianH(dd,0);
888  }
889  r = hessian(tAgs[2],activeVariablesYH.size(),&activeVariablesYH[0],&hessian_h[0]);
890  if(r<adloc_return_value) {
891  SETERRQ(PETSC_COMM_SELF,MOFEM_OPERATION_UNSUCCESSFUL,"ADOL-C function evaluation with error");
892  }
893  partialHSigma = VectorAdaptor(6,ublas::shallow_array_adaptor<double>(6,&gradientH[0]));
895  internalVariables.size(),
896  ublas::shallow_array_adaptor<double>(
897  internalVariables.size(),
898  &gradientH[6]
899  )
900  );
901  partialHFlux *= -1;
902  partial2HSigma.resize(6,6,false);
903  for(int ii = 0;ii<6;ii++) {
904  for(int jj = 0;jj<=ii;jj++) {
905  partial2HSigma(ii,jj) = hessianH(ii,jj);
906  }
907  }
908  partial2HFlux.resize(internalVariables.size(),internalVariables.size(),false);
909  for(unsigned int ii = 0;ii<internalVariables.size();ii++) {
910  for(unsigned int jj = 0;jj<=ii;jj++) {
911  partial2HFlux(ii,jj) = -hessianH(6+ii,6+jj);
912  }
913  }
914  partial2HSigmaFlux.resize(6,internalVariables.size(),false);
915  partial2HSigmaFlux.clear();
916  for(int ii = 0;ii<6;ii++) {
917  for(unsigned int jj = 0;jj<internalVariables.size();jj++) {
918  partial2HSigmaFlux(ii,jj) = -hessianH(6+jj,ii);
919  }
920  }
921  PetscFunctionReturn(0);
922 }

◆ pLayW()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::pLayW ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 762 of file SmallStrainPlasticity.cpp.

762  {
763  PetscFunctionBegin;
764  int r;
765  int adloc_return_value = 0;
766  r = ::function(tAgs[0],1,activeVariablesW.size(),&activeVariablesW[0],&w);
767  if(r<adloc_return_value) {
768  SETERRQ(PETSC_COMM_SELF,MOFEM_OPERATION_UNSUCCESSFUL,"ADOL-C function evaluation with error");
769  }
770  gradientW.resize(activeVariablesW.size(),false);
771  r = gradient(tAgs[0],activeVariablesW.size(),&activeVariablesW[0],&gradientW[0]);
772  if(r<adloc_return_value) {
773  SETERRQ(PETSC_COMM_SELF,MOFEM_OPERATION_UNSUCCESSFUL,"ADOL-C function evaluation with error");
774  }
775  hessianW.resize(activeVariablesW.size(),activeVariablesW.size(),false);
776  hessianW.clear();
777  vector<double*> hessian_w(activeVariablesW.size());
778  for(unsigned int dd = 0;dd<activeVariablesW.size();dd++) {
779  hessian_w[dd] = &hessianW(dd,0);
780  }
781  r = hessian(tAgs[0],activeVariablesW.size(),&activeVariablesW[0],&hessian_w[0]);
782  if(r<adloc_return_value) {
783  SETERRQ(PETSC_COMM_SELF,MOFEM_OPERATION_UNSUCCESSFUL,"ADOL-C function evaluation with error");
784  }
785  sTress = VectorAdaptor(6,ublas::shallow_array_adaptor<double>(6,&gradientW[0]));
787  internalVariables.size(),
788  ublas::shallow_array_adaptor<double>(
789  internalVariables.size(),
790  &gradientW[12]
791  )
792  );
793  C.resize(6,6,false);
794  for(int ii = 0;ii<6;ii++) {
795  for(int jj = 0;jj<=ii;jj++) {
796  C(ii,jj) = hessianW(ii,jj);
797  }
798  }
799  partialWStrainPlasticStrain.resize(6,6,false);
800  for(int ii = 0;ii<6;ii++) {
801  for(int jj = 0;jj<6;jj++) {
802  partialWStrainPlasticStrain(ii,jj) = hessianW(6+jj,ii);
803  }
804  }
805  D.resize(internalVariables.size(),internalVariables.size(),false);
806  for(unsigned int ii = 0;ii<internalVariables.size();ii++) {
807  for(unsigned int jj = 0;jj<=ii;jj++) {
808  D(ii,jj) = hessianW(12+ii,12+jj);
809  }
810  }
811  PetscFunctionReturn(0);
812 }

◆ pLayW_NoHessian()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::pLayW_NoHessian ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 813 of file SmallStrainPlasticity.cpp.

813  {
814  PetscFunctionBegin;
815  int r;
816  int adloc_return_value = 0;
817  r = ::function(tAgs[0],1,activeVariablesW.size(),&activeVariablesW[0],&w);
818  if(r<adloc_return_value) {
819  SETERRQ(PETSC_COMM_SELF,MOFEM_OPERATION_UNSUCCESSFUL,"ADOL-C function evaluation with error");
820  }
821  gradientW.resize(activeVariablesW.size(),false);
822  r = gradient(tAgs[0],activeVariablesW.size(),&activeVariablesW[0],&gradientW[0]);
823  if(r<adloc_return_value) {
824  SETERRQ(PETSC_COMM_SELF,MOFEM_OPERATION_UNSUCCESSFUL,"ADOL-C function evaluation with error");
825  }
826  sTress = VectorAdaptor(6,ublas::shallow_array_adaptor<double>(6,&gradientW[0]));
828  internalVariables.size(),
829  ublas::shallow_array_adaptor<double>(
830  internalVariables.size(),
831  &gradientW[12]
832  )
833  );
834  PetscFunctionReturn(0);
835 }

◆ pLayY()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::pLayY ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 837 of file SmallStrainPlasticity.cpp.

837  {
838  PetscFunctionBegin;
839  int r;
840  int adloc_return_value = 0;
841  r = ::function(tAgs[1],1,activeVariablesYH.size(),&activeVariablesYH[0],&y);
842  if(r<adloc_return_value) {
843  SETERRQ(PETSC_COMM_SELF,MOFEM_OPERATION_UNSUCCESSFUL,"ADOL-C function evaluation with error");
844  }
845  gradientY.resize(activeVariablesYH.size());
846  r = gradient(tAgs[1],activeVariablesYH.size(),&activeVariablesYH[0],&gradientY[0]);
847  if(r<adloc_return_value) {
848  SETERRQ(PETSC_COMM_SELF,MOFEM_OPERATION_UNSUCCESSFUL,"ADOL-C function evaluation with error");
849  }
850  partialYSigma = VectorAdaptor(6,ublas::shallow_array_adaptor<double>(6,&gradientY[0]));
852  internalVariables.size(),
853  ublas::shallow_array_adaptor<double>(
854  internalVariables.size(),
855  &gradientY[6]
856  )
857  );
858  PetscFunctionReturn(0);
859 }

◆ pLayY_NoGradient()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::pLayY_NoGradient ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 860 of file SmallStrainPlasticity.cpp.

860  {
861  PetscFunctionBegin;
862  int r;
863  int adloc_return_value = 0;
864  r = ::function(tAgs[1],1,activeVariablesYH.size(),&activeVariablesYH[0],&y);
865  if(r<adloc_return_value) {
866  SETERRQ(PETSC_COMM_SELF,MOFEM_OPERATION_UNSUCCESSFUL,"ADOL-C function evaluation with error");
867  }
868  PetscFunctionReturn(0);
869 }

◆ rEcordH()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::rEcordH ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 740 of file SmallStrainPlasticity.cpp.

740  {
741  PetscFunctionBegin;
742  trace_on(tAgs[2]);
743  {
744  //active variables
745  a_sTress.resize(6,false);
746  for(int ii = 0;ii<6;ii++) {
747  a_sTress[ii] <<= sTress[ii];
748  }
749  a_internalFluxes.resize(internalVariables.size(),false);
750  for(unsigned int ii = 0;ii<internalVariables.size();ii++) {
751  a_internalFluxes[ii] <<= internalFluxes[ii];
752  }
753  //evaluete functions
754  ierr = flowPotential(); CHKERRQ(ierr);
755  //return variables
756  a_h >>= h;
757  }
758  trace_off();
759  PetscFunctionReturn(0);
760 }

◆ rEcordW()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::rEcordW ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 692 of file SmallStrainPlasticity.cpp.

692  {
693  PetscFunctionBegin;
694  trace_on(tAgs[0]);
695  {
696  //active variables
697  a_sTrain.resize(6,false);
698  for(int ii = 0;ii<6;ii++) {
699  a_sTrain[ii] <<= sTrain[ii];
700  }
701  a_plasticStrain.resize(6,false);
702  for(int ii = 0;ii<6;ii++) {
703  a_plasticStrain[ii] <<= plasticStrain[ii];
704  }
705  a_internalVariables.resize(internalVariables.size(),false);
706  for(unsigned int ii = 0;ii<internalVariables.size();ii++) {
708  }
709  //evaluete functions
710  ierr = freeHemholtzEnergy(); CHKERRQ(ierr);
711  //return variables
712  a_w >>= w;
713  }
714  trace_off();
715  PetscFunctionReturn(0);
716 }

◆ rEcordY()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::rEcordY ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 718 of file SmallStrainPlasticity.cpp.

718  {
719  PetscFunctionBegin;
720  trace_on(tAgs[1]);
721  {
722  //active variables
723  a_sTress.resize(6,false);
724  for(int ii = 0;ii<6;ii++) {
725  a_sTress[ii] <<= sTress[ii];
726  }
727  a_internalFluxes.resize(internalVariables.size(),false);
728  for(unsigned int ii = 0;ii<internalVariables.size();ii++) {
729  a_internalFluxes[ii] <<= internalFluxes[ii];
730  }
731  //evaluete functions
732  ierr = yieldFunction(); CHKERRQ(ierr);
733  //return variables
734  a_y >>= y;
735  }
736  trace_off();
737  PetscFunctionReturn(0);
738 }

◆ setActiveVariablesW()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::setActiveVariablesW ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 665 of file SmallStrainPlasticity.cpp.

665  {
666  PetscFunctionBegin;
667  activeVariablesW.resize(12+internalVariables.size(),false);
668  for(int ii = 0;ii<6;ii++) {
669  activeVariablesW[ii] = sTrain[ii];
670  }
671  for(int ii = 0;ii<6;ii++) {
672  activeVariablesW[6+ii] = plasticStrain[ii];
673  }
674  for(unsigned int ii = 0;ii<internalVariables.size();ii++) {
675  activeVariablesW[12+ii] = internalVariables[ii];
676  }
677  PetscFunctionReturn(0);
678 }

◆ setActiveVariablesYH()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::setActiveVariablesYH ( )
virtual
Examples
SmallStrainPlasticity.hpp.

Definition at line 680 of file SmallStrainPlasticity.cpp.

680  {
681  PetscFunctionBegin;
682  activeVariablesYH.resize(6+internalVariables.size(),false);
683  for(int ii = 0;ii<6;ii++) {
684  activeVariablesYH[ii] = sTress[ii];
685  }
686  for(unsigned int ii = 0;ii<internalVariables.size();ii++) {
687  activeVariablesYH[6+ii] = internalFluxes[ii];
688  }
689  PetscFunctionReturn(0);
690 }

◆ snesCreate()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::snesCreate ( )
Examples
SmallStrainPlasticity.hpp.

Definition at line 1027 of file SmallStrainPlasticity.cpp.

1027  {
1028  PetscFunctionBegin;
1029  ierr = SNESCreate(PETSC_COMM_SELF,&sNes); CHKERRQ(ierr);
1030  ierr = SNESSetFunction(sNes,R,SmallStrainPlasticityfRes,(void*)this); CHKERRQ(ierr);
1031  ierr = SNESSetJacobian(sNes,A,A,SmallStrainPlasticityfJac,(void*)this); CHKERRQ(ierr);
1032  ierr = SNESGetKSP(sNes,&kSp); CHKERRQ(ierr);
1033  ierr = KSPGetPC(kSp,&pC); CHKERRQ(ierr);
1034  ierr = KSPSetType(kSp,KSPPREONLY); CHKERRQ(ierr);
1035  ierr = PCSetType(pC,PCLU); CHKERRQ(ierr);
1036  ierr = SNESSetTolerances(sNes,tOl,1e-12,0,20,1000); CHKERRQ(ierr);
1037  ierr = SNESGetLineSearch(sNes,&lineSearch); CHKERRQ(ierr); CHKERRQ(ierr);
1038  //ierr = SNESLineSearchSetType(lineSearch,SNESLINESEARCHBASIC); CHKERRQ(ierr);
1039  ierr = SNESLineSearchSetType(lineSearch,SNESLINESEARCHBT); CHKERRQ(ierr);
1040 
1041  PetscFunctionReturn(0);
1042 }

◆ snesDestroy()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::snesDestroy ( )
Examples
SmallStrainPlasticity.hpp.

Definition at line 1044 of file SmallStrainPlasticity.cpp.

1044  {
1045  PetscFunctionBegin;
1046  ierr = SNESDestroy(&sNes); CHKERRQ(ierr);
1047  PetscFunctionReturn(0);
1048 }

◆ solveColasetProjection()

PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::solveColasetProjection ( )
Examples
SmallStrainPlasticity.hpp.

Definition at line 1050 of file SmallStrainPlasticity.cpp.

1050  {
1051  PetscFunctionBegin;
1052  {
1053  double *array;
1054  ierr = VecGetArray(Chi,&array); CHKERRQ(ierr);
1055  for(int ii = 0;ii<6;ii++) {
1056  array[ii] = plasticStrain[ii];
1057  }
1058  int nb_internal_variables = internalVariables.size();
1059  for(unsigned int ii = 0;ii<nb_internal_variables;ii++) {
1060  array[6+ii] = internalVariables[ii];
1061  }
1062  deltaGamma = 0;
1063  ierr = VecRestoreArray(Chi,&array); CHKERRQ(ierr);
1064  }
1065  ierr = SNESSolve(sNes,PETSC_NULL,Chi); CHKERRQ(ierr);
1066  SNESConvergedReason reason;
1067  ierr = SNESGetConvergedReason(sNes,&reason); CHKERRQ(ierr);
1068  if(reason < 0) {
1069  int its;
1070  ierr = SNESGetIterationNumber(sNes,&its); CHKERRQ(ierr);
1071 // ierr = PetscPrintf(
1072 // PETSC_COMM_SELF,
1073 // "** WARNING Plasticity Closet Point Projection - number of Newton iterations = %D\n",
1074 // its
1075 // ); CHKERRQ(ierr);
1076  noalias(plasticStrain) = plasticStrain0;
1078  deltaGamma = 0;
1079  } else {
1080  double *array;
1081  ierr = VecGetArray(Chi,&array); CHKERRQ(ierr);
1082  for(int ii = 0;ii<6;ii++) {
1083  plasticStrain[ii] = array[ii];
1084  }
1085  int nb_internal_variables = internalVariables.size();
1086  for(unsigned int ii = 0;ii<nb_internal_variables;ii++) {
1087  internalVariables[ii] = array[6+ii];
1088  }
1089  deltaGamma = array[6+nb_internal_variables];
1090  ierr = VecRestoreArray(Chi,&array); CHKERRQ(ierr);
1091  }
1092  PetscFunctionReturn(0);
1093 }

◆ yieldFunction()

virtual PetscErrorCode SmallStrainPlasticity::ClosestPointProjection::yieldFunction ( )
inlinevirtual

Reimplemented in SmallStrainParaboloidalPlasticity, and SmallStrainJ2Plasticity.

Examples
SmallStrainPlasticity.hpp.

Definition at line 581 of file SmallStrainPlasticity.hpp.

581  {
582  PetscFunctionBegin;
583  SETERRQ(
584  PETSC_COMM_SELF,
586  "Not implemented (overload class and implement this function)"
587  );
588  PetscFunctionReturn(0);
589  }

Friends And Related Function Documentation

◆ SmallStrainPlasticityfJac

PetscErrorCode SmallStrainPlasticityfJac ( SNES  ,
Vec  ,
Mat  ,
Mat  ,
void *  ctx 
)
friend
Examples
SmallStrainPlasticity.hpp.

Definition at line 1207 of file SmallStrainPlasticity.cpp.

1207  {
1208  PetscFunctionBegin;
1211 
1212  {
1213  #if PETSC_VERSION_GE(3,6,0)
1214  const double *array;
1215  ierr = VecGetArrayRead(chi,&array); CHKERRQ(ierr);
1216  #else
1217  double *array;
1218  ierr = VecGetArray(chi,&array); CHKERRQ(ierr);
1219  #endif
1220  for(int ii = 0;ii<6;ii++) {
1221  cp->plasticStrain[ii] = array[ii];
1222  }
1223  for(unsigned int ii = 0;ii<cp->internalVariables.size();ii++) {
1224  cp->internalVariables[ii] = array[6+ii];
1225  }
1226  cp->deltaGamma = array[6+cp->internalVariables.size()];
1227  #if PETSC_VERSION_GE(3,6,0)
1228  ierr = VecRestoreArrayRead(chi,&array); CHKERRQ(ierr);
1229  #else
1230  ierr = VecRestoreArray(chi,&array); CHKERRQ(ierr);
1231  #endif
1232  }
1233  ierr = cp->evaluatePotentials(); CHKERRQ(ierr);
1234  ierr = cp->cAlculateA(); CHKERRQ(ierr);
1235  PetscFunctionReturn(0);
1236 }

◆ SmallStrainPlasticityfRes

PetscErrorCode SmallStrainPlasticityfRes ( SNES  ,
Vec  ,
Vec  ,
void *  ctx 
)
friend
Examples
SmallStrainPlasticity.hpp.

Definition at line 1154 of file SmallStrainPlasticity.cpp.

1154  {
1155  PetscFunctionBegin;
1158 
1159  //ierr = VecView(chi,PETSC_VIEWER_STDOUT_WORLD); CHKERRQ(ierr);
1160  {
1161  #if PETSC_VERSION_GE(3,6,0)
1162  const double *array;
1163  ierr = VecGetArrayRead(chi,&array); CHKERRQ(ierr);
1164  #else
1165  const double *array;
1166  ierr = VecGetArray(chi,&array); CHKERRQ(ierr);
1167  #endif
1168  for(int ii = 0;ii<6;ii++) {
1169  cp->plasticStrain[ii] = array[ii];
1170  }
1171  for(unsigned int ii = 0;ii<cp->internalVariables.size();ii++) {
1172  cp->internalVariables[ii] = array[6+ii];
1173  }
1174  cp->deltaGamma = array[6+cp->internalVariables.size()];
1175  #if PETSC_VERSION_GE(3,6,0)
1176  ierr = VecRestoreArrayRead(chi,&array); CHKERRQ(ierr);
1177  #else
1178  ierr = VecRestoreArray(chi,&array); CHKERRQ(ierr);
1179  #endif
1180  }
1181  ierr = cp->evaluatePotentials(); CHKERRQ(ierr);
1182  ierr = cp->cAlculateR(r); CHKERRQ(ierr);
1183  // ierr = VecView(r,PETSC_VIEWER_STDOUT_WORLD); CHKERRQ(ierr);
1184  // cout << "sTress " << cp->sTress << endl;
1185  // cout << "internalVariables " << cp->internalVariables << endl;
1186  // cout << "C " << cp->C << endl;
1187  // cout << "D " << cp->D << endl;
1188  // cout << "y " << cp->y << endl;
1189  // cout << "h " << cp->h << endl;
1190  // cout << "partialYSigma " << cp->partialYSigma << endl;
1191  // cout << "partialYFlux " << cp->partialYFlux << endl;
1192  // cout << "partialHSigma " << cp->partialHSigma << endl;
1193  // cout << "partialHFlux " << cp->partialHFlux << endl;
1194  // cout << "partial2HSigma " << cp->partial2HSigma << endl;
1195  // cout << "partial2HFlux " << cp->partial2HFlux << endl;
1196  // cout << "partial2HSigmaFlux " << cp->partial2HSigmaFlux << endl;
1197  // cout << "C " << cp->C << endl;
1198  // cout << "D " << cp->D << endl;
1199  // cout << "partialWStrainPlasticStrain " << cp->partialWStrainPlasticStrain << endl;
1200  // cout << "plasticStrain " << cp->plasticStrain << endl;
1201  // cout << "internalVariables " << cp->internalVariables << endl;
1202  // cout << "delta plasticStrain " << cp->plasticStrain-cp->plasticStrain0 << endl;
1203  // cout << "delta internalVariables " << cp->internalVariables-cp->internalVariables << endl;
1204  PetscFunctionReturn(0);
1205 }

Member Data Documentation

◆ A

Mat SmallStrainPlasticity::ClosestPointProjection::A
Examples
SmallStrainPlasticity.hpp.

Definition at line 540 of file SmallStrainPlasticity.hpp.

◆ a_h

adouble SmallStrainPlasticity::ClosestPointProjection::a_h

◆ a_internalFluxes

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

◆ a_internalVariables

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

◆ a_plasticStrain

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

◆ a_sTrain

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

◆ a_sTress

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

◆ a_w

adouble SmallStrainPlasticity::ClosestPointProjection::a_w

◆ a_y

adouble SmallStrainPlasticity::ClosestPointProjection::a_y

◆ activeVariablesW

VectorDouble SmallStrainPlasticity::ClosestPointProjection::activeVariablesW
Examples
SmallStrainPlasticity.hpp.

Definition at line 510 of file SmallStrainPlasticity.hpp.

◆ activeVariablesYH

VectorDouble SmallStrainPlasticity::ClosestPointProjection::activeVariablesYH
Examples
SmallStrainPlasticity.hpp.

Definition at line 513 of file SmallStrainPlasticity.hpp.

◆ C

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

Definition at line 496 of file SmallStrainPlasticity.hpp.

◆ Cep

MatrixDouble SmallStrainPlasticity::ClosestPointProjection::Cep
Examples
SmallStrainPlasticity.hpp.

Definition at line 566 of file SmallStrainPlasticity.hpp.

◆ Chi

Vec SmallStrainPlasticity::ClosestPointProjection::Chi
Examples
SmallStrainPlasticity.hpp.

Definition at line 542 of file SmallStrainPlasticity.hpp.

◆ Cp

MatrixDouble SmallStrainPlasticity::ClosestPointProjection::Cp
Examples
SmallStrainPlasticity.hpp.

Definition at line 566 of file SmallStrainPlasticity.hpp.

◆ D

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

Definition at line 496 of file SmallStrainPlasticity.hpp.

◆ dataA

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

Definition at line 541 of file SmallStrainPlasticity.hpp.

◆ dChi

Vec SmallStrainPlasticity::ClosestPointProjection::dChi
Examples
SmallStrainPlasticity.hpp.

Definition at line 542 of file SmallStrainPlasticity.hpp.

◆ deltaGamma

double SmallStrainPlasticity::ClosestPointProjection::deltaGamma
Examples
SmallStrainPlasticity.hpp.

Definition at line 490 of file SmallStrainPlasticity.hpp.

◆ Ep

MatrixDouble SmallStrainPlasticity::ClosestPointProjection::Ep
Examples
SmallStrainPlasticity.hpp.

Definition at line 566 of file SmallStrainPlasticity.hpp.

◆ gG

int SmallStrainPlasticity::ClosestPointProjection::gG
Examples
SmallStrainPlasticity.hpp.

Definition at line 492 of file SmallStrainPlasticity.hpp.

◆ gradientH

VectorDouble SmallStrainPlasticity::ClosestPointProjection::gradientH
Examples
SmallStrainPlasticity.hpp.

Definition at line 536 of file SmallStrainPlasticity.hpp.

◆ gradientW

VectorDouble SmallStrainPlasticity::ClosestPointProjection::gradientW
Examples
SmallStrainPlasticity.hpp.

Definition at line 527 of file SmallStrainPlasticity.hpp.

◆ gradientY

VectorDouble SmallStrainPlasticity::ClosestPointProjection::gradientY
Examples
SmallStrainPlasticity.hpp.

Definition at line 532 of file SmallStrainPlasticity.hpp.

◆ h

double SmallStrainPlasticity::ClosestPointProjection::h
Examples
SmallStrainPlasticity.hpp.

Definition at line 516 of file SmallStrainPlasticity.hpp.

◆ hessianH

MatrixDouble SmallStrainPlasticity::ClosestPointProjection::hessianH
Examples
SmallStrainPlasticity.hpp.

Definition at line 537 of file SmallStrainPlasticity.hpp.

◆ hessianW

MatrixDouble SmallStrainPlasticity::ClosestPointProjection::hessianW
Examples
SmallStrainPlasticity.hpp.

Definition at line 528 of file SmallStrainPlasticity.hpp.

◆ internalFluxes

VectorAdaptor SmallStrainPlasticity::ClosestPointProjection::internalFluxes
Examples
SmallStrainPlasticity.hpp.

Definition at line 495 of file SmallStrainPlasticity.hpp.

◆ internalVariables

VectorDouble SmallStrainPlasticity::ClosestPointProjection::internalVariables

◆ internalVariables0

VectorDouble SmallStrainPlasticity::ClosestPointProjection::internalVariables0
Examples
SmallStrainPlasticity.hpp.

Definition at line 487 of file SmallStrainPlasticity.hpp.

◆ kSp

KSP SmallStrainPlasticity::ClosestPointProjection::kSp
Examples
SmallStrainPlasticity.hpp.

Definition at line 557 of file SmallStrainPlasticity.hpp.

◆ lineSearch

SNESLineSearch SmallStrainPlasticity::ClosestPointProjection::lineSearch
Examples
SmallStrainPlasticity.hpp.

Definition at line 559 of file SmallStrainPlasticity.hpp.

◆ Lp

MatrixDouble SmallStrainPlasticity::ClosestPointProjection::Lp
Examples
SmallStrainPlasticity.hpp.

Definition at line 566 of file SmallStrainPlasticity.hpp.

◆ partial2HFlux

ublas::symmetric_matrix<double,ublas::lower> SmallStrainPlasticity::ClosestPointProjection::partial2HFlux
Examples
SmallStrainPlasticity.hpp.

Definition at line 503 of file SmallStrainPlasticity.hpp.

◆ partial2HSigma

ublas::symmetric_matrix<double,ublas::lower> SmallStrainPlasticity::ClosestPointProjection::partial2HSigma
Examples
SmallStrainPlasticity.hpp.

Definition at line 502 of file SmallStrainPlasticity.hpp.

◆ partial2HSigmaFlux

MatrixDouble SmallStrainPlasticity::ClosestPointProjection::partial2HSigmaFlux
Examples
SmallStrainPlasticity.hpp.

Definition at line 504 of file SmallStrainPlasticity.hpp.

◆ partialHFlux

VectorAdaptor SmallStrainPlasticity::ClosestPointProjection::partialHFlux
Examples
SmallStrainPlasticity.hpp.

Definition at line 501 of file SmallStrainPlasticity.hpp.

◆ partialHSigma

VectorAdaptor SmallStrainPlasticity::ClosestPointProjection::partialHSigma
Examples
SmallStrainPlasticity.hpp.

Definition at line 500 of file SmallStrainPlasticity.hpp.

◆ partialWStrainPlasticStrain

MatrixDouble SmallStrainPlasticity::ClosestPointProjection::partialWStrainPlasticStrain
Examples
SmallStrainPlasticity.hpp.

Definition at line 497 of file SmallStrainPlasticity.hpp.

◆ partialYFlux

VectorAdaptor SmallStrainPlasticity::ClosestPointProjection::partialYFlux
Examples
SmallStrainPlasticity.hpp.

Definition at line 499 of file SmallStrainPlasticity.hpp.

◆ partialYSigma

VectorAdaptor SmallStrainPlasticity::ClosestPointProjection::partialYSigma
Examples
SmallStrainPlasticity.hpp.

Definition at line 498 of file SmallStrainPlasticity.hpp.

◆ pC

PC SmallStrainPlasticity::ClosestPointProjection::pC
Examples
SmallStrainPlasticity.hpp.

Definition at line 558 of file SmallStrainPlasticity.hpp.

◆ plasticStrain

VectorDouble SmallStrainPlasticity::ClosestPointProjection::plasticStrain
Examples
SmallStrainPlasticity.hpp.

Definition at line 489 of file SmallStrainPlasticity.hpp.

◆ plasticStrain0

VectorDouble SmallStrainPlasticity::ClosestPointProjection::plasticStrain0
Examples
SmallStrainPlasticity.hpp.

Definition at line 488 of file SmallStrainPlasticity.hpp.

◆ R

Vec SmallStrainPlasticity::ClosestPointProjection::R
Examples
SmallStrainPlasticity.hpp.

Definition at line 542 of file SmallStrainPlasticity.hpp.

◆ sNes

SNES SmallStrainPlasticity::ClosestPointProjection::sNes
Examples
SmallStrainPlasticity.hpp.

Definition at line 556 of file SmallStrainPlasticity.hpp.

◆ sTrain

VectorDouble SmallStrainPlasticity::ClosestPointProjection::sTrain
Examples
SmallStrainPlasticity.hpp.

Definition at line 485 of file SmallStrainPlasticity.hpp.

◆ sTress

VectorAdaptor SmallStrainPlasticity::ClosestPointProjection::sTress
Examples
SmallStrainPlasticity.hpp.

Definition at line 494 of file SmallStrainPlasticity.hpp.

◆ tAgs

vector<int> SmallStrainPlasticity::ClosestPointProjection::tAgs
Examples
SmallStrainPlasticity.hpp.

Definition at line 506 of file SmallStrainPlasticity.hpp.

◆ tOl

double SmallStrainPlasticity::ClosestPointProjection::tOl
Examples
SmallStrainPlasticity.hpp.

Definition at line 491 of file SmallStrainPlasticity.hpp.

◆ w

double SmallStrainPlasticity::ClosestPointProjection::w
Examples
SmallStrainPlasticity.hpp.

Definition at line 516 of file SmallStrainPlasticity.hpp.

◆ y

double SmallStrainPlasticity::ClosestPointProjection::y
Examples
SmallStrainPlasticity.hpp.

Definition at line 516 of file SmallStrainPlasticity.hpp.

◆ Yp

VectorDouble SmallStrainPlasticity::ClosestPointProjection::Yp
Examples
SmallStrainPlasticity.hpp.

Definition at line 567 of file SmallStrainPlasticity.hpp.


The documentation for this struct was generated from the following files:
SmallStrainPlasticity::ClosestPointProjection::SmallStrainPlasticityfRes
friend PetscErrorCode SmallStrainPlasticityfRes(SNES, Vec, Vec, void *ctx)
Definition: SmallStrainPlasticity.cpp:1154
SmallStrainPlasticity::ClosestPointProjection::Cp
MatrixDouble Cp
Definition: SmallStrainPlasticity.hpp:566
SmallStrainPlasticity::ClosestPointProjection::evaluatePotentials
virtual PetscErrorCode evaluatePotentials()
Definition: SmallStrainPlasticity.cpp:945
SmallStrainPlasticity::ClosestPointProjection::plasticStrain0
VectorDouble plasticStrain0
Definition: SmallStrainPlasticity.hpp:488
SmallStrainPlasticity::ClosestPointProjection::a_plasticStrain
ublas::vector< adouble > a_plasticStrain
Definition: SmallStrainPlasticity.hpp:517
SmallStrainPlasticity::ClosestPointProjection::lineSearch
SNESLineSearch lineSearch
Definition: SmallStrainPlasticity.hpp:559
SmallStrainPlasticity::ClosestPointProjection::internalVariables0
VectorDouble internalVariables0
Definition: SmallStrainPlasticity.hpp:487
SmallStrainPlasticity::ClosestPointProjection::yieldFunction
virtual PetscErrorCode yieldFunction()
Definition: SmallStrainPlasticity.hpp:581
SmallStrainPlasticity::ClosestPointProjection::setActiveVariablesW
virtual PetscErrorCode setActiveVariablesW()
Definition: SmallStrainPlasticity.cpp:665
SmallStrainPlasticity::ClosestPointProjection::partialHFlux
VectorAdaptor partialHFlux
Definition: SmallStrainPlasticity.hpp:501
SmallStrainPlasticity::ClosestPointProjection::pLayH
virtual PetscErrorCode pLayH()
Definition: SmallStrainPlasticity.cpp:870
SmallStrainPlasticity::ClosestPointProjection::plasticStrain
VectorDouble plasticStrain
Definition: SmallStrainPlasticity.hpp:489
SmallStrainPlasticity::ClosestPointProjection::rEcordH
virtual PetscErrorCode rEcordH()
Definition: SmallStrainPlasticity.cpp:740
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
SmallStrainPlasticity::ClosestPointProjection::A
Mat A
Definition: SmallStrainPlasticity.hpp:540
SmallStrainPlasticity::ClosestPointProjection
Closest point projection algorithm.
Definition: SmallStrainPlasticity.hpp:482
SmallStrainPlasticity::ClosestPointProjection::partialHSigma
VectorAdaptor partialHSigma
Definition: SmallStrainPlasticity.hpp:500
SmallStrainPlasticity::ClosestPointProjection::rEcordY
virtual PetscErrorCode rEcordY()
Definition: SmallStrainPlasticity.cpp:718
SmallStrainPlasticity::ClosestPointProjection::Ep
MatrixDouble Ep
Definition: SmallStrainPlasticity.hpp:566
SmallStrainPlasticity::ClosestPointProjection::hessianH
MatrixDouble hessianH
Definition: SmallStrainPlasticity.hpp:537
SmallStrainPlasticity::ClosestPointProjection::h
double h
Definition: SmallStrainPlasticity.hpp:516
SmallStrainPlasticity::ClosestPointProjection::tAgs
vector< int > tAgs
Definition: SmallStrainPlasticity.hpp:506
SmallStrainPlasticity::ClosestPointProjection::hessianW
MatrixDouble hessianW
Definition: SmallStrainPlasticity.hpp:528
SmallStrainPlasticity::ClosestPointProjection::D
ublas::symmetric_matrix< double, ublas::lower > D
Definition: SmallStrainPlasticity.hpp:496
SmallStrainPlasticity::ClosestPointProjection::dataA
ublas::matrix< double, ublas::column_major > dataA
Definition: SmallStrainPlasticity.hpp:541
SmallStrainPlasticity::ClosestPointProjection::partialYSigma
VectorAdaptor partialYSigma
Definition: SmallStrainPlasticity.hpp:498
SmallStrainPlasticity::ClosestPointProjection::partial2HFlux
ublas::symmetric_matrix< double, ublas::lower > partial2HFlux
Definition: SmallStrainPlasticity.hpp:503
SmallStrainPlasticity::ClosestPointProjection::gG
int gG
Definition: SmallStrainPlasticity.hpp:492
sdf.r
int r
Definition: sdf.py:8
SmallStrainPlasticity::ClosestPointProjection::SmallStrainPlasticityfJac
friend PetscErrorCode SmallStrainPlasticityfJac(SNES, Vec, Mat, Mat, void *ctx)
Definition: SmallStrainPlasticity.cpp:1207
SmallStrainPlasticity::ClosestPointProjection::a_internalFluxes
ublas::vector< adouble > a_internalFluxes
Definition: SmallStrainPlasticity.hpp:518
SmallStrainPlasticity::ClosestPointProjection::internalFluxes
VectorAdaptor internalFluxes
Definition: SmallStrainPlasticity.hpp:495
SmallStrainPlasticity::ClosestPointProjection::sNes
SNES sNes
Definition: SmallStrainPlasticity.hpp:556
SmallStrainPlasticity::ClosestPointProjection::gradientY
VectorDouble gradientY
Definition: SmallStrainPlasticity.hpp:532
SmallStrainPlasticity::ClosestPointProjection::tOl
double tOl
Definition: SmallStrainPlasticity.hpp:491
SmallStrainPlasticity::ClosestPointProjection::partial2HSigmaFlux
MatrixDouble partial2HSigmaFlux
Definition: SmallStrainPlasticity.hpp:504
SmallStrainPlasticity::ClosestPointProjection::partialWStrainPlasticStrain
MatrixDouble partialWStrainPlasticStrain
Definition: SmallStrainPlasticity.hpp:497
SmallStrainPlasticity::ClosestPointProjection::Lp
MatrixDouble Lp
Definition: SmallStrainPlasticity.hpp:566
SmallStrainPlasticity::ClosestPointProjection::Chi
Vec Chi
Definition: SmallStrainPlasticity.hpp:542
SmallStrainPlasticity::ClosestPointProjection::y
double y
Definition: SmallStrainPlasticity.hpp:516
SmallStrainPlasticity::ClosestPointProjection::partialYFlux
VectorAdaptor partialYFlux
Definition: SmallStrainPlasticity.hpp:499
SmallStrainPlasticity::ClosestPointProjection::activeVariablesYH
VectorDouble activeVariablesYH
Definition: SmallStrainPlasticity.hpp:513
SmallStrainPlasticity::ClosestPointProjection::activeVariablesW
VectorDouble activeVariablesW
Definition: SmallStrainPlasticity.hpp:510
MOFEM_OPERATION_UNSUCCESSFUL
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition: definitions.h:34
SmallStrainPlasticity::ClosestPointProjection::rEcordW
virtual PetscErrorCode rEcordW()
Definition: SmallStrainPlasticity.cpp:692
SmallStrainPlasticity::ClosestPointProjection::a_internalVariables
ublas::vector< adouble > a_internalVariables
Definition: SmallStrainPlasticity.hpp:517
MoFEM::Types::VectorAdaptor
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
SmallStrainPlasticity::ClosestPointProjection::partial2HSigma
ublas::symmetric_matrix< double, ublas::lower > partial2HSigma
Definition: SmallStrainPlasticity.hpp:502
SmallStrainPlasticity::ClosestPointProjection::freeHemholtzEnergy
virtual PetscErrorCode freeHemholtzEnergy()
Definition: SmallStrainPlasticity.hpp:571
SmallStrainPlasticity::ClosestPointProjection::sTress
VectorAdaptor sTress
Definition: SmallStrainPlasticity.hpp:494
SmallStrainPlasticity::ClosestPointProjection::a_sTress
ublas::vector< adouble > a_sTress
Definition: SmallStrainPlasticity.hpp:518
SmallStrainPlasticity::ClosestPointProjection::w
double w
Definition: SmallStrainPlasticity.hpp:516
convert.n
n
Definition: convert.py:82
SmallStrainPlasticity::ClosestPointProjection::dChi
Vec dChi
Definition: SmallStrainPlasticity.hpp:542
SmallStrainPlasticity::ClosestPointProjection::internalVariables
VectorDouble internalVariables
Definition: SmallStrainPlasticity.hpp:486
SmallStrainPlasticity::ClosestPointProjection::a_w
adouble a_w
Definition: SmallStrainPlasticity.hpp:519
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
SmallStrainPlasticity::ClosestPointProjection::pLayW
virtual PetscErrorCode pLayW()
Definition: SmallStrainPlasticity.cpp:762
SmallStrainPlasticity::ClosestPointProjection::flowPotential
virtual PetscErrorCode flowPotential()
Definition: SmallStrainPlasticity.hpp:591
SmallStrainPlasticity::ClosestPointProjection::a_h
adouble a_h
Definition: SmallStrainPlasticity.hpp:519
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
SmallStrainPlasticity::ClosestPointProjection::kSp
KSP kSp
Definition: SmallStrainPlasticity.hpp:557
SmallStrainPlasticity::ClosestPointProjection::cAlculateR
virtual PetscErrorCode cAlculateR(Vec R)
Definition: SmallStrainPlasticity.cpp:962
SmallStrainPlasticity::ClosestPointProjection::setActiveVariablesYH
virtual PetscErrorCode setActiveVariablesYH()
Definition: SmallStrainPlasticity.cpp:680
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
SmallStrainPlasticity::ClosestPointProjection::a_sTrain
ublas::vector< adouble > a_sTrain
Definition: SmallStrainPlasticity.hpp:517
SmallStrainPlasticity::ClosestPointProjection::deltaGamma
double deltaGamma
Definition: SmallStrainPlasticity.hpp:490
SmallStrainPlasticity::ClosestPointProjection::gradientW
VectorDouble gradientW
Definition: SmallStrainPlasticity.hpp:527
SmallStrainPlasticity::ClosestPointProjection::pLayY
virtual PetscErrorCode pLayY()
Definition: SmallStrainPlasticity.cpp:837
SmallStrainPlasticity::ClosestPointProjection::pC
PC pC
Definition: SmallStrainPlasticity.hpp:558
SmallStrainPlasticity::ClosestPointProjection::gradientH
VectorDouble gradientH
Definition: SmallStrainPlasticity.hpp:536
SmallStrainPlasticity::ClosestPointProjection::cAlculateA
virtual PetscErrorCode cAlculateA()
Definition: SmallStrainPlasticity.cpp:978
SmallStrainPlasticity::ClosestPointProjection::a_y
adouble a_y
Definition: SmallStrainPlasticity.hpp:519
SmallStrainPlasticity::ClosestPointProjection::C
ublas::symmetric_matrix< double, ublas::lower > C
Definition: SmallStrainPlasticity.hpp:496
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
SmallStrainPlasticity::ClosestPointProjection::Cep
MatrixDouble Cep
Definition: SmallStrainPlasticity.hpp:566
SmallStrainPlasticity::ClosestPointProjection::Yp
VectorDouble Yp
Definition: SmallStrainPlasticity.hpp:567