/* Get list of nodes for this element and the number of nodes nodes will be in nds[1]... */ void ElementBase::GetNodes(int *numnds,int *nds) const { int i; *numnds = NumberNodes(); for(i=1;i<=*numnds;i++) nds[i] = nodes[i-1]; }
/* find next node (in ccw direction) from specified node Assumes nodes[NumberNodes()]=nodes[0] return 0 if gridNode not in this element */ int ElementBase::NextNode(int gridNode) { int i; for(i=0;i<NumberNodes();i++) { if(nodes[i]==gridNode) return nodes[i+1]; } return 0; }
/* find an edge (in ccw direction) with a pair of nodes Assumes nodes[NumberNodes()]=nodes[0] return edge number (0 based) or -1 if no such edge */ int ElementBase::FindEdge(int beginNode,int endNode) { int i; for(i=0;i<NumberNodes();i++) { if(nodes[i]==beginNode && nodes[i+1]==endNode) return i; } return -1; }
// Find Cartesion position from natural coordinates void ElementBase::GetPosition(Vector *xipos,Vector *xyzpos) { double fn[MaxElNd]; ShapeFunction(xipos,FALSE,&fn[1],NULL,NULL,NULL); int i; ZeroVector(xyzpos); for(i=1;i<=NumberNodes();i++) { xyzpos->x += fn[i]*nd[nodes[i-1]]->x; xyzpos->y += fn[i]*nd[nodes[i-1]]->y; xyzpos->z += fn[i]*nd[nodes[i-1]]->z; } }
// increment element of an interface stiffness matrix void Interface2D::IncrementStiffnessElements(double dStiff,double *fn, double xpxi,double ypxi,double dlxi,double Dn,double Dt) { int irow,jcol,nst=2*NumberNodes(); int n1,n2,term; // loop over upper half diagonal of stiffness matrix for(irow=1;irow<=nst;irow++) { for(jcol=irow;jcol<=nst;jcol++) { // get term type and shape function numbers if(IsEven(irow)) { n1=irow/2; if(IsEven(jcol)) { n2=jcol/2; term=3; // both even } else { n2=(jcol+1)/2; term=2; // even, odd } } else { n1=(irow+1)/2; if(IsEven(jcol)) { n2=jcol/2; term=2; // odd, even } else { n2=(jcol+1)/2; term=1; // both odd } } /* increment element depending on term and current weight type Note: signs attached to shape functions take care of minus signs for matrix elements with displacement on opposite surface */ switch(term) { case 1: // both odd se[irow][jcol]+=dStiff*fn[n1]*fn[n2]*(Dt*xpxi*xpxi + Dn*ypxi*ypxi)/dlxi; break; case 2: // one even and one odd se[irow][jcol]+=dStiff*fn[n1]*fn[n2]*(Dt-Dn)*xpxi*ypxi/dlxi; break; case 3: // both even se[irow][jcol]+=dStiff*fn[n1]*fn[n2]*(Dt*ypxi*ypxi + Dn*xpxi*xpxi)/dlxi; break; default: break; } } } }
/* Find node of this element closest to specified point Only checks nodes of this element. In distorted mesh, could be closer node in neighboring element, even if this point is in this element */ int ElementBase::NearestNode(double x,double y,int *secondChoice) { int ind,nearest,nextNearest; double dist,distMin,distNextMin,dltx,dlty; int k,numnds=NumberNodes(); // distance to first node nearest=nodes[0]; dltx=nd[nearest]->x-x; dlty=nd[nearest]->y-y; distMin=dltx*dltx+dlty*dlty; // distance to second node nextNearest=nodes[1]; dltx=nd[nextNearest]->x-x; dlty=nd[nextNearest]->y-y; distNextMin=dltx*dltx+dlty*dlty; if(distNextMin<distMin) { dist=distMin; distMin=distNextMin; distNextMin=distMin; nearest=nodes[1]; nextNearest=nodes[0]; } // check the rest for(k=2;k<numnds;k++) { ind=nodes[k]; dltx=nd[ind]->x-x; dlty=nd[ind]->y-y; dist=dltx*dltx+dlty*dlty; if(dist<distMin) { distNextMin=distMin; nextNearest=nearest; distMin=dist; nearest=ind; } else if(dist<distNextMin) { distNextMin=dist; nextNearest=ind; } } // return the node *secondChoice=nextNearest; return nearest; }
// Ignore zero shape functions in GIMP, but only need to check remote nodes void ElementBase::GimpCompact(int *numnds,int *nds,double *sfxn,double *xDeriv,double *yDeriv,double *zDeriv) const { int local=NumberNodes(); int i,found=local; for(i=local+1;i<=*numnds;i++) { if(sfxn[i]>0.0) { found++; nds[found]=nds[i]; sfxn[found]=sfxn[i]; if(xDeriv!=NULL) { xDeriv[found]=xDeriv[i]; yDeriv[found]=yDeriv[i]; if(zDeriv!=NULL) zDeriv[found]=zDeriv[i]; } } } *numnds=found; }
/* Find dimensionless coordinates by numerical methods input: pos is position in the element output: xipos is dimensionless position only used in MPM and only here if non-rectangular elements */ void ElementBase::GetXiPos(Vector *pos,Vector *xipos) const { double xt,yt,dxxi,dxeta,dyxi,dyeta; double deter,dxi,deta,dist; double gfn[MaxElNd],gdfnxi[MaxElNd],gdfnet[MaxElNd]; int numnds=NumberNodes(),i,j; // initial guess GetCentroid(xipos); // nodal coordinates Vector eNode[MaxElNd]; for(i=0;i<numnds;i++) { eNode[i].x=nd[nodes[i]]->x; eNode[i].y=nd[nodes[i]]->y; } /* solve for xipos using Newton-Rapheson (see FEA Notes) using shape functions and their derivatives */ for(j=1;j<=MAXITER;j++) { ShapeFunction(xipos,TRUE,gfn,gdfnxi,gdfnet,NULL,NULL,NULL,NULL); xt=-pos->x; yt=-pos->y; dxxi=0.; dxeta=0.; dyxi=0.; dyeta=0.; for(i=0;i<numnds;i++) { xt+=eNode[i].x*gfn[i]; yt+=eNode[i].y*gfn[i]; dxxi+=eNode[i].x*gdfnxi[i]; dxeta+=eNode[i].x*gdfnet[i]; dyxi+=eNode[i].y*gdfnxi[i]; dyeta+=eNode[i].y*gdfnet[i]; } deter=dxxi*dyeta-dxeta*dyxi; dxi=(-xt*dyeta+yt*dxeta)/deter; deta=(xt*dyxi-yt*dxxi)/deter; xipos->x+=dxi; xipos->y+=deta; dist=sqrt(dxi*dxi+deta*deta); if(dist<.001) break; } }
/* Calculate area of element (in mm^2 because nodes in mm) nodes is pointer to 0-based array NodalPoints */ double Quad2D::GetArea(void) const { short i,ns=NumberSides(),nn=NumberNodes()-1; double area; area=0.; for(i=0;i<ns-1;i++) { area+=nd[nodes[i]]->x*nd[nodes[i+ns]]->y -nd[nodes[i]]->y*nd[nodes[i+ns]]->x +nd[nodes[i+ns]]->x*nd[nodes[i+1]]->y -nd[nodes[i+ns]]->y*nd[nodes[i+1]]->x; } ns--; area+=nd[nodes[ns]]->x*nd[nodes[nn]]->y -nd[nodes[ns]]->y*nd[nodes[nn]]->x +nd[nodes[nn]]->x*nd[nodes[0]]->y -nd[nodes[nn]]->y*nd[nodes[0]]->x; return area/2.; }
/* Calculate Element forces, stresses, and strain energy */ void Interface2D::ForceStress(double *rm,int np,int nfree) { int numnds=NumberNodes(); int i,j,nst=2*numnds; int indg; int nameEl=ElementName(); double Fxy[(MaxElNd-1)*MxFree+1]; double dx,dy,len,Dn,Dt; double pi=3.141592653589793,dsx,dsy,xpxi,ypxi; MaterialBase *matl=theMaterials[material-1]; // Get stiffness matrix (and also load nodal coordinates (ce[]) and material props (pr.C[][])) Stiffness(np); // Load nodal displacements into re[] int ind=0; for(j=1;j<=numnds;j++) { indg=nfree*(nodes[j-1]-1); for(i=1;i<=nfree;i++) re[++ind]=rm[indg+i]; } // forces in cartensian frame for(i=1;i<=nst;i++) { Fxy[i]=0.; for(j=1;j<=nst;j++) { if(np!=AXI_SYM) Fxy[i]+=se[i][j]*re[j]; else Fxy[i]+=2.*pi*se[i][j]*re[j]; } } // transfer force to se[][] and get strain energy of interphase strainEnergy=0.; for(i=1;i<=nst;i++) { strainEnergy+=0.5*re[i]*Fxy[i]; se[i][7]=Fxy[i]; } for(i=1;i<=numnds;i++) { se[i][2]=0.; se[i][4]=0.; } // get normal and shear stresses from interface law Dn=matl->pr.C[1][1]; Dt=matl->pr.C[1][2]; if(nameEl==LINEAR_INTERFACE) { dx=ce[2].x-ce[1].x; dy=ce[2].y-ce[1].y; len=sqrt(dx*dx + dy*dy); // nodes 1 and 4 and then 2 and 3 InterfaceTraction(1,4,dx,dy,len,Dn,Dt); InterfaceTraction(2,3,dx,dy,len,Dn,Dt); } else { dx=(ce[3].x-ce[1].x)/2.; dy=(ce[3].y-ce[1].y)/2.; dsx=ce[1].x+ce[3].x-2.*ce[2].x; dsy=ce[1].y+ce[3].y-2.*ce[2].y; // nodes 1 and 6 xpxi=dx-dsx; ypxi=dy-dsy; len=sqrt(xpxi*xpxi + ypxi*ypxi); InterfaceTraction(1,6,xpxi,ypxi,len,Dn,Dt); // nodes 2 and 5 xpxi=dx; ypxi=dy; len=sqrt(xpxi*xpxi + ypxi*ypxi); InterfaceTraction(2,5,xpxi,ypxi,len,Dn,Dt); // nodes 3 and 4 xpxi=dx+dsx; ypxi=dy+dsy; len=sqrt(xpxi*xpxi + ypxi*ypxi); InterfaceTraction(3,4,xpxi,ypxi,len,Dn,Dt); } }
// return dimensionless location for material points void ElementBase::MPMPoints(short numPerElement,Vector *mpos) const { int i,j,k; double fxn[MaxElNd],row,zrow; if(NumberSides()==4) { switch(numPerElement) { case 4: // ENI or FNI - 2D only mpos[0].x=-.5; mpos[0].y=-.5; mpos[0].z=0.; mpos[1].x=.5; mpos[1].y=-.5; mpos[1].z=0.; mpos[2].x=-.5; mpos[2].y=.5; mpos[2].z=0.; mpos[3].x=.5; mpos[3].y=.5; mpos[3].z=0.; break; case 1: // CM of square or brick mpos[0].x=0.; mpos[0].y=0.; mpos[0].z=0.; break; case 8: // 3D box mpos[0].x=-.5; mpos[0].y=-.5; mpos[0].z=-.5; mpos[1].x=.5; mpos[1].y=-.5; mpos[1].z=-.5; mpos[2].x=-.5; mpos[2].y=.5; mpos[2].z=-.5; mpos[3].x=.5; mpos[3].y=.5; mpos[3].z=-.5; mpos[4].x=-.5; mpos[4].y=-.5; mpos[4].z=.5; mpos[5].x=.5; mpos[5].y=-.5; mpos[5].z=.5; mpos[6].x=-.5; mpos[6].y=.5; mpos[6].z=.5; mpos[7].x=.5; mpos[7].y=.5; mpos[7].z=.5; break; case 9: // 2D k=0; row = -2./3.; for(j=0;j<3;j++) { mpos[k].x=-2./3.; mpos[k].y=row; mpos[k].z=0.; mpos[k+1].x=0.; mpos[k+1].y=row; mpos[k+1].z=0.; mpos[k+2].x=2./3.; mpos[k+2].y=row; mpos[k+2].z=0.; k += 3; row += 2./3.; } break; case 16: // 2D k=0; row = -0.75; for(j=0;j<4;j++) { mpos[k].x=-0.75; mpos[k].y=row; mpos[k].z=0.; mpos[k+1].x=-.25; mpos[k+1].y=row; mpos[k+1].z=0.; mpos[k+2].x=.25; mpos[k+2].y=row; mpos[k+2].z=0.; mpos[k+3].x=.75; mpos[k+3].y=row; mpos[k+3].z=0.; k += 4; row += 0.5; } break; case 25: // 2D k=0; row = -0.8; for(j=0;j<5;j++) { mpos[k].x=-0.8; mpos[k].y=row; mpos[k].z=0.; mpos[k+1].x=-.4; mpos[k+1].y=row; mpos[k+1].z=0.; mpos[k+2].x=0.; mpos[k+2].y=row; mpos[k+2].z=0.; mpos[k+3].x=.4; mpos[k+3].y=row; mpos[k+3].z=0.; mpos[k+4].x=.8; mpos[k+4].y=row; mpos[k+4].z=0.; k += 5; row += 0.4; } break; case 27: // 3D k=0; zrow = -2./3.; for(i=0;i<3;i++) { row = -2./3.; for(j=0;j<3;j++) { mpos[k].x=-2./3.; mpos[k].y=row; mpos[k].z=zrow; mpos[k+1].x=0.; mpos[k+1].y=row; mpos[k+1].z=zrow; mpos[k+2].x=2./3.; mpos[k+2].y=row; mpos[k+2].z=zrow; k += 3; row += 2./3.; } zrow += 2./3.; } break; default: throw CommonException("Invalid number of material points per element.","ElementBase::MPMPoints"); break; } } // covert to x-y-z locations for(k=0;k<numPerElement;k++) { ShapeFunction(&mpos[k],FALSE,fxn,NULL,NULL,NULL); ZeroVector(&mpos[k]); for(j=0;j<NumberNodes();j++) { mpos[k].x+=nd[nodes[j]]->x*fxn[j]; mpos[k].y+=nd[nodes[j]]->y*fxn[j]; mpos[k].z+=nd[nodes[j]]->z*fxn[j]; } } }
/* Calculate Stiffness Matrix */ void CSTriangle::Stiffness(int np) { double detjac,asr,dv; double xiDeriv[MaxElNd],etaDeriv[MaxElNd],asbe[MaxElNd],sfxn[MaxElNd]; double bte[MxFree*MaxElNd][5],temp; double thck=thickness,deltaT; int numnds=NumberNodes(); int ind1,ind2,i,j,irow,jcol,nst=2*numnds; MaterialBase *matl=theMaterials[material-1]; Vector xi; // Load nodal coordinates (ce[]), temperature (te[]), and // material props (pr.C[][] and pr.alpha[]) GetProperties(np); // Zero upper hald element stiffness matrix (se[]) and reaction vector (re[]) for(irow=1;irow<=nst;irow++) { re[irow]=0.; for(jcol=irow;jcol<=nst;jcol++) se[irow][jcol]=0.; } /* Call shape routine to calculate element B (be,asbe) matrix and the determinant of the Jacobian - both at centriod */ xi.x=(ce[1].x+ce[2].x+ce[3].x)/3.; xi.y=(ce[1].y+ce[2].y+ce[3].y)/3.; ShapeFunction(&xi,BMATRIX,&sfxn[1],&xiDeriv[1],&etaDeriv[1],&ce[1], &detjac,&asr,&asbe[1]); /* Form matrix product BT E - exploit known sparcity of B matrix and only include multiplications by nonzero elements */ deltaT=0.; ind1=-1; if(np!=AXI_SYM) { dv=thck*detjac; for(i=1;i<=numnds;i++) { ind1=ind1+2; ind2=ind1+1; for(j=1;j<=3;j++) { bte[ind1][j]=dv*(xiDeriv[i]*matl->pr.C[1][j]+etaDeriv[i]*matl->pr.C[3][j]); bte[ind2][j]=dv*(etaDeriv[i]*matl->pr.C[2][j]+xiDeriv[i]*matl->pr.C[3][j]); } deltaT+=te[i]*sfxn[i]; } } else { dv=asr*detjac; for(i=1;i<=numnds;i++) { ind1=ind1+2; ind2=ind1+1; for(j=1;j<=4;j++) { bte[ind1][j]=dv*(xiDeriv[i]*matl->pr.C[1][j]+etaDeriv[i]*matl->pr.C[3][j] +asbe[i]*matl->pr.C[4][j]); bte[ind2][j]=dv*(etaDeriv[i]*matl->pr.C[2][j]+xiDeriv[i]*matl->pr.C[3][j]); } deltaT+=te[i]*sfxn[i]; } } /* Form stiffness matrix by getting BT E B and add in initial strains into element load vector */ for(irow=1;irow<=nst;irow++) { for(j=1;j<=3;j++) { re[irow]+=bte[irow][j]*matl->pr.alpha[j]*deltaT; } if(np!=AXI_SYM) { for(jcol=irow;jcol<=nst;jcol++) { if(IsEven(jcol)) { ind1=jcol/2; temp=bte[irow][2]*etaDeriv[ind1] +bte[irow][3]*xiDeriv[ind1]; } else { ind1=(jcol+1)/2; temp=bte[irow][1]*xiDeriv[ind1] +bte[irow][3]*etaDeriv[ind1]; } se[irow][jcol]+=temp; } } else { re[irow]+=bte[irow][4]*matl->pr.alpha[j]*deltaT; for(jcol=irow;jcol<=nst;jcol++) { if(IsEven(jcol)) { ind1=jcol/2; temp=bte[irow][2]*etaDeriv[ind1] +bte[irow][3]*xiDeriv[ind1]; } else { ind1=(jcol+1)/2; temp=+bte[irow][1]*xiDeriv[ind1] +bte[irow][3]*etaDeriv[ind1] +bte[irow][4]*asbe[ind1]; } se[irow][jcol]+=temp; } } } /* Fill in lower half of stiffness matrix */ for(irow=1;irow<=nst-1;irow++) { for(jcol=irow+1;jcol<=nst;jcol++) { se[jcol][irow]=se[irow][jcol]; } } }