CP2K 2.4 (Revision 12889)

dbcsr_index_operations.f90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------------!
00002 !   CP2K: A general program to perform molecular dynamics simulations         !
00003 !   Copyright (C) 2000 - 2013  CP2K developers group                          !
00004 !-----------------------------------------------------------------------------!
00005 
00006 ! *****************************************************************************
00015 MODULE dbcsr_index_operations
00016 
00017   USE array_types,                     ONLY: array_data,&
00018                                              array_size
00019   USE dbcsr_dist_operations,           ONLY: get_stored_canonical
00020   USE dbcsr_error_handling
00021   USE dbcsr_kinds,                     ONLY: int_4,&
00022                                              int_8
00023   USE dbcsr_methods,                   ONLY: &
00024        dbcsr_distribution, dbcsr_distribution_local_cols, &
00025        dbcsr_distribution_local_rows, dbcsr_distribution_ncols, &
00026        dbcsr_distribution_nlocal_cols, dbcsr_distribution_nlocal_rows, &
00027        dbcsr_distribution_nrows, dbcsr_distribution_thread_dist, &
00028        dbcsr_get_index_memory_type, dbcsr_nblkcols_total, &
00029        dbcsr_nblkrows_total, dbcsr_wm_use_mutable
00030   USE dbcsr_ptr_util,                  ONLY: default_resize_factor,&
00031                                              ensure_array_size,&
00032                                              memory_allocate,&
00033                                              memory_deallocate
00034   USE dbcsr_toollib,                   ONLY: joaat_hash,&
00035                                              sort,&
00036                                              swap
00037   USE dbcsr_types,                     ONLY: &
00038        dbcsr_distribution_obj, dbcsr_meta_size, dbcsr_num_slots, dbcsr_obj, &
00039        dbcsr_slot_blk_p, dbcsr_slot_col_i, dbcsr_slot_coo_l, &
00040        dbcsr_slot_home_prow, dbcsr_slot_home_vprow, &
00041        dbcsr_slot_nblkcols_local, dbcsr_slot_nblkcols_total, &
00042        dbcsr_slot_nblkrows_local, dbcsr_slot_nblkrows_total, &
00043        dbcsr_slot_nblks, dbcsr_slot_nfullcols_local, dbcsr_slot_nze, &
00044        dbcsr_slot_row_p, dbcsr_slot_size, dbcsr_slot_thr_c, dbcsr_type
00045 
00046   !$ USE OMP_LIB
00047 
00048   IMPLICIT NONE
00049 
00050   PRIVATE
00051 
00052   CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'dbcsr_index_operations'
00053 
00054 
00055   LOGICAL, PARAMETER :: careful_mod = .FALSE.
00056   LOGICAL, PARAMETER :: debug_mod = .FALSE.
00057 
00058 
00059 
00060 
00061   ! Index transformations
00062   PUBLIC :: dbcsr_make_index_canonical, make_index_triangular,&
00063             transpose_index_local, dbcsr_order_matrix_index
00064   PUBLIC :: dbcsr_make_index_local_row, dbcsr_has_local_row_index
00065   PUBLIC :: dbcsr_make_index_list
00066   ! Dense/Sparse
00067   PUBLIC :: make_dense_index, make_undense_index
00068   ! Working with DBCSR and linear indices
00069   PUBLIC :: dbcsr_make_dbcsr_index, dbcsr_sort_indices,&
00070             dbcsr_sort_many_indices, merge_index_arrays,&
00071             dbcsr_expand_row_index,&
00072             dbcsr_count_row_index, dbcsr_build_row_index,&
00073             dbcsr_index_prune_deleted,&
00074             dbcsr_index_compact
00075   PUBLIC :: dbcsr_index_checksum
00076   ! Index array manipulation
00077   PUBLIC :: dbcsr_addto_index_array, dbcsr_clearfrom_index_array,&
00078             dbcsr_repoint_index, dbcsr_make_index_exist
00079 
00080   INTERFACE dbcsr_count_row_index
00081      MODULE PROCEDURE dbcsr_count_row_index_copy,&
00082                       dbcsr_count_row_index_inplace
00083   END INTERFACE
00084 
00085   INTERFACE dbcsr_build_row_index
00086      MODULE PROCEDURE dbcsr_build_row_index_copy,&
00087                       dbcsr_build_row_index_inplace
00088   END INTERFACE
00089 
00090 CONTAINS
00091 
00092 ! *****************************************************************************
00102   SUBROUTINE make_index_canonical (new_row_p, new_col_i, new_blk_p,&
00103        old_row_p, old_col_i, old_blk_p, matrix)
00104     INTEGER, DIMENSION(:), INTENT(OUT)       :: new_row_p, new_col_i, 
00105                                                 new_blk_p
00106     INTEGER, DIMENSION(:), INTENT(IN)        :: old_row_p, old_col_i, 
00107                                                 old_blk_p
00108     TYPE(dbcsr_type), INTENT(IN)             :: matrix
00109 
00110     CHARACTER(len=*), PARAMETER :: routineN = 'make_index_canonical', 
00111       routineP = moduleN//':'//routineN
00112 
00113     INTEGER                                  :: blk, col, nblks, row, 
00114                                                 stored_col, stored_row
00115     INTEGER, ALLOCATABLE, DIMENSION(:)       :: row_i
00116     LOGICAL                                  :: tr
00117 
00118 !   ---------------------------------------------------------------------------
00119 
00120     nblks = SIZE(old_blk_p)
00121     ALLOCATE(row_i (nblks))
00122     IF (debug_mod) THEN
00123        WRITE(*,*)"old row_p", old_row_p
00124        WRITE(*,*)"old col_i", old_col_i
00125        WRITE(*,*)"old blk_p", old_blk_p
00126     ENDIF
00127     DO row = 1, SIZE(old_row_p)-1
00128        DO blk = old_row_p(row)+1, old_row_p(row+1)
00129           col = old_col_i(blk)
00130           stored_row = row
00131           stored_col = col
00132           tr = .FALSE.
00133           CALL get_stored_canonical (matrix, stored_row, stored_col, tr)
00134           IF (debug_mod) &
00135                WRITE(*,'(A,2(1X,I5),A,2(1X,I5),";",I7,1X,L1)')&
00136                routineN//" X->",row,col,"->",&
00137                stored_row, stored_col,blk,tr
00138           row_i(blk) = stored_row
00139           new_col_i(blk) = stored_col
00140           IF (.NOT. tr) THEN
00141              new_blk_p(blk) = old_blk_p(blk)
00142           ELSE
00143              new_blk_p(blk) = -old_blk_p(blk)
00144           ENDIF
00145        ENDDO
00146     ENDDO
00147     CALL dbcsr_sort_indices(nblks, row_i, new_col_i, blk_p=new_blk_p)
00148     ! Re-create the index
00149     CALL dbcsr_make_dbcsr_index (new_row_p, row_i, SIZE(new_row_p)-1, nblks)
00150     IF (debug_mod) THEN
00151        WRITE(*,*)"new row_p", new_row_p
00152        WRITE(*,*)"new row_i", row_i
00153        WRITE(*,*)"new col_i", new_col_i
00154        WRITE(*,*)"new blk_p", new_blk_p
00155     ENDIF
00156   END SUBROUTINE make_index_canonical
00157 
00158 
00159 ! *****************************************************************************
00169   SUBROUTINE make_index_triangular (new_row_p, new_col_i, new_blk_p,&
00170        old_row_p, old_col_i, old_blk_p, matrix)
00171     INTEGER, DIMENSION(:), INTENT(OUT)       :: new_row_p, new_col_i, 
00172                                                 new_blk_p
00173     INTEGER, DIMENSION(:), INTENT(IN)        :: old_row_p, old_col_i, 
00174                                                 old_blk_p
00175     TYPE(dbcsr_type), INTENT(IN)             :: matrix
00176 
00177     CHARACTER(len=*), PARAMETER :: routineN = 'make_index_triangular', 
00178       routineP = moduleN//':'//routineN
00179 
00180     INTEGER                                  :: blk, col, nblks, row, 
00181                                                 stored_col, stored_row
00182     INTEGER, ALLOCATABLE, DIMENSION(:)       :: row_i
00183     LOGICAL                                  :: tr
00184 
00185 !   ---------------------------------------------------------------------------
00186 
00187     nblks = SIZE(old_blk_p)
00188     ALLOCATE(row_i (nblks))
00189     IF (debug_mod) THEN
00190        WRITE(*,*)"old row_p", old_row_p
00191        WRITE(*,*)"old col_i", old_col_i
00192        WRITE(*,*)"old blk_p", old_blk_p
00193     ENDIF
00194     DO row = 1, SIZE(old_row_p)-1
00195        DO blk = old_row_p(row)+1, old_row_p(row+1)
00196           col = old_col_i(blk)
00197           stored_row = row
00198           stored_col = col
00199           tr = .FALSE.
00200           CALL get_stored_canonical (matrix, stored_row, stored_col, tr)
00201           IF (stored_row .GT. stored_col) THEN
00202              CALL swap (stored_row, stored_col)
00203              tr = .NOT. tr
00204           ENDIF
00205           IF (debug_mod) &
00206                WRITE(*,'(A,2(1X,I5),A,2(1X,I5),";",I7,1X,L1)')&
00207                routineN//" X->",row,col,"->",&
00208                stored_row, stored_col,blk,tr
00209           row_i(blk) = stored_row
00210           new_col_i(blk) = stored_col
00211           IF (.NOT. tr) THEN
00212              new_blk_p(blk) = old_blk_p(blk)
00213           ELSE
00214              new_blk_p(blk) = -old_blk_p(blk)
00215           ENDIF
00216        ENDDO
00217     ENDDO
00218     CALL dbcsr_sort_indices(nblks, row_i, new_col_i, blk_p=new_blk_p)
00219     ! Re-create the index
00220     CALL dbcsr_make_dbcsr_index (new_row_p, row_i, SIZE(new_row_p)-1, nblks)
00221     IF (debug_mod) THEN
00222        WRITE(*,*)"new row_p", new_row_p
00223        WRITE(*,*)"new row_i", row_i
00224        WRITE(*,*)"new col_i", new_col_i
00225        WRITE(*,*)"new blk_p", new_blk_p
00226     ENDIF
00227   END SUBROUTINE make_index_triangular
00228 
00229 ! *****************************************************************************
00232   PURE SUBROUTINE dbcsr_make_dbcsr_index (row_p, row_i, nrows, nblks)
00233     INTEGER, INTENT(in)                      :: nrows, nblks
00234     INTEGER, DIMENSION(1:nrows+1), 
00235       INTENT(out)                            :: row_p
00236     INTEGER, DIMENSION(1:nblks), INTENT(in)  :: row_i
00237 
00238     INTEGER                                  :: blk, row
00239 
00240 !
00241 !
00242 
00243     row_p(1) = 0
00244     row_p(nrows+1) = nblks
00245     row = 1
00246     blk = 1
00247     DO WHILE (row .LE. nrows)
00248        IF (blk .LE. nblks) THEN
00249           DO WHILE (row_i(blk) .EQ. row)
00250              blk = blk+1
00251              IF (blk .GT. nblks) THEN
00252                 row_p(row+1) = nblks-1
00253                 EXIT
00254              ENDIF
00255           ENDDO
00256        ENDIF
00257        row_p(row+1) = blk-1
00258        row = row+1
00259     ENDDO
00260   END SUBROUTINE dbcsr_make_dbcsr_index
00261 
00262 ! *****************************************************************************
00265   PURE SUBROUTINE dbcsr_expand_row_index (row_p, row_i, nrows, nblks)
00266     INTEGER, INTENT(in)                      :: nrows, nblks
00267     INTEGER, DIMENSION(1:nrows+1), 
00268       INTENT(IN)                             :: row_p
00269     INTEGER, DIMENSION(1:nblks), INTENT(out) :: row_i
00270 
00271     INTEGER                                  :: row
00272 
00273     DO row = 1, nrows
00274        row_i(row_p(row)+1 : row_p(row+1)) = row
00275     ENDDO
00276   END SUBROUTINE dbcsr_expand_row_index
00277 
00278 
00279 ! *****************************************************************************
00285   PURE SUBROUTINE dbcsr_count_row_index_inplace (rows, nrows)
00286     INTEGER, INTENT(IN)                      :: nrows
00287     INTEGER, DIMENSION(1:nrows+1), 
00288       INTENT(INOUT)                          :: rows
00289 
00290     INTEGER                                  :: row
00291 
00292     DO row = 1, nrows
00293        rows(row) = rows(row+1) - rows(row)
00294     ENDDO
00295     rows(nrows+1) = 0
00296   END SUBROUTINE dbcsr_count_row_index_inplace
00297 
00298 ! *****************************************************************************
00304   PURE SUBROUTINE dbcsr_count_row_index_copy (rows, counts, nrows)
00305     INTEGER, INTENT(IN)                      :: nrows
00306     INTEGER, DIMENSION(1:nrows), INTENT(OUT) :: counts
00307     INTEGER, DIMENSION(1:nrows+1), 
00308       INTENT(IN)                             :: rows
00309 
00310     INTEGER                                  :: row
00311 
00312     FORALL (row = 1:nrows)
00313        counts(row) = rows(row+1) - rows(row)
00314     END FORALL
00315   END SUBROUTINE dbcsr_count_row_index_copy
00316 
00317 ! *****************************************************************************
00323   PURE SUBROUTINE dbcsr_build_row_index_inplace (rows, nrows)
00324     INTEGER, INTENT(IN)                      :: nrows
00325     INTEGER, DIMENSION(1:nrows+1), 
00326       INTENT(INOUT)                          :: rows
00327 
00328     INTEGER                                  :: o, old_count, row
00329 
00330     old_count = rows(1)
00331     rows(1) = 0
00332     IF (nrows .GE. 1) THEN
00333        DO row = 2, nrows+1
00334           o = rows(row)
00335           rows(row) = rows(row-1) + old_count
00336           old_count = o
00337        ENDDO
00338     ENDIF
00339   END SUBROUTINE dbcsr_build_row_index_inplace
00340 
00341 ! *****************************************************************************
00348   PURE SUBROUTINE dbcsr_build_row_index_copy (counts, rows, nrows)
00349     INTEGER, INTENT(IN)                      :: nrows
00350     INTEGER, DIMENSION(1:nrows+1), 
00351       INTENT(OUT)                            :: rows
00352     INTEGER, DIMENSION(1:nrows), INTENT(IN)  :: counts
00353 
00354 !WTF?!rows(1) = 0
00355 !WTF?!IF (nrows .GE. 1) THEN
00356 !WTF?!   DO row = 2, nrows+1
00357 !WTF?!      rows(row) = rows(row-1) + counts(rows-1)
00358 !WTF?!   ENDDO
00359 !WTF?!ENDIF
00360 
00361     rows(1:nrows) = counts(1:nrows)
00362     CALL dbcsr_build_row_index_inplace (rows, nrows)
00363   END SUBROUTINE dbcsr_build_row_index_copy
00364 
00365 
00366 ! *****************************************************************************
00376   SUBROUTINE dbcsr_addto_index_array(matrix, slot, DATA, reservation, extra, error)
00377     TYPE(dbcsr_type), INTENT(INOUT)          :: matrix
00378     INTEGER, INTENT(IN)                      :: slot
00379     INTEGER, DIMENSION(:), INTENT(IN), 
00380       OPTIONAL                               :: DATA
00381     INTEGER, INTENT(IN), OPTIONAL            :: reservation, extra
00382     TYPE(dbcsr_error_type), INTENT(inout)    :: error
00383 
00384     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_addto_index_array', 
00385       routineP = moduleN//':'//routineN
00386 
00387     INTEGER                                  :: deplus, error_handler, space, 
00388                                                 ub, ub_new
00389 
00390 !   ---------------------------------------------------------------------------
00391 
00392     CALL dbcsr_error_set(routineN, error_handler, error)
00393     IF (debug_mod) THEN
00394         CALL dbcsr_assert (ASSOCIATED (matrix%index), dbcsr_fatal_level,&
00395              dbcsr_internal_error, routineP,&
00396              "Index must be preallocated.",__LINE__,error)
00397         CALL dbcsr_assert (UBOUND(matrix%index,1),'GE',dbcsr_num_slots,&
00398              dbcsr_failure_level, dbcsr_internal_error, routineP,&
00399              "Actual index size less than declared size",__LINE__,error)
00400         CALL dbcsr_assert(PRESENT(DATA), 'OR', PRESENT(reservation), dbcsr_fatal_level,&
00401              dbcsr_caller_error, routineP,&
00402              'Either an array or its size must be specified.',__LINE__,error)
00403         WRITE(*,*)routineP//' index', matrix%index(:dbcsr_num_slots)
00404     ENDIF
00405     IF (PRESENT (reservation)) THEN
00406        space = reservation
00407     ELSE
00408        space = SIZE(DATA)
00409     ENDIF
00410     IF (PRESENT (extra)) THEN
00411        deplus = extra
00412     ELSE
00413        deplus = 0
00414     ENDIF
00415     ub = UBOUND(matrix%index,1)
00417     IF (matrix%index(slot).EQ.0 .OR.&
00418          space.GT.matrix%index(slot+1)-matrix%index(slot)+1) THEN
00419        IF(debug_mod) WRITE(*,*)routineP//' Slot',slot,'not filled, adding at',&
00420             matrix%index(dbcsr_slot_size)+1,'sized',space
00421        matrix%index(slot) = matrix%index(dbcsr_slot_size)+1
00422        matrix%index(slot+1) = matrix%index(slot) + space - 1
00423        matrix%index(dbcsr_slot_size) = matrix%index(slot+1)
00424     ENDIF
00425     ! Shorten an index entry.
00426     IF (space .LT. matrix%index(slot+1) - matrix%index(slot)+1) THEN
00427        IF(debug_mod) WRITE(*,*)routineP//' Shortening index'
00428        matrix%index(slot+1) = matrix%index(slot) + space -1
00429        CALL dbcsr_repoint_index(matrix, slot)
00430     ENDIF
00431     ub_new = matrix%index(slot+1) + deplus
00432     IF(debug_mod) WRITE(*,*)routineP//' need',space,'at',matrix%index(slot),&
00433          'to',matrix%index(slot+1),'(',ub_new,')','have',ub
00434     IF (ub_new .GT. ub) THEN
00435        IF(debug_mod) WRITE(*,*)routineP//' Reallocating index to ubound', ub_new
00436        !CALL reallocate(matrix%index, 1, ub_new)
00437        CALL ensure_array_size(matrix%index, lb=1, ub=ub_new,&
00438             factor=default_resize_factor,&
00439             nocopy=.FALSE.,&
00440             memory_type=matrix%index_memory_type, error=error)
00441        CALL dbcsr_repoint_index(matrix)
00442     ENDIF
00443     IF(debug_mod) WRITE(*,*)routineP//' Adding slot',slot,'at',&
00444          matrix%index(slot),'sized',space
00445     CALL dbcsr_repoint_index(matrix, slot)
00446     IF (PRESENT(DATA)) &
00447          matrix%index(matrix%index(slot):matrix%index(slot+1)) = DATA(:)
00448     CALL dbcsr_error_stop(error_handler, error)
00449   END SUBROUTINE dbcsr_addto_index_array
00450 
00451 ! *****************************************************************************
00456   SUBROUTINE dbcsr_clearfrom_index_array(matrix, slot)
00457     TYPE(dbcsr_type), INTENT(INOUT)          :: matrix
00458     INTEGER, INTENT(IN)                      :: slot
00459 
00460     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_clearfrom_index_array', 
00461       routineP = moduleN//':'//routineN
00462 
00463     INTEGER                                  :: space
00464     INTEGER, DIMENSION(5)                    :: max_extents
00465     TYPE(dbcsr_error_type)                   :: error
00466 
00467 !   ---------------------------------------------------------------------------
00468 
00469     CALL dbcsr_assert (ASSOCIATED (matrix%index), dbcsr_fatal_level,&
00470          dbcsr_internal_error, routineP,&
00471          "Index must be preallocated.",__LINE__,error)
00472     CALL dbcsr_assert (UBOUND(matrix%index,1),'GE',dbcsr_num_slots,&
00473          dbcsr_failure_level, dbcsr_internal_error, routineP,&
00474          "Actual index size less than declared size",__LINE__,error)
00475     IF(debug_mod) WRITE(*,*)routineP//' index',&
00476          matrix%index(:dbcsr_num_slots)
00477     ! Clear index entry pointer
00478     matrix%index(slot) = 1
00479     matrix%index(slot+1) = 0
00480     CALL dbcsr_repoint_index(matrix, slot)
00481     ! Update the declared index size
00482     max_extents = (/ &
00483          matrix%index(dbcsr_slot_row_p+1),&
00484          matrix%index(dbcsr_slot_col_i+1),&
00485          matrix%index(dbcsr_slot_blk_p+1),&
00486          matrix%index(dbcsr_slot_thr_c+1),&
00487          matrix%index(dbcsr_slot_coo_l+1) /)
00488     space = MAX (MAXVAL (max_extents), dbcsr_num_slots)
00489     matrix%index (dbcsr_slot_size) = space
00490   END SUBROUTINE dbcsr_clearfrom_index_array
00491 
00492 
00493 
00494 ! *****************************************************************************
00499   SUBROUTINE dbcsr_repoint_index(m, slot)
00500     TYPE(dbcsr_type), INTENT(INOUT)          :: m
00501     INTEGER, INTENT(IN), OPTIONAL            :: slot
00502 
00503     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_repoint_index', 
00504       routineP = moduleN//':'//routineN
00505 
00506     INTEGER                                  :: s
00507     LOGICAL                                  :: all
00508 
00509 !   ---------------------------------------------------------------------------
00510 
00511     IF (m%nblks .NE. m%index(dbcsr_slot_nblks)) THEN
00512        m%nblks = m%index(dbcsr_slot_nblks)
00513        m%nze = m%index(dbcsr_slot_nze)
00514     ENDIF
00515     all = .TRUE.
00516     IF (PRESENT (slot)) THEN
00517        all = .FALSE.
00518        s = slot
00519     ELSE
00520        s = 0
00521     ENDIF
00522 
00523     IF (m%index(dbcsr_slot_row_p).GT.0.AND.all .OR.&
00524          s.EQ.dbcsr_slot_row_p) THEN
00525        IF (m%index(dbcsr_slot_row_p).GT.0) THEN
00526           m%row_p => m%index(m%index(dbcsr_slot_row_p):&
00527                &             m%index(dbcsr_slot_row_p+1))
00528        ELSE
00529           NULLIFY (m%row_p)
00530        ENDIF
00531     ENDIF
00532     IF (m%index(dbcsr_slot_col_i).GT.0.AND.all .OR.&
00533          s.EQ.dbcsr_slot_col_i) THEN
00534        IF (m%index(dbcsr_slot_col_i).GT.0) THEN
00535           m%col_i => m%index(m%index(dbcsr_slot_col_i):&
00536                &             m%index(dbcsr_slot_col_i+1))
00537        ELSE
00538           NULLIFY (m%col_i)
00539        ENDIF
00540     ENDIF
00541     IF (m%index(dbcsr_slot_blk_p).GT.0.AND.all .OR.&
00542          s.EQ.dbcsr_slot_blk_p) THEN
00543        IF (m%index(dbcsr_slot_blk_p).GT.0) THEN
00544           m%blk_p => m%index(m%index(dbcsr_slot_blk_p):&
00545                &             m%index(dbcsr_slot_blk_p+1))
00546        ELSE
00547           NULLIFY (m%blk_p)
00548        ENDIF
00549     ENDIF
00550     IF (m%index(dbcsr_slot_thr_c).GT.0.AND.all .OR.&
00551          s.EQ.dbcsr_slot_thr_c) THEN
00552        IF (m%index(dbcsr_slot_thr_c).GT.0) THEN
00553           m%thr_c => m%index(m%index(dbcsr_slot_thr_c):&
00554                &             m%index(dbcsr_slot_thr_c+1))
00555        ELSE
00556           NULLIFY (m%thr_c)
00557        ENDIF
00558     ENDIF
00559     IF (m%index(dbcsr_slot_coo_l).GT.0.AND.all .OR.&
00560          s.EQ.dbcsr_slot_coo_l) THEN
00561        IF (m%index(dbcsr_slot_coo_l).GT.0) THEN
00562           m%coo_l => m%index(m%index(dbcsr_slot_coo_l):&
00563                &             m%index(dbcsr_slot_coo_l+1))
00564        ELSE
00565           NULLIFY (m%coo_l)
00566        ENDIF
00567     ENDIF
00568     IF (all) THEN
00569        m%index(dbcsr_slot_nblks) = m%nblks
00570        m%index(dbcsr_slot_nze) = m%nze
00571     ENDIF
00572   END SUBROUTINE dbcsr_repoint_index
00573 
00574 
00575 ! *****************************************************************************
00577 
00581   SUBROUTINE dbcsr_make_index_exist(m, error)
00582     TYPE(dbcsr_type), INTENT(INOUT)          :: m
00583     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
00584 
00585     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_make_index_exist', 
00586       routineP = moduleN//':'//routineN
00587 
00588     INTEGER                                  :: error_handle
00589 
00590 !   ---------------------------------------------------------------------------
00591 
00592     CALL dbcsr_error_set (routineN, error_handle, error)
00593     CALL dbcsr_assert (ASSOCIATED (m%index),&
00594          dbcsr_fatal_level, dbcsr_internal_error, routineN,&
00595          "Index array does not yet exist.", __LINE__, error=error)
00596     IF (.NOT. ASSOCIATED (m%row_p)) THEN
00597        CALL dbcsr_addto_index_array (m, dbcsr_slot_row_p,&
00598             reservation=m%nblkrows_total+1, error=error)
00599        m%row_p(:) = 0
00600     ENDIF
00601     IF (.NOT. ASSOCIATED (m%col_i)) THEN
00602        CALL dbcsr_addto_index_array (m, dbcsr_slot_col_i,&
00603             reservation=0, error=error)
00604     ENDIF
00605     IF (.NOT. ASSOCIATED (m%blk_p)) THEN
00606        CALL dbcsr_addto_index_array (m, dbcsr_slot_blk_p,&
00607             reservation=0, error=error)
00608     ENDIF
00609     CALL dbcsr_repoint_index (m)
00610     CALL dbcsr_error_stop (error_handle, error)
00611   END SUBROUTINE dbcsr_make_index_exist
00612 
00613 
00614 
00615 ! *****************************************************************************
00629   SUBROUTINE dbcsr_sort_indices(n, row_i, col_i, blk_p, blk_d)
00630     INTEGER, INTENT(IN)                      :: n
00631     INTEGER, DIMENSION(1:n), INTENT(INOUT)   :: row_i, col_i
00632     INTEGER, DIMENSION(1:n), INTENT(INOUT), 
00633       OPTIONAL                               :: blk_p, blk_d
00634 
00635     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_sort_indices', 
00636       routineP = moduleN//':'//routineN
00637     INTEGER(KIND=int_8), PARAMETER           :: lmask8 = 4294967295_int_8
00638 
00639     INTEGER                                  :: error_handle, i, stat
00640     INTEGER(KIND=int_8), ALLOCATABLE, 
00641       DIMENSION(:)                           :: sort_keys
00642     INTEGER, ALLOCATABLE, DIMENSION(:)       :: buf, buf_d
00643     TYPE(dbcsr_error_type)                   :: error
00644 
00645 !   ---------------------------------------------------------------------------
00646 
00647     IF (SIZE (row_i) .EQ. 0) RETURN
00648 
00649     CALL dbcsr_error_set(routineN, error_handle, error)
00650 
00651     CALL dbcsr_assert(SIZE(row_i),'GE',n, dbcsr_failure_level,&
00652          dbcsr_caller_error, routineP, 'row_i too small',__LINE__,error)
00653     CALL dbcsr_assert(SIZE(col_i),'GE',n, dbcsr_failure_level,&
00654          dbcsr_caller_error, routineP, 'col_i too small',__LINE__,error)
00655     IF (PRESENT (blk_p)) CALL dbcsr_assert(SIZE(blk_p),'GE',n, dbcsr_failure_level,&
00656          dbcsr_caller_error, routineP, 'blk_p too small',__LINE__,error)
00657     IF (PRESENT (blk_p)) THEN
00658        ALLOCATE(buf(n), stat=stat)
00659        CALL dbcsr_assert(stat == 0, dbcsr_failure_level,&
00660             dbcsr_caller_error, routineP, 'buf',__LINE__,error)
00661        buf(1:n) = blk_p(1:n)
00662     ENDIF
00663     IF (PRESENT (blk_d)) THEN
00664        ALLOCATE(buf_d(n), stat=stat)
00665        CALL dbcsr_assert(stat == 0, dbcsr_failure_level,&
00666             dbcsr_caller_error, routineP, 'buf_d',__LINE__,error)
00667        buf_d(1:n) = blk_d(1:n)
00668     ENDIF
00672     ALLOCATE (sort_keys(n), stat=stat)
00673     CALL dbcsr_assert (stat, "EQ", 0, dbcsr_fatal_level, dbcsr_internal_error,&
00674          routineN, "Could not allocate memory for sorting buffer.", __LINE__,&
00675          error=error)
00676     sort_keys(:) = IOR (ISHFT(INT(row_i(:),int_8), 32), INT(col_i(:),int_8))
00677     IF (PRESENT (blk_p)) col_i(1:n) = (/ (i, i=1,n) /)
00679     CALL sort(sort_keys, n, col_i)
00680     ! Since blk_d is usually not present we can have two loops that
00681     ! are essentially the same.
00682     IF (PRESENT (blk_p)) THEN
00683        FORALL (i = 1:n)
00684           blk_p(i) = buf(col_i(i))
00685        END FORALL
00686        DEALLOCATE (buf)
00687     END IF
00688     IF (PRESENT (blk_d)) THEN
00689        FORALL (i = 1:n)
00690           blk_d(i) = buf_d(col_i(i))
00691        END FORALL
00692        DEALLOCATE (buf_d)
00693     ENDIF
00694     FORALL (i = 1:n)
00695        col_i(i) = INT (IAND(sort_keys(i), lmask8), int_4)
00696        row_i(i) = INT (ISHFT(sort_keys(i), -32), int_4)
00697     END FORALL
00698     DEALLOCATE (sort_keys)
00699     IF(debug_mod.AND.PRESENT(blk_p))&
00700          WRITE(*,*)routineP//' sort, blk_p =',blk_p
00701     CALL dbcsr_error_stop (error_handle, error)
00702 
00703   END SUBROUTINE dbcsr_sort_indices
00704 
00705 ! *****************************************************************************
00709   SUBROUTINE dbcsr_sort_many_indices(matrix)
00710     TYPE(dbcsr_type), INTENT(INOUT)          :: matrix
00711 
00712     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_sort_many_indices', 
00713       routineP = moduleN//':'//routineN
00714 
00715     INTEGER                                  :: i
00716 
00717 !   ---------------------------------------------------------------------------
00718 
00719     DO i = 1, SIZE (matrix%wms)
00720        IF (.NOT. dbcsr_wm_use_mutable (matrix%wms(i))) THEN
00721           IF (matrix%wms(i)%lastblk .GT. 0) THEN
00722              CALL dbcsr_sort_indices (matrix%wms(i)%lastblk,&
00723                   matrix%wms(i)%row_i,&
00724                   matrix%wms(i)%col_i, matrix%wms(i)%blk_p)
00725           ENDIF
00726        ENDIF
00727     ENDDO
00728   END SUBROUTINE dbcsr_sort_many_indices
00729 
00730 ! *****************************************************************************
00738   SUBROUTINE dbcsr_order_matrix_index(matrix, error)
00739     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
00740     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
00741 
00742     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_order_matrix_index', 
00743       routineP = moduleN//':'//routineN
00744 
00745     INTEGER                                  :: error_handle, first_blk, i, 
00746                                                 last_blk, ncols, row, stat
00747     INTEGER, ALLOCATABLE, DIMENSION(:)       :: buff, permutation
00748     INTEGER, DIMENSION(:), POINTER           :: blk_p, col_i, full_blk_p, 
00749                                                 full_col_i, full_row_p
00750 
00751 !   ---------------------------------------------------------------------------
00752 
00753     CALL dbcsr_error_set(routineN, error_handle, error)
00754 
00755     full_row_p => matrix%m%row_p
00756     full_col_i => matrix%m%col_i
00757     full_blk_p => matrix%m%blk_p
00758 !test!$OMP PARALLEL DEFAULT (none) &
00759 !test!$OMP PRIVATE (permutation, buff, stat, &
00760 !test!$             row, first_blk, last_blk, ncols, &
00761 !test!$             col_i, blk_p, i) &
00762 !test!$OMP SHARED (matrix, full_row_p, full_col_i, full_blk_p)
00763     ALLOCATE (permutation(dbcsr_nblkcols_total (matrix)), stat=stat)
00764     ALLOCATE (buff(dbcsr_nblkcols_total (matrix)), stat=stat)
00765     ! Goes through all rows and sorts the column in each row. It also
00766     ! reorders the blk_p correspondingly.
00767 !test!$OMP DO
00768     DO row = 1, dbcsr_nblkrows_total (matrix)
00769        first_blk = full_row_p(row)+1
00770        last_blk = full_row_p(row+1)
00771        ncols = last_blk - first_blk + 1
00772        col_i => full_col_i(first_blk : last_blk)
00773        blk_p => full_blk_p(first_blk : last_blk)
00774        permutation(1:ncols) = (/ (i, i=1,ncols) /)
00775        buff(1:ncols) = blk_p(1:ncols)
00776        CALL sort(col_i, ncols, permutation)
00777        DO i = 1, ncols
00778           blk_p(i) = buff(permutation(i))
00779        ENDDO
00780     ENDDO
00781 !test!OMP END DO NOWAIT
00782     DEALLOCATE (permutation, buff)
00783 !test!$OMP END PARALLEL
00784     CALL dbcsr_error_stop(error_handle, error)
00785   END SUBROUTINE dbcsr_order_matrix_index
00786 
00787 
00788 ! *****************************************************************************
00794   SUBROUTINE dbcsr_index_prune_deleted(matrix, error)
00795     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
00796     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
00797 
00798     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_index_prune_deleted', 
00799       routineP = moduleN//':'//routineN
00800 
00801     INTEGER                                  :: error_handle, nblks_max, 
00802                                                 new_blk, nrows, old_blk, row, 
00803                                                 stat
00804     INTEGER, ALLOCATABLE, DIMENSION(:)       :: new_blk_p, new_col_i, 
00805                                                 new_row_count
00806     INTEGER, DIMENSION(:), POINTER           :: old_blk_p, old_col_i, 
00807                                                 old_row_p
00808 
00809 !   ---------------------------------------------------------------------------
00810 
00811     CALL dbcsr_error_set(routineN, error_handle, error)
00812     !
00813     old_row_p => matrix%m%row_p
00814     old_col_i => matrix%m%col_i
00815     old_blk_p => matrix%m%blk_p
00816     !
00817     nrows = matrix%m%nblkrows_total
00818     nblks_max = old_row_p(nrows+1)
00819     ALLOCATE (new_row_count(nrows), stat=stat)
00820     CALL dbcsr_assert (stat, "EQ", 0, dbcsr_fatal_level, dbcsr_internal_error,&
00821          routineN, "Can not allocate memory.", __LINE__, error=error)
00822     ALLOCATE (new_col_i(nblks_max), stat=stat)
00823     CALL dbcsr_assert (stat, "EQ", 0, dbcsr_fatal_level, dbcsr_internal_error,&
00824          routineN, "Can not allocate memory.", __LINE__, error=error)
00825     ALLOCATE (new_blk_p(nblks_max), stat=stat)
00826     CALL dbcsr_assert (stat, "EQ", 0, dbcsr_fatal_level, dbcsr_internal_error,&
00827          routineN, "Can not allocate memory.", __LINE__, error=error)
00828     !
00829     ! Build up the new index from all non-deleted blocks in the
00830     ! existing index.
00831     new_blk = 0
00832     DO row = 1, nrows
00833        new_row_count(row) = 0
00834        DO old_blk = old_row_p(row)+1, old_row_p(row+1)
00835           IF (old_blk_p(old_blk) .GT. 0) THEN
00836              new_blk = new_blk + 1
00837              new_row_count(row) = new_row_count(row) + 1
00838              new_col_i(new_blk) = old_col_i(old_blk)
00839              new_blk_p(new_blk) = old_blk_p(old_blk)
00840           ENDIF
00841        ENDDO
00842     ENDDO
00843     !
00844     CALL dbcsr_clearfrom_index_array(matrix%m, dbcsr_slot_row_p)
00845     CALL dbcsr_clearfrom_index_array(matrix%m, dbcsr_slot_col_i)
00846     CALL dbcsr_clearfrom_index_array(matrix%m, dbcsr_slot_blk_p)
00847     CALL dbcsr_addto_index_array(matrix%m, dbcsr_slot_row_p,&
00848          reservation=nrows+1, extra=2*new_blk, error=error)
00849     old_row_p => matrix%m%row_p
00850     CALL dbcsr_build_row_index(counts=new_row_count, rows=old_row_p,&
00851          nrows=nrows)
00852     CALL dbcsr_addto_index_array(matrix%m, dbcsr_slot_col_i, DATA=new_col_i(1:new_blk),&
00853          error=error)
00854     CALL dbcsr_addto_index_array(matrix%m, dbcsr_slot_blk_p, DATA=new_blk_p(1:new_blk),&
00855          error=error)
00856     matrix%m%nblks = new_blk
00857     matrix%m%index(dbcsr_slot_nblks) = new_blk
00858     !
00859     DEALLOCATE (new_col_i, new_blk_p, new_row_count)
00860     !
00861     CALL dbcsr_error_stop(error_handle, error)
00862   END SUBROUTINE dbcsr_index_prune_deleted
00863 
00864 
00865 ! *****************************************************************************
00875   SUBROUTINE transpose_index_local (new_col_p, new_row_i, old_row_p,&
00876        old_col_i, new_blk_p, old_blk_p)
00877     INTEGER, DIMENSION(:), INTENT(OUT)       :: new_col_p, new_row_i
00878     INTEGER, DIMENSION(:), INTENT(IN)        :: old_row_p, old_col_i
00879     INTEGER, DIMENSION(:), INTENT(OUT), 
00880       OPTIONAL                               :: new_blk_p
00881     INTEGER, DIMENSION(:), INTENT(IN), 
00882       OPTIONAL                               :: old_blk_p
00883 
00884     CHARACTER(len=*), PARAMETER :: routineN = 'transpose_index_local', 
00885       routineP = moduleN//':'//routineN
00886 
00887     INTEGER                                  :: error_handle, nblks, 
00888                                                 ncols_new, nrows_old
00889     INTEGER, ALLOCATABLE, DIMENSION(:)       :: new_col_i
00890     LOGICAL                                  :: blks
00891     TYPE(dbcsr_error_type)                   :: error
00892 
00893 !   ---------------------------------------------------------------------------
00894 
00895     CALL dbcsr_error_set (routineN, error_handle, error)
00896     blks = PRESENT (new_blk_p) .AND. PRESENT (old_blk_p)
00897     nblks = SIZE (old_col_i)
00898     nrows_old = SIZE (old_row_p)-1
00899     ncols_new = SIZE (new_col_p)-1
00900     IF (blks) new_blk_p(:) = old_blk_p(:)
00901     ALLOCATE (new_col_i(nblks))
00902     CALL dbcsr_expand_row_index(old_row_p, new_row_i, nrows_old, nblks)
00903     new_col_i(:) = old_col_i(:)
00904     CALL dbcsr_sort_indices (nblks, new_col_i, new_row_i, new_blk_p)
00905     CALL dbcsr_make_dbcsr_index (new_col_p, new_col_i, ncols_new, nblks)
00906     DEALLOCATE(new_col_i)
00907     CALL dbcsr_error_stop (error_handle, error)
00908   END SUBROUTINE transpose_index_local
00909 
00910 
00911 ! *****************************************************************************
00918   SUBROUTINE dbcsr_make_index_canonical (matrix, cp2k)
00919     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
00920     LOGICAL, INTENT(IN), OPTIONAL            :: cp2k
00921 
00922     INTEGER                                  :: nb, nc, nr
00923     INTEGER, ALLOCATABLE, DIMENSION(:)       :: new_blk_p, new_col_i, 
00924                                                 new_row_p
00925     LOGICAL                                  :: rev
00926 
00927 !   ---------------------------------------------------------------------------
00928 
00929     rev = .FALSE.
00930     IF (PRESENT (cp2k)) rev = cp2k
00931     nr = SIZE(matrix%m%row_p)
00932     ALLOCATE (new_row_p (nr))
00933     nc = SIZE(matrix%m%col_i)
00934     ALLOCATE (new_col_i (nc))
00935     nb = SIZE(matrix%m%blk_p)
00936     ALLOCATE (new_blk_p (nb))
00937     IF (rev) THEN
00938        CALL make_index_triangular (new_row_p, new_col_i, new_blk_p,&
00939             matrix%m%row_p, matrix%m%col_i, matrix%m%blk_p, matrix%m)
00940     ELSE
00941        CALL make_index_canonical (new_row_p, new_col_i, new_blk_p,&
00942             matrix%m%row_p, matrix%m%col_i, matrix%m%blk_p, matrix%m)
00943     ENDIF
00944     matrix%m%row_p(:) = new_row_p
00945     matrix%m%col_i(:) = new_col_i
00946     matrix%m%blk_p(:) = new_blk_p
00947   END SUBROUTINE dbcsr_make_index_canonical
00948 
00949 
00950 ! *****************************************************************************
00965   SUBROUTINE make_dense_index (row_p, col_i, blk_p,&
00966        nblkrows_total, nblkcols_total, myblkrows, myblkcols,&
00967        row_blk_offsets, col_blk_offsets, meta, make_tr, error)
00968 
00969     !INTEGER, DIMENSION(:), INTENT(OUT)       :: row_p, col_i, blk_p
00970     INTEGER, INTENT(IN)                      :: nblkrows_total
00971     INTEGER, DIMENSION(:), INTENT(OUT)       :: blk_p, col_i
00972     INTEGER, DIMENSION(1:nblkrows_total+1), 
00973       INTENT(OUT)                            :: row_p
00974     INTEGER, INTENT(IN)                      :: nblkcols_total
00975     INTEGER, DIMENSION(:), INTENT(IN)        :: myblkrows, myblkcols, 
00976                                                 row_blk_offsets, 
00977                                                 col_blk_offsets
00978     INTEGER, DIMENSION(dbcsr_meta_size), 
00979       INTENT(INOUT)                          :: meta
00980     LOGICAL, INTENT(IN), OPTIONAL            :: make_tr
00981     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
00982 
00983     CHARACTER(len=*), PARAMETER :: routineN = 'make_dense_index', 
00984       routineP = moduleN//':'//routineN
00985 
00986     INTEGER                                  :: blk, c, col_l, mynblkcols, 
00987                                                 mynblkrows, nblks, nze, 
00988                                                 prev_row, row, row_l, 
00989                                                 sign_carrier, sz
00990 
00991 !   ---------------------------------------------------------------------------
00992 
00993     sign_carrier = 1
00994     IF (PRESENT (make_tr)) THEN
00995        IF (make_tr) sign_carrier = -1
00996     ENDIF
00997     mynblkrows = SIZE(myblkrows)
00998     mynblkcols = SIZE(myblkcols)
00999     meta(dbcsr_slot_nblkrows_local) = mynblkrows
01000     meta(dbcsr_slot_nblkcols_local) = mynblkcols
01001     nblks = mynblkrows * mynblkcols
01002     nze = 1
01003     IF (nblks .EQ. 0) THEN
01004        row_p(1:) = 0
01005     ELSE
01006        row_p(1) = 0
01007        !row_p(nrows+1) = nblks
01008        prev_row = 1
01009        blk = 0
01010        DO row_l = 1, mynblkrows
01011           row = myblkrows(row_l)
01012           row_p(prev_row+1:row) = blk
01013           DO col_l = 1, mynblkcols
01014              c = myblkcols(col_l)
01015              col_i(blk+col_l) = c
01016              sz = (row_blk_offsets(row+1)-row_blk_offsets(row))*&
01017                   (col_blk_offsets(c+1)-col_blk_offsets(c))
01018              IF (sz .GT. 0) THEN
01019                 blk_p(blk+col_l) = SIGN (nze, sign_carrier)
01020                 nze = nze + sz
01021              ELSE
01022                 blk_p(blk+col_l) = 0
01023              ENDIF
01024           ENDDO
01025           prev_row = row
01026           blk = blk + mynblkcols
01027        END DO
01028        CALL dbcsr_assert (blk, "EQ", nblks, dbcsr_fatal_level,&
01029             dbcsr_internal_error, routineN, "Block mismatch", __LINE__,&
01030             error=error)
01031        row_p(prev_row+1:nblkrows_total+1) = nblks
01032     ENDIF
01033     IF (debug_mod) THEN
01034        WRITE(*,*)routineN//" new index"
01035        WRITE(*,*)"row_p=",row_p
01036        WRITE(*,*)"col_i=",col_i
01037        WRITE(*,*)"blk_p=",blk_p
01038     ENDIF
01039     meta(dbcsr_slot_nblkrows_total) = nblkrows_total
01040     meta(dbcsr_slot_nblkcols_total) = nblkcols_total
01041   END SUBROUTINE make_dense_index
01042 
01043 ! *****************************************************************************
01052   SUBROUTINE make_undense_index (&
01053        row_p, col_i, blk_p,&
01054        distribution,  local_row_offsets, local_col_offsets,&
01055        meta)
01056     INTEGER, DIMENSION(:), INTENT(OUT)       :: row_p, col_i, blk_p
01057     TYPE(dbcsr_distribution_obj)             :: distribution
01058     INTEGER, DIMENSION(:), INTENT(IN)        :: local_row_offsets, 
01059                                                 local_col_offsets
01060     INTEGER, DIMENSION(dbcsr_meta_size), 
01061       INTENT(INOUT)                          :: meta
01062 
01063     CHARACTER(len=*), PARAMETER :: routineN = 'make_undense_index', 
01064       routineP = moduleN//':'//routineN
01065 
01066     INTEGER :: col, lr, lrow, nblkcols_local, nblkrows_local, nblkrows_total, 
01067       nfullcols_local, prev_row, row
01068     INTEGER, DIMENSION(:), POINTER           :: local_cols, local_rows
01069 
01070 !   ---------------------------------------------------------------------------
01071 
01072     local_cols => array_data (dbcsr_distribution_local_cols (distribution))
01073     local_rows => array_data (dbcsr_distribution_local_rows (distribution))
01074     meta(dbcsr_slot_nblkrows_total) = dbcsr_distribution_nrows (distribution)
01075     meta(dbcsr_slot_nblkcols_total) = dbcsr_distribution_ncols (distribution)
01076     meta(dbcsr_slot_nblkrows_local) = dbcsr_distribution_nlocal_rows (distribution)
01077     meta(dbcsr_slot_nblkcols_local) = dbcsr_distribution_nlocal_cols (distribution)
01078     nblkrows_total = meta(dbcsr_slot_nblkrows_total)
01079     nblkcols_local = meta(dbcsr_slot_nblkcols_local)
01080     nblkrows_local = meta(dbcsr_slot_nblkrows_local)
01081     nfullcols_local = meta(dbcsr_slot_nfullcols_local)
01082     ! Fill the row_p array.
01083     lr = 0
01084     row_p(1) = 0
01085     prev_row = 1
01086     DO lrow = 1, nblkrows_local
01087        row = local_rows(lrow)
01088        row_p(prev_row+1:row) = lr
01089        lr = lr + nblkcols_local
01090        row_p(row+1) = lr
01091        prev_row = row
01092     ENDDO
01093     row_p(prev_row+1:nblkrows_total+1) = lr
01094     !
01095     FORALL (row = 1 : nblkrows_local)
01096        FORALL (col = 1 : nblkcols_local)
01097           col_i(nblkcols_local*(row-1)+col) = local_cols(col)
01098           blk_p(nblkcols_local*(row-1)+col) = 1 + &
01099                (local_row_offsets(row)-1)*nfullcols_local&
01100                + (local_col_offsets(col)-1)*&
01101                (local_row_offsets(row+1)-local_row_offsets(row))
01102        END FORALL
01103     END FORALL
01104   END SUBROUTINE make_undense_index
01105 
01106 
01107 
01108 ! *****************************************************************************
01134   SUBROUTINE merge_index_arrays (new_row_i, new_col_i, new_blk_p, new_size,&
01135        old_row_i, old_col_i, old_blk_p, old_size,&
01136        add_ip, add_size, new_blk_d, old_blk_d,&
01137        added_size_offset, added_sizes, added_size, added_nblks, error)
01138     INTEGER, INTENT(IN)                      :: new_size
01139     INTEGER, DIMENSION(new_size), 
01140       INTENT(OUT)                            :: new_blk_p, new_col_i, 
01141                                                 new_row_i
01142     INTEGER, INTENT(IN)                      :: old_size
01143     INTEGER, DIMENSION(old_size), INTENT(IN) :: old_blk_p, old_col_i, 
01144                                                 old_row_i
01145     INTEGER, INTENT(IN)                      :: add_size
01146     INTEGER, DIMENSION(3, add_size), 
01147       INTENT(IN)                             :: add_ip
01148     INTEGER, DIMENSION(new_size), 
01149       INTENT(OUT), OPTIONAL                  :: new_blk_d
01150     INTEGER, DIMENSION(old_size), 
01151       INTENT(IN), OPTIONAL                   :: old_blk_d
01152     INTEGER, INTENT(IN), OPTIONAL            :: added_size_offset
01153     INTEGER, DIMENSION(:), INTENT(IN), 
01154       OPTIONAL                               :: added_sizes
01155     INTEGER, INTENT(OUT), OPTIONAL           :: added_size, added_nblks
01156     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
01157 
01158     CHARACTER(len=*), PARAMETER :: routineN = 'merge_index_arrays', 
01159       routineP = moduleN//':'//routineN
01160 
01161     INTEGER                                  :: add_blk, bp, i, 
01162                                                 merge_from_whom, new_blk, 
01163                                                 old_blk
01164     LOGICAL                                  :: multidata
01165 
01166 !   ---------------------------------------------------------------------------
01167 
01168     multidata = PRESENT (old_blk_d) .AND. PRESENT (new_blk_d)
01169     CALL dbcsr_assert(old_size+add_size .EQ. new_size, &
01170          dbcsr_warning_level, dbcsr_internal_error ,"merge_arrays",&
01171          "Mismatch of new and old size",__LINE__,error)
01172     CALL dbcsr_assert (PRESENT (added_size_offset), "EQV",&
01173          PRESENT (added_sizes), dbcsr_fatal_level, dbcsr_wrong_args_error,&
01174          routineN, "Must specify a set of arguments", __LINE__, error=error)
01175     CALL dbcsr_assert (PRESENT (added_sizes), "EQV", PRESENT (added_size),&
01176          dbcsr_fatal_level, dbcsr_wrong_args_error,&
01177          routineN, "Must specify a set of arguments", __LINE__, error=error)
01178     IF (debug_mod) THEN
01179        WRITE (*,*) " Old array", old_size
01180        DO i = 1, old_size
01181           WRITE(*,'(I7,2X,I7,2X,I7)')old_row_i(i),old_col_i(i),old_blk_p(i)
01182        ENDDO
01183        WRITE (*,*) " Add array", add_size
01184        DO i = 1, add_size
01185           WRITE(*,'(I7,2X,I7,2X,I7)')add_ip (1:3, i)
01186        ENDDO
01187     ENDIF
01188     IF (PRESENT (added_nblks)) added_nblks = 0
01189     IF (PRESENT (added_size)) THEN
01190        added_size = 0
01191        bp = added_size_offset
01192     ENDIF
01193     IF (add_size .GT. 0) THEN
01194        old_blk = 1
01195        add_blk = 1
01196        new_blk = 1
01197        IF (old_size .EQ. 0) THEN
01198           new_row_i(1:add_size) = add_ip(1, 1:add_size)
01199           new_col_i(1:add_size) = add_ip(2, 1:add_size)
01200           new_blk_p(1:add_size) = add_ip(3, 1:add_size)
01201           !IF (multidata) new_blk_d(1:add_size) = add_ip(4, 1:add_size)
01202           IF (PRESENT (added_nblks)) added_nblks = add_size
01203           IF (PRESENT (added_size)) added_size = SUM (added_sizes)
01204        ELSE
01205           DO WHILE (new_blk .LE. new_size)
01206              merge_from_whom = 0
01207              IF (old_blk .LE. old_size .AND. add_blk .LE. add_size) THEN
01208                 IF (add_ip(1, add_blk) .EQ. old_row_i(old_blk)&
01209                      .AND.add_ip(2, add_blk) .EQ. old_col_i(old_blk)) THEN
01210                    IF (debug_mod) THEN
01211                       WRITE(*,*)"Duplicate block! addblk",&
01212                            add_blk, "oldblk", old_blk
01213                    ENDIF
01214                 ENDIF
01215                 ! Rows come first
01216                 IF (add_ip(1, add_blk) .LT. old_row_i(old_blk)) THEN
01217                    merge_from_whom = 2
01218                 ELSEIF (add_ip(1, add_blk) .GT. old_row_i(old_blk)) THEN
01219                    merge_from_whom = 1
01220                 ELSE ! Same rows, so now come the columns
01221                    IF (add_ip(2, add_blk) .LT. old_col_i(old_blk)) THEN
01222                       ! Merges from the add array
01223                       merge_from_whom = 2
01224                    ELSEIF (add_ip(2, add_blk) .GT. old_col_i(old_blk)) THEN
01225                       ! Merges from the old array
01226                       merge_from_whom = 1
01227                    ELSE
01228                       ! Merge from old array and skip one in the new array
01229                       IF (debug_mod) THEN
01230                          WRITE(*,*)"Duplicate, keeping old",&
01231                          add_ip(1, add_blk), add_ip(2, add_blk)
01232                       ENDIF
01233                       merge_from_whom = 1
01234                       add_blk = add_blk + 1
01235                    ENDIF
01236                 ENDIF
01237              ELSE
01238                 IF (add_blk .LE. add_size) THEN
01239                    ! Merges from the add array
01240                    merge_from_whom = 2
01241                 ELSEIF (old_blk .LE. old_size) THEN
01242                    ! Merges from the old array
01243                    merge_from_whom = 1
01244                 ELSE
01245                    ! Hmmm, nothing to merge...
01246                    merge_from_whom = 0
01247                    !WRITE(*,*)"Ran out of data to merge"
01248                 ENDIF
01249              ENDIF
01250              SELECT CASE (merge_from_whom)
01251              CASE (2)
01252                 ! Merges from the add array
01253                 new_row_i(new_blk) = add_ip(1, add_blk)
01254                 new_col_i(new_blk) = add_ip(2, add_blk)
01255                 new_blk_p(new_blk) = add_ip(3, add_blk)
01256                 !IF (multidata) new_blk_d(new_blk) = add_ip(4, add_blk)
01257                 IF (PRESENT (added_nblks)) added_nblks = added_nblks + 1
01258                 IF (PRESENT (added_sizes)) THEN
01259                    new_blk_p(new_blk) = bp
01260                    bp = bp + added_sizes(add_blk)
01261                    added_size = added_size + added_sizes(add_blk)
01262                 ENDIF
01263                 add_blk = add_blk + 1
01264              CASE (1)
01265                 ! Merges from the old array
01266                 new_row_i(new_blk) = old_row_i(old_blk)
01267                 new_col_i(new_blk) = old_col_i(old_blk)
01268                 new_blk_p(new_blk) = old_blk_p(old_blk)
01269                 IF (multidata) new_blk_p(new_blk) = old_blk_d(old_blk)
01270                 old_blk = old_blk + 1
01271              CASE DEFAULT
01272                 !WRITE(*,*)"Nothing to merge"
01273              END SELECT
01274              new_blk = new_blk + 1
01275           ENDDO
01276        ENDIF
01277     ELSE
01278        new_row_i(1:old_size) = old_row_i(1:old_size)
01279        new_col_i(1:old_size) = old_col_i(1:old_size)
01280        new_blk_p(1:old_size) = old_blk_p(1:old_size)
01281        IF (multidata) new_blk_d(1:old_size) = old_blk_d(1:old_size)
01282     ENDIF
01283     IF (debug_mod) THEN
01284        WRITE (*,*) " New array"
01285        DO i = 1, new_size
01286           WRITE(*,'(4(2X,I7))')new_row_i(i),new_col_i(i),new_blk_p(i)
01287        ENDDO
01288     ENDIF
01289   END SUBROUTINE merge_index_arrays
01290 
01291 
01292 ! *****************************************************************************
01298   SUBROUTINE dbcsr_make_index_local_row (matrix, error)
01299     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
01300     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
01301 
01302     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_make_index_local_row', 
01303       routineP = moduleN//':'//routineN
01304 
01305     INTEGER                                  :: error_handle, lrow, 
01306                                                 nlocal_rows, ntotal_rows, 
01307                                                 prow, stat
01308     INTEGER, ALLOCATABLE, DIMENSION(:)       :: local_row_p
01309     INTEGER, DIMENSION(:), POINTER           :: local_rows
01310 
01311 !   ---------------------------------------------------------------------------
01312 
01313     CALL dbcsr_error_set (routineN, error_handle, error=error)
01314     CALL dbcsr_assert (ASSOCIATED (matrix%m%row_p),&
01315          dbcsr_fatal_level, dbcsr_wrong_args_error, routineN,&
01316          "The row index must be initialized.", __LINE__, error=error)
01317     CALL dbcsr_assert ("NOT", matrix%m%bcsc,&
01318          dbcsr_fatal_level, dbcsr_unimplemented_error_nr, routineN,&
01319          "Not support for BCSC yet.", __LINE__, error=error)
01320     !
01321     prow = matrix%m%index(dbcsr_slot_home_vprow)
01322     IF (prow .LT. 0) THEN
01323        prow = matrix%m%index(dbcsr_slot_home_prow)
01324     ENDIF
01325     nlocal_rows = matrix%m%nblkrows_local
01326     ALLOCATE (local_row_p(nlocal_rows+1), stat=stat)
01327     CALL dbcsr_assert (stat, "EQ", 0,&
01328          dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01329          "Could not allocate memory.", __LINE__, error=error)
01333     ntotal_rows = matrix%m%nblkrows_total
01334     CALL dbcsr_count_row_index (matrix%m%row_p, ntotal_rows)
01335     ! We first have to find the local rows for the given prow.
01336     local_rows => array_data (matrix%m%local_rows)
01337     CALL dbcsr_assert (SIZE(local_rows), "EQ", nlocal_rows,&
01338          dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01339          "Mismatch in the number of local rows.", __LINE__, error=error)
01341     FORALL (lrow = 1:nlocal_rows)
01342        local_row_p(lrow) = matrix%m%row_p(local_rows(lrow))
01343     END FORALL
01344     CALL dbcsr_assert (SUM (matrix%m%row_p(1:ntotal_rows)), "EQ",&
01345          SUM(local_row_p(1:nlocal_rows)),&
01346          dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01347          "Inconsistent row counts. Perhaps non-local rows contain data?.",&
01348          __LINE__, error=error)
01350     CALL dbcsr_build_row_index (local_row_p, nlocal_rows)
01352     CALL dbcsr_clearfrom_index_array(matrix%m, dbcsr_slot_row_p)
01353     CALL dbcsr_addto_index_array(matrix%m, dbcsr_slot_row_p, DATA=local_row_p,&
01354          error=error)
01356     matrix%m%local_indexing = .TRUE.
01357     IF (careful_mod) THEN
01358        CALL dbcsr_assert (array_size(matrix%m%local_rows), "EQ", nlocal_rows,&
01359             dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01360             "Inconsistent local row counts.",&
01361             __LINE__, error=error)
01362        CALL dbcsr_assert (array_size(matrix%m%global_rows), "EQ", ntotal_rows,&
01363             dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01364             "Inconsistent global row counts.",&
01365             __LINE__, error=error)
01366        IF (array_size (matrix%m%global_rows) .EQ. 0) THEN
01367           CALL dbcsr_assert (nlocal_rows, "EQ", 0,&
01368                dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01369                "Invalid number of local or global rows.",&
01370                __LINE__, error=error)
01371        ENDIF
01372     ENDIF
01373     DEALLOCATE (local_row_p)
01374     !
01375     CALL dbcsr_error_stop (error_handle, error)
01376   END SUBROUTINE dbcsr_make_index_local_row
01377 
01378 
01379 ! *****************************************************************************
01385   SUBROUTINE dbcsr_make_index_list (matrix, thread_redist, error)
01386     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
01387     LOGICAL, INTENT(IN)                      :: thread_redist
01388     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
01389 
01390     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_make_index_list', 
01391       routineP = moduleN//':'//routineN
01392     LOGICAL, PARAMETER                       :: dbg = .FALSE.
01393 
01394     INTEGER                                  :: blk, error_handle, ithread, 
01395                                                 my_cnt, nblks, nrows, nthreads
01396     INTEGER, ALLOCATABLE, DIMENSION(:, :)    :: blki
01397     INTEGER, DIMENSION(0)                    :: zero_len_array
01398     INTEGER, DIMENSION(:), POINTER           :: global_cols, local_rows, td, 
01399                                                 thr_c
01400 
01401 !   ---------------------------------------------------------------------------
01402 
01403     CALL dbcsr_error_set (routineN, error_handle, error=error)
01404     CALL dbcsr_assert (ASSOCIATED (matrix%m%row_p),&
01405          dbcsr_fatal_level, dbcsr_wrong_args_error, routineN,&
01406          "The row index must be initialized.", __LINE__, error=error)
01407     CALL dbcsr_assert ("NOT", matrix%m%list_indexing,&
01408          dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01409          "List index already exists?", __LINE__, error=error)
01410     CALL dbcsr_assert ("NOT", matrix%m%bcsc,&
01411          dbcsr_fatal_level, dbcsr_unimplemented_error_nr, routineN,&
01412          "Not support for BCSC yet.", __LINE__, error=error)
01413     CALL dbcsr_assert (matrix%m%nblks, "EQ", SIZE(matrix%m%col_i),&
01414          dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01415          "Block count mismatch.", __LINE__, error=error)
01416     CALL dbcsr_assert (matrix%m%nblks, "EQ", SIZE(matrix%m%blk_p),&
01417          dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01418          "Block count mismatch", __LINE__, error=error)
01419     !
01420     IF (matrix%m%local_indexing) THEN
01421        CALL dbcsr_assert (SIZE(matrix%m%row_p)-1, "EQ", matrix%m%nblkrows_local,&
01422             dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01423             "Local row index incorrectly sized.", __LINE__, error=error)
01424     ELSE
01425        CALL dbcsr_assert (SIZE(matrix%m%row_p)-1, "EQ", matrix%m%nblkrows_total,&
01426             dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01427             "Global row index incorrectly sized", __LINE__, error=error)
01428     ENDIF
01429     !
01430     matrix%m%list_indexing = .TRUE.
01431     !
01432     IF (matrix%m%local_indexing) THEN
01433        global_cols => array_data (matrix%m%global_cols)
01434        nrows = matrix%m%nblkrows_local
01435     ELSE
01436        nrows = matrix%m%nblkrows_total
01437     ENDIF
01438     !
01439     nblks = matrix%m%nblks
01440     ALLOCATE (blki(3, nblks))
01441     CALL dbcsr_expand_row_index (matrix%m%row_p, blki(1,:), nrows, nblks)
01442     IF (matrix%m%local_indexing) THEN
01443        ! If local indexing is enabled, then the rows but not the
01444        ! columns are already localized
01445        IF (dbg) THEN
01446           WRITE(*,*)routineN//" Making local columns"
01447           WRITE(*,'(10(1X,i7))')global_cols
01448           WRITE(*,*)'local'
01449           WRITE(*,'(10(1X,i7))')array_data(matrix%m%local_cols)
01450        ENDIF
01451        FORALL (blk = 1 : nblks)
01452           blki(2,blk) = global_cols(matrix%m%col_i(blk))
01453           blki(3,blk) = matrix%m%blk_p(blk)
01454        END FORALL
01455     ELSE
01456        IF (dbg) WRITE(*,*)routineN//" Not making local columns"
01457        FORALL (blk = 1 : nblks)
01458           blki(2,blk) = matrix%m%col_i(blk)
01459           blki(3,blk) = matrix%m%blk_p(blk)
01460        END FORALL
01461     ENDIF
01462     CALL dbcsr_clearfrom_index_array (matrix%m, dbcsr_slot_row_p)
01463     CALL dbcsr_clearfrom_index_array (matrix%m, dbcsr_slot_col_i)
01464     CALL dbcsr_clearfrom_index_array (matrix%m, dbcsr_slot_blk_p)
01465     nthreads = 0
01466     !$ nthreads = OMP_GET_MAX_THREADS ()
01467     !
01468     ! Reshuffle according to threads
01469     IF (nthreads .GT. 0 .AND. thread_redist) THEN
01470        td => array_data (dbcsr_distribution_thread_dist (dbcsr_distribution (matrix)))
01471        IF (matrix%m%local_indexing) THEN
01472           local_rows => array_data (matrix%m%local_rows)
01473        ENDIF
01474        !$OMP PARALLEL DEFAULT (none) &
01475        !$OMP PRIVATE (my_cnt, ithread, blk) &
01476        !$OMP SHARED (td, blki, nthreads, thr_c, nblks, matrix, error, local_rows)
01477        !
01478        ithread = 0
01479        !$ ithread = OMP_GET_THREAD_NUM ()
01480        !$OMP MASTER
01481        !$ nthreads = OMP_GET_NUM_THREADS ()
01482        CALL dbcsr_addto_index_array (matrix%m, dbcsr_slot_thr_c,&
01483             reservation=nthreads+1, extra=3*nblks, error=error)
01484        thr_c => matrix%m%thr_c
01485        CALL dbcsr_addto_index_array (matrix%m, dbcsr_slot_coo_l,&
01486             reservation=3*nblks, error=error)
01487        !$OMP END MASTER
01488        !$OMP BARRIER
01489        my_cnt = 0
01490        IF (matrix%m%local_indexing) THEN
01491           my_cnt = COUNT (td(local_rows(blki(1,:))) .EQ. ithread)
01492        ELSE
01493           my_cnt = COUNT (td(blki(1,:)) .EQ. ithread)
01494        ENDIF
01495        !DO blk = 1, nblks
01496        !   IF (td(blki(1, blk)) .EQ. ithread) my_cnt = my_cnt+1
01497        !ENDDO
01498        thr_c(ithread+1) = my_cnt
01499        !$OMP BARRIER
01500        !$OMP MASTER
01501        CALL dbcsr_build_row_index_inplace (thr_c, nthreads)
01502        !$OMP END MASTER
01503        !$OMP BARRIER
01504        my_cnt = (thr_c(ithread+1)+1) * 3 -2
01505        IF (matrix%m%local_indexing) THEN
01506           DO blk = 1, nblks
01507              IF (td(local_rows(blki(1, blk))) .EQ. ithread) THEN
01508                 matrix%m%coo_l(my_cnt:my_cnt+2) = blki(1:3, blk)
01509                 my_cnt = my_cnt + 3
01510              ENDIF
01511           ENDDO
01512        ELSE
01513           DO blk = 1, nblks
01514              IF (td(blki(1, blk)) .EQ. ithread) THEN
01515                 matrix%m%coo_l(my_cnt:my_cnt+2) = blki(1:3, blk)
01516                 my_cnt = my_cnt + 3
01517              ENDIF
01518           ENDDO
01519        ENDIF
01520        !$OMP END PARALLEL
01521     ELSE
01522        ! Small price to pay for avoiding infinite recursions.
01523        DO blk = 2, nblks
01524           IF (       blki(1, blk) .EQ. blki(1, blk-1)&
01525                .AND. blki(2, blk) .EQ. blki(2, blk-1)) &
01526                ! Weird assertion to give some idea of the two blocks.
01527                CALL dbcsr_assert (-blki(1,blk), "EQ", blki(2,blk),&
01528                dbcsr_fatal_level, dbcsr_caller_error, routineN,&
01529       "Should not have duplicate matrix blocks. (-row, col) is duplicated.",&
01530                __LINE__, error=error)
01531        ENDDO
01532        !
01533        IF (nblks>0) THEN
01534           CALL dbcsr_addto_index_array (matrix%m, dbcsr_slot_coo_l,&
01535                DATA=RESHAPE(blki, (/3*nblks/)), error=error)
01536        ELSE
01537           CALL dbcsr_addto_index_array (matrix%m, dbcsr_slot_coo_l,&
01538                DATA=zero_len_array, error=error)
01539        ENDIF
01540     ENDIF
01541 
01542     DEALLOCATE (blki)
01543     !
01544     CALL dbcsr_error_stop (error_handle, error)
01545   END SUBROUTINE dbcsr_make_index_list
01546 
01547 
01548 
01549 ! *****************************************************************************
01554   SUBROUTINE dbcsr_index_compact(matrix, error)
01555     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
01556     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
01557 
01558     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_index_compact', 
01559       routineP = moduleN//':'//routineN
01560 
01561     INTEGER                                  :: error_handle, new_size, 
01562                                                 size_blk_p, size_col_i, 
01563                                                 size_coo_l, size_row_p, 
01564                                                 size_thr_c
01565     INTEGER, ALLOCATABLE, DIMENSION(:)       :: blk_p, col_i, coo_l, meta, 
01566                                                 row_p, thr_c
01567     LOGICAL                                  :: compact, has_blk_p, 
01568                                                 has_col_i, has_coo_l, 
01569                                                 has_row_p, has_thr_c
01570 
01571 !   ---------------------------------------------------------------------------
01572 
01573     CALL dbcsr_error_set (routineN, error_handle, error=error)
01575     CALL dbcsr_repoint_index (matrix%m)
01576     ! Check that compaction is even needed.
01577     has_row_p = ASSOCIATED (matrix%m%row_p)
01578     IF (has_row_p) THEN
01579        size_row_p = SIZE(matrix%m%row_p)
01580     ELSE
01581        size_row_p = 0
01582     ENDIF
01583     has_col_i = ASSOCIATED (matrix%m%col_i)
01584     IF (has_col_i) THEN
01585        size_col_i = SIZE(matrix%m%col_i)
01586     ELSE
01587        size_col_i = 0
01588     ENDIF
01589     has_blk_p = ASSOCIATED (matrix%m%blk_p)
01590     IF (has_blk_p) THEN
01591        size_blk_p = SIZE(matrix%m%blk_p)
01592     ELSE
01593        size_blk_p = 0
01594     ENDIF
01595     has_thr_c = ASSOCIATED (matrix%m%thr_c)
01596     IF (has_thr_c) THEN
01597        size_thr_c = SIZE(matrix%m%thr_c)
01598     ELSE
01599        size_thr_c = 0
01600     ENDIF
01601     has_coo_l = ASSOCIATED (matrix%m%coo_l)
01602     IF (has_coo_l) THEN
01603        size_coo_l = SIZE(matrix%m%coo_l)
01604     ELSE
01605        size_coo_l = 0
01606     ENDIF
01607     !
01608     new_size = dbcsr_num_slots +&
01609          size_row_p + size_col_i + size_blk_p + size_thr_c + size_coo_l
01610     compact = new_size .LT. SIZE(matrix%m%index)
01611     IF (compact) THEN
01612        ! Store old index arrays.
01613        IF (has_blk_p) THEN
01614           ALLOCATE (row_p(size_row_p))
01615           row_p(:) = matrix%m%row_p(:)
01616        ENDIF
01617        IF (has_col_i) THEN
01618           ALLOCATE (col_i(size_col_i))
01619           col_i(:) = matrix%m%col_i(:)
01620        ENDIF
01621        IF (has_blk_p) THEN
01622           ALLOCATE (blk_p(size_blk_p))
01623           blk_p(:) = matrix%m%blk_p(:)
01624        ENDIF
01625        IF (has_thr_c) THEN
01626           ALLOCATE (thr_c(size_thr_c))
01627           thr_c(:) = matrix%m%thr_c(:)
01628        ENDIF
01629        IF (has_coo_l) THEN
01630           ALLOCATE (coo_l(size_coo_l))
01631           coo_l(:) = matrix%m%coo_l(:)
01632        ENDIF
01633        ALLOCATE (meta(dbcsr_num_slots))
01634        meta(:) = matrix%m%index(1:dbcsr_num_slots)
01635        ! Clear the index.
01636        CALL memory_deallocate (matrix%m%index,&
01637             dbcsr_get_index_memory_type (matrix), error=error)
01638        NULLIFY (matrix%m%index)
01639        CALL memory_allocate (matrix%m%index, new_size,&
01640             dbcsr_get_index_memory_type (matrix), error=error)
01641        !
01642        ! Now copy the old index arrays into the index. We must not
01643        ! copy the positions of the old pointers.
01644        matrix%m%index(1:dbcsr_meta_size) = meta(1:dbcsr_meta_size)
01645        matrix%m%index(dbcsr_meta_size+1:) = 0
01646        matrix%m%index(dbcsr_slot_size) = dbcsr_num_slots
01647        IF (has_thr_c) THEN
01648           CALL dbcsr_addto_index_array(matrix%m, dbcsr_slot_thr_c, thr_c,&
01649                error=error)
01650           DEALLOCATE (thr_c)
01651        ENDIF
01652        IF (has_row_p) THEN
01653           CALL dbcsr_addto_index_array(matrix%m, dbcsr_slot_row_p, row_p,&
01654                error=error)
01655           DEALLOCATE (row_p)
01656        ENDIF
01657        IF (has_col_i) THEN
01658           CALL dbcsr_addto_index_array(matrix%m, dbcsr_slot_col_i, col_i,&
01659                error=error)
01660           DEALLOCATE (col_i)
01661        ENDIF
01662        IF (has_blk_p) THEN
01663           CALL dbcsr_addto_index_array(matrix%m, dbcsr_slot_blk_p, blk_p,&
01664                error=error)
01665           DEALLOCATE (blk_p)
01666        ENDIF
01667        IF (has_coo_l) THEN
01668           CALL dbcsr_addto_index_array(matrix%m, dbcsr_slot_coo_l, coo_l,&
01669                error=error)
01670           DEALLOCATE (coo_l)
01671        ENDIF
01672        DEALLOCATE (meta)
01673        IF (careful_mod) THEN
01674           ! This is pretty strong but it should be true.
01675           CALL dbcsr_assert (matrix%m%index(dbcsr_slot_size), "EQ", new_size,&
01676                dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01677                "Unexpected index size.", __LINE__, error=error)
01678           CALL dbcsr_assert (SIZE(matrix%m%index), "EQ", new_size,&
01679                dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01680                "Unexpected index size.", __LINE__, error=error)
01681        ENDIF
01682        CALL dbcsr_repoint_index (matrix%m)
01683     ENDIF
01684     CALL dbcsr_error_stop (error_handle, error)
01685   END SUBROUTINE dbcsr_index_compact
01686 
01687 ! *****************************************************************************
01692   SUBROUTINE dbcsr_index_checksum(matrix, checksum, error)
01693     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
01694     INTEGER, INTENT(OUT)                     :: checksum
01695     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
01696 
01697     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_index_checksum', 
01698       routineP = moduleN//':'//routineN
01699 
01700     INTEGER                                  :: error_handle
01701 
01702 !   ---------------------------------------------------------------------------
01703 
01704     CALL dbcsr_error_set (routineN, error_handle, error=error)
01705     !
01706     checksum = joaat_hash ( (/ joaat_hash (matrix%m%row_p),&
01707                                joaat_hash (matrix%m%col_i),&
01708                                joaat_hash (matrix%m%blk_p) /) )
01709     !
01710     CALL dbcsr_error_stop (error_handle, error)
01711   END SUBROUTINE dbcsr_index_checksum
01712 
01713   FUNCTION dbcsr_has_local_row_index (matrix) RESULT (local_indexing)
01714     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
01715     LOGICAL                                  :: local_indexing
01716 
01717     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_has_local_row_index', 
01718       routineP = moduleN//':'//routineN
01719 
01720     local_indexing = matrix%m%local_indexing
01721   END FUNCTION dbcsr_has_local_row_index
01722 
01723 END MODULE dbcsr_index_operations