// allocate crack and material velocity fields needed for time step on real nodes
// tried critical sections when nodes changed, but it was slower
// can't use ghost nodes, because need to test all on real nodes
//
// This task only used if have cracks or in multimaterial mode
void InitVelocityFieldsTask::Execute(void)
{
	CommonException *initErr = NULL;
	
	int tp = fmobj->GetTotalNumberOfPatches();
	
#pragma omp parallel
	{
		int nds[maxShapeNodes];
		double fn[maxShapeNodes];
		
		int pn = GetPatchNumber();
		
		// do non-rigid and rigid contact materials in patch pn
		for(int block=FIRST_NONRIGID;block<=FIRST_RIGID_CONTACT;block++)
		{   // get material point (only in this patch)
			MPMBase *mpmptr = patches[pn]->GetFirstBlockPointer(block);
			
			while(mpmptr!=NULL)
			{	const MaterialBase *matID = theMaterials[mpmptr->MatID()];		// material object for this particle
				const int matfld = matID->GetField();                           // material velocity field
				
				// get nodes and shape function for material point p
				const ElementBase *elref = theElements[mpmptr->ElemID()];		// element containing this particle
				
				// don't actually need shape functions, but need to screen out zero shape function
				// like done in subsequent tasks, otherwise node numbers will not align correctly
				// only thing used from return are numnds and nds
				int numnds;
				elref->GetShapeFunctions(&numnds,fn,nds,mpmptr);
				
				// Only need to decipher crack velocity field if has cracks (firstCrack!=NULL)
				//      and if this material allows cracks.
#ifdef COMBINE_RIGID_MATERIALS
				bool decipherCVF = firstCrack!=NULL && block!=FIRST_RIGID_CONTACT;
#else
				bool decipherCVF = firstCrack!=NULL;
#endif
				
				// Check each node
				for(int i=1;i<=numnds;i++)
				{	// use real node in this loop
					NodalPoint *ndptr = nd[nds[i]];
					
					// always zero when no cracks (or when ignoring cracks)
					short vfld = 0;
					
					// If need, find vlocity field and for each field set location
					// (above or below crack) and crack number (1 based) or 0 for NO_CRACK
					if(decipherCVF)
					{	// in CRAMP, find crack crossing and appropriate velocity field
						CrackField cfld[2];
						cfld[0].loc = NO_CRACK;			// NO_CRACK=0, ABOVE_CRACK=1, or BELOW_CRACK=2
						cfld[1].loc = NO_CRACK;
						int cfound=0;
						Vector norm;					// track normal vector for crack plane
						
						CrackHeader *nextCrack = firstCrack;
						while(nextCrack!=NULL)
						{	vfld = nextCrack->CrackCross(mpmptr->pos.x,mpmptr->pos.y,ndptr->x,ndptr->y,&norm);
							if(vfld!=NO_CRACK)
							{	cfld[cfound].loc=vfld;
								cfld[cfound].norm=norm;
								
#ifdef IGNORE_CRACK_INTERACTIONS
								// appears to always be same crack, and stop when found one
								cfld[cfound].crackNum=1;
								break;
#endif
								
								// Get crack number (default code does not ignore interactions)
								cfld[cfound].crackNum=nextCrack->GetNumber();
								cfound++;
								
								// stop if found two because code can only handle two interacting cracks
								// It exits loop now to go ahead with the first two found, by physics may be off
								if(cfound>1) break;
							}
							nextCrack=(CrackHeader *)nextCrack->GetNextObject();
						}
						
						
						// find (and allocate if needed) the velocity field
						// Use vfld=0 if no cracks found
						if(cfound>0)
						{   // In parallel, this is critical code
#pragma omp critical
							{   try
								{   vfld = ndptr->AddCrackVelocityField(matfld,cfld);
								}
								catch(CommonException err)
								{   if(initErr==NULL)
									initErr = new CommonException(err);
								}
							}
						}
						
						// set material point velocity field for this node
						mpmptr->vfld[i] = vfld;
					}
					
					// make sure material velocity field is created too
					// (Note: when maxMaterialFields==1 (Singe Mat Mode), mvf[0] is always there
					//        so no need to create it here)
					if(maxMaterialFields>1 && ndptr->NeedsMatVelocityField(vfld,matfld))
					{   // If parallel, this is critical code
#pragma omp critical
						{   try
							{   ndptr->AddMatVelocityField(vfld,matfld);
							}
							catch(CommonException err)
							{   if(initErr==NULL)
								initErr = new CommonException(err);
							}
						}
						
					}
				}
				
				// next material point
				mpmptr = (MPMBase *)mpmptr->GetNextObject();
			}
		}
	}
		
	// was there an error?
	if(initErr!=NULL) throw *initErr;
	
	// copy crack and material fields on real nodes to ghost nodes
	if(tp>1)
	{   for(int pn=0;pn<tp;pn++)
		patches[pn]->InitializationReduction();
	}
}
Esempio n. 2
0
// Calculate J and K at crack tips
CustomTask *PropagateTask::StepCalculation(void)
{
    // if not needed, just exit
    if(!doPropCalcs) return nextTask;
	
	// particle extrapolation if needed
    if(doEnergyBalanceCalcs)
    {   totalPlastic = 0.;
        totalPotential = 0.;
        for(int p=0;p<nmpmsNR;p++)
        {   MPMBase *mpnt = mpm[p];
            
            // track total energies in J = N-m
            //	mp is g, stored energy is N/m^2 cm^3/g, vel is mm/sec
            // workEnergy in J =  1.0e-6*mp*mpm[p]->GetWorkEnergy()
            // plastic 1.0e-6*mp*mpm[p]->GetPlastEnergy()
            // external work 1.e-9*mpm[p]->GetExtWork()
            // kinetic energy 0.5e-9*mp*(vel.x*vel.x+vel.y*vel.y)
            
            // plastic energy per unit thickness (units of N) (only needed energy balance crack growth)
            double mp = mpnt->mp;
            totalPlastic += 1.0e-3*mp*mpnt->GetPlastEnergy()/mpnt->thickness();
            //totalPotential += 1.0e-3*(mp*mpnt->GetStrainEnergy()
            //                        + 0.5e-3*mp*(mpnt->vel.x*mpnt->vel.x+mpnt->vel.y*mpnt->vel.y)
            //                        - 1.e-3*mpnt->GetExtWork())/mpnt->thickness();
			throw "external work is no longer available";
        }
    }
    
    CrackHeader *nextCrack;
    CrackSegment *crkTip;
    double cSize;
    int i,inMat;
    Vector tipDir,growTo,grow;
    int tipElem;
	char isAlt[10];
    
    // loop over cracks
    nextCrack=firstCrack;
    while(nextCrack!=NULL)
    {	// each crack tip
        for(i=START_OF_CRACK;i<=END_OF_CRACK;i++)
        {   // find crack tip and direction
			nextCrack->CrackTipAndDirection(i,&crkTip,tipDir);
            
            // crack propagation
            inMat=crkTip->tipMatnum;
            if(inMat>0)
            {	// crack tip terms last two times propagated
                crkTip->potential[2]=crkTip->potential[1];
                crkTip->potential[1]=crkTip->potential[0];
                crkTip->potential[0]=totalPotential;
                crkTip->plastic[2]=crkTip->plastic[1];
                crkTip->plastic[1]=crkTip->plastic[0];
                crkTip->plastic[0]=totalPlastic;
                crkTip->clength[2]=crkTip->clength[1];
                crkTip->clength[1]=crkTip->clength[0];
                crkTip->clength[0]=nextCrack->Length();
            
                // see if it grows
				int shouldGo=theMaterials[inMat-1]->ShouldPropagate(crkTip,tipDir,nextCrack,fmobj->np,0);
				isAlt[0] = 0;
				if(shouldGo==GROWNOW)
				{	nextCrack->SetAllowAlternate(i,FALSE);
				}
				else if(nextCrack->GetAllowAlternate(i))
				{	shouldGo=theMaterials[inMat-1]->ShouldPropagate(crkTip,tipDir,nextCrack,fmobj->np,1);
					if(shouldGo==GROWNOW) strcpy(isAlt," (alt)");
				}
				if(shouldGo==GROWNOW)
                {   theResult=GROWNOW;
                    tipElem=crkTip->planeInElem-1;
                    if(fabs(tipDir.x)>fabs(tipDir.y))
                        cSize=theElements[tipElem]->xmax-theElements[tipElem]->xmin;
                    else
                        cSize=theElements[tipElem]->ymax-theElements[tipElem]->ymin;
                    grow.x=cellsPerPropagationStep*cSize*tipDir.x;
                    grow.y=cellsPerPropagationStep*cSize*tipDir.y;
					
					// adjust if crossing another crack
					double p = nextCrack->AdjustGrowForCrossing(&grow,crkTip);
					
					// crack number and tip
					archiver->IncrementPropagationCounter();
                    cout << "# propagation" << isAlt << " crack-tip " << nextCrack->GetNumber() << "-" << i;
					
					// summarize
					cout << " at t=" << 1000*mtime << " with J=Jtip+Jzone : " << 1000.*crkTip->Jint.z <<
							" = " << 1000.*crkTip->Jint.x << " + " << 1000.*(crkTip->Jint.z-crkTip->Jint.x) << endl;
                    
                    // if jump is .7 or more cells, make more than 1 segment
                    int iseg,numSegs = 1;
                    if(p*cellsPerPropagationStep>.7) numSegs= 2*(p*cellsPerPropagationStep+.25);
                    CrackSegment *newCrkTip;
                    for(iseg=1;iseg<=numSegs;iseg++)
                    {   growTo.x=crkTip->x+(double)iseg*grow.x/(double)numSegs;
                        growTo.y=crkTip->y+(double)iseg*grow.y/(double)numSegs;
                        if(fmobj->dflag[0]==4) growTo.y=0.;			// force cutting simulation to stay in cut plane at 0
                        newCrkTip=nextCrack->Propagate(growTo,(int)i,theMaterials[inMat-1]->tractionMat[0]);
                    }
                    crkTip = newCrkTip;
                    
					if(crkTip!=NULL)
					{	// check if crack speed is being controlled
						if(theMaterials[inMat-1]->ControlCrackSpeed(crkTip,propTime))
							nextPropTime=mtime+propTime;
						
						// crack tip heating (if activated)
						if(ConductionTask::active)
							conduction->StartCrackTipHeating(crkTip,grow,nextCrack->GetThickness());
					}
                }
                else
                {   // when no growth, restore previous growth results
                    crkTip->potential[0]=crkTip->potential[1];
                    crkTip->potential[1]=crkTip->potential[2];
                    crkTip->plastic[0]=crkTip->plastic[1];
                    crkTip->plastic[1]=crkTip->plastic[2];
                    crkTip->clength[0]=crkTip->clength[1];
                    crkTip->clength[1]=crkTip->clength[2];
                }
            }
        }

        // next crack
        nextCrack=(CrackHeader *)nextCrack->GetNextObject();
    }
    
    return nextTask;
}
Esempio n. 3
0
// Get mass matrix, find dimensionless particle locations,
//	and find grid momenta
void InitializationTask::Execute(void)
{
	CommonException *initErr = NULL;
	
	// Zero Mass Matrix and vectors
	warnings.BeginStep();
    
	int tp = fmobj->GetTotalNumberOfPatches();
#pragma omp parallel
	{
        // zero all nodal variables on real nodes
#pragma omp for
		for(int i=1;i<=nnodes;i++)
			nd[i]->InitializeForTimeStep();
		
        // zero ghost nodes in patch for this thread
        int pn = GetPatchNumber();
        patches[pn]->InitializeForTimeStep();

		// particle calculations get xipos for particles and if doing CPDI
        // precalculate CPDI info needed for subsequent shape functions
#pragma omp for nowait
		for(int p=0;p<nmpmsRC;p++)
        {   MPMBase *mpmptr = mpm[p];                                       // pointer
			const ElementBase *elref = theElements[mpmptr->ElemID()];		// element containing this particle
			try
			{	elref->GetShapeFunctionData(mpmptr);
			}
			catch(CommonException err)
			{	if(initErr==NULL)
				{
#pragma omp critical
					initErr = new CommonException(err);
				}
			}
		}
	}
	
	// was there an error?
	if(initErr!=NULL) throw *initErr;
    
	// allocate crack and material velocity fields needed for time step on real nodes
    // tried critical sections when nodes changed, but it was slower
    // can't use ghost nodes, because need to test all on real nodes
	if(firstCrack!=NULL || maxMaterialFields>1)
	{
#pragma omp parallel
        {
            int nds[maxShapeNodes];
            double fn[maxShapeNodes];
            
            //for(int pn=0;pn<tp;pn++) {
            int pn = GetPatchNumber();
		
            // do non-rigid and rigid contact materials in patch pn
            for(int block=FIRST_NONRIGID;block<=FIRST_RIGID_CONTACT;block++)
            {   // get material point (only in this patch)
                MPMBase *mpmptr = patches[pn]->GetFirstBlockPointer(block);
                
                while(mpmptr!=NULL)
                {	const MaterialBase *matID = theMaterials[mpmptr->MatID()];		// material object for this particle
                    const int matfld = matID->GetField();                           // material velocity field
                    
                    // get nodes and shape function for material point p
                    const ElementBase *elref = theElements[mpmptr->ElemID()];		// element containing this particle
                    
                    // don't actually need shape functions, but need to screen out zero shape function
                    // like done in subsequent tasks, otherwise node numbers will not align correctly
                    // only think used from return are numnds and nds
                    int numnds;
                    elref->GetShapeFunctions(&numnds,fn,nds,mpmptr);
                    
                    // Add particle property to each node in the element
                    for(int i=1;i<=numnds;i++)
                    {	// use real node in this loop
                        NodalPoint *ndptr = nd[nds[i]];
                        
                        // always zero when no cracks
                        short vfld = 0;
#ifdef COMBINE_RIGID_MATERIALS
                        // when combining rigid particles, extrapolate all to field 0 and later
                        // copy to other active fields
                        if(firstCrack!=NULL && block!=FIRST_RIGID_CONTACT)
#else
                        if(firstCrack!=NULL)
#endif
                        {	// in CRAMP, find crack crossing and appropriate velocity field
                            CrackField cfld[2];
                            cfld[0].loc = NO_CRACK;			// NO_CRACK, ABOVE_CRACK, or BELOW_CRACK
                            cfld[1].loc = NO_CRACK;
                            int cfound=0;
                            Vector norm;
                            
                            CrackHeader *nextCrack = firstCrack;
                            while(nextCrack!=NULL)
                            {	vfld = nextCrack->CrackCross(mpmptr->pos.x,mpmptr->pos.y,ndptr->x,ndptr->y,&norm);
                                if(vfld!=NO_CRACK)
                                {	cfld[cfound].loc=vfld;
                                    cfld[cfound].norm=norm;
#ifdef IGNORE_CRACK_INTERACTIONS
                                    cfld[cfound].crackNum=1;	// appears to always be same crack, and stop when found one
                                    break;
#else
                                    cfld[cfound].crackNum=nextCrack->GetNumber();
                                    cfound++;
                                    if(cfound>1) break;			// stop if found two, if there are more then two, physics will be off
#endif
                                }
                                nextCrack=(CrackHeader *)nextCrack->GetNextObject();
                            }
                                
                            
                            // find (and allocate if needed) the velocity field
                            // Use vfld=0 if no cracks found
                            if(cfound>0)
                            {   // In parallel, this is critical code
#pragma omp critical
                                {   try
                                    {   vfld = ndptr->AddCrackVelocityField(matfld,cfld);
                                    }
                                    catch(CommonException err)
                                    {   if(initErr==NULL)
                                            initErr = new CommonException(err);
                                    }
                                }
                            }
                            
                            // set material point velocity field for this node
                            mpmptr->vfld[i] = vfld;
                        }
                        
                        // make sure material velocity field is created too
                        if(maxMaterialFields>1 && ndptr->NeedsMatVelocityField(vfld,matfld))
                        {   // If parallel, this is critical code
#pragma omp critical
                            {   try
                                {   ndptr->AddMatVelocityField(vfld,matfld);
                                }
                                catch(CommonException err)
                                {   if(initErr==NULL)
                                        initErr = new CommonException(err);
                                }
                            }
                            
                        }
                    }
                    
                    // next material point
                    mpmptr = (MPMBase *)mpmptr->GetNextObject();
                }
            }
            //}           // end for loop when not in parallel
		}
    
        // was there an error?
        if(initErr!=NULL) throw *initErr;

		// copy crack and material fields on real nodes to ghost nodes
		if(tp>1)
        {   for(int pn=0;pn<tp;pn++)
				patches[pn]->InitializationReduction();
		}
	}
    	
    // Update forces applied to particles
	MatPtLoadBC::SetParticleFext(mtime);
	
	// remove contact conditions
	CrackNode::RemoveCrackNodes();
	MaterialInterfaceNode::RemoveInterfaceNodes();
	
    // turn off isothermal ramp when done and ramp step initialization
	thermal.CheckDone(mtime);
	
}	
Esempio n. 4
0
// Calculate J and K at crack tips
CustomTask *PropagateTask::StepCalculation(void)
{
    // if not needed, just exit
    if(!doPropCalcs) return nextTask;
	
    CrackHeader *nextCrack;
    CrackSegment *crkTip;
    double cSize;
    int i,inMat;
    Vector tipDir,growTo,grow;
    int tipElem;
	char isAlt[10];
    
    // loop over cracks
    nextCrack=firstCrack;
    while(nextCrack!=NULL)
    {	// each crack tip
        for(i=START_OF_CRACK;i<=END_OF_CRACK;i++)
        {   // find crack tip and direction
			nextCrack->CrackTipAndDirection(i,&crkTip,tipDir);
            
            // crack propagation
            inMat=crkTip->tipMatnum;
            if(inMat>0)
            {	// see if it grows
				int shouldGo=theMaterials[inMat-1]->ShouldPropagate(crkTip,tipDir,nextCrack,fmobj->np,0);
				isAlt[0] = 0;
				if(shouldGo==GROWNOW)
				{	nextCrack->SetAllowAlternate(i,false);
				}
				else if(nextCrack->GetAllowAlternate(i))
				{	shouldGo=theMaterials[inMat-1]->ShouldPropagate(crkTip,tipDir,nextCrack,fmobj->np,1);
					if(shouldGo==GROWNOW) strcpy(isAlt," (alt)");
				}
				if(shouldGo==GROWNOW)
                {   theResult=GROWNOW;
                    tipElem=crkTip->planeElemID();
                    if(fabs(tipDir.x)>fabs(tipDir.y))
                        cSize=theElements[tipElem]->xmax-theElements[tipElem]->xmin;
                    else
                        cSize=theElements[tipElem]->ymax-theElements[tipElem]->ymin;
                    grow.x=cellsPerPropagationStep*cSize*tipDir.x;
                    grow.y=cellsPerPropagationStep*cSize*tipDir.y;
					
					// adjust if crossing another crack - adjusts grow is need and returns resulting relative change (in p)
					double p = nextCrack->AdjustGrowForCrossing(&grow,crkTip,cSize,&tipDir);
					
					// crack number and tip
					archiver->IncrementPropagationCounter();
                    cout << "# propagation" << isAlt << " crack-tip " << nextCrack->GetNumber() << "-" << i;
					
					// summarize
					cout << " at t=" << mtime*UnitsController::Scaling(1.e3)
						<< " with J=Jtip+Jzone : " << crkTip->Jint.z*UnitsController::Scaling(1.e-3)
						<< " = " << crkTip->Jint.x*UnitsController::Scaling(1.e-3)
						<< " + " << (crkTip->Jint.z-crkTip->Jint.x)*UnitsController::Scaling(1.e-3) << endl;
                    
                    // if jump is .7 or more cells, divide propagation into multiple segments
                    int iseg,numSegs = 1;
                    if(p*cellsPerPropagationStep>.7) numSegs= (int)(2*(p*cellsPerPropagationStep+.25));
                    CrackSegment *newCrkTip = NULL;
                    for(iseg=1;iseg<=numSegs;iseg++)
                    {   growTo.x=crkTip->x+(double)iseg*grow.x/(double)numSegs;
                        growTo.y=crkTip->y+(double)iseg*grow.y/(double)numSegs;
                        newCrkTip=nextCrack->Propagate(growTo,(int)i,theMaterials[inMat-1]->tractionMat[0]);
                    }
                    crkTip = newCrkTip;
                    
					if(crkTip!=NULL)
					{	// check if crack speed is being controlled
						if(theMaterials[inMat-1]->ControlCrackSpeed(crkTip,propTime))
							nextPropTime=mtime+propTime;
						
						// crack tip heating (if activated)
						if(ConductionTask::active)
							conduction->StartCrackTipHeating(crkTip,grow,nextCrack->GetThickness());
					}
                }
            }
        }

        // next crack
        nextCrack=(CrackHeader *)nextCrack->GetNextObject();
    }
    
    return nextTask;
}