[5725] | 1 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 2 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 3 | !----------------------------------------------------------------------- |
---|
| 4 | ! CVS $Id$ |
---|
| 5 | ! CVS $Name$ |
---|
| 6 | !BOP ------------------------------------------------------------------- |
---|
| 7 | ! |
---|
| 8 | ! !MODULE: m_SparseMatrix -- Sparse Matrix Object |
---|
| 9 | ! |
---|
| 10 | ! !DESCRIPTION: |
---|
| 11 | ! The {\tt SparseMatrix} data type is MCT's object for storing sparse |
---|
| 12 | ! matrices. In MCT, intergrid interpolation is implemented as a sparse |
---|
| 13 | ! matrix-vector multiplication, with the {\tt AttrVect} type playing the |
---|
| 14 | ! roles of the input and output vectors. The interpolation matrices tend |
---|
| 15 | ! to be {\em extremely} sparse. For ${\bf x} \in \Re^{N_x}$, and |
---|
| 16 | ! ${\bf y} \in \Re^{N_y}$, the interpolation matrix {\bf M} used to effect |
---|
| 17 | ! ${\bf y} = {\bf M} {\bf x}$ will typically have ${\cal O}({N_y})$ |
---|
| 18 | ! non-zero elements. For that reason, the {\tt SparseMatrix} type |
---|
| 19 | ! stores {\em only} information about non-zero matrix elements, along |
---|
| 20 | ! with the number of rows and columns in the full matrix. The nonzero |
---|
| 21 | ! matrix elements are stored in {\tt AttrVect} form (see the module |
---|
| 22 | ! {\tt m\_AttrVect} for more details), and the set of attributes are |
---|
| 23 | ! listed below: |
---|
| 24 | ! |
---|
| 25 | !\begin{table}[htbp] |
---|
| 26 | !\begin{center} |
---|
| 27 | !\begin{tabular}{|l|l|l|} |
---|
| 28 | !\hline |
---|
| 29 | !{\bf Attribute Name} & {\bf Significance} & {\tt Type} \\ |
---|
| 30 | !\hline |
---|
| 31 | !{\tt grow} & Global Row Index & {\tt INTEGER} \\ |
---|
| 32 | !\hline |
---|
| 33 | !{\tt gcol} & Global Column Index & {\tt INTEGER} \\ |
---|
| 34 | !\hline |
---|
| 35 | !{\tt lrow} & Local Row Index & {\tt INTEGER} \\ |
---|
| 36 | !\hline |
---|
| 37 | !{\tt lcol} & Local Column Index & {\tt INTEGER} \\ |
---|
| 38 | !\hline |
---|
| 39 | !{\tt weight} & Matrix Element ${M_{ij}}$ & {\tt REAL} \\ |
---|
| 40 | !\hline |
---|
| 41 | !\end{tabular} |
---|
| 42 | !\end{center} |
---|
| 43 | !\end{table} |
---|
| 44 | ! |
---|
| 45 | ! The provision of both local and global column and row indices is |
---|
| 46 | ! made because this datatype can be used in either shared-memory or |
---|
| 47 | ! distributed-memory parallel matrix-vector products. |
---|
| 48 | ! |
---|
| 49 | ! This module contains the definition of the {\tt SparseMatrix} type, |
---|
| 50 | ! creation and destruction methods, a variety of accessor methods, |
---|
| 51 | ! routines for testing the suitability of the matrix for interpolation |
---|
| 52 | ! (i.e. the sum of each row is either zero or unity), and methods for |
---|
| 53 | ! sorting and permuting matrix entries. |
---|
| 54 | ! |
---|
| 55 | ! For better performance of the Matrix-Vector multiply on vector |
---|
| 56 | ! architectures, the {\tt SparseMatrix} object also contains arrays |
---|
| 57 | ! for holding the sparse matrix data in a more vector-friendly form. |
---|
| 58 | ! |
---|
| 59 | ! |
---|
| 60 | ! !INTERFACE: |
---|
| 61 | |
---|
| 62 | module m_SparseMatrix |
---|
| 63 | ! |
---|
| 64 | ! !USES: |
---|
| 65 | ! |
---|
| 66 | use m_realkinds, only : FP |
---|
| 67 | use m_AttrVect, only : AttrVect |
---|
| 68 | |
---|
| 69 | |
---|
| 70 | private ! except |
---|
| 71 | |
---|
| 72 | ! !PUBLIC TYPES: |
---|
| 73 | |
---|
| 74 | public :: SparseMatrix ! The class data structure |
---|
| 75 | |
---|
| 76 | Type SparseMatrix |
---|
| 77 | #ifdef SEQUENCE |
---|
| 78 | sequence |
---|
| 79 | #endif |
---|
| 80 | integer :: nrows |
---|
| 81 | integer :: ncols |
---|
| 82 | type(AttrVect) :: data |
---|
| 83 | logical :: vecinit ! additional data for the vectorized sMat |
---|
| 84 | integer,dimension(:),pointer :: row_s, row_e |
---|
| 85 | integer, dimension(:,:), pointer :: tcol |
---|
| 86 | real(FP), dimension(:,:), pointer :: twgt |
---|
| 87 | integer :: row_max, row_min |
---|
| 88 | integer :: tbl_end |
---|
| 89 | End Type SparseMatrix |
---|
| 90 | |
---|
| 91 | ! !PUBLIC MEMBER FUNCTIONS: |
---|
| 92 | |
---|
| 93 | public :: init ! Create a SparseMatrix |
---|
| 94 | public :: vecinit ! Initialize the vector parts |
---|
| 95 | public :: clean ! Destroy a SparseMatrix |
---|
| 96 | public :: lsize ! Local number of elements |
---|
| 97 | public :: indexIA ! Index integer attribute |
---|
| 98 | public :: indexRA ! Index real attribute |
---|
| 99 | public :: nRows ! Total number of rows |
---|
| 100 | public :: nCols ! Total number of columns |
---|
| 101 | |
---|
| 102 | public :: exportGlobalRowIndices ! Return global row indices |
---|
| 103 | ! for matrix elements |
---|
| 104 | public :: exportGlobalColumnIndices ! Return global column indices |
---|
| 105 | ! for matrix elements |
---|
| 106 | public :: exportLocalRowIndices ! Return local row indices |
---|
| 107 | ! for matrix elements |
---|
| 108 | public :: exportLocalColumnIndices ! Return local column indices |
---|
| 109 | ! for matrix elements |
---|
| 110 | public :: exportMatrixElements ! Return matrix elements |
---|
| 111 | |
---|
| 112 | public :: importGlobalRowIndices ! Set global row indices |
---|
| 113 | ! using |
---|
| 114 | public :: importGlobalColumnIndices ! Return global column indices |
---|
| 115 | ! for matrix elements |
---|
| 116 | public :: importLocalRowIndices ! Return local row indices |
---|
| 117 | ! for matrix elements |
---|
| 118 | public :: importLocalColumnIndices ! Return local column indices |
---|
| 119 | ! for matrix elements |
---|
| 120 | public :: importMatrixElements ! Return matrix elements |
---|
| 121 | public :: Copy ! Copy a SparseMatrix |
---|
| 122 | |
---|
| 123 | public :: GlobalNumElements ! Total number of nonzero elements |
---|
| 124 | public :: ComputeSparsity ! Fraction of matrix that is nonzero |
---|
| 125 | public :: local_row_range ! Local (on-process) row range |
---|
| 126 | public :: global_row_range ! Local (on-process) row range |
---|
| 127 | public :: local_col_range ! Local (on-process) column range |
---|
| 128 | public :: global_col_range ! Local (on-process) column range |
---|
| 129 | public :: CheckBounds ! Check row and column values |
---|
| 130 | ! for out-of-bounds values |
---|
| 131 | public :: row_sum ! Return SparseMatrix row sums |
---|
| 132 | public :: row_sum_check ! Check SparseMatrix row sums against |
---|
| 133 | ! input "valid" values |
---|
| 134 | public :: Sort ! Sort matrix entries to generate an |
---|
| 135 | ! index permutation (to be used by |
---|
| 136 | ! Permute() |
---|
| 137 | public :: Permute ! Permute matrix entries using index |
---|
| 138 | ! permutation gernerated by Sort() |
---|
| 139 | public :: SortPermute ! Sort/Permute matrix entries |
---|
| 140 | |
---|
| 141 | interface init ; module procedure init_ ; end interface |
---|
| 142 | interface vecinit ; module procedure vecinit_ ; end interface |
---|
| 143 | interface clean ; module procedure clean_ ; end interface |
---|
| 144 | interface lsize ; module procedure lsize_ ; end interface |
---|
| 145 | interface indexIA ; module procedure indexIA_ ; end interface |
---|
| 146 | interface indexRA ; module procedure indexRA_ ; end interface |
---|
| 147 | interface nRows ; module procedure nRows_ ; end interface |
---|
| 148 | interface nCols ; module procedure nCols_ ; end interface |
---|
| 149 | |
---|
| 150 | interface exportGlobalRowIndices ; module procedure & |
---|
| 151 | exportGlobalRowIndices_ |
---|
| 152 | end interface |
---|
| 153 | |
---|
| 154 | interface exportGlobalColumnIndices ; module procedure & |
---|
| 155 | exportGlobalColumnIndices_ |
---|
| 156 | end interface |
---|
| 157 | |
---|
| 158 | interface exportLocalRowIndices ; module procedure & |
---|
| 159 | exportLocalRowIndices_ |
---|
| 160 | end interface |
---|
| 161 | |
---|
| 162 | interface exportLocalColumnIndices ; module procedure & |
---|
| 163 | exportLocalColumnIndices_ |
---|
| 164 | end interface |
---|
| 165 | |
---|
| 166 | interface exportMatrixElements ; module procedure & |
---|
| 167 | exportMatrixElementsSP_, & |
---|
| 168 | exportMatrixElementsDP_ |
---|
| 169 | end interface |
---|
| 170 | |
---|
| 171 | interface importGlobalRowIndices ; module procedure & |
---|
| 172 | importGlobalRowIndices_ |
---|
| 173 | end interface |
---|
| 174 | |
---|
| 175 | interface importGlobalColumnIndices ; module procedure & |
---|
| 176 | importGlobalColumnIndices_ |
---|
| 177 | end interface |
---|
| 178 | |
---|
| 179 | interface importLocalRowIndices ; module procedure & |
---|
| 180 | importLocalRowIndices_ |
---|
| 181 | end interface |
---|
| 182 | |
---|
| 183 | interface importLocalColumnIndices ; module procedure & |
---|
| 184 | importLocalColumnIndices_ |
---|
| 185 | end interface |
---|
| 186 | |
---|
| 187 | interface importMatrixElements ; module procedure & |
---|
| 188 | importMatrixElementsSP_, & |
---|
| 189 | importMatrixElementsDP_ |
---|
| 190 | end interface |
---|
| 191 | |
---|
| 192 | interface Copy ; module procedure Copy_ ; end interface |
---|
| 193 | |
---|
| 194 | interface GlobalNumElements ; module procedure & |
---|
| 195 | GlobalNumElements_ |
---|
| 196 | end interface |
---|
| 197 | |
---|
| 198 | interface ComputeSparsity ; module procedure & |
---|
| 199 | ComputeSparsitySP_, & |
---|
| 200 | ComputeSparsityDP_ |
---|
| 201 | end interface |
---|
| 202 | |
---|
| 203 | interface local_row_range ; module procedure & |
---|
| 204 | local_row_range_ |
---|
| 205 | end interface |
---|
| 206 | |
---|
| 207 | interface global_row_range ; module procedure & |
---|
| 208 | global_row_range_ |
---|
| 209 | end interface |
---|
| 210 | |
---|
| 211 | interface local_col_range ; module procedure & |
---|
| 212 | local_col_range_ |
---|
| 213 | end interface |
---|
| 214 | |
---|
| 215 | interface global_col_range ; module procedure & |
---|
| 216 | global_col_range_ |
---|
| 217 | end interface |
---|
| 218 | |
---|
| 219 | interface CheckBounds; module procedure & |
---|
| 220 | CheckBounds_ |
---|
| 221 | end interface |
---|
| 222 | |
---|
| 223 | interface row_sum ; module procedure & |
---|
| 224 | row_sumSP_, & |
---|
| 225 | row_sumDP_ |
---|
| 226 | end interface |
---|
| 227 | |
---|
| 228 | interface row_sum_check ; module procedure & |
---|
| 229 | row_sum_checkSP_, & |
---|
| 230 | row_sum_checkDP_ |
---|
| 231 | end interface |
---|
| 232 | |
---|
| 233 | interface Sort ; module procedure Sort_ ; end interface |
---|
| 234 | interface Permute ; module procedure Permute_ ; end interface |
---|
| 235 | interface SortPermute ; module procedure SortPermute_ ; end interface |
---|
| 236 | |
---|
| 237 | ! !REVISION HISTORY: |
---|
| 238 | ! 19Sep00 - J.W. Larson <larson@mcs.anl.gov> - initial prototype |
---|
| 239 | ! 15Jan01 - J.W. Larson <larson@mcs.anl.gov> - added numerous APIs |
---|
| 240 | ! 25Feb01 - J.W. Larson <larson@mcs.anl.gov> - changed from row/column |
---|
| 241 | ! attributes to global and local row and column attributes |
---|
| 242 | ! 23Apr01 - J.W. Larson <larson@mcs.anl.gov> - added number of rows |
---|
| 243 | ! and columns to the SparseMatrix type. This means the |
---|
| 244 | ! SparseMatrix is no longer a straight AttrVect type. This |
---|
| 245 | ! also made necessary the addition of lsize(), indexIA(), |
---|
| 246 | ! and indexRA(). |
---|
| 247 | ! 29Oct03 - R. Jacob <jacob@mcs.anl.gov> - extend the SparseMatrix type |
---|
| 248 | ! to include mods from Fujitsu for a vector-friendly MatVecMul |
---|
| 249 | !EOP ___________________________________________________________________ |
---|
| 250 | |
---|
| 251 | character(len=*),parameter :: myname='MCT::m_SparseMatrix' |
---|
| 252 | |
---|
| 253 | ! SparseMatrix_iList components: |
---|
| 254 | character(len=*),parameter :: SparseMatrix_iList='grow:gcol:lrow:lcol' |
---|
| 255 | integer,parameter :: SparseMatrix_igrow=1 |
---|
| 256 | integer,parameter :: SparseMatrix_igcol=2 |
---|
| 257 | integer,parameter :: SparseMatrix_ilrow=3 |
---|
| 258 | integer,parameter :: SparseMatrix_ilcol=4 |
---|
| 259 | |
---|
| 260 | ! SparseMatrix_rList components: |
---|
| 261 | character(len=*),parameter :: SparseMatrix_rList='weight' |
---|
| 262 | integer,parameter :: SparseMatrix_iweight=1 |
---|
| 263 | |
---|
| 264 | contains |
---|
| 265 | |
---|
| 266 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 267 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 268 | !BOP ------------------------------------------------------------------- |
---|
| 269 | ! |
---|
| 270 | ! !IROUTINE: init_ - Initialize an Empty SparseMatrix |
---|
| 271 | ! |
---|
| 272 | ! !DESCRIPTION: This routine creates the storage space for the |
---|
| 273 | ! entries of a {\tt SparseMatrix}, and sets the number of rows and |
---|
| 274 | ! columns in it. The input {\tt INTEGER} arguments {\tt nrows} and |
---|
| 275 | ! {\tt ncols} specify the number of rows and columns respectively. |
---|
| 276 | ! The optional input argument {\tt lsize} specifies the number of |
---|
| 277 | ! nonzero entries in the {\tt SparseMatrix}. The initialized |
---|
| 278 | ! {\tt SparseMatrix} is returned in the output argument {\tt sMat}. |
---|
| 279 | ! |
---|
| 280 | ! {\bf N.B.}: This routine is allocating dynamical memory in the form |
---|
| 281 | ! of a {\tt SparseMatrix}. The user must deallocate this space when |
---|
| 282 | ! the {\tt SparseMatrix} is no longer needed by invoking the routine |
---|
| 283 | ! {\tt clean\_()}. |
---|
| 284 | ! |
---|
| 285 | ! !INTERFACE: |
---|
| 286 | |
---|
| 287 | subroutine init_(sMat, nrows, ncols, lsize) |
---|
| 288 | ! |
---|
| 289 | ! !USES: |
---|
| 290 | ! |
---|
| 291 | use m_AttrVect, only : AttrVect_init => init |
---|
| 292 | use m_die |
---|
| 293 | |
---|
| 294 | implicit none |
---|
| 295 | |
---|
| 296 | ! !INPUT PARAMETERS: |
---|
| 297 | |
---|
| 298 | integer, intent(in) :: nrows |
---|
| 299 | integer, intent(in) :: ncols |
---|
| 300 | integer, optional, intent(in) :: lsize |
---|
| 301 | |
---|
| 302 | ! !OUTPUT PARAMETERS: |
---|
| 303 | |
---|
| 304 | type(SparseMatrix), intent(out) :: sMat |
---|
| 305 | |
---|
| 306 | ! !REVISION HISTORY: |
---|
| 307 | ! 19Sep00 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
| 308 | ! 23Apr01 - Jay Larson <larson@mcs.anl.gov> - added arguments |
---|
| 309 | ! nrows and ncols--number of rows and columns in the |
---|
| 310 | ! SparseMatrix |
---|
| 311 | !EOP ___________________________________________________________________ |
---|
| 312 | ! |
---|
| 313 | character(len=*),parameter :: myname_=myname//'::init_' |
---|
| 314 | |
---|
| 315 | integer :: n |
---|
| 316 | |
---|
| 317 | ! if lsize is present, use it to set n; if not, set n=0 |
---|
| 318 | |
---|
| 319 | n = 0 |
---|
| 320 | if(present(lsize)) n=lsize |
---|
| 321 | |
---|
| 322 | ! Initialize number of rows and columns: |
---|
| 323 | |
---|
| 324 | sMat%nrows = nrows |
---|
| 325 | sMat%ncols = ncols |
---|
| 326 | |
---|
| 327 | ! Initialize sMat%data using AttrVect_init |
---|
| 328 | |
---|
| 329 | call AttrVect_init(sMat%data, SparseMatrix_iList, & |
---|
| 330 | SparseMatrix_rList, n) |
---|
| 331 | |
---|
| 332 | ! vecinit is off by default |
---|
| 333 | sMat%vecinit = .FALSE. |
---|
| 334 | |
---|
| 335 | end subroutine init_ |
---|
| 336 | |
---|
| 337 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 338 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 339 | !BOP ------------------------------------------------------------------- |
---|
| 340 | ! |
---|
| 341 | ! !IROUTINE: vecinit_ - Initialize vector parts of a SparseMatrix |
---|
| 342 | ! |
---|
| 343 | ! !DESCRIPTION: This routine creates the storage space for |
---|
| 344 | ! and intializes the vector parts of a {\tt SparseMatrix}. |
---|
| 345 | ! |
---|
| 346 | ! {\bf N.B.}: This routine assumes the locally indexed parts of a |
---|
| 347 | ! {\tt SparseMatrix} have been initialized. This is |
---|
| 348 | ! accomplished by either importing the values directly with |
---|
| 349 | ! {\tt importLocalRowIndices} and {\tt importLocalColIndices} or by |
---|
| 350 | ! importing the Global Row and Col Indices and making two calls to |
---|
| 351 | ! {\tt GlobalToLocalMatrix}. |
---|
| 352 | ! |
---|
| 353 | ! {\bf N.B.}: The vector portion can use a large amount of |
---|
| 354 | ! memory so it is highly recommended that this routine only |
---|
| 355 | ! be called on a {\tt SparseMatrix} that has been scattered |
---|
| 356 | ! or otherwise sized locally. |
---|
| 357 | ! |
---|
| 358 | ! !INTERFACE: |
---|
| 359 | |
---|
| 360 | subroutine vecinit_(sMat) |
---|
| 361 | ! |
---|
| 362 | ! !USES: |
---|
| 363 | ! |
---|
| 364 | use m_die |
---|
| 365 | use m_stdio |
---|
| 366 | |
---|
| 367 | implicit none |
---|
| 368 | |
---|
| 369 | ! !INPUT/OUTPUT PARAMETERS: |
---|
| 370 | |
---|
| 371 | type(SparseMatrix), intent(inout) :: sMat |
---|
| 372 | |
---|
| 373 | ! !REVISION HISTORY: |
---|
| 374 | ! 27Oct03 - R. Jacob <jacob@mcs.anl.gov> - initial version |
---|
| 375 | ! using code provided by Yoshi et. al. |
---|
| 376 | !EOP ___________________________________________________________________ |
---|
| 377 | ! |
---|
| 378 | character(len=*),parameter :: myname_=myname//'::vecinit_' |
---|
| 379 | |
---|
| 380 | integer :: irow,icol,iwgt |
---|
| 381 | integer :: num_elements |
---|
| 382 | integer :: row,col |
---|
| 383 | integer :: ier,l,n |
---|
| 384 | integer, dimension(:) , allocatable :: nr, rn |
---|
| 385 | |
---|
| 386 | if(sMat%vecinit) then |
---|
| 387 | write(stderr,'(2a)') myname_, & |
---|
| 388 | 'MCTERROR: sMat vector parts have already been initialized...Continuing' |
---|
| 389 | RETURN |
---|
| 390 | endif |
---|
| 391 | |
---|
| 392 | write(6,*) myname_,'Initializing vecMat' |
---|
| 393 | irow = indexIA_(sMat,'lrow',dieWith=myname_) |
---|
| 394 | icol = indexIA_(sMat,'lcol',dieWith=myname_) |
---|
| 395 | iwgt = indexRA_(sMat,'weight',dieWith=myname_) |
---|
| 396 | |
---|
| 397 | num_elements = lsize_(sMat) |
---|
| 398 | |
---|
| 399 | sMat%row_min = sMat%data%iAttr(irow,1) |
---|
| 400 | sMat%row_max = sMat%row_min |
---|
| 401 | do n=1,num_elements |
---|
| 402 | row = sMat%data%iAttr(irow,n) |
---|
| 403 | if ( row > sMat%row_max ) sMat%row_max = row |
---|
| 404 | if ( row < sMat%row_min ) sMat%row_min = row |
---|
| 405 | enddo |
---|
| 406 | |
---|
| 407 | allocate( nr(sMat%row_max), rn(num_elements), stat=ier) |
---|
| 408 | if(ier/=0) call die(myname_,'allocate(nr,rn)',ier) |
---|
| 409 | |
---|
| 410 | sMat%tbl_end = 0 |
---|
| 411 | nr(:) = 0 |
---|
| 412 | do n=1,num_elements |
---|
| 413 | row = sMat%data%iAttr(irow,n) |
---|
| 414 | nr(row) = nr(row)+1 |
---|
| 415 | rn(n) = nr(row) |
---|
| 416 | enddo |
---|
| 417 | sMat%tbl_end = maxval(rn) |
---|
| 418 | |
---|
| 419 | allocate( sMat%tcol(sMat%row_max,sMat%tbl_end), & |
---|
| 420 | sMat%twgt(sMat%row_max,sMat%tbl_end), stat=ier ) |
---|
| 421 | if(ier/=0) call die(myname_,'allocate(tcol,twgt)',ier) |
---|
| 422 | |
---|
| 423 | !CDIR COLLAPSE |
---|
| 424 | sMat%tcol(:,:) = -1 |
---|
| 425 | do n=1,num_elements |
---|
| 426 | row = sMat%data%iAttr(irow,n) |
---|
| 427 | sMat%tcol(row,rn(n)) = sMat%data%iAttr(icol,n) |
---|
| 428 | sMat%twgt(row,rn(n)) = sMat%data%rAttr(iwgt,n) |
---|
| 429 | enddo |
---|
| 430 | |
---|
| 431 | allocate( sMat%row_s(sMat%tbl_end) , sMat%row_e(sMat%tbl_end), & |
---|
| 432 | stat=ier ) |
---|
| 433 | if(ier/=0) call die(myname_,'allocate(row_s,row_e',ier) |
---|
| 434 | sMat%row_s = sMat%row_min |
---|
| 435 | sMat%row_e = sMat%row_max |
---|
| 436 | do l=1,sMat%tbl_end |
---|
| 437 | do n=sMat%row_min,sMat%row_max |
---|
| 438 | if (nr(n) >= l) then |
---|
| 439 | sMat%row_s(l) = n |
---|
| 440 | exit |
---|
| 441 | endif |
---|
| 442 | enddo |
---|
| 443 | do n = sMat%row_max,sMat%row_min,-1 |
---|
| 444 | if (nr(n) >= l) then |
---|
| 445 | sMat%row_e(l) = n |
---|
| 446 | exit |
---|
| 447 | endif |
---|
| 448 | enddo |
---|
| 449 | enddo |
---|
| 450 | |
---|
| 451 | deallocate(nr,rn, stat=ier) |
---|
| 452 | if(ier/=0) call die(myname_,'deallocate()',ier) |
---|
| 453 | |
---|
| 454 | sMat%vecinit = .TRUE. |
---|
| 455 | |
---|
| 456 | end subroutine vecinit_ |
---|
| 457 | |
---|
| 458 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 459 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 460 | !BOP ------------------------------------------------------------------- |
---|
| 461 | ! |
---|
| 462 | ! !IROUTINE: clean_ - Destroy a SparseMatrix. |
---|
| 463 | ! |
---|
| 464 | ! !DESCRIPTION: This routine deallocates dynamical memory held by the |
---|
| 465 | ! input {\tt SparseMatrix} argument {\tt sMat}. It also sets the number |
---|
| 466 | ! of rows and columns in the {\tt SparseMatrix} to zero. |
---|
| 467 | ! |
---|
| 468 | ! !INTERFACE: |
---|
| 469 | |
---|
| 470 | subroutine clean_(sMat,stat) |
---|
| 471 | ! |
---|
| 472 | ! !USES: |
---|
| 473 | ! |
---|
| 474 | use m_AttrVect,only : AttrVect_clean => clean |
---|
| 475 | use m_die |
---|
| 476 | |
---|
| 477 | implicit none |
---|
| 478 | |
---|
| 479 | ! !INPUT/OUTPTU PARAMETERS: |
---|
| 480 | |
---|
| 481 | type(SparseMatrix), intent(inout) :: sMat |
---|
| 482 | |
---|
| 483 | ! !OUTPUT PARAMETERS: |
---|
| 484 | |
---|
| 485 | integer, optional, intent(out) :: stat |
---|
| 486 | |
---|
| 487 | ! !REVISION HISTORY: |
---|
| 488 | ! 19Sep00 - J.W. Larson <larson@mcs.anl.gov> - initial prototype |
---|
| 489 | ! 23Apr00 - J.W. Larson <larson@mcs.anl.gov> - added changes to |
---|
| 490 | ! accomodate clearing nrows and ncols. |
---|
| 491 | ! 01Mar02 - E.T. Ong <eong@mcs.anl.gov> Added stat argument. |
---|
| 492 | ! 03Oct03 - R. Jacob <jacob@mcs.anl.gov> - clean vector parts |
---|
| 493 | !EOP ___________________________________________________________________ |
---|
| 494 | |
---|
| 495 | character(len=*),parameter :: myname_=myname//'::clean_' |
---|
| 496 | integer :: ier |
---|
| 497 | |
---|
| 498 | ! Deallocate memory held by sMat: |
---|
| 499 | |
---|
| 500 | if(present(stat)) then |
---|
| 501 | call AttrVect_clean(sMat%data,stat) |
---|
| 502 | else |
---|
| 503 | call AttrVect_clean(sMat%data) |
---|
| 504 | endif |
---|
| 505 | |
---|
| 506 | ! Set the number of rows and columns in sMat to zero: |
---|
| 507 | |
---|
| 508 | sMat%nrows = 0 |
---|
| 509 | sMat%ncols = 0 |
---|
| 510 | |
---|
| 511 | if(sMat%vecinit) then |
---|
| 512 | sMat%row_max = 0 |
---|
| 513 | sMat%row_min = 0 |
---|
| 514 | sMat%tbl_end = 0 |
---|
| 515 | deallocate(sMat%row_s,sMat%row_e,stat=ier) |
---|
| 516 | if(ier/=0) then |
---|
| 517 | if(present(stat)) then |
---|
| 518 | stat=ier |
---|
| 519 | else |
---|
| 520 | call warn(myname_,'deallocate(row_s,row_e)',ier) |
---|
| 521 | endif |
---|
| 522 | endif |
---|
| 523 | |
---|
| 524 | deallocate(sMat%tcol,sMat%twgt,stat=ier) |
---|
| 525 | if(ier/=0) then |
---|
| 526 | if(present(stat)) then |
---|
| 527 | stat=ier |
---|
| 528 | else |
---|
| 529 | call warn(myname_,'deallocate(tcol,twgt)',ier) |
---|
| 530 | endif |
---|
| 531 | endif |
---|
| 532 | sMat%vecinit = .FALSE. |
---|
| 533 | endif |
---|
| 534 | |
---|
| 535 | |
---|
| 536 | end subroutine clean_ |
---|
| 537 | |
---|
| 538 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 539 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 540 | !BOP ------------------------------------------------------------------- |
---|
| 541 | ! |
---|
| 542 | ! !IROUTINE: lsize_ - Local Number Non-zero Elements |
---|
| 543 | ! |
---|
| 544 | ! !DESCRIPTION: This {\tt INTEGER} function reports on-processor storage |
---|
| 545 | ! of the number of nonzero elements in the input {\tt SparseMatrix} |
---|
| 546 | ! argument {\tt sMat}. |
---|
| 547 | ! |
---|
| 548 | ! !INTERFACE: |
---|
| 549 | |
---|
| 550 | integer function lsize_(sMat) |
---|
| 551 | ! |
---|
| 552 | ! !USES: |
---|
| 553 | ! |
---|
| 554 | use m_AttrVect,only : AttrVect_lsize => lsize |
---|
| 555 | |
---|
| 556 | implicit none |
---|
| 557 | |
---|
| 558 | ! !INPUT PARAMETERS: |
---|
| 559 | |
---|
| 560 | type(SparseMatrix), intent(in) :: sMat |
---|
| 561 | |
---|
| 562 | ! !REVISION HISTORY: |
---|
| 563 | ! 23Apr00 - J.W. Larson <larson@mcs.anl.gov> - initial version. |
---|
| 564 | !EOP ___________________________________________________________________ |
---|
| 565 | |
---|
| 566 | character(len=*),parameter :: myname_=myname//'::lsize_' |
---|
| 567 | |
---|
| 568 | lsize_ = AttrVect_lsize(sMat%data) |
---|
| 569 | |
---|
| 570 | end function lsize_ |
---|
| 571 | |
---|
| 572 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 573 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 574 | !BOP ------------------------------------------------------------------- |
---|
| 575 | ! |
---|
| 576 | ! !IROUTINE: GlobalNumElements_ - Global Number of Non-zero Elements |
---|
| 577 | ! |
---|
| 578 | ! !DESCRIPTION: This routine computes the number of nonzero elements |
---|
| 579 | ! in a distributed {\tt SparseMatrix} variable {\tt sMat}. The input |
---|
| 580 | ! {\tt SparseMatrix} argument {\tt sMat} is examined on each process |
---|
| 581 | ! to determine the number of nonzero elements it holds, and this value |
---|
| 582 | ! is summed across the communicator associated with the input |
---|
| 583 | ! {\tt INTEGER} handle {\tt comm}, with the total returned {\em on each |
---|
| 584 | ! process on the communicator}. |
---|
| 585 | ! |
---|
| 586 | ! !INTERFACE: |
---|
| 587 | |
---|
| 588 | integer function GlobalNumElements_(sMat, comm) |
---|
| 589 | |
---|
| 590 | ! |
---|
| 591 | ! !USES: |
---|
| 592 | ! |
---|
| 593 | use m_die |
---|
| 594 | use m_mpif90 |
---|
| 595 | |
---|
| 596 | implicit none |
---|
| 597 | |
---|
| 598 | ! !INPUT PARAMETERS: |
---|
| 599 | |
---|
| 600 | type(SparseMatrix), intent(in) :: sMat |
---|
| 601 | integer, optional, intent(in) :: comm |
---|
| 602 | |
---|
| 603 | ! !REVISION HISTORY: |
---|
| 604 | ! 24Apr01 - Jay Larson <larson@mcs.anl.gov> - New routine. |
---|
| 605 | ! |
---|
| 606 | !EOP ___________________________________________________________________ |
---|
| 607 | ! |
---|
| 608 | character(len=*),parameter :: myname_=myname//':: GlobalNumElements_' |
---|
| 609 | |
---|
| 610 | integer :: MyNumElements, GNumElements, ierr |
---|
| 611 | |
---|
| 612 | ! Determine the number of locally held nonzero elements: |
---|
| 613 | |
---|
| 614 | MyNumElements = lsize_(sMat) |
---|
| 615 | |
---|
| 616 | call MPI_ALLREDUCE(MyNumElements, GNumElements, 1, MP_INTEGER, & |
---|
| 617 | MP_SUM, comm, ierr) |
---|
| 618 | if(ierr /= 0) then |
---|
| 619 | call MP_perr_die(myname_,"MPI_ALLREDUCE(MyNumElements...",ierr) |
---|
| 620 | endif |
---|
| 621 | |
---|
| 622 | GlobalNumElements_ = GNumElements |
---|
| 623 | |
---|
| 624 | end function GlobalNumElements_ |
---|
| 625 | |
---|
| 626 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 627 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 628 | !BOP ------------------------------------------------------------------- |
---|
| 629 | ! |
---|
| 630 | ! !IROUTINE: indexIA_ - Index an Integer Attribute |
---|
| 631 | ! |
---|
| 632 | ! !DESCRIPTION: This {\tt INTEGER} function reports the row index |
---|
| 633 | ! for a given {\tt INTEGER} attribute of the input {\tt SparseMatrix} |
---|
| 634 | ! argument {\tt sMat}. The attribute requested is represented by the |
---|
| 635 | ! input {\tt CHARACTER} variable {\tt attribute}. The list of integer |
---|
| 636 | ! attributes one can request is defined in the description block of the |
---|
| 637 | ! header of this module ({\tt m\_SparseMatrix}). |
---|
| 638 | ! |
---|
| 639 | ! Here is how {\tt indexIA\_} provides access to integer attribute data |
---|
| 640 | ! in a {\tt SparseMatrix} variable {\tt sMat}. Suppose we wish to access |
---|
| 641 | ! global row information. This attribute has associated with it the |
---|
| 642 | ! string tag {\tt grow}. The corresponding index returned ({\tt igrow}) |
---|
| 643 | ! is determined by invoking {\tt indexIA\_}: |
---|
| 644 | ! \begin{verbatim} |
---|
| 645 | ! igrow = indexIA_(sMat, 'grow') |
---|
| 646 | ! \end{verbatim} |
---|
| 647 | ! |
---|
| 648 | ! Access to the global row index data in {\tt sMat} is thus obtained by |
---|
| 649 | ! referencing {\tt sMat\%data\%iAttr(igrow,:)}. |
---|
| 650 | ! |
---|
| 651 | ! |
---|
| 652 | ! !INTERFACE: |
---|
| 653 | |
---|
| 654 | integer function indexIA_(sMat, item, perrWith, dieWith) |
---|
| 655 | ! |
---|
| 656 | ! !USES: |
---|
| 657 | ! |
---|
| 658 | use m_String, only : String |
---|
| 659 | use m_String, only : String_init => init |
---|
| 660 | use m_String, only : String_clean => clean |
---|
| 661 | use m_String, only : String_ToChar => ToChar |
---|
| 662 | |
---|
| 663 | use m_TraceBack, only : GenTraceBackString |
---|
| 664 | |
---|
| 665 | use m_AttrVect,only : AttrVect_indexIA => indexIA |
---|
| 666 | |
---|
| 667 | implicit none |
---|
| 668 | |
---|
| 669 | ! !INPUT PARAMETERS: |
---|
| 670 | |
---|
| 671 | type(SparseMatrix), intent(in) :: sMat |
---|
| 672 | character(len=*), intent(in) :: item |
---|
| 673 | character(len=*), optional, intent(in) :: perrWith |
---|
| 674 | character(len=*), optional, intent(in) :: dieWith |
---|
| 675 | |
---|
| 676 | ! !REVISION HISTORY: |
---|
| 677 | ! 23Apr00 - J.W. Larson <larson@mcs.anl.gov> - initial version. |
---|
| 678 | !EOP ___________________________________________________________________ |
---|
| 679 | |
---|
| 680 | character(len=*),parameter :: myname_=myname//'::indexIA_' |
---|
| 681 | type(String) :: myTrace |
---|
| 682 | |
---|
| 683 | ! Generate a traceback String |
---|
| 684 | |
---|
| 685 | if(present(dieWith)) then |
---|
| 686 | call GenTraceBackString(myTrace, dieWith, myname_) |
---|
| 687 | else |
---|
| 688 | if(present(perrWith)) then |
---|
| 689 | call GenTraceBackString(myTrace, perrWith, myname_) |
---|
| 690 | else |
---|
| 691 | call GenTraceBackString(myTrace, myname_) |
---|
| 692 | endif |
---|
| 693 | endif |
---|
| 694 | |
---|
| 695 | ! Call AttrVect_indexIA() accordingly: |
---|
| 696 | |
---|
| 697 | if( present(dieWith) .or. & |
---|
| 698 | ((.not. present(dieWith)) .and. (.not. present(perrWith))) ) then |
---|
| 699 | indexIA_ = AttrVect_indexIA(sMat%data, item, & |
---|
| 700 | dieWith=String_ToChar(myTrace)) |
---|
| 701 | else ! perrWith but no dieWith case |
---|
| 702 | indexIA_ = AttrVect_indexIA(sMat%data, item, & |
---|
| 703 | perrWith=String_ToChar(myTrace)) |
---|
| 704 | endif |
---|
| 705 | |
---|
| 706 | call String_clean(myTrace) |
---|
| 707 | |
---|
| 708 | end function indexIA_ |
---|
| 709 | |
---|
| 710 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 711 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 712 | !BOP ------------------------------------------------------------------- |
---|
| 713 | ! |
---|
| 714 | ! !IROUTINE: indexRA_ - Index a Real Attribute |
---|
| 715 | ! |
---|
| 716 | ! !DESCRIPTION: This {\tt INTEGER} function reports the row index |
---|
| 717 | ! for a given {\tt REAL} attribute of the input {\tt SparseMatrix} |
---|
| 718 | ! argument {\tt sMat}. The attribute requested is represented by the |
---|
| 719 | ! input {\tt CHARACTER} variable {\tt attribute}. The list of real |
---|
| 720 | ! attributes one can request is defined in the description block of the |
---|
| 721 | ! header of this module ({\tt m\_SparseMatrix}). |
---|
| 722 | ! |
---|
| 723 | ! Here is how {\tt indexRA\_} provides access to integer attribute data |
---|
| 724 | ! in a {\tt SparseMatrix} variable {\tt sMat}. Suppose we wish to access |
---|
| 725 | ! matrix element values. This attribute has associated with it the |
---|
| 726 | ! string tag {\tt weight}. The corresponding index returned ({\tt iweight}) |
---|
| 727 | ! is determined by invoking {\tt indexRA\_}: |
---|
| 728 | ! \begin{verbatim} |
---|
| 729 | ! iweight = indexRA_(sMat, 'weight') |
---|
| 730 | ! \end{verbatim} |
---|
| 731 | ! |
---|
| 732 | ! Access to the matrix element data in {\tt sMat} is thus obtained by |
---|
| 733 | ! referencing {\tt sMat\%data\%rAttr(iweight,:)}. |
---|
| 734 | ! |
---|
| 735 | ! !INTERFACE: |
---|
| 736 | |
---|
| 737 | integer function indexRA_(sMat, item, perrWith, dieWith) |
---|
| 738 | ! |
---|
| 739 | ! !USES: |
---|
| 740 | ! |
---|
| 741 | use m_String, only : String |
---|
| 742 | use m_String, only : String_init => init |
---|
| 743 | use m_String, only : String_clean => clean |
---|
| 744 | use m_String, only : String_ToChar => ToChar |
---|
| 745 | |
---|
| 746 | use m_TraceBack, only : GenTraceBackString |
---|
| 747 | |
---|
| 748 | use m_AttrVect,only : AttrVect_indexRA => indexRA |
---|
| 749 | |
---|
| 750 | implicit none |
---|
| 751 | |
---|
| 752 | ! !INPUT PARAMETERS: |
---|
| 753 | |
---|
| 754 | type(SparseMatrix), intent(in) :: sMat |
---|
| 755 | character(len=*), intent(in) :: item |
---|
| 756 | character(len=*), optional, intent(in) :: perrWith |
---|
| 757 | character(len=*), optional, intent(in) :: dieWith |
---|
| 758 | |
---|
| 759 | ! !REVISION HISTORY: |
---|
| 760 | ! 24Apr00 - J.W. Larson <larson@mcs.anl.gov> - initial version. |
---|
| 761 | !EOP ___________________________________________________________________ |
---|
| 762 | |
---|
| 763 | character(len=*),parameter :: myname_=myname//'::indexRA_' |
---|
| 764 | |
---|
| 765 | type(String) :: myTrace |
---|
| 766 | |
---|
| 767 | ! Generate a traceback String |
---|
| 768 | |
---|
| 769 | if(present(dieWith)) then ! append myname_ onto dieWith |
---|
| 770 | call GenTraceBackString(myTrace, dieWith, myname_) |
---|
| 771 | else |
---|
| 772 | if(present(perrWith)) then ! append myname_ onto perrwith |
---|
| 773 | call GenTraceBackString(myTrace, perrWith, myname_) |
---|
| 774 | else ! Start a TraceBack String |
---|
| 775 | call GenTraceBackString(myTrace, myname_) |
---|
| 776 | endif |
---|
| 777 | endif |
---|
| 778 | |
---|
| 779 | ! Call AttrVect_indexRA() accordingly: |
---|
| 780 | |
---|
| 781 | if( present(dieWith) .or. & |
---|
| 782 | ((.not. present(dieWith)) .and. (.not. present(perrWith))) ) then |
---|
| 783 | indexRA_ = AttrVect_indexRA(sMat%data, item, & |
---|
| 784 | dieWith=String_ToChar(myTrace)) |
---|
| 785 | else ! perrWith but no dieWith case |
---|
| 786 | indexRA_ = AttrVect_indexRA(sMat%data, item, & |
---|
| 787 | perrWith=String_ToChar(myTrace)) |
---|
| 788 | endif |
---|
| 789 | |
---|
| 790 | call String_clean(myTrace) |
---|
| 791 | |
---|
| 792 | end function indexRA_ |
---|
| 793 | |
---|
| 794 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 795 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 796 | !BOP ------------------------------------------------------------------- |
---|
| 797 | ! |
---|
| 798 | ! !IROUTINE: nRows_ - Return the Number of Rows |
---|
| 799 | ! |
---|
| 800 | ! !DESCRIPTION: This routine returns the {\em total} number of rows |
---|
| 801 | ! in the input {\tt SparseMatrix} argument {\tt sMat}. This number of |
---|
| 802 | ! rows is a constant, and not dependent on the decomposition of the |
---|
| 803 | ! {\tt SparseMatrix}. |
---|
| 804 | ! |
---|
| 805 | ! !INTERFACE: |
---|
| 806 | |
---|
| 807 | integer function nRows_(sMat) |
---|
| 808 | ! |
---|
| 809 | ! !USES: |
---|
| 810 | ! |
---|
| 811 | implicit none |
---|
| 812 | |
---|
| 813 | ! !INPUT PARAMETERS: |
---|
| 814 | |
---|
| 815 | type(SparseMatrix), intent(in) :: sMat |
---|
| 816 | |
---|
| 817 | ! !REVISION HISTORY: |
---|
| 818 | ! 19Apr01 - J.W. Larson <larson@mcs.anl.gov> - initial prototype |
---|
| 819 | !EOP ___________________________________________________________________ |
---|
| 820 | |
---|
| 821 | character(len=*),parameter :: myname_=myname//'::nRows_' |
---|
| 822 | |
---|
| 823 | nRows_ = sMat%nrows |
---|
| 824 | |
---|
| 825 | end function nRows_ |
---|
| 826 | |
---|
| 827 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 828 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 829 | !BOP ------------------------------------------------------------------- |
---|
| 830 | ! |
---|
| 831 | ! !IROUTINE: nCols_ - Return the Number of Columns |
---|
| 832 | ! |
---|
| 833 | ! !DESCRIPTION: This routine returns the {\em total} number of columns |
---|
| 834 | ! in the input {\tt SparseMatrix} argument {\tt sMat}. This number of |
---|
| 835 | ! columns is a constant, and not dependent on the decomposition of the |
---|
| 836 | ! {\tt SparseMatrix}. |
---|
| 837 | ! |
---|
| 838 | ! !INTERFACE: |
---|
| 839 | |
---|
| 840 | integer function nCols_(sMat) |
---|
| 841 | ! |
---|
| 842 | ! !USES: |
---|
| 843 | ! |
---|
| 844 | implicit none |
---|
| 845 | |
---|
| 846 | ! !INPUT PARAMETERS: |
---|
| 847 | |
---|
| 848 | type(SparseMatrix), intent(in) :: sMat |
---|
| 849 | |
---|
| 850 | ! !REVISION HISTORY: |
---|
| 851 | ! 19Apr01 - J.W. Larson <larson@mcs.anl.gov> - initial prototype |
---|
| 852 | !EOP ___________________________________________________________________ |
---|
| 853 | |
---|
| 854 | character(len=*),parameter :: myname_=myname//'::nCols_' |
---|
| 855 | |
---|
| 856 | nCols_ = sMat%ncols |
---|
| 857 | |
---|
| 858 | end function nCols_ |
---|
| 859 | |
---|
| 860 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 861 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 862 | !BOP ------------------------------------------------------------------- |
---|
| 863 | ! |
---|
| 864 | ! !IROUTINE: exportGlobalRowIndices_ - Return Global Row Indices |
---|
| 865 | ! |
---|
| 866 | ! !DESCRIPTION: |
---|
| 867 | ! This routine extracts from the input {\tt SparseMatrix} argument |
---|
| 868 | ! {\tt sMat} its global row indices, and returns them in the {\tt INTEGER} |
---|
| 869 | ! output array {\tt GlobalRows}, and its length in the output {\tt INTEGER} |
---|
| 870 | ! argument {\tt length}. |
---|
| 871 | ! |
---|
| 872 | ! {\bf N.B.:} The flexibility of this routine regarding the pointer |
---|
| 873 | ! association status of the output argument {\tt GlobalRows} means the |
---|
| 874 | ! user must invoke this routine with care. If the user wishes this |
---|
| 875 | ! routine to fill a pre-allocated array, then obviously this array |
---|
| 876 | ! must be allocated prior to calling this routine. If the user wishes |
---|
| 877 | ! that the routine {\em create} the output argument array {\tt GlobalRows}, |
---|
| 878 | ! then the user must ensure this pointer is not allocated (i.e. the user |
---|
| 879 | ! must nullify this pointer) at the time this routine is invoked. |
---|
| 880 | ! |
---|
| 881 | ! {\bf N.B.:} If the user has relied on this routine to allocate memory |
---|
| 882 | ! associated with the pointer {\tt GlobalRows}, then the user is responsible |
---|
| 883 | ! for deallocating this array once it is no longer needed. Failure to |
---|
| 884 | ! do so will result in a memory leak. |
---|
| 885 | ! |
---|
| 886 | ! !INTERFACE: |
---|
| 887 | |
---|
| 888 | subroutine exportGlobalRowIndices_(sMat, GlobalRows, length) |
---|
| 889 | ! |
---|
| 890 | ! !USES: |
---|
| 891 | ! |
---|
| 892 | use m_die |
---|
| 893 | use m_stdio |
---|
| 894 | |
---|
| 895 | use m_AttrVect, only : AttrVect_exportIAttr => exportIAttr |
---|
| 896 | |
---|
| 897 | implicit none |
---|
| 898 | |
---|
| 899 | ! !INPUT PARAMETERS: |
---|
| 900 | |
---|
| 901 | type(SparseMatrix), intent(in) :: sMat |
---|
| 902 | |
---|
| 903 | ! !OUTPUT PARAMETERS: |
---|
| 904 | |
---|
| 905 | integer, dimension(:), pointer :: GlobalRows |
---|
| 906 | integer, optional, intent(out) :: length |
---|
| 907 | |
---|
| 908 | ! !REVISION HISTORY: |
---|
| 909 | ! 7May02 - J.W. Larson <larson@mcs.anl.gov> - initial version. |
---|
| 910 | ! |
---|
| 911 | !EOP ___________________________________________________________________ |
---|
| 912 | |
---|
| 913 | character(len=*),parameter :: myname_=myname//'::exportGlobalRowIndices_' |
---|
| 914 | |
---|
| 915 | ! Export the data (inheritance from AttrVect) |
---|
| 916 | if(present(length)) then |
---|
| 917 | call AttrVect_exportIAttr(sMat%data, 'grow', GlobalRows, length) |
---|
| 918 | else |
---|
| 919 | call AttrVect_exportIAttr(sMat%data, 'grow', GlobalRows) |
---|
| 920 | endif |
---|
| 921 | |
---|
| 922 | end subroutine exportGlobalRowIndices_ |
---|
| 923 | |
---|
| 924 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 925 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 926 | !BOP ------------------------------------------------------------------- |
---|
| 927 | ! |
---|
| 928 | ! !IROUTINE: exportGlobalColumnIndices_ - Return Global Column Indices |
---|
| 929 | ! |
---|
| 930 | ! !DESCRIPTION: |
---|
| 931 | ! This routine extracts from the input {\tt SparseMatrix} argument |
---|
| 932 | ! {\tt sMat} its global column indices, and returns them in the {\tt INTEGER} |
---|
| 933 | ! output array {\tt GlobalColumns}, and its length in the output {\tt INTEGER} |
---|
| 934 | ! argument {\tt length}. |
---|
| 935 | ! |
---|
| 936 | ! {\bf N.B.:} The flexibility of this routine regarding the pointer |
---|
| 937 | ! association status of the output argument {\tt GlobalColumns} means the |
---|
| 938 | ! user must invoke this routine with care. If the user wishes this |
---|
| 939 | ! routine to fill a pre-allocated array, then obviously this array |
---|
| 940 | ! must be allocated prior to calling this routine. If the user wishes |
---|
| 941 | ! that the routine {\em create} the output argument array {\tt GlobalColumns}, |
---|
| 942 | ! then the user must ensure this pointer is not allocated (i.e. the user |
---|
| 943 | ! must nullify this pointer) at the time this routine is invoked. |
---|
| 944 | ! |
---|
| 945 | ! {\bf N.B.:} If the user has relied on this routine to allocate memory |
---|
| 946 | ! associated with the pointer {\tt GlobalColumns}, then the user is responsible |
---|
| 947 | ! for deallocating this array once it is no longer needed. Failure to |
---|
| 948 | ! do so will result in a memory leak. |
---|
| 949 | ! |
---|
| 950 | ! !INTERFACE: |
---|
| 951 | |
---|
| 952 | subroutine exportGlobalColumnIndices_(sMat, GlobalColumns, length) |
---|
| 953 | |
---|
| 954 | ! |
---|
| 955 | ! !USES: |
---|
| 956 | ! |
---|
| 957 | use m_die |
---|
| 958 | use m_stdio |
---|
| 959 | |
---|
| 960 | use m_AttrVect, only : AttrVect_exportIAttr => exportIAttr |
---|
| 961 | |
---|
| 962 | implicit none |
---|
| 963 | |
---|
| 964 | ! !INPUT PARAMETERS: |
---|
| 965 | |
---|
| 966 | type(SparseMatrix), intent(in) :: sMat |
---|
| 967 | |
---|
| 968 | ! !OUTPUT PARAMETERS: |
---|
| 969 | |
---|
| 970 | integer, dimension(:), pointer :: GlobalColumns |
---|
| 971 | integer, optional, intent(out) :: length |
---|
| 972 | |
---|
| 973 | ! !REVISION HISTORY: |
---|
| 974 | ! 7May02 - J.W. Larson <larson@mcs.anl.gov> - initial version. |
---|
| 975 | ! |
---|
| 976 | !EOP ___________________________________________________________________ |
---|
| 977 | |
---|
| 978 | character(len=*),parameter :: myname_=myname//'::exportGlobalColumnIndices_' |
---|
| 979 | |
---|
| 980 | ! Export the data (inheritance from AttrVect) |
---|
| 981 | if(present(length)) then |
---|
| 982 | call AttrVect_exportIAttr(sMat%data, 'gcol', GlobalColumns, length) |
---|
| 983 | else |
---|
| 984 | call AttrVect_exportIAttr(sMat%data, 'gcol', GlobalColumns) |
---|
| 985 | endif |
---|
| 986 | |
---|
| 987 | end subroutine exportGlobalColumnIndices_ |
---|
| 988 | |
---|
| 989 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 990 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 991 | !BOP ------------------------------------------------------------------- |
---|
| 992 | ! |
---|
| 993 | ! !IROUTINE: exportLocalRowIndices_ - Return Local Row Indices |
---|
| 994 | ! |
---|
| 995 | ! !DESCRIPTION: |
---|
| 996 | ! This routine extracts from the input {\tt SparseMatrix} argument |
---|
| 997 | ! {\tt sMat} its local row indices, and returns them in the {\tt INTEGER} |
---|
| 998 | ! output array {\tt LocalRows}, and its length in the output {\tt INTEGER} |
---|
| 999 | ! argument {\tt length}. |
---|
| 1000 | ! |
---|
| 1001 | ! {\bf N.B.:} The flexibility of this routine regarding the pointer |
---|
| 1002 | ! association status of the output argument {\tt LocalRows} means the |
---|
| 1003 | ! user must invoke this routine with care. If the user wishes this |
---|
| 1004 | ! routine to fill a pre-allocated array, then obviously this array |
---|
| 1005 | ! must be allocated prior to calling this routine. If the user wishes |
---|
| 1006 | ! that the routine {\em create} the output argument array {\tt LocalRows}, |
---|
| 1007 | ! then the user must ensure this pointer is not allocated (i.e. the user |
---|
| 1008 | ! must nullify this pointer) at the time this routine is invoked. |
---|
| 1009 | ! |
---|
| 1010 | ! {\bf N.B.:} If the user has relied on this routine to allocate memory |
---|
| 1011 | ! associated with the pointer {\tt LocalRows}, then the user is responsible |
---|
| 1012 | ! for deallocating this array once it is no longer needed. Failure to |
---|
| 1013 | ! do so will result in a memory leak. |
---|
| 1014 | ! |
---|
| 1015 | ! !INTERFACE: |
---|
| 1016 | |
---|
| 1017 | subroutine exportLocalRowIndices_(sMat, LocalRows, length) |
---|
| 1018 | ! |
---|
| 1019 | ! !USES: |
---|
| 1020 | ! |
---|
| 1021 | use m_die |
---|
| 1022 | use m_stdio |
---|
| 1023 | |
---|
| 1024 | use m_AttrVect, only : AttrVect_exportIAttr => exportIAttr |
---|
| 1025 | |
---|
| 1026 | implicit none |
---|
| 1027 | |
---|
| 1028 | ! !INPUT PARAMETERS: |
---|
| 1029 | |
---|
| 1030 | type(SparseMatrix), intent(in) :: sMat |
---|
| 1031 | |
---|
| 1032 | ! !OUTPUT PARAMETERS: |
---|
| 1033 | |
---|
| 1034 | integer, dimension(:), pointer :: LocalRows |
---|
| 1035 | integer, optional, intent(out) :: length |
---|
| 1036 | |
---|
| 1037 | ! !REVISION HISTORY: |
---|
| 1038 | ! 7May02 - J.W. Larson <larson@mcs.anl.gov> - initial version. |
---|
| 1039 | ! |
---|
| 1040 | !EOP ___________________________________________________________________ |
---|
| 1041 | |
---|
| 1042 | character(len=*),parameter :: myname_=myname//'::exportLocalRowIndices_' |
---|
| 1043 | |
---|
| 1044 | ! Export the data (inheritance from AttrVect) |
---|
| 1045 | if(present(length)) then |
---|
| 1046 | call AttrVect_exportIAttr(sMat%data, 'lrow', LocalRows, length) |
---|
| 1047 | else |
---|
| 1048 | call AttrVect_exportIAttr(sMat%data, 'lrow', LocalRows) |
---|
| 1049 | endif |
---|
| 1050 | |
---|
| 1051 | end subroutine exportLocalRowIndices_ |
---|
| 1052 | |
---|
| 1053 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1054 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1055 | !BOP ------------------------------------------------------------------- |
---|
| 1056 | ! |
---|
| 1057 | ! !IROUTINE: exportLocalColumnIndices_ - Return Local Column Indices |
---|
| 1058 | ! |
---|
| 1059 | ! !DESCRIPTION: |
---|
| 1060 | ! This routine extracts from the input {\tt SparseMatrix} argument |
---|
| 1061 | ! {\tt sMat} its local column indices, and returns them in the {\tt INTEGER} |
---|
| 1062 | ! output array {\tt LocalColumns}, and its length in the output {\tt INTEGER} |
---|
| 1063 | ! argument {\tt length}. |
---|
| 1064 | ! |
---|
| 1065 | ! {\bf N.B.:} The flexibility of this routine regarding the pointer |
---|
| 1066 | ! association status of the output argument {\tt LocalColumns} means the |
---|
| 1067 | ! user must invoke this routine with care. If the user wishes this |
---|
| 1068 | ! routine to fill a pre-allocated array, then obviously this array |
---|
| 1069 | ! must be allocated prior to calling this routine. If the user wishes |
---|
| 1070 | ! that the routine {\em create} the output argument array {\tt LocalColumns}, |
---|
| 1071 | ! then the user must ensure this pointer is not allocated (i.e. the user |
---|
| 1072 | ! must nullify this pointer) at the time this routine is invoked. |
---|
| 1073 | ! |
---|
| 1074 | ! {\bf N.B.:} If the user has relied on this routine to allocate memory |
---|
| 1075 | ! associated with the pointer {\tt LocalColumns}, then the user is responsible |
---|
| 1076 | ! for deallocating this array once it is no longer needed. Failure to |
---|
| 1077 | ! do so will result in a memory leak. |
---|
| 1078 | ! |
---|
| 1079 | ! !INTERFACE: |
---|
| 1080 | |
---|
| 1081 | subroutine exportLocalColumnIndices_(sMat, LocalColumns, length) |
---|
| 1082 | |
---|
| 1083 | ! |
---|
| 1084 | ! !USES: |
---|
| 1085 | ! |
---|
| 1086 | use m_die |
---|
| 1087 | use m_stdio |
---|
| 1088 | |
---|
| 1089 | use m_AttrVect, only : AttrVect_exportIAttr => exportIAttr |
---|
| 1090 | |
---|
| 1091 | implicit none |
---|
| 1092 | |
---|
| 1093 | ! !INPUT PARAMETERS: |
---|
| 1094 | |
---|
| 1095 | type(SparseMatrix), intent(in) :: sMat |
---|
| 1096 | |
---|
| 1097 | ! !OUTPUT PARAMETERS: |
---|
| 1098 | |
---|
| 1099 | integer, dimension(:), pointer :: LocalColumns |
---|
| 1100 | integer, optional, intent(out) :: length |
---|
| 1101 | |
---|
| 1102 | ! !REVISION HISTORY: |
---|
| 1103 | ! 7May02 - J.W. Larson <larson@mcs.anl.gov> - initial version. |
---|
| 1104 | ! |
---|
| 1105 | !EOP ___________________________________________________________________ |
---|
| 1106 | |
---|
| 1107 | character(len=*),parameter :: myname_=myname//'::exportLocalColumnIndices_' |
---|
| 1108 | |
---|
| 1109 | ! Export the data (inheritance from AttrVect) |
---|
| 1110 | if(present(length)) then |
---|
| 1111 | call AttrVect_exportIAttr(sMat%data, 'lcol', LocalColumns, length) |
---|
| 1112 | else |
---|
| 1113 | call AttrVect_exportIAttr(sMat%data, 'lcol', LocalColumns) |
---|
| 1114 | endif |
---|
| 1115 | |
---|
| 1116 | end subroutine exportLocalColumnIndices_ |
---|
| 1117 | |
---|
| 1118 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1119 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1120 | !BOP ------------------------------------------------------------------- |
---|
| 1121 | ! |
---|
| 1122 | ! !IROUTINE: exportMatrixElementsSP_ - Return Matrix Elements as Array |
---|
| 1123 | ! |
---|
| 1124 | ! !DESCRIPTION: |
---|
| 1125 | ! This routine extracts the matrix elements from the input {\tt SparseMatrix} |
---|
| 1126 | ! argument {\tt sMat}, and returns them in the {\tt REAL} output array |
---|
| 1127 | ! {\tt MatrixElements}, and its length in the output {\tt INTEGER} |
---|
| 1128 | ! argument {\tt length}. |
---|
| 1129 | ! |
---|
| 1130 | ! {\bf N.B.:} The flexibility of this routine regarding the pointer |
---|
| 1131 | ! association status of the output argument {\tt MatrixElements} means the |
---|
| 1132 | ! user must invoke this routine with care. If the user wishes this |
---|
| 1133 | ! routine to fill a pre-allocated array, then obviously this array |
---|
| 1134 | ! must be allocated prior to calling this routine. If the user wishes |
---|
| 1135 | ! that the routine {\em create} the output argument array {\tt MatrixElements}, |
---|
| 1136 | ! then the user must ensure this pointer is not allocated (i.e. the user |
---|
| 1137 | ! must nullify this pointer) at the time this routine is invoked. |
---|
| 1138 | ! |
---|
| 1139 | ! {\bf N.B.:} If the user has relied on this routine to allocate memory |
---|
| 1140 | ! associated with the pointer {\tt MatrixElements}, then the user is responsible |
---|
| 1141 | ! for deallocating this array once it is no longer needed. Failure to |
---|
| 1142 | ! do so will result in a memory leak. |
---|
| 1143 | ! |
---|
| 1144 | ! The native precision version is described here. A double precision version |
---|
| 1145 | ! is also available. |
---|
| 1146 | ! |
---|
| 1147 | ! !INTERFACE: |
---|
| 1148 | |
---|
| 1149 | subroutine exportMatrixelementsSP_(sMat, MatrixElements, length) |
---|
| 1150 | |
---|
| 1151 | ! |
---|
| 1152 | ! !USES: |
---|
| 1153 | ! |
---|
| 1154 | use m_die |
---|
| 1155 | use m_stdio |
---|
| 1156 | use m_realkinds, only : SP |
---|
| 1157 | |
---|
| 1158 | use m_AttrVect, only : AttrVect_exportRAttr => exportRAttr |
---|
| 1159 | |
---|
| 1160 | implicit none |
---|
| 1161 | |
---|
| 1162 | ! !INPUT PARAMETERS: |
---|
| 1163 | |
---|
| 1164 | type(SparseMatrix), intent(in) :: sMat |
---|
| 1165 | |
---|
| 1166 | ! !OUTPUT PARAMETERS: |
---|
| 1167 | |
---|
| 1168 | real(SP), dimension(:), pointer :: MatrixElements |
---|
| 1169 | integer, optional, intent(out) :: length |
---|
| 1170 | |
---|
| 1171 | ! !REVISION HISTORY: |
---|
| 1172 | ! 7May02 - J.W. Larson <larson@mcs.anl.gov> - initial version. |
---|
| 1173 | ! 6Jan04 - R. Jacob <jacob@mcs.anl.gov> - SP and DP versions |
---|
| 1174 | ! |
---|
| 1175 | !EOP ___________________________________________________________________ |
---|
| 1176 | |
---|
| 1177 | character(len=*),parameter :: myname_=myname//'::exportMatrixElementsSP_' |
---|
| 1178 | |
---|
| 1179 | ! Export the data (inheritance from AttrVect) |
---|
| 1180 | if(present(length)) then |
---|
| 1181 | call AttrVect_exportRAttr(sMat%data, 'weight', MatrixElements, length) |
---|
| 1182 | else |
---|
| 1183 | call AttrVect_exportRAttr(sMat%data, 'weight', MatrixElements) |
---|
| 1184 | endif |
---|
| 1185 | |
---|
| 1186 | end subroutine exportMatrixElementsSP_ |
---|
| 1187 | |
---|
| 1188 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1189 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1190 | ! ------------------------------------------------------------------- |
---|
| 1191 | ! |
---|
| 1192 | ! !IROUTINE: exportMatrixElementsDP_ - Return Matrix Elements as Array |
---|
| 1193 | ! |
---|
| 1194 | ! !DESCRIPTION: |
---|
| 1195 | ! Double precision version of exportMatrixElementsSP_ |
---|
| 1196 | ! |
---|
| 1197 | ! !INTERFACE: |
---|
| 1198 | |
---|
| 1199 | subroutine exportMatrixelementsDP_(sMat, MatrixElements, length) |
---|
| 1200 | |
---|
| 1201 | ! |
---|
| 1202 | ! !USES: |
---|
| 1203 | ! |
---|
| 1204 | use m_die |
---|
| 1205 | use m_stdio |
---|
| 1206 | use m_realkinds, only : DP |
---|
| 1207 | |
---|
| 1208 | use m_AttrVect, only : AttrVect_exportRAttr => exportRAttr |
---|
| 1209 | |
---|
| 1210 | implicit none |
---|
| 1211 | |
---|
| 1212 | ! !INPUT PARAMETERS: |
---|
| 1213 | |
---|
| 1214 | type(SparseMatrix), intent(in) :: sMat |
---|
| 1215 | |
---|
| 1216 | ! !OUTPUT PARAMETERS: |
---|
| 1217 | |
---|
| 1218 | real(DP), dimension(:), pointer :: MatrixElements |
---|
| 1219 | integer, optional, intent(out) :: length |
---|
| 1220 | |
---|
| 1221 | ! !REVISION HISTORY: |
---|
| 1222 | ! 7May02 - J.W. Larson <larson@mcs.anl.gov> - initial version. |
---|
| 1223 | ! |
---|
| 1224 | ! ___________________________________________________________________ |
---|
| 1225 | |
---|
| 1226 | character(len=*),parameter :: myname_=myname//'::exportMatrixElementsDP_' |
---|
| 1227 | |
---|
| 1228 | ! Export the data (inheritance from AttrVect) |
---|
| 1229 | if(present(length)) then |
---|
| 1230 | call AttrVect_exportRAttr(sMat%data, 'weight', MatrixElements, length) |
---|
| 1231 | else |
---|
| 1232 | call AttrVect_exportRAttr(sMat%data, 'weight', MatrixElements) |
---|
| 1233 | endif |
---|
| 1234 | |
---|
| 1235 | end subroutine exportMatrixElementsDP_ |
---|
| 1236 | |
---|
| 1237 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1238 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1239 | !BOP ------------------------------------------------------------------- |
---|
| 1240 | ! |
---|
| 1241 | ! !IROUTINE: importGlobalRowIndices_ - Set Global Row Indices of Elements |
---|
| 1242 | ! |
---|
| 1243 | ! !DESCRIPTION: |
---|
| 1244 | ! This routine imports global row index data into the {\tt SparseMatrix} |
---|
| 1245 | ! argument {\tt sMat}. The user provides the index data in the input |
---|
| 1246 | ! {\tt INTEGER} vector {\tt inVect}. The input {\tt INTEGER} argument |
---|
| 1247 | ! {\tt lsize} is used as a consistencey check to ensure the user is |
---|
| 1248 | ! sufficient space in the {\tt SparseMatrix} to store the data. |
---|
| 1249 | ! |
---|
| 1250 | ! !INTERFACE: |
---|
| 1251 | |
---|
| 1252 | subroutine importGlobalRowIndices_(sMat, inVect, lsize) |
---|
| 1253 | |
---|
| 1254 | ! |
---|
| 1255 | ! !USES: |
---|
| 1256 | ! |
---|
| 1257 | use m_die |
---|
| 1258 | use m_stdio |
---|
| 1259 | |
---|
| 1260 | use m_AttrVect, only : AttrVect_importIAttr => importIAttr |
---|
| 1261 | |
---|
| 1262 | implicit none |
---|
| 1263 | |
---|
| 1264 | ! !INPUT PARAMETERS: |
---|
| 1265 | |
---|
| 1266 | integer, dimension(:), pointer :: inVect |
---|
| 1267 | integer, intent(in) :: lsize |
---|
| 1268 | |
---|
| 1269 | ! !INPUT/OUTPUT PARAMETERS: |
---|
| 1270 | |
---|
| 1271 | type(SparseMatrix), intent(inout) :: sMat |
---|
| 1272 | |
---|
| 1273 | ! !REVISION HISTORY: |
---|
| 1274 | ! 7May02 - J.W. Larson <larson@mcs.anl.gov> - initial prototype. |
---|
| 1275 | ! |
---|
| 1276 | !EOP ___________________________________________________________________ |
---|
| 1277 | |
---|
| 1278 | character(len=*),parameter :: myname_=myname//'::importGlobalRowIndices_' |
---|
| 1279 | |
---|
| 1280 | ! Argument Check: |
---|
| 1281 | |
---|
| 1282 | if(lsize > lsize_(sMat)) then |
---|
| 1283 | write(stderr,*) myname_,':: ERROR, lsize > lsize_(sMat).', & |
---|
| 1284 | 'lsize = ',lsize,'lsize_(sMat) = ',lsize_(sMat) |
---|
| 1285 | call die(myname_) |
---|
| 1286 | endif |
---|
| 1287 | |
---|
| 1288 | ! Import the data (inheritance from AttrVect) |
---|
| 1289 | |
---|
| 1290 | call AttrVect_importIAttr(sMat%data, 'grow', inVect, lsize) |
---|
| 1291 | |
---|
| 1292 | end subroutine importGlobalRowIndices_ |
---|
| 1293 | |
---|
| 1294 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1295 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1296 | !BOP ------------------------------------------------------------------- |
---|
| 1297 | ! |
---|
| 1298 | ! !IROUTINE: importGlobalColumnIndices_ - Set Global Column Indices of Elements |
---|
| 1299 | ! |
---|
| 1300 | ! !DESCRIPTION: |
---|
| 1301 | ! This routine imports global column index data into the {\tt SparseMatrix} |
---|
| 1302 | ! argument {\tt sMat}. The user provides the index data in the input |
---|
| 1303 | ! {\tt INTEGER} vector {\tt inVect}. The input {\tt INTEGER} argument |
---|
| 1304 | ! {\tt lsize} is used as a consistencey check to ensure the user is |
---|
| 1305 | ! sufficient space in the {\tt SparseMatrix} to store the data. |
---|
| 1306 | ! |
---|
| 1307 | ! !INTERFACE: |
---|
| 1308 | |
---|
| 1309 | subroutine importGlobalColumnIndices_(sMat, inVect, lsize) |
---|
| 1310 | |
---|
| 1311 | ! |
---|
| 1312 | ! !USES: |
---|
| 1313 | ! |
---|
| 1314 | use m_die |
---|
| 1315 | use m_stdio |
---|
| 1316 | |
---|
| 1317 | use m_AttrVect, only : AttrVect_importIAttr => importIAttr |
---|
| 1318 | |
---|
| 1319 | implicit none |
---|
| 1320 | |
---|
| 1321 | ! !INPUT PARAMETERS: |
---|
| 1322 | |
---|
| 1323 | integer, dimension(:), pointer :: inVect |
---|
| 1324 | integer, intent(in) :: lsize |
---|
| 1325 | |
---|
| 1326 | ! !INPUT/OUTPUT PARAMETERS: |
---|
| 1327 | |
---|
| 1328 | type(SparseMatrix), intent(inout) :: sMat |
---|
| 1329 | |
---|
| 1330 | ! !REVISION HISTORY: |
---|
| 1331 | ! 7May02 - J.W. Larson <larson@mcs.anl.gov> - initial prototype. |
---|
| 1332 | ! |
---|
| 1333 | !EOP ___________________________________________________________________ |
---|
| 1334 | |
---|
| 1335 | character(len=*),parameter :: myname_=myname//'::importGlobalColumnIndices_' |
---|
| 1336 | |
---|
| 1337 | ! Argument Check: |
---|
| 1338 | |
---|
| 1339 | if(lsize > lsize_(sMat)) then |
---|
| 1340 | write(stderr,*) myname_,':: ERROR, lsize > lsize_(sMat).', & |
---|
| 1341 | 'lsize = ',lsize,'lsize_(sMat) = ',lsize_(sMat) |
---|
| 1342 | call die(myname_) |
---|
| 1343 | endif |
---|
| 1344 | |
---|
| 1345 | ! Import the data (inheritance from AttrVect) |
---|
| 1346 | |
---|
| 1347 | call AttrVect_importIAttr(sMat%data, 'gcol', inVect, lsize) |
---|
| 1348 | |
---|
| 1349 | end subroutine importGlobalColumnIndices_ |
---|
| 1350 | |
---|
| 1351 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1352 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1353 | !BOP ------------------------------------------------------------------- |
---|
| 1354 | ! |
---|
| 1355 | ! !IROUTINE: importLocalRowIndices_ - Set Local Row Indices of Elements |
---|
| 1356 | ! |
---|
| 1357 | ! !DESCRIPTION: |
---|
| 1358 | ! This routine imports local row index data into the {\tt SparseMatrix} |
---|
| 1359 | ! argument {\tt sMat}. The user provides the index data in the input |
---|
| 1360 | ! {\tt INTEGER} vector {\tt inVect}. The input {\tt INTEGER} argument |
---|
| 1361 | ! {\tt lsize} is used as a consistencey check to ensure the user is |
---|
| 1362 | ! sufficient space in the {\tt SparseMatrix} to store the data. |
---|
| 1363 | ! |
---|
| 1364 | ! !INTERFACE: |
---|
| 1365 | |
---|
| 1366 | subroutine importLocalRowIndices_(sMat, inVect, lsize) |
---|
| 1367 | |
---|
| 1368 | ! |
---|
| 1369 | ! !USES: |
---|
| 1370 | ! |
---|
| 1371 | use m_die |
---|
| 1372 | use m_stdio |
---|
| 1373 | |
---|
| 1374 | use m_AttrVect, only : AttrVect_importIAttr => importIAttr |
---|
| 1375 | |
---|
| 1376 | implicit none |
---|
| 1377 | |
---|
| 1378 | ! !INPUT PARAMETERS: |
---|
| 1379 | |
---|
| 1380 | integer, dimension(:), pointer :: inVect |
---|
| 1381 | integer, intent(in) :: lsize |
---|
| 1382 | |
---|
| 1383 | ! !INPUT/OUTPUT PARAMETERS: |
---|
| 1384 | |
---|
| 1385 | type(SparseMatrix), intent(inout) :: sMat |
---|
| 1386 | |
---|
| 1387 | ! !REVISION HISTORY: |
---|
| 1388 | ! 7May02 - J.W. Larson <larson@mcs.anl.gov> - initial prototype. |
---|
| 1389 | ! |
---|
| 1390 | !EOP ___________________________________________________________________ |
---|
| 1391 | |
---|
| 1392 | character(len=*),parameter :: myname_=myname//'::importLocalRowIndices_' |
---|
| 1393 | |
---|
| 1394 | ! Argument Check: |
---|
| 1395 | |
---|
| 1396 | if(lsize > lsize_(sMat)) then |
---|
| 1397 | write(stderr,*) myname_,':: ERROR, lsize > lsize_(sMat).', & |
---|
| 1398 | 'lsize = ',lsize,'lsize_(sMat) = ',lsize_(sMat) |
---|
| 1399 | call die(myname_) |
---|
| 1400 | endif |
---|
| 1401 | |
---|
| 1402 | ! Import the data (inheritance from AttrVect) |
---|
| 1403 | |
---|
| 1404 | call AttrVect_importIAttr(sMat%data, 'lrow', inVect, lsize) |
---|
| 1405 | |
---|
| 1406 | end subroutine importLocalRowIndices_ |
---|
| 1407 | |
---|
| 1408 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1409 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1410 | !BOP ------------------------------------------------------------------- |
---|
| 1411 | ! |
---|
| 1412 | ! !IROUTINE: importLocalColumnIndices_ - Set Local Column Indices of Elements |
---|
| 1413 | ! |
---|
| 1414 | ! !DESCRIPTION: |
---|
| 1415 | ! This routine imports local column index data into the {\tt SparseMatrix} |
---|
| 1416 | ! argument {\tt sMat}. The user provides the index data in the input |
---|
| 1417 | ! {\tt INTEGER} vector {\tt inVect}. The input {\tt INTEGER} argument |
---|
| 1418 | ! {\tt lsize} is used as a consistencey check to ensure the user is |
---|
| 1419 | ! sufficient space in the {\tt SparseMatrix} to store the data. |
---|
| 1420 | ! |
---|
| 1421 | ! !INTERFACE: |
---|
| 1422 | |
---|
| 1423 | subroutine importLocalColumnIndices_(sMat, inVect, lsize) |
---|
| 1424 | |
---|
| 1425 | ! |
---|
| 1426 | ! !USES: |
---|
| 1427 | ! |
---|
| 1428 | use m_die |
---|
| 1429 | use m_stdio |
---|
| 1430 | |
---|
| 1431 | use m_AttrVect, only : AttrVect_importIAttr => importIAttr |
---|
| 1432 | |
---|
| 1433 | implicit none |
---|
| 1434 | |
---|
| 1435 | ! !INPUT PARAMETERS: |
---|
| 1436 | |
---|
| 1437 | integer, dimension(:), pointer :: inVect |
---|
| 1438 | integer, intent(in) :: lsize |
---|
| 1439 | |
---|
| 1440 | ! !INPUT/OUTPUT PARAMETERS: |
---|
| 1441 | |
---|
| 1442 | type(SparseMatrix), intent(inout) :: sMat |
---|
| 1443 | |
---|
| 1444 | ! !REVISION HISTORY: |
---|
| 1445 | ! 7May02 - J.W. Larson <larson@mcs.anl.gov> - initial prototype. |
---|
| 1446 | ! |
---|
| 1447 | !EOP ___________________________________________________________________ |
---|
| 1448 | |
---|
| 1449 | character(len=*),parameter :: myname_=myname//'::importLocalColumnIndices_' |
---|
| 1450 | |
---|
| 1451 | ! Argument Check: |
---|
| 1452 | |
---|
| 1453 | if(lsize > lsize_(sMat)) then |
---|
| 1454 | write(stderr,*) myname_,':: ERROR, lsize > lsize_(sMat).', & |
---|
| 1455 | 'lsize = ',lsize,'lsize_(sMat) = ',lsize_(sMat) |
---|
| 1456 | call die(myname_) |
---|
| 1457 | endif |
---|
| 1458 | |
---|
| 1459 | ! Import the data (inheritance from AttrVect) |
---|
| 1460 | |
---|
| 1461 | call AttrVect_importIAttr(sMat%data, 'lcol', inVect, lsize) |
---|
| 1462 | |
---|
| 1463 | end subroutine importLocalColumnIndices_ |
---|
| 1464 | |
---|
| 1465 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1466 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1467 | !BOP ------------------------------------------------------------------- |
---|
| 1468 | ! |
---|
| 1469 | ! !IROUTINE: importMatrixElementsSP_ - Import Non-zero Matrix Elements |
---|
| 1470 | ! |
---|
| 1471 | ! !DESCRIPTION: |
---|
| 1472 | ! This routine imports matrix elements index data into the |
---|
| 1473 | ! {\tt SparseMatrix} argument {\tt sMat}. The user provides the index |
---|
| 1474 | ! data in the input {\tt REAL} vector {\tt inVect}. The input |
---|
| 1475 | ! {\tt INTEGER} argument {\tt lsize} is used as a consistencey check |
---|
| 1476 | ! to ensure the user is sufficient space in the {\tt SparseMatrix} |
---|
| 1477 | ! to store the data. |
---|
| 1478 | ! |
---|
| 1479 | ! !INTERFACE: |
---|
| 1480 | |
---|
| 1481 | subroutine importMatrixElementsSP_(sMat, inVect, lsize) |
---|
| 1482 | |
---|
| 1483 | ! |
---|
| 1484 | ! !USES: |
---|
| 1485 | ! |
---|
| 1486 | use m_die |
---|
| 1487 | use m_stdio |
---|
| 1488 | use m_realkinds, only : SP |
---|
| 1489 | |
---|
| 1490 | use m_AttrVect, only : AttrVect_importRAttr => importRAttr |
---|
| 1491 | |
---|
| 1492 | implicit none |
---|
| 1493 | |
---|
| 1494 | ! !INPUT PARAMETERS: |
---|
| 1495 | |
---|
| 1496 | real(SP), dimension(:), pointer :: inVect |
---|
| 1497 | integer, intent(in) :: lsize |
---|
| 1498 | |
---|
| 1499 | ! !INPUT/OUTPUT PARAMETERS: |
---|
| 1500 | |
---|
| 1501 | type(SparseMatrix), intent(inout) :: sMat |
---|
| 1502 | |
---|
| 1503 | ! !REVISION HISTORY: |
---|
| 1504 | ! 7May02 - J.W. Larson <larson@mcs.anl.gov> - initial prototype. |
---|
| 1505 | ! 6Jan04 - R. Jacob <jacob@mcs.anl.gov> - Make SP and DP versions. |
---|
| 1506 | ! |
---|
| 1507 | !EOP ___________________________________________________________________ |
---|
| 1508 | |
---|
| 1509 | character(len=*),parameter :: myname_=myname//'::importMatrixElementsSP_' |
---|
| 1510 | |
---|
| 1511 | ! Argument Check: |
---|
| 1512 | |
---|
| 1513 | if(lsize > lsize_(sMat)) then |
---|
| 1514 | write(stderr,*) myname_,':: ERROR, lsize > lsize_(sMat).', & |
---|
| 1515 | 'lsize = ',lsize,'lsize_(sMat) = ',lsize_(sMat) |
---|
| 1516 | call die(myname_) |
---|
| 1517 | endif |
---|
| 1518 | |
---|
| 1519 | ! Import the data (inheritance from AttrVect) |
---|
| 1520 | |
---|
| 1521 | call AttrVect_importRAttr(sMat%data, 'weight', inVect, lsize) |
---|
| 1522 | |
---|
| 1523 | end subroutine importMatrixElementsSP_ |
---|
| 1524 | |
---|
| 1525 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1526 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1527 | ! ------------------------------------------------------------------- |
---|
| 1528 | ! |
---|
| 1529 | ! !IROUTINE: importMatrixElementsDP_ - Import Non-zero Matrix Elements |
---|
| 1530 | ! |
---|
| 1531 | ! !DESCRIPTION: |
---|
| 1532 | ! Double precision version of importMatrixElementsSP_ |
---|
| 1533 | ! |
---|
| 1534 | ! !INTERFACE: |
---|
| 1535 | |
---|
| 1536 | subroutine importMatrixElementsDP_(sMat, inVect, lsize) |
---|
| 1537 | |
---|
| 1538 | ! |
---|
| 1539 | ! !USES: |
---|
| 1540 | ! |
---|
| 1541 | use m_die |
---|
| 1542 | use m_stdio |
---|
| 1543 | use m_realkinds, only : DP |
---|
| 1544 | |
---|
| 1545 | use m_AttrVect, only : AttrVect_importRAttr => importRAttr |
---|
| 1546 | |
---|
| 1547 | implicit none |
---|
| 1548 | |
---|
| 1549 | ! !INPUT PARAMETERS: |
---|
| 1550 | |
---|
| 1551 | real(DP), dimension(:), pointer :: inVect |
---|
| 1552 | integer, intent(in) :: lsize |
---|
| 1553 | |
---|
| 1554 | ! !INPUT/OUTPUT PARAMETERS: |
---|
| 1555 | |
---|
| 1556 | type(SparseMatrix), intent(inout) :: sMat |
---|
| 1557 | |
---|
| 1558 | ! !REVISION HISTORY: |
---|
| 1559 | ! 7May02 - J.W. Larson <larson@mcs.anl.gov> - initial prototype. |
---|
| 1560 | ! |
---|
| 1561 | ! ___________________________________________________________________ |
---|
| 1562 | |
---|
| 1563 | character(len=*),parameter :: myname_=myname//'::importMatrixElementsDP_' |
---|
| 1564 | |
---|
| 1565 | ! Argument Check: |
---|
| 1566 | |
---|
| 1567 | if(lsize > lsize_(sMat)) then |
---|
| 1568 | write(stderr,*) myname_,':: ERROR, lsize > lsize_(sMat).', & |
---|
| 1569 | 'lsize = ',lsize,'lsize_(sMat) = ',lsize_(sMat) |
---|
| 1570 | call die(myname_) |
---|
| 1571 | endif |
---|
| 1572 | |
---|
| 1573 | ! Import the data (inheritance from AttrVect) |
---|
| 1574 | |
---|
| 1575 | call AttrVect_importRAttr(sMat%data, 'weight', inVect, lsize) |
---|
| 1576 | |
---|
| 1577 | end subroutine importMatrixElementsDP_ |
---|
| 1578 | |
---|
| 1579 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1580 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1581 | !BOP ------------------------------------------------------------------- |
---|
| 1582 | ! |
---|
| 1583 | ! !IROUTINE: Copy_ - Create a Copy of an Input SparseMatrix |
---|
| 1584 | ! |
---|
| 1585 | ! !DESCRIPTION: |
---|
| 1586 | ! This routine creates a copy of the input {\tt SparseMatrix} argument |
---|
| 1587 | ! {\tt sMat}, returning it as the output {\tt SparseMatrix} argument |
---|
| 1588 | ! {\tt sMatCopy}. |
---|
| 1589 | ! |
---|
| 1590 | ! {\bf N.B.:} The output argument {\tt sMatCopy} represents allocated |
---|
| 1591 | ! memory the user must deallocate when it is no longer needed. The |
---|
| 1592 | ! MCT routine to use for this purpose is {\tt clean()} from this module. |
---|
| 1593 | ! |
---|
| 1594 | ! !INTERFACE: |
---|
| 1595 | |
---|
| 1596 | subroutine Copy_(sMat, sMatCopy) |
---|
| 1597 | |
---|
| 1598 | ! |
---|
| 1599 | ! !USES: |
---|
| 1600 | ! |
---|
| 1601 | use m_die |
---|
| 1602 | use m_stdio |
---|
| 1603 | |
---|
| 1604 | use m_AttrVect, only : AttrVect |
---|
| 1605 | use m_AttrVect, only : AttrVect_init => init |
---|
| 1606 | use m_AttrVect, only : AttrVect_lsize => lsize |
---|
| 1607 | use m_AttrVect, only : AttrVect_Copy => Copy |
---|
| 1608 | |
---|
| 1609 | implicit none |
---|
| 1610 | |
---|
| 1611 | ! !INPUT PARAMETERS: |
---|
| 1612 | |
---|
| 1613 | type(SparseMatrix), intent(in) :: sMat |
---|
| 1614 | |
---|
| 1615 | ! !OUTPUT PARAMETERS: |
---|
| 1616 | |
---|
| 1617 | type(SparseMatrix), intent(out) :: sMatCopy |
---|
| 1618 | |
---|
| 1619 | ! !REVISION HISTORY: |
---|
| 1620 | ! 27Sep02 - J.W. Larson <larson@mcs.anl.gov> - initial prototype. |
---|
| 1621 | ! |
---|
| 1622 | !EOP ___________________________________________________________________ |
---|
| 1623 | |
---|
| 1624 | character(len=*),parameter :: myname_=myname//'::Copy_' |
---|
| 1625 | |
---|
| 1626 | ! Step one: copy the integer components of sMat: |
---|
| 1627 | |
---|
| 1628 | sMatCopy%nrows = sMat%nrows |
---|
| 1629 | sMatCopy%ncols = sMat%ncols |
---|
| 1630 | |
---|
| 1631 | sMatCopy%vecinit = .FALSE. |
---|
| 1632 | |
---|
| 1633 | ! Step two: Initialize the AttrVect sMatCopy%data off of sMat: |
---|
| 1634 | |
---|
| 1635 | call AttrVect_init(sMatCopy%data, sMat%data, AttrVect_lsize(sMat%data)) |
---|
| 1636 | |
---|
| 1637 | ! Step three: Copy sMat%data to sMatCopy%data: |
---|
| 1638 | |
---|
| 1639 | call AttrVect_Copy(sMat%data, aVout=sMatCopy%data) |
---|
| 1640 | |
---|
| 1641 | if(sMat%vecinit) call vecinit_(sMatCopy) |
---|
| 1642 | |
---|
| 1643 | end subroutine Copy_ |
---|
| 1644 | |
---|
| 1645 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1646 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1647 | !BOP ------------------------------------------------------------------- |
---|
| 1648 | ! |
---|
| 1649 | ! !IROUTINE: local_row_range_ - Local Row Extent of Non-zero Elements |
---|
| 1650 | ! |
---|
| 1651 | ! !DESCRIPTION: This routine examines the input distributed |
---|
| 1652 | ! {\tt SparseMatrix} variable {\tt sMat}, and returns the range of local |
---|
| 1653 | ! row values having nonzero elements. The first local row with |
---|
| 1654 | ! nonzero elements is returned in the {\tt INTEGER} argument |
---|
| 1655 | ! {\tt start\_row}, the last row in {\tt end\_row}. |
---|
| 1656 | ! |
---|
| 1657 | ! !INTERFACE: |
---|
| 1658 | |
---|
| 1659 | subroutine local_row_range_(sMat, start_row, end_row) |
---|
| 1660 | ! |
---|
| 1661 | ! !USES: |
---|
| 1662 | ! |
---|
| 1663 | use m_die |
---|
| 1664 | |
---|
| 1665 | use m_AttrVect, only : AttrVect_lsize => lsize |
---|
| 1666 | use m_AttrVect, only : AttrVect_indexIA => indexIA |
---|
| 1667 | |
---|
| 1668 | implicit none |
---|
| 1669 | |
---|
| 1670 | ! !INPUT PARAMETERS: |
---|
| 1671 | |
---|
| 1672 | type(SparseMatrix), intent(in) :: sMat |
---|
| 1673 | |
---|
| 1674 | ! !OUTPUT PARAMETERS: |
---|
| 1675 | |
---|
| 1676 | integer, intent(out) :: start_row |
---|
| 1677 | integer, intent(out) :: end_row |
---|
| 1678 | |
---|
| 1679 | ! !REVISION HISTORY: |
---|
| 1680 | ! 15Jan01 - Jay Larson <larson@mcs.anl.gov> - API specification. |
---|
| 1681 | ! 25Feb01 - Jay Larson <larson@mcs.anl.gov> - Initial prototype. |
---|
| 1682 | ! 23Apr01 - Jay Larson <larson@mcs.anl.gov> - Modified to accomodate |
---|
| 1683 | ! changes to the SparseMatrix type. |
---|
| 1684 | !EOP ___________________________________________________________________ |
---|
| 1685 | ! |
---|
| 1686 | character(len=*),parameter :: myname_=myname//'::local_row_range_' |
---|
| 1687 | |
---|
| 1688 | integer :: i, ilrow, lsize |
---|
| 1689 | |
---|
| 1690 | ilrow = AttrVect_indexIA(sMat%data, 'lrow') |
---|
| 1691 | lsize = AttrVect_lsize(sMat%data) |
---|
| 1692 | |
---|
| 1693 | ! Initialize start_row and end_row: |
---|
| 1694 | |
---|
| 1695 | start_row = sMat%data%iAttr(ilrow,1) |
---|
| 1696 | end_row = sMat%data%iAttr(ilrow,1) |
---|
| 1697 | |
---|
| 1698 | do i=1,lsize |
---|
| 1699 | start_row = min(start_row, sMat%data%iAttr(ilrow,i)) |
---|
| 1700 | end_row = max(end_row, sMat%data%iAttr(ilrow,i)) |
---|
| 1701 | end do |
---|
| 1702 | |
---|
| 1703 | end subroutine local_row_range_ |
---|
| 1704 | |
---|
| 1705 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1706 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1707 | !BOP ------------------------------------------------------------------- |
---|
| 1708 | ! |
---|
| 1709 | ! !IROUTINE: global_row_range_ - Global Row Extent of Non-zero Elements |
---|
| 1710 | ! |
---|
| 1711 | ! !DESCRIPTION: This routine examines the input distributed |
---|
| 1712 | ! {\tt SparseMatrix} variable {\tt sMat}, and returns the range of |
---|
| 1713 | ! global row values having nonzero elements. The first local row with |
---|
| 1714 | ! nonzero elements is returned in the {\tt INTEGER} argument |
---|
| 1715 | ! {\tt start\_row}, the last row in {\tt end\_row}. |
---|
| 1716 | ! |
---|
| 1717 | ! !INTERFACE: |
---|
| 1718 | |
---|
| 1719 | subroutine global_row_range_(sMat, comm, start_row, end_row) |
---|
| 1720 | ! |
---|
| 1721 | ! !USES: |
---|
| 1722 | ! |
---|
| 1723 | use m_die |
---|
| 1724 | |
---|
| 1725 | use m_AttrVect, only : AttrVect_lsize => lsize |
---|
| 1726 | use m_AttrVect, only : AttrVect_indexIA => indexIA |
---|
| 1727 | |
---|
| 1728 | implicit none |
---|
| 1729 | |
---|
| 1730 | ! !INPUT PARAMETERS: |
---|
| 1731 | |
---|
| 1732 | type(SparseMatrix), intent(in) :: sMat |
---|
| 1733 | integer, intent(in) :: comm |
---|
| 1734 | |
---|
| 1735 | ! !OUTPUT PARAMETERS: |
---|
| 1736 | |
---|
| 1737 | integer, intent(out) :: start_row |
---|
| 1738 | integer, intent(out) :: end_row |
---|
| 1739 | |
---|
| 1740 | ! !REVISION HISTORY: |
---|
| 1741 | ! 15Jan01 - Jay Larson <larson@mcs.anl.gov> - API specification. |
---|
| 1742 | ! 25Feb01 - Jay Larson <larson@mcs.anl.gov> - Initial prototype. |
---|
| 1743 | ! 23Apr01 - Jay Larson <larson@mcs.anl.gov> - Modified to accomodate |
---|
| 1744 | ! changes to the SparseMatrix type. |
---|
| 1745 | !EOP ___________________________________________________________________ |
---|
| 1746 | ! |
---|
| 1747 | character(len=*),parameter :: myname_=myname//'::global_row_range_' |
---|
| 1748 | |
---|
| 1749 | integer :: i, igrow, lsize |
---|
| 1750 | |
---|
| 1751 | igrow = AttrVect_indexIA(sMat%data, 'grow', dieWith=myname_) |
---|
| 1752 | lsize = AttrVect_lsize(sMat%data) |
---|
| 1753 | |
---|
| 1754 | ! Initialize start_row and end_row: |
---|
| 1755 | |
---|
| 1756 | start_row = sMat%data%iAttr(igrow,1) |
---|
| 1757 | end_row = sMat%data%iAttr(igrow,1) |
---|
| 1758 | |
---|
| 1759 | do i=1,lsize |
---|
| 1760 | start_row = min(start_row, sMat%data%iAttr(igrow,i)) |
---|
| 1761 | end_row = max(end_row, sMat%data%iAttr(igrow,i)) |
---|
| 1762 | end do |
---|
| 1763 | |
---|
| 1764 | end subroutine global_row_range_ |
---|
| 1765 | |
---|
| 1766 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1767 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1768 | !BOP ------------------------------------------------------------------- |
---|
| 1769 | ! |
---|
| 1770 | ! !IROUTINE: local_col_range_ - Local Column Extent of Non-zero Elements |
---|
| 1771 | ! |
---|
| 1772 | ! !DESCRIPTION: This routine examines the input distributed |
---|
| 1773 | ! {\tt SparseMatrix} variable {\tt sMat}, and returns the range of |
---|
| 1774 | ! local column values having nonzero elements. The first local column |
---|
| 1775 | ! with nonzero elements is returned in the {\tt INTEGER} argument |
---|
| 1776 | ! {\tt start\_col}, the last column in {\tt end\_col}. |
---|
| 1777 | ! |
---|
| 1778 | ! !INTERFACE: |
---|
| 1779 | |
---|
| 1780 | subroutine local_col_range_(sMat, start_col, end_col) |
---|
| 1781 | ! |
---|
| 1782 | ! !USES: |
---|
| 1783 | ! |
---|
| 1784 | use m_die |
---|
| 1785 | |
---|
| 1786 | use m_AttrVect, only : AttrVect_lsize => lsize |
---|
| 1787 | use m_AttrVect, only : AttrVect_indexIA => indexIA |
---|
| 1788 | |
---|
| 1789 | implicit none |
---|
| 1790 | |
---|
| 1791 | ! !INPUT PARAMETERS: |
---|
| 1792 | |
---|
| 1793 | type(SparseMatrix), intent(in) :: sMat |
---|
| 1794 | |
---|
| 1795 | ! !OUTPUT PARAMETERS: |
---|
| 1796 | |
---|
| 1797 | integer, intent(out) :: start_col |
---|
| 1798 | integer, intent(out) :: end_col |
---|
| 1799 | |
---|
| 1800 | ! !REVISION HISTORY: |
---|
| 1801 | ! 15Jan01 - Jay Larson <larson@mcs.anl.gov> - API specification. |
---|
| 1802 | ! 25Feb01 - Jay Larson <larson@mcs.anl.gov> - Initial prototype. |
---|
| 1803 | ! 23Apr01 - Jay Larson <larson@mcs.anl.gov> - Modified to accomodate |
---|
| 1804 | ! changes to the SparseMatrix type. |
---|
| 1805 | !EOP ___________________________________________________________________ |
---|
| 1806 | ! |
---|
| 1807 | character(len=*),parameter :: myname_=myname//'::local_col_range_' |
---|
| 1808 | |
---|
| 1809 | integer :: i, ilcol, lsize |
---|
| 1810 | |
---|
| 1811 | ilcol = AttrVect_indexIA(sMat%data, 'lcol') |
---|
| 1812 | lsize = AttrVect_lsize(sMat%data) |
---|
| 1813 | |
---|
| 1814 | ! Initialize start_col and end_col: |
---|
| 1815 | |
---|
| 1816 | start_col = sMat%data%iAttr(ilcol,1) |
---|
| 1817 | end_col = sMat%data%iAttr(ilcol,1) |
---|
| 1818 | |
---|
| 1819 | do i=1,lsize |
---|
| 1820 | start_col = min(start_col, sMat%data%iAttr(ilcol,i)) |
---|
| 1821 | end_col = max(end_col, sMat%data%iAttr(ilcol,i)) |
---|
| 1822 | end do |
---|
| 1823 | |
---|
| 1824 | end subroutine local_col_range_ |
---|
| 1825 | |
---|
| 1826 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1827 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1828 | !BOP ------------------------------------------------------------------- |
---|
| 1829 | ! |
---|
| 1830 | ! !IROUTINE: global_col_range_ - Global Column Extent of Non-zero Elements |
---|
| 1831 | ! |
---|
| 1832 | ! !DESCRIPTION: This routine examines the input distributed |
---|
| 1833 | ! {\tt SparseMatrix} variable {\tt sMat}, and returns the range of |
---|
| 1834 | ! global column values having nonzero elements. The first global |
---|
| 1835 | ! column with nonzero elements is returned in the {\tt INTEGER} argument |
---|
| 1836 | ! {\tt start\_col}, the last column in {\tt end\_col}. |
---|
| 1837 | ! |
---|
| 1838 | ! !INTERFACE: |
---|
| 1839 | |
---|
| 1840 | subroutine global_col_range_(sMat, comm, start_col, end_col) |
---|
| 1841 | ! |
---|
| 1842 | ! !USES: |
---|
| 1843 | ! |
---|
| 1844 | use m_die |
---|
| 1845 | |
---|
| 1846 | use m_AttrVect, only : AttrVect_lsize => lsize |
---|
| 1847 | use m_AttrVect, only : AttrVect_indexIA => indexIA |
---|
| 1848 | |
---|
| 1849 | implicit none |
---|
| 1850 | |
---|
| 1851 | ! !INPUT PARAMETERS: |
---|
| 1852 | |
---|
| 1853 | type(SparseMatrix), intent(in) :: sMat |
---|
| 1854 | integer, intent(in) :: comm |
---|
| 1855 | |
---|
| 1856 | ! !OUTPUT PARAMETERS: |
---|
| 1857 | |
---|
| 1858 | integer, intent(out) :: start_col |
---|
| 1859 | integer, intent(out) :: end_col |
---|
| 1860 | |
---|
| 1861 | ! !REVISION HISTORY: |
---|
| 1862 | ! 15Jan01 - Jay Larson <larson@mcs.anl.gov> - API specification. |
---|
| 1863 | ! 25Feb01 - Jay Larson <larson@mcs.anl.gov> - Initial prototype. |
---|
| 1864 | ! 23Apr01 - Jay Larson <larson@mcs.anl.gov> - Modified to accomodate |
---|
| 1865 | ! changes to the SparseMatrix type. |
---|
| 1866 | !EOP ___________________________________________________________________ |
---|
| 1867 | ! |
---|
| 1868 | character(len=*),parameter :: myname_=myname//'::global_col_range_' |
---|
| 1869 | |
---|
| 1870 | integer :: i, igcol, lsize |
---|
| 1871 | |
---|
| 1872 | igcol = AttrVect_indexIA(sMat%data, 'gcol') |
---|
| 1873 | lsize = AttrVect_lsize(sMat%data) |
---|
| 1874 | |
---|
| 1875 | ! Initialize start_col and end_col: |
---|
| 1876 | |
---|
| 1877 | start_col = sMat%data%iAttr(igcol,1) |
---|
| 1878 | end_col = sMat%data%iAttr(igcol,1) |
---|
| 1879 | |
---|
| 1880 | do i=1,lsize |
---|
| 1881 | start_col = min(start_col, sMat%data%iAttr(igcol,i)) |
---|
| 1882 | end_col = max(end_col, sMat%data%iAttr(igcol,i)) |
---|
| 1883 | end do |
---|
| 1884 | |
---|
| 1885 | end subroutine global_col_range_ |
---|
| 1886 | |
---|
| 1887 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1888 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1889 | !BOP ------------------------------------------------------------------- |
---|
| 1890 | ! |
---|
| 1891 | ! !IROUTINE: ComputeSparsitySP_ - Compute Matrix Sparsity |
---|
| 1892 | ! |
---|
| 1893 | ! !DESCRIPTION: This routine computes the sparsity of a consolidated |
---|
| 1894 | ! (all on one process) or distributed {\tt SparseMatrix}. The input |
---|
| 1895 | ! {\tt SparseMatrix} argument {\tt sMat} is examined to determine the |
---|
| 1896 | ! number of nonzero elements it holds, and this value is divided by the |
---|
| 1897 | ! product of the number of rows and columns in {\tt sMat}. If the |
---|
| 1898 | ! optional input argument {\tt comm} is given, the number of elements is |
---|
| 1899 | ! reduced with a SUM operation and the sparsity computed from that. |
---|
| 1900 | ! The resulting value of {\tt sparsity} is computed identically on all |
---|
| 1901 | ! processors. |
---|
| 1902 | ! |
---|
| 1903 | ! The rows and columns in an {\tt sMat} are always the global values so |
---|
| 1904 | ! there is no need to do a reduction on those. |
---|
| 1905 | ! |
---|
| 1906 | ! The calculation is always done in double precision even with the argument is |
---|
| 1907 | ! single precision. |
---|
| 1908 | ! |
---|
| 1909 | ! !INTERFACE: |
---|
| 1910 | |
---|
| 1911 | subroutine ComputeSparsitySP_(sMat, sparsity, comm) |
---|
| 1912 | |
---|
| 1913 | ! |
---|
| 1914 | ! !USES: |
---|
| 1915 | ! |
---|
| 1916 | use m_die |
---|
| 1917 | use m_mpif90 |
---|
| 1918 | use m_realkinds, only : SP, FP, DP |
---|
| 1919 | |
---|
| 1920 | use m_AttrVect, only : AttrVect_lsize => lsize |
---|
| 1921 | |
---|
| 1922 | implicit none |
---|
| 1923 | |
---|
| 1924 | ! !INPUT PARAMETERS: |
---|
| 1925 | |
---|
| 1926 | type(SparseMatrix), intent(in) :: sMat |
---|
| 1927 | integer, optional, intent(in) :: comm |
---|
| 1928 | |
---|
| 1929 | ! !OUTPUT PARAMETERS: |
---|
| 1930 | |
---|
| 1931 | real(SP), intent(out) :: sparsity |
---|
| 1932 | |
---|
| 1933 | ! !REVISION HISTORY: |
---|
| 1934 | ! 04Aug20 - Robert Jacob <jacob@anl.gov> - new algorithm |
---|
| 1935 | ! 23Apr01 - Jay Larson <larson@mcs.anl.gov> - New routine. |
---|
| 1936 | ! |
---|
| 1937 | !EOP ___________________________________________________________________ |
---|
| 1938 | ! |
---|
| 1939 | character(len=*),parameter :: myname_=myname//'::ComputeSparsitySP_' |
---|
| 1940 | |
---|
| 1941 | integer :: num_elements, num_rows, num_cols |
---|
| 1942 | real(DP) :: Lnum_elements, Lnum_rows, Lnum_cols |
---|
| 1943 | real(DP) :: tot_elements |
---|
| 1944 | integer :: ierr |
---|
| 1945 | |
---|
| 1946 | ! Extract number of nonzero elements and convert to real |
---|
| 1947 | |
---|
| 1948 | num_elements = lsize_(sMat) |
---|
| 1949 | Lnum_elements = REAL(num_elements,DP) |
---|
| 1950 | |
---|
| 1951 | ! If a communicator handle is present, sum up the |
---|
| 1952 | ! distributed num_elements to all processes. If not, |
---|
| 1953 | ! use the value of lnum_elements |
---|
| 1954 | |
---|
| 1955 | if(present(comm)) then |
---|
| 1956 | call MPI_ALLREDUCE(Lnum_elements, tot_elements, 1, MP_REAL8, & |
---|
| 1957 | MP_SUM, comm, ierr) |
---|
| 1958 | if(ierr /= 0) then |
---|
| 1959 | call MP_perr_die(myname_,"MPI_ALLREDUCE(MySparsity...",ierr) |
---|
| 1960 | endif |
---|
| 1961 | else |
---|
| 1962 | tot_elements = Lnum_elements |
---|
| 1963 | endif |
---|
| 1964 | |
---|
| 1965 | ! Extract number of rows and convert to real |
---|
| 1966 | |
---|
| 1967 | num_rows = nRows_(sMat) |
---|
| 1968 | Lnum_rows = REAL(num_rows,DP) |
---|
| 1969 | |
---|
| 1970 | ! Extract number of columns and convert to real |
---|
| 1971 | |
---|
| 1972 | num_cols = nCols_(sMat) |
---|
| 1973 | Lnum_cols = REAL(num_cols,DP) |
---|
| 1974 | |
---|
| 1975 | ! Compute sparsity |
---|
| 1976 | |
---|
| 1977 | sparsity = tot_elements / (Lnum_rows * Lnum_cols) |
---|
| 1978 | |
---|
| 1979 | end subroutine ComputeSparsitySP_ |
---|
| 1980 | |
---|
| 1981 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 1982 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 1983 | ! ---------------------------------------------------------------------- |
---|
| 1984 | ! |
---|
| 1985 | ! !IROUTINE: ComputeSparsityDP_ - Compute Matrix Sparsity |
---|
| 1986 | ! |
---|
| 1987 | ! !DESCRIPTION: |
---|
| 1988 | ! Double precision version of ComputeSparsitySP_ |
---|
| 1989 | ! |
---|
| 1990 | ! !INTERFACE: |
---|
| 1991 | |
---|
| 1992 | subroutine ComputeSparsityDP_(sMat, sparsity, comm) |
---|
| 1993 | |
---|
| 1994 | ! |
---|
| 1995 | ! !USES: |
---|
| 1996 | ! |
---|
| 1997 | use m_die |
---|
| 1998 | use m_mpif90 |
---|
| 1999 | use m_realkinds, only : DP, FP |
---|
| 2000 | |
---|
| 2001 | use m_AttrVect, only : AttrVect_lsize => lsize |
---|
| 2002 | |
---|
| 2003 | implicit none |
---|
| 2004 | |
---|
| 2005 | ! !INPUT PARAMETERS: |
---|
| 2006 | |
---|
| 2007 | type(SparseMatrix), intent(in) :: sMat |
---|
| 2008 | integer, optional, intent(in) :: comm |
---|
| 2009 | |
---|
| 2010 | ! !OUTPUT PARAMETERS: |
---|
| 2011 | |
---|
| 2012 | real(DP), intent(out) :: sparsity |
---|
| 2013 | |
---|
| 2014 | ! !REVISION HISTORY: |
---|
| 2015 | ! 04Aug20 - Robert Jacob <jacob@anl.gov> - new algorithm |
---|
| 2016 | ! 23Apr01 - Jay Larson <larson@mcs.anl.gov> - New routine. |
---|
| 2017 | ! |
---|
| 2018 | ! ______________________________________________________________________ |
---|
| 2019 | ! |
---|
| 2020 | character(len=*),parameter :: myname_=myname//'::ComputeSparsityDP_' |
---|
| 2021 | |
---|
| 2022 | integer :: num_elements, num_rows, num_cols |
---|
| 2023 | real(FP) :: Lnum_elements, Lnum_rows, Lnum_cols |
---|
| 2024 | real(FP) :: tot_elements |
---|
| 2025 | integer :: ierr |
---|
| 2026 | |
---|
| 2027 | ! Extract number of nonzero elements and convert to real |
---|
| 2028 | |
---|
| 2029 | num_elements = lsize_(sMat) |
---|
| 2030 | Lnum_elements = REAL(num_elements,FP) |
---|
| 2031 | |
---|
| 2032 | ! If a communicator handle is present, sum up the |
---|
| 2033 | ! distributed num_elements to all processes. If not, |
---|
| 2034 | ! use the value of lnum_elements |
---|
| 2035 | |
---|
| 2036 | if(present(comm)) then |
---|
| 2037 | call MPI_ALLREDUCE(Lnum_elements, tot_elements, 1, MP_REAL8, & |
---|
| 2038 | MP_SUM, comm, ierr) |
---|
| 2039 | if(ierr /= 0) then |
---|
| 2040 | call MP_perr_die(myname_,"MPI_ALLREDUCE(MySparsity...",ierr) |
---|
| 2041 | endif |
---|
| 2042 | else |
---|
| 2043 | tot_elements = Lnum_elements |
---|
| 2044 | endif |
---|
| 2045 | |
---|
| 2046 | ! Extract number of rows and convert to real |
---|
| 2047 | |
---|
| 2048 | num_rows = nRows_(sMat) |
---|
| 2049 | Lnum_rows = REAL(num_rows,FP) |
---|
| 2050 | |
---|
| 2051 | ! Extract number of columns and convert to real |
---|
| 2052 | |
---|
| 2053 | num_cols = nCols_(sMat) |
---|
| 2054 | Lnum_cols = REAL(num_cols,FP) |
---|
| 2055 | |
---|
| 2056 | ! Compute sparsity |
---|
| 2057 | |
---|
| 2058 | sparsity = tot_elements / (Lnum_rows * Lnum_cols) |
---|
| 2059 | |
---|
| 2060 | end subroutine ComputeSparsityDP_ |
---|
| 2061 | |
---|
| 2062 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 2063 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 2064 | !BOP ------------------------------------------------------------------- |
---|
| 2065 | ! |
---|
| 2066 | ! !IROUTINE: CheckBounds_ - Check for Out-of-Bounds Row/Column Values |
---|
| 2067 | ! |
---|
| 2068 | ! !DESCRIPTION: This routine examines the input distributed |
---|
| 2069 | ! {\tt SparseMatrix} variable {\tt sMat}, and examines the global row |
---|
| 2070 | ! and column index for each element, comparing them with the known |
---|
| 2071 | ! maximum values for each (as returned by the routines {\tt nRows\_()} |
---|
| 2072 | ! and {\tt nCols\_()}, respectively). If global row or column entries |
---|
| 2073 | ! are non-positive, or greater than the defined maximum values, this |
---|
| 2074 | ! routine stops execution with an error message. If no out-of-bounds |
---|
| 2075 | ! values are detected, the output {\tt INTEGER} status {\tt ierror} is |
---|
| 2076 | ! set to zero. |
---|
| 2077 | ! |
---|
| 2078 | ! !INTERFACE: |
---|
| 2079 | |
---|
| 2080 | subroutine CheckBounds_(sMat, ierror) |
---|
| 2081 | ! |
---|
| 2082 | ! !USES: |
---|
| 2083 | ! |
---|
| 2084 | use m_die |
---|
| 2085 | |
---|
| 2086 | use m_AttrVect, only : AttrVect_lsize => lsize |
---|
| 2087 | use m_AttrVect, only : AttrVect_indexIA => indexIA |
---|
| 2088 | |
---|
| 2089 | implicit none |
---|
| 2090 | |
---|
| 2091 | ! !INPUT PARAMETERS: |
---|
| 2092 | |
---|
| 2093 | type(SparseMatrix), intent(in) :: sMat |
---|
| 2094 | |
---|
| 2095 | ! !OUTPUT PARAMETERS: |
---|
| 2096 | |
---|
| 2097 | integer, intent(out) :: ierror |
---|
| 2098 | |
---|
| 2099 | ! !REVISION HISTORY: |
---|
| 2100 | ! 24Apr01 - Jay Larson <larson@mcs.anl.gov> - Initial prototype. |
---|
| 2101 | !EOP ___________________________________________________________________ |
---|
| 2102 | ! |
---|
| 2103 | character(len=*),parameter :: myname_=myname//'::CheckBounds_' |
---|
| 2104 | |
---|
| 2105 | integer :: MaxRow, MaxCol, NumElements |
---|
| 2106 | integer :: igrow, igcol |
---|
| 2107 | integer :: i |
---|
| 2108 | |
---|
| 2109 | ! Initially, set ierror to zero (success): |
---|
| 2110 | |
---|
| 2111 | ierror = 0 |
---|
| 2112 | |
---|
| 2113 | ! Query sMat to find the number of rows and columns: |
---|
| 2114 | |
---|
| 2115 | MaxRow = nRows_(sMat) |
---|
| 2116 | MaxCol = nCols_(sMat) |
---|
| 2117 | |
---|
| 2118 | ! Query sMat for the number of nonzero elements: |
---|
| 2119 | |
---|
| 2120 | NumElements = lsize_(sMat) |
---|
| 2121 | |
---|
| 2122 | ! Query sMat to index global row and column storage indices: |
---|
| 2123 | |
---|
| 2124 | igrow = indexIA_(sMat=sMat,item='grow',dieWith=myname_) |
---|
| 2125 | igcol = indexIA_(sMat=sMat,item='gcol',dieWith=myname_) |
---|
| 2126 | |
---|
| 2127 | ! Scan the entries of sMat for row or column elements that |
---|
| 2128 | ! are out-of-bounds. Here, out-of-bounds means: 1) non- |
---|
| 2129 | ! positive row or column indices; 2) row or column indices |
---|
| 2130 | ! exceeding the stated number of rows or columns. |
---|
| 2131 | |
---|
| 2132 | do i=1,NumElements |
---|
| 2133 | |
---|
| 2134 | ! Row index out of bounds? |
---|
| 2135 | |
---|
| 2136 | if((sMat%data%iAttr(igrow,i) > MaxRow) .or. & |
---|
| 2137 | (sMat%data%iAttr(igrow,i) <= 0)) then |
---|
| 2138 | ierror = 1 |
---|
| 2139 | call die(myname_,"Row index out of bounds",ierror) |
---|
| 2140 | endif |
---|
| 2141 | |
---|
| 2142 | ! Column index out of bounds? |
---|
| 2143 | |
---|
| 2144 | if((sMat%data%iAttr(igcol,i) > MaxCol) .or. & |
---|
| 2145 | (sMat%data%iAttr(igcol,i) <= 0)) then |
---|
| 2146 | ierror = 2 |
---|
| 2147 | call die(myname_,"Column index out of bounds",ierror) |
---|
| 2148 | endif |
---|
| 2149 | |
---|
| 2150 | end do |
---|
| 2151 | |
---|
| 2152 | end subroutine CheckBounds_ |
---|
| 2153 | |
---|
| 2154 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 2155 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 2156 | !BOP ------------------------------------------------------------------- |
---|
| 2157 | ! |
---|
| 2158 | ! !IROUTINE: row_sumSP_ - Sum Elements in Each Row |
---|
| 2159 | ! |
---|
| 2160 | ! !DESCRIPTION: |
---|
| 2161 | ! Given an input {\tt SparseMatrix} argument {\tt sMat}, {\tt row\_sum\_()} |
---|
| 2162 | ! returns the number of the rows {\tt num\_rows} in the sparse matrix and |
---|
| 2163 | ! the sum of the elements in each row in the array {\tt sums}. The input |
---|
| 2164 | ! argument {\tt comm} is the Fortran 90 MPI communicator handle used to |
---|
| 2165 | ! determine the number of rows and perform the sums. The output arguments |
---|
| 2166 | ! {\tt num\_rows} and {\tt sums} are valid on all processes. |
---|
| 2167 | ! |
---|
| 2168 | ! {\bf N.B.: } This routine allocates an array {\tt sums}. The user is |
---|
| 2169 | ! responsible for deallocating this array when it is no longer needed. |
---|
| 2170 | ! Failure to do so will cause a memory leak. |
---|
| 2171 | ! |
---|
| 2172 | ! !INTERFACE: |
---|
| 2173 | |
---|
| 2174 | subroutine row_sumSP_(sMat, num_rows, sums, comm) |
---|
| 2175 | |
---|
| 2176 | ! |
---|
| 2177 | ! !USES: |
---|
| 2178 | ! |
---|
| 2179 | use m_die |
---|
| 2180 | use m_mpif90 |
---|
| 2181 | use m_realkinds, only : SP, FP |
---|
| 2182 | |
---|
| 2183 | use m_AttrVect, only : AttrVect_lsize => lsize |
---|
| 2184 | use m_AttrVect, only : AttrVect_indexIA => indexIA |
---|
| 2185 | use m_AttrVect, only : AttrVect_indexRA => indexRA |
---|
| 2186 | |
---|
| 2187 | implicit none |
---|
| 2188 | |
---|
| 2189 | ! !INPUT PARAMETERS: |
---|
| 2190 | |
---|
| 2191 | type(SparseMatrix), intent(in) :: sMat |
---|
| 2192 | integer, intent(in) :: comm |
---|
| 2193 | |
---|
| 2194 | ! !OUTPUT PARAMETERS: |
---|
| 2195 | |
---|
| 2196 | integer, intent(out) :: num_rows |
---|
| 2197 | real(SP), dimension(:), pointer :: sums |
---|
| 2198 | |
---|
| 2199 | |
---|
| 2200 | |
---|
| 2201 | ! !REVISION HISTORY: |
---|
| 2202 | ! 15Jan01 - Jay Larson <larson@mcs.anl.gov> - API specification. |
---|
| 2203 | ! 25Jan01 - Jay Larson <larson@mcs.anl.gov> - Prototype code. |
---|
| 2204 | ! 23Apr01 - Jay Larson <larson@mcs.anl.gov> - Modified to accomodate |
---|
| 2205 | ! changes to the SparseMatrix type. |
---|
| 2206 | ! 18May01 - R. Jacob <jacob@mcs.anl.gov> - Use MP_TYPE function |
---|
| 2207 | ! to set type in the mpi_allreduce |
---|
| 2208 | !EOP ___________________________________________________________________ |
---|
| 2209 | ! |
---|
| 2210 | character(len=*),parameter :: myname_=myname//'::row_sumSP_' |
---|
| 2211 | |
---|
| 2212 | integer :: i, igrow, ierr, iwgt, lsize, myID |
---|
| 2213 | integer :: start_row, end_row |
---|
| 2214 | integer :: mp_Type_lsums |
---|
| 2215 | real(FP), dimension(:), allocatable :: lsums |
---|
| 2216 | real(FP), dimension(:), allocatable :: gsums |
---|
| 2217 | |
---|
| 2218 | ! Determine local rank |
---|
| 2219 | |
---|
| 2220 | call MP_COMM_RANK(comm, myID, ierr) |
---|
| 2221 | |
---|
| 2222 | ! Determine on each process the row of global row indices: |
---|
| 2223 | |
---|
| 2224 | call global_row_range_(sMat, comm, start_row, end_row) |
---|
| 2225 | |
---|
| 2226 | ! Determine across the communicator the _maximum_ value of |
---|
| 2227 | ! end_row, which will be assigned to num_rows on each process: |
---|
| 2228 | |
---|
| 2229 | call MPI_ALLREDUCE(end_row, num_rows, 1, MP_INTEGER, MP_MAX, & |
---|
| 2230 | comm, ierr) |
---|
| 2231 | if(ierr /= 0) then |
---|
| 2232 | call MP_perr_die(myname_,"MPI_ALLREDUCE(end_row...",ierr) |
---|
| 2233 | endif |
---|
| 2234 | |
---|
| 2235 | ! Allocate storage for the sums on each process. |
---|
| 2236 | |
---|
| 2237 | allocate(lsums(num_rows), gsums(num_rows), sums(num_rows), stat=ierr) |
---|
| 2238 | |
---|
| 2239 | if(ierr /= 0) then |
---|
| 2240 | call die(myname_,"allocate(lsums(...",ierr) |
---|
| 2241 | endif |
---|
| 2242 | |
---|
| 2243 | ! Compute the local entries to lsum(1:num_rows) for each process: |
---|
| 2244 | |
---|
| 2245 | lsize = AttrVect_lsize(sMat%data) |
---|
| 2246 | igrow = AttrVect_indexIA(aV=sMat%data,item='grow',dieWith=myname_) |
---|
| 2247 | iwgt = AttrVect_indexRA(aV=sMat%data,item='weight',dieWith=myname_) |
---|
| 2248 | |
---|
| 2249 | lsums = 0._FP |
---|
| 2250 | do i=1,lsize |
---|
| 2251 | lsums(sMat%data%iAttr(igrow,i)) = lsums(sMat%data%iAttr(igrow,i)) + & |
---|
| 2252 | sMat%data%rAttr(iwgt,i) |
---|
| 2253 | end do |
---|
| 2254 | |
---|
| 2255 | ! Compute the global sum of the entries of lsums so that all |
---|
| 2256 | ! processes own the global sums. |
---|
| 2257 | |
---|
| 2258 | mp_Type_lsums=MP_Type(lsums) |
---|
| 2259 | call MPI_ALLREDUCE(lsums, gsums, num_rows, mp_Type_lsums, MP_SUM, comm, ierr) |
---|
| 2260 | if(ierr /= 0) then |
---|
| 2261 | call MP_perr_die(myname_,"MPI_ALLREDUCE(lsums...",ierr) |
---|
| 2262 | endif |
---|
| 2263 | |
---|
| 2264 | ! Copy our temporary array gsums into the output pointer sums |
---|
| 2265 | ! This was done so that lsums and gsums have the same precision (FP) |
---|
| 2266 | ! Precision conversion occurs here from FP to (SP or DP) |
---|
| 2267 | |
---|
| 2268 | sums = gsums |
---|
| 2269 | |
---|
| 2270 | ! Clean up... |
---|
| 2271 | |
---|
| 2272 | deallocate(lsums, gsums, stat=ierr) |
---|
| 2273 | if(ierr /= 0) then |
---|
| 2274 | call die(myname_,"deallocate(lsums...",ierr) |
---|
| 2275 | endif |
---|
| 2276 | |
---|
| 2277 | end subroutine row_sumSP_ |
---|
| 2278 | |
---|
| 2279 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 2280 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 2281 | ! ---------------------------------------------------------------------- |
---|
| 2282 | ! |
---|
| 2283 | ! !IROUTINE: row_sumDP_ - Sum Elements in Each Row |
---|
| 2284 | ! |
---|
| 2285 | ! !DESCRIPTION: |
---|
| 2286 | ! Double precision version of row_sumSP_ |
---|
| 2287 | ! |
---|
| 2288 | ! {\bf N.B.: } This routine allocates an array {\tt sums}. The user is |
---|
| 2289 | ! responsible for deallocating this array when it is no longer needed. |
---|
| 2290 | ! Failure to do so will cause a memory leak. |
---|
| 2291 | ! |
---|
| 2292 | ! !INTERFACE: |
---|
| 2293 | |
---|
| 2294 | subroutine row_sumDP_(sMat, num_rows, sums, comm) |
---|
| 2295 | |
---|
| 2296 | ! |
---|
| 2297 | ! !USES: |
---|
| 2298 | ! |
---|
| 2299 | use m_die |
---|
| 2300 | use m_mpif90 |
---|
| 2301 | |
---|
| 2302 | use m_realkinds, only : DP, FP |
---|
| 2303 | |
---|
| 2304 | use m_AttrVect, only : AttrVect_lsize => lsize |
---|
| 2305 | use m_AttrVect, only : AttrVect_indexIA => indexIA |
---|
| 2306 | use m_AttrVect, only : AttrVect_indexRA => indexRA |
---|
| 2307 | |
---|
| 2308 | implicit none |
---|
| 2309 | |
---|
| 2310 | ! !INPUT PARAMETERS: |
---|
| 2311 | |
---|
| 2312 | type(SparseMatrix), intent(in) :: sMat |
---|
| 2313 | integer, intent(in) :: comm |
---|
| 2314 | |
---|
| 2315 | ! !OUTPUT PARAMETERS: |
---|
| 2316 | |
---|
| 2317 | integer, intent(out) :: num_rows |
---|
| 2318 | real(DP), dimension(:), pointer :: sums |
---|
| 2319 | |
---|
| 2320 | |
---|
| 2321 | |
---|
| 2322 | ! !REVISION HISTORY: |
---|
| 2323 | ! 15Jan01 - Jay Larson <larson@mcs.anl.gov> - API specification. |
---|
| 2324 | ! 25Jan01 - Jay Larson <larson@mcs.anl.gov> - Prototype code. |
---|
| 2325 | ! 23Apr01 - Jay Larson <larson@mcs.anl.gov> - Modified to accomodate |
---|
| 2326 | ! changes to the SparseMatrix type. |
---|
| 2327 | ! 18May01 - R. Jacob <jacob@mcs.anl.gov> - Use MP_TYPE function |
---|
| 2328 | ! to set type in the mpi_allreduce |
---|
| 2329 | ! ______________________________________________________________________ |
---|
| 2330 | ! |
---|
| 2331 | character(len=*),parameter :: myname_=myname//'::row_sumDP_' |
---|
| 2332 | |
---|
| 2333 | integer :: i, igrow, ierr, iwgt, lsize, myID |
---|
| 2334 | integer :: start_row, end_row |
---|
| 2335 | integer :: mp_Type_lsums |
---|
| 2336 | real(FP), dimension(:), allocatable :: lsums |
---|
| 2337 | real(FP), dimension(:), allocatable :: gsums |
---|
| 2338 | |
---|
| 2339 | ! Determine local rank |
---|
| 2340 | |
---|
| 2341 | call MP_COMM_RANK(comm, myID, ierr) |
---|
| 2342 | |
---|
| 2343 | ! Determine on each process the row of global row indices: |
---|
| 2344 | |
---|
| 2345 | call global_row_range_(sMat, comm, start_row, end_row) |
---|
| 2346 | |
---|
| 2347 | ! Determine across the communicator the _maximum_ value of |
---|
| 2348 | ! end_row, which will be assigned to num_rows on each process: |
---|
| 2349 | |
---|
| 2350 | call MPI_ALLREDUCE(end_row, num_rows, 1, MP_INTEGER, MP_MAX, & |
---|
| 2351 | comm, ierr) |
---|
| 2352 | if(ierr /= 0) then |
---|
| 2353 | call MP_perr_die(myname_,"MPI_ALLREDUCE(end_row...",ierr) |
---|
| 2354 | endif |
---|
| 2355 | |
---|
| 2356 | ! Allocate storage for the sums on each process. |
---|
| 2357 | |
---|
| 2358 | allocate(lsums(num_rows), gsums(num_rows), sums(num_rows), stat=ierr) |
---|
| 2359 | |
---|
| 2360 | if(ierr /= 0) then |
---|
| 2361 | call die(myname_,"allocate(lsums(...",ierr) |
---|
| 2362 | endif |
---|
| 2363 | |
---|
| 2364 | ! Compute the local entries to lsum(1:num_rows) for each process: |
---|
| 2365 | |
---|
| 2366 | lsize = AttrVect_lsize(sMat%data) |
---|
| 2367 | igrow = AttrVect_indexIA(aV=sMat%data,item='grow',dieWith=myname_) |
---|
| 2368 | iwgt = AttrVect_indexRA(aV=sMat%data,item='weight',dieWith=myname_) |
---|
| 2369 | |
---|
| 2370 | lsums = 0._FP |
---|
| 2371 | do i=1,lsize |
---|
| 2372 | lsums(sMat%data%iAttr(igrow,i)) = lsums(sMat%data%iAttr(igrow,i)) + & |
---|
| 2373 | sMat%data%rAttr(iwgt,i) |
---|
| 2374 | end do |
---|
| 2375 | |
---|
| 2376 | ! Compute the global sum of the entries of lsums so that all |
---|
| 2377 | ! processes own the global sums. |
---|
| 2378 | |
---|
| 2379 | mp_Type_lsums=MP_Type(lsums) |
---|
| 2380 | call MPI_ALLREDUCE(lsums, gsums, num_rows, mp_Type_lsums, MP_SUM, comm, ierr) |
---|
| 2381 | if(ierr /= 0) then |
---|
| 2382 | call MP_perr_die(myname_,"MPI_ALLREDUCE(lsums...",ierr) |
---|
| 2383 | endif |
---|
| 2384 | |
---|
| 2385 | ! Copy our temporary array gsums into the output pointer sums |
---|
| 2386 | ! This was done so that lsums and gsums have the same precision (FP) |
---|
| 2387 | ! Precision conversion occurs here from FP to (SP or DP) |
---|
| 2388 | |
---|
| 2389 | sums = gsums |
---|
| 2390 | |
---|
| 2391 | ! Clean up... |
---|
| 2392 | |
---|
| 2393 | deallocate(lsums, gsums, stat=ierr) |
---|
| 2394 | if(ierr /= 0) then |
---|
| 2395 | call die(myname_,"deallocate(lsums...",ierr) |
---|
| 2396 | endif |
---|
| 2397 | |
---|
| 2398 | end subroutine row_sumDP_ |
---|
| 2399 | |
---|
| 2400 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 2401 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 2402 | !BOP ------------------------------------------------------------------- |
---|
| 2403 | ! |
---|
| 2404 | ! !IROUTINE: row_sum_checkSP_ - Check Row Sums vs. Valid Values |
---|
| 2405 | ! |
---|
| 2406 | ! !DESCRIPTION: The routine {\tt row\_sum\_check()} sums the rows of |
---|
| 2407 | ! the input distributed (across the communicator identified by {\tt comm}) |
---|
| 2408 | ! {\tt SparseMatrix} variable {\tt sMat}. It then compares these sums |
---|
| 2409 | ! with the {\tt num\_valid} input "valid" values stored in the array |
---|
| 2410 | ! {\tt valid\_sums}. If all of the sums are within the absolute tolerence |
---|
| 2411 | ! specified by the input argument {\tt abs\_tol} of any of the valid values, |
---|
| 2412 | ! the output {\tt LOGICAL} flag {\tt valid} is set to {\tt .TRUE}. |
---|
| 2413 | ! Otherwise, this flag is returned with value {\tt .FALSE}. |
---|
| 2414 | ! |
---|
| 2415 | ! !INTERFACE: |
---|
| 2416 | |
---|
| 2417 | subroutine row_sum_checkSP_(sMat, comm, num_valid, valid_sums, abs_tol, valid) |
---|
| 2418 | |
---|
| 2419 | ! |
---|
| 2420 | ! !USES: |
---|
| 2421 | ! |
---|
| 2422 | use m_die |
---|
| 2423 | use m_realkinds, only : SP, FP |
---|
| 2424 | |
---|
| 2425 | implicit none |
---|
| 2426 | |
---|
| 2427 | ! !INPUT PARAMETERS: |
---|
| 2428 | |
---|
| 2429 | type(SparseMatrix), intent(in) :: sMat |
---|
| 2430 | integer, intent(in) :: comm |
---|
| 2431 | integer, intent(in) :: num_valid |
---|
| 2432 | real(SP), intent(in) :: valid_sums(num_valid) |
---|
| 2433 | real(SP), intent(in) :: abs_tol |
---|
| 2434 | |
---|
| 2435 | ! !OUTPUT PARAMETERS: |
---|
| 2436 | |
---|
| 2437 | logical, intent(out) :: valid |
---|
| 2438 | |
---|
| 2439 | ! !REVISION HISTORY: |
---|
| 2440 | ! 15Jan01 - Jay Larson <larson@mcs.anl.gov> - API specification. |
---|
| 2441 | ! 25Feb01 - Jay Larson <larson@mcs.anl.gov> - Prototype code. |
---|
| 2442 | ! 06Jan03 - R. Jacob <jacob@mcs.anl.gov> - create DP and SP versions |
---|
| 2443 | !EOP ___________________________________________________________________ |
---|
| 2444 | ! |
---|
| 2445 | character(len=*),parameter :: myname_=myname//'::row_sum_checkSP_' |
---|
| 2446 | |
---|
| 2447 | integer :: i, j, num_invalid, num_rows |
---|
| 2448 | real(FP), dimension(:), pointer :: sums |
---|
| 2449 | |
---|
| 2450 | ! Compute row sums: |
---|
| 2451 | |
---|
| 2452 | call row_sum(sMat, num_rows, sums, comm) |
---|
| 2453 | |
---|
| 2454 | ! Initialize for the scanning loop (assume the matrix row |
---|
| 2455 | ! sums are valid): |
---|
| 2456 | |
---|
| 2457 | valid = .TRUE. |
---|
| 2458 | i = 1 |
---|
| 2459 | |
---|
| 2460 | SCAN_LOOP: do |
---|
| 2461 | |
---|
| 2462 | ! Count the number of elements in valid_sums(:) that |
---|
| 2463 | ! are separated from sums(i) by more than abs_tol |
---|
| 2464 | |
---|
| 2465 | num_invalid = 0 |
---|
| 2466 | |
---|
| 2467 | do j=1,num_valid |
---|
| 2468 | if(abs(sums(i) - valid_sums(j)) > abs_tol) then |
---|
| 2469 | num_invalid = num_invalid + 1 |
---|
| 2470 | endif |
---|
| 2471 | end do |
---|
| 2472 | |
---|
| 2473 | ! If num_invalid = num_valid, then we have failed to |
---|
| 2474 | ! find a valid sum value within abs_tol of sums(i). This |
---|
| 2475 | ! one failure is enough to halt the process. |
---|
| 2476 | |
---|
| 2477 | if(num_invalid == num_valid) then |
---|
| 2478 | valid = .FALSE. |
---|
| 2479 | EXIT |
---|
| 2480 | endif |
---|
| 2481 | |
---|
| 2482 | ! Prepare index i for the next element of sums(:) |
---|
| 2483 | |
---|
| 2484 | i = i + 1 |
---|
| 2485 | if( i > num_rows) EXIT |
---|
| 2486 | |
---|
| 2487 | end do SCAN_LOOP |
---|
| 2488 | |
---|
| 2489 | end subroutine row_sum_checkSP_ |
---|
| 2490 | |
---|
| 2491 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 2492 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 2493 | ! ---------------------------------------------------------------------- |
---|
| 2494 | ! |
---|
| 2495 | ! !IROUTINE: row_sum_checkDP_ - Check Row Sums vs. Valid Values |
---|
| 2496 | ! |
---|
| 2497 | ! !DESCRIPTION: |
---|
| 2498 | ! Double precision version of row_sum_checkSP |
---|
| 2499 | ! |
---|
| 2500 | ! !INTERFACE: |
---|
| 2501 | |
---|
| 2502 | subroutine row_sum_checkDP_(sMat, comm, num_valid, valid_sums, abs_tol, valid) |
---|
| 2503 | |
---|
| 2504 | ! |
---|
| 2505 | ! !USES: |
---|
| 2506 | ! |
---|
| 2507 | use m_die |
---|
| 2508 | use m_realkinds, only : DP, FP |
---|
| 2509 | |
---|
| 2510 | implicit none |
---|
| 2511 | |
---|
| 2512 | ! !INPUT PARAMETERS: |
---|
| 2513 | |
---|
| 2514 | type(SparseMatrix), intent(in) :: sMat |
---|
| 2515 | integer, intent(in) :: comm |
---|
| 2516 | integer, intent(in) :: num_valid |
---|
| 2517 | real(DP), intent(in) :: valid_sums(num_valid) |
---|
| 2518 | real(DP), intent(in) :: abs_tol |
---|
| 2519 | |
---|
| 2520 | ! !OUTPUT PARAMETERS: |
---|
| 2521 | |
---|
| 2522 | logical, intent(out) :: valid |
---|
| 2523 | |
---|
| 2524 | ! !REVISION HISTORY: |
---|
| 2525 | ! 15Jan01 - Jay Larson <larson@mcs.anl.gov> - API specification. |
---|
| 2526 | ! 25Feb01 - Jay Larson <larson@mcs.anl.gov> - Prototype code. |
---|
| 2527 | ! 06Jan03 - R. Jacob <jacob@mcs.anl.gov> - create DP and SP versions |
---|
| 2528 | ! ______________________________________________________________________ |
---|
| 2529 | ! |
---|
| 2530 | character(len=*),parameter :: myname_=myname//'::row_sum_checkDP_' |
---|
| 2531 | |
---|
| 2532 | integer :: i, j, num_invalid, num_rows |
---|
| 2533 | real(FP), dimension(:), pointer :: sums |
---|
| 2534 | |
---|
| 2535 | ! Compute row sums: |
---|
| 2536 | |
---|
| 2537 | call row_sum(sMat, num_rows, sums, comm) |
---|
| 2538 | |
---|
| 2539 | ! Initialize for the scanning loop (assume the matrix row |
---|
| 2540 | ! sums are valid): |
---|
| 2541 | |
---|
| 2542 | valid = .TRUE. |
---|
| 2543 | i = 1 |
---|
| 2544 | |
---|
| 2545 | SCAN_LOOP: do |
---|
| 2546 | |
---|
| 2547 | ! Count the number of elements in valid_sums(:) that |
---|
| 2548 | ! are separated from sums(i) by more than abs_tol |
---|
| 2549 | |
---|
| 2550 | num_invalid = 0 |
---|
| 2551 | |
---|
| 2552 | do j=1,num_valid |
---|
| 2553 | if(abs(sums(i) - valid_sums(j)) > abs_tol) then |
---|
| 2554 | num_invalid = num_invalid + 1 |
---|
| 2555 | endif |
---|
| 2556 | end do |
---|
| 2557 | |
---|
| 2558 | ! If num_invalid = num_valid, then we have failed to |
---|
| 2559 | ! find a valid sum value within abs_tol of sums(i). This |
---|
| 2560 | ! one failure is enough to halt the process. |
---|
| 2561 | |
---|
| 2562 | if(num_invalid == num_valid) then |
---|
| 2563 | valid = .FALSE. |
---|
| 2564 | EXIT |
---|
| 2565 | endif |
---|
| 2566 | |
---|
| 2567 | ! Prepare index i for the next element of sums(:) |
---|
| 2568 | |
---|
| 2569 | i = i + 1 |
---|
| 2570 | if( i > num_rows) EXIT |
---|
| 2571 | |
---|
| 2572 | end do SCAN_LOOP |
---|
| 2573 | |
---|
| 2574 | end subroutine row_sum_checkDP_ |
---|
| 2575 | |
---|
| 2576 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 2577 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 2578 | !BOP ------------------------------------------------------------------- |
---|
| 2579 | ! |
---|
| 2580 | ! !IROUTINE: Sort_ - Generate Index Permutation |
---|
| 2581 | ! |
---|
| 2582 | ! !DESCRIPTION: |
---|
| 2583 | ! The subroutine {\tt Sort\_()} uses a list of sorting keys defined by |
---|
| 2584 | ! the input {\tt List} argument {\tt key\_list}, searches for the appropriate |
---|
| 2585 | ! integer or real attributes referenced by the items in {\tt key\_list} |
---|
| 2586 | ! ( that is, it identifies the appropriate entries in {sMat\%data\%iList} |
---|
| 2587 | ! and {\tt sMat\%data\%rList}), and then uses these keys to generate an index |
---|
| 2588 | ! permutation {\tt perm} that will put the nonzero matrix entries of stored |
---|
| 2589 | ! in {\tt sMat\%data} in lexicographic order as defined by {\tt key\_ist} |
---|
| 2590 | ! (the ordering in {\tt key\_list} being from left to right. The optional |
---|
| 2591 | ! {\tt LOGICAL} array input argument {\tt descend} specifies whether or |
---|
| 2592 | ! not to sort by each key in {\em descending} order or {\em ascending} |
---|
| 2593 | ! order. Entries in {\tt descend} that have value {\tt .TRUE.} correspond |
---|
| 2594 | ! to a sort by the corresponding key in descending order. If the argument |
---|
| 2595 | ! {\tt descend} is not present, the sort is performed for all keys in |
---|
| 2596 | ! ascending order. |
---|
| 2597 | ! |
---|
| 2598 | ! !INTERFACE: |
---|
| 2599 | |
---|
| 2600 | subroutine Sort_(sMat, key_list, perm, descend) |
---|
| 2601 | |
---|
| 2602 | ! |
---|
| 2603 | ! !USES: |
---|
| 2604 | ! |
---|
| 2605 | use m_die , only : die |
---|
| 2606 | use m_stdio , only : stderr |
---|
| 2607 | |
---|
| 2608 | use m_List , only : List |
---|
| 2609 | |
---|
| 2610 | use m_AttrVect, only: AttrVect_Sort => Sort |
---|
| 2611 | |
---|
| 2612 | implicit none |
---|
| 2613 | ! |
---|
| 2614 | ! !INPUT PARAMETERS: |
---|
| 2615 | |
---|
| 2616 | type(SparseMatrix), intent(in) :: sMat |
---|
| 2617 | type(List), intent(in) :: key_list |
---|
| 2618 | logical, dimension(:), optional, intent(in) :: descend |
---|
| 2619 | ! |
---|
| 2620 | ! !OUTPUT PARAMETERS: |
---|
| 2621 | |
---|
| 2622 | integer, dimension(:), pointer :: perm |
---|
| 2623 | |
---|
| 2624 | |
---|
| 2625 | ! !REVISION HISTORY: |
---|
| 2626 | ! 24Apr01 - J.W. Larson <larson@mcs.anl.gov> - initial prototype |
---|
| 2627 | !EOP ___________________________________________________________________ |
---|
| 2628 | |
---|
| 2629 | character(len=*),parameter :: myname_=myname//'::Sort_' |
---|
| 2630 | |
---|
| 2631 | if(present(descend)) then |
---|
| 2632 | call AttrVect_Sort(sMat%data, key_list, perm, descend) |
---|
| 2633 | else |
---|
| 2634 | call AttrVect_Sort(sMat%data, key_list, perm) |
---|
| 2635 | endif |
---|
| 2636 | |
---|
| 2637 | end Subroutine Sort_ |
---|
| 2638 | |
---|
| 2639 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 2640 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 2641 | !BOP ------------------------------------------------------------------- |
---|
| 2642 | ! |
---|
| 2643 | ! !IROUTINE: Permute_ - Permute Matrix Elements using Supplied Index Permutation |
---|
| 2644 | ! |
---|
| 2645 | ! !DESCRIPTION: |
---|
| 2646 | ! The subroutine {\tt Permute\_()} uses an input index permutation |
---|
| 2647 | ! {\tt perm} to re-order the entries of the {\tt SparseMatrix} argument |
---|
| 2648 | ! {\tt sMat}. The index permutation {\tt perm} is generated using the |
---|
| 2649 | ! routine {\tt Sort\_()} (in this module). |
---|
| 2650 | ! |
---|
| 2651 | ! !INTERFACE: |
---|
| 2652 | |
---|
| 2653 | subroutine Permute_(sMat, perm) |
---|
| 2654 | |
---|
| 2655 | ! |
---|
| 2656 | ! !USES: |
---|
| 2657 | ! |
---|
| 2658 | use m_die , only : die |
---|
| 2659 | use m_stdio , only : stderr |
---|
| 2660 | |
---|
| 2661 | use m_AttrVect, only: AttrVect_Permute => Permute |
---|
| 2662 | |
---|
| 2663 | implicit none |
---|
| 2664 | ! |
---|
| 2665 | ! !INPUT PARAMETERS: |
---|
| 2666 | |
---|
| 2667 | |
---|
| 2668 | integer, dimension(:), pointer :: perm |
---|
| 2669 | ! |
---|
| 2670 | ! !INPUT/OUTPUT PARAMETERS: |
---|
| 2671 | |
---|
| 2672 | type(SparseMatrix), intent(inout) :: sMat |
---|
| 2673 | |
---|
| 2674 | |
---|
| 2675 | ! !REVISION HISTORY: |
---|
| 2676 | ! 24Apr01 - J.W. Larson <larson@mcs.anl.gov> - initial prototype |
---|
| 2677 | !EOP ___________________________________________________________________ |
---|
| 2678 | |
---|
| 2679 | character(len=*),parameter :: myname_=myname//'::Permute_' |
---|
| 2680 | |
---|
| 2681 | call AttrVect_Permute(sMat%data, perm) |
---|
| 2682 | |
---|
| 2683 | end Subroutine Permute_ |
---|
| 2684 | |
---|
| 2685 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 2686 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 2687 | !BOP ------------------------------------------------------------------- |
---|
| 2688 | ! |
---|
| 2689 | ! !IROUTINE: SortPermute_ - Sort and Permute Matrix Elements |
---|
| 2690 | ! |
---|
| 2691 | ! !DESCRIPTION: |
---|
| 2692 | ! The subroutine {\tt SortPermute\_()} uses a list of sorting keys defined |
---|
| 2693 | ! by the input {\tt List} argument {\tt key\_list}, searches for the |
---|
| 2694 | ! appropriate integer or real attributes referenced by the items in |
---|
| 2695 | ! {\tt key\_ist} ( that is, it identifies the appropriate entries in |
---|
| 2696 | ! {sMat\%data\%iList} and {\tt sMat\%data\%rList}), and then uses these |
---|
| 2697 | ! keys to generate an index permutation that will put the nonzero matrix |
---|
| 2698 | ! entries of stored in {\tt sMat\%data} in lexicographic order as defined |
---|
| 2699 | ! by {\tt key\_list} (the ordering in {\tt key\_list} being from left to |
---|
| 2700 | ! right. The optional {\tt LOGICAL} array input argument {\tt descend} |
---|
| 2701 | ! specifies whether or not to sort by each key in {\em descending} order |
---|
| 2702 | ! or {\em ascending} order. Entries in {\tt descend} that have value |
---|
| 2703 | ! {\tt .TRUE.} correspond to a sort by the corresponding key in descending |
---|
| 2704 | ! order. If the argument {\tt descend} is not present, the sort is |
---|
| 2705 | ! performed for all keys in ascending order. |
---|
| 2706 | ! |
---|
| 2707 | ! Once this index permutation is created, it is applied to re-order the |
---|
| 2708 | ! entries of the {\tt SparseMatrix} argument {\tt sMat} accordingly. |
---|
| 2709 | ! |
---|
| 2710 | ! !INTERFACE: |
---|
| 2711 | |
---|
| 2712 | subroutine SortPermute_(sMat, key_list, descend) |
---|
| 2713 | |
---|
| 2714 | ! |
---|
| 2715 | ! !USES: |
---|
| 2716 | ! |
---|
| 2717 | use m_die , only : die |
---|
| 2718 | use m_stdio , only : stderr |
---|
| 2719 | |
---|
| 2720 | use m_List , only : List |
---|
| 2721 | |
---|
| 2722 | implicit none |
---|
| 2723 | ! |
---|
| 2724 | ! !INPUT PARAMETERS: |
---|
| 2725 | |
---|
| 2726 | type(List), intent(in) :: key_list |
---|
| 2727 | logical, dimension(:), optional, intent(in) :: descend |
---|
| 2728 | ! |
---|
| 2729 | ! !INPUT/OUTPUT PARAMETERS: |
---|
| 2730 | |
---|
| 2731 | type(SparseMatrix), intent(inout) :: sMat |
---|
| 2732 | |
---|
| 2733 | ! !REVISION HISTORY: |
---|
| 2734 | ! 24Apr01 - J.W. Larson <larson@mcs.anl.gov> - initial prototype |
---|
| 2735 | !EOP ___________________________________________________________________ |
---|
| 2736 | |
---|
| 2737 | character(len=*),parameter :: myname_=myname//'::SortPermute_' |
---|
| 2738 | |
---|
| 2739 | integer :: ier |
---|
| 2740 | integer, dimension(:), pointer :: perm |
---|
| 2741 | |
---|
| 2742 | ! Create index permutation perm(:) |
---|
| 2743 | |
---|
| 2744 | if(present(descend)) then |
---|
| 2745 | call Sort_(sMat, key_list, perm, descend) |
---|
| 2746 | else |
---|
| 2747 | call Sort_(sMat, key_list, perm) |
---|
| 2748 | endif |
---|
| 2749 | |
---|
| 2750 | ! Apply index permutation perm(:) to re-order sMat: |
---|
| 2751 | |
---|
| 2752 | call Permute_(sMat, perm) |
---|
| 2753 | |
---|
| 2754 | ! Clean up |
---|
| 2755 | |
---|
| 2756 | deallocate(perm, stat=ier) |
---|
| 2757 | if(ier/=0) call die(myname_, "deallocate(perm)", ier) |
---|
| 2758 | |
---|
| 2759 | end subroutine SortPermute_ |
---|
| 2760 | |
---|
| 2761 | end module m_SparseMatrix |
---|
| 2762 | |
---|
| 2763 | |
---|
| 2764 | |
---|