Changeset 17016


Ignore:
Timestamp:
12/11/13 12:29:25 (11 years ago)
Author:
insa
Message:

NEW: added new finite element for FS

Location:
issm/trunk/src
Files:
5 edited

Legend:

Unmodified
Added
Removed
  • issm/trunk/src/c/classes/Elements/PentaRef.cpp

    r16560 r17016  
    2222#define NUMNODESP2xP1 12
    2323#define NUMNODESP2    15
     24#define NUMNODESP2xP4 30
     25#define NUMNODESP1xP3 12
    2426
    2527/*Object constructors and destructor*/
     
    13331335                        basis[14]=2.*gauss->coord1*gauss->coord2*zeta*(zeta+1.);
    13341336                        return;
     1337                       
     1338                        /*
     1339                                                MODIFICATION
     1340                        */
     1341                case P2xP4Enum :
     1342                        /*Corner nodes*/
     1343                        basis[ 0]=(2./3.)*gauss->coord1*(2.*gauss->coord1-1.)*(zeta-1.)*(zeta-0.5 )*(zeta)*(zeta+0.5);
     1344                        basis[ 1]=(2./3.)*gauss->coord2*(2.*gauss->coord2-1.)*(zeta-1.)*(zeta-0.5 )*(zeta)*(zeta+0.5);
     1345                        basis[ 2]=(2./3.)*gauss->coord3*(2.*gauss->coord3-1.)*(zeta-1.)*(zeta-0.5 )*(zeta)*(zeta+0.5);
     1346                        basis[ 3]=(2./3.)*gauss->coord1*(2.*gauss->coord1-1.)*(zeta-0.5)*(zeta)*(zeta+0.5)*(zeta +1.);
     1347                        basis[ 4]=(2./3.)*gauss->coord2*(2.*gauss->coord2-1.)*(zeta-0.5)*(zeta)*(zeta+0.5)*(zeta +1.);
     1348                        basis[ 5]=(2./3.)*gauss->coord3*(2.*gauss->coord3-1.)*(zeta-0.5)*(zeta)*(zeta+0.5)*(zeta +1.);
     1349               
     1350                        /*mid-sides of quads*/
     1351                        basis[ 6]=4.*gauss->coord1*(zeta-1.)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
     1352                        basis[ 7]=4.*gauss->coord2*(zeta-1.)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
     1353                        basis[ 8]=4.*gauss->coord3*(zeta-1.)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
     1354                       
     1355                        /*mid-sides of triangles*/
     1356                        basis[ 9]=(2./3.)*gauss->coord2*gauss->coord3*(zeta-1.)*(zeta-0.5)*zeta*(zeta+0.5);
     1357                        basis[10]=(2./3.)*gauss->coord1*gauss->coord3*(zeta-1.)*(zeta-0.5)*zeta*(zeta+0.5);
     1358                        basis[11]=(2./3.)*gauss->coord1*gauss->coord2*(zeta-1.)*(zeta-0.5)*zeta*(zeta+0.5);
     1359                        basis[12]=(2./3.)*gauss->coord2*gauss->coord3*(zeta-0.5)*zeta*(zeta+0.5)*(zeta+1.);
     1360                        basis[13]=(2./3.)*gauss->coord1*gauss->coord3*(zeta-0.5)*zeta*(zeta+0.5)*(zeta+1.);
     1361                        basis[14]=(2./3.)*gauss->coord1*gauss->coord2*(zeta-0.5)*zeta*(zeta+0.5)*(zeta+1.);
     1362               
     1363                        /*quarter-sides of quads*/
     1364                        basis[15]=-(2./3.)*gauss->coord1*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
     1365                        basis[16]=-(2./3.)*gauss->coord2*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
     1366                        basis[17]=-(2./3.)*gauss->coord3*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
     1367                        basis[18]=-(2./3.)*gauss->coord1*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
     1368                        basis[19]=-(2./3.)*gauss->coord2*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
     1369                        basis[20]=-(2./3.)*gauss->coord3*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
     1370
     1371                        /* mid-sides of interior triangles*/
     1372                        basis[21]=-(8./3.)*gauss->coord2*gauss->coord3*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
     1373                        basis[22]=-(8./3.)*gauss->coord1*gauss->coord3*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
     1374                        basis[23]=-(8./3.)*gauss->coord1*gauss->coord2*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
     1375                        basis[24]=16.*gauss->coord2*gauss->coord3*(zeta-1.0)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
     1376                        basis[25]=16.*gauss->coord1*gauss->coord3*(zeta-1.0)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
     1377                        basis[26]=16.*gauss->coord1*gauss->coord2*(zeta-1.0)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
     1378                        basis[27]=-(8./3.)*gauss->coord2*gauss->coord3*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
     1379                        basis[28]=-(8./3.)*gauss->coord1*gauss->coord3*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
     1380                        basis[29]=-(8./3.)*gauss->coord1*gauss->coord2*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
     1381                        return;
     1382                       
     1383                       
     1384                case P1xP3Enum :
     1385                        /*Corner nodes*/
     1386                        basis[ 0]=-(9.)/(16.)*gauss->coord1*(zeta-1)*(zeta-1./3.)*(zeta+1./3.);
     1387                        basis[ 1]=-(9.)/(16.)*gauss->coord2*(zeta-1)*(zeta-1./3.)*(zeta+1./3.);
     1388                        basis[ 2]=-(9.)/(16.)*gauss->coord3*(zeta-1)*(zeta-1./3.)*(zeta+1./3.);
     1389                        basis[ 3]=(9.)/(16.)*gauss->coord1*(zeta-1./3.)*(zeta+1./3.)*(zeta+1.);
     1390                        basis[ 4]=(9.)/(16.)*gauss->coord2*(zeta-1./3.)*(zeta+1./3.)*(zeta+1.);
     1391                        basis[ 5]=(9.)/(16.)*gauss->coord3*(zeta-1./3.)*(zeta+1./3.)*(zeta+1.);
     1392               
     1393                        /*quarter-sides of quads*/
     1394                        basis[ 6]=(27.)/(16.)*gauss->coord1*(zeta-1)*(zeta-1./3.)*(zeta+1.);
     1395                        basis[ 7]=(27.)/(16.)*gauss->coord2*(zeta-1)*(zeta-1./3.)*(zeta+1.);
     1396                        basis[ 8]=(27.)/(16.)*gauss->coord3*(zeta-1)*(zeta-1./3.)*(zeta+1.);
     1397                        basis[ 9]=-(27.)/(16.)*gauss->coord1*(zeta-1)*(zeta+1./3.)*(zeta+1.);
     1398                        basis[10]=-(27.)/(16.)*gauss->coord2*(zeta-1)*(zeta+1./3.)*(zeta+1.);
     1399                        basis[11]=-(27.)/(16.)*gauss->coord3*(zeta-1)*(zeta+1./3.)*(zeta+1.);
     1400                        return;
     1401                       
    13351402                default:
    13361403                        _error_("Element type "<<EnumToStringx(this->element_type)<<" not supported yet");
     
    17311798                        dbasis[NUMNODESP2*2+14] = 2.*gauss->coord1*gauss->coord2*(2.*zeta+1.);
    17321799                        return;
     1800                       
     1801                        /*
     1802                        Modification
     1803                        */
     1804                case P2xP4Enum :
     1805                        /*Nodal function 1*/
     1806                        dbasis[NUMNODESP2xP4*0+0 ] = (-2* gauss->coord1 + 0.5 ) *(2./3.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta+0.5) ;
     1807                        dbasis[NUMNODESP2xP4*1+0 ] = (-((2.*SQRT3)/(3.))*gauss->coord1 + (SQRT3/6.) )*(2./3.)*(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta+0.5) ;
     1808                        dbasis[NUMNODESP2xP4*2+0 ] = (2./3.)* gauss->coord1 *(2.* gauss->coord1 -1)*( (2.*zeta-1)*(zeta -0.5)*(zeta +0.5) + 2.* zeta *zeta *(zeta -1.));
     1809      /*Nodal function 2*/
     1810                        dbasis[NUMNODESP2xP4*0+1 ] = (2.*gauss->coord2 - 0.5 ) *(2./3.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta+0.5) ;
     1811                        dbasis[NUMNODESP2xP4*1+1 ] = (-((2.*SQRT3)/(3.))*gauss->coord2 + (SQRT3/6.) )*(2./3.)*(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta+0.5) ;
     1812      dbasis[NUMNODESP2xP4*2+1 ] = (2./3.)*gauss->coord2*(2.*gauss->coord2 -1.)*((2.*zeta-1.)*(zeta -0.5)*(zeta +0.5) + 2. * zeta *zeta*(zeta -1.));
     1813      /*Nodal function 3*/
     1814                        dbasis[NUMNODESP2xP4*0+2 ] = 0. ;
     1815                        dbasis[NUMNODESP2xP4*1+2 ] = (((4.*SQRT3)/(3.))*gauss->coord3 - (SQRT3)/(3.))*(2./3.)*(zeta -1.)*(zeta-0.5)*(zeta)*(zeta+0.5);
     1816      dbasis[NUMNODESP2xP4*2+2 ] = (2./3.)*gauss->coord3*(2.* gauss->coord3 -1.)*( (2.*zeta-1.)*(zeta -0.5)*(zeta +0.5) + 2.*zeta *zeta*(zeta -1.));
     1817      /*Nodal function 4*/
     1818                        dbasis[NUMNODESP2xP4*0+3 ] = (-2.* gauss->coord1 + 0.5 ) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
     1819                        dbasis[NUMNODESP2xP4*1+3 ] = (-((2.*SQRT3)/(3.)) *gauss->coord1 + (SQRT3)/(6.) ) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
     1820      dbasis[NUMNODESP2xP4*2+3 ] = (2./3.)*gauss->coord1*(2.*gauss->coord1 -1.)*( (2.*zeta+1.)*(zeta -0.5)*(zeta +0.5) + 2.*zeta *zeta*(zeta +1.));
     1821      /*Nodal function 5*/
     1822                        dbasis[NUMNODESP2xP4*0+4 ] = (2*gauss->coord2 - 0.5 ) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
     1823                        dbasis[NUMNODESP2xP4*1+4 ] = -((2.*SQRT3)/(3.)) *gauss->coord2 + (SQRT3/(6.))*(2./3.)*(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
     1824      dbasis[NUMNODESP2xP4*2+4 ] = (2./3.)*gauss->coord2 *(2.*gauss->coord2 -1.)*( (2.*zeta+1.)*(zeta -0.5)*(zeta +0.5) + 2.*zeta *zeta*(zeta +1.));
     1825      /*Nodal function 6*/
     1826                        dbasis[NUMNODESP2xP4*0+5 ] = 0. ;
     1827                        dbasis[NUMNODESP2xP4*1+5 ] = (((4.*SQRT3)/(3.))*gauss->coord3 - SQRT3/3. ) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
     1828      dbasis[NUMNODESP2xP4*2+5 ] = (2./3.)*gauss->coord3 *(2.*gauss->coord3 -1.)*( (2.*zeta+1.)*(zeta -0.5)*(zeta +0.5) + 2.*zeta *zeta*(zeta +1));
     1829      /*Nodal function 7*/
     1830                        dbasis[NUMNODESP2xP4*0+6 ] =  -2. *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. ) ;
     1831                        dbasis[NUMNODESP2xP4*1+6 ] =  (-4.*SQRT3)/(6.) *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1.)  ;
     1832      dbasis[NUMNODESP2xP4*2+6 ] =  4.*gauss->coord1*( 4.*zeta *zeta*zeta - (5./2.)*zeta );
     1833      /*Nodal function 8*/
     1834                        dbasis[NUMNODESP2xP4*0+7 ] =  2.*(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. ) ;
     1835                        dbasis[NUMNODESP2xP4*1+7 ] =  ((-4.*SQRT3)/(6.)) *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. )  ;
     1836      dbasis[NUMNODESP2xP4*2+7 ] =  4.* gauss->coord2*( 4.*zeta *zeta*zeta - (5./2.)*zeta );
     1837      /*Nodal function 9*/
     1838                        dbasis[NUMNODESP2xP4*0+8 ] = 0. ;
     1839                        dbasis[NUMNODESP2xP4*1+8 ] = ((4.*SQRT3)/(3.)) *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. )  ;
     1840                        dbasis[NUMNODESP2xP4*2+8 ] =  4.*gauss->coord3*( 4.*zeta *zeta*zeta - (5./2.)*zeta );
     1841      /*Nodal function 10*/
     1842                        dbasis[NUMNODESP2xP4*0+9 ] = (4./3.)*(gauss->coord3 )*(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta+0.5) ;
     1843                        dbasis[NUMNODESP2xP4*1+9 ] =  4.* (((- SQRT3)/6.)* gauss->coord3+ (SQRT3/3.) *gauss->coord2) *(2./3.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta+0.5);
     1844                        dbasis[NUMNODESP2xP4*2+9 ] = 4.* gauss->coord2 * gauss->coord3 *(2./3.)*((2.*zeta-1.)*(zeta -0.5)*(zeta +0.5) + 2.* zeta*zeta*(zeta -1.));
     1845      /*Nodal function 11*/
     1846                        dbasis[NUMNODESP2xP4*0+10] = -(4./3.)* gauss->coord3*(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta+0.5) ;
     1847                        dbasis[NUMNODESP2xP4*1+10] =  4.* ((SQRT3/3.)*gauss->coord1- (SQRT3/6.)*gauss->coord3) *(2./3.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta+0.5);
     1848                        dbasis[NUMNODESP2xP4*2+10] = 4.* gauss->coord3*gauss->coord1 *(2./3.)*( (2*zeta-1.)*(zeta -0.5)*(zeta +0.5) + 2* zeta * zeta*(zeta -1));
     1849      /*Nodal function 12*/
     1850                        dbasis[NUMNODESP2xP4*0+11] = 4.* (gauss->coord1*0.5- gauss->coord2 *0.5 )*(2./3.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta+0.5) ;
     1851                        dbasis[NUMNODESP2xP4*1+11] =  4.* (- (SQRT3/6.) *gauss->coord2 - (SQRT3/6.) *gauss->coord1) *(2./3.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta+0.5);
     1852                        dbasis[NUMNODESP2xP4*2+11] = 4.* gauss->coord1*gauss->coord2 *(2./3.) *( (2.*zeta-1)*(zeta -0.5)*(zeta +0.5) + 2* zeta* zeta*(zeta -1));
     1853      /*Nodal function 13*/
     1854                        dbasis[NUMNODESP2xP4*0+12] = (4./3.) * gauss->coord3 *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1855                        dbasis[NUMNODESP2xP4*1+12] =  4.* ( (- SQRT3/6. ) * gauss->coord3 + (SQRT3/3.) *gauss->coord2) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
     1856                        dbasis[NUMNODESP2xP4*2+12] = 4.* gauss->coord2*gauss->coord3 *(2./3.)*( (2.*zeta+1.)*(zeta -0.5)*(zeta +0.5) + 2*zeta*zeta*(zeta +1.));
     1857      /*Nodal function 14*/
     1858                        dbasis[NUMNODESP2xP4*0+13] = (- 4./3.) *gauss->coord3*(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1859                        dbasis[NUMNODESP2xP4*1+13] =  4.* ((SQRT3/3.)*gauss->coord1- (SQRT3/6.)*gauss->coord3) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
     1860                        dbasis[NUMNODESP2xP4*2+13] = 4.* gauss->coord3*gauss->coord1 *(2./3.)*( (2.*zeta+1.)*(zeta -0.5)*(zeta +0.5) + 2.*zeta *zeta*(zeta +1.));
     1861      /*Nodal function 15*/
     1862                        dbasis[NUMNODESP2xP4*0+14] = 4.* ( gauss->coord1*0.5- gauss->coord2*0.5) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1863                        dbasis[NUMNODESP2xP4*1+14] =  4.* (- (SQRT3/6.)*gauss->coord2 - (SQRT3/6.)*gauss->coord1) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
     1864                        dbasis[NUMNODESP2xP4*2+14] = 4.* gauss->coord1*gauss->coord2 *(2./3.)*( (2.*zeta+1.)*(zeta -0.5)*(zeta +0.5) + 2* zeta *zeta*(zeta +1.));
     1865      /*Nodal function 16*/
     1866                        dbasis[NUMNODESP2xP4*0+15] =  (8./6.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1867                        dbasis[NUMNODESP2xP4*1+15] =  (8./3.) *(SQRT3/6.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1868                        dbasis[NUMNODESP2xP4*2+15] =  -(8./3.) *gauss->coord1 *((2.*zeta -1. )*(zeta-0.5)*(zeta +1.) +  zeta*(zeta -1.)*( 2.*zeta + 0.5));
     1869      /*Nodal function 17*/
     1870                        dbasis[NUMNODESP2xP4*0+16] = - (8./6.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1871                        dbasis[NUMNODESP2xP4*1+16] =  (8./3.) * (SQRT3/6.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1872                        dbasis[NUMNODESP2xP4*2+16] =  -(8./3.)*gauss->coord2 *((2.*zeta -1. )*(zeta-0.5)*(zeta +1.) +  zeta *(zeta -1.)*( 2.*zeta + 0.5));
     1873      /*Nodal function 18*/
     1874                        dbasis[NUMNODESP2xP4*0+17] = 0. ;
     1875                        dbasis[NUMNODESP2xP4*1+17] =  (8./3.) *(SQRT3/3.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1876                        dbasis[NUMNODESP2xP4*2+17] =  -(8./3.)*gauss->coord3 *((2.* zeta -1. ) *(zeta-0.5)* (zeta +1.) +  zeta *(zeta -1.)*( 2.*zeta + 0.5));
     1877      /*Nodal function 19*/
     1878                        dbasis[NUMNODESP2xP4*0+18] =  (4./3.) *(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1879                        dbasis[NUMNODESP2xP4*1+18] =  (8./3.)* (SQRT3/6.)*(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1880                        dbasis[NUMNODESP2xP4*2+18] =  -(8./3.)*gauss->coord1*((2.*zeta -1. ) *(zeta+0.5)* (zeta +1.) +  zeta* (zeta -1.)*( 2.*zeta + 3./2.));
     1881      /*Nodal function 20*/
     1882                        dbasis[NUMNODESP2xP4*0+19] =  (8./6.) *(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1883                        dbasis[NUMNODESP2xP4*1+19] =  (8./3.)*(SQRT3/6) *(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1884                        dbasis[NUMNODESP2xP4*2+19] =  -(8./3.)*gauss->coord2*((2.*zeta -1. )*(zeta+0.5)*(zeta +1.) +  zeta*(zeta -1.)*( 2.*zeta + 3./2.));
     1885      /*Nodal function 21*/
     1886                        dbasis[NUMNODESP2xP4*0+20] = 0 ;
     1887                        dbasis[NUMNODESP2xP4*1+20] =  -(8./3.)* (SQRT3/3) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1888                        dbasis[NUMNODESP2xP4*2+20] =  -(8./3.)*gauss->coord3 *(2. *zeta -1. )*(zeta+0.5)*(zeta +1.) +  zeta*(zeta -1.)*( 2.*zeta + 3./2.);
     1889      /*Nodal function 22*/
     1890                        dbasis[NUMNODESP2xP4*0+21] = -(32./6.) *gauss->coord3 *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. );
     1891                        dbasis[NUMNODESP2xP4*1+21] = -(32./3.)*(- (SQRT3/6.)*gauss->coord3+ (SQRT3/3.)*gauss->coord2 ) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1892                        dbasis[NUMNODESP2xP4*2+21] = -(32./3.)*gauss->coord2 *gauss->coord3*((2.*zeta -1. )*(zeta-0.5)*(zeta +1.) +  zeta*(zeta -1.)*( 2.*zeta + 0.5));
     1893      /*Nodal function 23*/
     1894                        dbasis[NUMNODESP2xP4*0+22] = (32./6.) *gauss->coord3 *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. );
     1895                        dbasis[NUMNODESP2xP4*1+22] = -(32./3.)*( (SQRT3/3.)*gauss->coord1- (SQRT3/6.)*gauss->coord3 )*(zeta -1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1896                        dbasis[NUMNODESP2xP4*2+22] = -(32./3.) *gauss->coord1*gauss->coord3*((2.*zeta -1. )*(zeta-0.5)*(zeta +1.) +  zeta*(zeta -1.)*( 2.*zeta + 0.5));
     1897      /*Nodal function 24*/
     1898                        dbasis[NUMNODESP2xP4*0+23] = -(32./3.)*( gauss->coord1*0.5- gauss->coord2*0.5) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. );
     1899                        dbasis[NUMNODESP2xP4*1+23] = -(32./3.)*(- (SQRT3/6.)*gauss->coord2- (SQRT3/6.)*gauss->coord1 ) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1900                        dbasis[NUMNODESP2xP4*2+23] = -(32./3.)*gauss->coord1*gauss->coord2*((2.*zeta -1. )* (zeta-0.5) *(zeta +1.) +  zeta* (zeta -1.)*( 2.*zeta + 0.5));
     1901      /*Nodal function 25*/
     1902                        dbasis[NUMNODESP2xP4*0+24] = 8. *gauss->coord3 *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. );
     1903                        dbasis[NUMNODESP2xP4*1+24] = 16.*(- (SQRT3/6.)*gauss->coord3 + (SQRT3/3.) *gauss->coord2) *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. ) ;
     1904                        dbasis[NUMNODESP2xP4*2+24] = 16.* gauss->coord2 * gauss->coord3* ( 4.* zeta*zeta*zeta - (5./2.) *zeta );
     1905      /*Nodal function 26*/
     1906                        dbasis[NUMNODESP2xP4*0+25] = 8. *gauss->coord3*(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. );
     1907                        dbasis[NUMNODESP2xP4*1+25] = -16.*((SQRT3/3.)*gauss->coord1- (SQRT3/6.)*gauss->coord3 )*(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. ) ;
     1908                        dbasis[NUMNODESP2xP4*2+25] = 16. * gauss->coord1 * gauss->coord3 *( 4. *zeta*zeta*zeta - (5./2.)* zeta );
     1909      /*Nodal function 27*/
     1910                        dbasis[NUMNODESP2xP4*0+26] = 16.*( gauss->coord1*0.5-gauss->coord2*0.5) *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. );
     1911                        dbasis[NUMNODESP2xP4*1+26] = 16.*(- (SQRT3/6.)*gauss->coord2- (SQRT3/6.)*gauss->coord1 )*(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. ) ;
     1912                        dbasis[NUMNODESP2xP4*2+26] = 16. *gauss->coord1 *gauss->coord2 *( 4.* zeta*zeta*zeta - (5./2.)*zeta );
     1913      /*Nodal function 28*/
     1914                        dbasis[NUMNODESP2xP4*0+27] = -(32./6.)* gauss->coord3*(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. );
     1915                        dbasis[NUMNODESP2xP4*1+27] = -(32./3.)*(- (SQRT3/6.)*gauss->coord3+ (SQRT3/3.)*gauss->coord2 ) *(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1916                        dbasis[NUMNODESP2xP4*2+27] = -(32./3.)* gauss->coord2*gauss->coord3* ((2.*zeta -1. )*(zeta+0.5)*(zeta +1.) +  zeta*(zeta -1.)*( 2.*zeta + (3./2.)));
     1917      /*Nodal function 29*/
     1918                        dbasis[NUMNODESP2xP4*0+28] = (32./6.) *gauss->coord3 *(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. );
     1919                        dbasis[NUMNODESP2xP4*1+28] = -(32./3.)*((SQRT3/3.)*gauss->coord1- (SQRT3/6.)*gauss->coord3 ) *(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1920                        dbasis[NUMNODESP2xP4*2+28] = -(32./3.)* gauss->coord1*gauss->coord3* ((2.*zeta -1. )*(zeta+0.5)*(zeta +1.) +  zeta*(zeta -1.)*( 2.*zeta + (3./2.)));
     1921      /*Nodal function 30*/
     1922                        dbasis[NUMNODESP2xP4*0+29] = -(32./3.)*( gauss->coord1*0.5- gauss->coord2*0.5)*(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. );
     1923                        dbasis[NUMNODESP2xP4*1+29] = -(32./3.)*(- (SQRT3/6.)*gauss->coord2- (SQRT3/6.)*gauss->coord1 ) *(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1924                        dbasis[NUMNODESP2xP4*2+29] = -(32./3.)* gauss->coord1*gauss->coord2 *((2.*zeta -1. )*(zeta+0.5)*(zeta +1) +   zeta*(zeta -1.)*( 2.*zeta + (3./2.)));
     1925                        return;
     1926                       
     1927                case P1xP3Enum :
     1928                        /*Nodal function 1*/
     1929            dbasis[NUMNODESP1xP3*0+0 ] =  (9./32.)*(zeta-1)*(zeta-1./3.)*(zeta+1./3.);
     1930            dbasis[NUMNODESP1xP3*1+0 ] = ((3.*SQRT3)/32.)*(zeta-1)*(zeta-1./3.)*(zeta+1./3.);
     1931            dbasis[NUMNODESP1xP3*2+0 ] =- (9./16.)* gauss->coord1 *( 2. *zeta *( zeta -1. ) + ( zeta - (1./3.) )*( zeta + (1./3.) ));
     1932                        /*Nodal function 2*/
     1933            dbasis[NUMNODESP1xP3*0+1 ] = - (9./32.)*(zeta-1)*(zeta-1./3.)*(zeta+1./3.);
     1934            dbasis[NUMNODESP1xP3*1+1 ] = ((3.*SQRT3)/32.) *(zeta-1)*(zeta-1./3.)*(zeta+1./3.);
     1935            dbasis[NUMNODESP1xP3*2+1 ] =- (9./16.)*gauss->coord2 *( 2.* zeta* ( zeta -1. ) + ( zeta - (1./3.) )*( zeta + (1./3.) ));
     1936                  /*Nodal function 3*/
     1937            dbasis[NUMNODESP1xP3*0+2 ] =  0.;
     1938            dbasis[NUMNODESP1xP3*1+2 ] = - ((3.*SQRT3)/32.)*(zeta-1)*(zeta-1./3.)*(zeta+1./3.);
     1939            dbasis[NUMNODESP1xP3*2+2 ] = - (9./16.)* gauss->coord3* ( 2. *zeta *( zeta -1. ) + ( zeta - (1./3.) )*( zeta + (1./3.) ));
     1940                        /*Nodal function 4*/     
     1941            dbasis[NUMNODESP1xP3*0+3 ] = -  (9./32.)*(zeta-1./3.)*(zeta+1./3.)*(zeta+1.);
     1942            dbasis[NUMNODESP1xP3*1+3 ] =  -((3.*SQRT3)/32.) *(zeta-1./3.)*(zeta+1./3.)*(zeta+1.);
     1943            dbasis[NUMNODESP1xP3*2+3 ] = (9./16.)* gauss->coord1*( 2.* zeta *( zeta +1. ) + ( zeta - (1./3.) )*( zeta + (1./3.) ));
     1944                        /*Nodal function 5*/   
     1945            dbasis[NUMNODESP1xP3*0+4 ] =   (9./32.)* (zeta-1./3.)*(zeta+1./3.)*(zeta+1.);
     1946            dbasis[NUMNODESP1xP3*1+4 ] = - ((3.*SQRT3)/32.) *(zeta-1./3.)*(zeta+1./3.)*(zeta+1.);
     1947            dbasis[NUMNODESP1xP3*2+4 ] = (9./16.)* gauss->coord2* ( 2.* zeta *( zeta +1. ) + ( zeta - (1./3.) )*( zeta + (1./3.) ));
     1948                        /*Nodal function 6*/   
     1949            dbasis[NUMNODESP1xP3*0+5 ] =  0.;
     1950            dbasis[NUMNODESP1xP3*1+5 ] =  ((3.*SQRT3)/16.)  *(zeta-1./3.)*(zeta+1./3.)*(zeta+1.);
     1951            dbasis[NUMNODESP1xP3*2+5 ] =  (9./16.)* gauss->coord3 *( 2.* zeta * ( zeta  + 1. ) + ( zeta - (1./3.) )*( zeta + (1./3.) ));
     1952                        /*Nodal function 7*/   
     1953            dbasis[NUMNODESP1xP3*0+6 ] = -  (27./32.) *(zeta-1)*(zeta-1./3.)*(zeta+1.);
     1954            dbasis[NUMNODESP1xP3*1+6 ] = -  (27./32.) *(zeta-1)*(zeta-1./3.)*(zeta+1.);
     1955            dbasis[NUMNODESP1xP3*2+6 ] =  gauss->coord1*(27./16.)*( 2.* zeta *( zeta - (1./3.)) + ( zeta - 1. )*( zeta + 1. ));
     1956                        /*Nodal function 8*/   
     1957            dbasis[NUMNODESP1xP3*0+7 ] =  (27./32.) *(zeta-1)*(zeta-1./3.)*(zeta+1.);
     1958            dbasis[NUMNODESP1xP3*1+7 ] = -(27./32.) *(zeta-1)*(zeta-1./3.)*(zeta+1.);
     1959            dbasis[NUMNODESP1xP3*2+7 ] =  gauss->coord2*(27./16.)*( 2.* zeta *( zeta - (1./3.)) + ( zeta - 1. )*( zeta + 1. ));
     1960                        /*Nodal function 9*/   
     1961            dbasis[NUMNODESP1xP3*0+8 ] = 0.;
     1962            dbasis[NUMNODESP1xP3*1+8 ] =  ((9.*SQRT3)/16.) *(zeta-1.)*(zeta-1./3.)*(zeta+1.);
     1963            dbasis[NUMNODESP1xP3*2+8 ] =  gauss->coord3*(27./16.)*( 2. *zeta *( zeta - (1./3.)) + ( zeta - 1. )*( zeta + 1. ));
     1964                        /*Nodal function 10*/   
     1965            dbasis[NUMNODESP1xP3*0+9 ] = - (27./32.) *(zeta-1.)*(zeta+1./3.)*(zeta+1.);
     1966            dbasis[NUMNODESP1xP3*1+9 ] = - ((9.*SQRT3)/32.) *(zeta-1.)*(zeta+1./3.)*(zeta+1.);
     1967            dbasis[NUMNODESP1xP3*2+9 ] =  gauss->coord1 *(27./16.)*( 2* zeta *( zeta + (1./3.)) + ( zeta - 1. )*( zeta + 1. ));
     1968                        /*Nodal function 11*/   
     1969            dbasis[NUMNODESP1xP3*0+10] =  (27./32.) *(zeta-1)*(zeta+1./3.)*(zeta+1);
     1970            dbasis[NUMNODESP1xP3*1+10] = - ((9.*SQRT3)/32.)  *(zeta-1.)*(zeta+1./3.)*(zeta+1);
     1971            dbasis[NUMNODESP1xP3*2+10] =  gauss->coord2 *(27./16.) * 2.* zeta *( zeta + (1./3.) + ( zeta - 1. )*( zeta + 1. ));
     1972                        /*Nodal function 12*/   
     1973            dbasis[NUMNODESP1xP3*0+11] = 0.;
     1974            dbasis[NUMNODESP1xP3*1+11] =  ((9.*SQRT3)/16.) *(zeta-1.)*(zeta+1./3.)*(zeta+1);
     1975            dbasis[NUMNODESP1xP3*2+11] =  gauss->coord3 *(27./16.)* 2.* zeta *( zeta + (1./3.) + ( zeta - 1. )*( zeta + 1. ));
     1976                        return;
     1977                       
    17331978                default:
    17341979                        _error_("Element type "<<EnumToStringx(this->element_type)<<" not supported yet");
     
    19732218                case MINIEnum:              return NUMNODESP1b+NUMNODESP1;
    19742219                case TaylorHoodEnum:        return NUMNODESP2+NUMNODESP1;
     2220                case P2xP4Enum:                         return NUMNODESP2xP4;
     2221                case P1xP3Enum:                         return NUMNODESP1xP3;
    19752222                default:       _error_("Element type "<<EnumToStringx(this->element_type)<<" not supported yet");
    19762223        }
  • issm/trunk/src/c/shared/Enum/EnumDefinitions.h

    r16560 r17016  
    524524        MINIcondensedEnum,
    525525        TaylorHoodEnum,
     526        P2xP4Enum,
     527        P1xP3Enum,
    526528        /*}}}*/
    527529        /*Results{{{*/
  • issm/trunk/src/c/shared/Enum/EnumToStringx.cpp

    r16560 r17016  
    516516                case MINIcondensedEnum : return "MINIcondensed";
    517517                case TaylorHoodEnum : return "TaylorHood";
     518                case P2xP4Enum : return "P2xP4";
     519                case P1xP3Enum : return "P1xP3";
    518520                case SaveResultsEnum : return "SaveResults";
    519521                case BoolExternalResultEnum : return "BoolExternalResult";
  • issm/trunk/src/c/shared/Enum/StringToEnumx.cpp

    r16560 r17016  
    528528              else if (strcmp(name,"MINIcondensed")==0) return MINIcondensedEnum;
    529529              else if (strcmp(name,"TaylorHood")==0) return TaylorHoodEnum;
     530              else if (strcmp(name,"P2xP4")==0) return P2xP4Enum;
     531              else if (strcmp(name,"P1xP3")==0) return P1xP3Enum;
    530532              else if (strcmp(name,"SaveResults")==0) return SaveResultsEnum;
    531533              else if (strcmp(name,"BoolExternalResult")==0) return BoolExternalResultEnum;
  • issm/trunk/src/m/enum/EnumDefinitions.py

    r16560 r17016  
    508508def MINIcondensedEnum(): return StringToEnum("MINIcondensed")[0]
    509509def TaylorHoodEnum(): return StringToEnum("TaylorHood")[0]
     510def P2xP4Enum(): return StringToEnum("P2xP4")[0]
     511def P1xP3Enum(): return StringToEnum("P1xP3")[0]
    510512def SaveResultsEnum(): return StringToEnum("SaveResults")[0]
    511513def BoolExternalResultEnum(): return StringToEnum("BoolExternalResult")[0]
Note: See TracChangeset for help on using the changeset viewer.