Actual source code: matstashspace.c

  1: #include <petsc/private/matimpl.h>

  3: /* Get new PetscMatStashSpace into the existing space */
  4: PetscErrorCode PetscMatStashSpaceGet(PetscInt bs2, PetscInt n, PetscMatStashSpace *space)
  5: {
  6:   PetscMatStashSpace a;

  8:   PetscFunctionBegin;
  9:   if (!n) PetscFunctionReturn(PETSC_SUCCESS);

 11:   PetscCall(PetscMalloc(sizeof(struct _MatStashSpace), &a));
 12:   PetscCall(PetscMalloc3(n * bs2, &a->space_head, n, &a->idx, n, &a->idy));

 14:   a->val              = a->space_head;
 15:   a->local_remaining  = n;
 16:   a->local_used       = 0;
 17:   a->total_space_size = 0;
 18:   a->next             = NULL;

 20:   if (*space) {
 21:     (*space)->next      = a;
 22:     a->total_space_size = (*space)->total_space_size;
 23:   }
 24:   a->total_space_size += n;
 25:   *space = a;
 26:   PetscFunctionReturn(PETSC_SUCCESS);
 27: }

 29: /* Copy the values in space into arrays val, idx and idy. Then destroy space */
 30: PetscErrorCode PetscMatStashSpaceContiguous(PetscInt bs2, PetscMatStashSpace *space, PetscScalar *val, PetscInt *idx, PetscInt *idy)
 31: {
 32:   PetscMatStashSpace a;

 34:   PetscFunctionBegin;
 35:   while (*space) {
 36:     a = (*space)->next;
 37:     PetscCall(PetscArraycpy(val, (*space)->val, (*space)->local_used * bs2));
 38:     val += bs2 * (*space)->local_used;
 39:     PetscCall(PetscArraycpy(idx, (*space)->idx, (*space)->local_used));
 40:     idx += (*space)->local_used;
 41:     PetscCall(PetscArraycpy(idy, (*space)->idy, (*space)->local_used));
 42:     idy += (*space)->local_used;

 44:     PetscCall(PetscFree3((*space)->space_head, (*space)->idx, (*space)->idy));
 45:     PetscCall(PetscFree(*space));
 46:     *space = a;
 47:   }
 48:   PetscFunctionReturn(PETSC_SUCCESS);
 49: }

 51: PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace *space)
 52: {
 53:   PetscMatStashSpace a;

 55:   PetscFunctionBegin;
 56:   while (*space) {
 57:     a = (*space)->next;
 58:     PetscCall(PetscFree3((*space)->space_head, (*space)->idx, (*space)->idy));
 59:     PetscCall(PetscFree(*space));
 60:     *space = a;
 61:   }
 62:   *space = NULL;
 63:   PetscFunctionReturn(PETSC_SUCCESS);
 64: }