mutils.F90 3.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139
  1. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  2. ! Math and Computer Science Division, Argonne National Laboratory !
  3. !-----------------------------------------------------------------------
  4. ! CVS mutils.F90,v 1.8 2005-11-18 23:15:38 rloy Exp
  5. ! CVS MCT_2_8_0
  6. !BOP -------------------------------------------------------------------
  7. !
  8. ! !MODULE: mutils -- utilities for the sequential climate example
  9. !
  10. ! !DESCRIPTION:
  11. !
  12. ! !INTERFACE:
  13. !
  14. module mutils
  15. ! module of utilties for the sequential climate example
  16. !
  17. implicit none
  18. private
  19. ! except
  20. ! !PUBLIC TYPES:
  21. public get_index
  22. contains
  23. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  24. ! Math and Computer Science Division, Argonne National Laboratory !
  25. !BOP -------------------------------------------------------------------
  26. !
  27. ! !IROUTINE: get_index - get local index array and size
  28. ! for 3 standard decompositions of a grid.
  29. !
  30. ! !DESCRIPTION:
  31. ! The routine get_index will return a local index array and size that can
  32. ! be passed to a GSMap_init routine for three possible decompositions:
  33. ! R - by row or latitude
  34. ! C - by column or longitude
  35. ! RC - row and column or checkerboard
  36. ! choice is determined by the value of ldecomp.
  37. !
  38. ! !INTERFACE:
  39. subroutine get_index(ldecomp,nprocs,myproc,gnx,gny,gridbuf)
  40. ! !INPUT PARAMETERS:
  41. !
  42. character(len=*),intent(inout) :: ldecomp ! decomp choice
  43. integer,intent(in) :: nprocs ! total number of MPI processes
  44. integer,intent(in) :: myproc ! my rank in local communicator
  45. integer,intent(in) :: gnx ! total points in X direction
  46. integer,intent(in) :: gny ! total points in Y direction
  47. ! !OUTPUT PARAMETERS:
  48. !
  49. integer,dimension(:),pointer :: gridbuf ! local index array
  50. !
  51. !EOP ___________________________________________________________________
  52. integer :: npesx,npesy,ng,ny,n,i,j,nx,ig,jg,nseg,factor
  53. ! default decomp is R
  54. if((trim(ldecomp) .ne. 'R') .and. (ldecomp .ne. 'C') .and. (ldecomp .ne. 'RC')) then
  55. ldecomp = 'R'
  56. endif
  57. ! A 'by-row' or 'by-latitude' decomposition
  58. if(trim(ldecomp) .eq. 'R') then
  59. npesx=1
  60. npesy=nprocs
  61. nx=gnx
  62. ny=gny/npesy
  63. allocate(gridbuf(nx*ny))
  64. n=0
  65. do j=1,ny
  66. do i=1,nx
  67. n=n+1
  68. ig=i
  69. jg = j + myProc*ny
  70. ng =(jg-1)*gnx + ig
  71. gridbuf(n)=ng
  72. enddo
  73. enddo
  74. ! A 'by-column' or 'by-longitude' decomposition
  75. else if (ldecomp .eq. 'C') then
  76. npesx=nprocs
  77. npesy=1
  78. nx=gnx/npesx
  79. ny=gny
  80. allocate(gridbuf(nx*ny))
  81. n=0
  82. do j=1,ny
  83. do i=1,nx
  84. n=n+1
  85. ig=i + myProc*nx
  86. jg= j
  87. ng=(jg-1)*gnx + ig
  88. gridbuf(n)=ng
  89. enddo
  90. enddo
  91. ! A 'row-columen' or 'checkerboard' decomposition
  92. else if (ldecomp .eq. 'RC') then
  93. ! find the closest square
  94. factor=1
  95. do i=2,INT(sqrt(FLOAT(nprocs)))
  96. if ( (nprocs/i) * i .eq. nprocs) then
  97. factor = i
  98. endif
  99. enddo
  100. npesx=factor
  101. npesy=nprocs/factor
  102. nx=gnx/npesx
  103. ny=gny/npesy
  104. ! write(6,*) 'RC',factor,npesy,nx,ny
  105. allocate(gridbuf(nx*ny))
  106. n=0
  107. do j=1,ny
  108. do i=1,nx
  109. n=n+1
  110. ig=mod(myProc,npesx)*nx+i
  111. jg=(myProc/npesx)*ny+j
  112. ng=(jg-1)*gnx + ig
  113. gridbuf(n)=ng
  114. enddo
  115. enddo
  116. endif
  117. end subroutine get_index
  118. end module mutils