37 #include "qrm_common.h"
91 integer :: porder(:), parent(:), rc(:), stair(:)
94 integer ::
i, pnt, svar, f, ss
95 integer,
allocatable :: work(:)
98 character(len=*),
parameter :: name=
'qrm_compress_data'
104 if(rc(
i) .ne. 0) adata%nnodes = adata%nnodes+1
109 __qrm_check_ret(name,
'qrm_aalloc',9999)
120 if (rc(svar) .ne. 0)
then
121 adata%cp_ptr(pnt) =
i
126 adata%cp_ptr(adata%nnodes+1) = n+1
135 call
qrm_aalloc(adata%childptr, adata%nnodes+1)
137 __qrm_check_ret(name,
'qrm_aalloc',9999)
142 svar = porder(adata%cp_ptr(
i))
144 adata%rc(
i) = rc(svar)
146 adata%stair(
i) = stair(svar)
151 adata%parent(
i) = work(f)
152 adata%childptr(work(f)+1) = adata%childptr(work(f)+1)+1
156 adata%childptr(1) = 1
157 do i=2, adata%nnodes+1
158 adata%childptr(
i) = adata%childptr(
i-1)+adata%childptr(
i)
163 work(1:adata%nnodes)=0
167 adata%child(adata%childptr(f)+work(f)) =
i
173 adata%icperm(adata%cperm(
i)) =
i
178 __qrm_check_ret(name,
'qrm_adealloc',9999)
185 if(err_act .eq. qrm_abort_)
then
Generic interface for the qrm_adealloc_i, qrm_adealloc_2i, qrm_adealloc_s, qrm_adealloc_2s, qrm_adealloc_3s, qrm_adealloc_d, qrm_adealloc_2d, qrm_adealloc_3d, qrm_adealloc_c, qrm_adealloc_2c, qrm_adealloc_3c, qrm_adealloc_z, qrm_adealloc_2z, qrm_adealloc_3z, routines.
This module contains the interfaces of all non-typed routines.
subroutine qrm_err_act_restore(err_act)
Restores the value of the qrm_err_act variable.
This module contains the definition of the analysis data type.
subroutine qrm_err_check()
This subroutine checks the errors stack. If something is found all the entries in the stack are poppe...
The main data type for the analysis phase.
Generic interface for the qrm_aalloc_i, qrm_aalloc_2i, qrm_aalloc_s, qrm_aalloc_2s, qrm_aalloc_3s, qrm_aalloc_d, qrm_aalloc_2d, qrm_aalloc_3d, qrm_aalloc_c, qrm_aalloc_2c, qrm_aalloc_3c, qrm_aalloc_z, qrm_aalloc_2z, qrm_aalloc_3z, routines.
Generic interface for the ::qrm_compress_data routine.
This module implements the memory handling routines. Pretty mucch allocations and deallocations...
subroutine qrm_compress_data(adata, porder, parent, rc, stair, n)
This routine compresses the results of a number of operations in the analysis phase. Basically, the input data is of size n and the output of size adatannodes (which is the number of nodes in the elimination tree).
subroutine qrm_err_act_save(err_act)
Saves a copy of the qrm_err_act variable.