srcmodel.F90 7.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248
  1. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  2. ! Math and Computer Science Division, Argonne National Laboratory !
  3. !-----------------------------------------------------------------------
  4. ! CVS srcmodel.F90,v 1.8 2005-11-18 23:15:38 rloy Exp
  5. ! CVS MCT_2_8_0
  6. !BOP -------------------------------------------------------------------
  7. !
  8. ! !MODULE: srcmodel -- generic model for unit tester
  9. !
  10. ! !DESCRIPTION:
  11. ! init run and finalize methods for source model
  12. !
  13. module srcmodel
  14. !
  15. ! !USES:
  16. !
  17. ! Get the things needed from MCT by "Use,only" with renaming:
  18. !
  19. !---Domain Decomposition Descriptor DataType and associated methods
  20. use m_GlobalSegMap,only: GlobalSegMap
  21. use m_GlobalSegMap,only: GlobalSegMap_init => init
  22. use m_GlobalSegMap,only: GlobalSegMap_lsize => lsize
  23. use m_GlobalSegMap,only: GlobalSegMap_clean => clean
  24. !---Field Storage DataType and associated methods
  25. use m_AttrVect,only : AttrVect
  26. use m_AttrVect,only : AttrVect_init => init
  27. use m_AttrVect,only : AttrVect_lsize => lsize
  28. use m_AttrVect,only : AttrVect_clean => clean
  29. use m_AttrVect,only : AttrVect_copy => copy
  30. use m_AttrVect,only : AttrVect_zero => zero
  31. use m_AttrVect,only : AttrVect_indxR => indexRA
  32. use m_AttrVect,only : AttrVect_importRAttr => importRAttr
  33. use m_AttrVectComms,only : AttrVect_scatter => scatter
  34. ! Get things from MPEU
  35. use m_inpak90 ! Resource files
  36. use m_stdio ! I/O utils
  37. use m_ioutil
  38. ! Get utilities for this program.
  39. use mutils
  40. implicit none
  41. private
  42. ! except
  43. ! !PUBLIC MEMBER FUNCTIONS:
  44. public srcinit
  45. public srcrun
  46. public srcfin
  47. ! private module variables
  48. character(len=*), parameter :: modelname='srcmodel.F90'
  49. integer :: rank
  50. real, dimension(:), pointer :: avdata
  51. !EOP -------------------------------------------------------------------
  52. contains
  53. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  54. ! Math and Computer Science Division, Argonne National Laboratory !
  55. !BOP -------------------------------------------------------------------
  56. !
  57. ! !IROUTINE: srcinit - Source model initialization
  58. subroutine srcinit(GSMap,IMPORT,EXPORT,comm,compid)
  59. ! !INPUT PARAMETERS:
  60. type(GlobalSegMap),intent(inout) :: GSMap ! decomposition
  61. type(AttrVect),intent(inout) :: IMPORT,EXPORT ! state data
  62. integer,intent(in) :: comm ! MPI communicator
  63. integer,intent(in) :: compid ! component ID
  64. !
  65. !EOP ___________________________________________________________________
  66. ! local variables
  67. ! parameters for this model
  68. integer :: nxa ! number of points in x-direction
  69. integer :: nya ! number of points in y-direction
  70. integer :: i,j,k,mdev,fx,fy
  71. integer :: nprocs, root, ier,fileno
  72. ! GlobalSegMap variables
  73. integer,dimension(:),pointer :: lindex
  74. ! AttrVect variables
  75. integer :: avsize
  76. type(AttrVect) :: GlobalD ! Av to hold global data
  77. real,dimension(:),pointer :: rootdata
  78. character*2 :: ldecomp
  79. call MPI_COMM_RANK(comm,rank, ier)
  80. call MPI_COMM_SIZE(comm,nprocs,ier)
  81. if(rank==0) then
  82. write(6,*) modelname, ' init start'
  83. write(6,*) modelname,' MyID ', compid
  84. write(6,*) modelname,' Num procs ', nprocs
  85. endif
  86. ! Get configuration
  87. call i90_LoadF('src.rc',ier)
  88. call i90_label('nx:',ier)
  89. nxa=i90_gint(ier)
  90. call i90_label('ny:',ier)
  91. nya=i90_gint(ier)
  92. if(rank==0) write(6,*) modelname, ' x,y ', nxa,nya
  93. call i90_label('decomp:',ier)
  94. call i90_Gtoken(ldecomp, ier)
  95. if(rank==0) write(6,*) modelname, ' decomp ', ldecomp
  96. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  97. ! Initialize a Global Segment Map
  98. call get_index(ldecomp,nprocs,rank,nxa,nya,lindex)
  99. call GlobalSegMap_init(GSMap,lindex,comm,compid,gsize=nxa*nya)
  100. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  101. if(rank==0) write(6,*) modelname, ' GSMap ',GSMap%ngseg,GSMap%gsize
  102. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  103. ! Initialize import and export Attribute vectors
  104. ! size is the number of grid points on this processor
  105. avsize = GlobalSegMap_lsize(GSMap,comm)
  106. if(rank==0) write(6,*) modelname, ' localsize ', avsize
  107. ! Initialize the IMPORT Av by scattering from a root Av
  108. ! with real data.
  109. ! Read in data from root and scatter to nodes
  110. if(rank==0) then
  111. call AttrVect_init(GlobalD,rList="field1:field2",lsize=nxa*nya)
  112. mdev=luavail()
  113. open(mdev, file="TS1.dat",status="old")
  114. read(mdev,*) fx,fy
  115. do i=1,nxa*nya
  116. read(mdev,*) GlobalD%rAttr(1,i)
  117. enddo
  118. write(6,*) modelname,'Global init ',GlobalD%rAttr(1,1),GlobalD%rAttr(1,8000)
  119. endif
  120. ! this scatter will create IMPORT if it hasn't already been initialized
  121. call AttrVect_scatter(GlobalD,IMPORT,GSMap,0,comm,ier)
  122. ! initialize EXPORT Av with two real attributes.
  123. call AttrVect_init(EXPORT,rList="field3:field4",lsize=avsize)
  124. call AttrVect_zero(EXPORT)
  125. if(rank==0) then
  126. write(6,*) modelname, rank,' IMPORT field1', IMPORT%rAttr(1,1)
  127. write(6,*) modelname, rank,' IMPORt field2', IMPORT%rAttr(2,1)
  128. write(6,*) modelname, rank,' EXPORT field3', EXPORT%rAttr(1,1)
  129. write(6,*) modelname, rank,' EXPORT field4', EXPORT%rAttr(2,1)
  130. endif
  131. ! allocate buffer for use in run method
  132. allocate(avdata(avsize),stat=ier)
  133. if(rank==0) write(6,*) modelname, ' init done'
  134. end subroutine srcinit
  135. !!! END OF INIT !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  136. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  137. ! RUN PHASE
  138. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  139. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  140. ! Math and Computer Science Division, Argonne National Laboratory !
  141. !BOP -------------------------------------------------------------------
  142. !
  143. ! !IROUTINE: srcrun - Source model run method
  144. subroutine srcrun(IMPORT,EXPORT)
  145. ! !INPUT PARAMETERS:
  146. type(AttrVect),intent(inout) :: IMPORT,EXPORT ! Input and Output states
  147. !EOP -------------------------------------------------------------------
  148. ! local variables
  149. integer :: avsize,ier,i
  150. ! Nothing to do with IMPORT
  151. ! Fill EXPORT with data
  152. if(rank==0) write(6,*) modelname, ' run start'
  153. ! Use Av copy to copy input data from field1 in Imp to field3 in EXPORT
  154. call AttrVect_copy(IMPORT,EXPORT,rList='field1',TrList='field3')
  155. ! Use import to load data in second field
  156. avdata=30.0
  157. call AttrVect_importRAttr(EXPORT,"field4",avdata)
  158. if(rank==0) write(6,*) modelname, ' In field1', IMPORT%rAttr(1,1)
  159. if(rank==0) write(6,*) modelname, ' In field2', IMPORT%rAttr(2,1)
  160. if(rank==0) write(6,*) modelname, ' Out field3', EXPORT%rAttr(1,1)
  161. if(rank==0) write(6,*) modelname, ' Out field4', EXPORT%rAttr(2,1)
  162. if(rank==0) write(6,*) modelname, ' run done'
  163. end subroutine srcrun
  164. !!! END OF RUN !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  165. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  166. ! FINALIZE PHASE
  167. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  168. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  169. ! Math and Computer Science Division, Argonne National Laboratory !
  170. !BOP -------------------------------------------------------------------
  171. !
  172. ! !IROUTINE: srcfin - Source model finalize method
  173. subroutine srcfin(IMPORT,EXPORT,GSMap)
  174. ! !INPUT PARAMETERS:
  175. type(AttrVect),intent(inout) :: IMPORT,EXPORT ! imp,exp states
  176. type(GlobalSegMap),intent(inout) :: GSMap
  177. !EOP -------------------------------------------------------------------
  178. ! clean up
  179. call AttrVect_clean(IMPORT)
  180. call AttrVect_clean(EXPORT)
  181. call GlobalSegMap_clean(GSMap)
  182. deallocate(avdata)
  183. if(rank==0) write(6,*) modelname,' fin done'
  184. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  185. endsubroutine srcfin
  186. end module srcmodel