224 subroutine initialize_struct_slepc (this, matrix_size, matrix_type, block_size)
226 class(slepcmatrix) :: this
227 integer,
intent(in) :: matrix_size, matrix_type, block_size
229 petscint :: non_zero_guess, petsc_row, to_insert, one = 1, dummy_int, matrix_size_petsc
230 petscerrorcode :: ierr
231 real(wp) :: mem_usage_begin, mem_usage_end, ratio, dummy_real
232 integer :: ido, total_vals, per_elm, l_update, c_update, indicies, total_elems
234 integer,
allocatable :: mem_use_total(:)
235 petscint,
allocatable :: test_insertion_ij(:)
236 petscscalar,
allocatable :: test_insertion_val(:)
241 call this % slepc_cache % construct
242 call this % destroy_PETSC_mat
244 this % matrix_dimension_PETSC = this % matrix_dimension
245 write (stdout,
"('Mat_size = ',i8)") this % matrix_dimension_PETSC
247 call this % create_PETSC_mat(this % matrix_dimension_PETSC, matrix_type)
249 matrix_size_petsc = matrix_size
250 this % local_size_PETSC = petsc_decide
253 call petscsplitownership(int(grid % gcomm, kind(petsc_comm_world)), &
254 this % local_size_PETSC, this % matrix_dimension_PETSC, ierr)
257 this % local_size = this % local_size_PETSC
259 write (stdout,
"('Local_size = ',i8)") this % local_size
260 write (stdout,
"('THRESHOLD = ',es16.8)") this % threshold
263 call mpi_mod_scan_sum(this % local_size, this % end_row, grid % gcomm)
265 this % start_row = this % end_row - this % local_size
268 this % start_row_PETSC = this % start_row
269 this % end_row_PETSC = this % end_row
270 write (stdout,
"('start,end = ',2i8)") this % start_row, this % end_row
273 if (this % matrix_type ==
mat_dense)
then
274 call petscmemorygetcurrentusage(mem_usage_begin, ierr)
275 call matcreatedense(grid % gcomm, petsc_decide, petsc_decide, matrix_size_petsc, matrix_size_petsc, &
276 petsc_null_scalar_array, this % A, ierr)
277 call petscmemorygetcurrentusage(mem_usage_end, ierr)
279 write (stdout,
"('Matrix uses ',f8.3,' GB of memory')") (mem_usage_end - mem_usage_begin) / 1024**3
280 this % mem_track = mem_usage_end - mem_usage_begin
282 call master_memory % track_memory(this % mem_track, 1, 0,
'PETSC MATRIX')
287 if (
allocated(this % diagonal_nonzero))
deallocate(this % diagonal_nonzero)
288 if (
allocated(this % offdiagonal_nonzero))
deallocate(this % offdiagonal_nonzero)
291 allocate(this % diagonal_nonzero(this % local_size), stat = ierr)
292 allocate(this % offdiagonal_nonzero(this % local_size), stat = ierr)
294 this % diagonal_nonzero = 0
295 this % offdiagonal_nonzero = 0
299 this % diagonal_start = this % start_row_PETSC
300 this % diagonal_end = min(this % start_row_PETSC + this % local_size, this % matrix_dimension)
302 do ido = 1, this % local_size
303 total_vals = this % matrix_dimension - this % start_row - ido + 1
304 this % diagonal_nonzero(ido) = this % diagonal_end - this % diagonal_start + 1 - ido
305 this % offdiagonal_nonzero(ido) = total_vals - this % diagonal_nonzero(ido)
308 if ((this % start_row + ido - 1) >= block_size)
then
311 ratio = real(this % diagonal_nonzero(ido)) / real(total_vals)
313 total_vals = real(total_vals) * 0.1
317 this % diagonal_nonzero(ido) = max(1, int(real(total_vals) * ratio))
318 this % offdiagonal_nonzero(ido) = int(real(total_vals) * (1.0 - ratio))
320 total_elems = total_elems + this % diagonal_nonzero(ido) + this % offdiagonal_nonzero(ido)
324 call petscmemorygetcurrentusage(mem_usage_begin, ierr)
326 this % continuum_counter = 0
327 this % L2_counter = 0
329 call matcreatesbaij(grid % gcomm, one, petsc_decide, petsc_decide, matrix_size_petsc, matrix_size_petsc, &
330 petsc_default_integer, this % diagonal_nonzero, petsc_default_integer, this % offdiagonal_nonzero, &
348 call mpi_xermsg(
'SLEPCMatrix_module',
'initialize_struct_SLEPC',
'Failed to create matrix', int(ierr), 1)
352 call matsetoption(this % A, mat_ignore_off_proc_entries, petsc_true, ierr)
353 call matsetoption(this % A, mat_no_off_proc_entries, petsc_true, ierr)
354 call matsetoption(this % A, mat_ignore_lower_triangular, petsc_true, ierr)
355 call matsetoption(this % A, mat_new_nonzero_allocation_err, petsc_false, ierr)
357 deallocate(this % diagonal_nonzero, this % offdiagonal_nonzero)
359 allocate(test_insertion_ij(this % matrix_dimension), test_insertion_val(this % matrix_dimension))
361 do ido = 1, this % matrix_dimension
362 test_insertion_ij(ido) = ido - 1
363 test_insertion_val(ido) = 1.0_wp
367 call master_timer % start_timer(
"Insertion time")
368 do ido = 1, this % local_size
369 total_vals = this % matrix_dimension - this % start_row - ido + 1
371 dummy_int = this % start_row + ido
372 if ((dummy_int - 1) < block_size)
then
373 to_insert =
size(test_insertion_ij(dummy_int:this % matrix_dimension))
376 call matsetvaluesblocked(this % A, one, [dummy_int - one], to_insert, &
377 test_insertion_ij(dummy_int : this % matrix_dimension), &
378 test_insertion_val(dummy_int : this % matrix_dimension), &
384 call matassemblybegin(this % A, mat_flush_assembly, ierr)
385 call matassemblyend(this % A, mat_flush_assembly, ierr)
388 call matzeroentries(this % A, ierr)
391 deallocate(test_insertion_ij, test_insertion_val)
393 write (stdout,
"('PETSC matrix creation completed!')")
395 call master_timer % stop_timer(
"Insertion time")
397 call petscmemorygetcurrentusage(mem_usage_end, ierr)
399 per_elm = (storage_size(test_insertion_ij) + storage_size(test_insertion_val)) / 8
401 write (stdout,
"('Matrix uses ',f0.0,' B of memory')") mem_usage_end - mem_usage_begin
402 write (stdout,
"('Estimated uses ',i0,' B of memory')") total_elems * per_elm
404 this % mem_track = mem_usage_end - mem_usage_begin
406 allocate(mem_use_total(grid % gprocs))
408 call mpi_mod_allgather(total_elems, mem_use_total, grid % gcomm)
410 total_elems = sum(mem_use_total)
411 total_elems = total_elems / grid % gprocs
412 this % mem_track = total_elems
414 deallocate(mem_use_total)
416 write (stdout,
"('Average uses ',i0,' B of memory')") total_elems*per_elm
418 call master_memory % track_memory(per_elm, total_elems, 0,
'PETSC MATRIX')
449 logical function compress_cache_to_csr_format (this, matrix_cache, csr_matrix)
450 class(slepcmatrix) :: this
451 type(csrformat),
allocatable,
intent(out) :: csr_matrix(:)
452 class(matrixcache) :: matrix_cache
454 integer :: ido, jdo, start_idx, end_idx, num_elems, row_count, row_idx, first_row
455 integer :: current_chunk, current_row, next_row, last_chunk, i, j
458 compress_cache_to_csr_format = .false.
461 if (this % matrix_type ==
mat_dense)
return
464 call matrix_cache % sort
468 do ido = 1, matrix_cache % get_size()
469 call matrix_cache % get_from_cache(ido, i, j, coeff)
470 if (ido == 1) first_row = i
471 if (i /= current_row)
then
472 row_count = row_count + 1
477 call master_timer % start_timer(
"CSR conversion")
479 num_elems = matrix_cache % get_size()
481 if (row_count == 0)
then
482 call master_timer % stop_timer(
"CSR conversion")
486 if (
allocated(csr_matrix))
deallocate(csr_matrix)
489 allocate(csr_matrix(row_count))
491 csr_matrix(:) % num_cols = 0
499 current_row = first_row
502 do while (ido <= num_elems)
506 do jdo = ido, num_elems
507 call matrix_cache % get_from_cache(jdo, i, j, coeff)
510 if (i == current_row)
then
511 if (abs(coeff) > this % threshold) row_count = row_count + 1
518 csr_matrix(row_idx) % row = current_row
519 csr_matrix(row_idx) % num_cols = row_count
522 if (row_count > 0)
then
524 call csr_matrix(row_idx) % construct(row_count)
527 do jdo = 1, row_count
529 call matrix_cache % get_from_cache(jdo + start_idx - 1, i, j, coeff)
531 csr_matrix(row_idx) % col(jdo) = j
532 csr_matrix(row_idx) % coeff(jdo) = coeff
533 current_chunk = matrix_cache % get_chunk_idx(jdo + start_idx - 1)
535 if (last_chunk /= current_chunk)
then
536 call matrix_cache % matrix_arrays(last_chunk) % destroy
537 last_chunk = current_chunk
544 start_idx = end_idx + 1
545 row_idx = row_idx + 1
546 current_row = next_row
549 compress_cache_to_csr_format = .true.
552 call matrix_cache % destroy
553 call matrix_cache % construct
554 call master_timer % stop_timer(
"CSR conversion")