void
PeridigmNS::InterfaceAwareDamageModel::initialize(const double dt,
                                                   const int numOwnedPoints,
                                                   const int* ownedIDs,
                                                   const int* neighborhoodList,
                                                   PeridigmNS::DataManager& dataManager) const
{
  double *damage, *bondDamage, *criticalStretch;
  dataManager.getData(m_damageFieldId, PeridigmField::STEP_NP1)->ExtractView(&damage);
  dataManager.getData(m_bondDamageFieldId, PeridigmField::STEP_NP1)->ExtractView(&bondDamage);
  dataManager.getData(m_criticalStretchFieldId, PeridigmField::STEP_NONE)->ExtractView(&criticalStretch);

//  int totalListSize = numOwnedPoints;

  // Initialize damage to zero
  int neighborhoodListIndex = 0;
  int bondIndex = 0;
  for(int iID=0 ; iID<numOwnedPoints ; ++iID){
    int nodeID = ownedIDs[iID];
    damage[nodeID] = 0.0;
    criticalStretch[nodeID] = m_criticalStretch;
    int numNeighbors = neighborhoodList[neighborhoodListIndex++];
    //	  totalListSize += numNeighbors;
    neighborhoodListIndex += numNeighbors;
    for(int iNID=0 ; iNID<numNeighbors ; ++iNID){
      bondDamage[bondIndex++] = 0.0;
    }
  }
}
void PeridigmNS::ElasticPlasticMaterial::initialize(const double dt,
                                                    const int numOwnedPoints,
                                                    const int* ownedIDs,
                                                    const int* neighborhoodList,
                                                    PeridigmNS::DataManager& dataManager)
{
  // Extract pointers to the underlying data
  double *xOverlap,  *cellVolumeOverlap, *weightedVolume;
  dataManager.getData(m_modelCoordinatesFieldId, PeridigmField::STEP_NONE)->ExtractView(&xOverlap);
  dataManager.getData(m_volumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&cellVolumeOverlap);
  dataManager.getData(m_weightedVolumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&weightedVolume);

  MATERIAL_EVALUATION::computeWeightedVolume(xOverlap,cellVolumeOverlap,weightedVolume,numOwnedPoints,neighborhoodList,m_horizon);
}
예제 #3
0
void
PeridigmNS::LCMMaterial::computeForce(const double dt,
                                        const int numOwnedPoints,
                                        const int* ownedIDs,
                                        const int* neighborhoodList,
                                        PeridigmNS::DataManager& dataManager) const
{
  // Zero out the forces
  dataManager.getData(m_forceDensityFieldId, PeridigmField::STEP_NP1)->PutScalar(0.0);
}
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());
  }
}
void
PeridigmNS::ElasticPlasticMaterial::computeForce(const double dt,
                                                 const int numOwnedPoints,
                                                 const int* ownedIDs,
                                                 const int* neighborhoodList,
                                                 PeridigmNS::DataManager& dataManager) const
{
  double *x, *y, *volume, *dilatation, *weightedVolume, *bondDamage, *edpN, *edpNP1, *lambdaN, *lambdaNP1, *force;
  dataManager.getData(m_modelCoordinatesFieldId, PeridigmField::STEP_NONE)->ExtractView(&x);
  dataManager.getData(m_coordinatesFieldId, PeridigmField::STEP_NP1)->ExtractView(&y);
  dataManager.getData(m_volumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&volume);
  dataManager.getData(m_dilatationFieldId, PeridigmField::STEP_NP1)->ExtractView(&dilatation);
  dataManager.getData(m_weightedVolumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&weightedVolume);
  dataManager.getData(m_bondDamageFieldId, PeridigmField::STEP_NP1)->ExtractView(&bondDamage);
  dataManager.getData(m_deviatoricPlasticExtensionFieldId, PeridigmField::STEP_N)->ExtractView(&edpN);
  dataManager.getData(m_deviatoricPlasticExtensionFieldId, PeridigmField::STEP_NP1)->ExtractView(&edpNP1);
  dataManager.getData(m_lambdaFieldId, PeridigmField::STEP_N)->ExtractView(&lambdaN);
  dataManager.getData(m_lambdaFieldId, PeridigmField::STEP_NP1)->ExtractView(&lambdaNP1);
  dataManager.getData(m_forceDensityFieldId, PeridigmField::STEP_NP1)->ExtractView(&force);

  // Zero out the force
  dataManager.getData(m_forceDensityFieldId, PeridigmField::STEP_NP1)->PutScalar(0.0);

  MATERIAL_EVALUATION::computeDilatation(x,y,weightedVolume,volume,bondDamage,dilatation,neighborhoodList,numOwnedPoints,m_horizon);
  MATERIAL_EVALUATION::computeInternalForceIsotropicElasticPlastic
     (
       x,
       y,
       weightedVolume,
       volume,
       dilatation,
       bondDamage,
       edpN,
       edpNP1,
       lambdaN,
       lambdaNP1,
       force,
       neighborhoodList,
       numOwnedPoints,
       m_bulkModulus,
       m_shearModulus,
       m_horizon,
       m_yieldStress,
       m_isPlanarProblem,
       m_thickness
    );
}
예제 #6
0
void
PeridigmNS::Pals_Model::computeStoredElasticEnergyDensity(const double dt,
                                                               const int numOwnedPoints,
                                                               const int* ownedIDs,
                                                               const int* neighborhoodList,
                                                               PeridigmNS::DataManager& dataManager) const
{
	// namespace PALS
	using namespace MATERIAL_EVALUATION::PALS;

	double K=m_bulkModulus;
	double MU=m_shearModulus;
	double horizon=m_horizon;

  // This function is intended to be called from a compute class.
  // The compute class should have already created the Stored_Elastic_Energy_Density field id.
  int storedElasticEnergyDensityFieldId = PeridigmNS::FieldManager::self().getFieldId("Stored_Elastic_Energy_Density");

	// Extract pointers to the underlying data
	double *xOverlap, *yOverlap, *volumeOverlap;
	double *omega_constants, *sigma_constants;
	double *dilatation, *energy_density;
	vector<const double *> omega_multipliers(num_lagrange_multipliers), sigma_multipliers(num_lagrange_multipliers);
	dataManager.getData(m_modelCoordinatesFieldId, PeridigmField::STEP_NONE)->ExtractView(&xOverlap);
	dataManager.getData(m_coordinatesFieldId, PeridigmField::STEP_NP1)->ExtractView(&yOverlap);
	dataManager.getData(m_volumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&volumeOverlap);
	dataManager.getData(m_dilatationFieldId, PeridigmField::STEP_NP1)->ExtractView(&dilatation);
	dataManager.getData(m_dilatationNormalizationFieldId, PeridigmField::STEP_NONE)->ExtractView(&omega_constants);
	dataManager.getData(m_deviatoricNormalizationFieldId, PeridigmField::STEP_NONE)->ExtractView(&sigma_constants);
  dataManager.getData(storedElasticEnergyDensityFieldId , PeridigmField::STEP_NONE)->ExtractView(&energy_density);
	for(int i=0;i<num_lagrange_multipliers;i++){
		double *dil, *dev;
		dataManager.getData(m_dilatationLagrangeMultiplersFieldIds[i],PeridigmField::STEP_NONE)->ExtractView(&dil);
		dataManager.getData(m_deviatoricLagrangeMultiplersFieldIds[i],PeridigmField::STEP_NONE)->ExtractView(&dev);
		omega_multipliers[i]=dil;
		sigma_multipliers[i]=dev;
	}

	double bond[3];
	const double *xOwned = xOverlap;
	const double *yOwned = yOverlap;
	const double *oc=omega_constants;
	const double *sc=sigma_constants;
	const double *theta=dilatation;
	const int *neighPtr = neighborhoodList;
	double cell_volume;
	for(int q=0; q<numOwnedPoints; q++, xOwned+=3, yOwned+=3, oc++, sc++, theta++){

		// Collect computed Lagrange multipliers for this point 'q'
		double sigma_X[NUM_LAGRANGE_MULTIPLIERS];
		for(int i=0;i<NUM_LAGRANGE_MULTIPLIERS;i++){
			sigma_X[i]=sigma_multipliers[i][q];
		}

		int numNeigh = *neighPtr; neighPtr++;
		const double *X = xOwned;
		const double *Y = yOwned;
		double sum=0.0;
		for(int n=0;n<numNeigh;n++,neighPtr++){
			int localId = *neighPtr;
			cell_volume = volumeOverlap[localId];
			const double *XP = &xOverlap[3*localId];
			const double *YP = &yOverlap[3*localId];
			bond[0]=XP[0]-X[0];
			bond[1]=XP[1]-X[1];
			bond[2]=XP[2]-X[2];
			double a = bond[0];
			double b = bond[1];
			double c = bond[2];
			double xi = sqrt(a*a+b*b+c*c);
			a = YP[0]-Y[0];
			b = YP[1]-Y[1];
			c = YP[2]-Y[2];
			double dY = sqrt(a*a+b*b+c*c);
			double e = dY-xi;
      double eps=e-xi*(*theta)/3.0;
			pals_influence<deviatoric_influence> SIGMA(m_SIGMA_0,*sc,sigma_X);
      double sigma = SIGMA(bond,horizon);
      sum += sigma * eps * eps * cell_volume;
		}
		double d=*theta;
		/*
		 * Deviatoric term does not include factor of 1/2
		 */
		energy_density[q] = 0.5 * K * d * d + MU * sum;
	}
}
예제 #7
0
void
PeridigmNS::Pals_Model::computeForce(const double dt,
                                          const int numOwnedPoints,
                                          const int* ownedIDs,
                                          const int* neighborhoodList,
                                          PeridigmNS::DataManager& dataManager) const
{
	// Zero out the forces
	dataManager.getData(m_forceDensityFieldId, PeridigmField::STEP_NP1)->PutScalar(0.0);

	// Extract pointers to the underlying data
	double *x, *y, *cellVolume, *pals_pressure, *force, *weightedVolume;
	double *omega_constants, *sigma_constants;
	double *dilatation;
	vector<const double *> omega_multipliers(num_lagrange_multipliers), sigma_multipliers(num_lagrange_multipliers);
	dataManager.getData(m_modelCoordinatesFieldId, PeridigmField::STEP_NONE)->ExtractView(&x);
	dataManager.getData(m_coordinatesFieldId, PeridigmField::STEP_NP1)->ExtractView(&y);
	dataManager.getData(m_volumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&cellVolume);
	dataManager.getData(m_weightedVolumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&weightedVolume);
	dataManager.getData(m_dilatationFieldId, PeridigmField::STEP_NP1)->ExtractView(&dilatation);
	dataManager.getData(m_palsPressureFieldId, PeridigmField::STEP_NP1)->ExtractView(&pals_pressure);
	dataManager.getData(m_forceDensityFieldId, PeridigmField::STEP_NP1)->ExtractView(&force);
	dataManager.getData(m_dilatationNormalizationFieldId, PeridigmField::STEP_NONE)->ExtractView(&omega_constants);
	dataManager.getData(m_deviatoricNormalizationFieldId, PeridigmField::STEP_NONE)->ExtractView(&sigma_constants);

	for(int i=0;i<num_lagrange_multipliers;i++){
		double *dil, *dev;
		dataManager.getData(m_dilatationLagrangeMultiplersFieldIds[i],PeridigmField::STEP_NONE)->ExtractView(&dil);
		dataManager.getData(m_deviatoricLagrangeMultiplersFieldIds[i],PeridigmField::STEP_NONE)->ExtractView(&dev);
		omega_multipliers[i]=dil;
		sigma_multipliers[i]=dev;
	}

	// namespace PALS
	using namespace MATERIAL_EVALUATION::PALS;
	computeDilatationAndPalsPressure
	(
		x,
		y,
		cellVolume,
		omega_multipliers,
		omega_constants,
		sigma_multipliers,
		sigma_constants,
		weightedVolume,
		dilatation,
		pals_pressure,
		neighborhoodList,
		numOwnedPoints,
		m_bulkModulus,
		m_shearModulus,
		m_horizon,
		m_OMEGA_0,
		m_SIGMA_0
	);

	computeInternalForcePals
	(
		x,
		y,
		cellVolume,
		omega_multipliers,
		omega_constants,
		sigma_multipliers,
		sigma_constants,
		dilatation,
		pals_pressure,
		force,
		neighborhoodList,
		numOwnedPoints,
		m_bulkModulus,
		m_shearModulus,
		m_horizon,
		m_OMEGA_0,
		m_SIGMA_0
	);


}
예제 #8
0
void
PeridigmNS::Pals_Model::
initialize(const double dt,
           const int numOwnedPoints,
           const int* ownedIDs,
           const int* neighborhoodList,
           PeridigmNS::DataManager& dataManager)
{

	// Extract pointers to the underlying data
	double *xOverlap,  *cellVolumeOverlap, *weightedVolume;
	double *omega_constants, *sigma_constants;
	double *normalized_weighted_volume;
	dataManager.getData(m_modelCoordinatesFieldId, PeridigmField::STEP_NONE)->ExtractView(&xOverlap);
	dataManager.getData(m_volumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&cellVolumeOverlap);
	dataManager.getData(m_weightedVolumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&weightedVolume);
	dataManager.getData(m_normalizedWeightedVolumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&normalized_weighted_volume);
	dataManager.getData(m_dilatationNormalizationFieldId, PeridigmField::STEP_NONE)->ExtractView(&omega_constants);
	dataManager.getData(m_deviatoricNormalizationFieldId, PeridigmField::STEP_NONE)->ExtractView(&sigma_constants);

	// namespace PALS
	using namespace MATERIAL_EVALUATION::PALS;
	{
		vector<double *> omega_multipliers(num_lagrange_multipliers), sigma_multipliers(num_lagrange_multipliers);
		for(int i=0;i<num_lagrange_multipliers;i++){
			double *dil, *dev;
			dataManager.getData(m_dilatationLagrangeMultiplersFieldIds[i],PeridigmField::STEP_NONE)->ExtractView(&dil);
			dataManager.getData(m_deviatoricLagrangeMultiplersFieldIds[i],PeridigmField::STEP_NONE)->ExtractView(&dev);
			omega_multipliers[i]=dil;
			sigma_multipliers[i]=dev;
		}

		compute_lagrange_multipliers
		(
			xOverlap,
			cellVolumeOverlap,
			numOwnedPoints,
			neighborhoodList,
			m_horizon,
			omega_multipliers,
			omega_constants,
			sigma_multipliers,
			sigma_constants,
			m_OMEGA_0,
			m_SIGMA_0
		);
	}


	{
		vector<const double *> sigma_multipliers(num_lagrange_multipliers);
		for(int i=0;i<num_lagrange_multipliers;i++){
			double *dev;
			dataManager.getData(m_deviatoricLagrangeMultiplersFieldIds[i],PeridigmField::STEP_NONE)->ExtractView(&dev);
			sigma_multipliers[i]=dev;
		}
		/*
		 * DEBUGGING
		 * Compute normalized weighted volume: all values should be 3.0
		computeNormalizedWeightedVolume
		(
			xOverlap,
			cellVolumeOverlap,
			omega_constants,
			bondDamage,
			normalized_weighted_volume,
			numOwnedPoints,
			neighborhoodList,
			m_horizon,
			m_SIGMA_0
		);
		*/
		/*
		 * Weighted volume with influence function $\sigma$ rather than $\omega$
		 */
		computeWeightedVolume
		(
			xOverlap,
			cellVolumeOverlap,
			sigma_multipliers,
			sigma_constants,
			weightedVolume,
			numOwnedPoints,
			neighborhoodList,
			m_horizon,
			m_SIGMA_0
		);
	}
}
void
PeridigmNS::InterfaceAwareDamageModel::computeDamage(const double dt,
                                                      const int numOwnedPoints,
                                                      const int* ownedIDs,
                                                      const int* neighborhoodList,
                                                      PeridigmNS::DataManager& dataManager) const
{
  double *x, *y, *damage, *bondDamage, *deltaTemperature, *criticalStretch;
  dataManager.getData(m_modelCoordinatesFieldId, PeridigmField::STEP_NONE)->ExtractView(&x);
  dataManager.getData(m_coordinatesFieldId, PeridigmField::STEP_NP1)->ExtractView(&y);
  dataManager.getData(m_damageFieldId, PeridigmField::STEP_NP1)->ExtractView(&damage);
  dataManager.getData(m_bondDamageFieldId, PeridigmField::STEP_NP1)->ExtractView(&bondDamage);
  dataManager.getData(m_criticalStretchFieldId, PeridigmField::STEP_NONE)->ExtractView(&criticalStretch);
  deltaTemperature = NULL;
  if(m_applyThermalStrains)
    dataManager.getData(m_deltaTemperatureFieldId, PeridigmField::STEP_NP1)->ExtractView(&deltaTemperature);

  double trialDamage(0.0);
  int neighborhoodListIndex(0), bondIndex(0);
  int nodeId, numNeighbors, neighborID, iID, iNID;
  double nodeInitialX[3], nodeCurrentX[3], initialDistance, currentDistance, relativeExtension, totalDamage;

  // Update the bond damage
  // Break bonds if the extension is greater than the critical extension
  for(iID=0 ; iID<numOwnedPoints ; ++iID){
    nodeId = ownedIDs[iID];
    nodeInitialX[0] = x[nodeId*3];
    nodeInitialX[1] = x[nodeId*3+1];
    nodeInitialX[2] = x[nodeId*3+2];
    nodeCurrentX[0] = y[nodeId*3];
    nodeCurrentX[1] = y[nodeId*3+1];
    nodeCurrentX[2] = y[nodeId*3+2];
    numNeighbors = neighborhoodList[neighborhoodListIndex++];

    for(iNID=0 ; iNID<numNeighbors ; ++iNID){
      neighborID = neighborhoodList[neighborhoodListIndex++];
      const double neighborCriticalStretch = criticalStretch[neighborID];
      // Skip neighbors that do not allow bond breaking or bonds that are already broken
      if(bondDamage[bondIndex]!=0.0 || neighborCriticalStretch == 0.0){bondIndex+=1;continue;}

      TEUCHOS_TEST_FOR_EXCEPTION(neighborCriticalStretch==0.0,std::logic_error,"Error: Neighbor critical stretch cannot be zero at this point.");
      const double minCriticalStretch = (neighborCriticalStretch < m_criticalStretch) ? neighborCriticalStretch : m_criticalStretch;

      initialDistance = 
          distance(nodeInitialX[0], nodeInitialX[1], nodeInitialX[2],
            x[neighborID*3], x[neighborID*3+1], x[neighborID*3+2]);
      currentDistance = 
          distance(nodeCurrentX[0], nodeCurrentX[1], nodeCurrentX[2],
            y[neighborID*3], y[neighborID*3+1], y[neighborID*3+2]);
      if(m_applyThermalStrains)
        currentDistance -= m_alpha*deltaTemperature[nodeId]*initialDistance;
      relativeExtension = (currentDistance - initialDistance)/initialDistance;
      trialDamage = 0.0;
      if(relativeExtension > minCriticalStretch)
      {
        trialDamage = 1.0;
      }
      if(trialDamage > bondDamage[bondIndex]){
        bondDamage[bondIndex] = trialDamage;
      }
      bondIndex += 1;
    }
  }

  TEUCHOS_TEST_FOR_EXCEPTION(m_bcManager==Teuchos::null,std::logic_error,"Error: the bc manager pointer should have been set by here.");
  Teuchos::RCP< std::map< std::string, std::vector<int> > > nodeSetMap = m_bcManager->getNodeSets();

  //  Update the element damage (percent of bonds broken)
  neighborhoodListIndex = 0;
  bondIndex = 0;
  for(iID=0 ; iID<numOwnedPoints ; ++iID){
	nodeId = ownedIDs[iID];
	numNeighbors = neighborhoodList[neighborhoodListIndex++];
    neighborhoodListIndex += numNeighbors;
	totalDamage = 0.0;
	for(iNID=0 ; iNID<numNeighbors ; ++iNID){
	  totalDamage += bondDamage[bondIndex++];
	}

	if(totalDamage >= numNeighbors-2) // This would imply rank deficiency and would lead to problems in CG
	{
	  const string deficientSetName = "RANK_DEFICIENT_NODES";
	  TEUCHOS_TEST_FOR_EXCEPTION(nodeSetMap->find(deficientSetName)==nodeSetMap->end(),std::logic_error,"Error: The placeholder nodeset for rank deficient nodes is missing.");
	  bool nodeAlreadyRegistered = false;
	  vector<int> * deficientSet = &nodeSetMap->find(deficientSetName)->second;
	  for(unsigned i=0;i<deficientSet->size();++i)
		if((*deficientSet)[i]==nodeId)
		  nodeAlreadyRegistered = true;

	  if(!nodeAlreadyRegistered){
		cout << "Warning: Potentially rank deficient node detected (Node with 2 or less bonds). " << endl
			 << "Node " << nodeId + 1 << " will be removed from the linear system." << endl;
		deficientSet->push_back(nodeId);

		// if the node is removed from the linear system break all its bonds
		bondIndex -= numNeighbors;
		for(iNID=0 ; iNID<numNeighbors ; ++iNID){
		  bondDamage[bondIndex++] = 1.0;
		}
		totalDamage = numNeighbors;
	  }
	}

	if(numNeighbors > 0)
	  totalDamage /= numNeighbors;
	else
	  totalDamage = 0.0;
 	damage[nodeId] = totalDamage;
  }
}
void
PeridigmNS::ShortRangeForceContactModel::computeForce(const double dt,
                                                      const int numOwnedPoints,
                                                      const int* ownedIDs,
                                                      const int* contactNeighborhoodList,
                                                      PeridigmNS::DataManager& dataManager) const
{
  // Zero out the forces
  dataManager.getData(m_contactForceDensityFieldId, PeridigmField::STEP_NP1)->PutScalar(0.0);

  double *cellVolume, *y, *contactForce, *velocity;
  dataManager.getData(m_volumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&cellVolume);
  dataManager.getData(m_coordinatesFieldId, PeridigmField::STEP_NP1)->ExtractView(&y);
  dataManager.getData(m_velocityFieldId, PeridigmField::STEP_NP1)->ExtractView(&velocity);
  dataManager.getData(m_contactForceDensityFieldId, PeridigmField::STEP_NP1)->ExtractView(&contactForce);

  int neighborhoodListIndex(0), numNeighbors, nodeID, neighborID, iID, iNID;
  double nodeCurrentX[3], nodeCurrentV[3], nodeVolume, currentDistance, c, temp, neighborVolume;
  double normal[3], currentDotNormal, currentDotNeighbor, nodeCurrentVperp[3], nodeNeighborVperp[3], Vcm[3], nodeCurrentVrel[3], nodeNeighborVrel[3];
  double normCurrentVrel, normNeighborVrel, currentNormalForce[3], neighborNormalForce[3], normCurrentNormalForce, normNeighborNormalForce, currentFrictionForce[3], neighborFrictionForce[3];

  double currentDistanceSquared;
  double contactRadiusSquared = m_contactRadius*m_contactRadius;

  const double pi = boost::math::constants::pi<double>();

  for(iID=0 ; iID<numOwnedPoints ; ++iID){
    numNeighbors = contactNeighborhoodList[neighborhoodListIndex++];
    if(numNeighbors > 0){
      nodeID = ownedIDs[iID];
      nodeCurrentX[0] = y[nodeID*3];
      nodeCurrentX[1] = y[nodeID*3+1];
      nodeCurrentX[2] = y[nodeID*3+2];
      nodeCurrentV[0] = velocity[nodeID*3];
      nodeCurrentV[1] = velocity[nodeID*3+1];
      nodeCurrentV[2] = velocity[nodeID*3+2];
      nodeVolume = cellVolume[nodeID];
      for(iNID=0 ; iNID<numNeighbors ; ++iNID){
        neighborID = contactNeighborhoodList[neighborhoodListIndex++];
        TEUCHOS_TEST_FOR_EXCEPT_MSG(neighborID < 0, "Invalid neighbor list\n");
	currentDistanceSquared =  distanceSquared(nodeCurrentX[0], nodeCurrentX[1], nodeCurrentX[2],
						  y[neighborID*3], y[neighborID*3+1], y[neighborID*3+2]);
        if(currentDistanceSquared < contactRadiusSquared){
	  currentDistance = distance(nodeCurrentX[0], nodeCurrentX[1], nodeCurrentX[2],
				     y[neighborID*3], y[neighborID*3+1], y[neighborID*3+2]);
          c = 9.0*m_springConstant/(pi*m_horizon*m_horizon*m_horizon*m_horizon);	// half value (of 18) due to force being applied to both nodes
          temp = c*(m_contactRadius - currentDistance)/m_horizon;
          neighborVolume = cellVolume[neighborID];
          
          if (m_frictionCoefficient != 0.0){

            // calculate the perpendicular velocity of the current node wrt the vector between the nodes 

            normal[0] = (y[neighborID*3] - nodeCurrentX[0])/currentDistance;
            normal[1] = (y[neighborID*3+1] - nodeCurrentX[1])/currentDistance;
            normal[2] = (y[neighborID*3+2] - nodeCurrentX[2])/currentDistance;

            currentDotNormal = nodeCurrentV[0]*normal[0] + 
              nodeCurrentV[1]*normal[1] + 
              nodeCurrentV[2]*normal[2];

            currentDotNeighbor = velocity[neighborID*3]*normal[0] + 
              velocity[neighborID*3+1]*normal[1] + 
              velocity[neighborID*3+2]*normal[2]; 

            nodeCurrentVperp[0] = nodeCurrentV[0] - currentDotNormal*normal[0];
            nodeCurrentVperp[1] = nodeCurrentV[1] - currentDotNormal*normal[1];
            nodeCurrentVperp[2] = nodeCurrentV[2] - currentDotNormal*normal[2];

            nodeNeighborVperp[0] = velocity[neighborID*3] - currentDotNeighbor*normal[0];
            nodeNeighborVperp[1] = velocity[neighborID*3+1] - currentDotNeighbor*normal[1];
            nodeNeighborVperp[2] = velocity[neighborID*3+2] - currentDotNeighbor*normal[2];

            // calculate frame of reference for the perpendicular velocities

            Vcm[0] = 0.5*(nodeCurrentVperp[0] + nodeNeighborVperp[0]);
            Vcm[1] = 0.5*(nodeCurrentVperp[1] + nodeNeighborVperp[1]);
            Vcm[2] = 0.5*(nodeCurrentVperp[2] + nodeNeighborVperp[2]);
          
            // calculate the relative velocity of the current node wrt the neighboring node and vice versa

            nodeCurrentVrel[0] = nodeCurrentVperp[0] - Vcm[0];
            nodeCurrentVrel[1] = nodeCurrentVperp[1] - Vcm[1];
            nodeCurrentVrel[2] = nodeCurrentVperp[2] - Vcm[2];

            nodeNeighborVrel[0] = nodeNeighborVperp[0] - Vcm[0];
            nodeNeighborVrel[1] = nodeNeighborVperp[1] - Vcm[1];
            nodeNeighborVrel[2] = nodeNeighborVperp[2] - Vcm[2];

            normCurrentVrel = sqrt(nodeCurrentVrel[0]*nodeCurrentVrel[0] + 
                                   nodeCurrentVrel[1]*nodeCurrentVrel[1] + 
                                   nodeCurrentVrel[2]*nodeCurrentVrel[2]); 

            normNeighborVrel = sqrt(nodeNeighborVrel[0]*nodeNeighborVrel[0] + 
                                    nodeNeighborVrel[1]*nodeNeighborVrel[1] + 
                                    nodeNeighborVrel[2]*nodeNeighborVrel[2]);         
            
            // calculate the normal forces

            currentNormalForce[0] = -(temp*neighborVolume*(y[neighborID*3]   - nodeCurrentX[0])/currentDistance);
            currentNormalForce[1] = -(temp*neighborVolume*(y[neighborID*3+1] - nodeCurrentX[1])/currentDistance);
            currentNormalForce[2] = -(temp*neighborVolume*(y[neighborID*3+2] - nodeCurrentX[2])/currentDistance);

            neighborNormalForce[0] = (temp*nodeVolume*(y[neighborID*3]   - nodeCurrentX[0])/currentDistance);
            neighborNormalForce[1] = (temp*nodeVolume*(y[neighborID*3+1] - nodeCurrentX[1])/currentDistance);
            neighborNormalForce[2] = (temp*nodeVolume*(y[neighborID*3+2] - nodeCurrentX[2])/currentDistance);

            normCurrentNormalForce = sqrt(currentNormalForce[0]*currentNormalForce[0] + 
                                          currentNormalForce[1]*currentNormalForce[1] + 
                                          currentNormalForce[2]*currentNormalForce[2]);

            normNeighborNormalForce = sqrt(neighborNormalForce[0]*neighborNormalForce[0] + 
                                           neighborNormalForce[1]*neighborNormalForce[1] + 
                                           neighborNormalForce[2]*neighborNormalForce[2]);
            
            // calculate the friction forces

            if (normCurrentVrel != 0.0) {
              currentFrictionForce[0] = -m_frictionCoefficient*normCurrentNormalForce*nodeCurrentVrel[0]/normCurrentVrel;
              currentFrictionForce[1] = -m_frictionCoefficient*normCurrentNormalForce*nodeCurrentVrel[1]/normCurrentVrel;
              currentFrictionForce[2] = -m_frictionCoefficient*normCurrentNormalForce*nodeCurrentVrel[2]/normCurrentVrel;
            }
            else {
              currentFrictionForce[0] = 0.0;
              currentFrictionForce[1] = 0.0;
              currentFrictionForce[2] = 0.0;
            }

            if (normNeighborVrel != 0.0) {
              neighborFrictionForce[0] = -m_frictionCoefficient*normNeighborNormalForce*nodeNeighborVrel[0]/normNeighborVrel;
              neighborFrictionForce[1] = -m_frictionCoefficient*normNeighborNormalForce*nodeNeighborVrel[1]/normNeighborVrel;
              neighborFrictionForce[2] = -m_frictionCoefficient*normNeighborNormalForce*nodeNeighborVrel[2]/normNeighborVrel;
            }
            else {
              neighborFrictionForce[0] = 0.0;
              neighborFrictionForce[1] = 0.0;
              neighborFrictionForce[2] = 0.0;
            }

            // compute total contributions to force density

            contactForce[nodeID*3]       += currentNormalForce[0] + currentFrictionForce[0];
            contactForce[nodeID*3+1]     += currentNormalForce[1] + currentFrictionForce[1];
            contactForce[nodeID*3+2]     += currentNormalForce[2] + currentFrictionForce[2];
            contactForce[neighborID*3]   += neighborNormalForce[0] + neighborFrictionForce[0];
            contactForce[neighborID*3+1] += neighborNormalForce[1] + neighborFrictionForce[1];
            contactForce[neighborID*3+2] += neighborNormalForce[2] + neighborFrictionForce[2];
          }
          else {          

            // compute contributions to force density (Normal Force Only)

            contactForce[nodeID*3]       -= temp*neighborVolume*(y[neighborID*3]   - nodeCurrentX[0])/currentDistance;
            contactForce[nodeID*3+1]     -= temp*neighborVolume*(y[neighborID*3+1] - nodeCurrentX[1])/currentDistance;
            contactForce[nodeID*3+2]     -= temp*neighborVolume*(y[neighborID*3+2] - nodeCurrentX[2])/currentDistance;
            contactForce[neighborID*3]   += temp*nodeVolume*(y[neighborID*3]   - nodeCurrentX[0])/currentDistance;
            contactForce[neighborID*3+1] += temp*nodeVolume*(y[neighborID*3+1] - nodeCurrentX[1])/currentDistance;
            contactForce[neighborID*3+2] += temp*nodeVolume*(y[neighborID*3+2] - nodeCurrentX[2])/currentDistance;

          }
        }
      }
    }
  }
}
void
PeridigmNS::CriticalStretchDamageModel::computeDamage(const double dt,
                                                      const int numOwnedPoints,
                                                      const int* ownedIDs,
                                                      const int* neighborhoodList,
                                                      PeridigmNS::DataManager& dataManager) const
{
  double *x, *y, *damage, *bondDamageN, *bondDamageNP1, *deltaTemperature;
  dataManager.getData(m_modelCoordinatesFieldId, PeridigmField::STEP_NONE)->ExtractView(&x);
  dataManager.getData(m_coordinatesFieldId, PeridigmField::STEP_NP1)->ExtractView(&y);
  dataManager.getData(m_damageFieldId, PeridigmField::STEP_NP1)->ExtractView(&damage);
  dataManager.getData(m_bondDamageFieldId, PeridigmField::STEP_N)->ExtractView(&bondDamageN);
  dataManager.getData(m_bondDamageFieldId, PeridigmField::STEP_NP1)->ExtractView(&bondDamageNP1);
  deltaTemperature = NULL;
  if(m_applyThermalStrains)
    dataManager.getData(m_deltaTemperatureFieldId, PeridigmField::STEP_NP1)->ExtractView(&deltaTemperature);

  double trialDamage(0.0);
  int neighborhoodListIndex(0), bondIndex(0);
  int nodeId, numNeighbors, neighborID, iID, iNID;
  double nodeInitialX[3], nodeCurrentX[3], initialDistance, currentDistance, relativeExtension, totalDamage;

  // Set the bond damage to the previous value
  dataManager.getData(m_bondDamageFieldId, PeridigmField::STEP_NP1) =
    dataManager.getData(m_bondDamageFieldId, PeridigmField::STEP_N);

  // Update the bond damage
  // Break bonds if the extension is greater than the critical extension

  for(iID=0 ; iID<numOwnedPoints ; ++iID){
	nodeId = ownedIDs[iID];
	nodeInitialX[0] = x[nodeId*3];
	nodeInitialX[1] = x[nodeId*3+1];
	nodeInitialX[2] = x[nodeId*3+2];
	nodeCurrentX[0] = y[nodeId*3];
	nodeCurrentX[1] = y[nodeId*3+1];
	nodeCurrentX[2] = y[nodeId*3+2];
	numNeighbors = neighborhoodList[neighborhoodListIndex++];
	for(iNID=0 ; iNID<numNeighbors ; ++iNID){
	  neighborID = neighborhoodList[neighborhoodListIndex++];
      initialDistance = 
        distance(nodeInitialX[0], nodeInitialX[1], nodeInitialX[2],
                 x[neighborID*3], x[neighborID*3+1], x[neighborID*3+2]);
      currentDistance = 
        distance(nodeCurrentX[0], nodeCurrentX[1], nodeCurrentX[2],
                 y[neighborID*3], y[neighborID*3+1], y[neighborID*3+2]);
      if(m_applyThermalStrains)
        currentDistance -= m_alpha*deltaTemperature[nodeId]*initialDistance;
      relativeExtension = (currentDistance - initialDistance)/initialDistance;
      trialDamage = 0.0;
      if(relativeExtension > m_criticalStretch)
        trialDamage = 1.0;
      if(trialDamage > bondDamageNP1[bondIndex]){
        bondDamageNP1[bondIndex] = trialDamage;
      }
      bondIndex += 1;
    }
  }

  //  Update the element damage (percent of bonds broken)

  neighborhoodListIndex = 0;
  bondIndex = 0;
  for(iID=0 ; iID<numOwnedPoints ; ++iID){
	nodeId = ownedIDs[iID];
	numNeighbors = neighborhoodList[neighborhoodListIndex++];
    neighborhoodListIndex += numNeighbors;
	totalDamage = 0.0;
	for(iNID=0 ; iNID<numNeighbors ; ++iNID){
	  totalDamage += bondDamageNP1[bondIndex++];
	}
	if(numNeighbors > 0)
	  totalDamage /= numNeighbors;
	else
	  totalDamage = 0.0;
 	damage[nodeId] = totalDamage;
  }
}
void
PeridigmNS::VectorPoissonMaterial::computeForce(const double dt,
						const int numOwnedPoints,
						const int* ownedIDs,
						const int* neighborhoodList,
						PeridigmNS::DataManager& dataManager) const
{
  // Zero out the forces
  dataManager.getData(m_forceDensityFieldId, PeridigmField::STEP_NP1)->PutScalar(0.0);

  // Extract pointers to the underlying data
  double *volume, *x, *y, *f;

  dataManager.getData(m_volumeFieldId, PeridigmField::STEP_NONE)->ExtractView(&volume);
  dataManager.getData(m_modelCoordinatesFieldId, PeridigmField::STEP_NONE)->ExtractView(&x);
  dataManager.getData(m_coordinatesFieldId, PeridigmField::STEP_NP1)->ExtractView(&y);
  dataManager.getData(m_forceDensityFieldId, PeridigmField::STEP_NP1)->ExtractView(&f);

  int neighborhoodListIndex(0);
  int numNeighbors, neighborID, iID, iNID;
  double nodeInitialPosition[3], initialDistance;
  double nodeU, neighborU, kernel, nodeVolume, neighborVolume, temp, nodeForce, neighborForce;

  const double pi = boost::math::constants::pi<double>();

  for(iID=0 ; iID<numOwnedPoints ; ++iID){

    nodeVolume = volume[iID];
    nodeInitialPosition[0] = x[iID*3];
    nodeInitialPosition[1] = x[iID*3+1];
    nodeInitialPosition[2] = x[iID*3+2];

    numNeighbors = neighborhoodList[neighborhoodListIndex++];
    for(iNID=0 ; iNID<numNeighbors ; ++iNID){

      neighborID = neighborhoodList[neighborhoodListIndex++];
      neighborVolume = volume[neighborID];
      initialDistance = distance(nodeInitialPosition[0], nodeInitialPosition[1], nodeInitialPosition[2],
				 x[neighborID*3], x[neighborID*3+1], x[neighborID*3+2]);
      kernel = 3.0/(pi*m_horizon*m_horizon*m_horizon*m_horizon*initialDistance);

      // We are solving a Poisson equation
      // The function is lives in three-dimensional space and has a one-dimensional scalar output
      // Because the code infrastructure assumes both 3D input and output, we'll just solve three
      // instances of the problem at once.

      for(int eqn=0 ; eqn<3 ; ++eqn){

	nodeU = y[iID*3+eqn] - x[iID*3+eqn];
	neighborU = y[neighborID*3+eqn] - x[neighborID*3+eqn];
	temp = (neighborU - nodeU)*kernel;
	nodeForce = temp*neighborVolume;
	neighborForce = -temp*nodeVolume;
	TEUCHOS_TEST_FOR_EXCEPT_MSG(!boost::math::isfinite(nodeForce), "**** NaN detected in VectorPoissonMaterial::computeForce().\n");
	TEUCHOS_TEST_FOR_EXCEPT_MSG(!boost::math::isfinite(neighborForce), "**** NaN detected in VectorPoissonMaterial::computeForce().\n");
	f[iID*3+eqn] += nodeForce;
	f[neighborID*3+eqn] += neighborForce;

      }
    }
  }
}