1 2 #include <petsc-private/matimpl.h> 3 4 /* Get new PetscMatStashSpace into the existing space */ 5 #undef __FUNCT__ 6 #define __FUNCT__ "PetscMatStashSpaceGet" 7 PetscErrorCode PetscMatStashSpaceGet(PetscInt bs2,PetscInt n,PetscMatStashSpace *space) 8 { 9 PetscMatStashSpace a; 10 PetscErrorCode ierr; 11 12 PetscFunctionBegin; 13 if (!n) PetscFunctionReturn(0); 14 15 ierr = PetscMalloc(sizeof(struct _MatStashSpace),&a);CHKERRQ(ierr); 16 ierr = PetscMalloc3(n*bs2,&(a->space_head),n,&a->idx,n,&a->idy);CHKERRQ(ierr); 17 18 a->val = a->space_head; 19 a->local_remaining = n; 20 a->local_used = 0; 21 a->total_space_size = 0; 22 a->next = NULL; 23 24 if (*space) { 25 (*space)->next = a; 26 a->total_space_size = (*space)->total_space_size; 27 } 28 a->total_space_size += n; 29 *space = a; 30 PetscFunctionReturn(0); 31 } 32 33 /* Copy the values in space into arrays val, idx and idy. Then destroy space */ 34 #undef __FUNCT__ 35 #define __FUNCT__ "PetscMatStashSpaceContiguous" 36 PetscErrorCode PetscMatStashSpaceContiguous(PetscInt bs2,PetscMatStashSpace *space,PetscScalar *val,PetscInt *idx,PetscInt *idy) 37 { 38 PetscMatStashSpace a; 39 PetscErrorCode ierr; 40 41 PetscFunctionBegin; 42 while ((*space) != NULL) { 43 a = (*space)->next; 44 ierr = PetscMemcpy(val,(*space)->val,((*space)->local_used*bs2)*sizeof(PetscScalar));CHKERRQ(ierr); 45 val += bs2*(*space)->local_used; 46 ierr = PetscMemcpy(idx,(*space)->idx,((*space)->local_used)*sizeof(PetscInt));CHKERRQ(ierr); 47 idx += (*space)->local_used; 48 ierr = PetscMemcpy(idy,(*space)->idy,((*space)->local_used)*sizeof(PetscInt));CHKERRQ(ierr); 49 idy += (*space)->local_used; 50 51 ierr = PetscFree3((*space)->space_head,(*space)->idx,(*space)->idy);CHKERRQ(ierr); 52 ierr = PetscFree(*space);CHKERRQ(ierr); 53 *space = a; 54 } 55 PetscFunctionReturn(0); 56 } 57 58 #undef __FUNCT__ 59 #define __FUNCT__ "PetscMatStashSpaceDestroy" 60 PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace *space) 61 { 62 PetscMatStashSpace a; 63 PetscErrorCode ierr; 64 65 PetscFunctionBegin; 66 while (*space) { 67 a = (*space)->next; 68 ierr = PetscFree3((*space)->space_head,(*space)->idx,(*space)->idy);CHKERRQ(ierr); 69 ierr = PetscFree((*space));CHKERRQ(ierr); 70 *space = a; 71 } 72 *space = NULL; 73 PetscFunctionReturn(0); 74 } 75