Changeset 5852


Ignore:
Timestamp:
09/16/10 15:34:45 (15 years ago)
Author:
seroussi
Message:

Added ElementMatrix for Coupling

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

Legend:

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

    r5851 r5852  
    21412141
    21422142        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;
    21442146        delete tria->matice; delete tria;
    21452147}
  • issm/trunk/src/c/objects/Elements/Tria.cpp

    r5849 r5852  
    28602860/*}}}*/
    28612861/*FUNCTION Tria::CreateKMatrixCouplingMacAyealPattynFriction {{{1*/
    2862 void  Tria::CreateKMatrixCouplingMacAyealPattynFriction(Mat Kgg){
    2863 
    2864         /* local declarations */
     2862ElementMatrix* Tria::CreateKMatrixCouplingMacAyealPattynFriction(void){
     2863
     2864        const int numdof   = 2 *NUMVERTICES;
    28652865        int       i,j;
     2866        int     ig;
    28662867        int       analysis_type;
    2867 
    2868         /* node data: */
    2869         const int numdof   = 2 *NUMVERTICES;
    28702868        double    xyz_list[NUMVERTICES][3];
    2871         int*      doflistm=NULL;
    2872         int*      doflistp=NULL;
    28732869        int       numberofdofspernode=2;
    2874 
    2875         /* gaussian points: */
    2876         int     ig;
    28772870        GaussTria *gauss=NULL;
    2878 
    2879         /* matrices: */
    28802871        double L[2][numdof];
    28812872        double DL[2][2]={{ 0,0 },{0,0}}; //for basal drag
    28822873        double DL_scalar;
    2883 
    2884         /* local element matrices: */
    28852874        double Ke_gg[numdof][numdof]={0.0};
    28862875        double Ke_gg_gaussian[numdof][numdof]; //stiffness matrix contribution from drag
    2887        
    28882876        double Jdet;
    2889        
    2890         /*slope: */
    28912877        double  slope[2]={0.0,0.0};
    28922878        double  slope_magnitude;
    2893 
    2894         /*friction: */
    28952879        Friction *friction = NULL;
    28962880        double    alpha2;
    2897 
    28982881        double MAXSLOPE=.06; // 6 %
    28992882        double MOUNTAINKEXPONENT=10;
    2900 
    2901         /*inputs: */
    29022883        int  drag_type;
    29032884
    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);
    29052894        parameters->FindParam(&analysis_type,AnalysisTypeEnum);
    2906 
    2907         /*retrieve inputs :*/
    29082895        inputs->GetParameterValue(&drag_type,DragTypeEnum);
    29092896        Input* surface_input=inputs->GetInput(SurfaceEnum); ISSMASSERT(surface_input);
     
    29112898        Input* vy_input=inputs->GetInput(VyEnum);           ISSMASSERT(vy_input);
    29122899        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         }
    29232900
    29242901        /*build friction object, used later on: */
     
    29622939                                        &Ke_gg_gaussian[0][0],0);
    29632940
    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];
    29712947
    29722948        /*Clean up and return*/
    29732949        delete gauss;
    29742950        delete friction;
    2975         xfree((void**)&doflistm);
    2976         xfree((void**)&doflistp);
     2951        return Ke;
    29772952}       
    29782953/*}}}*/
  • issm/trunk/src/c/objects/Elements/Tria.h

    r5849 r5852  
    129129                ElementMatrix* CreateKMatrixDiagnosticMacAyealFriction(void);
    130130                void    CreateKMatrixDiagnosticMacAyealFriction(Mat Kgg,Mat Kff, Mat Kfs);
    131                 void      CreateKMatrixCouplingMacAyealPattynFriction(Mat Kgg);
     131                ElementMatrix* CreateKMatrixCouplingMacAyealPattynFriction(void);
    132132                ElementMatrix* CreateKMatrixDiagnosticPattynFriction(void);
    133133                ElementMatrix* CreateKMatrixDiagnosticHutter(void);
Note: See TracChangeset for help on using the changeset viewer.