xref: /petsc/src/mat/utils/freespace.c (revision 30cb48eedbb25d9f81c36be9f48fe7f4d174065c)
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 
50*30cb48eeSHong Zhang /*
51*30cb48eeSHong Zhang   Copy a linket list obtained from matrix symbolic ilu or lu factorization into a contiguous array that enables
52*30cb48eeSHong Zhang   an efficient matrix solve.
53*30cb48eeSHong Zhang 
54*30cb48eeSHong Zhang    Input Parameters:
55*30cb48eeSHong Zhang +  head - linked list of column indices obtained from matrix symbolic ilu or lu factorization
56*30cb48eeSHong Zhang .  space - an allocated int arry with length nnz of factored matrix.
57*30cb48eeSHong Zhang .  n - order of the matrix
58*30cb48eeSHong Zhang .  bi - row pointer of factored matrix with length 2n+2. First n+1 entries are set based on the traditional layout of L and U matrices
59*30cb48eeSHong Zhang .
60*30cb48eeSHong Zhang -  bdiag - int array holding the number of nonzeros in each row of L matrix, excluding diagonals.
61*30cb48eeSHong Zhang 
62*30cb48eeSHong Zhang    Output Parameter:
63*30cb48eeSHong Zhang .  space - column indices are copied into this int array with contiguous layout of L and U
64*30cb48eeSHong Zhang    See MatILUFactorSymbolic_SeqAIJ_ilu0_newdatastruct() for detailed description.
65*30cb48eeSHong Zhang */
667a48dd6fSHong Zhang #undef __FUNCT__
6712b5cbf3SHong Zhang #define __FUNCT__ "PetscFreeSpaceContiguous_newdatastruct"
6812b5cbf3SHong Zhang PetscErrorCode PetscFreeSpaceContiguous_newdatastruct(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *bi,PetscInt *bdiag)
6912b5cbf3SHong Zhang {
7012b5cbf3SHong Zhang   PetscFreeSpaceList a;
7112b5cbf3SHong Zhang   PetscErrorCode     ierr;
7212b5cbf3SHong Zhang   PetscInt           row,nnz,*bj,*array,total;
7312b5cbf3SHong Zhang   PetscInt           nnzL,nnzU;
7412b5cbf3SHong Zhang 
7512b5cbf3SHong Zhang   PetscFunctionBegin;
7612b5cbf3SHong Zhang   bi[2*n+1] = bi[n];
77*30cb48eeSHong Zhang   row       = 0;
7812b5cbf3SHong Zhang   total     = 0;
7912b5cbf3SHong Zhang   nnzL  = bdiag[0];
8012b5cbf3SHong Zhang   while ((*head)!=NULL) {
8112b5cbf3SHong Zhang     total += (*head)->local_used;
8212b5cbf3SHong Zhang     array  = (*head)->array_head;
8312b5cbf3SHong Zhang 
84*30cb48eeSHong Zhang     while (bi[row+1] <= total && row < n){
85*30cb48eeSHong Zhang       /* copy array entries into bj for this row */
86*30cb48eeSHong Zhang       nnz  = bi[row+1] - bi[row];
87*30cb48eeSHong Zhang       /* set bi[row] for new datastruct */
88*30cb48eeSHong Zhang       if (row == 0 ){
89*30cb48eeSHong Zhang         bi[row] = 0;
9012b5cbf3SHong Zhang       } else {
91*30cb48eeSHong Zhang         bi[row] = bi[row-1] + nnzL; /* nnzL of previous row */
9212b5cbf3SHong Zhang       }
9312b5cbf3SHong Zhang 
9412b5cbf3SHong Zhang       /* L part */
95*30cb48eeSHong Zhang       nnzL = bdiag[row];
96*30cb48eeSHong Zhang       bj   = space+bi[row];
9712b5cbf3SHong Zhang       ierr = PetscMemcpy(bj,array,nnzL*sizeof(PetscInt));CHKERRQ(ierr);
9812b5cbf3SHong Zhang 
9912b5cbf3SHong Zhang       /* diagonal entry */
100*30cb48eeSHong Zhang       bdiag[row]        = bi[2*n-row+1]-1;
101*30cb48eeSHong Zhang       space[bdiag[row]] = row;
10212b5cbf3SHong Zhang 
10312b5cbf3SHong Zhang       /* U part */
10412b5cbf3SHong Zhang       nnzU        = nnz - nnzL;
105*30cb48eeSHong Zhang       bi[2*n-row] = bi[2*n-row+1] - nnzU;
10612b5cbf3SHong Zhang       nnzU --;      /* exclude diagonal */
107*30cb48eeSHong Zhang       bj   = space + bi[2*n-(row)];
10812b5cbf3SHong Zhang       ierr = PetscMemcpy(bj,array+nnzL+1,nnzU*sizeof(PetscInt));CHKERRQ(ierr);
10912b5cbf3SHong Zhang 
11012b5cbf3SHong Zhang       array += nnz;
11112b5cbf3SHong Zhang       row++;
11212b5cbf3SHong Zhang     }
11312b5cbf3SHong Zhang 
11412b5cbf3SHong Zhang     a     = (*head)->more_space;
11512b5cbf3SHong Zhang     ierr  = PetscFree((*head)->array_head);CHKERRQ(ierr);
11612b5cbf3SHong Zhang     ierr  = PetscFree(*head);CHKERRQ(ierr);
11712b5cbf3SHong Zhang     *head = a;
11812b5cbf3SHong Zhang   }
11912b5cbf3SHong Zhang   bi[n] = bi[n-1] + nnzL;
12012b5cbf3SHong Zhang   if (bi[n] != bi[n+1]) SETERRQ2(1,"bi[n] %d != bi[n+1] %d",bi[n],bi[n+1]);
12112b5cbf3SHong Zhang   PetscFunctionReturn(0);
12212b5cbf3SHong Zhang }
12312b5cbf3SHong Zhang 
12412b5cbf3SHong Zhang #undef __FUNCT__
125a1a86e44SBarry Smith #define __FUNCT__ "PetscFreeSpaceDestroy"
126a1a86e44SBarry Smith PetscErrorCode PetscFreeSpaceDestroy(PetscFreeSpaceList head)
1277a48dd6fSHong Zhang {
128a1a86e44SBarry Smith   PetscFreeSpaceList a;
1297a48dd6fSHong Zhang   PetscErrorCode     ierr;
1307a48dd6fSHong Zhang 
1317a48dd6fSHong Zhang   PetscFunctionBegin;
1327a48dd6fSHong Zhang   while ((head)!=NULL) {
1337a48dd6fSHong Zhang     a    = (head)->more_space;
1347a48dd6fSHong Zhang     ierr = PetscFree((head)->array_head);CHKERRQ(ierr);
1357a48dd6fSHong Zhang     ierr = PetscFree(head);CHKERRQ(ierr);
1367a48dd6fSHong Zhang     head = a;
1377a48dd6fSHong Zhang   }
1387a48dd6fSHong Zhang   PetscFunctionReturn(0);
1397a48dd6fSHong Zhang }
140