Changeset 352


Ignore:
Timestamp:
05/12/09 11:07:01 (16 years ago)
Author:
seroussi
Message:

KMatrixBaseVert not needed

Location:
issm/trunk/src/c/objects
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • issm/trunk/src/c/objects/Tria.cpp

    r336 r352  
    20542054
    20552055#undef __FUNCT__
    2056 #define __FUNCT__ "Tria::CreateKMatrixDiagnosticBaseVert"
    2057 void  Tria::CreateKMatrixDiagnosticBaseVert(Mat Kgg,void* vinputs,int analysis_type){
    2058 
    2059         int i,j;
    2060 
    2061         /* node data: */
    2062         const int    numgrids=3;
    2063         const int    NDOF1=1;
    2064         const int    numdofs=NDOF1*numgrids;
    2065         double       xyz_list[numgrids][3];
    2066         int          doflist[numdofs];
    2067         int          numberofdofspernode;
    2068        
    2069         /* gaussian points: */
    2070         int     num_gauss,ig;
    2071         double* first_gauss_area_coord  =  NULL;
    2072         double* second_gauss_area_coord =  NULL;
    2073         double* third_gauss_area_coord  =  NULL;
    2074         double* gauss_weights           =  NULL;
    2075         double  gauss_weight;
    2076         double  gauss_l1l2l3[3];
    2077 
    2078 
    2079         /* matrices: */
    2080         double L[1][3];
    2081         double DL_scalar;
    2082 
    2083         double Ke_gg[3][3]; //stiffness matrix
    2084         double Ke_gg_gaussian[3][3]; //stiffness matrix evaluated at the gaussian point.
    2085         double Jdet;
    2086 
    2087         ParameterInputs* inputs=NULL;
    2088 
    2089         /*recover pointers: */
    2090         inputs=(ParameterInputs*)vinputs;
    2091 
    2092         /* Get node coordinates and dof list: */
    2093         GetElementNodeData( &xyz_list[0][0], nodes, numgrids);
    2094         GetDofList(&doflist[0],&numberofdofspernode);
    2095 
    2096         /* Set Ke_gg to 0: */
    2097         for(i=0;i<numdofs;i++) for(j=0;j<numdofs;j++) Ke_gg[i][j]=0.0;
    2098 
    2099         /* Get gaussian points and weights (make this a statically initialized list of points? fstd): */
    2100         GaussTria( &num_gauss, &first_gauss_area_coord, &second_gauss_area_coord, &third_gauss_area_coord, &gauss_weights, 2);
    2101 
    2102         /* Start  looping on the number of gaussian points: */
    2103         for (ig=0; ig<num_gauss; ig++){
    2104                 /*Pick up the gaussian point: */
    2105                 gauss_weight=*(gauss_weights+ig);
    2106                 gauss_l1l2l3[0]=*(first_gauss_area_coord+ig);
    2107                 gauss_l1l2l3[1]=*(second_gauss_area_coord+ig);
    2108                 gauss_l1l2l3[2]=*(third_gauss_area_coord+ig);
    2109 
    2110                 /* Get Jacobian determinant: */
    2111                 GetJacobianDeterminant(&Jdet, &xyz_list[0][0],gauss_l1l2l3);
    2112        
    2113                 //Get L matrix if viscous basal drag present:
    2114                 GetL(&L[0][0], &xyz_list[0][0], gauss_l1l2l3,NDOF1);
    2115 
    2116                 DL_scalar=gauss_weight*Jdet;
    2117                
    2118                 /*  Do the triple producte tL*D*L: */
    2119                 TripleMultiply( &L[0][0],1,3,1,
    2120                                         &DL_scalar,1,1,0,
    2121                                         &L[0][0],1,3,0,
    2122                                         &Ke_gg_gaussian[0][0],0);
    2123 
    2124                
    2125                 /* Add the Ke_gg_gaussian, and optionally Ke_gg_drag_gaussian onto Ke_gg: */
    2126                 for( i=0; i<numdofs; i++) for(j=0;j<numdofs;j++) Ke_gg[i][j]+=Ke_gg_gaussian[i][j];
    2127         } //for (ig=0; ig<num_gauss; ig++
    2128 
    2129         /*Add Ke_gg to global matrix Kgg: */
    2130         MatSetValues(Kgg,numdofs,doflist,numdofs,doflist,(const double*)Ke_gg,ADD_VALUES);
    2131                
    2132         cleanup_and_return:
    2133         xfree((void**)&first_gauss_area_coord);
    2134         xfree((void**)&second_gauss_area_coord);
    2135         xfree((void**)&third_gauss_area_coord);
    2136         xfree((void**)&gauss_weights);
    2137 }
    2138 
    2139 #undef __FUNCT__
    21402056#define __FUNCT__ "Tria::CreateKMatrixDiagnosticSurfaceVert"
    21412057void  Tria::CreateKMatrixDiagnosticSurfaceVert(Mat Kgg,void* vinputs,int analysis_type){
  • issm/trunk/src/c/objects/Tria.h

    r301 r352  
    7272                void  CreateKMatrixDiagnosticHoriz(Mat Kgg,void* inputs,int analysis_type);
    7373                void  CreateKMatrixDiagnosticHorizFriction(Mat Kgg,void* inputs,int analysis_type);
    74                 void  CreateKMatrixDiagnosticBaseVert(Mat Kgg,void* inputs,int analysis_type);
    7574                void  CreateKMatrixDiagnosticSurfaceVert(Mat Kgg,void* inputs,int analysis_type);
    7675                void  GetParameterValue(double* pp, double* plist, double* gauss_l1l2l3);
Note: See TracChangeset for help on using the changeset viewer.