Changeset 25442
- Timestamp:
- 08/21/20 11:35:36 (5 years ago)
- Location:
- issm/trunk-jpl/src/c/classes
- Files:
-
- 8 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/c/classes/Elements/Element.h
r25406 r25442 303 303 virtual Gauss* NewGauss(int point1,IssmDouble fraction1,IssmDouble fraction2,bool mainlyfloating,int order)=0; 304 304 virtual Gauss* NewGauss(int point1,IssmDouble fraction1,IssmDouble fraction2,int order)=0; 305 virtual Gauss* NewGauss(IssmDouble fraction1,IssmDouble fraction2,int order)=0; 305 306 virtual Gauss* NewGaussBase(int order)=0; 306 307 virtual Gauss* NewGaussLine(int vertex1,int vertex2,int order)=0; -
issm/trunk-jpl/src/c/classes/Elements/Penta.h
r25406 r25442 136 136 Gauss* NewGauss(int point1,IssmDouble fraction1,IssmDouble fraction2,bool mainlyfloating,int order); 137 137 Gauss* NewGauss(int point1,IssmDouble fraction1,IssmDouble fraction2,int order){_error_("not implemented yet");}; 138 Gauss* NewGauss(IssmDouble fraction1,IssmDouble fraction2,int order){_error_("not implemented yet");}; 138 139 Gauss* NewGaussBase(int order); 139 140 Gauss* NewGaussLine(int vertex1,int vertex2,int order); -
issm/trunk-jpl/src/c/classes/Elements/Seg.h
r25406 r25442 108 108 Gauss* NewGauss(int point1,IssmDouble fraction1,IssmDouble fraction2,bool mainlyfloating,int order){_error_("not implemented yet");}; 109 109 Gauss* NewGauss(int point1,IssmDouble fraction1,IssmDouble fraction2,int order){_error_("not implemented yet");}; 110 Gauss* NewGauss(IssmDouble fraction1,IssmDouble fraction2,int order){_error_("not implemented yet");}; 110 111 Gauss* NewGaussBase(int order){_error_("not implemented yet");}; 111 112 Gauss* NewGaussLine(int vertex1,int vertex2,int order){_error_("not implemented yet");}; -
issm/trunk-jpl/src/c/classes/Elements/Tetra.h
r25406 r25442 114 114 Gauss* NewGauss(int point1,IssmDouble fraction1,IssmDouble fraction2,bool mainlyfloating,int order){_error_("not implemented yet");}; 115 115 Gauss* NewGauss(int point1,IssmDouble fraction1,IssmDouble fraction2,int order){_error_("not implemented yet");}; 116 Gauss* NewGauss(IssmDouble fraction1,IssmDouble fraction2,int order){_error_("not implemented yet");}; 116 117 Gauss* NewGaussBase(int order); 117 118 Gauss* NewGaussLine(int vertex1,int vertex2,int order){_error_("not implemented yet");}; -
issm/trunk-jpl/src/c/classes/Elements/Tria.cpp
r25436 r25442 3547 3547 3548 3548 return new GaussTria(point1,fraction1,fraction2,order); 3549 } 3550 /*}}}*/ 3551 Gauss* Tria::NewGauss(IssmDouble fraction1,IssmDouble fraction2,int order){/*{{{*/ 3552 3553 return new GaussTria(fraction1,fraction2,order); 3549 3554 } 3550 3555 /*}}}*/ -
issm/trunk-jpl/src/c/classes/Elements/Tria.h
r25406 r25442 208 208 Gauss* NewGauss(int point1,IssmDouble fraction1,IssmDouble fraction2,bool mainlyfloating,int order); 209 209 Gauss* NewGauss(int point1,IssmDouble fraction1,IssmDouble fraction2,int order); 210 Gauss* NewGauss(IssmDouble fraction1,IssmDouble fraction2,int order); 210 211 Gauss* NewGauss(IssmDouble* xyz_list, IssmDouble* xyz_list_front,int order_horiz,int order_vert); 211 212 Gauss* NewGaussBase(int order); -
issm/trunk-jpl/src/c/classes/gauss/GaussTria.cpp
r25441 r25442 349 349 } 350 350 /*}}}*/ 351 GaussTria::GaussTria(IssmDouble r1,IssmDouble r2,int order){/*{{{*/ 352 353 /* 354 * ^ 355 * ------------------ 356 * 1|\ | 357 * | \ | 358 * | \ | 359 * | \ | 360 * | \ | 361 * | \ | 362 * | +(x,y) \ | 363 * | \| 364 * +---------------+--> 365 * 0 1 366 * 367 */ 368 int ig; 369 IssmDouble x,y; 370 IssmDouble xy_list[3][2]; 371 372 /*Double number of gauss points*/ 373 GaussTria *gauss1 = NULL; //blue 374 GaussTria *gauss2 = NULL; //green 375 GaussTria *gauss3 = NULL; //red 376 gauss1=new GaussTria(order); 377 gauss2=new GaussTria(order); 378 gauss3=new GaussTria(order); 379 380 this->numgauss = gauss1->numgauss + gauss2->numgauss + gauss3->numgauss; 381 this->coords1=xNew<IssmDouble>(this->numgauss); 382 this->coords2=xNew<IssmDouble>(this->numgauss); 383 this->coords3=xNew<IssmDouble>(this->numgauss); 384 this->weights=xNew<IssmDouble>(this->numgauss); 385 386 for(ig=0;ig<gauss1->numgauss;ig++){ // Add the first triangle gauss points (BLUE) 387 this->coords1[ig]=gauss1->coords1[ig]; 388 this->coords2[ig]=gauss1->coords2[ig]; 389 this->coords3[ig]=gauss1->coords3[ig]; 390 this->weights[ig]=gauss1->weights[ig]*r1*r2; 391 } 392 for(ig=0;ig<gauss2->numgauss;ig++){ // Add the second triangle gauss points (GREEN) 393 this->coords1[gauss1->numgauss+ig]=gauss2->coords1[ig]; 394 this->coords2[gauss1->numgauss+ig]=gauss2->coords2[ig]; 395 this->coords3[gauss1->numgauss+ig]=gauss2->coords3[ig]; 396 this->weights[gauss1->numgauss+ig]=gauss2->weights[ig]*r1*(1-r2); 397 } 398 for(ig=0;ig<gauss3->numgauss;ig++){ // Add the second triangle gauss points (RED) 399 this->coords1[gauss1->numgauss+gauss2->numgauss+ig]=gauss3->coords1[ig]; 400 this->coords2[gauss1->numgauss+gauss2->numgauss+ig]=gauss3->coords2[ig]; 401 this->coords3[gauss1->numgauss+gauss2->numgauss+ig]=gauss3->coords3[ig]; 402 this->weights[gauss1->numgauss+gauss2->numgauss+ig]=gauss3->weights[ig]*(1-r1); 403 } 404 405 /*Delete gauss points*/ 406 delete gauss1; 407 delete gauss2; 408 delete gauss3; 409 410 /*Initialize static fields as undefined*/ 411 ig = -1; 412 weight=UNDEF; 413 coord1=UNDEF; 414 coord2=UNDEF; 415 coord3=UNDEF; 416 } 417 /*}}}*/ 351 418 GaussTria::~GaussTria(){/*{{{*/ 352 419 xDelete<IssmDouble>(weights); -
issm/trunk-jpl/src/c/classes/gauss/GaussTria.h
r25438 r25442 33 33 GaussTria(int index,IssmDouble r1, IssmDouble r2,bool maintlyfloating,int order); 34 34 GaussTria(int index,IssmDouble r1, IssmDouble r2,int order); 35 GaussTria(IssmDouble r1, IssmDouble r2,int order); 35 36 GaussTria(IssmDouble area_coordinates[2][3],int order); 36 37 ~GaussTria();
Note:
See TracChangeset
for help on using the changeset viewer.