Eötvös Quantum Utilities  v4.9.146
Providing the Horsepowers in the Quantum Realm
getDiagInv.F90
Go to the documentation of this file.
1 !======================================================================
2 ! Copyright (C) 2017 Peter Rakyta, Ph.D.
3 !
4 ! This program is free software: you can redistribute it and/or modify
5 ! it under the terms of the GNU General Public License as published by
6 ! the Free Software Foundation, either version 3 of the License, or
7 ! (at your option) any later version.
8 !
9 ! This program is distributed in the hope that it will be useful,
10 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
11 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 ! GNU General Public License for more details.
13 !
14 ! You should have received a copy of the GNU General Public License
15 ! along with this program. If not, see http://www.gnu.org/licenses/.
16 !
17 !======================================================================
18 ! DESCRIPTION:
22 ! MODULE getPartialInv
24 MODULE getdiaginv
25 
26 #ifndef MAXNRHS
27 #define MAXNRHS 100
28 #endif
29 
30 #ifndef MAXMEM
31 #define MAXMEM 16000000 !In kilobytes
32 #endif
33 
34 !****************************************************************
35  CONTAINS
36 !****************************************************************
37 
38 !!
47 #ifdef MIC
48 !dir$ attributes offload:mic :: dgetPartialInv
49 #endif
50  subroutine dgetdiaginv(sizeA, nonzerosA, a, ia, ja, invA, error)
51  use mkl_pardiso
52  use mkl_service
53 
54  IMPLICIT NONE
55 
56 
57 ! Internal solver memory pointer
58  TYPE(mkl_pardiso_handle), pointer :: PT(:) => null()
59  integer*4, pointer :: idum(:) => null(), pia(:) => null(), pja(:) => null()
60  real(KIND=8), pointer :: ddum(:) => null(), pa(:) => null(), pb(:,:) => null(), px(:,:) => null()
61  INTEGER*4 maxfct, mnum, mtype, phase, nrhs, error, msglvl, i, nnz, numIterations, idxIterations, offset
62  INTEGER*4 iparm(64)
63 
64 
65 ! The matrix is given in the CSR format A = (a, ia, ja)
66  INTEGER*4 sizeA, nonzerosA !matrix A is of size sizeA x sizeA
67  INTEGER*4, TARGET :: ia(sizeA+1), ja(nonzerosA)
68  real*8, TARGET :: a(nonzerosa)
69 
70  INTEGER*4, ALLOCATABLE :: perm(:)
71 
72  real*8 :: inva(sizea)
73 
74  INTEGER*4 idx, jdx
75  integer maxmem
76 
77  DATA nrhs /1/, maxfct /1/, mnum /1/
78  allocate( pt(64), idum(0), ddum(0) )
79 
80  do i = 1, 64
81  pt(i)%DUMMY = 0
82  end do
83  pa => a
84  pia => ia
85  pja => ja
86 ! initialize error flag
87  error = 0
88 
89 #ifdef DEBUG
90  write(0,*)
91  write(0,*) 'The input matrix: (getDiagInv)'
92  idx = 0
93  do jdx = 1, nonzerosa
94  if (jdx >= ia(idx+1)) then
95  idx=idx+1
96  end if
97 
98  if (jdx < ia(idx+1)) then
99  !write(0,*) 'row', idx, 'col', ja(jdx), 'element of a = ', a(jdx)
100  end if
101  end do
102  write(0,*) 'The size of the diagonal vector to be calculated:', sizea
103 #endif
104 
105 !
106 ! .. Setup Pardiso control parameters und initialize the solvers
107 ! internal adress pointers. This is only necessary for the FIRST
108 ! call of the PARDISO solver.
109 !
110 ! choose the matrix type for pardiso factorization
111  mtype = 1 ! structurally symmetric real matrix
112 
113 ! .. PARDISO license check and initialize solver
114  call pardisoinit(pt, mtype, iparm)
115 
116 
117 
118 
119 
120 ! Reordering and Factorization, This step also allocates
121 ! all memory that is necessary for the factorization
122 #ifdef DEBUG
123  msglvl = 1 ! with statistical information
124  iparm(27) = 1 ! perform matrix check
125 #else
126  msglvl = 0 ! without statistical information
127  iparm(27) = 0 ! perform matrix check
128 #endif
129  iparm(6) = 0 ! The solution is written into x
130  iparm(28) = 0 ! double real/complex
131 
132 #ifdef MIC
133  iparm(2) = 3 ! parallel OpenMP version of the nested dissection
134  iparm(24) = 1 ! Intel MKL PARDISO uses a two-level factorization algorithm
135  iparm(25) = 0 ! Intel MKL PARDISO uses a parallel algorithm for the solve step.
136  iparm(34) = 0 !
137  iparm(60) = 0 ! in-core (IC) mode
138 #else
139  iparm(2) = 2 !
140  iparm(24) = 0 ! Intel MKL PARDISO uses a single-level factorization algorithm
141  iparm(25) = 0 ! Intel MKL PARDISO uses a parallel algorithm for the solve step.
142  iparm(34) = 0 !
143  iparm(60) = 0 ! in-core (IC) mode
144 #endif
145 
146  phase = 11 ! analysis
147 
148  ! setting the appropriate number of right hand sides.
149 #ifdef MIC
150  nrhs = maxnrhs
151  if (nrhs.GT.sizea) then
152  nrhs = sizea
153  end if
154 #else
155  nrhs = 1
156 #endif
157 
158 #ifdef DEBUG
159  write(0,*) 'number of right sides: ', nrhs
160 #endif
161 
162  CALL pardiso (pt, maxfct, mnum, mtype, phase, sizea, pa, pia, pja, &
163  idum, nrhs, iparm, msglvl, ddum, ddum, error)
164  maxmem = (iparm(16)+iparm(17))
165 
166 #ifdef DEBUG
167  write(0,*) 'maxmem: ', maxmem, ' kbytes'
168  write(0,*) 'maxmem total: ', maxmem+sizea*nrhs*16/1000, ' kbytes'
169 #endif
170 
171  IF ((maxmem+sizea*nrhs*16/1000).GT.maxmem) then
172  nrhs = floor(real(maxmem-maxmem)/real(sizea*16)*1000)
173  write(0,*) 'number of right sides: ', nrhs
174  IF (nrhs.EQ.0) THEN
175  error = -2
176  RETURN
177  ENDIF
178  ENDIF
179 
180  phase = 22 ! factorization
181  CALL pardiso (pt, maxfct, mnum, mtype, phase, sizea, pa, pia, pja, &
182  idum, nrhs, iparm, msglvl, ddum, ddum, error)
183 
184 #ifdef DEBUG
185  write(0,*) 'Reordering and factorization completed ... '
186 #endif
187  IF (error .NE. 0) THEN
188  write(0,*) 'The following ERROR was detected at the factorization step: ', error
189  RETURN
190  ENDIF
191 
192 
193 
194 ! Set parameters to calculate only the diagonal components of the solution
195  iparm(31) = 3
196  phase = 33 ! solver
197  allocate( perm(sizea) )
198 
199 ! Calculating the chunk of the solution to be calculated at once
200  if (modulo(sizea, nrhs).EQ.0) then
201  numiterations = sizea/nrhs
202  else
203  numiterations = floor(real(sizea)/real(nrhs)) + 1
204  end if
205 
206 #ifdef DEBUG
207  write(0,*) 'Number of iterations to obtain the right solutions ', numiterations
208 #endif
209 
210 
211 ! Allocating the right hand side
212  allocate( pb(sizea, nrhs) )
213  allocate( px(sizea, nrhs) )
214 
215  do idxiterations = 1, numiterations
216 
217  offset = (idxiterations-1)*nrhs
218 
219  if ((sizea-offset).LT.nrhs) then
220  nrhs = sizea-offset
221  deallocate(pb, px)
222  allocate( pb(sizea, nrhs) )
223  allocate( px(sizea, nrhs) )
224  end if
225 
226 ! identify matrix component to be calculated
227  perm = 0
228  do idx = 1, nrhs
229  perm(offset+idx) = 1
230  end do
231 
232 ! Creating the right hand side
233  pb = 0.d0
234  px = 0.d0
235  do idx = 1, nrhs
236  pb(offset+idx, idx) = 1.d0
237  end do
238 ! Solver
239  CALL pardiso (pt, maxfct, mnum, mtype, phase, sizea, pa, pia, pja, &
240  perm, nrhs, iparm, msglvl, pb, px, error)
241 
242  IF (error .NE. 0) THEN
243  exit
244  ENDIF
245 
246 
247  do idx = 1, nrhs
248  inva(offset+idx) = px(offset+idx, idx)
249  end do
250 
251  end do
252 
253  deallocate(pb, px, perm)
254 
255 
256  IF (error .NE. 0) THEN
257  write(0,*) 'The following ERROR was detected at the solution step: ', error
258  RETURN
259  ENDIF
260 
261 
262 #ifdef DEBUG
263  do idx = 1,sizea
264  write(0,*) idx, 'th element of a = ', inva(idx)
265  end do
266 #endif
267 
268 
269 
270 !
271 ! Termination and release of memory
272  phase = -1 ! release internal memory
273  CALL pardiso (pt, maxfct, mnum, mtype, phase, sizea, ddum, &
274  idum, idum, &
275  idum, nrhs, iparm, msglvl, ddum, ddum, error)
276 
277  deallocate(idum, ddum, pt)
278 
279  IF (error .NE. 0) THEN
280  write(0,*) 'The following ERROR was detected at releasing the memory: ', error
281  RETURN
282  ENDIF
283 
284 
285  end
286 
287 
288 
289 ! DESCRIPTION:
299 #ifdef MIC
300 !dir$ attributes offload:mic :: zgetDiagInv
301 #endif
302  subroutine zgetdiaginv(sizeA, nonzerosA, a, ia, ja, real_invA, imag_invA, error)
304  use mkl_pardiso
305 
306  IMPLICIT NONE
307 
308 ! Internal solver memory pointer
309  TYPE(mkl_pardiso_handle), pointer :: PT(:) => null()
310  integer*4, pointer :: idum(:) => null(), pia(:) => null(), pja(:) => null()
311  complex(KIND=8), pointer :: ddum(:) => null(), pa(:) => null(), pb(:,:) => null(), px(:,:) => null()
312  INTEGER*4 maxfct, mnum, mtype, phase, nrhs, error, msglvl, nnz, numIterations, idxIterations, offset
313  INTEGER*4 iparm(64)
314 
315 
316 ! The matrix is given in the CSR format A = (a, ia, ja)
317  INTEGER*4 sizeA, nonzerosA !matrix A is of size sizeA x sizeA
318  INTEGER*4, TARGET :: ia(sizeA+1), ja(nonzerosA)
319  COMPLEX*16, TARGET :: a(nonzerosA)
320 
321  INTEGER*4, ALLOCATABLE :: perm(:)
322 
323  real*8 :: real_inva(sizea)
324  real*8 :: imag_inva(sizea)
325 
326  INTEGER*4 idx, jdx, maxmem
327 
328 
329 
330  DATA maxfct /1/, mnum /1/
331  allocate( pt(64), idum(0), ddum(0) )
332 
333  do idx = 1, 64
334  pt(idx)%DUMMY = 0
335  end do
336  pa => a
337  pia => ia
338  pja => ja
339 ! initialize error flag
340  error = 0
341 
342 
343 #ifdef DEBUG
344  write(0,*)
345  write(0,*) 'The input matrix: (getPartialInv)'
346  idx = 0
347  do jdx = 1, nonzerosa
348  if (jdx >= ia(idx+1)) then
349  idx=idx+1
350  end if
351 
352  if (jdx < ia(idx+1)) then
353  !write(0,*) 'row', idx, 'col', ja(jdx), 'element of a = ', a(jdx)
354  end if
355  end do
356  write(0,*) 'The size of the diagonal vector to be calculated:', sizea
357 #endif
358 
359 !
360 ! .. Setup Pardiso control parameters und initialize the solvers
361 ! internal adress pointers. This is only necessary for the FIRST
362 ! call of the PARDISO solver.
363 !
364 ! choose the matrix type for pardiso factorization
365  mtype = 3 ! structurally symmetric complex matrix
366 
367 ! .. PARDISO license check and initialize solver
368  call pardisoinit(pt, mtype, iparm)
369 
370 
371 
372 
373 ! Reordering and Factorization, This step also allocates
374 ! all memory that is necessary for the factorization
375 #ifdef DEBUG
376  msglvl = 1 ! with statistical information
377  iparm(27) = 1 ! perform matrix check
378 #else
379  msglvl = 0 ! without statistical information
380  iparm(27) = 0 ! perform matrix check
381 #endif
382  iparm(6) = 0 ! The solution is written into x
383  iparm(28) = 0 ! double real/complex
384 
385 #ifdef MIC
386  iparm(2) = 3 ! parallel OpenMP version of the nested dissection
387  iparm(24) = 1 ! Intel MKL PARDISO uses a two-level factorization algorithm
388  iparm(25) = 0 ! Intel MKL PARDISO uses a parallel algorithm for the solve step.
389  iparm(34) = 0 !
390  iparm(60) = 0 ! in-core (IC) mode
391 #else
392  iparm(2) = 2 ! parallel OpenMP version of the nested dissection
393  iparm(24) = 0 ! Intel MKL PARDISO uses a two-level factorization algorithm
394  iparm(25) = 0 ! Intel MKL PARDISO uses a parallel algorithm for the solve step.
395  iparm(34) = 0 !
396  iparm(60) = 0 ! in-core (IC) mode
397 #endif
398 
399  phase = 11 ! Analysis
400 
401  ! setting the appropriate number of right hand sides.
402 #ifdef MIC
403  nrhs = maxnrhs
404  if (nrhs.GT.sizea) then
405  nrhs = sizea
406  end if
407 #else
408  nrhs = 1
409 #endif
410 
411 #ifdef DEBUG
412  write(0,*) 'number of right sides: ', nrhs
413 #endif
414 
415 
416 
417  CALL pardiso (pt, maxfct, mnum, mtype, phase, sizea, pa, pia, pja, &
418  idum, nrhs, iparm, msglvl, ddum, ddum, error)
419  maxmem = (iparm(16)+iparm(17))
420 
421 #ifdef DEBUG
422  write(0,*) 'maxmem: ', maxmem, ' kbytes'
423  write(0,*) 'maxmem total: ', maxmem+sizea*nrhs*32/1000, ' kbytes'
424 #endif
425 
426  IF ((maxmem+sizea*nrhs*32/1000).GT.maxmem) then
427  nrhs = floor(real(maxmem-maxmem)/real(sizea*32)*1000)
428  write(0,*) 'number of right sides: ', nrhs
429  IF (nrhs.EQ.0) THEN
430  error = -2
431  RETURN
432  ENDIF
433  ENDIF
434 
435  phase = 22 ! factorization
436  CALL pardiso (pt, maxfct, mnum, mtype, phase, sizea, pa, pia, pja, &
437  idum, nrhs, iparm, msglvl, ddum, ddum, error)
438 
439 #ifdef DEBUG
440  write(0,*) 'Reordering and factorization completed ... '
441 #endif
442  IF (error .NE. 0) THEN
443  write(0,*) 'The following ERROR was detected at the factorization step: ', error
444  RETURN
445  ENDIF
446 
447 
448 
449 ! Set parameters to calculate only the diagonal components of the solution
450  iparm(31) = 3
451  phase = 33 ! solver
452  allocate( perm(sizea) )
453 
454 ! Calculating the chunk of the solution to be calculated at once
455  if (modulo(sizea, nrhs).EQ.0) then
456  numiterations = sizea/nrhs
457  else
458  numiterations = floor(real(sizea)/real(nrhs)) + 1
459  end if
460 
461 #ifdef DEBUG
462  write(0,*) 'Number of iterations to obtain the right solutions ', numiterations
463 #endif
464 
465 ! Allocating the right hand side
466  allocate( pb(sizea, nrhs) )
467  allocate( px(sizea, nrhs) )
468 
469  do idxiterations = 1, numiterations
470 
471  offset = (idxiterations-1)*nrhs
472 
473  if ((sizea-offset).LT.nrhs) then
474  nrhs = sizea-offset
475  deallocate(pb, px)
476  allocate( pb(sizea, nrhs) )
477  allocate( px(sizea, nrhs) )
478  end if
479 
480 ! identify matrix component to be calculated
481  perm = 0
482  do idx = 1, nrhs
483  perm(offset+idx) = 1
484  end do
485 
486 ! Creating the right hand side
487  pb = (0.d0, 0.d0)
488  px = (0.d0, 0.d0)
489  do idx = 1, nrhs
490  pb(offset+idx, idx) = (1.d0, 0.d0)
491  end do
492 ! Solver
493  CALL pardiso (pt, maxfct, mnum, mtype, phase, sizea, pa, pia, pja, &
494  perm, nrhs, iparm, msglvl, pb, px, error)
495 
496  IF (error .NE. 0) THEN
497  exit
498  ENDIF
499 
500  do idx = 1, nrhs
501  real_inva(offset+idx) = real(px(offset+idx, idx))
502  imag_inva(offset+idx) = aimag(px(offset+idx, idx))
503  end do
504 
505  end do
506 
507  deallocate(pb, px, perm)
508 
509 
510  IF (error .NE. 0) THEN
511  write(0,*) 'The following ERROR was detected at the solution step: ', error
512  RETURN
513  ENDIF
514 
515 !
516 ! Termination and release of memory
517  phase = -1 ! release internal memory
518  CALL pardiso (pt, maxfct, mnum, mtype, phase, sizea, ddum, &
519  idum, idum, &
520  idum, nrhs, iparm, msglvl, ddum, ddum, error)
521 
522  deallocate(idum, ddum, pt)
523 
524  IF (error .NE. 0) THEN
525  write(0,*) 'The following ERROR was detected at releasing the memory: ', error
526  RETURN
527  ENDIF
528 
529  end
530 
531 
532 
533 end module
534 
535 
Module to calculate the diagonal of the inverse of a real/complex sparse matrix via the PARDISO libra...
Definition: getDiagInv.F90:24
subroutine zgetdiaginv(sizeA, nonzerosA, a, ia, ja, real_invA, imag_invA, error)
Get the diagonal part of the inverse of double valued complex matrices.
Definition: getDiagInv.F90:303
subroutine dgetdiaginv(sizeA, nonzerosA, a, ia, ja, invA, error)
Get the diagonal part of the inverse of double valued real matrices.
Definition: getDiagInv.F90:51