Ignore:
Timestamp:
04/16/12 14:57:18 (13 years ago)
Author:
Mathieu Morlighem
Message:

merged trunk-jpl and trunk for revision 11994M

Location:
issm/trunk
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • issm/trunk

  • issm/trunk/src/c/modules/SystemMatricesx/SystemMatricesx.cpp

    r11527 r11995  
    1010#include "../../EnumDefinitions/EnumDefinitions.h"
    1111
    12 void SystemMatricesx(Mat* pKff, Mat* pKfs, Vec* ppf, Vec* pdf, double* pkmax,Elements* elements,Nodes* nodes, Vertices* vertices,Loads* loads,Materials* materials, Parameters* parameters,bool kflag,bool pflag,bool penalty_kflag,bool penalty_pflag){
     12void SystemMatricesx(Matrix** pKff, Matrix** pKfs, Vector** ppf, Vector** pdf, double* pkmax,Elements* elements,Nodes* nodes, Vertices* vertices,Loads* loads,Materials* materials, Parameters* parameters,bool kflag,bool pflag,bool penalty_kflag,bool penalty_pflag){
    1313       
    1414        /*intermediary: */
     
    2121       
    2222        /*output: */
    23         Mat    Kff  = NULL;
    24         Mat    Kfs  = NULL;
    25         Vec    pf   = NULL;
    26         Vec    df=NULL;
     23        Matrix*    Kff  = NULL;
     24        Matrix*    Kfs  = NULL;
     25        Vector*    pf   = NULL;
     26        Vector*    df=NULL;
    2727        double kmax = 0;
    2828
     
    4949        if(kflag){
    5050
    51                 Kff=NewMat(fsize,fsize,connectivity,numberofdofspernode);
    52                 Kfs=NewMat(fsize,ssize,connectivity,numberofdofspernode);
    53                 df=NewVec(fsize);
     51                Kff=new Matrix(fsize,fsize,connectivity,numberofdofspernode);
     52                Kfs=new Matrix(fsize,ssize,connectivity,numberofdofspernode);
     53                df=new Vector(fsize);
    5454
    5555                /*Fill stiffness matrix from elements: */
     
    6666
    6767                /*Assemble matrix and doftypes and compress matrix to save memory: */
    68                 MatAssemblyBegin(Kff,MAT_FINAL_ASSEMBLY);
    69                 MatAssemblyEnd(Kff,MAT_FINAL_ASSEMBLY);
    70                 #if _PETSC_MAJOR_ == 2
    71                 MatCompress(Kff);
    72                 #endif
    73 
    74                 MatAssemblyBegin(Kfs,MAT_FINAL_ASSEMBLY);
    75                 MatAssemblyEnd(Kfs,MAT_FINAL_ASSEMBLY);
    76                 #if _PETSC_MAJOR_ == 2
    77                 MatCompress(Kfs);
    78                 #endif
    79                 VecAssemblyBegin(df);
    80                 VecAssemblyEnd(df);
    81                
     68                Kff->Assemble();
     69                Kfs->Assemble();
     70                df->Assemble();
    8271        }
    8372       
    8473        if(pflag){
    8574
    86                 pf=NewVec(fsize);
     75                pf=new Vector(fsize);
    8776
    8877                /*Fill right hand side vector, from elements: */
     
    9786                        if (load->InAnalysis(configuration_type)) load->CreatePVector(pf);
    9887                }
    99 
    100                 VecAssemblyBegin(pf);
    101                 VecAssemblyEnd(pf);
     88                pf->Assemble();
    10289        }
    10390
    10491        /*Now, figure out maximum value of K_gg, so that we can penalize it correctly: */
    105         MatNorm(Kff,NORM_INFINITY,&kmax);
     92        kmax=Kff->Norm(NORM_INF);
    10693
    10794        /*Now, deal with penalties*/
     
    115102
    116103                /*Assemble matrix and compress matrix to save memory: */
    117                 MatAssemblyBegin(Kff,MAT_FINAL_ASSEMBLY);
    118                 MatAssemblyEnd(Kff,MAT_FINAL_ASSEMBLY);
    119                 #if _PETSC_MAJOR_ == 2
    120                 MatCompress(Kff);
    121                 #endif
    122 
    123                 MatAssemblyBegin(Kfs,MAT_FINAL_ASSEMBLY);
    124                 MatAssemblyEnd(Kfs,MAT_FINAL_ASSEMBLY);
    125                 #if _PETSC_MAJOR_ == 2
    126                 MatCompress(Kfs);
    127                 #endif
     104                Kff->Assemble();
     105                Kfs->Assemble();
    128106        }
    129107
     
    137115                }
    138116
    139                 VecAssemblyBegin(pf);
    140                 VecAssemblyEnd(pf);
     117                pf->Assemble();
    141118        }
    142119
    143120        /*Assign output pointers: */
    144121        if(pKff) *pKff=Kff;
    145         else      MatFree(&Kff);
     122        else      xdelete(&Kff);
    146123        if(pKfs) *pKfs=Kfs;
    147         else      MatFree(&Kfs);
     124        else      xdelete(&Kfs);
    148125        if(ppf)  *ppf=pf;
    149         else      VecFree(&pf);
     126        else      xdelete(&pf);
    150127        if(pdf)  *pdf=df;
    151         else      VecFree(&df);
     128        else      xdelete(&df);
    152129        if(pkmax) *pkmax=kmax;
    153130}
Note: See TracChangeset for help on using the changeset viewer.