Changeset 5848
- Timestamp:
- 09/16/10 13:50:55 (15 years ago)
- Location:
- issm/trunk/src/c/objects/Elements
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk/src/c/objects/Elements/Penta.cpp
r5847 r5848 2189 2189 2190 2190 int approximation; 2191 ElementMatrix* Ke=NULL; 2192 2191 2193 inputs->GetParameterValue(&approximation,ApproximationEnum); 2192 2194 2193 2195 switch(approximation){ 2194 2196 case MacAyealApproximationEnum: 2195 CreateKMatrixDiagnosticMacAyeal2d( Kgg); 2197 Ke=CreateKMatrixDiagnosticMacAyeal2d(); 2198 if(Ke) Ke->AddToGlobal(Kgg,NULL,NULL); 2199 delete Ke; 2196 2200 break; 2197 2201 case PattynApproximationEnum: … … 2291 2295 } 2292 2296 /*FUNCTION Penta::CreateKMatrixDiagnosticMacAyeal2d{{{1*/ 2293 void Penta::CreateKMatrixDiagnosticMacAyeal2d(Mat Kgg){2297 ElementMatrix* Penta::CreateKMatrixDiagnosticMacAyeal2d(void){ 2294 2298 2295 /*If on water, skip stiffness: */ 2296 if(IsOnWater()) return; 2297 2298 /*Figure out if this pentaelem is collapsed. If so, then bailout, except if it is at the 2299 /*Figure out if this penta is collapsed. If so, then bailout, except if it is at the 2299 2300 bedrock, in which case we spawn a tria element using the 3 first grids, and use it to build 2300 2301 the stiffness matrix. */ 2301 2302 if (!IsOnBed()){ 2303 /*This element should be collapsed, but this element is not on the bedrock, therefore all its 2304 * dofs have already been frozen! Do nothing: */ 2305 return; 2306 } 2307 else{ 2308 2309 /*This element should be collapsed into a tria element at its base. Create this tria element, 2310 *and use its CreateKMatrix functionality to fill the global stiffness matrix: */ 2311 2312 /*Depth Averaging B*/ 2313 this->InputDepthAverageAtBase(RheologyBEnum,RheologyBbarEnum,MaterialsEnum); 2314 2315 /*Call Tria function*/ 2316 Tria* tria=(Tria*)SpawnTria(0,1,2); //grids 0, 1 and 2 make the new tria. 2317 tria->CreateKMatrix(Kgg,NULL,NULL); 2318 delete tria->matice; delete tria; 2319 2320 /*Delete B averaged*/ 2321 this->matice->inputs->DeleteInput(RheologyBbarEnum); 2322 2323 return; 2324 } 2302 if (!IsOnBed()) return NULL; 2303 2304 /*Depth Averaging B*/ 2305 this->InputDepthAverageAtBase(RheologyBEnum,RheologyBbarEnum,MaterialsEnum); 2306 2307 /*Call Tria function*/ 2308 Tria* tria=(Tria*)SpawnTria(0,1,2); //grids 0, 1 and 2 make the new tria. 2309 ElementMatrix* Ke=tria->CreateKMatrixDiagnosticMacAyeal(); 2310 delete tria->matice; delete tria; 2311 2312 /*Delete B averaged*/ 2313 this->matice->inputs->DeleteInput(RheologyBbarEnum); 2314 2315 /*clean up and return*/ 2316 return Ke; 2325 2317 } 2326 2318 /*}}}*/ … … 5408 5400 int i; 5409 5401 5410 const int numdofpervertex=1; 5411 const int numdof=numdofpervertex*NUMVERTICES; 5402 const int numdof=NDOF1*NUMVERTICES; 5412 5403 int* doflist=NULL; 5413 5404 double values[numdof]; -
issm/trunk/src/c/objects/Elements/Penta.h
r5847 r5848 127 127 void CreateKMatrixDiagnosticHoriz( Mat Kgg); 128 128 void CreateKMatrixDiagnosticHutter( Mat Kgg); 129 void CreateKMatrixDiagnosticMacAyeal2d(Mat Kgg);129 ElementMatrix* CreateKMatrixDiagnosticMacAyeal2d(void); 130 130 void CreateKMatrixDiagnosticMacAyeal3d(Mat Kgg); 131 131 void CreateKMatrixDiagnosticPattyn( Mat Kgg);
Note:
See TracChangeset
for help on using the changeset viewer.