Changeset 5924
- Timestamp:
- 09/21/10 14:33:07 (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
r5923 r5924 3002 3002 void Penta::CreatePVectorBalancedthickness( Vec pg){ 3003 3003 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; 3012 3005 3013 3006 /*Depth Averaging Vx and Vy*/ … … 3015 3008 this->InputDepthAverageAtBase(VyEnum,VyAverageEnum); 3016 3009 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(); 3020 3013 delete tria->matice; delete tria; 3014 pe->AddToGlobal(pg,NULL); 3015 delete pe; 3021 3016 3022 3017 /*Delete Vx and Vy averaged*/ … … 3024 3019 this->inputs->DeleteInput(VyAverageEnum); 3025 3020 3021 /*Clean up and return*/ 3026 3022 return; 3027 3023 } … … 3030 3026 void Penta::CreatePVectorBalancedvelocities( Vec pg){ 3031 3027 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; 3040 3029 3041 3030 /*Depth Averaging Vx and Vy*/ … … 3043 3032 this->InputDepthAverageAtBase(VyEnum,VyAverageEnum); 3044 3033 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(); 3048 3037 delete tria->matice; delete tria; 3038 pe->AddToGlobal(pg,NULL); 3039 delete pe; 3049 3040 3050 3041 /*Delete Vx and Vy averaged*/ … … 3052 3043 this->inputs->DeleteInput(VyAverageEnum); 3053 3044 3045 /*Clean up and return*/ 3054 3046 return; 3055 3047 } … … 3729 3721 void Penta::CreatePVectorPrognostic( Vec pg){ 3730 3722 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; 3739 3724 3740 3725 /*Depth Averaging Vx and Vy*/ … … 3742 3727 this->InputDepthAverageAtBase(VyEnum,VyAverageEnum); 3743 3728 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(); 3747 3732 delete tria->matice; delete tria; 3733 pe->AddToGlobal(pg,NULL); 3734 delete pe; 3748 3735 3749 3736 /*Delete Vx and Vy averaged*/ … … 3751 3738 this->inputs->DeleteInput(VyAverageEnum); 3752 3739 3740 /*Clean up and return*/ 3753 3741 return; 3754 3742 } -
issm/trunk/src/c/objects/Elements/Tria.cpp
r5922 r5924 752 752 case DiagnosticHorizAnalysisEnum: 753 753 pe=CreatePVectorDiagnosticMacAyeal(); 754 pe->AddToGlobal(pg,pf);755 delete pe;756 754 break; 757 755 case AdjointHorizAnalysisEnum: 758 756 pe=CreatePVectorAdjointHoriz(); 759 pe->AddToGlobal(pg,pf);760 delete pe;761 757 break; 762 758 case DiagnosticHutterAnalysisEnum: 763 759 pe=CreatePVectorDiagnosticHutter(); 764 pe->AddToGlobal(pg,pf);765 delete pe;766 760 break; 767 761 case BedSlopeXAnalysisEnum: case SurfaceSlopeXAnalysisEnum: case BedSlopeYAnalysisEnum: case SurfaceSlopeYAnalysisEnum: 768 762 pe=CreatePVectorSlope(); 769 pe->AddToGlobal(pg,pf);770 delete pe;771 763 break; 772 764 case PrognosticAnalysisEnum: 773 CreatePVectorPrognostic(pg);765 pe=CreatePVectorPrognostic(); 774 766 break; 775 767 case BalancedthicknessAnalysisEnum: 776 CreatePVectorBalancedthickness(pg);768 pe=CreatePVectorBalancedthickness(); 777 769 break; 778 770 case AdjointBalancedthicknessAnalysisEnum: 779 CreatePVectorAdjointBalancedthickness(pg);771 pe=CreatePVectorAdjointBalancedthickness(); 780 772 break; 781 773 case BalancedvelocitiesAnalysisEnum: 782 CreatePVectorBalancedvelocities( pg);774 pe=CreatePVectorBalancedvelocities(); 783 775 break; 784 776 default: … … 786 778 } 787 779 780 /*Add to global Vector*/ 781 if(pe){ 782 pe->AddToGlobal(pg,pf); 783 delete pe; 784 } 788 785 } 789 786 /*}}}*/ … … 3451 3448 /*}}}*/ 3452 3449 /*FUNCTION Tria::CreatePVectorBalancedthickness{{{1*/ 3453 void Tria::CreatePVectorBalancedthickness(Vec pg){3450 ElementVector* Tria::CreatePVectorBalancedthickness(void){ 3454 3451 3455 3452 switch(GetElementType()){ 3456 3453 case P1Enum: 3457 CreatePVectorBalancedthickness_CG( pg);3454 return CreatePVectorBalancedthickness_CG(); 3458 3455 break; 3459 3456 case P1DGEnum: 3460 CreatePVectorBalancedthickness_DG( pg); 3461 break; 3457 return CreatePVectorBalancedthickness_DG(); 3462 3458 default: 3463 3459 ISSMERROR("Element type %s not supported yet",EnumToString(GetElementType())); … … 3466 3462 /*}}}*/ 3467 3463 /*FUNCTION Tria::CreatePVectorBalancedthickness_CG{{{1*/ 3468 void Tria::CreatePVectorBalancedthickness_CG(Vec pg){3464 ElementVector* Tria::CreatePVectorBalancedthickness_CG(void){ 3469 3465 3470 3466 /*Constants*/ … … 3477 3473 double dhdt_g,melting_g,accumulation_g,Jdettria; 3478 3474 double L[NUMVERTICES]; 3479 double pe_g[NUMVERTICES] ={0.0};3480 3475 GaussTria* gauss=NULL; 3481 3476 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*/ 3483 3482 GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES); 3484 GetDofList(&doflist,NoneApproximationEnum,GsetEnum);3485 3486 /*retrieve inputs :*/3487 3483 Input* accumulation_input=inputs->GetInput(AccumulationRateEnum); ISSMASSERT(accumulation_input); 3488 3484 Input* melting_input=inputs->GetInput(MeltingRateEnum); ISSMASSERT(melting_input); … … 3502 3498 GetL(&L[0], &xyz_list[0][0], gauss,NDOF1); 3503 3499 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 } 3508 3502 3509 3503 /*Clean up and return*/ 3510 3504 delete gauss; 3511 xfree((void**)&doflist);3505 return pe; 3512 3506 } 3513 3507 /*}}}*/ 3514 3508 /*FUNCTION Tria::CreatePVectorBalancedthickness_DG {{{1*/ 3515 void Tria::CreatePVectorBalancedthickness_DG(Vec pg){3509 ElementVector* Tria::CreatePVectorBalancedthickness_DG(void){ 3516 3510 3517 3511 /*Constants*/ … … 3524 3518 double melting_g,accumulation_g,dhdt_g,Jdettria; 3525 3519 double L[NUMVERTICES]; 3526 double pe_g[NUMVERTICES] ={0.0};3527 3520 GaussTria* gauss=NULL; 3528 3521 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*/ 3530 3527 GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES); 3531 GetDofList(&doflist,NoneApproximationEnum,GsetEnum);3532 3533 /*Retrieve all inputs we will be needing: */3534 3528 Input* accumulation_input=inputs->GetInput(AccumulationRateEnum); ISSMASSERT(accumulation_input); 3535 3529 Input* melting_input=inputs->GetInput(MeltingRateEnum); ISSMASSERT(melting_input); … … 3549 3543 GetL(&L[0], &xyz_list[0][0], gauss,NDOF1); 3550 3544 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 } 3555 3547 3556 3548 /*Clean up and return*/ 3557 3549 delete gauss; 3558 xfree((void**)&doflist);3550 return pe; 3559 3551 } 3560 3552 /*}}}*/ 3561 3553 /*FUNCTION Tria::CreatePVectorBalancedvelocities {{{1*/ 3562 void Tria::CreatePVectorBalancedvelocities(Vec pg){3554 ElementVector* Tria::CreatePVectorBalancedvelocities(void){ 3563 3555 3564 3556 /*Constants*/ … … 3571 3563 double Jdettria,accumulation_g,melting_g; 3572 3564 double L[NUMVERTICES]; 3573 double pe_g[NUMVERTICES]={0.0};3574 3565 GaussTria* gauss=NULL; 3575 3566 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*/ 3576 3572 GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES); 3577 GetDofList(&doflist,NoneApproximationEnum,GsetEnum);3578 3579 /*Retrieve all inputs we will be needing: */3580 3573 Input* accumulation_input=inputs->GetInput(AccumulationRateEnum); ISSMASSERT(accumulation_input); 3581 3574 Input* melting_input=inputs->GetInput(MeltingRateEnum); ISSMASSERT(melting_input); … … 3593 3586 melting_input->GetParameterValue(&melting_g,gauss); 3594 3587 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 } 3600 3590 3601 3591 /*Clean up and return*/ 3602 3592 delete gauss; 3603 xfree((void**)&doflist);3593 return pe; 3604 3594 } 3605 3595 /*}}}*/ … … 3742 3732 /*}}}*/ 3743 3733 /*FUNCTION Tria::CreatePVectorAdjointBalancedthickness{{{1*/ 3744 void Tria::CreatePVectorAdjointBalancedthickness(Vec p_g){3734 ElementVector* Tria::CreatePVectorAdjointBalancedthickness(void){ 3745 3735 3746 3736 /*Constants*/ … … 3755 3745 double l1l2l3[3]; 3756 3746 double pe_g_gaussian[numdof]; 3757 double pe_g[numdof] ={0.0};3758 3747 GaussTria* gauss=NULL; 3759 3748 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*/ 3761 3754 GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES); 3762 GetDofList(&doflist,NoneApproximationEnum,GsetEnum);3763 3764 /*Retrieve all Inputs and parameters: */3765 3755 Input* thickness_input =inputs->GetInput(ThicknessEnum); ISSMASSERT(thickness_input); 3766 3756 Input* thicknessobs_input=inputs->GetInput(ThicknessObsEnum);ISSMASSERT(thicknessobs_input); … … 3780 3770 weights_input->GetParameterValue(&weight, gauss); 3781 3771 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 } 3791 3774 3792 3775 /*Clean up and return*/ 3793 3776 delete gauss; 3794 xfree((void**)&doflist);3777 return pe; 3795 3778 } 3796 3779 /*}}}*/ … … 4256 4239 /*}}}*/ 4257 4240 /*FUNCTION Tria::CreatePVectorPrognostic{{{1*/ 4258 void Tria::CreatePVectorPrognostic(Vec pg){4241 ElementVector* Tria::CreatePVectorPrognostic(void){ 4259 4242 4260 4243 switch(GetElementType()){ 4261 4244 case P1Enum: 4262 CreatePVectorPrognostic_CG( pg); 4263 break; 4245 return CreatePVectorPrognostic_CG(); 4264 4246 case P1DGEnum: 4265 CreatePVectorPrognostic_DG( pg); 4266 break; 4247 return CreatePVectorPrognostic_DG(); 4267 4248 default: 4268 4249 ISSMERROR("Element type %s not supported yet",EnumToString(GetElementType())); … … 4271 4252 /*}}}*/ 4272 4253 /*FUNCTION Tria::CreatePVectorPrognostic_CG {{{1*/ 4273 void Tria::CreatePVectorPrognostic_CG(Vec pg){4254 ElementVector* Tria::CreatePVectorPrognostic_CG(void){ 4274 4255 4275 4256 /* local declarations */ … … 4281 4262 int* doflist=NULL; 4282 4263 int numberofdofspernode=1; 4283 4284 /* gaussian points: */4285 4264 int ig; 4286 4265 GaussTria* gauss=NULL; 4287 4288 /* matrix */4289 double pe_g[NUMVERTICES]={0.0};4290 4266 double L[NUMVERTICES]; 4291 4267 double Jdettria; 4292 4293 /*input parameters for structural analysis (diagnostic): */4294 4268 double accumulation_g; 4295 4269 double melting_g; … … 4297 4271 double dt; 4298 4272 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); 4300 4279 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: */4307 4280 Input* accumulation_input=inputs->GetInput(AccumulationRateEnum); ISSMASSERT(accumulation_input); 4308 4281 Input* melting_input=inputs->GetInput(MeltingRateEnum); ISSMASSERT(melting_input); … … 4315 4288 gauss->GaussPoint(ig); 4316 4289 4317 /* Get Jacobian determinant: */4318 4290 GetJacobianDeterminant2d(&Jdettria, &xyz_list[0][0],gauss); 4319 4320 /*Get L matrix: */4321 4291 GetL(&L[0], &xyz_list[0][0], gauss,numberofdofspernode); 4322 4292 4323 /* Get accumulation, melting and thickness at gauss point */4324 4293 accumulation_input->GetParameterValue(&accumulation_g,gauss); 4325 4294 melting_input->GetParameterValue(&melting_g,gauss); 4326 4295 thickness_input->GetParameterValue(&thickness_g,gauss); 4327 4296 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 } 4335 4299 4336 4300 /*Clean up and return*/ 4337 4301 delete gauss; 4338 xfree((void**)&doflist); 4339 4302 return pe; 4340 4303 } 4341 4304 /*}}}*/ 4342 4305 /*FUNCTION Tria::CreatePVectorPrognostic_DG {{{1*/ 4343 void Tria::CreatePVectorPrognostic_DG(Vec pg){4306 ElementVector* Tria::CreatePVectorPrognostic_DG(void){ 4344 4307 4345 4308 /* local declarations */ 4309 const int numdof=NDOF1*NUMVERTICES; 4346 4310 int i,j; 4347 4348 /* node data: */ 4349 const int numdof=NDOF1*NUMVERTICES; 4311 int ig; 4350 4312 double xyz_list[NUMVERTICES][3]; 4351 4313 int* doflist=NULL; 4352 4314 int numberofdofspernode=1; 4353 4354 /* gaussian points: */4355 int ig;4356 4315 GaussTria* gauss=NULL; 4357 4358 /* matrix */4359 double pe_g[NUMVERTICES]={0.0};4360 4316 double L[NUMVERTICES]; 4361 4317 double Jdettria; 4362 4363 /*input parameters for structural analysis (diagnostic): */4364 4318 double accumulation_g; 4365 4319 double melting_g; … … 4367 4321 double dt; 4368 4322 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*/ 4370 4328 this->parameters->FindParam(&dt,DtEnum); 4371 4372 /* Get node coordinates and dof list: */4373 4329 GetVerticesCoordinates(&xyz_list[0][0], nodes, NUMVERTICES); 4374 GetDofList(&doflist,NoneApproximationEnum,GsetEnum);4375 4376 /*Retrieve all inputs we will be needing: */4377 4330 Input* accumulation_input=inputs->GetInput(AccumulationRateEnum); ISSMASSERT(accumulation_input); 4378 4331 Input* melting_input=inputs->GetInput(MeltingRateEnum); ISSMASSERT(melting_input); … … 4385 4338 gauss->GaussPoint(ig); 4386 4339 4387 /* Get Jacobian determinant: */4388 4340 GetJacobianDeterminant2d(&Jdettria, &xyz_list[0][0],gauss); 4389 4390 /*Get L matrix: */4391 4341 GetL(&L[0], &xyz_list[0][0], gauss,numberofdofspernode); 4392 4342 4393 /* Get accumulation, melting and thickness at gauss point */4394 4343 accumulation_input->GetParameterValue(&accumulation_g,gauss); 4395 4344 melting_input->GetParameterValue(&melting_g,gauss); 4396 4345 thickness_input->GetParameterValue(&thickness_g,gauss); 4397 4346 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 } 4405 4349 4406 4350 /*Clean up and return*/ 4407 4351 delete gauss; 4408 xfree((void**)&doflist);4352 return pe; 4409 4353 } 4410 4354 /*}}}*/ -
issm/trunk/src/c/objects/Elements/Tria.h
r5922 r5924 138 138 ElementMatrix* CreateKMatrixSlope(void); 139 139 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); 144 144 void CreatePVectorDiagnosticBaseVert(Vec pg); 145 145 ElementVector* CreatePVectorDiagnosticMacAyeal(void); 146 146 ElementVector* CreatePVectorAdjointHoriz(void); 147 147 void CreatePVectorAdjointStokes(Vec pg); 148 void CreatePVectorAdjointBalancedthickness(Vec pg);148 ElementVector* CreatePVectorAdjointBalancedthickness(void); 149 149 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); 153 153 ElementVector* CreatePVectorSlope(void); 154 154 void CreatePVectorThermalSheet( Vec pg);
Note:
See TracChangeset
for help on using the changeset viewer.