• Main Page
  • Related Pages
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

LamentStress_Def.hpp

Go to the documentation of this file.
00001 //*****************************************************************//
00002 //    Albany 2.0:  Copyright 2012 Sandia Corporation               //
00003 //    This Software is released under the BSD license detailed     //
00004 //    in the file "license.txt" in the top-level Albany directory  //
00005 //*****************************************************************//
00006 
00007 #include <Intrepid_MiniTensor.h>
00008 #include "Teuchos_TestForException.hpp"
00009 #include "Phalanx_DataLayout.hpp"
00010 #include "Intrepid_FunctionSpaceTools.hpp"
00011 #include "QCAD_MaterialDatabase.hpp"
00012 
00013 using namespace std;
00014 
00015 namespace LCM {
00016 
00017 template<typename EvalT, typename Traits>
00018 LamentStress<EvalT, Traits>::
00019 LamentStress(Teuchos::ParameterList& p) :
00020   defGradField(p.get<std::string>("DefGrad Name"),
00021                p.get<Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout")),
00022   stressField(p.get<std::string>("Stress Name"),
00023               p.get<Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout") ),
00024   lamentMaterialModel(Teuchos::RCP<lament::Material<ScalarT> >())
00025 {
00026   // Pull out numQPs and numDims from a Layout
00027   tensor_dl = p.get< Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout");
00028   std::vector<PHX::DataLayout::size_type> dims;
00029   tensor_dl->dimensions(dims);
00030   numQPs  = dims[1];
00031   numDims = dims[2];
00032 
00033   TEUCHOS_TEST_FOR_EXCEPTION(this->numDims != 3, Teuchos::Exceptions::InvalidParameter, " LAMENT materials enabled only for three-dimensional analyses.");
00034 
00035   defGradName = p.get<std::string>("DefGrad Name")+"_old";
00036   this->addDependentField(defGradField);
00037 
00038   stressName = p.get<std::string>("Stress Name")+"_old";
00039   this->addEvaluatedField(stressField);
00040 
00041   this->setName("LamentStress"+PHX::TypeString<EvalT>::value);
00042 
00043   // Default to getting material info form base input file (possibley overwritten later)
00044   lamentMaterialModelName = p.get<string>("Lame Material Model", "Elastic");
00045   std::cout << "Material Model Name : " << lamentMaterialModelName << std::endl;
00046   Teuchos::ParameterList& lamentMaterialParameters = p.sublist("Lame Material Parameters");
00047 
00048   // Code to allow material data to come from materials.xml data file
00049   int haveMatDB = p.get<bool>("Have MatDB", false);
00050 
00051   std::string ebName = p.get<std::string>("Element Block Name", "Missing");
00052 
00053   // Check for material database file
00054   if (haveMatDB) {
00055     // Check if material database will be supplying the data
00056     bool dataFromDatabase = lamentMaterialParameters.get<bool>("Material Dependent Data Source",false);
00057 
00058     // If so, overwrite material model and data from database file
00059     if (dataFromDatabase) {
00060        Teuchos::RCP<QCAD::MaterialDatabase> materialDB = p.get< Teuchos::RCP<QCAD::MaterialDatabase> >("MaterialDB");
00061 
00062        lamentMaterialModelName = materialDB->getElementBlockParam<std::string>(ebName, "Lame Material Model");
00063        lamentMaterialParameters = materialDB->getElementBlockSublist(ebName, "Lame Material Parameters");
00064      }
00065   }
00066 
00067   // Initialize the LAMENT material model
00068   // This assumes that there is a single material model associated with this
00069   // evaluator and that the material properties are constant (read directly
00070   // from input deck parameter list)
00071   lamentMaterialModel = LameUtils::constructLamentMaterialModel<ScalarT>(lamentMaterialModelName, lamentMaterialParameters);
00072 
00073   // Get a list of the LAMENT material model state variable names
00074   lamentMaterialModelStateVariableNames = LameUtils::getStateVariableNames(lamentMaterialModelName, lamentMaterialParameters);
00075 
00076   // Declare the state variables as evaluated fields (type is always double)
00077   Teuchos::RCP<PHX::DataLayout> dataLayout = p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout");
00078   for(unsigned int i=0 ; i<lamentMaterialModelStateVariableNames.size() ; ++i){
00079     PHX::MDField<ScalarT,Cell,QuadPoint,Dim,Dim> lamentMaterialModelStateVariableField(lamentMaterialModelStateVariableNames[i], dataLayout);
00080     this->addEvaluatedField(lamentMaterialModelStateVariableField);
00081     lamentMaterialModelStateVariableFields.push_back(lamentMaterialModelStateVariableField);
00082   }
00083 }
00084 
00085 template<typename EvalT, typename Traits>
00086 void LamentStress<EvalT, Traits>::
00087 postRegistrationSetup(typename Traits::SetupData d,
00088                       PHX::FieldManager<Traits>& fm)
00089 {
00090   this->utils.setFieldData(defGradField,fm);
00091   this->utils.setFieldData(stressField,fm);
00092   for(unsigned int i=0 ; i<lamentMaterialModelStateVariableFields.size() ; ++i)
00093     this->utils.setFieldData(lamentMaterialModelStateVariableFields[i],fm);
00094 }
00095 
00096 template<typename EvalT, typename Traits>
00097 void LamentStress<EvalT, Traits>::
00098 evaluateFields(typename Traits::EvalData workset)
00099 {
00100   Teuchos::RCP<lament::matParams<ScalarT> > matp = Teuchos::rcp(new lament::matParams<ScalarT>());
00101 
00102   // Get the old state data
00103   Albany::MDArray oldDefGrad = (*workset.stateArrayPtr)[defGradName];
00104   Albany::MDArray oldStress = (*workset.stateArrayPtr)[stressName];
00105 
00106   int numStateVariables = (int)(this->lamentMaterialModelStateVariableNames.size());
00107 
00108   // \todo Get actual time step for calls to LAMENT materials.
00109   double deltaT = 1.0;
00110 
00111   vector<ScalarT> strainRate(6);              // symmetric tensor
00112   vector<ScalarT> spin(3);                     // skew-symmetric tensor
00113   vector<ScalarT> defGrad(9);              // symmetric tensor
00114   vector<ScalarT> leftStretch(6);              // symmetric tensor
00115   vector<ScalarT> rotation(9);                 // full tensor
00116   vector<double> stressOld(6);                // symmetric tensor
00117   vector<ScalarT> stressNew(6);               // symmetric tensor
00118   vector<double> stateOld(numStateVariables); // a single scalar for each state variable
00119   vector<double> stateNew(numStateVariables); // a single scalar for each state variable
00120 
00121   // \todo Set up scratch space for material models using getNumScratchVars() and setScratchPtr().
00122 
00123   // Create the matParams structure, which is passed to Lament
00124   matp->nelements = 1;
00125   matp->dt = deltaT;
00126   matp->time = 0.0;
00127   matp->strain_rate = &strainRate[0];
00128   matp->spin = &spin[0];
00129   matp->deformation_gradient = &defGrad[0];
00130   matp->left_stretch = &leftStretch[0];
00131   matp->rotation = &rotation[0];
00132   matp->state_old = &stateOld[0];
00133   matp->state_new = &stateNew[0];
00134   matp->stress_old = &stressOld[0];
00135   matp->stress_new = &stressNew[0];
00136 //   matp->dt_mat = std::numeric_limits<double>::max();
00137   
00138   // matParams that still need to be added:
00139   // matp->temp_old  (temperature)
00140   // matp->temp_new
00141   // matp->sound_speed_old
00142   // matp->sound_speed_new
00143   // matp->volume
00144   // scratch pointer
00145   // function pointers (lots to be done here)
00146 
00147   for (int cell=0; cell < (int)workset.numCells; ++cell) {
00148     for (int qp=0; qp < (int)numQPs; ++qp) {
00149 
00150       // std::cout << "QP: " << qp << std::endl;
00151 
00152       // Fill the following entries in matParams for call to LAMENT
00153       //
00154       // nelements     - number of elements 
00155       // dt            - time step, this one is tough because Albany does not currently have a concept of time step for implicit integration
00156       // time          - current time, again Albany does not currently have a concept of time for implicit integration
00157       // strain_rate   - what Sierra calls the rate of deformation, it is the symmetric part of the velocity gradient
00158       // spin          - anti-symmetric part of the velocity gradient
00159       // left_stretch  - found as V in the polar decomposition of the deformation gradient F = VR
00160       // rotation      - found as R in the polar decomposition of the deformation gradient F = VR
00161       // state_old     - material state data for previous time step (material dependent, none for lament::Elastic)
00162       // state_new     - material state data for current time step (material dependent, none for lament::Elastic)
00163       // stress_old    - stress at previous time step
00164       // stress_new    - stress at current time step, filled by material model
00165       //
00166       // The total deformation gradient is available as field data
00167       // 
00168       // The velocity gradient is not available but can be computed at the logarithm of the incremental deformation gradient divided by deltaT
00169       // The incremental deformation gradient is computed as F_new F_old^-1
00170 
00171       // JTO:  here is how I think this will go (of course the first two lines won't work as is...)
00172       // Intrepid::Tensor<RealType> F = newDefGrad;
00173       // Intrepid::Tensor<RealType> Fn = oldDefGrad;
00174       // Intrepid::Tensor<RealType> f = F*Intrepid::inverse(Fn);
00175       // Intrepid::Tensor<RealType> V;
00176       // Intrepid::Tensor<RealType> R;
00177       // boost::tie(V,R) = Intrepid::polar_left(F);
00178       // Intrepid::Tensor<RealType> Vinc;
00179       // Intrepid::Tensor<RealType> Rinc;
00180       // Intrepid::Tensor<RealType> logVinc;
00181       // boost::tie(Vinc,Rinc,logVinc) = Intrepid::polar_left_logV(f)
00182       // Intrepid::Tensor<RealType> logRinc = Intrepid::log_rotation(Rinc);
00183       // Intrepid::Tensor<RealType> logf = Intrepid::bch(logVinc,logRinc);
00184       // Intrepid::Tensor<RealType> L = (1.0/deltaT)*logf;
00185       // Intrepid::Tensor<RealType> D = Intrepid::sym(L);
00186       // Intrepid::Tensor<RealType> W = Intrepid::skew(L);
00187       // and then fill data into the vectors below
00188 
00189       // new deformation gradient (the current deformation gradient as computed in the current configuration)
00190       Intrepid::Tensor<ScalarT> Fnew( 3, &defGradField(cell,qp,0,0) );
00191 
00192       // old deformation gradient (deformation gradient at previous load step)
00193       Intrepid::Tensor<ScalarT> Fold( oldDefGrad(cell,qp,0,0), oldDefGrad(cell,qp,0,1), oldDefGrad(cell,qp,0,2),
00194                                     oldDefGrad(cell,qp,1,0), oldDefGrad(cell,qp,1,1), oldDefGrad(cell,qp,1,2),
00195                                     oldDefGrad(cell,qp,2,0), oldDefGrad(cell,qp,2,1), oldDefGrad(cell,qp,2,2) );
00196 
00197       // incremental deformation gradient
00198       Intrepid::Tensor<ScalarT> Finc = Fnew * Intrepid::inverse(Fold);
00199 
00200       
00201       // DEBUGGING //
00202       //if(cell==0 && qp==0){
00203       // std::cout << "Fnew(0,0) " << Fnew(0,0) << endl;
00204       // std::cout << "Fnew(1,0) " << Fnew(1,0) << endl;
00205       // std::cout << "Fnew(2,0) " << Fnew(2,0) << endl;
00206       // std::cout << "Fnew(0,1) " << Fnew(0,1) << endl;
00207       // std::cout << "Fnew(1,1) " << Fnew(1,1) << endl;
00208       // std::cout << "Fnew(2,1) " << Fnew(2,1) << endl;
00209       // std::cout << "Fnew(0,2) " << Fnew(0,2) << endl;
00210       // std::cout << "Fnew(1,2) " << Fnew(1,2) << endl;
00211       // std::cout << "Fnew(2,2) " << Fnew(2,2) << endl;
00212         //}
00213       // END DEBUGGING //
00214 
00215       // left stretch V, and rotation R, from left polar decomposition of new deformation gradient
00216       Intrepid::Tensor<ScalarT> V(3), R(3), U(3);
00217       boost::tie(V,R) = Intrepid::polar_left(Fnew);
00218       //V = R * U * transpose(R);
00219       
00220       // DEBUGGING //
00221       //if(cell==0 && qp==0){
00222       // std::cout << "U(0,0) " << U(0,0) << endl;
00223       // std::cout << "U(1,0) " << U(1,0) << endl;
00224       // std::cout << "U(2,0) " << U(2,0) << endl;
00225       // std::cout << "U(0,1) " << U(0,1) << endl;
00226       // std::cout << "U(1,1) " << U(1,1) << endl;
00227       // std::cout << "U(2,1) " << U(2,1) << endl;
00228       // std::cout << "U(0,2) " << U(0,2) << endl;
00229       // std::cout << "U(1,2) " << U(1,2) << endl;
00230       // std::cout << "U(2,2) " << U(2,2) << endl;
00231       // std::cout << "========\n";
00232       // std::cout << "V(0,0) " << V(0,0) << endl;
00233       // std::cout << "V(1,0) " << V(1,0) << endl;
00234       // std::cout << "V(2,0) " << V(2,0) << endl;
00235       // std::cout << "V(0,1) " << V(0,1) << endl;
00236       // std::cout << "V(1,1) " << V(1,1) << endl;
00237       // std::cout << "V(2,1) " << V(2,1) << endl;
00238       // std::cout << "V(0,2) " << V(0,2) << endl;
00239       // std::cout << "V(1,2) " << V(1,2) << endl;
00240       // std::cout << "V(2,2) " << V(2,2) << endl;
00241       // std::cout << "========\n";
00242       // std::cout << "R(0,0) " << R(0,0) << endl;
00243       // std::cout << "R(1,0) " << R(1,0) << endl;
00244       // std::cout << "R(2,0) " << R(2,0) << endl;
00245       // std::cout << "R(0,1) " << R(0,1) << endl;
00246       // std::cout << "R(1,1) " << R(1,1) << endl;
00247       // std::cout << "R(2,1) " << R(2,1) << endl;
00248       // std::cout << "R(0,2) " << R(0,2) << endl;
00249       // std::cout << "R(1,2) " << R(1,2) << endl;
00250       // std::cout << "R(2,2) " << R(2,2) << endl;
00251         //}
00252       // END DEBUGGING //
00253 
00254       // incremental left stretch Vinc, incremental rotation Rinc, and log of incremental left stretch, logVinc
00255       
00256       Intrepid::Tensor<ScalarT> Uinc(3), Vinc(3), Rinc(3), logVinc(3);
00257       //boost::tie(Vinc,Rinc,logVinc) = Intrepid::polar_left_logV(Finc);
00258       boost::tie(Vinc,Rinc) = Intrepid::polar_left(Finc);
00259       //Vinc = Rinc * Uinc * transpose(Rinc);
00260       logVinc = Intrepid::log(Vinc);
00261 
00262       // log of incremental rotation
00263       Intrepid::Tensor<ScalarT> logRinc = Intrepid::log_rotation(Rinc);
00264 
00265       // log of incremental deformation gradient
00266       Intrepid::Tensor<ScalarT> logFinc = Intrepid::bch(logVinc, logRinc);
00267 
00268       // velocity gradient
00269       Intrepid::Tensor<ScalarT> L = (1.0/deltaT)*logFinc;
00270 
00271       // strain rate (a.k.a rate of deformation)
00272       Intrepid::Tensor<ScalarT> D = Intrepid::sym(L);
00273 
00274       // spin
00275       Intrepid::Tensor<ScalarT> W = Intrepid::skew(L);
00276 
00277       // load everything into the Lament data structure
00278 
00279       strainRate[0] = ( D(0,0) );
00280       strainRate[1] = ( D(1,1) );
00281       strainRate[2] = ( D(2,2) );
00282       strainRate[3] = ( D(0,1) );
00283       strainRate[4] = ( D(1,2) );
00284       strainRate[5] = ( D(2,0) );
00285 
00286       spin[0] = ( W(0,1) );
00287       spin[1] = ( W(1,2) );
00288       spin[2] = ( W(2,0) );
00289 
00290       leftStretch[0] = ( V(0,0) );
00291       leftStretch[1] = ( V(1,1) );
00292       leftStretch[2] = ( V(2,2) );
00293       leftStretch[3] = ( V(0,1) );
00294       leftStretch[4] = ( V(1,2) );
00295       leftStretch[5] = ( V(2,0) );
00296 
00297       rotation[0] = ( R(0,0) );
00298       rotation[1] = ( R(1,1) );
00299       rotation[2] = ( R(2,2) );
00300       rotation[3] = ( R(0,1) );
00301       rotation[4] = ( R(1,2) );
00302       rotation[5] = ( R(2,0) );
00303       rotation[6] = ( R(1,0) );
00304       rotation[7] = ( R(2,1) );
00305       rotation[8] = ( R(0,2) );
00306 
00307       defGrad[0] = ( Fnew(0,0) );
00308       defGrad[1] = ( Fnew(1,1) );
00309       defGrad[2] = ( Fnew(2,2) );
00310       defGrad[3] = ( Fnew(0,1) );
00311       defGrad[4] = ( Fnew(1,2) );
00312       defGrad[5] = ( Fnew(2,0) );
00313       defGrad[6] = ( Fnew(1,0) );
00314       defGrad[7] = ( Fnew(2,1) );
00315       defGrad[8] = ( Fnew(0,2) );
00316 
00317       stressOld[0] = oldStress(cell,qp,0,0);
00318       stressOld[1] = oldStress(cell,qp,1,1);
00319       stressOld[2] = oldStress(cell,qp,2,2);
00320       stressOld[3] = oldStress(cell,qp,0,1);
00321       stressOld[4] = oldStress(cell,qp,1,2);
00322       stressOld[5] = oldStress(cell,qp,2,0);
00323 
00324       // copy data from the state manager to the LAMENT data structure
00325       for(int iVar=0 ; iVar<numStateVariables ; iVar++){
00326         const std::string& variableName = this->lamentMaterialModelStateVariableNames[iVar]+"_old";
00327         Albany::MDArray stateVar = (*workset.stateArrayPtr)[variableName];
00328         stateOld[iVar] = stateVar(cell,qp);
00329       }
00330 
00331       // Make a call to the LAMENT material model to initialize the load step
00332       this->lamentMaterialModel->loadStepInit(matp.get());
00333 
00334       // Get the stress from the LAMENT material
00335 
00336       // std::cout << "about to call lament->getStress()" << std::endl;
00337 
00338       this->lamentMaterialModel->getStress(matp.get());
00339 
00340       // std::cout << "after calling lament->getStress() 2" << std::endl;
00341 
00342       // rotate to get the Cauchy Stress
00343       Intrepid::Tensor<ScalarT> lameStress( stressNew[0], stressNew[3], stressNew[5],
00344                                           stressNew[3], stressNew[1], stressNew[4],
00345                                           stressNew[5], stressNew[4], stressNew[2] );
00346       Intrepid::Tensor<ScalarT> cauchy = R * lameStress * transpose(R);
00347 
00348       // DEBUGGING //
00349       //if(cell==0 && qp==0){
00350   // std::cout << "check strainRate[0] " << strainRate[0] << endl;
00351   // std::cout << "check strainRate[1] " << strainRate[1] << endl;
00352   // std::cout << "check strainRate[2] " << strainRate[2] << endl;
00353   // std::cout << "check strainRate[3] " << strainRate[3] << endl;
00354   // std::cout << "check strainRate[4] " << strainRate[4] << endl;
00355   // std::cout << "check strainRate[5] " << strainRate[5] << endl;
00356         //}
00357       // END DEBUGGING //
00358 
00359       // Copy the new stress into the stress field
00360       for (int i(0); i < 3; ++i)
00361         for (int j(0); j < 3; ++j)
00362           stressField(cell,qp,i,j) = cauchy(i,j);
00363 
00364       // stressField(cell,qp,0,0) = stressNew[0];
00365       // stressField(cell,qp,1,1) = stressNew[1];
00366       // stressField(cell,qp,2,2) = stressNew[2];
00367       // stressField(cell,qp,0,1) = stressNew[3];
00368       // stressField(cell,qp,1,2) = stressNew[4];
00369       // stressField(cell,qp,2,0) = stressNew[5];
00370       // stressField(cell,qp,1,0) = stressNew[3]; 
00371       // stressField(cell,qp,2,1) = stressNew[4]; 
00372       // stressField(cell,qp,0,2) = stressNew[5];
00373 
00374       // copy state_new data from the LAMENT data structure to the corresponding state variable field
00375       for(int iVar=0 ; iVar<numStateVariables ; iVar++)
00376   this->lamentMaterialModelStateVariableFields[iVar](cell,qp) = stateNew[iVar];
00377 
00378       // DEBUGGING //
00379       //if(cell==0 && qp==0){
00380       //   std::cout << "stress(0,0) " << this->stressField(cell,qp,0,0) << endl;
00381       //   std::cout << "stress(1,1) " << this->stressField(cell,qp,1,1) << endl;
00382       //   std::cout << "stress(2,2) " << this->stressField(cell,qp,2,2) << endl;
00383       //   std::cout << "stress(0,1) " << this->stressField(cell,qp,0,1) << endl;
00384       //   std::cout << "stress(1,2) " << this->stressField(cell,qp,1,2) << endl;
00385       //   std::cout << "stress(0,2) " << this->stressField(cell,qp,0,2) << endl;
00386       //   std::cout << "stress(1,0) " << this->stressField(cell,qp,1,0) << endl;
00387       //   std::cout << "stress(2,1) " << this->stressField(cell,qp,2,1) << endl;
00388       //   std::cout << "stress(2,0) " << this->stressField(cell,qp,2,0) << endl;
00389       //   //}
00390       // // END DEBUGGING //
00391 
00392     }
00393   }
00394 }
00395 
00396 }
00397 

Generated on Wed Mar 26 2014 18:36:39 for Albany: a Trilinos-based PDE code by  doxygen 1.7.1