xref: /petsc/src/mat/utils/freespace.c (revision 12b5cbf38412732d5e37f9780c28bef77aed689b)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
37c4f633dSBarry Smith #include "../src/mat/utils/freespace.h"
470f19b1fSKris Buschelman 
570f19b1fSKris Buschelman #undef __FUNCT__
6a1a86e44SBarry Smith #define __FUNCT__ "PetscFreeSpaceGet"
7a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceGet(PetscInt n,PetscFreeSpaceList *list)
82e111b49SBarry Smith {
9a1a86e44SBarry Smith   PetscFreeSpaceList a;
10dfbe8321SBarry Smith   PetscErrorCode     ierr;
1170f19b1fSKris Buschelman 
1270f19b1fSKris Buschelman   PetscFunctionBegin;
13a1a86e44SBarry Smith   ierr = PetscMalloc(sizeof(struct _Space),&a);CHKERRQ(ierr);
142e111b49SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&(a->array_head));CHKERRQ(ierr);
1570f19b1fSKris Buschelman   a->array            = a->array_head;
162e111b49SBarry Smith   a->local_remaining  = n;
1770f19b1fSKris Buschelman   a->local_used       = 0;
1870f19b1fSKris Buschelman   a->total_array_size = 0;
1970f19b1fSKris Buschelman   a->more_space       = NULL;
2070f19b1fSKris Buschelman 
2170f19b1fSKris Buschelman   if (*list) {
2270f19b1fSKris Buschelman     (*list)->more_space = a;
2370f19b1fSKris Buschelman     a->total_array_size = (*list)->total_array_size;
2470f19b1fSKris Buschelman   }
2570f19b1fSKris Buschelman 
262e111b49SBarry Smith   a->total_array_size += n;
2770f19b1fSKris Buschelman   *list               =  a;
2870f19b1fSKris Buschelman   PetscFunctionReturn(0);
2970f19b1fSKris Buschelman }
3070f19b1fSKris Buschelman 
3170f19b1fSKris Buschelman #undef __FUNCT__
32a1a86e44SBarry Smith #define __FUNCT__ "PetscFreeSpaceContiguous"
33a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceContiguous(PetscFreeSpaceList *head,PetscInt *space)
3438baddfdSBarry Smith {
35a1a86e44SBarry Smith   PetscFreeSpaceList a;
36dfbe8321SBarry Smith   PetscErrorCode     ierr;
3770f19b1fSKris Buschelman 
3870f19b1fSKris Buschelman   PetscFunctionBegin;
3970f19b1fSKris Buschelman   while ((*head)!=NULL) {
4070f19b1fSKris Buschelman     a     =  (*head)->more_space;
412e111b49SBarry Smith     ierr  =  PetscMemcpy(space,(*head)->array_head,((*head)->local_used)*sizeof(PetscInt));CHKERRQ(ierr);
4270f19b1fSKris Buschelman     space += (*head)->local_used;
4370f19b1fSKris Buschelman     ierr  =  PetscFree((*head)->array_head);CHKERRQ(ierr);
4470f19b1fSKris Buschelman     ierr  =  PetscFree(*head);CHKERRQ(ierr);
4570f19b1fSKris Buschelman     *head =  a;
4670f19b1fSKris Buschelman   }
4770f19b1fSKris Buschelman   PetscFunctionReturn(0);
4870f19b1fSKris Buschelman }
497a48dd6fSHong Zhang 
507a48dd6fSHong Zhang #undef __FUNCT__
51*12b5cbf3SHong Zhang #define __FUNCT__ "PetscFreeSpaceContiguous_newdatastruct"
52*12b5cbf3SHong Zhang PetscErrorCode PetscFreeSpaceContiguous_newdatastruct(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *bi,PetscInt *bdiag)
53*12b5cbf3SHong Zhang {
54*12b5cbf3SHong Zhang   PetscFreeSpaceList a;
55*12b5cbf3SHong Zhang   PetscErrorCode     ierr;
56*12b5cbf3SHong Zhang   PetscInt           row,nnz,*bj,*array,total;
57*12b5cbf3SHong Zhang   PetscInt           nnzL,nnzU;
58*12b5cbf3SHong Zhang 
59*12b5cbf3SHong Zhang   PetscFunctionBegin;
60*12b5cbf3SHong Zhang   bi[2*n+1] = bi[n];
61*12b5cbf3SHong Zhang   row   = 1;
62*12b5cbf3SHong Zhang   total = 0;
63*12b5cbf3SHong Zhang   nnzL  = bdiag[0];
64*12b5cbf3SHong Zhang   while ((*head)!=NULL) {
65*12b5cbf3SHong Zhang     total += (*head)->local_used;
66*12b5cbf3SHong Zhang     array  = (*head)->array_head;
67*12b5cbf3SHong Zhang 
68*12b5cbf3SHong Zhang     while (bi[row] <= total && row <=n){
69*12b5cbf3SHong Zhang       /* copy array entries into bj for row-1 */
70*12b5cbf3SHong Zhang       nnz  = bi[row] - bi[row-1];
71*12b5cbf3SHong Zhang       /* set bi[row-1] for new datastruct */
72*12b5cbf3SHong Zhang       if (row -1 <= 1 ){
73*12b5cbf3SHong Zhang         bi[row -1] = 0;
74*12b5cbf3SHong Zhang       } else {
75*12b5cbf3SHong Zhang         bi[row-1] = bi[row-2] + nnzL; /* nnzL of previous row */
76*12b5cbf3SHong Zhang       }
77*12b5cbf3SHong Zhang 
78*12b5cbf3SHong Zhang       /* L part */
79*12b5cbf3SHong Zhang       nnzL = bdiag[row-1];
80*12b5cbf3SHong Zhang       bj   = space+bi[row-1];
81*12b5cbf3SHong Zhang       ierr = PetscMemcpy(bj,array,nnzL*sizeof(PetscInt));CHKERRQ(ierr);
82*12b5cbf3SHong Zhang 
83*12b5cbf3SHong Zhang       /* diagonal entry */
84*12b5cbf3SHong Zhang       bdiag[row-1]        = bi[2*n-(row-1)+1]-1;
85*12b5cbf3SHong Zhang       space[bdiag[row-1]] = row-1;
86*12b5cbf3SHong Zhang 
87*12b5cbf3SHong Zhang       /* U part */
88*12b5cbf3SHong Zhang       nnzU = nnz - nnzL;
89*12b5cbf3SHong Zhang       bi[2*n-(row-1)] = bi[2*n-(row-1)+1] - nnzU;
90*12b5cbf3SHong Zhang       nnzU --; /* exclude diagonal */
91*12b5cbf3SHong Zhang       bj  = space + bi[2*n-(row-1)];
92*12b5cbf3SHong Zhang       ierr = PetscMemcpy(bj,array+nnzL+1,nnzU*sizeof(PetscInt));CHKERRQ(ierr);
93*12b5cbf3SHong Zhang 
94*12b5cbf3SHong Zhang       array += nnz;
95*12b5cbf3SHong Zhang       row++;
96*12b5cbf3SHong Zhang     }
97*12b5cbf3SHong Zhang 
98*12b5cbf3SHong Zhang     a     =  (*head)->more_space;
99*12b5cbf3SHong Zhang     ierr  =  PetscFree((*head)->array_head);CHKERRQ(ierr);
100*12b5cbf3SHong Zhang     ierr  =  PetscFree(*head);CHKERRQ(ierr);
101*12b5cbf3SHong Zhang     *head =  a;
102*12b5cbf3SHong Zhang   }
103*12b5cbf3SHong Zhang   bi[n] = bi[n-1] + nnzL;
104*12b5cbf3SHong Zhang   if (bi[n] != bi[n+1]) SETERRQ2(1,"bi[n] %d != bi[n+1] %d",bi[n],bi[n+1]);
105*12b5cbf3SHong Zhang   PetscFunctionReturn(0);
106*12b5cbf3SHong Zhang }
107*12b5cbf3SHong Zhang 
108*12b5cbf3SHong Zhang #undef __FUNCT__
109a1a86e44SBarry Smith #define __FUNCT__ "PetscFreeSpaceDestroy"
110a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceDestroy(PetscFreeSpaceList head)
1117a48dd6fSHong Zhang {
112a1a86e44SBarry Smith   PetscFreeSpaceList a;
1137a48dd6fSHong Zhang   PetscErrorCode     ierr;
1147a48dd6fSHong Zhang 
1157a48dd6fSHong Zhang   PetscFunctionBegin;
1167a48dd6fSHong Zhang   while ((head)!=NULL) {
1177a48dd6fSHong Zhang     a    = (head)->more_space;
1187a48dd6fSHong Zhang     ierr = PetscFree((head)->array_head);CHKERRQ(ierr);
1197a48dd6fSHong Zhang     ierr = PetscFree(head);CHKERRQ(ierr);
1207a48dd6fSHong Zhang     head = a;
1217a48dd6fSHong Zhang   }
1227a48dd6fSHong Zhang   PetscFunctionReturn(0);
1237a48dd6fSHong Zhang }
124