Changeset 5924


Ignore:
Timestamp:
09/21/10 14:33:07 (15 years ago)
Author:
Mathieu Morlighem
Message:

more ElementVector

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

Legend:

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

    r5923 r5924  
    30023002void Penta::CreatePVectorBalancedthickness( Vec pg){
    30033003
    3004         /*Collapsed formulation: */
    3005         Tria*  tria=NULL;
    3006 
    3007         /*If on water, skip: */
    3008         if(IsOnWater())return;
    3009 
    3010         /*Is this element on the bed? :*/
    3011         if(!IsOnBed())return;
     3004        if (!IsOnBed() || IsOnWater()) return;
    30123005
    30133006        /*Depth Averaging Vx and Vy*/
     
    30153008        this->InputDepthAverageAtBase(VyEnum,VyAverageEnum);
    30163009
    3017         /*Spawn Tria element from the base of the Penta: */
    3018         tria=(Tria*)SpawnTria(0,1,2); //grids 0, 1 and 2 make the new tria.
    3019         tria->CreatePVector(pg,NULL);
     3010        /*Call Tria function*/
     3011        Tria* tria=(Tria*)SpawnTria(0,1,2); //grids 0, 1 and 2 make the new tria.
     3012        ElementVector* pe=tria->CreatePVectorBalancedthickness();
    30203013        delete tria->matice; delete tria;
     3014        pe->AddToGlobal(pg,NULL);
     3015        delete pe;
    30213016
    30223017        /*Delete Vx and Vy averaged*/
     
    30243019        this->inputs->DeleteInput(VyAverageEnum);
    30253020
     3021        /*Clean up and return*/
    30263022        return;
    30273023}
     
    30303026void Penta::CreatePVectorBalancedvelocities( Vec pg){
    30313027
    3032         /*Collapsed formulation: */
    3033         Tria*  tria=NULL;
    3034 
    3035         /*If on water, skip: */
    3036         if(IsOnWater())return;
    3037 
    3038         /*Is this element on the bed? :*/
    3039         if(!IsOnBed())return;
     3028        if (!IsOnBed() || IsOnWater()) return;
    30403029
    30413030        /*Depth Averaging Vx and Vy*/
     
    30433032        this->InputDepthAverageAtBase(VyEnum,VyAverageEnum);
    30443033
    3045         /*Spawn Tria element from the base of the Penta: */
    3046         tria=(Tria*)SpawnTria(0,1,2); //grids 0, 1 and 2 make the new tria.
    3047         tria->CreatePVector(pg,NULL);
     3034        /*Call Tria function*/
     3035        Tria* tria=(Tria*)SpawnTria(0,1,2); //grids 0, 1 and 2 make the new tria.
     3036        ElementVector* pe=tria->CreatePVectorBalancedvelocities();
    30483037        delete tria->matice; delete tria;
     3038        pe->AddToGlobal(pg,NULL);
     3039        delete pe;
    30493040
    30503041        /*Delete Vx and Vy averaged*/
     
    30523043        this->inputs->DeleteInput(VyAverageEnum);
    30533044
     3045        /*Clean up and return*/
    30543046        return;
    30553047}
     
    37293721void Penta::CreatePVectorPrognostic( Vec pg){
    37303722
    3731         /*Collapsed formulation: */
    3732         Tria*  tria=NULL;
    3733        
    3734         /*If on water, skip: */
    3735         if(IsOnWater())return;
    3736 
    3737         /*Is this element on the bed? :*/
    3738         if(!IsOnBed())return;
     3723        if (!IsOnBed() || IsOnWater()) return;
    37393724
    37403725        /*Depth Averaging Vx and Vy*/
     
    37423727        this->InputDepthAverageAtBase(VyEnum,VyAverageEnum);
    37433728
    3744         /*Spawn Tria element from the base of the Penta: */
    3745         tria=(Tria*)SpawnTria(0,1,2); //grids 0, 1 and 2 make the new tria.
    3746         tria->CreatePVector(pg,NULL);
     3729        /*Call Tria function*/
     3730        Tria* tria=(Tria*)SpawnTria(0,1,2); //grids 0, 1 and 2 make the new tria.
     3731        ElementVector* pe=tria->CreatePVectorPrognostic();
    37473732        delete tria->matice; delete tria;
     3733        pe->AddToGlobal(pg,NULL);
     3734        delete pe;
    37483735
    37493736        /*Delete Vx and Vy averaged*/
     
    37513738        this->inputs->DeleteInput(VyAverageEnum);
    37523739
     3740        /*Clean up and return*/
    37533741        return;
    37543742}
  • issm/trunk/src/c/objects/Elements/Tria.cpp

    r5922 r5924  
    752752                case DiagnosticHorizAnalysisEnum:
    753753                        pe=CreatePVectorDiagnosticMacAyeal();
    754                         pe->AddToGlobal(pg,pf);
    755                         delete pe;
    756754                        break;
    757755                case AdjointHorizAnalysisEnum:
    758756                        pe=CreatePVectorAdjointHoriz();
    759                         pe->AddToGlobal(pg,pf);
    760                         delete pe;
    761757                        break;
    762758                case DiagnosticHutterAnalysisEnum:
    763759                        pe=CreatePVectorDiagnosticHutter();
    764                         pe->AddToGlobal(pg,pf);
    765                         delete pe;
    766760                        break;
    767761                case BedSlopeXAnalysisEnum: case SurfaceSlopeXAnalysisEnum: case BedSlopeYAnalysisEnum: case SurfaceSlopeYAnalysisEnum:
    768762                        pe=CreatePVectorSlope();
    769                         pe->AddToGlobal(pg,pf);
    770                         delete pe;
    771763                        break;
    772764                case PrognosticAnalysisEnum:
    773                         CreatePVectorPrognostic(pg);
     765                        pe=CreatePVectorPrognostic();
    774766                        break;
    775767                case BalancedthicknessAnalysisEnum:
    776                         CreatePVectorBalancedthickness(pg);
     768                        pe=CreatePVectorBalancedthickness();
    777769                        break;
    778770                case AdjointBalancedthicknessAnalysisEnum:
    779                         CreatePVectorAdjointBalancedthickness(pg);
     771                        pe=CreatePVectorAdjointBalancedthickness();
    780772                        break;
    781773                case BalancedvelocitiesAnalysisEnum:
    782                         CreatePVectorBalancedvelocities( pg);
     774                        pe=CreatePVectorBalancedvelocities();
    783775                        break;
    784776                default:
     
    786778        }
    787779
     780        /*Add to global Vector*/
     781        if(pe){
     782                pe->AddToGlobal(pg,pf);
     783                delete pe;
     784        }
    788785}
    789786/*}}}*/
     
    34513448/*}}}*/
    34523449/*FUNCTION Tria::CreatePVectorBalancedthickness{{{1*/
    3453 void  Tria::CreatePVectorBalancedthickness(Vec pg){
     3450ElementVector* Tria::CreatePVectorBalancedthickness(void){
    34543451
    34553452        switch(GetElementType()){
    34563453                case P1Enum:
    3457                         CreatePVectorBalancedthickness_CG( pg);
     3454                        return CreatePVectorBalancedthickness_CG();
    34583455                        break;
    34593456                case P1DGEnum:
    3460                         CreatePVectorBalancedthickness_DG( pg);
    3461                         break;
     3457                        return CreatePVectorBalancedthickness_DG();
    34623458                default:
    34633459                        ISSMERROR("Element type %s not supported yet",EnumToString(GetElementType()));
     
    34663462/*}}}*/
    34673463/*FUNCTION Tria::CreatePVectorBalancedthickness_CG{{{1*/
    3468 void  Tria::CreatePVectorBalancedthickness_CG(Vec pg ){
     3464ElementVector* Tria::CreatePVectorBalancedthickness_CG(void){
    34693465
    34703466        /*Constants*/
     
    34773473        double     dhdt_g,melting_g,accumulation_g,Jdettria;
    34783474        double     L[NUMVERTICES];
    3479         double     pe_g[NUMVERTICES]                        ={0.0};
    34803475        GaussTria* gauss=NULL;
    34813476
    3482         /* Get node coordinates and dof list: */
     3477        /*Initialize Element vector and return if necessary*/
     3478        if(IsOnWater()) return NULL;
     3479        ElementVector* pe=NewElementVector(nodes,NUMVERTICES,this->parameters);
     3480
     3481        /*Retrieve all inputs and parameters*/
    34833482        GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES);
    3484         GetDofList(&doflist,NoneApproximationEnum,GsetEnum);
    3485 
    3486         /*retrieve inputs :*/
    34873483        Input* accumulation_input=inputs->GetInput(AccumulationRateEnum); ISSMASSERT(accumulation_input);
    34883484        Input* melting_input=inputs->GetInput(MeltingRateEnum);           ISSMASSERT(melting_input);
     
    35023498                GetL(&L[0], &xyz_list[0][0], gauss,NDOF1);
    35033499
    3504                 for( i=0; i<numdof; i++) pe_g[i]+=Jdettria*gauss->weight*(accumulation_g-melting_g-dhdt_g)*L[i];
    3505         }
    3506 
    3507         VecSetValues(pg,numdof,doflist,(const double*)pe_g,ADD_VALUES);
     3500                for(i=0;i<numdof;i++) pe->values[i]+=Jdettria*gauss->weight*(accumulation_g-melting_g-dhdt_g)*L[i];
     3501        }
    35083502
    35093503        /*Clean up and return*/
    35103504        delete gauss;
    3511         xfree((void**)&doflist);
     3505        return pe;
    35123506}
    35133507/*}}}*/
    35143508/*FUNCTION Tria::CreatePVectorBalancedthickness_DG {{{1*/
    3515 void  Tria::CreatePVectorBalancedthickness_DG(Vec pg){
     3509ElementVector* Tria::CreatePVectorBalancedthickness_DG(void){
    35163510
    35173511        /*Constants*/
     
    35243518        double     melting_g,accumulation_g,dhdt_g,Jdettria;
    35253519        double     L[NUMVERTICES];
    3526         double     pe_g[NUMVERTICES]                       ={0.0};
    35273520        GaussTria* gauss=NULL;
    35283521
    3529         /* Get node coordinates and dof list: */
     3522        /*Initialize Element vector and return if necessary*/
     3523        if(IsOnWater()) return NULL;
     3524        ElementVector* pe=NewElementVector(nodes,NUMVERTICES,this->parameters);
     3525
     3526        /*Retrieve all inputs and parameters*/
    35303527        GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES);
    3531         GetDofList(&doflist,NoneApproximationEnum,GsetEnum);
    3532 
    3533         /*Retrieve all inputs we will be needing: */
    35343528        Input* accumulation_input=inputs->GetInput(AccumulationRateEnum); ISSMASSERT(accumulation_input);
    35353529        Input* melting_input=inputs->GetInput(MeltingRateEnum);           ISSMASSERT(melting_input);
     
    35493543                GetL(&L[0], &xyz_list[0][0], gauss,NDOF1);
    35503544
    3551                 for(i=0;i<numdof;i++) pe_g[i]+=Jdettria*gauss->weight*(accumulation_g-melting_g-dhdt_g)*L[i];
    3552         }
    3553 
    3554         VecSetValues(pg,numdof,doflist,(const double*)pe_g,ADD_VALUES);
     3545                for(i=0;i<numdof;i++) pe->values[i]+=Jdettria*gauss->weight*(accumulation_g-melting_g-dhdt_g)*L[i];
     3546        }
    35553547
    35563548        /*Clean up and return*/
    35573549        delete gauss;
    3558         xfree((void**)&doflist);
     3550        return pe;
    35593551}
    35603552/*}}}*/
    35613553/*FUNCTION Tria::CreatePVectorBalancedvelocities {{{1*/
    3562 void  Tria::CreatePVectorBalancedvelocities(Vec pg){
     3554ElementVector* Tria::CreatePVectorBalancedvelocities(void){
    35633555
    35643556        /*Constants*/
     
    35713563        double     Jdettria,accumulation_g,melting_g;
    35723564        double     L[NUMVERTICES];
    3573         double     pe_g[NUMVERTICES]={0.0};
    35743565        GaussTria* gauss=NULL;
    35753566
     3567        /*Initialize Element vector and return if necessary*/
     3568        if(IsOnWater()) return NULL;
     3569        ElementVector* pe=NewElementVector(nodes,NUMVERTICES,this->parameters);
     3570
     3571        /*Retrieve all inputs and parameters*/
    35763572        GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES);
    3577         GetDofList(&doflist,NoneApproximationEnum,GsetEnum);
    3578 
    3579         /*Retrieve all inputs we will be needing: */
    35803573        Input* accumulation_input=inputs->GetInput(AccumulationRateEnum); ISSMASSERT(accumulation_input);
    35813574        Input* melting_input=inputs->GetInput(MeltingRateEnum);           ISSMASSERT(melting_input);
     
    35933586                melting_input->GetParameterValue(&melting_g,gauss);
    35943587
    3595                 for( i=0; i<numdof; i++) pe_g[i]+=Jdettria*gauss->weight*(accumulation_g-melting_g)*L[i];
    3596 
    3597         }
    3598 
    3599         VecSetValues(pg,numdof,doflist,(const double*)pe_g,ADD_VALUES);
     3588                for(i=0;i<numdof;i++) pe->values[i]+=Jdettria*gauss->weight*(accumulation_g-melting_g)*L[i];
     3589        }
    36003590
    36013591        /*Clean up and return*/
    36023592        delete gauss;
    3603         xfree((void**)&doflist);
     3593        return pe;
    36043594}
    36053595/*}}}*/
     
    37423732/*}}}*/
    37433733/*FUNCTION Tria::CreatePVectorAdjointBalancedthickness{{{1*/
    3744 void Tria::CreatePVectorAdjointBalancedthickness(Vec p_g){
     3734ElementVector* Tria::CreatePVectorAdjointBalancedthickness(void){
    37453735
    37463736        /*Constants*/
     
    37553745        double      l1l2l3[3];
    37563746        double      pe_g_gaussian[numdof];
    3757         double      pe_g[numdof]                 ={0.0};
    37583747        GaussTria* gauss=NULL;
    37593748
    3760         /* Get node coordinates and dof list: */
     3749        /*Initialize Element vector and return if necessary*/
     3750        if(IsOnWater()) return NULL;
     3751        ElementVector* pe=NewElementVector(nodes,NUMVERTICES,this->parameters);
     3752
     3753        /*Retrieve all inputs and parameters*/
    37613754        GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES);
    3762         GetDofList(&doflist,NoneApproximationEnum,GsetEnum);
    3763 
    3764         /*Retrieve all Inputs and parameters: */
    37653755        Input* thickness_input   =inputs->GetInput(ThicknessEnum);   ISSMASSERT(thickness_input);
    37663756        Input* thicknessobs_input=inputs->GetInput(ThicknessObsEnum);ISSMASSERT(thicknessobs_input);
     
    37803770                weights_input->GetParameterValue(&weight, gauss);
    37813771
    3782                 for (i=0;i<NUMVERTICES;i++){
    3783                         pe_g_gaussian[i]=(thicknessobs-thickness)*weight*Jdet*gauss->weight*l1l2l3[i];
    3784                 }
    3785 
    3786                 /*Add pe_g_gaussian vector to pe_g: */
    3787                 for(i=0; i<numdof;i++) pe_g[i]+=pe_g_gaussian[i];
    3788         }
    3789 
    3790         VecSetValues(p_g,numdof,doflist,(const double*)pe_g,ADD_VALUES);
     3772                for(i=0;i<numdof;i++) pe->values[i]+=(thicknessobs-thickness)*weight*Jdet*gauss->weight*l1l2l3[i];
     3773        }
    37913774
    37923775        /*Clean up and return*/
    37933776        delete gauss;
    3794         xfree((void**)&doflist);
     3777        return pe;
    37953778}
    37963779/*}}}*/
     
    42564239/*}}}*/
    42574240/*FUNCTION Tria::CreatePVectorPrognostic{{{1*/
    4258 void  Tria::CreatePVectorPrognostic(Vec pg){
     4241ElementVector* Tria::CreatePVectorPrognostic(void){
    42594242
    42604243        switch(GetElementType()){
    42614244                case P1Enum:
    4262                         CreatePVectorPrognostic_CG( pg);
    4263                         break;
     4245                        return CreatePVectorPrognostic_CG();
    42644246                case P1DGEnum:
    4265                         CreatePVectorPrognostic_DG( pg);
    4266                         break;
     4247                        return CreatePVectorPrognostic_DG();
    42674248                default:
    42684249                        ISSMERROR("Element type %s not supported yet",EnumToString(GetElementType()));
     
    42714252/*}}}*/
    42724253/*FUNCTION Tria::CreatePVectorPrognostic_CG {{{1*/
    4273 void  Tria::CreatePVectorPrognostic_CG(Vec pg){
     4254ElementVector* Tria::CreatePVectorPrognostic_CG(void){
    42744255
    42754256        /* local declarations */
     
    42814262        int*         doflist=NULL;
    42824263        int          numberofdofspernode=1;
    4283 
    4284         /* gaussian points: */
    42854264        int     ig;
    42864265        GaussTria* gauss=NULL;
    4287 
    4288         /* matrix */
    4289         double pe_g[NUMVERTICES]={0.0};
    42904266        double L[NUMVERTICES];
    42914267        double Jdettria;
    4292 
    4293         /*input parameters for structural analysis (diagnostic): */
    42944268        double  accumulation_g;
    42954269        double  melting_g;
     
    42974271        double  dt;
    42984272
    4299         /*retrieve some parameters: */
     4273        /*Initialize Element vector and return if necessary*/
     4274        if(IsOnWater()) return NULL;
     4275        ElementVector* pe=NewElementVector(nodes,NUMVERTICES,this->parameters);
     4276
     4277        /*Retrieve all inputs and parameters*/
     4278        GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES);
    43004279        this->parameters->FindParam(&dt,DtEnum);
    4301 
    4302         /* Get node coordinates and dof list: */
    4303         GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES);
    4304         GetDofList(&doflist,NoneApproximationEnum,GsetEnum);
    4305 
    4306         /*Retrieve all inputs we will be needing: */
    43074280        Input* accumulation_input=inputs->GetInput(AccumulationRateEnum); ISSMASSERT(accumulation_input);
    43084281        Input* melting_input=inputs->GetInput(MeltingRateEnum);           ISSMASSERT(melting_input);
     
    43154288                gauss->GaussPoint(ig);
    43164289
    4317                 /* Get Jacobian determinant: */
    43184290                GetJacobianDeterminant2d(&Jdettria, &xyz_list[0][0],gauss);
    4319 
    4320                 /*Get L matrix: */
    43214291                GetL(&L[0], &xyz_list[0][0], gauss,numberofdofspernode);
    43224292
    4323                 /* Get accumulation, melting and thickness at gauss point */
    43244293                accumulation_input->GetParameterValue(&accumulation_g,gauss);
    43254294                melting_input->GetParameterValue(&melting_g,gauss);
    43264295                thickness_input->GetParameterValue(&thickness_g,gauss);
    43274296
    4328                 /* Add value into pe_g: */
    4329                 for( i=0; i<numdof; i++) pe_g[i]+=Jdettria*gauss->weight*(thickness_g+dt*(accumulation_g-melting_g))*L[i];
    4330 
    4331         }
    4332 
    4333         /*Add pe_g to global matrix Kgg: */
    4334         VecSetValues(pg,numdof,doflist,(const double*)pe_g,ADD_VALUES);
     4297                for(i=0;i<numdof;i++) pe->values[i]+=Jdettria*gauss->weight*(thickness_g+dt*(accumulation_g-melting_g))*L[i];
     4298        }
    43354299
    43364300        /*Clean up and return*/
    43374301        delete gauss;
    4338         xfree((void**)&doflist);
    4339 
     4302        return pe;
    43404303}
    43414304/*}}}*/
    43424305/*FUNCTION Tria::CreatePVectorPrognostic_DG {{{1*/
    4343 void  Tria::CreatePVectorPrognostic_DG(Vec pg){
     4306ElementVector* Tria::CreatePVectorPrognostic_DG(void){
    43444307
    43454308        /* local declarations */
     4309        const int    numdof=NDOF1*NUMVERTICES;
    43464310        int             i,j;
    4347 
    4348         /* node data: */
    4349         const int    numdof=NDOF1*NUMVERTICES;
     4311        int     ig;
    43504312        double       xyz_list[NUMVERTICES][3];
    43514313        int*         doflist=NULL;
    43524314        int          numberofdofspernode=1;
    4353 
    4354         /* gaussian points: */
    4355         int     ig;
    43564315        GaussTria* gauss=NULL;
    4357 
    4358         /* matrix */
    4359         double pe_g[NUMVERTICES]={0.0};
    43604316        double L[NUMVERTICES];
    43614317        double Jdettria;
    4362 
    4363         /*input parameters for structural analysis (diagnostic): */
    43644318        double  accumulation_g;
    43654319        double  melting_g;
     
    43674321        double  dt;
    43684322
    4369         /*retrieve some parameters: */
     4323        /*Initialize Element vector and return if necessary*/
     4324        if(IsOnWater()) return NULL;
     4325        ElementVector* pe=NewElementVector(nodes,NUMVERTICES,this->parameters);
     4326
     4327        /*Retrieve all inputs and parameters*/
    43704328        this->parameters->FindParam(&dt,DtEnum);
    4371 
    4372         /* Get node coordinates and dof list: */
    43734329        GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES);
    4374         GetDofList(&doflist,NoneApproximationEnum,GsetEnum);
    4375 
    4376         /*Retrieve all inputs we will be needing: */
    43774330        Input* accumulation_input=inputs->GetInput(AccumulationRateEnum); ISSMASSERT(accumulation_input);
    43784331        Input* melting_input=inputs->GetInput(MeltingRateEnum);           ISSMASSERT(melting_input);
     
    43854338                gauss->GaussPoint(ig);
    43864339
    4387                 /* Get Jacobian determinant: */
    43884340                GetJacobianDeterminant2d(&Jdettria, &xyz_list[0][0],gauss);
    4389 
    4390                 /*Get L matrix: */
    43914341                GetL(&L[0], &xyz_list[0][0], gauss,numberofdofspernode);
    43924342
    4393                 /* Get accumulation, melting and thickness at gauss point */
    43944343                accumulation_input->GetParameterValue(&accumulation_g,gauss);
    43954344                melting_input->GetParameterValue(&melting_g,gauss);
    43964345                thickness_input->GetParameterValue(&thickness_g,gauss);
    43974346
    4398                 /* Add value into pe_g: */
    4399                 for( i=0; i<numdof; i++) pe_g[i]+=Jdettria*gauss->weight*(thickness_g+dt*(accumulation_g-melting_g))*L[i];
    4400 
    4401         }
    4402 
    4403         /*Add pe_g to global matrix Kgg: */
    4404         VecSetValues(pg,numdof,doflist,(const double*)pe_g,ADD_VALUES);
     4347                for(i=0;i<numdof;i++) pe->values[i]+=Jdettria*gauss->weight*(thickness_g+dt*(accumulation_g-melting_g))*L[i];
     4348        }
    44054349
    44064350        /*Clean up and return*/
    44074351        delete gauss;
    4408         xfree((void**)&doflist);
     4352        return pe;
    44094353}
    44104354/*}}}*/
  • issm/trunk/src/c/objects/Elements/Tria.h

    r5922 r5924  
    138138                ElementMatrix* CreateKMatrixSlope(void);
    139139                ElementMatrix* CreateKMatrixThermal(void);
    140                 void      CreatePVectorBalancedthickness(Vec pg);
    141                 void      CreatePVectorBalancedthickness_DG(Vec pg);
    142                 void      CreatePVectorBalancedthickness_CG(Vec pg);
    143                 void      CreatePVectorBalancedvelocities(Vec pg);
     140                ElementVector* CreatePVectorBalancedthickness(void);
     141                ElementVector* CreatePVectorBalancedthickness_DG(void);
     142                ElementVector* CreatePVectorBalancedthickness_CG(void);
     143                ElementVector* CreatePVectorBalancedvelocities(void);
    144144                void      CreatePVectorDiagnosticBaseVert(Vec pg);
    145145                ElementVector* CreatePVectorDiagnosticMacAyeal(void);
    146146                ElementVector* CreatePVectorAdjointHoriz(void);
    147147                void      CreatePVectorAdjointStokes(Vec pg);
    148                 void      CreatePVectorAdjointBalancedthickness(Vec pg);
     148                ElementVector* CreatePVectorAdjointBalancedthickness(void);
    149149                ElementVector* CreatePVectorDiagnosticHutter(void);
    150                 void      CreatePVectorPrognostic(Vec pg);
    151                 void      CreatePVectorPrognostic_CG(Vec pg);
    152                 void      CreatePVectorPrognostic_DG(Vec pg);
     150                ElementVector* CreatePVectorPrognostic(void);
     151                ElementVector* CreatePVectorPrognostic_CG(void);
     152                ElementVector* CreatePVectorPrognostic_DG(void);
    153153                ElementVector* CreatePVectorSlope(void);
    154154                void      CreatePVectorThermalSheet( Vec pg);
Note: See TracChangeset for help on using the changeset viewer.