model.F90 6.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198
  1. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  2. ! Math and Computer Science Division, Argonne National Laboratory !
  3. !-----------------------------------------------------------------------
  4. ! CVS model.F90,v 1.8 2004-04-23 20:56:23 jacob Exp
  5. ! CVS MCT_2_8_0
  6. !BOP -------------------------------------------------------------------
  7. !
  8. ! !ROUTINE: model -- generic model for unit tester
  9. !
  10. ! !DESCRIPTION:
  11. ! A generic model subroutine to test functionality of MCT.
  12. !
  13. ! !INTERFACE:
  14. !
  15. subroutine model (comm,ncomps,compid)
  16. !
  17. ! !USES:
  18. !
  19. ! Get the things needed from MCT by "Use,only" with renaming:
  20. !
  21. !---Component Model Registry
  22. use m_MCTWorld,only: MCTWorld_init => init
  23. use m_MCTWorld,only: MCTWorld_clean => clean
  24. !---Domain Decomposition Descriptor DataType and associated methods
  25. use m_GlobalSegMap,only: GlobalSegMap
  26. use m_GlobalSegMap,only: GlobalSegMap_init => init
  27. use m_GlobalSegMap,only: GlobalSegMap_lsize => lsize
  28. use m_GlobalSegMap,only: GlobalSegMap_clean => clean
  29. use m_GlobalSegMap,only: GlobalSegMap_Ordpnts => OrderedPoints
  30. !---Field Storage DataType and associated methods
  31. use m_AttrVect,only : AttrVect
  32. use m_AttrVect,only : AttrVect_init => init
  33. use m_AttrVect,only : AttrVect_clean => clean
  34. use m_AttrVect,only : AttrVect_indxR => indexRA
  35. use m_AttrVect,only : AttrVect_importRAttr => importRAttr
  36. !---Intercomponent communications scheduler
  37. use m_Router,only: Router
  38. use m_Router,only: Router_init => init
  39. use m_Router,only: Router_clean => clean
  40. !---Intercomponent transfer
  41. use m_Transfer,only : MCT_Send => send
  42. use m_Transfer,only : MCT_Recv => recv
  43. !---Stored Grid data
  44. implicit none
  45. include "mpif.h"
  46. ! !INPUT PARAMETERS:
  47. integer,intent(in) :: comm ! MPI communicator for this component
  48. integer,intent(in) :: ncomps ! total number of models in coupled system
  49. integer,intent(in) :: compid ! the integer id of this model
  50. !
  51. !EOP ___________________________________________________________________
  52. ! local variables
  53. ! parameters for this model
  54. character(len=*), parameter :: modelname='model.F90'
  55. integer,parameter :: nxa = 128 ! number of points in x-direction
  56. integer,parameter :: nya = 64 ! number of points in y-direction
  57. integer :: i,j,k
  58. ! note decleration of instances of MCT defined types.
  59. ! MPI variables
  60. integer :: rank, nprocs, root, CplID, ierr
  61. ! Grid variables
  62. integer :: localsize
  63. ! GlobalSegMap variables
  64. type(GlobalSegMap) :: GSMap ! MCT defined type
  65. integer,dimension(1) :: start,length
  66. integer, dimension(:), pointer :: points
  67. ! AttrVect variables
  68. type(AttrVect) :: AV ! MCT defined type
  69. real, dimension(:), pointer :: avdata
  70. integer :: avsize
  71. ! Router variables
  72. type(Router) :: Rout ! MCT defined type
  73. ! _____________________________________________________________________
  74. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  75. ! INITIALIZATION PHASE
  76. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  77. ! Get local rank and size
  78. call MPI_COMM_RANK (comm,rank, ierr)
  79. call MPI_COMM_SIZE(comm,nprocs,ierr)
  80. root = 0
  81. if(rank==0) write(6,*) modelname,' MyID ', compid
  82. if(rank==0) write(6,*) modelname,' Num procs ', nprocs
  83. ! Initialize MCTworld
  84. call MCTWorld_init(ncomps,MPI_COMM_WORLD,comm,compid)
  85. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  86. ! Initialize a Global Segment Map
  87. ! set up a 1-d decomposition.
  88. ! there is just 1 segment per processor
  89. localsize = nxa*nya / nprocs
  90. ! we'll use the distributed init of GSMap so
  91. ! initialize start and length arrays for this processor
  92. start(1) = (rank*localsize) + 1
  93. length(1) = localsize
  94. ! initialize the GSMap
  95. call GlobalSegMap_init(GSMap,start,length,root,comm,compid)
  96. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  97. ! Use a GSMap function:
  98. ! return the points local to this processor
  99. ! in their assumed order.
  100. call GlobalSegMap_Ordpnts(GSMap,rank,points)
  101. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  102. ! Initialize an Attribute vector
  103. ! size is the number of grid point on this processor
  104. avsize = GlobalSegMap_lsize(GSMap,comm)
  105. if(rank==0) write(6,*) modelname, ' localsize ', avsize
  106. ! initialize Av with two real attributes.
  107. call AttrVect_init(AV,rList="field1:field2",lsize=avsize)
  108. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  109. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  110. ! Initialize a router to the coupler component.
  111. !
  112. ! Need to know the integer ID of the coupler.
  113. CplID = 2
  114. call Router_init(CplID,GSMap,comm,Rout)
  115. ! create an array used in RUN
  116. allocate(avdata(avsize),stat=ierr)
  117. !!! END OF INIT !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  118. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  119. ! RUN PHASE
  120. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  121. do j=1,10 ! "timestep" loop
  122. ! model calculations
  123. ! load data into aV
  124. ! load the first field using "import" method.
  125. ! First field will be a constant real number.
  126. avdata=30.0
  127. call AttrVect_importRAttr(AV,"field1",avdata)
  128. ! Load the second field using direct access
  129. ! Second field will be the indicies of each grid point
  130. ! in the grid point numbering scheme.
  131. do i=1,avsize
  132. AV%rAttr(AttrVect_indxR(AV,"field2"),i) = points(i)
  133. enddo
  134. ! Send the data
  135. ! this is a synchronization point between the coupler and
  136. ! this model.
  137. if(rank==0) write(6,*) modelname,' sending data step ',j
  138. call MCT_Send(AV,Rout)
  139. ! more model calculations
  140. enddo
  141. !!! END OF RUN !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  142. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  143. ! FINALIZE PHASE
  144. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  145. ! clean up
  146. call Router_clean(Rout)
  147. call AttrVect_clean(AV)
  148. call GlobalSegMap_clean(GSMap)
  149. call MCTWorld_clean()
  150. if(rank==0) write(6,*) modelname,' done'
  151. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  152. end subroutine model