Changeset 17058


Ignore:
Timestamp:
01/04/14 03:52:50 (11 years ago)
Author:
insa
Message:

Fixed several errors in the derivatives of the P2xP4 Element.

File:
1 edited

Legend:

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

    r17056 r17058  
    13731373                        basis[ 5]=(2./3.)*gauss->coord3*(2.*gauss->coord3-1.)*(zeta-0.5)*(zeta)*(zeta+0.5)*(zeta +1.);
    13741374                        /*mid-sides of quads*/
    1375                         basis[ 6]=4.*gauss->coord1*(zeta-1.)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
    1376                         basis[ 7]=4.*gauss->coord2*(zeta-1.)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
    1377                         basis[ 8]=4.*gauss->coord3*(zeta-1.)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
     1375                        basis[ 6]=4.*gauss->coord1*(2.*gauss->coord1-1)*(zeta-1.)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
     1376                        basis[ 7]=4.*gauss->coord2*(2.*gauss->coord2-1)*(zeta-1.)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
     1377                        basis[ 8]=4.*gauss->coord3*(2.*gauss->coord3-1)*(zeta-1.)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
    13781378                        /*mid-sides of triangles*/
    1379                         basis[ 9]=(2./3.)*gauss->coord2*gauss->coord3*(zeta-1.)*(zeta-0.5)*zeta*(zeta+0.5);
    1380                         basis[10]=(2./3.)*gauss->coord1*gauss->coord3*(zeta-1.)*(zeta-0.5)*zeta*(zeta+0.5);
    1381                         basis[11]=(2./3.)*gauss->coord1*gauss->coord2*(zeta-1.)*(zeta-0.5)*zeta*(zeta+0.5);
    1382                         basis[12]=(2./3.)*gauss->coord2*gauss->coord3*(zeta-0.5)*zeta*(zeta+0.5)*(zeta+1.);
    1383                         basis[13]=(2./3.)*gauss->coord1*gauss->coord3*(zeta-0.5)*zeta*(zeta+0.5)*(zeta+1.);
    1384                         basis[14]=(2./3.)*gauss->coord1*gauss->coord2*(zeta-0.5)*zeta*(zeta+0.5)*(zeta+1.);
     1379                        basis[ 9]=(8./3.)*gauss->coord2*gauss->coord3*(zeta-1.)*(zeta-0.5)*zeta*(zeta+0.5);
     1380                        basis[10]=(8./3.)*gauss->coord1*gauss->coord3*(zeta-1.)*(zeta-0.5)*zeta*(zeta+0.5);
     1381                        basis[11]=(8./3.)*gauss->coord1*gauss->coord2*(zeta-1.)*(zeta-0.5)*zeta*(zeta+0.5);
     1382                        basis[12]=(8./3.)*gauss->coord2*gauss->coord3*(zeta-0.5)*zeta*(zeta+0.5)*(zeta+1.);
     1383                        basis[13]=(8./3.)*gauss->coord1*gauss->coord3*(zeta-0.5)*zeta*(zeta+0.5)*(zeta+1.);
     1384                        basis[14]=(8./3.)*gauss->coord1*gauss->coord2*(zeta-0.5)*zeta*(zeta+0.5)*(zeta+1.);
    13851385                        /*quarter-sides of quads*/
    1386                         basis[15]=-(2./3.)*gauss->coord1*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
    1387                         basis[16]=-(2./3.)*gauss->coord2*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
    1388                         basis[17]=-(2./3.)*gauss->coord3*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
    1389                         basis[18]=-(2./3.)*gauss->coord1*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
    1390                         basis[19]=-(2./3.)*gauss->coord2*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
    1391                         basis[20]=-(2./3.)*gauss->coord3*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
     1386                        basis[15]=-(8./3.)*gauss->coord1*(2.*gauss->coord1-1.)*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
     1387                        basis[16]=-(8./3.)*gauss->coord2*(2.*gauss->coord2-1.)*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
     1388                        basis[17]=-(8./3.)*gauss->coord3*(2.*gauss->coord3-1.)*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
     1389                        basis[18]=-(8./3.)*gauss->coord1*(2.*gauss->coord1-1.)*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
     1390                        basis[19]=-(8./3.)*gauss->coord2*(2.*gauss->coord2-1.)*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
     1391                        basis[20]=-(8./3.)*gauss->coord3*(2.*gauss->coord3-1.)*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
    13921392                        /* mid-sides of interior triangles*/
    1393                         basis[21]=-(8./3.)*gauss->coord2*gauss->coord3*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
    1394                         basis[22]=-(8./3.)*gauss->coord1*gauss->coord3*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
    1395                         basis[23]=-(8./3.)*gauss->coord1*gauss->coord2*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
     1393                        basis[21]=-(32./3.)*gauss->coord2*gauss->coord3*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
     1394                        basis[22]=-(32./3.)*gauss->coord1*gauss->coord3*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
     1395                        basis[23]=-(32./3.)*gauss->coord1*gauss->coord2*(zeta-1.0)*(zeta-0.5)*zeta*(zeta+1.);
    13961396                        basis[24]=16.*gauss->coord2*gauss->coord3*(zeta-1.0)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
    13971397                        basis[25]=16.*gauss->coord1*gauss->coord3*(zeta-1.0)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
    13981398                        basis[26]=16.*gauss->coord1*gauss->coord2*(zeta-1.0)*(zeta-0.5)*(zeta+0.5)*(zeta+1.);
    1399                         basis[27]=-(8./3.)*gauss->coord2*gauss->coord3*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
    1400                         basis[28]=-(8./3.)*gauss->coord1*gauss->coord3*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
    1401                         basis[29]=-(8./3.)*gauss->coord1*gauss->coord2*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
    1402                         return;
    1403                 case P1xP3Enum :
     1399                        basis[27]=-(32./3.)*gauss->coord2*gauss->coord3*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
     1400                        basis[28]=-(32./3.)*gauss->coord1*gauss->coord3*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
     1401                        basis[29]=-(32./3.)*gauss->coord1*gauss->coord2*(zeta-1.0)*zeta*(zeta+0.5)*(zeta+1.);
     1402                        return;
     1403                                        case P1xP3Enum :
    14041404                        /*Corner nodes*/
    14051405                        basis[ 0]=-(9.)/(16.)*gauss->coord1*(zeta-1)*(zeta-1./3.)*(zeta+1./3.);
     
    18391839                        /*Nodal function 5*/
    18401840                        dbasis[NUMNODESP2xP4*0+4 ] = (2*gauss->coord2 - 0.5 ) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
    1841                         dbasis[NUMNODESP2xP4*1+4 ] = -((2.*SQRT3)/(3.)) *gauss->coord2 + (SQRT3/(6.))*(2./3.)*(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
     1841                        dbasis[NUMNODESP2xP4*1+4 ] = (-((2.*SQRT3)/(3.)) *gauss->coord2 + (SQRT3/6.))*(2./3.)*(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
    18421842                        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.));
    18431843                        /*Nodal function 6*/
     
    18461846                        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));
    18471847                        /*Nodal function 7*/
    1848                         dbasis[NUMNODESP2xP4*0+6 ] =  -2. *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. ) ;
    1849                         dbasis[NUMNODESP2xP4*1+6 ] =  (-4.*SQRT3)/(6.) *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1.)  ;
    1850                         dbasis[NUMNODESP2xP4*2+6 ] =  4.*gauss->coord1*( 4.*zeta *zeta*zeta - (5./2.)*zeta );
     1848                        dbasis[NUMNODESP2xP4*0+6 ] =  4.*(-2.* gauss->coord1 + 0.5 ) *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. ) ;
     1849                        dbasis[NUMNODESP2xP4*1+6 ] =  4.*(-((2.*SQRT3)/(3.)) *gauss->coord1 + (SQRT3)/(6.) )*(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1.)  ;
     1850                        dbasis[NUMNODESP2xP4*2+6 ] =  4.*gauss->coord1*(2.*gauss->coord1-1)*( 4.*zeta *zeta*zeta - (5./2.)*zeta );
    18511851                        /*Nodal function 8*/
    1852                         dbasis[NUMNODESP2xP4*0+7 ] =  2.*(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. ) ;
    1853                         dbasis[NUMNODESP2xP4*1+7 ] =  ((-4.*SQRT3)/(6.)) *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. )  ;
    1854                         dbasis[NUMNODESP2xP4*2+7 ] =  4.* gauss->coord2*( 4.*zeta *zeta*zeta - (5./2.)*zeta );
     1852                        dbasis[NUMNODESP2xP4*0+7 ] =  4.*(2*gauss->coord2 - 0.5 )*(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. ) ;
     1853                        dbasis[NUMNODESP2xP4*1+7 ] =  4.*(-((2.*SQRT3)/(3.)) *gauss->coord2 + (SQRT3)/(6.)) *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. )  ;
     1854                        dbasis[NUMNODESP2xP4*2+7 ] =  4.*gauss->coord2*(2.*gauss->coord2-1)*( 4.*zeta *zeta*zeta - (5./2.)*zeta );
    18551855                        /*Nodal function 9*/
    18561856                        dbasis[NUMNODESP2xP4*0+8 ] = 0. ;
    1857                         dbasis[NUMNODESP2xP4*1+8 ] = ((4.*SQRT3)/(3.)) *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. )  ;
    1858                         dbasis[NUMNODESP2xP4*2+8 ] =  4.*gauss->coord3*( 4.*zeta *zeta*zeta - (5./2.)*zeta );
     1857                        dbasis[NUMNODESP2xP4*1+8 ] = 4.*(((4.*SQRT3)/(3.))*gauss->coord3 - SQRT3/3. ) *(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. )  ;
     1858                        dbasis[NUMNODESP2xP4*2+8 ] =  4.*gauss->coord3*(2.*gauss->coord3-1)*( 4.*zeta *zeta*zeta - (5./2.)*zeta );
    18591859                        /*Nodal function 10*/
    18601860                        dbasis[NUMNODESP2xP4*0+9 ] = (4./3.)*(gauss->coord3 )*(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta+0.5) ;
     
    18751875                        /*Nodal function 14*/
    18761876                        dbasis[NUMNODESP2xP4*0+13] = (- 4./3.) *gauss->coord3*(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
    1877                         dbasis[NUMNODESP2xP4*1+13] =  4.* ((SQRT3/3.)*gauss->coord1- (SQRT3/6.)*gauss->coord3) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
     1877                        dbasis[NUMNODESP2xP4*1+13] = 4.* ((SQRT3/3.)*gauss->coord1- (SQRT3/6.)*gauss->coord3) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
    18781878                        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.));
    18791879                        /*Nodal function 15*/
    18801880                        dbasis[NUMNODESP2xP4*0+14] = 4.* ( gauss->coord1*0.5- gauss->coord2*0.5) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
    1881                         dbasis[NUMNODESP2xP4*1+14] =  4.* (- (SQRT3/6.)*gauss->coord2 - (SQRT3/6.)*gauss->coord1) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
     1881                        dbasis[NUMNODESP2xP4*1+14] = 4.* (- (SQRT3/6.)*gauss->coord2 - (SQRT3/6.)*gauss->coord1) *(2./3.) *(zeta - 0.5)*(zeta)*(zeta+0.5)*(zeta +1. );
    18821882                        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.));
    18831883                        /*Nodal function 16*/
    1884                         dbasis[NUMNODESP2xP4*0+15] =  (8./6.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
    1885                         dbasis[NUMNODESP2xP4*1+15] =  (8./3.) *(SQRT3/6.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
    1886                         dbasis[NUMNODESP2xP4*2+15] =  -(8./3.) *gauss->coord1 *((2.*zeta -1. )*(zeta-0.5)*(zeta +1.) +  zeta*(zeta -1.)*( 2.*zeta + 0.5));
     1884                        dbasis[NUMNODESP2xP4*0+15] =  -(8./3.) *(-2.* gauss->coord1 + 0.5 )*(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1885                        dbasis[NUMNODESP2xP4*1+15] =  -(8./3.) *(-((2.*SQRT3)/(3.)) *gauss->coord1 + (SQRT3)/(6.) ) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1886                        dbasis[NUMNODESP2xP4*2+15] =  -(8./3.) *gauss->coord1*(2.*gauss->coord1-1) *((2.*zeta -1. )*(zeta-0.5)*(zeta +1.) +  zeta*(zeta -1.)*( 2.*zeta + 0.5));
    18871887                        /*Nodal function 17*/
    1888                         dbasis[NUMNODESP2xP4*0+16] = - (8./6.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
    1889                         dbasis[NUMNODESP2xP4*1+16] =  (8./3.) * (SQRT3/6.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
    1890                         dbasis[NUMNODESP2xP4*2+16] =  -(8./3.)*gauss->coord2 *((2.*zeta -1. )*(zeta-0.5)*(zeta +1.) +  zeta *(zeta -1.)*( 2.*zeta + 0.5));
     1888                        dbasis[NUMNODESP2xP4*0+16] = - (8./3.) *(2*gauss->coord2 - 0.5 )*(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1889                        dbasis[NUMNODESP2xP4*1+16] = - (8./3.)*(-((2.*SQRT3)/(3.)) *gauss->coord2 + (SQRT3/6.))*(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1890                        dbasis[NUMNODESP2xP4*2+16] = - (8./3.)*gauss->coord2*(2.*gauss->coord2-1) *((2.*zeta -1. )*(zeta-0.5)*(zeta +1.) +  zeta *(zeta -1.)*( 2.*zeta + 0.5));
    18911891                        /*Nodal function 18*/
    18921892                        dbasis[NUMNODESP2xP4*0+17] = 0. ;
    1893                         dbasis[NUMNODESP2xP4*1+17] =  (8./3.) *(SQRT3/3.) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
    1894                         dbasis[NUMNODESP2xP4*2+17] =  -(8./3.)*gauss->coord3 *((2.* zeta -1. ) *(zeta-0.5)* (zeta +1.) +  zeta *(zeta -1.)*( 2.*zeta + 0.5));
     1893                        dbasis[NUMNODESP2xP4*1+17] =  -(8./3.) *(((4.*SQRT3)/(3.))*gauss->coord3 - SQRT3/3. )*(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
     1894                        dbasis[NUMNODESP2xP4*2+17] =  -(8./3.)*gauss->coord3*(2*gauss->coord3-1) *((2.* zeta -1. ) *(zeta-0.5)* (zeta +1.) +  zeta *(zeta -1.)*( 2.*zeta + 0.5));
    18951895                        /*Nodal function 19*/
    1896                         dbasis[NUMNODESP2xP4*0+18] =  (4./3.) *(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
    1897                         dbasis[NUMNODESP2xP4*1+18] =  (8./3.)* (SQRT3/6.)*(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
    1898                         dbasis[NUMNODESP2xP4*2+18] =  -(8./3.)*gauss->coord1*((2.*zeta -1. ) *(zeta+0.5)* (zeta +1.) +  zeta* (zeta -1.)*( 2.*zeta + 3./2.));
     1896                        dbasis[NUMNODESP2xP4*0+18] =  -(8./3.)*(-2.* gauss->coord1 + 0.5 ) *(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1897                        dbasis[NUMNODESP2xP4*1+18] =  -(8./3.)* (-((2.*SQRT3)/(3.)) *gauss->coord1 + (SQRT3)/(6.) )*(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1898                        dbasis[NUMNODESP2xP4*2+18] =  -(8./3.)*gauss->coord1*(2.*gauss->coord1-1) *((2.*zeta -1. ) *(zeta+0.5)* (zeta +1.) +  zeta* (zeta -1.)*( 2.*zeta + 3./2.));
    18991899                        /*Nodal function 20*/
    1900                         dbasis[NUMNODESP2xP4*0+19] =  (8./6.) *(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
    1901                         dbasis[NUMNODESP2xP4*1+19] =  (8./3.)*(SQRT3/6) *(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
    1902                         dbasis[NUMNODESP2xP4*2+19] =  -(8./3.)*gauss->coord2*((2.*zeta -1. )*(zeta+0.5)*(zeta +1.) +  zeta*(zeta -1.)*( 2.*zeta + 3./2.));
     1900                        dbasis[NUMNODESP2xP4*0+19] =  -(8./3.)*(2*gauss->coord2 - 0.5 )*(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1901                        dbasis[NUMNODESP2xP4*1+19] =  -(8./3.)*(-((2.*SQRT3)/(3.)) *gauss->coord2 + (SQRT3/6.)) *(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
     1902                        dbasis[NUMNODESP2xP4*2+19] =  -(8./3.)*gauss->coord2*(2.*gauss->coord2-1)*((2.*zeta -1. )*(zeta+0.5)*(zeta +1.) +  zeta*(zeta -1.)*( 2.*zeta + 3./2.));
    19031903                        /*Nodal function 21*/
    19041904                        dbasis[NUMNODESP2xP4*0+20] = 0 ;
    1905                         dbasis[NUMNODESP2xP4*1+20] =  -(8./3.)* (SQRT3/3) *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. ) ;
    1906                         dbasis[NUMNODESP2xP4*2+20] =  -(8./3.)*gauss->coord3 *(2. *zeta -1. )*(zeta+0.5)*(zeta +1.) +  zeta*(zeta -1.)*( 2.*zeta + 3./2.);
     1905                        dbasis[NUMNODESP2xP4*1+20] =  -(8./3.)* (((4.*SQRT3)/(3.))*gauss->coord3 - SQRT3/3. ) *(zeta - 1.)*(zeta + 0.5)*(zeta)*(zeta +1. ) ;
     1906                        dbasis[NUMNODESP2xP4*2+20] =  -(8./3.)*gauss->coord3*(2*gauss->coord3-1) *((2. *zeta -1. )*(zeta+0.5)*(zeta +1.) +  zeta*(zeta -1.)*( 2.*zeta + 3./2.));
    19071907                        /*Nodal function 22*/
    19081908                        dbasis[NUMNODESP2xP4*0+21] = -(32./6.) *gauss->coord3 *(zeta - 1.)*(zeta - 0.5)*(zeta)*(zeta +1. );
     
    19221922                        dbasis[NUMNODESP2xP4*2+24] = 16.* gauss->coord2 * gauss->coord3* ( 4.* zeta*zeta*zeta - (5./2.) *zeta );
    19231923                        /*Nodal function 26*/
    1924                         dbasis[NUMNODESP2xP4*0+25] = 8. *gauss->coord3*(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. );
    1925                         dbasis[NUMNODESP2xP4*1+25] = -16.*((SQRT3/3.)*gauss->coord1- (SQRT3/6.)*gauss->coord3 )*(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. ) ;
     1924                        dbasis[NUMNODESP2xP4*0+25] = -8. *gauss->coord3*(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. );
     1925                        dbasis[NUMNODESP2xP4*1+25] = 16.*((SQRT3/3.)*gauss->coord1- (SQRT3/6.)*gauss->coord3 )*(zeta - 1.)*(zeta - 0.5)*(zeta+0.5)*(zeta +1. ) ;
    19261926                        dbasis[NUMNODESP2xP4*2+25] = 16. * gauss->coord1 * gauss->coord3 *( 4. *zeta*zeta*zeta - (5./2.)* zeta );
    19271927                        /*Nodal function 27*/
     
    19411941                        dbasis[NUMNODESP2xP4*1+29] = -(32./3.)*(- (SQRT3/6.)*gauss->coord2- (SQRT3/6.)*gauss->coord1 ) *(zeta - 1.)*(zeta)*(zeta+0.5)*(zeta +1. ) ;
    19421942                        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.)));
    1943                         return;
     1943                        return;                 return;
    19441944                case P1xP3Enum :
    19451945                        /*Nodal function 1*/
Note: See TracChangeset for help on using the changeset viewer.