coupler.F90 6.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214
  1. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  2. ! Math and Computer Science Division, Argonne National Laboratory !
  3. !-----------------------------------------------------------------------
  4. ! CVS coupler.F90,v 1.6 2006-10-17 21:46:35 jacob Exp
  5. ! CVS MCT_2_8_0
  6. !BOP -------------------------------------------------------------------
  7. !
  8. ! !ROUTINE: coupler -- coupler for sequential model example
  9. !
  10. ! !DESCRIPTION:
  11. ! A coupler subroutine for sequential climate model example.
  12. !
  13. ! !INTERFACE:
  14. !
  15. module coupler
  16. !
  17. ! !USES:
  18. !
  19. ! Get the things needed from MCT by "Use,only" with renaming:
  20. !
  21. !---Domain Decomposition Descriptor DataType and associated methods
  22. use m_GlobalSegMap,only: GlobalSegMap
  23. !---Field Storage DataType and associated methods
  24. use m_AttrVect,only : AttrVect
  25. !---Sparse Matrix DataType and associated methods
  26. use m_SparseMatrix, only : SparseMatrix
  27. use m_SparseMatrix, only : SparseMatrix_clean => clean
  28. use m_SparseMatrix, only : SparseMatrix_init => init
  29. use m_SparseMatrix, only : SparseMatrix_importGRowInd => &
  30. importGlobalRowIndices
  31. use m_SparseMatrix, only : SparseMatrix_importGColInd => &
  32. importGlobalColumnIndices
  33. use m_SparseMatrix, only : SparseMatrix_importMatrixElts => &
  34. importMatrixElements
  35. use m_SparseMatrixPlus, only : SparseMatrixPlus
  36. use m_SparseMatrixPlus, only : SparseMatrixPlus_init => init
  37. use m_SparseMatrixPlus, only : SparseMatrixPlus_clean => clean
  38. use m_SparseMatrixPlus, only : Xonly ! Decompose matrix by row
  39. !---Matrix-Vector multiply methods
  40. use m_MatAttrVectMul, only: MCT_MatVecMul => sMatAvMult
  41. !---MPEU I/O utilities
  42. use m_stdio
  43. use m_ioutil
  44. implicit none
  45. private
  46. ! !PUBLIC MEMBER FUNCTIONS:
  47. public cplinit
  48. public cplrun
  49. public cplfin
  50. ! !PRIVATE DATA MEMBERS
  51. type(SparseMatrixPlus) :: Src2DstMatPlus ! the mapping weights
  52. character(len=*), parameter :: cplname='coupler.F90'
  53. integer :: rank
  54. !EOP ___________________________________________________________________
  55. contains
  56. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  57. ! Math and Computer Science Division, Argonne National Laboratory !
  58. !BOP -------------------------------------------------------------------
  59. !
  60. ! !IROUTINE: cplinit - initialize the coupler
  61. !
  62. ! !INTERFACE:
  63. subroutine cplinit(SrcGSMap,DstGSMap,comm,compid)
  64. ! !INPUT PARAMETERS:
  65. type(GlobalSegMap),intent(in) :: SrcGSMap,DstGSMap ! GSmaps for source and dst
  66. integer,intent(in) :: comm ! local MPI communicator
  67. integer,intent(in) :: compid ! coupler's component ID
  68. !
  69. !EOP ___________________________________________________________________
  70. ! Local variables
  71. character(len=100),parameter :: &
  72. RemapMatrixFile='../../data/t42_to_popx1_c_mat.asc'
  73. ! Loop indicies
  74. integer :: i,j,k,n
  75. ! MPI variables
  76. integer :: nprocs, root, ierr
  77. ! SparseMatrix variables
  78. integer :: mdev
  79. integer :: num_elements, nRows, nColumns
  80. integer, dimension(2) :: src_dims, dst_dims
  81. integer, dimension(:), pointer :: rows, columns
  82. real, dimension(:), pointer :: weights
  83. ! SparseMatrix elements on root
  84. type(SparseMatrix) :: sMat
  85. ! _____________________________________________________________________
  86. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  87. ! INITIALIZATION PHASE
  88. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  89. ! LOCAL RANK AND SIZE
  90. call MPI_COMM_RANK(comm,rank,ierr)
  91. call MPI_COMM_SIZE(comm,nprocs,ierr)
  92. root = 0
  93. if(rank==0) write(6,*) cplname,' init start'
  94. if(rank==0) write(6,*) cplname,' MyID ', compid
  95. if(rank==0) write(6,*) cplname,' Num procs ', nprocs
  96. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  97. ! Read matrix weights for interpolation from a file.
  98. if (rank == root) then
  99. mdev = luavail()
  100. open(mdev, file=trim(RemapMatrixFile), status="old")
  101. read(mdev,*) num_elements
  102. read(mdev,*) src_dims(1), src_dims(2)
  103. read(mdev,*) dst_dims(1), dst_dims(2)
  104. allocate(rows(num_elements), columns(num_elements), &
  105. weights(num_elements), stat=ierr)
  106. do n=1, num_elements
  107. read(mdev,*) rows(n), columns(n), weights(n)
  108. end do
  109. close(mdev)
  110. ! Initialize a Sparsematrix
  111. nRows = dst_dims(1) * dst_dims(2)
  112. nColumns = src_dims(1) * src_dims(2)
  113. call SparseMatrix_init(sMat,nRows,nColumns,num_elements)
  114. call SparseMatrix_importGRowInd(sMat, rows, size(rows))
  115. call SparseMatrix_importGColInd(sMat, columns, size(columns))
  116. call SparseMatrix_importMatrixElts(sMat, weights, size(weights))
  117. deallocate(rows, columns, weights, stat=ierr)
  118. endif
  119. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  120. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  121. ! Build a SparseMatrixPlus for doing the interpolation
  122. ! Specify matrix decomposition to be by row.
  123. ! following the atmosphere's decomposition.
  124. call SparseMatrixPlus_init(Src2DstMatPlus, sMat, SrcGSMap, DstGSMap, &
  125. Xonly, root, comm, compid)
  126. ! no longer need the matrix defined on root
  127. if(rank==0) call SparseMatrix_clean(sMat)
  128. if(rank==0) write(6,*) cplname, ' init done'
  129. !!! END OF INIT !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  130. end subroutine cplinit
  131. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  132. ! RUN PHASE
  133. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  134. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  135. ! Math and Computer Science Division, Argonne National Laboratory !
  136. !BOP -------------------------------------------------------------------
  137. !
  138. ! !IROUTINE: cplrun - coupler's run method
  139. subroutine cplrun(IMPORT,EXPORT)
  140. ! !INPUT PARAMETERS:
  141. type(AttrVect),intent(in) :: IMPORT
  142. type(AttrVect),intent(out) :: EXPORT
  143. !EOP -------------------------------------------------------------------
  144. if(rank==0) write(6,*) cplname, ' run start'
  145. ! Interpolate by doing a parallel sparsematrix-attrvect multiply
  146. ! Note: this will interpolate all fields with the same names
  147. call MCT_MatVecMul(IMPORT, Src2DstMatPlus, EXPORT)
  148. ! possibly do more calculations
  149. if(rank==0) write(6,*) cplname, ' run done'
  150. !!! END OF RUN !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  151. end subroutine cplrun
  152. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  153. ! FINALIZE PHASE
  154. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  155. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  156. ! Math and Computer Science Division, Argonne National Laboratory !
  157. !BOP -------------------------------------------------------------------
  158. !
  159. ! !IROUTINE: cplfin - coupler's finalize method
  160. subroutine cplfin
  161. !
  162. !EOP -------------------------------------------------------------------
  163. call SparseMatrixPlus_clean(Src2DstMatPlus)
  164. if(rank==0) write(6,*) cplname, " done"
  165. end subroutine cplfin
  166. end module coupler