void PeridigmNS::ElasticPlasticMaterial::computeAutomaticDifferentiationJacobian(const double dt, const int numOwnedPoints, const int* ownedIDs, const int* neighborhoodList, PeridigmNS::DataManager& dataManager, PeridigmNS::SerialMatrix& jacobian, PeridigmNS::Material::JacobianType jacobianType) const { // Compute contributions to the tangent matrix on an element-by-element basis // To reduce memory re-allocation, use static variable to store Fad types for // current coordinates (independent variables). static vector<Sacado::Fad::DFad<double> > y_AD; // Loop over all points. int neighborhoodListIndex = 0; for(int iID=0 ; iID<numOwnedPoints ; ++iID){ // Create a temporary neighborhood consisting of a single point and its neighbors. int numNeighbors = neighborhoodList[neighborhoodListIndex++]; int numEntries = numNeighbors+1; int numDof = 3*numEntries; vector<int> tempMyGlobalIDs(numEntries); // Put the node at the center of the neighborhood at the beginning of the list. tempMyGlobalIDs[0] = dataManager.getOwnedScalarPointMap()->GID(iID); vector<int> tempNeighborhoodList(numEntries); tempNeighborhoodList[0] = numNeighbors; for(int iNID=0 ; iNID<numNeighbors ; ++iNID){ int neighborID = neighborhoodList[neighborhoodListIndex++]; tempMyGlobalIDs[iNID+1] = dataManager.getOverlapScalarPointMap()->GID(neighborID); tempNeighborhoodList[iNID+1] = iNID+1; } Epetra_SerialComm serialComm; Teuchos::RCP<Epetra_BlockMap> tempOneDimensionalMap = Teuchos::rcp(new Epetra_BlockMap(numEntries, numEntries, &tempMyGlobalIDs[0], 1, 0, serialComm)); Teuchos::RCP<Epetra_BlockMap> tempThreeDimensionalMap = Teuchos::rcp(new Epetra_BlockMap(numEntries, numEntries, &tempMyGlobalIDs[0], 3, 0, serialComm)); Teuchos::RCP<Epetra_BlockMap> tempBondMap = Teuchos::rcp(new Epetra_BlockMap(1, 1, &tempMyGlobalIDs[0], numNeighbors, 0, serialComm)); // Create a temporary DataManager containing data for this point and its neighborhood. PeridigmNS::DataManager tempDataManager; tempDataManager.setMaps(Teuchos::RCP<const Epetra_BlockMap>(), tempOneDimensionalMap, Teuchos::RCP<const Epetra_BlockMap>(), tempThreeDimensionalMap, tempBondMap); // The temporary data manager will have the same fields and data as the real data manager. vector<int> fieldIds = dataManager.getFieldIds(); tempDataManager.allocateData(fieldIds); tempDataManager.copyLocallyOwnedDataFromDataManager(dataManager); // Set up numOwnedPoints and ownedIDs. // There is only one owned ID, and it has local ID zero in the tempDataManager. int tempNumOwnedPoints = 1; vector<int> tempOwnedIDs(tempNumOwnedPoints); tempOwnedIDs[0] = 0; // Use the scratchMatrix as sub-matrix for storing tangent values prior to loading them into the global tangent matrix. // Resize scratchMatrix if necessary if(scratchMatrix.Dimension() < numDof) scratchMatrix.Resize(numDof); // Create a list of global indices for the rows/columns in the scratch matrix. vector<int> globalIndices(numDof); for(int i=0 ; i<numEntries ; ++i){ int globalID = tempOneDimensionalMap->GID(i); for(int j=0 ; j<3 ; ++j) globalIndices[3*i+j] = 3*globalID+j; } // Extract pointers to the underlying data in the constitutiveData array. double *x, *y, *cellVolume, *weightedVolume, *damage, *bondDamage, *edpN, *lambdaN; tempDataManager.getData(m_modelCoordinatesFieldId, PeridigmField::STEP_NONE)->ExtractView(&x); tempDataManager.getData(m_coordinatesFieldId, PeridigmField::STEP_NP1)->ExtractView(&y); tempDataManager.getData(m_volumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&cellVolume); tempDataManager.getData(m_weightedVolumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&weightedVolume); tempDataManager.getData(m_damageFieldId, PeridigmField::STEP_NP1)->ExtractView(&damage); tempDataManager.getData(m_bondDamageFieldId, PeridigmField::STEP_NP1)->ExtractView(&bondDamage); tempDataManager.getData(m_deviatoricPlasticExtensionFieldId, PeridigmField::STEP_N)->ExtractView(&edpN); tempDataManager.getData(m_lambdaFieldId, PeridigmField::STEP_N)->ExtractView(&lambdaN); // Create arrays of Fad objects for the current coordinates, dilatation, and force density // Modify the existing vector of Fad objects for the current coordinates if((int)y_AD.size() < numDof) y_AD.resize(numDof); for(int i=0 ; i<numDof ; ++i){ y_AD[i].diff(i, numDof); y_AD[i].val() = y[i]; } // Create vectors of empty AD types for the dependent variables vector<Sacado::Fad::DFad<double> > dilatation_AD(numEntries); vector<Sacado::Fad::DFad<double> > lambdaNP1_AD(numEntries); int numBonds = tempDataManager.getData(m_deviatoricPlasticExtensionFieldId, PeridigmField::STEP_N)->MyLength(); vector<Sacado::Fad::DFad<double> > edpNP1(numBonds); vector<Sacado::Fad::DFad<double> > force_AD(numDof); // Evaluate the constitutive model using the AD types MATERIAL_EVALUATION::computeDilatation(x,&y_AD[0],weightedVolume,cellVolume,bondDamage,&dilatation_AD[0],&tempNeighborhoodList[0],tempNumOwnedPoints,m_horizon); MATERIAL_EVALUATION::computeInternalForceIsotropicElasticPlastic ( x, &y_AD[0], weightedVolume, cellVolume, &dilatation_AD[0], bondDamage, edpN, &edpNP1[0], lambdaN, &lambdaNP1_AD[0], &force_AD[0], &tempNeighborhoodList[0], tempNumOwnedPoints, m_bulkModulus, m_shearModulus, m_horizon, m_yieldStress, m_isPlanarProblem, m_thickness); // Load derivative values into scratch matrix // Multiply by volume along the way to convert force density to force for(int row=0 ; row<numDof ; ++row){ for(int col=0 ; col<numDof ; ++col){ scratchMatrix(row, col) = force_AD[row].dx(col) * cellVolume[row/3]; } } // Sum the values into the global tangent matrix (this is expensive). jacobian.addValues((int)globalIndices.size(), &globalIndices[0], scratchMatrix.Data()); } }