CP2K 2.4 (Revision 12889)

dbcsr_methods.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_methods
00016   USE array_types,                     ONLY: array_data,&
00017                                              array_exists,&
00018                                              array_hold,&
00019                                              array_i1d_obj,&
00020                                              array_new,&
00021                                              array_nullify,&
00022                                              array_release,&
00023                                              array_size
00024   USE btree_I8_k_cp2d_v,               ONLY: btree_destroy_c => btree_delete,&
00025                                              btree_new_c => btree_new
00026   USE btree_I8_k_dp2d_v,               ONLY: btree_destroy_d => btree_delete,&
00027                                              btree_new_d => btree_new
00028   USE btree_I8_k_sp2d_v,               ONLY: btree_destroy_s => btree_delete,&
00029                                              btree_new_s => btree_new
00030   USE btree_I8_k_zp2d_v,               ONLY: btree_destroy_z => btree_delete,&
00031                                              btree_new_z => btree_new
00032   USE dbcsr_block_buffers,             ONLY: dbcsr_buffers_2d_needed,&
00033                                              dbcsr_buffers_flush,&
00034                                              dbcsr_buffers_init,&
00035                                              dbcsr_buffers_new,&
00036                                              dbcsr_buffers_release
00037   USE dbcsr_config,                    ONLY: comm_thread_load,&
00038                                              use_subcommunicators
00039   USE dbcsr_data_methods,              ONLY: dbcsr_data_get_size,&
00040                                              dbcsr_data_get_size_referenced,&
00041                                              dbcsr_data_get_type,&
00042                                              dbcsr_data_hold,&
00043                                              dbcsr_data_init,&
00044                                              dbcsr_data_release,&
00045                                              dbcsr_data_set_size_referenced
00046   USE dbcsr_error_handling,            ONLY: &
00047        dbcsr_assert, dbcsr_caller_error, dbcsr_error_set, dbcsr_error_stop, &
00048        dbcsr_error_type, dbcsr_failure_level, dbcsr_fatal_level, &
00049        dbcsr_internal_error, dbcsr_warning_level, dbcsr_wrong_args_error
00050   USE dbcsr_kinds,                     ONLY: default_string_length,&
00051                                              real_8,&
00052                                              sp
00053   USE dbcsr_message_passing,           ONLY: mp_cart_create,&
00054                                              mp_cart_sub,&
00055                                              mp_comm_free,&
00056                                              mp_sum
00057   USE dbcsr_ptr_util,                  ONLY: memory_deallocate
00058   USE dbcsr_toollib,                   ONLY: sort
00059   USE dbcsr_types,                     ONLY: &
00060        dbcsr_1d_array_obj, dbcsr_2d_array_obj, dbcsr_2d_array_type, &
00061        dbcsr_array_type, dbcsr_block_buffer_obj, dbcsr_data_obj, &
00062        dbcsr_distribution_obj, dbcsr_imagedistribution_obj, &
00063        dbcsr_imagedistribution_type, dbcsr_magic_number, &
00064        dbcsr_memory_default, dbcsr_mp_obj, dbcsr_mutable_obj, dbcsr_obj, &
00065        dbcsr_type, dbcsr_type_antihermitian, dbcsr_type_antisymmetric, &
00066        dbcsr_type_complex_4, dbcsr_type_complex_8, dbcsr_type_hermitian, &
00067        dbcsr_type_invalid, dbcsr_type_no_symmetry, dbcsr_type_real_4, &
00068        dbcsr_type_real_8, dbcsr_type_symmetric, dbcsr_work_type
00069   USE min_heap,                        ONLY: heap_fill,&
00070                                              heap_get_first,&
00071                                              heap_new,&
00072                                              heap_release,&
00073                                              heap_reset_first,&
00074                                              heap_t
00075 
00076   !$ USE OMP_LIB
00077   IMPLICIT NONE
00078 
00079   PRIVATE
00080 
00081   CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'dbcsr_methods'
00082 
00083   INTEGER, PUBLIC, SAVE :: dbcsr_matrix_counter = 111111
00084 
00085   PUBLIC :: dbcsr_init, dbcsr_init_p
00086   PUBLIC :: dbcsr_hold, dbcsr_release, dbcsr_release_p, dbcsr_destroy
00087   PUBLIC :: dbcsr_valid_index, dbcsr_is_initialized
00088   PUBLIC :: dbcsr_mp_new, dbcsr_mp_hold, dbcsr_mp_release
00089   PUBLIC :: dbcsr_mp_pgrid, dbcsr_mp_numnodes, dbcsr_mp_mynode, dbcsr_mp_group,&
00090             dbcsr_mp_new_transposed, dbcsr_mp_nprows, dbcsr_mp_npcols,&
00091             dbcsr_mp_myprow, dbcsr_mp_mypcol, dbcsr_mp_pgrid_equal,&
00092             dbcsr_mp_get_coordinates,&
00093             dbcsr_mp_my_row_group, dbcsr_mp_my_col_group,&
00094             dbcsr_mp_has_subgroups, dbcsr_mp_get_process,&
00095             dbcsr_mp_grid_setup, dbcsr_mp_grid_remove,&
00096             dbcsr_mp_init, dbcsr_mp_exists, dbcsr_mp_active
00097   PUBLIC :: dbcsr_distribution_new, dbcsr_distribution_hold,&
00098             dbcsr_distribution_release, dbcsr_distribution_init
00099   PUBLIC :: dbcsr_distribution_mp, dbcsr_distribution_processor,&
00100             dbcsr_distribution_nrows, dbcsr_distribution_ncols,&
00101             dbcsr_distribution_row_dist, dbcsr_distribution_col_dist,&
00102             dbcsr_distribution_nlocal_rows, dbcsr_distribution_nlocal_cols,&
00103             dbcsr_distribution_local_rows, dbcsr_distribution_local_cols
00104   PUBLIC :: dbcsr_distribution_thread_dist, dbcsr_distribution_has_threads,&
00105             dbcsr_distribution_make_threads, dbcsr_distribution_no_threads,&
00106             dbcsr_distribution_num_threads
00107   PUBLIC :: dbcsr_distribution_add_row_map, dbcsr_distribution_add_col_map,&
00108             dbcsr_distribution_del_row_map, dbcsr_distribution_del_col_map
00109   PUBLIC :: dbcsr_release_locals, dbcsr_dist_release_locals
00110   PUBLIC :: dbcsr_get_info, dbcsr_distribution,&
00111             dbcsr_get_matrix_type, dbcsr_get_data_type, dbcsr_get_replication_type,&
00112             dbcsr_row_block_sizes, dbcsr_col_block_sizes,&
00113             dbcsr_nblkrows_total, dbcsr_nblkcols_total, dbcsr_nfullrows_total,&
00114             dbcsr_nfullcols_total, dbcsr_nblkcols_local, dbcsr_nblkrows_local,&
00115             dbcsr_nfullrows_local, dbcsr_nfullcols_local,&
00116             dbcsr_max_row_size, dbcsr_max_col_size,&
00117             dbcsr_get_occupation,&
00118             dbcsr_get_index_memory_type, dbcsr_get_data_memory_type, &
00119             dbcsr_is_real, dbcsr_is_complex,&
00120             dbcsr_name, dbcsr_may_be_dense,&
00121             dbcsr_use_mutable, dbcsr_wm_use_mutable, dbcsr_has_symmetry,&
00122             dbcsr_get_data_size
00123   PUBLIC :: dbcsr_get_data_size_used, dbcsr_get_data_size_referenced
00124   PUBLIC :: dbcsr_col_block_offsets, dbcsr_row_block_offsets
00125   PUBLIC :: dbcsr_data_area, dbcsr_switch_data_area
00126   PUBLIC :: dbcsr_get_num_blocks, dbcsr_get_nze
00127 
00128 
00129   PUBLIC :: dbcsr_blk_row_size, dbcsr_blk_column_size,&
00130             dbcsr_blk_row_offset, dbcsr_blk_col_offset
00131 
00132   PUBLIC :: dbcsr_allocate_matrix_array, dbcsr_array_new, dbcsr_array_put,&
00133             dbcsr_array_get, dbcsr_array_destroy, dbcsr_array_hold,&
00134             dbcsr_array_release
00135   PUBLIC :: dbcsr_destroy_array
00136   PUBLIC :: dbcsr_image_dist_init, dbcsr_image_dist_hold, dbcsr_image_dist_release
00137 
00138   PUBLIC :: dbcsr_mutable_init, dbcsr_mutable_new, dbcsr_mutable_destroy,&
00139             dbcsr_mutable_release, dbcsr_mutable_hold,&
00140             dbcsr_mutable_instantiated
00141 
00142   PUBLIC :: dbcsr_modify_lock, dbcsr_modify_unlock
00143 
00144 
00145 
00146   INTERFACE dbcsr_init
00147      MODULE PROCEDURE dbcsr_init_type, dbcsr_init_obj
00148   END INTERFACE
00149 
00150   INTERFACE dbcsr_init_p
00151      MODULE PROCEDURE dbcsr_init_obj_p
00152   END INTERFACE
00153 
00154   INTERFACE dbcsr_valid_index
00155      MODULE PROCEDURE dbcsr_valid_index_type, dbcsr_valid_index_obj
00156   END INTERFACE
00157 
00158   INTERFACE dbcsr_is_initialized
00159      MODULE PROCEDURE dbcsr_is_initialized_type, dbcsr_is_initialized_obj
00160   END INTERFACE
00161 
00162   !INTERFACE dbcsr_uses_special_memory
00163   !   MODULE PROCEDURE uses_special_memory_matrix, uses_special_memory_area
00164   !END INTERFACE
00165 
00166   INTERFACE dbcsr_get_data_size
00167      !MODULE PROCEDURE get_data_size_area, get_data_size_matrix
00168      MODULE PROCEDURE get_data_size_matrix
00169   END INTERFACE
00170 
00171   INTERFACE dbcsr_get_data_size_used
00172      !MODULE PROCEDURE get_data_size_ref_area
00173      MODULE PROCEDURE dbcsr_get_data_size_used
00174   END INTERFACE
00175 
00176   INTERFACE dbcsr_get_data_size_referenced
00177      !MODULE PROCEDURE get_data_size_ref_area, get_data_size_ref_matrix
00178      MODULE PROCEDURE get_data_size_ref_matrix
00179   END INTERFACE
00180 
00181   INTERFACE dbcsr_get_data_type
00182      !MODULE PROCEDURE dbcsr_get_data_type_obj, data_get_data_type
00183      MODULE PROCEDURE dbcsr_get_data_type_obj
00184   END INTERFACE
00185 
00186   INTERFACE dbcsr_is_real
00187      MODULE PROCEDURE dbcsr_is_real_obj
00188   END INTERFACE
00189 
00190   INTERFACE dbcsr_is_complex
00191      MODULE PROCEDURE dbcsr_is_complex_obj
00192   END INTERFACE
00193 
00194   INTERFACE dbcsr_blk_row_size
00195      MODULE PROCEDURE dbcsr_blk_row_size_type, dbcsr_blk_row_size_obj
00196   END INTERFACE
00197 
00198   INTERFACE dbcsr_blk_column_size
00199      MODULE PROCEDURE dbcsr_blk_column_size_type, dbcsr_blk_column_size_obj
00200   END INTERFACE
00201 
00202   ! For the 1-D and 2-D arrays
00203 
00204   INTERFACE dbcsr_allocate_matrix_array
00205      MODULE PROCEDURE array_init_1d, array_init_2d
00206   END INTERFACE
00207 
00208   INTERFACE dbcsr_array_new
00209      MODULE PROCEDURE array_new_1d, array_new_2d
00210   END INTERFACE
00211 
00212   INTERFACE dbcsr_array_put
00213      MODULE PROCEDURE array_put_1d, array_put_2d
00214   END INTERFACE
00215 
00216   INTERFACE dbcsr_array_get
00217      MODULE PROCEDURE array_get_1d, array_get_2d
00218   END INTERFACE
00219 
00220   INTERFACE dbcsr_array_destroy
00221      MODULE PROCEDURE array_destroy_1d, array_destroy_2d
00222   END INTERFACE
00223 
00224   INTERFACE dbcsr_array_hold
00225      MODULE PROCEDURE array_hold_1d, array_hold_2d
00226   END INTERFACE
00227 
00228   INTERFACE dbcsr_array_release
00229      MODULE PROCEDURE array_release_1d, array_release_2d
00230   END INTERFACE
00231 
00232   INTERFACE dbcsr_destroy_array
00233      MODULE PROCEDURE dbcsr_destroy_2d_array, dbcsr_destroy_1d_array
00234   END INTERFACE
00235 
00236 
00237 
00238 CONTAINS
00239 
00240 
00241 ! *****************************************************************************
00246   SUBROUTINE dbcsr_init_type (matrix)
00247     TYPE(dbcsr_type), INTENT(OUT)            :: matrix
00248 
00249     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_init_type', 
00250       routineP = moduleN//':'//routineN
00251 
00252 !TYPE(dbcsr_error_type), INTENT(INOUT)   :: error
00253 !   ---------------------------------------------------------------------------
00254 
00255     matrix%initialized = dbcsr_magic_number
00256     matrix%valid = .FALSE.
00257     matrix%refcount = 0
00258     ! Nullifies all pointers.
00259     NULLIFY (matrix%index, matrix%row_p, matrix%col_i,&
00260          matrix%blk_p, matrix%thr_c, matrix%coo_l)
00261     CALL dbcsr_data_init (matrix%data_area)
00262     CALL dbcsr_distribution_init (matrix%dist)
00263     CALL array_nullify (matrix%row_blk_size)
00264     CALL array_nullify (matrix%col_blk_size)
00265     CALL array_nullify (matrix%row_blk_offset)
00266     CALL array_nullify (matrix%col_blk_offset)
00267     CALL array_nullify (matrix%local_rows)
00268     CALL array_nullify (matrix%global_rows)
00269     CALL array_nullify (matrix%local_cols)
00270     CALL array_nullify (matrix%global_cols)
00271     NULLIFY (matrix%wms)
00272     NULLIFY (matrix%predistributed)
00273   END SUBROUTINE dbcsr_init_type
00274 
00275 ! *****************************************************************************
00280   SUBROUTINE dbcsr_init_obj (matrix)
00281     TYPE(dbcsr_obj), INTENT(OUT)             :: matrix
00282 
00283     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_init_obj', 
00284       routineP = moduleN//':'//routineN
00285 
00286 !TYPE(dbcsr_error_type), INTENT(INOUT)   :: error
00287 !   ---------------------------------------------------------------------------
00288 
00289     CALL dbcsr_init(matrix%m)
00290     matrix%m%initialized = 0
00291   END SUBROUTINE dbcsr_init_obj
00292 
00293 
00294 ! *****************************************************************************
00299   SUBROUTINE dbcsr_init_obj_p (matrix, error)
00300     TYPE(dbcsr_obj), POINTER                 :: matrix
00301     TYPE(dbcsr_error_type), INTENT(inout)    :: error
00302 
00303     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_init_obj_p', 
00304       routineP = moduleN//':'//routineN
00305 
00306 !   ---------------------------------------------------------------------------
00307 
00308     IF(ASSOCIATED(matrix)) CALL dbcsr_release_p(matrix, error)
00309     ALLOCATE(matrix)
00310     CALL dbcsr_init(matrix%m)
00311     matrix%m%initialized = 0
00312   END SUBROUTINE dbcsr_init_obj_p
00313 
00314 ! *****************************************************************************
00319   PURE FUNCTION dbcsr_valid_index_type (matrix) RESULT (valid_index)
00320     TYPE(dbcsr_type), INTENT(IN)             :: matrix
00321     LOGICAL                                  :: valid_index
00322 
00323 !   ---------------------------------------------------------------------------
00324 !valid_index = .FALSE.
00325 !IF (ASSOCIATED (matrix%row_p)) THEN
00326 !   valid_index = SIZE (matrix%row_p) .GT. 0
00327 !ENDIF
00328 
00329     valid_index = matrix%valid
00330   END FUNCTION dbcsr_valid_index_type
00331 
00332 ! *****************************************************************************
00337   PURE FUNCTION dbcsr_valid_index_obj (matrix) RESULT (valid_index)
00338     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
00339     LOGICAL                                  :: valid_index
00340 
00341 !   ---------------------------------------------------------------------------
00342 
00343     valid_index = dbcsr_valid_index_type (matrix%m)
00344     !valid_index = .FALSE.
00345     !IF (matrix%m%initialized .EQ. dbcsr_magic_number) THEN
00346     !   IF (ASSOCIATED (matrix%m%row_p)) THEN
00347     !      valid_index = SIZE (matrix%m%row_p) .GT. 0
00348     !   ENDIF
00349     !ELSE
00350     !   valid_index = .FALSE.
00351     !ENDIF
00352   END FUNCTION dbcsr_valid_index_obj
00353 
00354 ! *****************************************************************************
00359   PURE FUNCTION dbcsr_is_initialized_type (matrix) RESULT (initialized)
00360     TYPE(dbcsr_type), INTENT(IN)             :: matrix
00361     LOGICAL                                  :: initialized
00362 
00363 !   ---------------------------------------------------------------------------
00364 
00365     initialized = matrix%initialized .EQ. dbcsr_magic_number
00366   END FUNCTION dbcsr_is_initialized_type
00367 
00368 ! *****************************************************************************
00373   PURE FUNCTION dbcsr_is_initialized_obj (matrix) RESULT (initialized)
00374     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
00375     LOGICAL                                  :: initialized
00376 
00377 !   ---------------------------------------------------------------------------
00378 
00379     initialized = dbcsr_is_initialized_type (matrix%m)
00380   END FUNCTION dbcsr_is_initialized_obj
00381 
00382 
00383 ! *****************************************************************************
00387   SUBROUTINE dbcsr_hold (matrix)
00388     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
00389 
00390     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_hold', 
00391       routineP = moduleN//':'//routineN
00392 
00393     TYPE(dbcsr_error_type)                   :: error
00394 
00395 !   ---------------------------------------------------------------------------
00396 
00397     CALL dbcsr_assert (matrix%m%initialized .EQ. dbcsr_magic_number,&
00398          dbcsr_failure_level, dbcsr_caller_error, routineN,&
00399          "Matrix not initialized",__LINE__,error)
00400     matrix%m%refcount = matrix%m%refcount + 1
00401   END SUBROUTINE dbcsr_hold
00402 
00403 
00404 ! *****************************************************************************
00410   RECURSIVE SUBROUTINE dbcsr_release (matrix)
00411     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
00412 
00413     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_release', 
00414       routineP = moduleN//':'//routineN
00415 
00416     TYPE(dbcsr_error_type)                   :: error
00417 
00418 !   ---------------------------------------------------------------------------
00419 
00420     CALL dbcsr_assert (dbcsr_is_initialized (matrix),&
00421          dbcsr_warning_level, dbcsr_caller_error, routineN,&
00422          "Matrix not initialized",__LINE__,error)
00423     IF (matrix%m%initialized .EQ. dbcsr_magic_number) THEN
00424        matrix%m%refcount = matrix%m%refcount - 1
00425        IF (matrix%m%refcount .EQ. 0) THEN
00426           CALL dbcsr_destroy (matrix,error)
00427        ENDIF
00428     ENDIF
00429   END SUBROUTINE dbcsr_release
00430 
00431 
00432 ! *****************************************************************************
00438   SUBROUTINE dbcsr_release_p (matrix,error)
00439     TYPE(dbcsr_obj), POINTER                 :: matrix
00440     TYPE(dbcsr_error_type), INTENT(inout)    :: error
00441 
00442     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_release_p', 
00443       routineP = moduleN//':'//routineN
00444 
00445 !   ---------------------------------------------------------------------------
00446 
00447     IF(ASSOCIATED(matrix)) THEN
00448        CALL dbcsr_assert (dbcsr_is_initialized (matrix),&
00449             dbcsr_warning_level, dbcsr_caller_error, routineN,&
00450             "Matrix not initialized",__LINE__,error)
00451        IF (matrix%m%initialized .EQ. dbcsr_magic_number) THEN
00452           matrix%m%refcount = matrix%m%refcount - 1
00453           IF (matrix%m%refcount .EQ. 0) THEN
00454              CALL dbcsr_destroy (matrix, error)
00455              DEALLOCATE(matrix)
00456           ENDIF
00457        ENDIF
00458     ENDIF
00459   END SUBROUTINE dbcsr_release_p
00460 
00461 ! *****************************************************************************
00467   RECURSIVE SUBROUTINE dbcsr_destroy(matrix, error, force)
00468     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
00469     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
00470     LOGICAL, INTENT(IN), OPTIONAL            :: force
00471 
00472     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_destroy', 
00473       routineP = moduleN//':'//routineN
00474 
00475     INTEGER                                  :: error_handle
00476     LOGICAL                                  :: force_all
00477 
00478 !   ---------------------------------------------------------------------------
00479 
00480     CALL dbcsr_error_set (routineN, error_handle, error=error)
00481     CALL dbcsr_assert (dbcsr_is_initialized (matrix),&
00482          dbcsr_warning_level, dbcsr_caller_error, routineN,&
00483          "Can not destroy uninitialized matrix object.",__LINE__,error)
00484     force_all = .FALSE.
00485     IF (PRESENT (force)) force_all = force
00486     IF (dbcsr_is_initialized (matrix)) THEN
00487        CALL dbcsr_assert (force_all .OR. matrix%m%refcount.EQ.0,&
00488             dbcsr_warning_level, dbcsr_caller_error,&
00489             routineN, "You should not destroy referenced matrix.",__LINE__,error)
00490        CALL dbcsr_assert (.NOT.force_all .OR. matrix%m%refcount.LE.1,&
00491             dbcsr_warning_level, dbcsr_caller_error,&
00492             routineN, "You should not destroy referenced matrix.",__LINE__,error)
00493        IF (force_all .OR. matrix%m%refcount .EQ. 0) THEN
00494           CALL dbcsr_assert (.NOT. ASSOCIATED (matrix%m%wms), dbcsr_warning_level,&
00495                dbcsr_caller_error, routineN, "Destroying unfinalized matrix",__LINE__,error)
00496           IF (dbcsr_buffers_2d_needed) THEN
00497              CALL dbcsr_buffers_flush (matrix%m%buffers, error=error)
00498              CALL dbcsr_buffers_release (matrix%m%buffers, error=error)
00499           ENDIF
00500           CALL memory_deallocate(matrix%m%index, matrix%m%index_memory_type,&
00501                error=error)
00502           !CALL cp_error_init(my_error)
00503           IF (ASSOCIATED (matrix%m%predistributed)) THEN
00504              CALL dbcsr_assert (ASSOCIATED (matrix%m%predistributed%mats),&
00505                   dbcsr_fatal_level, dbcsr_internal_error, routineN,&
00506                   "Predistributed data should exist if allocated",__LINE__,error)
00507              CALL dbcsr_destroy_array (matrix%m%predistributed,error)
00508              DEALLOCATE (matrix%m%predistributed)
00509           ENDIF
00510           NULLIFY (matrix%m%predistributed)
00511           CALL dbcsr_data_release (matrix%m%data_area)
00512           CALL array_release (matrix%m%row_blk_size)
00513           CALL array_release (matrix%m%col_blk_size)
00514           CALL array_release (matrix%m%row_blk_offset)
00515           CALL array_release (matrix%m%col_blk_offset)
00516           CALL dbcsr_distribution_release(matrix%m%dist)
00517           IF (matrix%m%local_indexing) THEN
00518              CALL dbcsr_assert (array_exists (matrix%m%local_rows),&
00519                   dbcsr_fatal_level, dbcsr_internal_error, routineN,&
00520                   "Global row mapping should exist", __LINE__, error=error)
00521              CALL dbcsr_assert (array_exists (matrix%m%global_rows),&
00522                   dbcsr_fatal_level, dbcsr_internal_error, routineN,&
00523                   "Global row mapping should exist", __LINE__, error=error)
00524              CALL dbcsr_assert (array_exists (matrix%m%local_cols),&
00525                   dbcsr_fatal_level, dbcsr_internal_error, routineN,&
00526                   "Global column mapping should exist", __LINE__, error=error)
00527              CALL dbcsr_assert (array_exists (matrix%m%global_cols),&
00528                   dbcsr_fatal_level, dbcsr_internal_error, routineN,&
00529                   "Global column mapping should exist", __LINE__, error=error)
00530           ENDIF
00531           CALL dbcsr_release_locals (matrix, error)
00532           matrix%m%valid = .FALSE.
00533 !$        CALL OMP_DESTROY_LOCK (matrix%m%modification_lock)
00534           CALL dbcsr_init (matrix%m)
00535        ENDIF
00536     ENDIF
00537     CALL dbcsr_error_stop (error_handle, error=error)
00538   END SUBROUTINE dbcsr_destroy
00539 
00540   SUBROUTINE dbcsr_release_locals (matrix, error)
00541     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
00542     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
00543 
00544     IF (matrix%m%has_local_rows) &
00545          CALL array_release (matrix%m%local_rows)
00546     IF (matrix%m%has_global_rows) &
00547          CALL array_release (matrix%m%global_rows)
00548     IF (matrix%m%has_local_cols) &
00549          CALL array_release (matrix%m%local_cols)
00550     IF (matrix%m%has_global_cols) &
00551          CALL array_release (matrix%m%global_cols)
00552     matrix%m%has_local_rows  = .FALSE.
00553     matrix%m%has_global_rows = .FALSE.
00554     matrix%m%has_local_cols  = .FALSE.
00555     matrix%m%has_global_cols = .FALSE.
00556   END SUBROUTINE dbcsr_release_locals
00557 
00558 !  SUBROUTINE dbcsr_release_vlocals (matrix, error)
00559 !    TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
00560 !    TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
00561 !
00562 !    IF (matrix%m%has_local_vrows) &
00563 !         CALL array_release (matrix%m%local_vrows)
00564 !    IF (matrix%m%has_global_vrows) &
00565 !         CALL array_release (matrix%m%global_vrows)
00566 !    IF (matrix%m%has_local_vcols) &
00567 !         CALL array_release (matrix%m%local_vcols)
00568 !    IF (matrix%m%has_global_vcols) &
00569 !         CALL array_release (matrix%m%global_vcols)
00570 !    matrix%m%has_local_vrows  = .FALSE.
00571 !    matrix%m%has_global_vrows = .FALSE.
00572 !    matrix%m%has_local_vcols  = .FALSE.
00573 !    matrix%m%has_global_vcols = .FALSE.
00574 !  END SUBROUTINE dbcsr_release_vlocals
00575 !
00576 
00577 
00578 ! *****************************************************************************
00581   SUBROUTINE dbcsr_mp_init (mp_env)
00582     TYPE(dbcsr_mp_obj), INTENT(OUT)          :: mp_env
00583 
00584     NULLIFY (mp_env%mp)
00585   END SUBROUTINE dbcsr_mp_init
00586 
00587 ! *****************************************************************************
00590   FUNCTION dbcsr_mp_exists (mp_env) RESULT (exists)
00591     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00592     LOGICAL                                  :: exists
00593 
00594     exists = ASSOCIATED (mp_env%mp)
00595   END FUNCTION dbcsr_mp_exists
00596 
00597 ! *****************************************************************************
00600   FUNCTION dbcsr_mp_active (mp_env) RESULT (active)
00601     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00602     LOGICAL                                  :: active
00603 
00604     active = ASSOCIATED (mp_env%mp)
00605   END FUNCTION dbcsr_mp_active
00606 
00607 
00608 ! *****************************************************************************
00618   SUBROUTINE dbcsr_mp_new(mp_env, pgrid, mp_group, mynode, numnodes, myprow,&
00619        mypcol)
00620     TYPE(dbcsr_mp_obj), INTENT(OUT)          :: mp_env
00621     INTEGER, DIMENSION(0:, 0:), INTENT(IN)   :: pgrid
00622     INTEGER, INTENT(IN)                      :: mp_group, mynode
00623     INTEGER, INTENT(IN), OPTIONAL            :: numnodes, myprow, mypcol
00624 
00625     INTEGER                                  :: pcol, prow
00626 
00627 !   ---------------------------------------------------------------------------
00628 
00629      ALLOCATE(mp_env%mp)
00630      mp_env%mp%refcount = 1
00631      ALLOCATE (mp_env%mp%pgrid (0:SIZE(pgrid, 1)-1, 0:SIZE(pgrid, 2)-1 ))
00632      mp_env%mp%pgrid(:,:) = pgrid(:,:)
00633      mp_env%mp%mynode = mynode
00634      mp_env%mp%mp_group = mp_group
00635      IF (PRESENT (numnodes)) THEN
00636         mp_env%mp%numnodes = numnodes
00637      ELSE
00638         mp_env%mp%numnodes = SIZE (pgrid)
00639      ENDIF
00640      IF (PRESENT (myprow) .AND. PRESENT (mypcol)) THEN
00641         mp_env%mp%myprow = myprow
00642         mp_env%mp%mypcol = mypcol
00643      ELSE
00644         mp_env%mp%myprow = -33777
00645         mp_env%mp%mypcol = -33777
00646         column_loop: DO pcol = LBOUND (pgrid, 2), UBOUND (pgrid, 2)
00647            row_loop: DO prow = LBOUND (pgrid, 1), UBOUND (pgrid, 1)
00648               test_position: IF (pgrid (prow, pcol) .EQ. mynode) THEN
00649                  mp_env%mp%myprow = prow
00650                  mp_env%mp%mypcol = pcol
00651                  EXIT column_loop
00652               ENDIF test_position
00653            ENDDO row_loop
00654         ENDDO column_loop
00655      ENDIF
00656      mp_env%mp%subgroups_defined = .FALSE.
00657 !!!! KG workaround in place.
00658 !!!! This will be the replacement:
00659 !    ALLOCATE(mp_env%mp)
00660 !    mp_env%mp%refcount = 1
00661 !    ndims = 2
00662 !    dims(1:2) = (/ SIZE (pgrid, 1), SIZE (pgrid, 2) /)
00663 !    CALL mp_cart_create (mp_group, ndims,&
00664 !         dims, my_pos,&
00665 !         mp_env%mp%mp_group)
00666 !    CALL mp_environ (mp_env%mp%numnodes, mp_env%mp%mynode, mp_env%mp%mp_group)
00667 !    mp_env%mp%myprow = my_pos(1)
00668 !    mp_env%mp%mypcol = my_pos(2)
00669 !    ALLOCATE (mp_env%mp%pgrid (0:SIZE(pgrid, 1)-1, 0:SIZE(pgrid, 2)-1 ))
00670 !    column_loop: DO pcol = 0, SIZE (mp_env%mp%pgrid, 2)-1
00671 !       row_loop: DO prow = 0, SIZE (mp_env%mp%pgrid, 1)-1
00672 !          CALL mp_cart_rank (mp_env%mp%mp_group, (/ prow, pcol /),&
00673 !               mp_env%mp%pgrid (prow, pcol))
00674 !       ENDDO row_loop
00675 !    ENDDO column_loop
00676 !    mp_env%mp%subgroups_defined = .FALSE.
00677   END SUBROUTINE dbcsr_mp_new
00678 
00679 ! *****************************************************************************
00684   SUBROUTINE dbcsr_mp_grid_setup(mp_env, force)
00685     TYPE(dbcsr_mp_obj), INTENT(INOUT)        :: mp_env
00686     LOGICAL, INTENT(IN), OPTIONAL            :: force
00687 
00688     INTEGER                                  :: ndims, tmp_group
00689     INTEGER, DIMENSION(2)                    :: dims, my_pos
00690     LOGICAL                                  :: should_make
00691     LOGICAL, DIMENSION(2)                    :: remain
00692     TYPE(dbcsr_error_type)                   :: error
00693 
00694 !   ---------------------------------------------------------------------------
00695 ! Intel MPI Workaround
00696 !CALL dbcsr_assert (..NOT. mp_env%mp%subgroups_defined,&
00697 !     dbcsr_fatal_level, dbcsr_internal_error, "mp_grid_setup",&
00698 !     "Subcommunicators should be turned off")
00699 
00700     IF (PRESENT (force)) THEN
00701        should_make = use_subcommunicators .OR. force
00702     ELSE
00703        should_make = use_subcommunicators
00704     ENDIF
00705 
00706     IF (.NOT. mp_env%mp%subgroups_defined .AND. should_make) THEN
00707        ! KG workaround.
00708        ! This will be deleted (replaced by code in mp_new).
00709        ndims = 2
00710        dims(1:2) = (/ SIZE (mp_env%mp%pgrid, 1), SIZE (mp_env%mp%pgrid, 2) /)
00711        CALL mp_cart_create (mp_env%mp%mp_group, ndims,&
00712             dims, my_pos,&
00713             tmp_group)
00714        CALL dbcsr_assert (my_pos(1) .EQ. mp_env%mp%myprow,&
00715             dbcsr_fatal_level, dbcsr_internal_error, "mp_grid_setup",&
00716             "Got different MPI process grid",__LINE__,error)
00717        CALL dbcsr_assert (my_pos(2) .EQ. mp_env%mp%mypcol,&
00718             dbcsr_fatal_level, dbcsr_internal_error, "mp_grid_setup",&
00719             "Got different MPI process grid",__LINE__,error)
00720        !
00721        remain = (/.FALSE., .TRUE./)
00722        CALL mp_cart_sub (tmp_group, remain, mp_env%mp%prow_group)
00723        remain = (/.TRUE., .FALSE./)
00724        CALL mp_cart_sub (tmp_group, remain, mp_env%mp%pcol_group)
00725        CALL mp_comm_free (tmp_group)
00726        mp_env%mp%subgroups_defined = .TRUE.
00727     ENDIF
00728   END SUBROUTINE dbcsr_mp_grid_setup
00729 
00730 ! *****************************************************************************
00734   SUBROUTINE dbcsr_mp_grid_remove (mp_env)
00735     TYPE(dbcsr_mp_obj), INTENT(INOUT)        :: mp_env
00736 
00737     IF (mp_env%mp%subgroups_defined) THEN
00738        CALL mp_comm_free (mp_env%mp%prow_group)
00739        CALL mp_comm_free (mp_env%mp%pcol_group)
00740     ENDIF
00741   END SUBROUTINE dbcsr_mp_grid_remove
00742 
00743 ! *****************************************************************************
00747   SUBROUTINE dbcsr_modify_lock (matrix)
00748     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
00749 
00750     !$ CALL OMP_SET_LOCK (matrix%m%modification_lock)
00751   END SUBROUTINE dbcsr_modify_lock
00752 
00753 
00754 ! *****************************************************************************
00758   SUBROUTINE dbcsr_modify_unlock (matrix)
00759     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
00760 
00761     !$ CALL OMP_UNSET_LOCK (matrix%m%modification_lock)
00762   END SUBROUTINE dbcsr_modify_unlock
00763 
00764 ! *****************************************************************************
00768   SUBROUTINE dbcsr_read_lock (matrix)
00769     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
00770 
00771     !$ CALL OMP_SET_LOCK (matrix%m%modification_lock)
00772   END SUBROUTINE dbcsr_read_lock
00773 
00774 ! *****************************************************************************
00778   SUBROUTINE dbcsr_read_unlock (matrix)
00779     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
00780 
00781     !$ CALL OMP_UNSET_LOCK (matrix%m%modification_lock)
00782   END SUBROUTINE dbcsr_read_unlock
00783 
00784 
00785 ! *****************************************************************************
00789   PURE SUBROUTINE dbcsr_mp_hold(mp_env)
00790     TYPE(dbcsr_mp_obj), INTENT(INOUT)        :: mp_env
00791 
00792 !   ---------------------------------------------------------------------------
00793 
00794     mp_env%mp%refcount = mp_env%mp%refcount+1
00795   END SUBROUTINE dbcsr_mp_hold
00796 
00797 ! *****************************************************************************
00801   SUBROUTINE dbcsr_mp_release(mp_env)
00802     TYPE(dbcsr_mp_obj), INTENT(INOUT)        :: mp_env
00803 
00804 !   ---------------------------------------------------------------------------
00805 
00806     IF (ASSOCIATED (mp_env%mp)) THEN
00807        mp_env%mp%refcount = mp_env%mp%refcount - 1
00808        IF (mp_env%mp%refcount .LE. 0) THEN
00809           CALL dbcsr_mp_grid_remove (mp_env)
00810           ! KG workaround
00811           !CALL mp_comm_free (mp_env%mp%mp_group)
00812           DEALLOCATE (mp_env%mp%pgrid)
00813           DEALLOCATE (mp_env%mp)
00814           NULLIFY (mp_env%mp)
00815        ENDIF
00816     ENDIF
00817   END SUBROUTINE dbcsr_mp_release
00818 
00819   PURE FUNCTION dbcsr_mp_get_process(mp_env, prow, pcol) RESULT (process)
00820     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00821     INTEGER, INTENT(IN)                      :: prow, pcol
00822     INTEGER                                  :: process
00823 
00824     process = mp_env%mp%pgrid(prow, pcol)
00825   END FUNCTION dbcsr_mp_get_process
00826 
00827   PURE SUBROUTINE dbcsr_mp_get_coordinates(mp_env, node, prow, pcol)
00828     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00829     INTEGER, INTENT(IN)                      :: node
00830     INTEGER, INTENT(OUT)                     :: prow, pcol
00831 
00832     INTEGER                                  :: tmp
00833 
00834     pcol = MOD (node, SIZE(mp_env%mp%pgrid,2))
00835     prow = node / SIZE(mp_env%mp%pgrid,2)
00836     IF (mp_env%mp%pgrid(prow, pcol) .NE. node) THEN
00837        tmp = pcol ; pcol = prow ; prow = tmp
00838     ENDIF
00839   END SUBROUTINE dbcsr_mp_get_coordinates
00840 
00841   FUNCTION dbcsr_mp_pgrid(mp_env) RESULT (pgrid)
00842     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00843     INTEGER, DIMENSION(:, :), POINTER        :: pgrid
00844 
00845     pgrid => mp_env%mp%pgrid
00846   END FUNCTION dbcsr_mp_pgrid
00847   PURE FUNCTION dbcsr_mp_numnodes(mp_env) RESULT (numnodes)
00848     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00849     INTEGER                                  :: numnodes
00850 
00851     numnodes = mp_env%mp%numnodes
00852   END FUNCTION dbcsr_mp_numnodes
00853   PURE FUNCTION dbcsr_mp_mynode(mp_env) RESULT (mynode)
00854     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00855     INTEGER                                  :: mynode
00856 
00857     mynode = mp_env%mp%mynode
00858   END FUNCTION dbcsr_mp_mynode
00859   PURE FUNCTION dbcsr_mp_group(mp_env) RESULT (mp_group)
00860     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00861     INTEGER                                  :: mp_group
00862 
00863     mp_group = mp_env%mp%mp_group
00864   END FUNCTION dbcsr_mp_group
00865   PURE FUNCTION dbcsr_mp_nprows(mp_env) RESULT (nprows)
00866     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00867     INTEGER                                  :: nprows
00868 
00869     nprows = SIZE (mp_env%mp%pgrid, 1)
00870   END FUNCTION dbcsr_mp_nprows
00871   PURE FUNCTION dbcsr_mp_npcols(mp_env) RESULT (npcols)
00872     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00873     INTEGER                                  :: npcols
00874 
00875     npcols = SIZE (mp_env%mp%pgrid, 2)
00876   END FUNCTION dbcsr_mp_npcols
00877   PURE FUNCTION dbcsr_mp_myprow(mp_env) RESULT (myprow)
00878     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00879     INTEGER                                  :: myprow
00880 
00881     myprow = mp_env%mp%myprow
00882   END FUNCTION dbcsr_mp_myprow
00883   PURE FUNCTION dbcsr_mp_mypcol(mp_env) RESULT (mypcol)
00884     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00885     INTEGER                                  :: mypcol
00886 
00887     mypcol = mp_env%mp%mypcol
00888   END FUNCTION dbcsr_mp_mypcol
00889   PURE FUNCTION dbcsr_mp_has_subgroups(mp_env) RESULT (has_subgroups)
00890     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00891     LOGICAL                                  :: has_subgroups
00892 
00893     has_subgroups = mp_env%mp%subgroups_defined
00894   END FUNCTION dbcsr_mp_has_subgroups
00895   PURE FUNCTION dbcsr_mp_my_row_group(mp_env) RESULT (row_group)
00896     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00897     INTEGER                                  :: row_group
00898 
00899     row_group = mp_env%mp%prow_group
00900   END FUNCTION dbcsr_mp_my_row_group
00901   PURE FUNCTION dbcsr_mp_my_col_group(mp_env) RESULT (col_group)
00902     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00903     INTEGER                                  :: col_group
00904 
00905     col_group = mp_env%mp%pcol_group
00906   END FUNCTION dbcsr_mp_my_col_group
00907   FUNCTION dbcsr_mp_pgrid_equal(mp_env1, mp_env2) RESULT (equal_pgrid)
00908     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env1, mp_env2
00909     LOGICAL                                  :: equal_pgrid
00910 
00911     INTEGER, DIMENSION(:, :), POINTER        :: pgrid1, pgrid2
00912 
00913     IF (dbcsr_mp_nprows (mp_env1) .EQ. dbcsr_mp_nprows (mp_env2)&
00914          .AND. dbcsr_mp_nprows (mp_env1) .EQ. dbcsr_mp_nprows (mp_env2)) THEN
00915        pgrid1 => dbcsr_mp_pgrid (mp_env1)
00916        pgrid2 => dbcsr_mp_pgrid (mp_env2)
00917        IF (ALL (pgrid1 .EQ. pgrid2)) THEN
00918           equal_pgrid = .TRUE.
00919        ELSE
00920           equal_pgrid = .FALSE.
00921        ENDIF
00922     ELSE
00923        equal_pgrid = .FALSE.
00924     ENDIF
00925   END FUNCTION dbcsr_mp_pgrid_equal
00926 
00927 ! *****************************************************************************
00932   SUBROUTINE dbcsr_mp_new_transposed(mp_t, mp)
00933     TYPE(dbcsr_mp_obj), INTENT(OUT)          :: mp_t
00934     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp
00935 
00936 !   ---------------------------------------------------------------------------
00937 
00938     CALL dbcsr_mp_new (mp_t, TRANSPOSE (dbcsr_mp_pgrid (mp)),&
00939          dbcsr_mp_group (mp),&
00940          dbcsr_mp_mynode (mp), dbcsr_mp_numnodes (mp),&
00941          dbcsr_mp_mypcol (mp), dbcsr_mp_myprow (mp))
00942   END SUBROUTINE dbcsr_mp_new_transposed
00943 
00944 
00945 ! *****************************************************************************
00951   SUBROUTINE dbcsr_distribution_new(dist, mp_env, row_dist, col_dist,&
00952        local_rows, local_cols)
00953     TYPE(dbcsr_distribution_obj), 
00954       INTENT(OUT)                            :: dist
00955     TYPE(dbcsr_mp_obj), INTENT(IN)           :: mp_env
00956     TYPE(array_i1d_obj), INTENT(IN)          :: row_dist, col_dist
00957     TYPE(array_i1d_obj), INTENT(IN), 
00958       OPTIONAL                               :: local_rows, local_cols
00959 
00960     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_distribution_new', 
00961       routineP = moduleN//':'//routineN
00962 
00963     INTEGER                                  :: i, mypcoor, seq
00964     INTEGER, DIMENSION(:), POINTER           :: dd, ld
00965     TYPE(dbcsr_error_type)                   :: error
00966 
00967 !   ---------------------------------------------------------------------------
00968 
00969     ALLOCATE (dist%d)
00970     dist%d%refcount = 1
00971     dist%d%row_dist = row_dist
00972     CALL array_hold (dist%d%row_dist)
00973     dist%d%col_dist = col_dist
00974     CALL array_hold (dist%d%col_dist)
00975     dist%d%mp_env = mp_env
00976     CALL dbcsr_mp_hold (dist%d%mp_env)
00977     ! Verify given process row distribution.
00978     dd => array_data (row_dist)
00979     CALL dbcsr_assert (MAXVAL (dd), "LT", dbcsr_mp_nprows (mp_env),&
00980          dbcsr_fatal_level, dbcsr_wrong_args_error, routineN,&
00981          "A process row is too big for process grid",&
00982          __LINE__, error=error)
00983     ! Verify given process column distribution.
00984     dd => array_data (col_dist)
00985     CALL dbcsr_assert (MAXVAL (dd), "LT", dbcsr_mp_npcols (mp_env),&
00986          dbcsr_fatal_level, dbcsr_wrong_args_error, routineN,&
00987          "A process column is too big for process grid",&
00988          __LINE__, error=error)
00989     IF (PRESENT (local_rows)) THEN
00990        dist%d%local_rows = local_rows
00991        CALL array_hold (dist%d%local_rows)
00992     ELSE
00993        dd => array_data (row_dist)
00994        mypcoor = dbcsr_mp_myprow (mp_env)
00995        i = COUNT (dd.EQ.mypcoor)
00996        NULLIFY (ld)
00997        ALLOCATE (ld(i))
00998        seq = 1
00999        DO i = 1, array_size (row_dist)
01000           IF (dd(i) .EQ. mypcoor) THEN
01001              ld(seq) = i
01002              seq = seq+1
01003           ENDIF
01004        ENDDO
01005        CALL array_new (dist%d%local_rows, ld, gift=.TRUE.)
01006     ENDIF
01007     IF (PRESENT (local_cols)) THEN
01008        dist%d%local_cols = local_cols
01009        CALL array_hold (dist%d%local_cols)
01010     ELSE
01011        dd => array_data (col_dist)
01012        mypcoor = dbcsr_mp_mypcol (mp_env)
01013        i = COUNT (dd.EQ.mypcoor)
01014        NULLIFY (ld)
01015        ALLOCATE (ld(i))
01016        seq = 1
01017        DO i = 1, array_size (col_dist)
01018           IF (dd(i) .EQ. mypcoor) THEN
01019              ld(seq) = i
01020              seq = seq+1
01021           ENDIF
01022        ENDDO
01023        CALL array_new (dist%d%local_cols, ld, gift=.TRUE.)
01024     ENDIF
01025     dist%d%num_threads = 1
01026 !$  dist%d%num_threads = OMP_GET_MAX_THREADS()
01027     dist%d%has_thread_dist = .FALSE.
01028     CALL array_nullify (dist%d%thread_dist)
01029     CALL array_nullify (dist%d%row_map)
01030     CALL array_nullify (dist%d%col_map)
01031     NULLIFY (dist%d%other_l_rows)
01032     NULLIFY (dist%d%other_l_cols)
01033     dist%d%has_other_l_rows = .FALSE.
01034     dist%d%has_other_l_cols = .FALSE.
01035     CALL array_nullify (dist%d%global_row_map)
01036     CALL array_nullify (dist%d%global_col_map)
01037     dist%d%has_global_row_map = .FALSE.
01038     dist%d%has_global_col_map = .FALSE.
01039   END SUBROUTINE dbcsr_distribution_new
01040 
01041 ! *****************************************************************************
01045   SUBROUTINE dbcsr_distribution_hold(dist)
01046     TYPE(dbcsr_distribution_obj), 
01047       INTENT(INOUT)                          :: dist
01048 
01049 !   ---------------------------------------------------------------------------
01050 
01051     dist%d%refcount = dist%d%refcount + 1
01052   END SUBROUTINE dbcsr_distribution_hold
01053 
01054 ! *****************************************************************************
01058   SUBROUTINE dbcsr_distribution_release(dist)
01059     TYPE(dbcsr_distribution_obj), 
01060       INTENT(INOUT)                          :: dist
01061 
01062     TYPE(dbcsr_error_type)                   :: error
01063 
01064 !   ---------------------------------------------------------------------------
01065 
01066     IF (ASSOCIATED (dist%d)) THEN
01067        dist%d%refcount = dist%d%refcount - 1
01068        IF (dist%d%refcount .EQ. 0) THEN
01069           CALL array_release (dist%d%row_dist)
01070           CALL array_release (dist%d%col_dist)
01071           CALL array_release (dist%d%local_rows)
01072           CALL array_release (dist%d%local_cols)
01073           CALL dbcsr_mp_release (dist%d%mp_env)
01074           IF (dist%d%has_thread_dist) &
01075                CALL array_release (dist%d%thread_dist)
01076           CALL array_release (dist%d%row_map)
01077           CALL array_release (dist%d%col_map)
01078           CALL dbcsr_dist_release_locals (dist, error)
01079           DEALLOCATE (dist%d)
01080           CALL dbcsr_distribution_init (dist)
01081        ENDIF
01082     ENDIF
01083   END SUBROUTINE dbcsr_distribution_release
01084 
01085   SUBROUTINE dbcsr_dist_release_locals (dist, error)
01086     TYPE(dbcsr_distribution_obj), 
01087       INTENT(INOUT)                          :: dist
01088     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
01089 
01090     INTEGER                                  :: i
01091 
01092     IF (dist%d%has_other_l_rows) THEN
01093        DO i = LBOUND (dist%d%other_l_rows,1), UBOUND (dist%d%other_l_rows,1)
01094           CALL array_release (dist%d%other_l_rows(i))
01095        END DO
01096        DEALLOCATE (dist%d%other_l_rows)
01097        NULLIFY (dist%d%other_l_rows)
01098     ENDIF
01099     IF (dist%d%has_other_l_cols) THEN
01100        DO i = LBOUND (dist%d%other_l_cols,1), UBOUND (dist%d%other_l_cols,1)
01101           CALL array_release (dist%d%other_l_cols(i))
01102        END DO
01103        DEALLOCATE (dist%d%other_l_cols)
01104        NULLIFY (dist%d%other_l_cols)
01105     ENDIF
01106     IF (dist%d%has_global_row_map) THEN
01107        CALL array_release (dist%d%global_row_map)
01108     ENDIF
01109     IF (dist%d%has_global_col_map) THEN
01110        CALL array_release (dist%d%global_col_map)
01111     ENDIF
01112     dist%d%has_other_l_rows = .FALSE.
01113     dist%d%has_other_l_cols = .FALSE.
01114     dist%d%has_global_row_map = .FALSE.
01115     dist%d%has_global_col_map = .FALSE.
01116   END SUBROUTINE dbcsr_dist_release_locals
01117 
01118 
01119 
01120   SUBROUTINE dbcsr_distribution_init (dist)
01121     TYPE(dbcsr_distribution_obj), 
01122       INTENT(OUT)                            :: dist
01123 
01124     NULLIFY (dist%d)
01125   END SUBROUTINE dbcsr_distribution_init
01126 
01127   SUBROUTINE dbcsr_distribution_add_row_map (dist, row_map)
01128     TYPE(dbcsr_distribution_obj), 
01129       INTENT(INOUT)                          :: dist
01130     TYPE(array_i1d_obj), INTENT(IN)          :: row_map
01131 
01132     dist%d%row_map = row_map
01133     CALL array_hold (dist%d%row_map)
01134   END SUBROUTINE dbcsr_distribution_add_row_map
01135   SUBROUTINE dbcsr_distribution_add_col_map (dist, col_map)
01136     TYPE(dbcsr_distribution_obj), 
01137       INTENT(INOUT)                          :: dist
01138     TYPE(array_i1d_obj), INTENT(IN)          :: col_map
01139 
01140     dist%d%col_map = col_map
01141     CALL array_hold (dist%d%col_map)
01142   END SUBROUTINE dbcsr_distribution_add_col_map
01143   SUBROUTINE dbcsr_distribution_del_row_map (dist)
01144     TYPE(dbcsr_distribution_obj), 
01145       INTENT(INOUT)                          :: dist
01146 
01147     CALL array_release (dist%d%row_map)
01148     CALL array_nullify (dist%d%row_map)
01149   END SUBROUTINE dbcsr_distribution_del_row_map
01150   SUBROUTINE dbcsr_distribution_del_col_map (dist)
01151     TYPE(dbcsr_distribution_obj), 
01152       INTENT(INOUT)                          :: dist
01153 
01154     CALL array_release (dist%d%col_map)
01155     CALL array_nullify (dist%d%col_map)
01156   END SUBROUTINE dbcsr_distribution_del_col_map
01157 
01158   FUNCTION dbcsr_distribution_mp(dist) RESULT (mp_env)
01159     TYPE(dbcsr_distribution_obj), INTENT(IN) :: dist
01160     TYPE(dbcsr_mp_obj)                       :: mp_env
01161 
01162 !   ---------------------------------------------------------------------------
01163 
01164     mp_env = dist%d%mp_env
01165   END FUNCTION dbcsr_distribution_mp
01166   PURE FUNCTION dbcsr_distribution_nrows(dist) RESULT (nrows)
01167     TYPE(dbcsr_distribution_obj), INTENT(IN) :: dist
01168     INTEGER                                  :: nrows
01169 
01170     nrows = array_size (dist%d%row_dist)
01171   END FUNCTION dbcsr_distribution_nrows
01172   PURE FUNCTION dbcsr_distribution_ncols(dist) RESULT (ncols)
01173     TYPE(dbcsr_distribution_obj), INTENT(IN) :: dist
01174     INTEGER                                  :: ncols
01175 
01176     ncols = array_size (dist%d%col_dist)
01177   END FUNCTION dbcsr_distribution_ncols
01178   FUNCTION dbcsr_distribution_row_dist(dist) RESULT (row_dist)
01179     TYPE(dbcsr_distribution_obj), INTENT(IN) :: dist
01180     TYPE(array_i1d_obj)                      :: row_dist
01181 
01182 !   ---------------------------------------------------------------------------
01183 
01184     row_dist = dist%d%row_dist
01185   END FUNCTION dbcsr_distribution_row_dist
01186   FUNCTION dbcsr_distribution_col_dist(dist) RESULT (col_dist)
01187     TYPE(dbcsr_distribution_obj), INTENT(IN) :: dist
01188     TYPE(array_i1d_obj)                      :: col_dist
01189 
01190 !   ---------------------------------------------------------------------------
01191 
01192     col_dist = dist%d%col_dist
01193   END FUNCTION dbcsr_distribution_col_dist
01194 
01195   PURE FUNCTION dbcsr_distribution_nlocal_rows(dist) RESULT (nlocalrows)
01196     TYPE(dbcsr_distribution_obj), INTENT(IN) :: dist
01197     INTEGER                                  :: nlocalrows
01198 
01199     nlocalrows = array_size (dist%d%local_rows)
01200   END FUNCTION dbcsr_distribution_nlocal_rows
01201   PURE FUNCTION dbcsr_distribution_nlocal_cols(dist) RESULT (nlocalcols)
01202     TYPE(dbcsr_distribution_obj), INTENT(IN) :: dist
01203     INTEGER                                  :: nlocalcols
01204 
01205     nlocalcols = array_size (dist%d%local_cols)
01206   END FUNCTION dbcsr_distribution_nlocal_cols
01207   FUNCTION dbcsr_distribution_local_rows(dist) RESULT (local_rows)
01208     TYPE(dbcsr_distribution_obj), INTENT(IN) :: dist
01209     TYPE(array_i1d_obj)                      :: local_rows
01210 
01211     local_rows = dist%d%local_rows
01212   END FUNCTION dbcsr_distribution_local_rows
01213   FUNCTION dbcsr_distribution_local_cols(dist) RESULT (local_cols)
01214     TYPE(dbcsr_distribution_obj), INTENT(IN) :: dist
01215     TYPE(array_i1d_obj)                      :: local_cols
01216 
01217     local_cols = dist%d%local_cols
01218   END FUNCTION dbcsr_distribution_local_cols
01219   !
01220   PURE FUNCTION dbcsr_distribution_processor(dist, row, col)&
01221        RESULT (processor)
01222     TYPE(dbcsr_distribution_obj), INTENT(IN) :: dist
01223     INTEGER, INTENT(IN)                      :: row, col
01224     INTEGER                                  :: processor
01225 
01226     INTEGER                                  :: c, r
01227 
01228     IF (ASSOCIATED (dist%d%row_map%low)) THEN ! instead of array_exists
01229        r = dist%d%row_map%low%data(row)
01230     ELSE
01231        r = row
01232     ENDIF
01233     IF (ASSOCIATED (dist%d%col_map%low)) THEN ! instead of array_exists
01234        c = dist%d%col_map%low%data(col)
01235     ELSE
01236        c = col
01237     ENDIF
01238     processor = dist%d%mp_env%mp%pgrid(dist%d%row_dist%low%data(r),&
01239          dist%d%col_dist%low%data(c))
01240   END FUNCTION dbcsr_distribution_processor
01241 
01242   FUNCTION dbcsr_distribution_thread_dist(dist) RESULT (thread_dist)
01243     TYPE(dbcsr_distribution_obj), INTENT(IN) :: dist
01244     TYPE(array_i1d_obj)                      :: thread_dist
01245 
01246 !   ---------------------------------------------------------------------------
01247 
01248     thread_dist = dist%d%thread_dist
01249   END FUNCTION dbcsr_distribution_thread_dist
01250 
01251 
01252   PURE FUNCTION dbcsr_distribution_has_threads(dist) RESULT (has_thread_dist)
01253     TYPE(dbcsr_distribution_obj), INTENT(IN) :: dist
01254     LOGICAL                                  :: has_thread_dist
01255 
01256 !   ---------------------------------------------------------------------------
01257 
01258     has_thread_dist = dist%d%has_thread_dist
01259   END FUNCTION dbcsr_distribution_has_threads
01260 
01261   PURE FUNCTION dbcsr_distribution_num_threads(dist) RESULT (num_threads)
01262     TYPE(dbcsr_distribution_obj), INTENT(IN) :: dist
01263     INTEGER                                  :: num_threads
01264 
01265 !   ---------------------------------------------------------------------------
01266 
01267     num_threads = dist%d%num_threads
01268   END FUNCTION dbcsr_distribution_num_threads
01269 
01270 ! *****************************************************************************
01275   SUBROUTINE dbcsr_distribution_make_threads(dist, row_sizes)
01276     TYPE(dbcsr_distribution_obj), 
01277       INTENT(INOUT), TARGET                  :: dist
01278     INTEGER, DIMENSION(:), INTENT(IN), 
01279       OPTIONAL                               :: row_sizes
01280 
01281     TYPE(dbcsr_distribution_obj), POINTER    :: dist_p
01282 
01283 !   ---------------------------------------------------------------------------
01284 
01285     dist_p => dist
01286     !$ IF (.NOT. OMP_IN_PARALLEL ()) THEN
01287        !$OMP PARALLEL DEFAULT(NONE) &
01288        !$OMP          SHARED(dist_p,row_sizes)
01289        !$    CALL make_threads (dist_p, row_sizes=row_sizes)
01290        !$OMP END PARALLEL
01291     !$ ELSE
01292        CALL make_threads (dist_p, row_sizes=row_sizes)
01293        !$OMP BARRIER
01294     !$ ENDIF
01295   END SUBROUTINE dbcsr_distribution_make_threads
01296 
01297 ! *****************************************************************************
01308   SUBROUTINE make_threads(dist, row_sizes)
01309     TYPE(dbcsr_distribution_obj), POINTER    :: dist
01310     INTEGER, DIMENSION(:), INTENT(IN), 
01311       OPTIONAL                               :: row_sizes
01312 
01313     CHARACTER(len=*), PARAMETER :: routineN = 'make_threads', 
01314       routineP = moduleN//':'//routineN
01315 
01316     INTEGER :: block_size, block_size0, cur_block, group_size, i, last_row, 
01317       nlrows, nrows, nthreads, row, t, t_cnt
01318     INTEGER, ALLOCATABLE, DIMENSION(:)       :: itemp1, itemp2, reorder, 
01319                                                 sorted_row_sizes
01320     INTEGER, DIMENSION(:), POINTER           :: lrows, td
01321     LOGICAL                                  :: assigned, found, heap_error
01322     REAL(kind=sp)                            :: load_fraction, rn, soft_thr
01323     TYPE(dbcsr_error_type)                   :: error
01324     TYPE(heap_t)                             :: t_heap
01325 
01326 !   ---------------------------------------------------------------------------
01327 
01328     nthreads = 1
01329     !$ nthreads = OMP_GET_NUM_THREADS () ;
01330 
01331     !$  CALL dbcsr_assert (dist%d%num_threads, "EQ", nthreads,&
01332     !$     dbcsr_fatal_level, dbcsr_internal_error, routineN,&
01333     !$     "Thread number has changed", __LINE__, error=error)
01334     nrows = dbcsr_distribution_nrows (dist)
01335     nlrows = dbcsr_distribution_nlocal_rows (dist)
01336     lrows => array_data (dbcsr_distribution_local_rows (dist))
01337 
01338     !$OMP BARRIER
01339     !$OMP MASTER
01340 
01341     IF (comm_thread_load .NE. 0) THEN
01342        load_fraction = REAL(comm_thread_load)/100.0
01343     ELSE
01344        load_fraction = 1.0
01345     END IF
01346 
01347     IF (.NOT. dist%d%has_thread_dist) THEN
01348        dist%d%num_threads = nthreads
01349        group_size = 0 ; cur_block = 0
01350 
01351        ALLOCATE (td(nrows))
01352        dist%d%has_thread_dist = .TRUE.
01353        CALL array_new (dist%d%thread_dist, td, gift=.TRUE.)
01354        td => array_data (dist%d%thread_dist)
01355 
01356        IF (PRESENT (row_sizes)) THEN
01357           ! The goal is to distribute rows to threads as equally as
01358           ! possible. The row sizes are first sorted. Each group of
01359           ! equally sized rows (group_size rows of size cur_block) is
01360           ! distributed to threads (keeping consecutive rows
01361           ! together). The group is divided into equally-sized blocks
01362           ! (block_size0, block_size).  Leftover rows (those that can
01363           ! not be equally distributed to threads) are then assigned
01364           ! to threads so that each thread's commulative load attempts
01365           ! to be equal. This distribution is achieved using a heap.
01366           !
01367           ! The heap is used to distribute "leftover"rows to threads.
01368           ! Lefotver rows are those of the same size that can not be
01369           ! evenly distributed among all threads.
01370           CALL heap_new (t_heap, nthreads-1, heap_error)
01371           ! We do not want thread 0 to be in the heap.
01372           ALLOCATE(itemp1(1:nthreads-1))
01373           ALLOCATE(itemp2(1:nthreads-1))
01374           DO i=1,nthreads-1
01375              itemp1(i)=i
01376              itemp2(i)=0
01377           ENDDO
01378           CALL heap_fill (t_heap,itemp1,itemp2,heap_error)
01379           DEALLOCATE(itemp1,itemp2)
01380           ALLOCATE (sorted_row_sizes (nrows))
01381           ALLOCATE (reorder (nrows))
01382           sorted_row_sizes(:) = row_sizes(:)
01383           CALL sort (sorted_row_sizes, nrows, reorder)
01384 
01385           row = 1
01386           DO WHILE ( row .LE. nrows)
01387              cur_block = sorted_row_sizes(nrows-row+1)
01388              assigned = .FALSE.
01389              group_size = 0
01390 
01391              last_row = nrows-row+1
01392              DO i = last_row, 1, -1
01393                 IF ( cur_block == sorted_row_sizes(i) ) THEN
01394                    group_size = group_size + 1
01395                    row = row + 1
01396                 ELSE
01397                    EXIT
01398                 END IF
01399              END DO
01400 
01401              soft_thr = load_fraction + nthreads - 1
01402              block_size0 = load_fraction*(group_size/soft_thr)
01403              block_size = group_size/soft_thr
01404 
01405              IF (.NOT. (block_size0 .EQ. 0 .AND. block_size .EQ. 0)) THEN
01406                 !blocks for master thread
01407                 td(reorder(last_row:last_row-block_size0+1:-1)) = 0
01408 
01409                 !Other threads
01410                 DO t=1, nthreads-1
01411                    td(reorder(last_row-block_size0-(t-1)*block_size:&
01412                         last_row-block_size0-(t)*block_size+1:-1)) = t
01413                 END DO
01414              END IF
01415 
01416              !Leftover bocks
01417              DO i=last_row-block_size0-(nthreads-1)*block_size, last_row+1-group_size, -1
01418                 CALL heap_get_first (t_heap, t, t_cnt, found,heap_error)
01419                 t_cnt = t_cnt + cur_block
01420                 CALL heap_reset_first (t_heap, t_cnt, heap_error)
01421                 td(reorder(i)) = t
01422              END DO
01423 
01424           END DO
01425           CALL heap_release (t_heap, heap_error)
01426           DEALLOCATE (sorted_row_sizes)
01427           DEALLOCATE (reorder)
01428        ELSE
01429           DO t = 1, nrows
01430              IF (.FALSE.) THEN
01431                 td(t) = MOD(t-1, nthreads)
01432              ELSE
01433                 CALL RANDOM_NUMBER (rn)
01434                 ! Makes sure the numbers are in the proper integer range.
01435                 td(t) = MOD (INT (rn*REAL(nthreads)), nthreads)
01436              ENDIF
01437           END DO
01438        ENDIF
01439     ENDIF
01440     !$OMP END MASTER
01441   END SUBROUTINE make_threads
01442 
01444   SUBROUTINE dbcsr_distribution_no_threads(dist)
01445     TYPE(dbcsr_distribution_obj), 
01446       INTENT(INOUT)                          :: dist
01447 
01448 !$OMP MASTER
01449     CALL array_release (dist%d%thread_dist)
01450     dist%d%has_thread_dist = .FALSE.
01451 !$OMP END MASTER
01452   END SUBROUTINE dbcsr_distribution_no_threads
01453 
01454 
01455 ! Pertaining to the dbcsr matrix.
01456 
01457   FUNCTION dbcsr_nblkrows_total(matrix) RESULT (nblkrows_total)
01458     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01459     INTEGER                                  :: nblkrows_total
01460 
01461     nblkrows_total = matrix%m%nblkrows_total
01462   END FUNCTION dbcsr_nblkrows_total
01463 
01464   FUNCTION dbcsr_nblkcols_total(matrix) RESULT (nblkcols_total)
01465     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01466     INTEGER                                  :: nblkcols_total
01467 
01468     nblkcols_total = matrix%m%nblkcols_total
01469   END FUNCTION dbcsr_nblkcols_total
01470   FUNCTION dbcsr_nfullrows_total(matrix) RESULT (nfullrows_total)
01471     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01472     INTEGER                                  :: nfullrows_total
01473 
01474     nfullrows_total = matrix%m%nfullrows_total
01475   END FUNCTION dbcsr_nfullrows_total
01476   FUNCTION dbcsr_nfullcols_total(matrix) RESULT (nfullcols_total)
01477     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01478     INTEGER                                  :: nfullcols_total
01479 
01480     nfullcols_total = matrix%m%nfullcols_total
01481   END FUNCTION dbcsr_nfullcols_total
01482   FUNCTION dbcsr_nblkrows_local(matrix) RESULT (nblkrows_local)
01483     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01484     INTEGER                                  :: nblkrows_local
01485 
01486     nblkrows_local = matrix%m%nblkrows_local
01487   END FUNCTION dbcsr_nblkrows_local
01488   FUNCTION dbcsr_nblkcols_local(matrix) RESULT (nblkcols_local)
01489     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01490     INTEGER                                  :: nblkcols_local
01491 
01492     nblkcols_local = matrix%m%nblkcols_local
01493   END FUNCTION dbcsr_nblkcols_local
01494   FUNCTION dbcsr_nfullrows_local(matrix) RESULT (nfullrows_local)
01495     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01496     INTEGER                                  :: nfullrows_local
01497 
01498     nfullrows_local = matrix%m%nfullrows_local
01499   END FUNCTION dbcsr_nfullrows_local
01500   FUNCTION dbcsr_nfullcols_local(matrix) RESULT (nfullcols_local)
01501     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01502     INTEGER                                  :: nfullcols_local
01503 
01504     nfullcols_local = matrix%m%nfullcols_local
01505   END FUNCTION dbcsr_nfullcols_local
01506   FUNCTION dbcsr_max_row_size(matrix) RESULT (max_row_size)
01507     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01508     INTEGER                                  :: max_row_size
01509 
01510     max_row_size = matrix%m%max_rbs
01511   END FUNCTION dbcsr_max_row_size
01512   FUNCTION dbcsr_max_col_size(matrix) RESULT (max_col_size)
01513     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01514     INTEGER                                  :: max_col_size
01515 
01516     max_col_size = matrix%m%max_cbs
01517   END FUNCTION dbcsr_max_col_size
01518 
01519   FUNCTION dbcsr_distribution (matrix) RESULT (distribution)
01520     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01521     TYPE(dbcsr_distribution_obj)             :: distribution
01522 
01523     distribution = matrix%m%dist
01524   END FUNCTION dbcsr_distribution
01525 
01526   FUNCTION dbcsr_name (matrix) RESULT (name)
01527     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01528     CHARACTER(len=default_string_length)     :: name
01529 
01530     name = matrix%m%name
01531   END FUNCTION dbcsr_name
01532 
01533 ! *****************************************************************************
01539   PURE FUNCTION dbcsr_wm_use_mutable (wm) RESULT (use_mutable)
01540     TYPE(dbcsr_work_type), INTENT(IN)        :: wm
01541     LOGICAL                                  :: use_mutable
01542 
01543 !   ---------------------------------------------------------------------------
01544 
01545     use_mutable = dbcsr_mutable_instantiated (wm%mutable)
01546   END FUNCTION dbcsr_wm_use_mutable
01547 
01548 ! *****************************************************************************
01554   PURE FUNCTION dbcsr_use_mutable (matrix) RESULT (use_mutable)
01555     TYPE(dbcsr_type), INTENT(IN)             :: matrix
01556     LOGICAL                                  :: use_mutable
01557 
01558 !   ---------------------------------------------------------------------------
01559 
01560     use_mutable = matrix%work_mutable
01561   END FUNCTION dbcsr_use_mutable
01562 
01563 
01564   FUNCTION dbcsr_row_block_sizes (matrix) RESULT (row_blk_sizes)
01565     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01566     TYPE(array_i1d_obj)                      :: row_blk_sizes
01567 
01568     row_blk_sizes = matrix%m%row_blk_size
01569   END FUNCTION dbcsr_row_block_sizes
01570 
01571   FUNCTION dbcsr_col_block_sizes (matrix) RESULT (col_blk_sizes)
01572     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01573     TYPE(array_i1d_obj)                      :: col_blk_sizes
01574 
01575     col_blk_sizes = matrix%m%col_blk_size
01576   END FUNCTION dbcsr_col_block_sizes
01577 
01578   FUNCTION dbcsr_col_block_offsets (matrix) RESULT (col_blk_offsets)
01579     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01580     TYPE(array_i1d_obj)                      :: col_blk_offsets
01581 
01582     col_blk_offsets = matrix%m%col_blk_offset
01583   END FUNCTION dbcsr_col_block_offsets
01584 
01585   FUNCTION dbcsr_row_block_offsets (matrix) RESULT (row_blk_offsets)
01586     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01587     TYPE(array_i1d_obj)                      :: row_blk_offsets
01588 
01589     row_blk_offsets = matrix%m%row_blk_offset
01590   END FUNCTION dbcsr_row_block_offsets
01591 
01592 ! *****************************************************************************
01618   SUBROUTINE dbcsr_get_info(matrix, nblkrows_total, nblkcols_total,&
01619        nfullrows_total, nfullcols_total,&
01620        nblkrows_local, nblkcols_local,&
01621        nfullrows_local, nfullcols_local,&
01622        my_prow, my_pcol,&
01623        local_rows, local_cols, proc_row_dist, proc_col_dist,&
01624        row_blk_size, col_blk_size, row_blk_offset, col_blk_offset, distribution, name, data_area,&
01625        matrix_type, data_type)
01626     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01627     INTEGER, INTENT(OUT), OPTIONAL :: nblkrows_total, nblkcols_total, 
01628       nfullrows_total, nfullcols_total, nblkrows_local, nblkcols_local, 
01629       nfullrows_local, nfullcols_local, my_prow, my_pcol
01630     INTEGER, DIMENSION(:), OPTIONAL, POINTER :: local_rows, local_cols, 
01631                                                 proc_row_dist, proc_col_dist
01632     TYPE(array_i1d_obj), INTENT(OUT), 
01633       OPTIONAL                               :: row_blk_size, col_blk_size, 
01634                                                 row_blk_offset, col_blk_offset
01635     TYPE(dbcsr_distribution_obj), 
01636       INTENT(OUT), OPTIONAL                  :: distribution
01637     CHARACTER(len=*), INTENT(OUT), OPTIONAL  :: name
01638     TYPE(dbcsr_data_obj), INTENT(OUT), 
01639       OPTIONAL                               :: data_area
01640     CHARACTER, OPTIONAL                      :: matrix_type
01641     INTEGER, OPTIONAL                        :: data_type
01642 
01643     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_get_info', 
01644       routineP = moduleN//':'//routineN
01645 
01646     TYPE(dbcsr_error_type)                   :: error
01647 
01648 !   ---------------------------------------------------------------------------
01649 
01650     CALL dbcsr_assert (matrix%m%initialized.EQ.dbcsr_magic_number,&
01651          dbcsr_fatal_level, dbcsr_caller_error, routineP,&
01652          "Matrix not initialized",__LINE__,error)
01653     !vw avoid massive priting of warnings
01654     !CALL cp_assert (matrix%m%valid, cp_warning_level, cp_caller_error,&
01655     !     routineP,"Invalid matrix")
01656     IF (PRESENT (nblkrows_total)) nblkrows_total = matrix%m%nblkrows_total
01657     IF (PRESENT (nblkcols_total)) nblkcols_total = matrix%m%nblkcols_total
01658     IF (PRESENT (nfullrows_total)) nfullrows_total = matrix%m%nfullrows_total
01659     IF (PRESENT (nfullcols_total)) nfullcols_total = matrix%m%nfullcols_total
01660     IF (PRESENT (nblkrows_local)) nblkrows_local = matrix%m%nblkrows_local
01661     IF (PRESENT (nblkcols_local)) nblkcols_local = matrix%m%nblkcols_local
01662     IF (PRESENT (nfullrows_local)) nfullrows_local = matrix%m%nfullrows_local
01663     IF (PRESENT (nfullcols_local)) nfullcols_local = matrix%m%nfullcols_local
01664     IF (PRESENT (row_blk_size)) row_blk_size = matrix%m%row_blk_size
01665     IF (PRESENT (col_blk_size)) col_blk_size = matrix%m%col_blk_size
01666     IF (PRESENT (row_blk_offset)) row_blk_offset = matrix%m%row_blk_offset
01667     IF (PRESENT (col_blk_offset)) col_blk_offset = matrix%m%col_blk_offset
01668     IF (PRESENT (distribution)) distribution = matrix%m%dist
01669     IF (PRESENT (name)) name = matrix%m%name
01670     IF (PRESENT (data_area)) data_area = matrix%m%data_area
01671     IF (PRESENT (matrix_type)) THEN
01672        matrix_type = dbcsr_get_matrix_type (matrix)
01673        CALL dbcsr_assert (matrix_type.NE.dbcsr_type_invalid, dbcsr_failure_level, dbcsr_caller_error,&
01674             routineN, "Incorrect symmetry",__LINE__,error)
01675     ENDIF
01676     IF (PRESENT (data_type)) data_type = matrix%m%data_type
01677     IF (PRESENT (local_rows)) &
01678          local_rows => array_data (dbcsr_distribution_local_rows (matrix%m%dist))
01679     IF (PRESENT (local_cols)) &
01680          local_cols => array_data (dbcsr_distribution_local_cols (matrix%m%dist))
01681     IF (PRESENT (proc_row_dist)) &
01682          proc_row_dist => array_data (dbcsr_distribution_row_dist (matrix%m%dist))
01683     IF (PRESENT (proc_col_dist)) &
01684          proc_col_dist => array_data (dbcsr_distribution_col_dist (matrix%m%dist))
01685     IF (PRESENT (my_prow)) &
01686        my_prow = dbcsr_mp_myprow (dbcsr_distribution_mp (matrix%m%dist))
01687     IF (PRESENT (my_pcol)) &
01688        my_pcol = dbcsr_mp_mypcol (dbcsr_distribution_mp (matrix%m%dist))
01689   END SUBROUTINE dbcsr_get_info
01690 
01691 ! *****************************************************************************
01697   FUNCTION dbcsr_may_be_dense (matrix, occ_thresh) RESULT (may_be_dense)
01698     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01699     REAL(real_8), INTENT(in)                 :: occ_thresh
01700     LOGICAL                                  :: may_be_dense
01701 
01702     REAL(real_8)                             :: occ
01703 
01704 !   ---------------------------------------------------------------------------
01705 
01706     occ = dbcsr_get_occupation(matrix)
01707     may_be_dense =  .NOT.( occ .LT. occ_thresh )
01708     ! make sure every proc sees the same
01709     CALL mp_sum(may_be_dense, dbcsr_mp_group (dbcsr_distribution_mp (matrix%m%dist)))
01710   END FUNCTION dbcsr_may_be_dense
01711 
01712 ! *****************************************************************************
01720   PURE FUNCTION dbcsr_blk_row_size_type (matrix, row) RESULT (row_size)
01721     TYPE(dbcsr_type), INTENT(IN)             :: matrix
01722     INTEGER, INTENT(IN)                      :: row
01723     INTEGER                                  :: row_size
01724 
01725     row_size = matrix%row_blk_size%low%data(row)
01726   END FUNCTION dbcsr_blk_row_size_type
01727   PURE FUNCTION dbcsr_blk_row_size_obj (matrix, row) RESULT (row_size)
01728     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01729     INTEGER, INTENT(IN)                      :: row
01730     INTEGER                                  :: row_size
01731 
01732     row_size = matrix%m%row_blk_size%low%data(row)
01733   END FUNCTION dbcsr_blk_row_size_obj
01734   PURE FUNCTION dbcsr_blk_row_offset (matrix, row) RESULT (row_offset)
01735     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01736     INTEGER, INTENT(IN)                      :: row
01737     INTEGER                                  :: row_offset
01738 
01739     row_offset = matrix%m%row_blk_offset%low%data(row)
01740   END FUNCTION dbcsr_blk_row_offset
01741 
01742 
01743 ! *****************************************************************************
01751   PURE FUNCTION dbcsr_blk_column_size_type (matrix, column) RESULT (column_size)
01752     TYPE(dbcsr_type), INTENT(IN)             :: matrix
01753     INTEGER, INTENT(IN)                      :: column
01754     INTEGER                                  :: column_size
01755 
01756     column_size = matrix%col_blk_size%low%data(column)
01757   END FUNCTION dbcsr_blk_column_size_type
01758   PURE FUNCTION dbcsr_blk_column_size_obj (matrix, column) RESULT (column_size)
01759     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01760     INTEGER, INTENT(IN)                      :: column
01761     INTEGER                                  :: column_size
01762 
01763     column_size = matrix%m%col_blk_size%low%data(column)
01764   END FUNCTION dbcsr_blk_column_size_obj
01765   PURE FUNCTION dbcsr_blk_col_offset (matrix, col) RESULT (col_offset)
01766     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01767     INTEGER, INTENT(IN)                      :: col
01768     INTEGER                                  :: col_offset
01769 
01770     col_offset = matrix%m%col_blk_offset%low%data(col)
01771   END FUNCTION dbcsr_blk_col_offset
01772 
01773 ! *****************************************************************************
01778   FUNCTION dbcsr_data_area (matrix) RESULT (data_area)
01779     TYPE(dbcsr_type), INTENT(IN)             :: matrix
01780     TYPE(dbcsr_data_obj)                     :: data_area
01781 
01782     data_area = matrix%data_area
01783   END FUNCTION dbcsr_data_area
01784 
01785 ! *****************************************************************************
01792   SUBROUTINE dbcsr_switch_data_area (matrix, data_area, previous_data_area, &
01793        error)
01794     TYPE(dbcsr_obj), INTENT(INOUT)           :: matrix
01795     TYPE(dbcsr_data_obj), INTENT(IN)         :: data_area
01796     TYPE(dbcsr_data_obj), INTENT(OUT), 
01797       OPTIONAL                               :: previous_data_area
01798     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
01799 
01800     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_switch_data_area', 
01801       routineP = moduleN//':'//routineN
01802 
01803     INTEGER                                  :: error_handler
01804     TYPE(dbcsr_block_buffer_obj)             :: buffers
01805 
01806 !   ---------------------------------------------------------------------------
01807 
01808     CALL dbcsr_error_set (routineN, error_handler, error)
01809     IF (dbcsr_buffers_2d_needed) THEN
01810        buffers = matrix%m%buffers
01811        CALL dbcsr_buffers_flush (buffers, error=error)
01812        CALL dbcsr_buffers_release (buffers, error=error)
01813        CALL dbcsr_buffers_init (buffers)
01814        CALL dbcsr_buffers_new (buffers, data_area, error=error)
01815        matrix%m%buffers = buffers
01816     ENDIF
01817     IF (PRESENT (previous_data_area)) THEN
01818        previous_data_area = matrix%m%data_area
01819     ELSE
01820        CALL dbcsr_data_release (matrix%m%data_area)
01821     ENDIF
01822     matrix%m%data_area = data_area
01823     CALL dbcsr_data_hold (matrix%m%data_area)
01824     CALL dbcsr_error_stop (error_handler, error)
01825   END SUBROUTINE dbcsr_switch_data_area
01826 
01827 
01828 ! *****************************************************************************
01834   PURE FUNCTION dbcsr_get_matrix_type (matrix) RESULT (matrix_type)
01835     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01836     CHARACTER                                :: matrix_type
01837 
01838     matrix_type = dbcsr_type_invalid
01839     IF (matrix%m%symmetry) THEN
01840        IF ( (.NOT.matrix%m%negate_real).AND. matrix%m%negate_imaginary ) THEN
01841           matrix_type = dbcsr_type_hermitian
01842        ELSEIF ( matrix%m%negate_real .AND. (.NOT.matrix%m%negate_imaginary) ) THEN
01843           matrix_type = dbcsr_type_antihermitian
01844        ELSEIF (matrix%m%negate_real .AND. matrix%m%negate_imaginary) THEN
01845           matrix_type = dbcsr_type_antisymmetric
01846        ELSEIF ( (.NOT.matrix%m%negate_real) .AND. (.NOT.matrix%m%negate_imaginary)) THEN
01847           matrix_type = dbcsr_type_symmetric
01848        ENDIF
01849     ELSE
01850        matrix_type = dbcsr_type_no_symmetry
01851     ENDIF
01852   END FUNCTION dbcsr_get_matrix_type
01853 
01854 ! *****************************************************************************
01859   PURE FUNCTION dbcsr_has_symmetry (matrix) RESULT (has_symmetry)
01860     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01861     LOGICAL                                  :: has_symmetry
01862 
01863     has_symmetry = matrix%m%symmetry
01864   END FUNCTION dbcsr_has_symmetry
01865 
01866 ! *****************************************************************************
01877 
01878 ! *****************************************************************************
01884   PURE FUNCTION dbcsr_get_replication_type (matrix) RESULT (repl_type)
01885     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01886     CHARACTER                                :: repl_type
01887 
01888     repl_type = matrix%m%replication_type
01889   END FUNCTION dbcsr_get_replication_type
01890 
01891 ! *****************************************************************************
01897   PURE FUNCTION dbcsr_get_data_type_obj (matrix) RESULT (data_type)
01898     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01899     INTEGER                                  :: data_type
01900 
01901     data_type = matrix%m%data_type
01902   END FUNCTION dbcsr_get_data_type_obj
01903 
01904 ! *****************************************************************************
01909   PURE FUNCTION data_get_data_type (area) RESULT (data_type)
01910     TYPE(dbcsr_data_obj), INTENT(IN)         :: area
01911     INTEGER                                  :: data_type
01912 
01913     data_type = dbcsr_data_get_type (area)
01914   END FUNCTION data_get_data_type
01915 
01916 
01917 ! *****************************************************************************
01922   PURE FUNCTION dbcsr_is_real_obj (matrix) RESULT (is_real)
01923     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01924     LOGICAL                                  :: is_real
01925 
01926     IF(matrix%m%data_type.EQ.dbcsr_type_real_4.OR.&
01927        matrix%m%data_type.EQ.dbcsr_type_real_8) THEN
01928        is_real = .TRUE.
01929     ELSE
01930        is_real = .FALSE.
01931     ENDIF
01932   END FUNCTION dbcsr_is_real_obj
01933 
01934 ! *****************************************************************************
01939   PURE FUNCTION dbcsr_is_complex_obj (matrix) RESULT (is_complex)
01940     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01941     LOGICAL                                  :: is_complex
01942 
01943     IF(matrix%m%data_type.EQ.dbcsr_type_complex_4.OR.&
01944        matrix%m%data_type.EQ.dbcsr_type_complex_8) THEN
01945        is_complex = .TRUE.
01946     ELSE
01947        is_complex = .FALSE.
01948     ENDIF
01949   END FUNCTION dbcsr_is_complex_obj
01950 
01951 
01952 ! *****************************************************************************
01960   PURE FUNCTION dbcsr_get_data_memory_type (matrix) &
01961        RESULT (memory_type)
01962     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01963     INTEGER                                  :: memory_type
01964 
01965     memory_type = matrix%m%data_memory_type
01966   END FUNCTION dbcsr_get_data_memory_type
01967 
01968 ! *****************************************************************************
01973   PURE FUNCTION dbcsr_get_index_memory_type (matrix) RESULT (memory_type)
01974     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01975     INTEGER                                  :: memory_type
01976 
01977     memory_type = matrix%m%index_memory_type
01978   END FUNCTION dbcsr_get_index_memory_type
01979 
01980 
01981 ! *****************************************************************************
01987   PURE FUNCTION uses_special_memory_matrix (matrix) RESULT (uses_special)
01988     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
01989     LOGICAL                                  :: uses_special
01990 
01991     uses_special = matrix%m%data_memory_type .NE. dbcsr_memory_default
01992   END FUNCTION uses_special_memory_matrix
01993 
01994 
01995 ! *****************************************************************************
02001   PURE FUNCTION uses_special_memory_area (area) RESULT (uses_special)
02002     TYPE(dbcsr_data_obj), INTENT(IN)         :: area
02003     LOGICAL                                  :: uses_special
02004 
02005     IF (ASSOCIATED (area%d)) THEN
02006        uses_special = area%d%memory_type .NE. dbcsr_memory_default
02007     ELSE
02008        uses_special = .FALSE.
02009     ENDIF
02010   END FUNCTION uses_special_memory_area
02011 
02012 
02013 ! *****************************************************************************
02018   FUNCTION get_data_size_area (area) RESULT (data_size)
02019     TYPE(dbcsr_data_obj), INTENT(IN)         :: area
02020     INTEGER                                  :: data_size
02021 
02022     CHARACTER(len=*), PARAMETER :: routineN = 'get_data_size_area', 
02023       routineP = moduleN//':'//routineN
02024 
02025     data_size = dbcsr_data_get_size (area)
02026   END FUNCTION get_data_size_area
02027 
02028 ! *****************************************************************************
02033   PURE FUNCTION data_data_get_type (area) RESULT (data_type)
02034     TYPE(dbcsr_data_obj), INTENT(IN)         :: area
02035     INTEGER                                  :: data_type
02036 
02037     data_type = dbcsr_data_get_type (area)
02038   END FUNCTION data_data_get_type
02039 
02040 
02041 ! *****************************************************************************
02046   FUNCTION get_data_size_matrix (matrix) RESULT (data_size)
02047     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
02048     INTEGER                                  :: data_size
02049 
02050     CHARACTER(len=*), PARAMETER :: routineN = 'get_data_size_matrix', 
02051       routineP = moduleN//':'//routineN
02052 
02053     TYPE(dbcsr_error_type)                   :: error
02054 
02055     data_size = 0
02056     IF (ASSOCIATED (matrix%m%data_area%d)) THEN
02057        SELECT CASE (matrix%m%data_area%d%data_type)
02058           CASE (dbcsr_type_real_8)
02059              IF (ASSOCIATED (matrix%m%data_area%d%r_dp))&
02060                   data_size = SIZE (matrix%m%data_area%d%r_dp)
02061           CASE (dbcsr_type_real_4)
02062              IF (ASSOCIATED (matrix%m%data_area%d%r_sp))&
02063                   data_size = SIZE (matrix%m%data_area%d%r_sp)
02064           CASE (dbcsr_type_complex_8)
02065              IF (ASSOCIATED (matrix%m%data_area%d%c_dp))&
02066                   data_size = SIZE (matrix%m%data_area%d%c_dp)
02067           CASE (dbcsr_type_complex_4)
02068              IF (ASSOCIATED (matrix%m%data_area%d%c_sp))&
02069                   data_size = SIZE (matrix%m%data_area%d%c_sp)
02070           CASE default
02071              CALL dbcsr_assert (.FALSE., dbcsr_failure_level, dbcsr_caller_error,&
02072                   routineN, "Incorrect data type",__LINE__,error)
02073           END SELECT
02074     ELSE
02075        CALL dbcsr_assert (.FALSE., dbcsr_warning_level, dbcsr_caller_error, routineN,&
02076             "Uninitialized data area",__LINE__,error)
02077        data_size = 0
02078     ENDIF
02079   END FUNCTION get_data_size_matrix
02080 
02081 ! *****************************************************************************
02086   PURE FUNCTION get_data_size_ref_matrix (matrix) RESULT (data_size_referenced)
02087     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
02088     INTEGER                                  :: data_size_referenced
02089 
02090     data_size_referenced = dbcsr_data_get_size_referenced (matrix%m%data_area)
02091   END FUNCTION get_data_size_ref_matrix
02092 
02093 ! *****************************************************************************
02098   PURE FUNCTION get_data_size_ref_area (area) RESULT (data_size_referenced)
02099     TYPE(dbcsr_data_obj), INTENT(IN)         :: area
02100     INTEGER                                  :: data_size_referenced
02101 
02102     data_size_referenced = dbcsr_data_get_size_referenced (area)
02103   END FUNCTION get_data_size_ref_area
02104 
02105 ! *****************************************************************************
02110   PURE SUBROUTINE set_data_size_ref_area (data_area, referenced_size)
02111     TYPE(dbcsr_data_obj), INTENT(INOUT)      :: data_area
02112     INTEGER, INTENT(IN)                      :: referenced_size
02113 
02114     CALL dbcsr_data_set_size_referenced (data_area, referenced_size)
02115   END SUBROUTINE set_data_size_ref_area
02116 
02117 
02118 ! *****************************************************************************
02124   FUNCTION dbcsr_get_data_size_used (matrix, error) RESULT (data_size)
02125     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
02126     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
02127     INTEGER                                  :: data_size
02128 
02129     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_get_data_size_used', 
02130       routineP = moduleN//':'//routineN
02131 
02132     INTEGER                                  :: blk, blk_p, col, 
02133                                                 error_handle, nze, row
02134     INTEGER, DIMENSION(:), POINTER           :: col_blk_sizes, row_blk_sizes
02135 
02136 !type(dbcsr_iterator_type) :: iter
02137 !   ---------------------------------------------------------------------------
02138 
02139     CALL dbcsr_error_set (routineN, error_handle, error)
02140     row_blk_sizes => array_data (dbcsr_row_block_sizes (matrix))
02141     col_blk_sizes => array_data (dbcsr_col_block_sizes (matrix))
02142     data_size = 0
02143     !$OMP DO
02144     DO row = 1, matrix%m%nblkrows_total
02145        DO blk = matrix%m%row_p(row)+1, matrix%m%row_p(row+1)
02146           col = matrix%m%col_i(blk)
02147           blk_p = ABS (matrix%m%blk_p(blk))
02148           IF (matrix%m%blk_p(blk) .NE. 0) THEN
02149              nze = row_blk_sizes(row) * col_blk_sizes(col)
02150              data_size = data_size + nze
02151           ENDIF
02152        ENDDO
02153     ENDDO
02154     !$OMP END DO
02155     CALL dbcsr_error_stop (error_handle, error)
02156   END FUNCTION dbcsr_get_data_size_used
02157 
02158 
02159 
02160 ! *****************************************************************************
02165   PURE FUNCTION dbcsr_get_num_blocks (matrix) RESULT (num_blocks)
02166     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
02167     INTEGER                                  :: num_blocks
02168 
02169     num_blocks = matrix%m%nblks
02170   END FUNCTION dbcsr_get_num_blocks
02171 
02172 ! *****************************************************************************
02177   PURE FUNCTION dbcsr_get_nze (matrix) RESULT (num_nze)
02178     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
02179     INTEGER                                  :: num_nze
02180 
02181     num_nze = matrix%m%nze
02182   END FUNCTION dbcsr_get_nze
02183 
02184 ! *****************************************************************************
02188   FUNCTION dbcsr_get_occupation (matrix) RESULT (occupation)
02189     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
02190     REAL(KIND=real_8)                        :: occupation
02191 
02192     INTEGER                                  :: nfullcols, nfullrows
02193     INTEGER, DIMENSION(:), POINTER           :: row_blk_size
02194 
02195     occupation = REAL(matrix%m%nze,real_8)
02196     CALL mp_sum(occupation,dbcsr_mp_group(dbcsr_distribution_mp(matrix%m%dist)))
02197 
02198     nfullrows = dbcsr_nfullrows_total(matrix)
02199     nfullcols = dbcsr_nfullcols_total(matrix)
02200 
02201     row_blk_size => array_data (matrix%m%row_blk_size)
02202 
02203     IF(nfullrows.NE.0.AND.nfullcols.NE.0) THEN
02204        IF(dbcsr_has_symmetry (matrix)) THEN
02205           occupation = 2.0_real_8 * occupation / &
02206                ( REAL(nfullrows,real_8) * REAL(nfullrows+1,real_8) + 
02207                SUM( REAL(row_blk_size,real_8) * REAL(row_blk_size-1,real_8) ) )
02208        ELSE
02209           occupation = occupation / ( REAL(nfullrows,real_8) * REAL(nfullcols,real_8) )
02210        ENDIF
02211     ELSE
02212        occupation = 0.0_real_8
02213     ENDIF
02214   END FUNCTION dbcsr_get_occupation
02215 
02216 
02217 ! *****************************************************************************
02218 ! Arrays
02219 ! *****************************************************************************
02220 
02221 
02222 ! *****************************************************************************
02227   SUBROUTINE array_init_1d (array, n)
02228     TYPE(dbcsr_1d_array_obj), INTENT(OUT)    :: array
02229     INTEGER, INTENT(IN), OPTIONAL            :: n
02230 
02231     NULLIFY (array%s, array%refcount)
02232     IF (PRESENT (n)) CALL dbcsr_array_new (array, n)
02233   END SUBROUTINE array_init_1d
02234 
02235 ! *****************************************************************************
02240   SUBROUTINE array_init_2d (array, n)
02241     TYPE(dbcsr_2d_array_obj), INTENT(OUT)    :: array
02242     INTEGER, DIMENSION(2), INTENT(IN), 
02243       OPTIONAL                               :: n
02244 
02245     NULLIFY (array%s, array%refcount)
02246     IF (PRESENT (n)) CALL dbcsr_array_new (array, n)
02247   END SUBROUTINE array_init_2d
02248 
02249 ! *****************************************************************************
02254   SUBROUTINE array_new_1d (array, n)
02255     TYPE(dbcsr_1d_array_obj), INTENT(INOUT)  :: array
02256     INTEGER, INTENT(IN)                      :: n
02257 
02258     CHARACTER(len=*), PARAMETER :: routineN = 'array_new_1d', 
02259       routineP = moduleN//':'//routineN
02260 
02261     TYPE(dbcsr_error_type)                   :: error
02262 
02263 !   ---------------------------------------------------------------------------
02264 
02265     CALL dbcsr_assert (.NOT. ASSOCIATED (array%s), dbcsr_warning_level,&
02266          dbcsr_caller_error, routineN,&
02267          "Array was already allocated, memory leakage will occur",__LINE__,error)
02268     ALLOCATE (array%s(n))
02269     ALLOCATE (array%refcount)
02270     array%refcount = 1
02271   END SUBROUTINE array_new_1d
02272 
02273 ! *****************************************************************************
02278   SUBROUTINE array_new_2d (array, n)
02279     TYPE(dbcsr_2d_array_obj), INTENT(INOUT)  :: array
02280     INTEGER, DIMENSION(2), INTENT(IN)        :: n
02281 
02282     CHARACTER(len=*), PARAMETER :: routineN = 'array_new_2d', 
02283       routineP = moduleN//':'//routineN
02284 
02285     TYPE(dbcsr_error_type)                   :: error
02286 
02287 !   ---------------------------------------------------------------------------
02288 
02289     CALL dbcsr_assert (.NOT. ASSOCIATED (array%s), dbcsr_warning_level,&
02290          dbcsr_caller_error, routineN,&
02291          "Array was already allocated, memory leakage will occur",__LINE__,error)
02292     ALLOCATE (array%s(n(1), n(2)))
02293     ALLOCATE (array%refcount)
02294     array%refcount = 1
02295   END SUBROUTINE array_new_2d
02296 
02297 ! *****************************************************************************
02305   SUBROUTINE array_put_1d (array, position, matrix)
02306     TYPE(dbcsr_1d_array_obj), INTENT(INOUT)  :: array
02307     INTEGER, INTENT(IN)                      :: position
02308     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
02309 
02310     CHARACTER(len=*), PARAMETER :: routineN = 'array_put_1d', 
02311       routineP = moduleN//':'//routineN
02312 
02313     TYPE(dbcsr_error_type)                   :: error
02314 
02315 !   ---------------------------------------------------------------------------
02316 
02317     CALL dbcsr_assert (ASSOCIATED (array%s), dbcsr_fatal_level, dbcsr_caller_error,&
02318          routineN, "Can not insert into unallocated array",__LINE__,error)
02319     CALL dbcsr_assert(position .GE. LBOUND (array%s, 1)&
02320          .AND. position .LE. UBOUND (array%s, 1),&
02321          dbcsr_failure_level, dbcsr_caller_error, routineN,&
02322          "Insertion position out of bounds.",__LINE__,error)
02323     array%s(position) = matrix
02324     CALL dbcsr_hold (array%s(position))
02325   END SUBROUTINE array_put_1d
02326 
02327 ! *****************************************************************************
02335   SUBROUTINE array_put_2d (array, position, matrix)
02336     TYPE(dbcsr_2d_array_obj), INTENT(INOUT)  :: array
02337     INTEGER, DIMENSION(2), INTENT(IN)        :: position
02338     TYPE(dbcsr_obj), INTENT(IN)              :: matrix
02339 
02340     CHARACTER(len=*), PARAMETER :: routineN = 'array_put_2d', 
02341       routineP = moduleN//':'//routineN
02342 
02343     TYPE(dbcsr_error_type)                   :: error
02344 
02345 !   ---------------------------------------------------------------------------
02346 
02347     CALL dbcsr_assert (ASSOCIATED (array%s), dbcsr_fatal_level, dbcsr_caller_error,&
02348          routineN, "Can not insert into unallocated array",__LINE__,error)
02349     CALL dbcsr_assert(position(1) .GE. LBOUND (array%s, 1)&
02350          .AND. position(1) .LE. UBOUND (array%s, 1)&
02351          .AND. position(2) .GE. LBOUND (array%s, 2)&
02352          .AND. position(2) .LE. UBOUND (array%s, 2),&
02353          dbcsr_failure_level, dbcsr_caller_error, routineN,&
02354          "Insertion position out of bounds.",__LINE__,error)
02355     array%s(position(1), position(2)) = matrix
02356     CALL dbcsr_hold (array%s(position(1), position(2)))
02357   END SUBROUTINE array_put_2d
02358 
02359 ! *****************************************************************************
02365   FUNCTION array_get_1d (array, position) RESULT (matrix)
02366     TYPE(dbcsr_1d_array_obj), INTENT(IN)     :: array
02367     INTEGER, INTENT(IN)                      :: position
02368     TYPE(dbcsr_obj)                          :: matrix
02369 
02370     CHARACTER(len=*), PARAMETER :: routineN = 'array_get_1d', 
02371       routineP = moduleN//':'//routineN
02372 
02373     TYPE(dbcsr_error_type)                   :: error
02374 
02375 !   ---------------------------------------------------------------------------
02376 
02377     CALL dbcsr_assert (ASSOCIATED (array%s), dbcsr_fatal_level, dbcsr_caller_error,&
02378          routineN, "Can not read from an unallocated array",__LINE__,error)
02379     CALL dbcsr_assert (position .GE. LBOUND (array%s,1)&
02380          .AND. position .LE. UBOUND (array%s,1),&
02381          dbcsr_failure_level, dbcsr_caller_error, routineN,&
02382          "Read position out of bounds.",__LINE__,error)
02383     matrix = array%s(position)
02384   END FUNCTION array_get_1d
02385 
02386 ! *****************************************************************************
02392   FUNCTION array_get_2d (array, position) RESULT (matrix)
02393     TYPE(dbcsr_2d_array_obj), INTENT(IN)     :: array
02394     INTEGER, DIMENSION(2), INTENT(IN)        :: position
02395     TYPE(dbcsr_obj)                          :: matrix
02396 
02397     CHARACTER(len=*), PARAMETER :: routineN = 'array_get_2d', 
02398       routineP = moduleN//':'//routineN
02399 
02400     TYPE(dbcsr_error_type)                   :: error
02401 
02402 !   ---------------------------------------------------------------------------
02403 
02404     CALL dbcsr_assert (ASSOCIATED (array%s), dbcsr_fatal_level, dbcsr_caller_error,&
02405          routineN, "Can not read from an unallocated array",__LINE__,error)
02406     CALL dbcsr_assert(position(1) .GE. LBOUND (array%s, 1)&
02407          .AND. position(1) .LE. UBOUND (array%s, 1)&
02408          .AND. position(2) .GE. LBOUND (array%s, 2)&
02409          .AND. position(2) .LE. UBOUND (array%s, 2),&
02410          dbcsr_failure_level, dbcsr_caller_error, routineN,&
02411          "Read position out of bounds.",__LINE__,error)
02412     matrix = array%s(position(1), position(2))
02413   END FUNCTION array_get_2d
02414 
02415 ! *****************************************************************************
02419   SUBROUTINE array_destroy_1d (array)
02420     TYPE(dbcsr_1d_array_obj), INTENT(INOUT)  :: array
02421 
02422     CHARACTER(len=*), PARAMETER :: routineN = 'array_destroy_1d', 
02423       routineP = moduleN//':'//routineN
02424 
02425     INTEGER                                  :: i
02426     TYPE(dbcsr_error_type)                   :: error
02427 
02428 !   ---------------------------------------------------------------------------
02429 
02430     CALL dbcsr_assert (ASSOCIATED (array%s), dbcsr_fatal_level, dbcsr_caller_error,&
02431          routineN, "Can not destroy an unallocated array",__LINE__,error)
02432     CALL dbcsr_assert (ASSOCIATED (array%refcount),&
02433          dbcsr_failure_level, dbcsr_caller_error, routineN,&
02434          "Reference count does not exist.",__LINE__,error)
02435     CALL dbcsr_assert (array%refcount .EQ. 0, dbcsr_warning_level, dbcsr_caller_error,&
02436          routineN, "References are still held",__LINE__,error)
02437     DO i = LBOUND (array%s, 1), UBOUND (array%s, 1)
02438        CALL dbcsr_release (array%s(i))
02439     ENDDO
02440     DEALLOCATE (array%s)
02441     DEALLOCATE (array%refcount)
02442     CALL dbcsr_allocate_matrix_array (array)
02443   END SUBROUTINE array_destroy_1d
02444 
02445 ! *****************************************************************************
02449   SUBROUTINE array_destroy_2d (array)
02450     TYPE(dbcsr_2d_array_obj), INTENT(INOUT)  :: array
02451 
02452     CHARACTER(len=*), PARAMETER :: routineN = 'array_destroy_2d', 
02453       routineP = moduleN//':'//routineN
02454 
02455     INTEGER                                  :: i, j
02456     TYPE(dbcsr_error_type)                   :: error
02457 
02458 !   ---------------------------------------------------------------------------
02459 
02460     CALL dbcsr_assert (ASSOCIATED (array%s), dbcsr_fatal_level, dbcsr_caller_error,&
02461          routineN, "Can not destroy an unallocated array",__LINE__,error)
02462     CALL dbcsr_assert (ASSOCIATED (array%refcount),&
02463          dbcsr_failure_level, dbcsr_caller_error, routineN,&
02464          "Reference count does not exist.",__LINE__,error)
02465     CALL dbcsr_assert (array%refcount .EQ. 0, dbcsr_warning_level, dbcsr_caller_error,&
02466          routineN, "References are still held",__LINE__,error)
02467     DO i = LBOUND (array%s, 1), UBOUND (array%s, 1)
02468        DO j = LBOUND (array%s, 2), UBOUND (array%s, 2)
02469           CALL dbcsr_release (array%s(i, j))
02470        ENDDO
02471     ENDDO
02472     DEALLOCATE (array%s)
02473     DEALLOCATE (array%refcount)
02474     CALL dbcsr_allocate_matrix_array (array)
02475   END SUBROUTINE array_destroy_2d
02476 
02477 
02478 ! *****************************************************************************
02482   SUBROUTINE array_hold_1d (array)
02483     TYPE(dbcsr_1d_array_obj), INTENT(INOUT)  :: array
02484 
02485     CHARACTER(len=*), PARAMETER :: routineN = 'array_hold_1d', 
02486       routineP = moduleN//':'//routineN
02487 
02488     TYPE(dbcsr_error_type)                   :: error
02489 
02490 !   ---------------------------------------------------------------------------
02491 
02492     CALL dbcsr_assert (ASSOCIATED (array%refcount),&
02493          dbcsr_failure_level, dbcsr_caller_error, routineN,&
02494          "Reference count does not exist.",__LINE__,error)
02495     array%refcount = array%refcount + 1
02496   END SUBROUTINE array_hold_1d
02497 
02498 ! *****************************************************************************
02502   SUBROUTINE array_hold_2d (array)
02503     TYPE(dbcsr_2d_array_obj), INTENT(INOUT)  :: array
02504 
02505     CHARACTER(len=*), PARAMETER :: routineN = 'array_hold_2d', 
02506       routineP = moduleN//':'//routineN
02507 
02508     TYPE(dbcsr_error_type)                   :: error
02509 
02510 !   ---------------------------------------------------------------------------
02511 
02512     CALL dbcsr_assert (ASSOCIATED (array%refcount),&
02513          dbcsr_failure_level, dbcsr_caller_error, routineN,&
02514          "Reference count does not exist.",__LINE__,error)
02515     array%refcount = array%refcount + 1
02516   END SUBROUTINE array_hold_2d
02517 
02518 
02519 ! *****************************************************************************
02525   SUBROUTINE array_release_1d (array)
02526     TYPE(dbcsr_1d_array_obj), INTENT(INOUT)  :: array
02527 
02528     CHARACTER(len=*), PARAMETER :: routineN = 'array_release_1d', 
02529       routineP = moduleN//':'//routineN
02530 
02531     TYPE(dbcsr_error_type)                   :: error
02532 
02533 !   ---------------------------------------------------------------------------
02534 
02535     CALL dbcsr_assert (ASSOCIATED (array%refcount),&
02536          dbcsr_failure_level, dbcsr_caller_error, routineN,&
02537          "Reference count does not exist.",__LINE__,error)
02538     array%refcount = array%refcount - 1
02539     IF (array%refcount .EQ. 0) THEN
02540        CALL dbcsr_array_destroy (array)
02541     ENDIF
02542   END SUBROUTINE array_release_1d
02543 
02544 ! *****************************************************************************
02550   SUBROUTINE array_release_2d (array)
02551     TYPE(dbcsr_2d_array_obj), INTENT(INOUT)  :: array
02552 
02553     CHARACTER(len=*), PARAMETER :: routineN = 'array_release_2d', 
02554       routineP = moduleN//':'//routineN
02555 
02556     TYPE(dbcsr_error_type)                   :: error
02557 
02558 !   ---------------------------------------------------------------------------
02559 
02560     CALL dbcsr_assert (ASSOCIATED (array%refcount),&
02561          dbcsr_failure_level, dbcsr_caller_error, routineN,&
02562          "Reference count does not exist.",__LINE__,error)
02563     array%refcount = array%refcount - 1
02564     IF (array%refcount .EQ. 0) THEN
02565        CALL dbcsr_array_destroy (array)
02566     ENDIF
02567   END SUBROUTINE array_release_2d
02568 
02569 ! *****************************************************************************
02574   SUBROUTINE dbcsr_destroy_1d_array(marray, error)
02575     TYPE(dbcsr_array_type), INTENT(INOUT)    :: marray
02576     TYPE(dbcsr_error_type), INTENT(inout)    :: error
02577 
02578     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_destroy_1d_array', 
02579       routineP = moduleN//':'//routineN
02580 
02581     INTEGER                                  :: i
02582 
02583 !   ---------------------------------------------------------------------------
02584 
02585     DO i = LBOUND(marray%mats,1), UBOUND (marray%mats,1)
02586        CALL dbcsr_destroy (marray%mats(i), error=error, force=.TRUE.)
02587     ENDDO
02588     DEALLOCATE (marray%mats)
02589     !CALL dbcsr_destroy_image_dist(marray%image_dist)
02590   END SUBROUTINE dbcsr_destroy_1d_array
02591 
02592 
02593 ! *****************************************************************************
02598   SUBROUTINE dbcsr_destroy_2d_array(marray,error)
02599     TYPE(dbcsr_2d_array_type), INTENT(INOUT) :: marray
02600     TYPE(dbcsr_error_type), INTENT(inout)    :: error
02601 
02602     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_destroy_2d_array', 
02603       routineP = moduleN//':'//routineN
02604 
02605     INTEGER                                  :: col, row
02606 
02607 !   ---------------------------------------------------------------------------
02608 
02609     DO row = LBOUND(marray%mats,1), UBOUND (marray%mats,1)
02610        DO col = LBOUND(marray%mats,2), UBOUND (marray%mats,2)
02611           CALL dbcsr_destroy (marray%mats(row, col), error=error, force=.TRUE.)
02612        ENDDO
02613     ENDDO
02614     CALL dbcsr_image_dist_release (marray%image_dist, error)
02615     DEALLOCATE (marray%mats)
02616     NULLIFY (marray%mats)
02617   END SUBROUTINE dbcsr_destroy_2d_array
02618 
02619 ! *****************************************************************************
02623   SUBROUTINE dbcsr_image_dist_release (imgdist, error)
02624     TYPE(dbcsr_imagedistribution_obj), 
02625       INTENT(INOUT)                          :: imgdist
02626     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
02627 
02628     IF (ASSOCIATED (imgdist%i)) THEN
02629        imgdist%i%refcount = imgdist%i%refcount - 1
02630        IF (imgdist%i%refcount .EQ. 0) THEN
02631           CALL dbcsr_destroy_image_dist (imgdist%i)
02632           DEALLOCATE (imgdist%i)
02633           NULLIFY (imgdist%i)
02634        ENDIF
02635     ENDIF
02636   END SUBROUTINE dbcsr_image_dist_release
02637 
02638 ! *****************************************************************************
02641   SUBROUTINE dbcsr_image_dist_hold (imgdist, error)
02642     TYPE(dbcsr_imagedistribution_obj), 
02643       INTENT(INOUT)                          :: imgdist
02644     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
02645 
02646     imgdist%i%refcount = imgdist%i%refcount + 1
02647   END SUBROUTINE dbcsr_image_dist_hold
02648 
02649 ! *****************************************************************************
02653   SUBROUTINE dbcsr_image_dist_init (imgdist, error)
02654     TYPE(dbcsr_imagedistribution_obj), 
02655       INTENT(OUT)                            :: imgdist
02656     TYPE(dbcsr_error_type), INTENT(INOUT)    :: error
02657 
02658     NULLIFY (imgdist%i)
02659   END SUBROUTINE dbcsr_image_dist_init
02660 
02661 ! *****************************************************************************
02666   SUBROUTINE dbcsr_destroy_image_dist(imgdist)
02667     TYPE(dbcsr_imagedistribution_type), 
02668       INTENT(INOUT)                          :: imgdist
02669 
02670     INTEGER                                  :: i
02671 
02672 !   ---------------------------------------------------------------------------
02673 
02674     CALL array_release (imgdist%row_image)
02675     CALL array_release (imgdist%col_image)
02676     CALL dbcsr_distribution_release (imgdist%main)
02677     !
02678     CALL array_release (imgdist%vrow_dist)
02679     CALL array_release (imgdist%vcol_dist)
02680     !
02681     IF (imgdist%has_other_vl_rows) THEN
02682        DO i = LBOUND (imgdist%other_vl_rows,1), UBOUND (imgdist%other_vl_rows,1)
02683           CALL array_release (imgdist%other_vl_rows(i))
02684        END DO
02685        DEALLOCATE (imgdist%other_vl_rows)
02686        NULLIFY (imgdist%other_vl_rows)
02687        imgdist%has_other_vl_rows = .FALSE.
02688     ENDIF
02689     !
02690     IF (imgdist%has_other_vl_cols) THEN
02691        DO i = LBOUND (imgdist%other_vl_cols,1), UBOUND (imgdist%other_vl_cols,1)
02692           CALL array_release (imgdist%other_vl_cols(i))
02693        END DO
02694        DEALLOCATE (imgdist%other_vl_cols)
02695        NULLIFY (imgdist%other_vl_cols)
02696        imgdist%has_other_vl_cols = .FALSE.
02697     ENDIF
02698     !
02699     IF (imgdist%has_global_vrow_map) THEN
02700        CALL array_release (imgdist%global_vrow_map)
02701     ENDIF
02702     IF (imgdist%has_global_vcol_map) THEN
02703        CALL array_release (imgdist%global_vcol_map)
02704     ENDIF
02705   END SUBROUTINE dbcsr_destroy_image_dist
02706 
02707 ! *****************************************************************************
02708 ! Mutable data
02709 ! *****************************************************************************
02710 
02711 ! *****************************************************************************
02715   SUBROUTINE dbcsr_mutable_init (mutable)
02716     TYPE(dbcsr_mutable_obj), INTENT(OUT)     :: mutable
02717 
02718     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_mutable_init', 
02719       routineP = moduleN//':'//routineN
02720 
02721 !   ---------------------------------------------------------------------------
02722 
02723     NULLIFY (mutable%m)
02724   END SUBROUTINE dbcsr_mutable_init
02725 
02726 ! *****************************************************************************
02730   SUBROUTINE dbcsr_mutable_destroy (mutable)
02731     TYPE(dbcsr_mutable_obj), INTENT(INOUT)   :: mutable
02732 
02733     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_mutable_destroy', 
02734       routineP = moduleN//':'//routineN
02735 
02736 !   ---------------------------------------------------------------------------
02737 
02738     IF (ASSOCIATED (mutable%m)) THEN
02739        !CALL dbcsr_assert (mutable%m%refcount .EQ. 0, dbcsr_warning_level,&
02740        !     dbcsr_caller_error, routineN, "Destroying with non-0 reference count")
02741        CALL btree_destroy_s (mutable%m%btree_s)
02742        CALL btree_destroy_d (mutable%m%btree_d)
02743        CALL btree_destroy_c (mutable%m%btree_c)
02744        CALL btree_destroy_z (mutable%m%btree_z)
02745        DEALLOCATE (mutable%m)
02746     ENDIF
02747     NULLIFY (mutable%m)
02748   END SUBROUTINE dbcsr_mutable_destroy
02749 
02750 
02751 ! *****************************************************************************
02755   SUBROUTINE dbcsr_mutable_hold (mutable)
02756     TYPE(dbcsr_mutable_obj), INTENT(INOUT)   :: mutable
02757 
02758     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_mutable_hold', 
02759       routineP = moduleN//':'//routineN
02760 
02761     TYPE(dbcsr_error_type)                   :: error
02762 
02763 !   ---------------------------------------------------------------------------
02764 
02765     CALL dbcsr_assert (ASSOCIATED (mutable%m), dbcsr_fatal_level,&
02766          dbcsr_caller_error, routineN, "Mutable data area not instantiated",__LINE__,error)
02767     mutable%m%refcount = mutable%m%refcount + 1
02768   END SUBROUTINE dbcsr_mutable_hold
02769 
02770 ! *****************************************************************************
02776   SUBROUTINE dbcsr_mutable_release (mutable)
02777     TYPE(dbcsr_mutable_obj), INTENT(INOUT)   :: mutable
02778 
02779     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_mutable_release', 
02780       routineP = moduleN//':'//routineN
02781 
02782     TYPE(dbcsr_error_type)                   :: error
02783 
02784 !   ---------------------------------------------------------------------------
02785 
02786     CALL dbcsr_assert (ASSOCIATED (mutable%m), dbcsr_fatal_level,&
02787          dbcsr_caller_error, routineN, "Mutable data area not instantiated",__LINE__,error)
02788     mutable%m%refcount = mutable%m%refcount - 1
02789     IF (mutable%m%refcount .EQ. 0) THEN
02790        CALL dbcsr_mutable_destroy (mutable)
02791     ENDIF
02792   END SUBROUTINE dbcsr_mutable_release
02793 
02794 ! *****************************************************************************
02801   SUBROUTINE dbcsr_mutable_new (mutable, data_type)
02802     TYPE(dbcsr_mutable_obj), INTENT(INOUT)   :: mutable
02803     INTEGER, INTENT(IN)                      :: data_type
02804 
02805     CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_mutable_new', 
02806       routineP = moduleN//':'//routineN
02807 
02808     TYPE(dbcsr_error_type)                   :: error
02809 
02810 !   ---------------------------------------------------------------------------
02811 
02812     CALL dbcsr_assert (.NOT.ASSOCIATED (mutable%m), dbcsr_fatal_level,&
02813          dbcsr_caller_error, routineN, "Mutable data area already instantiated",__LINE__,error)
02814     CALL dbcsr_assert (data_type.EQ.dbcsr_type_real_4&
02815          .OR. data_type.EQ.dbcsr_type_real_8&
02816          .OR. data_type.EQ.dbcsr_type_complex_4&
02817          .OR. data_type.EQ.dbcsr_type_complex_8, dbcsr_fatal_level,&
02818          dbcsr_wrong_args_error, routineN, "Invalid data type",__LINE__,error)
02819     ALLOCATE (mutable%m)
02820     mutable%m%refcount = 1
02821     mutable%m%data_type = data_type
02822     CALL btree_new_s (mutable%m%btree_s)
02823     CALL btree_new_d (mutable%m%btree_d)
02824     CALL btree_new_c (mutable%m%btree_c)
02825     CALL btree_new_z (mutable%m%btree_z)
02826   END SUBROUTINE dbcsr_mutable_new
02827 
02828 ! *****************************************************************************
02835   PURE FUNCTION dbcsr_mutable_instantiated (mutable) RESULT (instantiated)
02836     TYPE(dbcsr_mutable_obj), INTENT(IN)      :: mutable
02837     LOGICAL                                  :: instantiated
02838 
02839 !   ---------------------------------------------------------------------------
02840 
02841     instantiated = ASSOCIATED (mutable%m)
02842   END FUNCTION dbcsr_mutable_instantiated
02843 
02844 
02845 !! *****************************************************************************
02846 !!> \brief Returns the entire data for a matrix.
02847 !!> \par Warning
02848 !!>      This routine should only be used within DBCSR code.
02849 !!> \param[in] area       data area
02850 !!> \param[in] coersion   force datatype
02851 !!> \param[out] data      pointer to data
02852 !! *****************************************************************************
02853 !  SUBROUTINE get_data_m_[nametype1] (matrix, DATA)
02854 !    TYPE(dbcsr_obj), INTENT(IN)              :: matrix
02855 !    [type1], DIMENSION(:), POINTER :: DATA
02856 !
02857 !    CALL get_data_[nametype1] (matrix%m%data_area, DATA)
02858 !  END SUBROUTINE get_data_m_[nametype1]
02859 
02860 
02861 END MODULE dbcsr_methods