Changeset 5877
- Timestamp:
- 09/17/10 17:18:54 (15 years ago)
- Location:
- issm/trunk/src/c/objects/Elements
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk/src/c/objects/Elements/Penta.cpp
r5876 r5877 709 709 break; 710 710 case DiagnosticVertAnalysisEnum: 711 CreateKMatrixDiagnosticVert( Kgg); 711 Ke=CreateKMatrixDiagnosticVert(); 712 if(Ke) Ke->AddToGlobal(Kgg,NULL,NULL); 713 delete Ke; 712 714 break; 713 715 case BedSlopeXAnalysisEnum: case SurfaceSlopeXAnalysisEnum: case BedSlopeYAnalysisEnum: case SurfaceSlopeYAnalysisEnum: … … 2281 2283 ElementMatrix* Penta::CreateKMatrixDiagnosticMacAyeal3d(void){ 2282 2284 2283 2284 2285 /*compute all stiffness matrices for this element*/ 2285 2286 ElementMatrix* Ke1=CreateKMatrixDiagnosticMacAyeal3dViscous(); … … 2719 2720 /*}}}*/ 2720 2721 /*FUNCTION Penta::CreateKMatrixDiagnosticVert {{{1*/ 2721 void Penta::CreateKMatrixDiagnosticVert( Mat Kgg){ 2722 ElementMatrix* Penta::CreateKMatrixDiagnosticVert(void){ 2723 2724 /*compute all stiffness matrices for this element*/ 2725 ElementMatrix* Ke1=CreateKMatrixDiagnosticVertVolume(); 2726 ElementMatrix* Ke2=CreateKMatrixDiagnosticVertVolume(); 2727 ElementMatrix* Ke =new ElementMatrix(Ke1,Ke2); 2728 2729 /*clean-up and return*/ 2730 delete Ke1; 2731 delete Ke2; 2732 return Ke; 2733 } 2734 /*}}}*/ 2735 /*FUNCTION Penta::CreateKMatrixDiagnosticVertVolume {{{1*/ 2736 ElementMatrix* Penta::CreateKMatrixDiagnosticVertVolume(void){ 2722 2737 2723 2738 /* local declarations */ 2724 int i,j;2725 2726 /* node data: */2727 2739 const int numdof=NDOF1*NUMVERTICES; 2740 int i,j,ig; 2728 2741 double xyz_list[NUMVERTICES][3]; 2729 int* doflist=NULL;2730 2731 /* 3d gaussian points: */2732 int ig;2733 2742 GaussPenta *gauss=NULL; 2734 2735 /* matrices: */2736 2743 double Ke_gg[numdof][numdof]={0.0}; 2737 double Ke_gg_gaussian[numdof][numdof];2738 2744 double Jdet; 2739 2745 double B[NDOF1][NUMVERTICES]; … … 2741 2747 double DL_scalar; 2742 2748 2743 /*Collapsed formulation: */ 2744 Tria* tria=NULL; 2745 2746 /*If on water, skip stiffness: */ 2747 if(IsOnWater())return; 2748 2749 /*If this element is on the surface, we have a dynamic boundary condition that applies, as a stiffness 2750 * matrix: */ 2751 if(IsOnSurface()){ 2752 tria=(Tria*)SpawnTria(3,4,5); //nodes 3,4 and 5 are on the surface 2753 tria->CreateKMatrixDiagnosticSurfaceVert(Kgg); 2754 delete tria->matice; delete tria; 2755 } 2756 2757 /*Now, onto the formulation for the vertical velocity: */ 2758 2759 /* Get node coordinates and dof list: */ 2749 /*Initialize Element matrix and return if necessary*/ 2750 if(IsOnWater()) return NULL; 2751 ElementMatrix* Ke=this->NewElementMatrix(NoneApproximationEnum); 2752 2753 /*Retrieve all inputs and parameters*/ 2760 2754 GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES); 2761 GetDofList(&doflist,NoneApproximationEnum,GsetEnum);2762 2755 2763 2756 /* Start looping on the number of gaussian points: */ … … 2767 2760 gauss->GaussPoint(ig); 2768 2761 2769 /*Get B and Bprime matrices: */2770 2762 GetBVert(&B[0][0], &xyz_list[0][0], gauss); 2771 2763 GetBprimeVert(&Bprime[0][0], &xyz_list[0][0], gauss); 2772 2764 2773 /* Get Jacobian determinant: */2774 2765 GetJacobianDeterminant(&Jdet, &xyz_list[0][0],gauss); 2775 2766 DL_scalar=gauss->weight*Jdet; 2776 2767 2777 /* Do the triple product tB*D*Bprime: */2778 2768 TripleMultiply( &B[0][0],1,NUMVERTICES,1, 2779 2769 &DL_scalar,1,1,0, 2780 2770 &Bprime[0][0],1,NUMVERTICES,0, 2781 &Ke_gg_gaussian[0][0],0); 2782 2783 /* Add the Ke_gg_gaussian, and optionally Ke_gg_gaussian onto Ke_gg: */ 2784 for( i=0; i<numdof; i++) for(j=0;j<numdof;j++) Ke_gg[i][j]+=Ke_gg_gaussian[i][j]; 2771 &Ke_gg[0][0],0); 2772 2773 for(i=0;i<numdof;i++) for(j=0;j<numdof;j++) Ke->values[i*numdof+j]+=Ke_gg[i][j]; 2785 2774 } 2786 2775 2787 /*Add Ke_gg to global matrix Kgg: */ 2788 MatSetValues(Kgg,numdof,doflist,numdof,doflist,(const double*)Ke_gg,ADD_VALUES); 2789 2790 xfree((void**)&doflist); 2776 /*Clean up and return*/ 2791 2777 delete gauss; 2778 return Ke; 2779 } 2780 /*}}}*/ 2781 /*FUNCTION Penta::CreateKMatrixDiagnosticVertSurface {{{1*/ 2782 ElementMatrix* Penta::CreateKMatrixDiagnosticVertSurface(void){ 2783 2784 if (!IsOnSurface() || IsOnWater()) return NULL; 2785 2786 /*Call Tria function*/ 2787 Tria* tria=(Tria*)SpawnTria(3,4,5); //nodes 3,4 and 5 are on the surface 2788 ElementMatrix* Ke=tria->CreateKMatrixDiagnosticVertSurface(); 2789 delete tria->matice; delete tria; 2790 2791 /*clean up and return*/ 2792 return Ke; 2792 2793 } 2793 2794 /*}}}*/ -
issm/trunk/src/c/objects/Elements/Penta.h
r5876 r5877 141 141 ElementMatrix* CreateKMatrixDiagnosticStokesViscous(void); 142 142 ElementMatrix* CreateKMatrixDiagnosticStokesFriction(void); 143 void CreateKMatrixDiagnosticVert( Mat Kgg); 143 ElementMatrix* CreateKMatrixDiagnosticVert(void); 144 ElementMatrix* CreateKMatrixDiagnosticVertVolume(void); 145 ElementMatrix* CreateKMatrixDiagnosticVertSurface(void); 144 146 ElementMatrix* CreateKMatrixMelting(void); 145 147 ElementMatrix* CreateKMatrixPrognostic(void); -
issm/trunk/src/c/objects/Elements/Tria.cpp
r5875 r5877 3045 3045 } 3046 3046 /*}}}*/ 3047 /*FUNCTION Tria::CreateKMatrixDiagnostic SurfaceVert{{{1*/3048 void Tria::CreateKMatrixDiagnosticSurfaceVert(Mat Kgg){3047 /*FUNCTION Tria::CreateKMatrixDiagnosticVertSurface {{{1*/ 3048 ElementMatrix* Tria::CreateKMatrixDiagnosticVertSurface(void){ 3049 3049 3050 3050 /*Constants*/ … … 3059 3059 double L[3]; 3060 3060 GaussTria *gauss=NULL; 3061 3062 /* local element matrices: */ 3063 double Ke_gg[numdof][numdof]={0.0}; //local element stiffness matrix 3064 double Ke_gg_gaussian[numdof][numdof]; //stiffness matrix evaluated at the gaussian point. 3065 3066 /* Get node coordinates and dof list: */ 3061 double Ke_g[numdof][numdof]; //stiffness matrix evaluated at the gaussian point. 3062 3063 /*Initialize Element matrix and return if necessary*/ 3064 if(IsOnWater()) return NULL; 3065 ElementMatrix* Ke=this->NewElementMatrix(NoneApproximationEnum); 3066 3067 /*Retrieve all inputs and parameters*/ 3067 3068 GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES); 3068 GetDofList(&doflist,NoneApproximationEnum,GsetEnum);3069 3069 SurfaceNormal(&surface_normal[0],xyz_list); 3070 3070 … … 3085 3085 &DL_scalar,1,1,0, 3086 3086 L,1,3,0, 3087 &Ke_gg_gaussian[0][0],0); 3088 3089 for( i=0; i<numdof; i++) for(j=0;j<numdof;j++) Ke_gg[i][j]+=Ke_gg_gaussian[i][j]; 3090 } 3091 3092 MatSetValues(Kgg,numdof,doflist,numdof,doflist,(const double*)Ke_gg,ADD_VALUES); 3087 &Ke_g[0][0],0); 3088 3089 for(i=0;i<numdof;i++) for(j=0;j<numdof;j++) Ke->values[i*numdof+j]+=Ke_g[i][j]; 3090 } 3093 3091 3094 3092 /*Clean up and return*/ 3095 3093 delete gauss; 3096 xfree((void**)&doflist);3094 return Ke; 3097 3095 } 3098 3096 /*}}}*/ -
issm/trunk/src/c/objects/Elements/Tria.h
r5854 r5877 132 132 ElementMatrix* CreateKMatrixDiagnosticPattynFriction(void); 133 133 ElementMatrix* CreateKMatrixDiagnosticHutter(void); 134 void CreateKMatrixDiagnosticSurfaceVert(Mat Kgg);134 ElementMatrix* CreateKMatrixDiagnosticVertSurface(void); 135 135 ElementMatrix* CreateKMatrixMelting(void); 136 136 ElementMatrix* CreateKMatrixPrognostic(void);
Note:
See TracChangeset
for help on using the changeset viewer.