Changeset 5877


Ignore:
Timestamp:
09/17/10 17:18:54 (15 years ago)
Author:
Mathieu Morlighem
Message:

ElementMatrix in Vert

Location:
issm/trunk/src/c/objects/Elements
Files:
4 edited

Legend:

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

    r5876 r5877  
    709709                        break;
    710710                case DiagnosticVertAnalysisEnum:
    711                         CreateKMatrixDiagnosticVert( Kgg);
     711                        Ke=CreateKMatrixDiagnosticVert();
     712                        if(Ke) Ke->AddToGlobal(Kgg,NULL,NULL);
     713                        delete Ke;
    712714                        break;
    713715                case BedSlopeXAnalysisEnum: case SurfaceSlopeXAnalysisEnum: case BedSlopeYAnalysisEnum: case SurfaceSlopeYAnalysisEnum:
     
    22812283ElementMatrix* Penta::CreateKMatrixDiagnosticMacAyeal3d(void){
    22822284
    2283 
    22842285        /*compute all stiffness matrices for this element*/
    22852286        ElementMatrix* Ke1=CreateKMatrixDiagnosticMacAyeal3dViscous();
     
    27192720/*}}}*/
    27202721/*FUNCTION Penta::CreateKMatrixDiagnosticVert {{{1*/
    2721 void Penta::CreateKMatrixDiagnosticVert( Mat Kgg){
     2722ElementMatrix* 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*/
     2736ElementMatrix* Penta::CreateKMatrixDiagnosticVertVolume(void){
    27222737
    27232738        /* local declarations */
    2724         int             i,j;
    2725 
    2726         /* node data: */
    27272739        const int    numdof=NDOF1*NUMVERTICES;
     2740        int             i,j,ig;
    27282741        double       xyz_list[NUMVERTICES][3];
    2729         int*         doflist=NULL;
    2730 
    2731         /* 3d gaussian points: */
    2732         int     ig;
    27332742        GaussPenta *gauss=NULL;
    2734 
    2735         /* matrices: */
    27362743        double  Ke_gg[numdof][numdof]={0.0};
    2737         double  Ke_gg_gaussian[numdof][numdof];
    27382744        double  Jdet;
    27392745        double  B[NDOF1][NUMVERTICES];
     
    27412747        double  DL_scalar;
    27422748
    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*/
    27602754        GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES);
    2761         GetDofList(&doflist,NoneApproximationEnum,GsetEnum);
    27622755
    27632756        /* Start  looping on the number of gaussian points: */
     
    27672760                gauss->GaussPoint(ig);
    27682761
    2769                 /*Get B and Bprime matrices: */
    27702762                GetBVert(&B[0][0], &xyz_list[0][0], gauss);
    27712763                GetBprimeVert(&Bprime[0][0], &xyz_list[0][0], gauss);
    27722764
    2773                 /* Get Jacobian determinant: */
    27742765                GetJacobianDeterminant(&Jdet, &xyz_list[0][0],gauss);
    27752766                DL_scalar=gauss->weight*Jdet;
    27762767
    2777                 /*  Do the triple product tB*D*Bprime: */
    27782768                TripleMultiply( &B[0][0],1,NUMVERTICES,1,
    27792769                                        &DL_scalar,1,1,0,
    27802770                                        &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];
    27852774        }
    27862775
    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*/
    27912777        delete gauss;
     2778        return Ke;
     2779}
     2780/*}}}*/
     2781/*FUNCTION Penta::CreateKMatrixDiagnosticVertSurface {{{1*/
     2782ElementMatrix* 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;
    27922793}
    27932794/*}}}*/
  • issm/trunk/src/c/objects/Elements/Penta.h

    r5876 r5877  
    141141                ElementMatrix* CreateKMatrixDiagnosticStokesViscous(void);
    142142                ElementMatrix* CreateKMatrixDiagnosticStokesFriction(void);
    143                 void      CreateKMatrixDiagnosticVert( Mat Kgg);
     143                ElementMatrix* CreateKMatrixDiagnosticVert(void);
     144                ElementMatrix* CreateKMatrixDiagnosticVertVolume(void);
     145                ElementMatrix* CreateKMatrixDiagnosticVertSurface(void);
    144146                ElementMatrix* CreateKMatrixMelting(void);
    145147                ElementMatrix* CreateKMatrixPrognostic(void);
  • issm/trunk/src/c/objects/Elements/Tria.cpp

    r5875 r5877  
    30453045}
    30463046/*}}}*/
    3047 /*FUNCTION Tria::CreateKMatrixDiagnosticSurfaceVert {{{1*/
    3048 void  Tria::CreateKMatrixDiagnosticSurfaceVert(Mat Kgg){
     3047/*FUNCTION Tria::CreateKMatrixDiagnosticVertSurface {{{1*/
     3048ElementMatrix* Tria::CreateKMatrixDiagnosticVertSurface(void){
    30493049
    30503050        /*Constants*/
     
    30593059        double    L[3];
    30603060        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*/
    30673068        GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES);
    3068         GetDofList(&doflist,NoneApproximationEnum,GsetEnum);
    30693069        SurfaceNormal(&surface_normal[0],xyz_list);
    30703070
     
    30853085                                        &DL_scalar,1,1,0,
    30863086                                        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        }
    30933091
    30943092        /*Clean up and return*/
    30953093        delete gauss;
    3096         xfree((void**)&doflist);
     3094        return Ke;
    30973095}
    30983096/*}}}*/
  • issm/trunk/src/c/objects/Elements/Tria.h

    r5854 r5877  
    132132                ElementMatrix* CreateKMatrixDiagnosticPattynFriction(void);
    133133                ElementMatrix* CreateKMatrixDiagnosticHutter(void);
    134                 void      CreateKMatrixDiagnosticSurfaceVert(Mat Kgg);
     134                ElementMatrix* CreateKMatrixDiagnosticVertSurface(void);
    135135                ElementMatrix* CreateKMatrixMelting(void);
    136136                ElementMatrix* CreateKMatrixPrognostic(void);
Note: See TracChangeset for help on using the changeset viewer.