Changeset 5895
- Timestamp:
- 09/20/10 07:39:10 (15 years ago)
- Location:
- issm/trunk/src/c
- Files:
-
- 7 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk/src/c/modules/Reduceloadfromgtofx/Reduceloadfromgtofx.cpp
r5787 r5895 25 25 int verbose; 26 26 bool kffpartition=false; 27 bool fromlocalsize=true; 27 28 28 29 parameters->FindParam(&verbose,VerboseEnum); … … 47 48 /*pf = pf - Kfs * y_s;*/ 48 49 MatGetLocalSize(Kfs,&Kfsm,&Kfsn); 49 Kfsy_s=NewVec FromLocalSize(Kfsm);50 Kfsy_s=NewVec(Kfsm,fromlocalsize); 50 51 if (flag_ys0){ 51 52 -
issm/trunk/src/c/modules/Reduceloadx/Reduceloadx.cpp
r5772 r5895 18 18 int Kfsm,Kfsn; 19 19 PetscScalar a; 20 bool fromlocalsize=true; 20 21 21 22 if(pf){ … … 23 24 /*pf = pf - Kfs * y_s;*/ 24 25 MatGetLocalSize(Kfs,&Kfsm,&Kfsn); 25 Kfsy_s=NewVec FromLocalSize(Kfsm);26 Kfsy_s=NewVec(Kfsm,fromlocalsize); 26 27 if (flag_ys0){ 27 28 -
issm/trunk/src/c/modules/SystemMatricesx/SystemMatricesx.cpp
r5782 r5895 53 53 if(kflag){ 54 54 55 if(!buildkff)Kgg=NewMat(gsize,gsize, NULL,&connectivity,&numberofdofspernode);55 if(!buildkff)Kgg=NewMat(gsize,gsize,connectivity,numberofdofspernode); 56 56 else{ 57 Kff=NewMat(fsize,fsize, NULL,&connectivity,&numberofdofspernode);58 Kfs=NewMat(fsize,ssize, NULL,&connectivity,&numberofdofspernode);57 Kff=NewMat(fsize,fsize,connectivity,numberofdofspernode); 58 Kfs=NewMat(fsize,ssize,connectivity,numberofdofspernode); 59 59 } 60 60 -
issm/trunk/src/c/objects/Params/PetscMatParam.cpp
r5103 r5895 149 149 serial_mat=(double*)xmalloc(M*N*sizeof(double)); 150 150 memcpy(serial_mat,marshalled_dataset,M*N*sizeof(double));marshalled_dataset+=(M*N*sizeof(double)); 151 value=NewMat(M,N, &sparsity,NULL,NULL);151 value=NewMat(M,N,sparsity); 152 152 idxm=(int*)xmalloc(M*sizeof(int)); 153 153 idxn=(int*)xmalloc(N*sizeof(int)); -
issm/trunk/src/c/toolkits/petsc/patches/MatInvert.cpp
r3775 r5895 32 32 33 33 /*Create identitiy matrix: */ 34 identity=NewMat(M,N, &sparsity,NULL,NULL);34 identity=NewMat(M,N,sparsity); 35 35 diagonal=NewVec(M); 36 36 VecSet(diagonal,1.0); -
issm/trunk/src/c/toolkits/petsc/patches/NewMat.cpp
r5890 r5895 77 77 int d_nz,o_nz; 78 78 int nnz; 79 MatType type; 80 79 81 80 82 /*Determine local sizes: */ … … 86 88 o_nz=(int)connectivity*numberofdofspernode/2; 87 89 88 MatCreateMPIAIJ(MPI_COMM_WORLD,m,n,M,N,d_nz,NULL,o_nz,NULL,&outmatrix); 90 MatCreate(MPI_COMM_WORLD,&outmatrix); 91 MatSetSizes(outmatrix,m,n,M,N); 92 MatSetFromOptions(outmatrix); 93 94 /*preallocation according to type: */ 95 MatGetType(outmatrix,&type); 96 if(strcmp(type,"mpiaij")==0){ 97 MatMPIAIJSetPreallocation(outmatrix,d_nz,NULL,o_nz,NULL); 98 } 89 99 90 100 return outmatrix; -
issm/trunk/src/c/toolkits/petsc/patches/petscpatches.h
r5890 r5895 22 22 23 23 int PetscDetermineLocalSize(int global_size); 24 Vec NewVec(int size,bool fromlocalsize );24 Vec NewVec(int size,bool fromlocalsize=false); 25 25 Mat NewMat(int M,int N); 26 26 Mat NewMat(int M,int N,double sparsity);
Note:
See TracChangeset
for help on using the changeset viewer.