int main(int ac, char* av[]) { KokkosGuard kokkos(ac, av); typedef PHX::MDField<PHAL::AlbanyTraits::Residual::ScalarT>::size_type size_type; typedef PHAL::AlbanyTraits::Residual Residual; typedef PHAL::AlbanyTraits::Residual::ScalarT ScalarT; typedef PHAL::AlbanyTraits Traits; std::cout.precision(15); // // Create a command line processor and parse command line options // Teuchos::CommandLineProcessor command_line_processor; command_line_processor.setDocString( "Material Point Simulator.\n" "For testing material models in LCM.\n"); std::string input_file = "materials.xml"; command_line_processor.setOption("input", &input_file, "Input File Name"); std::string timing_file = "timing.csv"; command_line_processor.setOption("timing", &timing_file, "Timing File Name"); int workset_size = 1; command_line_processor.setOption("wsize", &workset_size, "Workset Size"); int num_pts = 1; command_line_processor.setOption( "npoints", &num_pts, "Number of Gaussian Points"); size_t memlimit = 1024; // 1GB heap limit by default command_line_processor.setOption( "memlimit", &memlimit, "Heap memory limit in MB for CUDA kernels"); // Throw a warning and not error for unrecognized options command_line_processor.recogniseAllOptions(true); // Don't throw exceptions for errors command_line_processor.throwExceptions(false); // Parse command line Teuchos::CommandLineProcessor::EParseCommandLineReturn parse_return = command_line_processor.parse(ac, av); std::ofstream tout(timing_file.c_str()); if (parse_return == Teuchos::CommandLineProcessor::PARSE_HELP_PRINTED) { return 0; } if (parse_return != Teuchos::CommandLineProcessor::PARSE_SUCCESSFUL) { return 1; } util::TimeMonitor& tmonitor = util::PerformanceContext::instance().timeMonitor(); Teuchos::RCP<Teuchos::Time> total_time = tmonitor["MPS: Total Time"]; Teuchos::RCP<Teuchos::Time> compute_time = tmonitor["MPS: Compute Time"]; // // Process material.xml file // Read into materialDB and get material model name // // A mpi object must be instantiated before using the comm to read // material file Teuchos::GlobalMPISession mpi_session(&ac, &av); Teuchos::RCP<const Teuchos_Comm> commT = Albany::createTeuchosCommFromMpiComm(Albany_MPI_COMM_WORLD); Teuchos::RCP<Albany::MaterialDatabase> material_db; material_db = Teuchos::rcp(new Albany::MaterialDatabase(input_file, commT)); // Get the name of the material model to be used (and make sure there is one) std::string element_block_name = "Block0"; std::string material_model_name; material_model_name = material_db->getElementBlockSublist(element_block_name, "Material Model") .get<std::string>("Model Name"); TEUCHOS_TEST_FOR_EXCEPTION( material_model_name.length() == 0, std::logic_error, "A material model must be defined for block: " + element_block_name); // // Preloading stage setup // set up evaluators, create field and state managers // // Set up the data layout // const int workset_size = 1; const int num_dims = 3; const int num_vertices = 8; const int num_nodes = 8; const Teuchos::RCP<Albany::Layouts> dl = Teuchos::rcp(new Albany::Layouts( workset_size, num_vertices, num_nodes, num_pts, num_dims)); // create field name strings LCM::FieldNameMap field_name_map(false); Teuchos::RCP<std::map<std::string, std::string>> fnm = field_name_map.getMap(); //--------------------------------------------------------------------------- // Deformation gradient // initially set the deformation gradient to the identity Teuchos::ArrayRCP<ScalarT> def_grad(workset_size * num_pts * 9); for (int i = 0; i < workset_size; ++i) { for (int j = 0; j < num_pts; ++j) { int base = i * num_pts * 9 + j * 9; for (int k = 0; k < 9; ++k) def_grad[base + k] = 0.0; def_grad[base + 0] = 1.0; def_grad[base + 4] = 1.0; def_grad[base + 8] = 1.0; } } // SetField evaluator, which will be used to manually assign a value // to the def_grad field Teuchos::ParameterList setDefGradP("SetFieldDefGrad"); setDefGradP.set<std::string>("Evaluated Field Name", "F"); setDefGradP.set<Teuchos::RCP<PHX::DataLayout>>( "Evaluated Field Data Layout", dl->qp_tensor); setDefGradP.set<Teuchos::ArrayRCP<ScalarT>>("Field Values", def_grad); auto setFieldDefGrad = Teuchos::rcp(new LCM::SetField<Residual, Traits>(setDefGradP)); //--------------------------------------------------------------------------- // Det(deformation gradient) Teuchos::ArrayRCP<ScalarT> detdefgrad(workset_size * num_pts); for (int i = 0; i < workset_size * num_pts; ++i) detdefgrad[i] = 1.0; // SetField evaluator, which will be used to manually assign a value // to the detdefgrad field Teuchos::ParameterList setDetDefGradP("SetFieldDetDefGrad"); setDetDefGradP.set<std::string>("Evaluated Field Name", "J"); setDetDefGradP.set<Teuchos::RCP<PHX::DataLayout>>( "Evaluated Field Data Layout", dl->qp_scalar); setDetDefGradP.set<Teuchos::ArrayRCP<ScalarT>>("Field Values", detdefgrad); auto setFieldDetDefGrad = Teuchos::rcp(new LCM::SetField<Residual, Traits>(setDetDefGradP)); //--------------------------------------------------------------------------- // Small strain tensor // initially set the strain tensor to zeros Teuchos::ArrayRCP<ScalarT> strain(workset_size * num_pts * 9); for (int i = 0; i < workset_size; ++i) { for (int j = 0; j < num_pts; ++j) { int base = i * num_pts * 9 + j * 9; for (int k = 0; k < 9; ++k) strain[base + k] = 0.0; } } // SetField evaluator, which will be used to manually assign a value // to the strain field Teuchos::ParameterList setStrainP("SetFieldStrain"); setStrainP.set<std::string>("Evaluated Field Name", "Strain"); setStrainP.set<Teuchos::RCP<PHX::DataLayout>>( "Evaluated Field Data Layout", dl->qp_tensor); setStrainP.set<Teuchos::ArrayRCP<ScalarT>>("Field Values", strain); auto setFieldStrain = Teuchos::rcp(new LCM::SetField<Residual, Traits>(setStrainP)); //--------------------------------------------------------------------------- // Instantiate a field manager PHX::FieldManager<Traits> fieldManager; // Instantiate a field manager for States PHX::FieldManager<Traits> stateFieldManager; // Register the evaluators with the field manager fieldManager.registerEvaluator<Residual>(setFieldDefGrad); fieldManager.registerEvaluator<Residual>(setFieldDetDefGrad); fieldManager.registerEvaluator<Residual>(setFieldStrain); // Register the evaluators with the state field manager stateFieldManager.registerEvaluator<Residual>(setFieldDefGrad); stateFieldManager.registerEvaluator<Residual>(setFieldDetDefGrad); stateFieldManager.registerEvaluator<Residual>(setFieldStrain); // Instantiate a state manager Albany::StateManager stateMgr; // extract the Material ParameterList for use below std::string matName = material_db->getElementBlockParam<std::string>( element_block_name, "material"); Teuchos::ParameterList& paramList = material_db->getElementBlockSublist(element_block_name, matName); Teuchos::ParameterList& mpsParams = paramList.sublist("Material Point Simulator"); // Get loading parameters from .xml file std::string load_case = mpsParams.get<std::string>("Loading Case Name", "uniaxial"); int number_steps = mpsParams.get<int>("Number of Steps", 10); double step_size = mpsParams.get<double>("Step Size", 1.0e-2); std::cout << "Loading parameters:" << "\n number of steps: " << number_steps << "\n step_size : " << step_size << std::endl; // determine if temperature is being used bool have_temperature = mpsParams.get<bool>("Use Temperature", false); std::cout << "have_temp: " << have_temperature << std::endl; //--------------------------------------------------------------------------- // Temperature (optional) if (have_temperature) { Teuchos::ArrayRCP<ScalarT> temperature(workset_size); ScalarT temp = mpsParams.get<double>("Temperature", 1.0); for (int i = 0; i < workset_size * num_pts; ++i) temperature[0] = temp; // SetField evaluator, which will be used to manually assign a value // to the detdefgrad field Teuchos::ParameterList setTempP("SetFieldTemperature"); setTempP.set<std::string>("Evaluated Field Name", "Temperature"); setTempP.set<Teuchos::RCP<PHX::DataLayout>>( "Evaluated Field Data Layout", dl->qp_scalar); setTempP.set<Teuchos::ArrayRCP<ScalarT>>("Field Values", temperature); auto setFieldTemperature = Teuchos::rcp(new LCM::SetField<Residual, Traits>(setTempP)); fieldManager.registerEvaluator<Residual>(setFieldTemperature); stateFieldManager.registerEvaluator<Residual>(setFieldTemperature); } //--------------------------------------------------------------------------- // Time step Teuchos::ArrayRCP<ScalarT> delta_time(1); delta_time[0] = step_size; Teuchos::ParameterList setDTP("SetFieldTimeStep"); setDTP.set<std::string>("Evaluated Field Name", "Delta Time"); setDTP.set<Teuchos::RCP<PHX::DataLayout>>( "Evaluated Field Data Layout", dl->workset_scalar); setDTP.set<Teuchos::ArrayRCP<ScalarT>>("Field Values", delta_time); auto setFieldDT = Teuchos::rcp(new LCM::SetField<Residual, Traits>(setDTP)); fieldManager.registerEvaluator<Residual>(setFieldDT); stateFieldManager.registerEvaluator<Residual>(setFieldDT); // check if the material wants the tangent to be computed bool check_stability; check_stability = mpsParams.get<bool>("Check Stability", false); paramList.set<bool>("Compute Tangent", check_stability); std::cout << "Check stability = " << check_stability << std::endl; //--------------------------------------------------------------------------- // std::cout << "// Constitutive Model Parameters" //<< std::endl; Teuchos::ParameterList cmpPL; paramList.set<Teuchos::RCP<std::map<std::string, std::string>>>( "Name Map", fnm); cmpPL.set<Teuchos::ParameterList*>("Material Parameters", ¶mList); if (have_temperature) { cmpPL.set<std::string>("Temperature Name", "Temperature"); paramList.set<bool>("Have Temperature", true); } auto CMP = Teuchos::rcp( new LCM::ConstitutiveModelParameters<Residual, Traits>(cmpPL, dl)); fieldManager.registerEvaluator<Residual>(CMP); stateFieldManager.registerEvaluator<Residual>(CMP); //--------------------------------------------------------------------------- // std::cout << "// Constitutive Model Interface Evaluator" // << std::endl; Teuchos::ParameterList cmiPL; cmiPL.set<Teuchos::ParameterList*>("Material Parameters", ¶mList); if (have_temperature) { cmiPL.set<std::string>("Temperature Name", "Temperature"); } Teuchos::RCP<LCM::ConstitutiveModelInterface<Residual, Traits>> CMI = Teuchos::rcp( new LCM::ConstitutiveModelInterface<Residual, Traits>(cmiPL, dl)); fieldManager.registerEvaluator<Residual>(CMI); stateFieldManager.registerEvaluator<Residual>(CMI); // Set the evaluated fields as required for (std::vector<Teuchos::RCP<PHX::FieldTag>>::const_iterator it = CMI->evaluatedFields().begin(); it != CMI->evaluatedFields().end(); ++it) { fieldManager.requireField<Residual>(**it); } // register state variables Teuchos::RCP<Teuchos::ParameterList> p; Teuchos::RCP<PHX::Evaluator<Traits>> ev; for (int sv(0); sv < CMI->getNumStateVars(); ++sv) { CMI->fillStateVariableStruct(sv); p = stateMgr.registerStateVariable( CMI->getName(), CMI->getLayout(), dl->dummy, element_block_name, CMI->getInitType(), CMI->getInitValue(), CMI->getStateFlag(), CMI->getOutputFlag()); ev = Teuchos::rcp(new PHAL::SaveStateField<Residual, Traits>(*p)); fieldManager.registerEvaluator<Residual>(ev); stateFieldManager.registerEvaluator<Residual>(ev); } //--------------------------------------------------------------------------- if (check_stability) { std::string parametrization_type = mpsParams.get<std::string>("Parametrization Type", "Spherical"); double parametrization_interval = mpsParams.get<double>("Parametrization Interval", 0.05); std::cout << "Bifurcation Check in Material Point Simulator:" << std::endl; std::cout << "Parametrization Type: " << parametrization_type << std::endl; Teuchos::ParameterList bcPL; bcPL.set<Teuchos::ParameterList*>("Material Parameters", ¶mList); bcPL.set<std::string>("Parametrization Type Name", parametrization_type); bcPL.set<double>("Parametrization Interval Name", parametrization_interval); bcPL.set<std::string>("Material Tangent Name", "Material Tangent"); bcPL.set<std::string>("Ellipticity Flag Name", "Ellipticity_Flag"); bcPL.set<std::string>("Bifurcation Direction Name", "Direction"); bcPL.set<std::string>("Min detA Name", "Min detA"); Teuchos::RCP<LCM::BifurcationCheck<Residual, Traits>> BC = Teuchos::rcp(new LCM::BifurcationCheck<Residual, Traits>(bcPL, dl)); fieldManager.registerEvaluator<Residual>(BC); stateFieldManager.registerEvaluator<Residual>(BC); // register the ellipticity flag p = stateMgr.registerStateVariable( "Ellipticity_Flag", dl->qp_scalar, dl->dummy, element_block_name, "scalar", 0.0, false, true); ev = Teuchos::rcp(new PHAL::SaveStateField<Residual, Traits>(*p)); fieldManager.registerEvaluator<Residual>(ev); stateFieldManager.registerEvaluator<Residual>(ev); // register the direction p = stateMgr.registerStateVariable( "Direction", dl->qp_vector, dl->dummy, element_block_name, "scalar", 0.0, false, true); ev = Teuchos::rcp(new PHAL::SaveStateField<Residual, Traits>(*p)); fieldManager.registerEvaluator<Residual>(ev); stateFieldManager.registerEvaluator<Residual>(ev); // register min(det(A)) p = stateMgr.registerStateVariable( "Min detA", dl->qp_scalar, dl->dummy, element_block_name, "scalar", 0.0, false, true); ev = Teuchos::rcp(new PHAL::SaveStateField<Residual, Traits>(*p)); fieldManager.registerEvaluator<Residual>(ev); stateFieldManager.registerEvaluator<Residual>(ev); } //--------------------------------------------------------------------------- // std::cout << "// register deformation gradient" // << std::endl; p = stateMgr.registerStateVariable( "F", dl->qp_tensor, dl->dummy, element_block_name, "identity", 1.0, true, true); ev = Teuchos::rcp(new PHAL::SaveStateField<Residual, Traits>(*p)); fieldManager.registerEvaluator<Residual>(ev); stateFieldManager.registerEvaluator<Residual>(ev); //--------------------------------------------------------------------------- // std::cout << "// register small strain tensor" // << std::endl; p = stateMgr.registerStateVariable( "Strain", dl->qp_tensor, dl->dummy, element_block_name, "scalar", 0.0, false, true); ev = Teuchos::rcp(new PHAL::SaveStateField<Residual, Traits>(*p)); fieldManager.registerEvaluator<Residual>(ev); stateFieldManager.registerEvaluator<Residual>(ev); //--------------------------------------------------------------------------- // Traits::SetupData setupData = "Test String"; // std::cout << "Calling postRegistrationSetup" << std::endl; fieldManager.postRegistrationSetup(setupData); // std::cout << "// set the required fields for the state manager" //<< std::endl; Teuchos::RCP<PHX::DataLayout> dummy = Teuchos::rcp(new PHX::MDALayout<Dummy>(0)); std::vector<std::string> responseIDs = stateMgr.getResidResponseIDsToRequire(element_block_name); std::vector<std::string>::const_iterator it; for (it = responseIDs.begin(); it != responseIDs.end(); it++) { const std::string& responseID = *it; PHX::Tag<PHAL::AlbanyTraits::Residual::ScalarT> res_response_tag( responseID, dummy); stateFieldManager.requireField<PHAL::AlbanyTraits::Residual>( res_response_tag); } stateFieldManager.postRegistrationSetup(""); // std::cout << "Process using 'dot -Tpng -O <name>'\n"; fieldManager.writeGraphvizFile<Residual>("FM", true, true); stateFieldManager.writeGraphvizFile<Residual>("SFM", true, true); //--------------------------------------------------------------------------- // grab the output file name // std::string output_file = mpsParams.get<std::string>("Output File Name", "output.exo"); //--------------------------------------------------------------------------- // Create discretization, as required by the StateManager // Teuchos::RCP<Teuchos::ParameterList> discretizationParameterList = Teuchos::rcp(new Teuchos::ParameterList("Discretization")); discretizationParameterList->set<int>("1D Elements", workset_size); discretizationParameterList->set<int>("2D Elements", 1); discretizationParameterList->set<int>("3D Elements", 1); discretizationParameterList->set<std::string>("Method", "STK3D"); discretizationParameterList->set<int>("Number Of Time Derivatives", 0); discretizationParameterList->set<std::string>( "Exodus Output File Name", output_file); discretizationParameterList->set<int>("Workset Size", workset_size); Teuchos::RCP<Tpetra_Map> mapT = Teuchos::rcp(new Tpetra_Map( workset_size * num_dims * num_nodes, 0, commT, Tpetra::LocallyReplicated)); Teuchos::RCP<Tpetra_Vector> solution_vectorT = Teuchos::rcp(new Tpetra_Vector(mapT)); int numberOfEquations = 3; Albany::AbstractFieldContainer::FieldContainerRequirements req; Teuchos::RCP<Albany::AbstractSTKMeshStruct> stkMeshStruct = Teuchos::rcp(new Albany::TmplSTKMeshStruct<3>( discretizationParameterList, Teuchos::null, commT)); stkMeshStruct->setFieldAndBulkData( commT, discretizationParameterList, numberOfEquations, req, stateMgr.getStateInfoStruct(), stkMeshStruct->getMeshSpecs()[0]->worksetSize); Teuchos::RCP<Albany::AbstractDiscretization> discretization = Teuchos::rcp(new Albany::STKDiscretization( discretizationParameterList, stkMeshStruct, commT)); //--------------------------------------------------------------------------- // Associate the discretization with the StateManager // stateMgr.setupStateArrays(discretization); //--------------------------------------------------------------------------- // Create a workset // PHAL::Workset workset; workset.numCells = workset_size; workset.stateArrayPtr = &stateMgr.getStateArray(Albany::StateManager::ELEM, 0); // create MDFields PHX::MDField<ScalarT, Cell, QuadPoint, Dim, Dim> stressField( "Cauchy_Stress", dl->qp_tensor); // construct the final deformation gradient based on the loading case std::vector<ScalarT> F_vector(9, 0.0); if (load_case == "uniaxial") { F_vector[0] = 1.0 + number_steps * step_size; F_vector[4] = 1.0; F_vector[8] = 1.0; } else if (load_case == "simple-shear") { F_vector[0] = 1.0; F_vector[1] = number_steps * step_size; F_vector[4] = 1.0; F_vector[8] = 1.0; } else if (load_case == "hydrostatic") { F_vector[0] = 1.0 + number_steps * step_size; F_vector[4] = 1.0 + number_steps * step_size; F_vector[8] = 1.0 + number_steps * step_size; } else if (load_case == "general") { F_vector = mpsParams.get<Teuchos::Array<double>>("Deformation Gradient Components") .toVector(); } else { TEUCHOS_TEST_FOR_EXCEPTION( true, std::runtime_error, "Improper Loading Case in Material Point Simulator block"); } minitensor::Tensor<ScalarT> F_tensor(3, &F_vector[0]); minitensor::Tensor<ScalarT> log_F_tensor = minitensor::log(F_tensor); std::cout << "F\n" << F_tensor << std::endl; // std::cout << "log F\n" << log_F_tensor << std::endl; // // Setup loading scenario and instantiate evaluatFields // PHX::MDField<ScalarT, Cell, QuadPoint> minDetA("Min detA", dl->qp_scalar); PHX::MDField<ScalarT, Cell, QuadPoint, Dim> direction( "Direction", dl->qp_vector); // Bifurcation check parameters double mu_0 = 0; double mu_k = 0; int bifurcationTime_rough = number_steps; bool bifurcation_flag = false; for (int istep(0); istep <= number_steps; ++istep) { util::TimeGuard total_time_guard(total_time); // std::cout << "****** in MPS step " << istep << " ****** " << std::endl; // alpha \in [0,1] double alpha = double(istep) / number_steps; // std::cout << "alpha: " << alpha << std::endl; minitensor::Tensor<ScalarT> scaled_log_F_tensor = alpha * log_F_tensor; minitensor::Tensor<ScalarT> current_F = minitensor::exp(scaled_log_F_tensor); // std::cout << "scaled log F\n" << scaled_log_F_tensor << std::endl; // std::cout << "current F\n" << current_F << std::endl; for (int i = 0; i < 3; ++i) { for (int j = 0; j < 3; ++j) { def_grad[3 * i + j] = current_F(i, j); } } // jacobian detdefgrad[0] = minitensor::det(current_F); // small strain tensor minitensor::Tensor<ScalarT> current_strain; current_strain = 0.5 * (current_F + minitensor::transpose(current_F)) - minitensor::eye<ScalarT>(3); for (int i = 0; i < 3; ++i) { for (int j = 0; j < 3; ++j) { strain[3 * i + j] = current_strain(i, j); } } // std::cout << "current strain\n" << current_strain << std::endl; // Call the evaluators, evaluateFields() is the function that // computes stress based on deformation gradient compute_time->start(); fieldManager.preEvaluate<Residual>(workset); fieldManager.evaluateFields<Residual>(workset); fieldManager.postEvaluate<Residual>(workset); compute_time->stop(); stateFieldManager.getFieldData<Residual>(stressField); // Call the state field manager // std::cout << "+++ calling the stateFieldManager\n"; compute_time->start(); stateFieldManager.preEvaluate<Residual>(workset); stateFieldManager.evaluateFields<Residual>(workset); stateFieldManager.postEvaluate<Residual>(workset); compute_time->stop(); stateMgr.updateStates(); // output to the exodus file // Don't include this in timing data... total_time->stop(); discretization->writeSolutionT( *solution_vectorT, Teuchos::as<double>(istep)); // if check for bifurcation, adaptive step total_time->start(); if (check_stability) { // get current minDet(A) stateFieldManager.getFieldData<Residual>(minDetA); if (istep == 0) { mu_0 = minDetA(0, 0); } if (minDetA(0, 0) <= 0 && !bifurcation_flag) { mu_k = minDetA(0, 0); bifurcationTime_rough = istep; bifurcation_flag = true; // adaptive step begin std::cout << "\nAdaptive step begin - step " << istep << std::endl; // initialization for adaptive step double tol = 1E-8; double alpha_local = 1.0; double alpha_local_step = 0.5; int k = 1; int maxIteration = 50; // small strain tensor minitensor::Tensor<ScalarT> current_strain; // iteration begin while (((mu_k <= 0) || (std::abs(mu_k / mu_0) > tol))) { alpha = double(bifurcationTime_rough - 1 + alpha_local) / number_steps; minitensor::Tensor<ScalarT> scaled_log_F_tensor = alpha * log_F_tensor; minitensor::Tensor<ScalarT> current_F = minitensor::exp(scaled_log_F_tensor); for (int i = 0; i < 3; ++i) { for (int j = 0; j < 3; ++j) { def_grad[3 * i + j] = current_F(i, j); } } // jacobian detdefgrad[0] = minitensor::det(current_F); current_strain = 0.5 * (current_F + minitensor::transpose(current_F)) - minitensor::eye<ScalarT>(3); for (int i = 0; i < 3; ++i) { for (int j = 0; j < 3; ++j) { strain[3 * i + j] = current_strain(i, j); } } // Call the evaluators, evaluateFields() is the function that // computes stress based on deformation gradient fieldManager.preEvaluate<Residual>(workset); fieldManager.evaluateFields<Residual>(workset); fieldManager.postEvaluate<Residual>(workset); // Call the state field manager // std::cout << "+++ calling the stateFieldManager\n"; stateFieldManager.preEvaluate<Residual>(workset); stateFieldManager.evaluateFields<Residual>(workset); stateFieldManager.postEvaluate<Residual>(workset); stateFieldManager.getFieldData<Residual>(minDetA); stateFieldManager.getFieldData<Residual>(direction); mu_k = minDetA(0, 0); if (mu_k > 0) { alpha_local += alpha_local_step; } else { alpha_local -= alpha_local_step; } alpha_local_step /= 2; k = k + 1; if (k >= maxIteration) { std::cout << "Adaptive step for bifurcation check not converging after " << k << " iterations" << std::endl; break; } } // adaptive step iteration end } // end adaptive step } // end check bifurcation stateMgr.updateStates(); // if (bifurcation_flag) { // break the loading step after adaptive time step loop break; } // } // end loading steps // Summarize with AlbanyUtil performance monitors if (tout) { util::PerformanceContext::instance().timeMonitor().summarize(tout); tout.close(); } }
void LamentStress<EvalT, Traits>:: evaluateFields(typename Traits::EvalData workset) { Teuchos::RCP<lament::matParams<ScalarT>> matp = Teuchos::rcp(new lament::matParams<ScalarT>()); // Get the old state data Albany::MDArray oldDefGrad = (*workset.stateArrayPtr)[defGradName]; Albany::MDArray oldStress = (*workset.stateArrayPtr)[stressName]; int numStateVariables = (int)(this->lamentMaterialModelStateVariableNames.size()); // \todo Get actual time step for calls to LAMENT materials. double deltaT = 1.0; vector<ScalarT> strainRate(6); // symmetric tensor vector<ScalarT> spin(3); // skew-symmetric tensor vector<ScalarT> defGrad(9); // symmetric tensor vector<ScalarT> leftStretch(6); // symmetric tensor vector<ScalarT> rotation(9); // full tensor vector<double> stressOld(6); // symmetric tensor vector<ScalarT> stressNew(6); // symmetric tensor vector<double> stateOld(numStateVariables); // a single scalar for each state variable vector<double> stateNew(numStateVariables); // a single scalar for each state variable // \todo Set up scratch space for material models using getNumScratchVars() and setScratchPtr(). // Create the matParams structure, which is passed to Lament matp->nelements = 1; matp->dt = deltaT; matp->time = 0.0; matp->strain_rate = &strainRate[0]; matp->spin = &spin[0]; matp->deformation_gradient = &defGrad[0]; matp->left_stretch = &leftStretch[0]; matp->rotation = &rotation[0]; matp->state_old = &stateOld[0]; matp->state_new = &stateNew[0]; matp->stress_old = &stressOld[0]; matp->stress_new = &stressNew[0]; // matp->dt_mat = std::numeric_limits<double>::max(); // matParams that still need to be added: // matp->temp_old (temperature) // matp->temp_new // matp->sound_speed_old // matp->sound_speed_new // matp->volume // scratch pointer // function pointers (lots to be done here) for (int cell=0; cell < (int)workset.numCells; ++cell) { for (int qp=0; qp < (int)numQPs; ++qp) { // std::cout << "QP: " << qp << std::endl; // Fill the following entries in matParams for call to LAMENT // // nelements - number of elements // dt - time step, this one is tough because Albany does not currently have a concept of time step for implicit integration // time - current time, again Albany does not currently have a concept of time for implicit integration // strain_rate - what Sierra calls the rate of deformation, it is the symmetric part of the velocity gradient // spin - anti-symmetric part of the velocity gradient // left_stretch - found as V in the polar decomposition of the deformation gradient F = VR // rotation - found as R in the polar decomposition of the deformation gradient F = VR // state_old - material state data for previous time step (material dependent, none for lament::Elastic) // state_new - material state data for current time step (material dependent, none for lament::Elastic) // stress_old - stress at previous time step // stress_new - stress at current time step, filled by material model // // The total deformation gradient is available as field data // // The velocity gradient is not available but can be computed at the logarithm of the incremental deformation gradient divided by deltaT // The incremental deformation gradient is computed as F_new F_old^-1 // JTO: here is how I think this will go (of course the first two lines won't work as is...) // Intrepid2::Tensor<RealType> F = newDefGrad; // Intrepid2::Tensor<RealType> Fn = oldDefGrad; // Intrepid2::Tensor<RealType> f = F*Intrepid2::inverse(Fn); // Intrepid2::Tensor<RealType> V; // Intrepid2::Tensor<RealType> R; // boost::tie(V,R) = Intrepid2::polar_left(F); // Intrepid2::Tensor<RealType> Vinc; // Intrepid2::Tensor<RealType> Rinc; // Intrepid2::Tensor<RealType> logVinc; // boost::tie(Vinc,Rinc,logVinc) = Intrepid2::polar_left_logV(f) // Intrepid2::Tensor<RealType> logRinc = Intrepid2::log_rotation(Rinc); // Intrepid2::Tensor<RealType> logf = Intrepid2::bch(logVinc,logRinc); // Intrepid2::Tensor<RealType> L = (1.0/deltaT)*logf; // Intrepid2::Tensor<RealType> D = Intrepid2::sym(L); // Intrepid2::Tensor<RealType> W = Intrepid2::skew(L); // and then fill data into the vectors below // new deformation gradient (the current deformation gradient as computed in the current configuration) Intrepid2::Tensor<ScalarT> Fnew( 3, defGradField,cell,qp,0,0); // old deformation gradient (deformation gradient at previous load step) Intrepid2::Tensor<ScalarT> Fold( oldDefGrad(cell,qp,0,0), oldDefGrad(cell,qp,0,1), oldDefGrad(cell,qp,0,2), oldDefGrad(cell,qp,1,0), oldDefGrad(cell,qp,1,1), oldDefGrad(cell,qp,1,2), oldDefGrad(cell,qp,2,0), oldDefGrad(cell,qp,2,1), oldDefGrad(cell,qp,2,2) ); // incremental deformation gradient Intrepid2::Tensor<ScalarT> Finc = Fnew * Intrepid2::inverse(Fold); // DEBUGGING // //if(cell==0 && qp==0){ // std::cout << "Fnew(0,0) " << Fnew(0,0) << endl; // std::cout << "Fnew(1,0) " << Fnew(1,0) << endl; // std::cout << "Fnew(2,0) " << Fnew(2,0) << endl; // std::cout << "Fnew(0,1) " << Fnew(0,1) << endl; // std::cout << "Fnew(1,1) " << Fnew(1,1) << endl; // std::cout << "Fnew(2,1) " << Fnew(2,1) << endl; // std::cout << "Fnew(0,2) " << Fnew(0,2) << endl; // std::cout << "Fnew(1,2) " << Fnew(1,2) << endl; // std::cout << "Fnew(2,2) " << Fnew(2,2) << endl; //} // END DEBUGGING // // left stretch V, and rotation R, from left polar decomposition of new deformation gradient Intrepid2::Tensor<ScalarT> V(3), R(3), U(3); boost::tie(V,R) = Intrepid2::polar_left(Fnew); //V = R * U * transpose(R); // DEBUGGING // //if(cell==0 && qp==0){ // std::cout << "U(0,0) " << U(0,0) << endl; // std::cout << "U(1,0) " << U(1,0) << endl; // std::cout << "U(2,0) " << U(2,0) << endl; // std::cout << "U(0,1) " << U(0,1) << endl; // std::cout << "U(1,1) " << U(1,1) << endl; // std::cout << "U(2,1) " << U(2,1) << endl; // std::cout << "U(0,2) " << U(0,2) << endl; // std::cout << "U(1,2) " << U(1,2) << endl; // std::cout << "U(2,2) " << U(2,2) << endl; // std::cout << "========\n"; // std::cout << "V(0,0) " << V(0,0) << endl; // std::cout << "V(1,0) " << V(1,0) << endl; // std::cout << "V(2,0) " << V(2,0) << endl; // std::cout << "V(0,1) " << V(0,1) << endl; // std::cout << "V(1,1) " << V(1,1) << endl; // std::cout << "V(2,1) " << V(2,1) << endl; // std::cout << "V(0,2) " << V(0,2) << endl; // std::cout << "V(1,2) " << V(1,2) << endl; // std::cout << "V(2,2) " << V(2,2) << endl; // std::cout << "========\n"; // std::cout << "R(0,0) " << R(0,0) << endl; // std::cout << "R(1,0) " << R(1,0) << endl; // std::cout << "R(2,0) " << R(2,0) << endl; // std::cout << "R(0,1) " << R(0,1) << endl; // std::cout << "R(1,1) " << R(1,1) << endl; // std::cout << "R(2,1) " << R(2,1) << endl; // std::cout << "R(0,2) " << R(0,2) << endl; // std::cout << "R(1,2) " << R(1,2) << endl; // std::cout << "R(2,2) " << R(2,2) << endl; //} // END DEBUGGING // // incremental left stretch Vinc, incremental rotation Rinc, and log of incremental left stretch, logVinc Intrepid2::Tensor<ScalarT> Uinc(3), Vinc(3), Rinc(3), logVinc(3); //boost::tie(Vinc,Rinc,logVinc) = Intrepid2::polar_left_logV(Finc); boost::tie(Vinc,Rinc) = Intrepid2::polar_left(Finc); //Vinc = Rinc * Uinc * transpose(Rinc); logVinc = Intrepid2::log(Vinc); // log of incremental rotation Intrepid2::Tensor<ScalarT> logRinc = Intrepid2::log_rotation(Rinc); // log of incremental deformation gradient Intrepid2::Tensor<ScalarT> logFinc = Intrepid2::bch(logVinc, logRinc); // velocity gradient Intrepid2::Tensor<ScalarT> L = (1.0/deltaT)*logFinc; // strain rate (a.k.a rate of deformation) Intrepid2::Tensor<ScalarT> D = Intrepid2::sym(L); // spin Intrepid2::Tensor<ScalarT> W = Intrepid2::skew(L); // load everything into the Lament data structure strainRate[0] = ( D(0,0) ); strainRate[1] = ( D(1,1) ); strainRate[2] = ( D(2,2) ); strainRate[3] = ( D(0,1) ); strainRate[4] = ( D(1,2) ); strainRate[5] = ( D(2,0) ); spin[0] = ( W(0,1) ); spin[1] = ( W(1,2) ); spin[2] = ( W(2,0) ); leftStretch[0] = ( V(0,0) ); leftStretch[1] = ( V(1,1) ); leftStretch[2] = ( V(2,2) ); leftStretch[3] = ( V(0,1) ); leftStretch[4] = ( V(1,2) ); leftStretch[5] = ( V(2,0) ); rotation[0] = ( R(0,0) ); rotation[1] = ( R(1,1) ); rotation[2] = ( R(2,2) ); rotation[3] = ( R(0,1) ); rotation[4] = ( R(1,2) ); rotation[5] = ( R(2,0) ); rotation[6] = ( R(1,0) ); rotation[7] = ( R(2,1) ); rotation[8] = ( R(0,2) ); defGrad[0] = ( Fnew(0,0) ); defGrad[1] = ( Fnew(1,1) ); defGrad[2] = ( Fnew(2,2) ); defGrad[3] = ( Fnew(0,1) ); defGrad[4] = ( Fnew(1,2) ); defGrad[5] = ( Fnew(2,0) ); defGrad[6] = ( Fnew(1,0) ); defGrad[7] = ( Fnew(2,1) ); defGrad[8] = ( Fnew(0,2) ); stressOld[0] = oldStress(cell,qp,0,0); stressOld[1] = oldStress(cell,qp,1,1); stressOld[2] = oldStress(cell,qp,2,2); stressOld[3] = oldStress(cell,qp,0,1); stressOld[4] = oldStress(cell,qp,1,2); stressOld[5] = oldStress(cell,qp,2,0); // copy data from the state manager to the LAMENT data structure for(int iVar=0 ; iVar<numStateVariables ; iVar++){ const std::string& variableName = this->lamentMaterialModelStateVariableNames[iVar]+"_old"; Albany::MDArray stateVar = (*workset.stateArrayPtr)[variableName]; stateOld[iVar] = stateVar(cell,qp); } // Make a call to the LAMENT material model to initialize the load step this->lamentMaterialModel->loadStepInit(matp.get()); // Get the stress from the LAMENT material // std::cout << "about to call lament->getStress()" << std::endl; this->lamentMaterialModel->getStress(matp.get()); // std::cout << "after calling lament->getStress() 2" << std::endl; // rotate to get the Cauchy Stress Intrepid2::Tensor<ScalarT> lameStress( stressNew[0], stressNew[3], stressNew[5], stressNew[3], stressNew[1], stressNew[4], stressNew[5], stressNew[4], stressNew[2] ); Intrepid2::Tensor<ScalarT> cauchy = R * lameStress * transpose(R); // DEBUGGING // //if(cell==0 && qp==0){ // std::cout << "check strainRate[0] " << strainRate[0] << endl; // std::cout << "check strainRate[1] " << strainRate[1] << endl; // std::cout << "check strainRate[2] " << strainRate[2] << endl; // std::cout << "check strainRate[3] " << strainRate[3] << endl; // std::cout << "check strainRate[4] " << strainRate[4] << endl; // std::cout << "check strainRate[5] " << strainRate[5] << endl; //} // END DEBUGGING // // Copy the new stress into the stress field for (int i(0); i < 3; ++i) for (int j(0); j < 3; ++j) stressField(cell,qp,i,j) = cauchy(i,j); // stressField(cell,qp,0,0) = stressNew[0]; // stressField(cell,qp,1,1) = stressNew[1]; // stressField(cell,qp,2,2) = stressNew[2]; // stressField(cell,qp,0,1) = stressNew[3]; // stressField(cell,qp,1,2) = stressNew[4]; // stressField(cell,qp,2,0) = stressNew[5]; // stressField(cell,qp,1,0) = stressNew[3]; // stressField(cell,qp,2,1) = stressNew[4]; // stressField(cell,qp,0,2) = stressNew[5]; // copy state_new data from the LAMENT data structure to the corresponding state variable field for(int iVar=0 ; iVar<numStateVariables ; iVar++) this->lamentMaterialModelStateVariableFields[iVar](cell,qp) = stateNew[iVar]; // DEBUGGING // //if(cell==0 && qp==0){ // std::cout << "stress(0,0) " << this->stressField(cell,qp,0,0) << endl; // std::cout << "stress(1,1) " << this->stressField(cell,qp,1,1) << endl; // std::cout << "stress(2,2) " << this->stressField(cell,qp,2,2) << endl; // std::cout << "stress(0,1) " << this->stressField(cell,qp,0,1) << endl; // std::cout << "stress(1,2) " << this->stressField(cell,qp,1,2) << endl; // std::cout << "stress(0,2) " << this->stressField(cell,qp,0,2) << endl; // std::cout << "stress(1,0) " << this->stressField(cell,qp,1,0) << endl; // std::cout << "stress(2,1) " << this->stressField(cell,qp,2,1) << endl; // std::cout << "stress(2,0) " << this->stressField(cell,qp,2,0) << endl; // //} // // END DEBUGGING // } } }