35 #include "qrm_common.h"
46 integer,
parameter :: qrm_busy_ = -huge(0)
47 integer,
parameter :: qrm_ready_ = 0
48 integer,
parameter :: qrm_activable_ = 1
49 integer,
parameter :: qrm_active_ = 2
50 integer,
parameter :: qrm_factorized_ = 3
51 integer,
parameter :: qrm_free_ = 4
52 integer,
parameter :: qrm_done_ = huge(0)
65 integer,
allocatable :: rows(:)
67 integer,
allocatable :: cols(:)
70 integer,
allocatable :: aiptr(:)
73 integer,
allocatable :: ajcn(:)
77 real(kind(1.d0)),
allocatable :: aval(:)
84 integer,
allocatable :: colmap(:)
90 integer,
allocatable :: rowmap(:)
94 integer,
allocatable :: stair(:)
97 real(kind(1.d0)),
allocatable :: front(:,:)
100 real(kind(1.d0)),
allocatable :: r(:)
103 real(kind(1.d0)),
allocatable :: h(:)
107 real(kind(1.d0)),
allocatable :: cb(:)
110 real(kind(1.d0)),
allocatable :: tau(:)
113 real(kind(1.d0)),
allocatable :: t(:,:)
116 real(kind(1.d0)),
allocatable :: b(:,:)
127 integer,
allocatable :: ptable(:)
156 integer(kind=8) :: rsize, hsize
159 integer(kind=omp_lock_kind) :: lock
178 integer(kind=omp_lock_kind) :: lock
181 logical :: ok=.false.
207 character(len=*),
parameter :: name=
'qrm_front_destroy'
228 call omp_destroy_lock(qrm_front%lock)
230 __qrm_check_ret(name,
'qrm_adelloc',9999)
240 if(err_act .eq. qrm_abort_)
then
262 character(len=*),
parameter :: name=
'qrm_fdata_destroy'
267 if(
allocated(qrm_fdata%front_list))
then
268 do i=1, qrm_fdata%nfronts
271 deallocate(qrm_fdata%front_list)
273 __qrm_check_ret(name,
'qrm_front_destroy',9999)
276 qrm_fdata%nfronts = 0
277 qrm_fdata%ok = .false.
287 if(err_act .eq. qrm_abort_)
then
This module contains all the error management routines and data.
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.
subroutine qrm_err_act_restore(err_act)
Restores the value of the qrm_err_act variable.
subroutine qrm_err_check()
This subroutine checks the errors stack. If something is found all the entries in the stack are poppe...
The data structure meant to store all the results of the factorization phase.
This module contains the definition of all the data related to the factorization phase.
subroutine dqrm_front_destroy(qrm_front)
Frees a qrm_front_type instance.
This type defines a data structure containing all the data related to a front.
This module implements the memory handling routines. Pretty mucch allocations and deallocations...
subroutine dqrm_fdata_destroy(qrm_fdata)
Destroys a dqrm_fdata_type instance.
subroutine qrm_err_act_save(err_act)
Saves a copy of the qrm_err_act variable.