Changeset 5852
- Timestamp:
- 09/16/10 15:34:45 (15 years ago)
- Location:
- issm/trunk/src/c/objects/Elements
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk/src/c/objects/Elements/Penta.cpp
r5851 r5852 2141 2141 2142 2142 Tria* tria=(Tria*)SpawnTria(0,1,2); //grids 0, 1 and 2 make the new tria. 2143 tria->CreateKMatrixCouplingMacAyealPattynFriction(Kgg); 2143 ElementMatrix* Ke=tria->CreateKMatrixCouplingMacAyealPattynFriction(); 2144 if(Ke) Ke->AddToGlobal(Kgg,NULL,NULL); 2145 delete Ke; 2144 2146 delete tria->matice; delete tria; 2145 2147 } -
issm/trunk/src/c/objects/Elements/Tria.cpp
r5849 r5852 2860 2860 /*}}}*/ 2861 2861 /*FUNCTION Tria::CreateKMatrixCouplingMacAyealPattynFriction {{{1*/ 2862 void Tria::CreateKMatrixCouplingMacAyealPattynFriction(Mat Kgg){2863 2864 /* local declarations */2862 ElementMatrix* Tria::CreateKMatrixCouplingMacAyealPattynFriction(void){ 2863 2864 const int numdof = 2 *NUMVERTICES; 2865 2865 int i,j; 2866 int ig; 2866 2867 int analysis_type; 2867 2868 /* node data: */2869 const int numdof = 2 *NUMVERTICES;2870 2868 double xyz_list[NUMVERTICES][3]; 2871 int* doflistm=NULL;2872 int* doflistp=NULL;2873 2869 int numberofdofspernode=2; 2874 2875 /* gaussian points: */2876 int ig;2877 2870 GaussTria *gauss=NULL; 2878 2879 /* matrices: */2880 2871 double L[2][numdof]; 2881 2872 double DL[2][2]={{ 0,0 },{0,0}}; //for basal drag 2882 2873 double DL_scalar; 2883 2884 /* local element matrices: */2885 2874 double Ke_gg[numdof][numdof]={0.0}; 2886 2875 double Ke_gg_gaussian[numdof][numdof]; //stiffness matrix contribution from drag 2887 2888 2876 double Jdet; 2889 2890 /*slope: */2891 2877 double slope[2]={0.0,0.0}; 2892 2878 double slope_magnitude; 2893 2894 /*friction: */2895 2879 Friction *friction = NULL; 2896 2880 double alpha2; 2897 2898 2881 double MAXSLOPE=.06; // 6 % 2899 2882 double MOUNTAINKEXPONENT=10; 2900 2901 /*inputs: */2902 2883 int drag_type; 2903 2884 2904 /*retrive parameters: */ 2885 /*Initialize Element matrix and return if necessary*/ 2886 if(IsOnWater() || IsOnShelf()) return NULL; 2887 ElementMatrix* Ke1=this->NewElementMatrix(MacAyealApproximationEnum); 2888 ElementMatrix* Ke2=this->NewElementMatrix(PattynApproximationEnum); 2889 ElementMatrix* Ke=new ElementMatrix(Ke1,Ke2); 2890 delete Ke1; delete Ke2; 2891 2892 /*retrieve inputs :*/ 2893 GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES); 2905 2894 parameters->FindParam(&analysis_type,AnalysisTypeEnum); 2906 2907 /*retrieve inputs :*/2908 2895 inputs->GetParameterValue(&drag_type,DragTypeEnum); 2909 2896 Input* surface_input=inputs->GetInput(SurfaceEnum); ISSMASSERT(surface_input); … … 2911 2898 Input* vy_input=inputs->GetInput(VyEnum); ISSMASSERT(vy_input); 2912 2899 Input* vz_input=inputs->GetInput(VzEnum); ISSMASSERT(vz_input); 2913 2914 /* Get node coordinates and dof list: */2915 GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES);2916 GetDofList(&doflistm,MacAyealApproximationEnum,GsetEnum);2917 GetDofList(&doflistp,PattynApproximationEnum,GsetEnum);2918 2919 if (IsOnShelf()){2920 /*no friction, do nothing*/2921 return;2922 }2923 2900 2924 2901 /*build friction object, used later on: */ … … 2962 2939 &Ke_gg_gaussian[0][0],0); 2963 2940 2964 for( i=0; i<numdof; i++) for(j=0;j<numdof;j++) Ke_gg[i][j]+=Ke_gg_gaussian[i][j]; 2965 2966 } 2967 2968 /*Add Ke_gg to global matrix Kgg: */ 2969 MatSetValues(Kgg,numdof,doflistm,numdof,doflistp,(const double*)Ke_gg,ADD_VALUES); 2970 MatSetValues(Kgg,numdof,doflistp,numdof,doflistm,(const double*)Ke_gg,ADD_VALUES); 2941 for(i=0;i<numdof;i++) for(j=0;j<numdof;j++) Ke_gg[i][j]+=Ke_gg_gaussian[i][j]; 2942 } 2943 2944 2945 for(i=0;i<numdof;i++) for(j=0;j<numdof;j++) Ke->values[i*numdof+(numdof+j)]+=Ke_gg[i][j]; 2946 for(i=0;i<numdof;i++) for(j=0;j<numdof;j++) Ke->values[(i+numdof)*numdof+j]+=Ke_gg[i][j]; 2971 2947 2972 2948 /*Clean up and return*/ 2973 2949 delete gauss; 2974 2950 delete friction; 2975 xfree((void**)&doflistm); 2976 xfree((void**)&doflistp); 2951 return Ke; 2977 2952 } 2978 2953 /*}}}*/ -
issm/trunk/src/c/objects/Elements/Tria.h
r5849 r5852 129 129 ElementMatrix* CreateKMatrixDiagnosticMacAyealFriction(void); 130 130 void CreateKMatrixDiagnosticMacAyealFriction(Mat Kgg,Mat Kff, Mat Kfs); 131 void CreateKMatrixCouplingMacAyealPattynFriction(Mat Kgg);131 ElementMatrix* CreateKMatrixCouplingMacAyealPattynFriction(void); 132 132 ElementMatrix* CreateKMatrixDiagnosticPattynFriction(void); 133 133 ElementMatrix* CreateKMatrixDiagnosticHutter(void);
Note:
See TracChangeset
for help on using the changeset viewer.