twocmp.seq.F90 6.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204
  1. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  2. ! Math and Computer Science Division, Argonne National Laboratory !
  3. !-----------------------------------------------------------------------
  4. ! CVS twocmp.seq.F90,v 1.6 2006-07-25 17:09:42 jacob Exp
  5. ! CVS MCT_2_8_0
  6. !BOP -------------------------------------------------------------------
  7. !
  8. ! !ROUTINE: twocomponent.sequential
  9. !
  10. !
  11. ! !DESCRIPTION: Provide a simple example of using MCT to connect
  12. ! two components executing in sequence in a single executable.
  13. !
  14. ! Data is passed between models by using input/output arguments
  15. ! in the run method. Compare with twocmp.seqNB.F90
  16. !
  17. ! !INTERFACE:
  18. !
  19. program twoseq
  20. !
  21. ! !USES:
  22. !
  23. !--- Get only the things needed from MCT
  24. use m_MCTWorld,only: MCTWorld_init => init
  25. use m_GlobalSegMap,only: GlobalSegMap
  26. use m_GlobalSegMap,only: MCT_GSMap_init => init
  27. use m_GlobalSegMap,only: MCT_GSMap_lsize => lsize
  28. use m_AttrVect,only : AttrVect
  29. use m_AttrVect,only : MCT_AtrVt_init => init
  30. use m_AttrVect,only : MCT_AtrVt_zero => zero
  31. use m_AttrVect,only : MCT_AtrVt_lsize => lsize
  32. use m_AttrVect,only : MCT_AtrVt_indexRA => indexRA
  33. use m_AttrVect,only : MCT_AtrVt_importRA => importRAttr
  34. use m_Rearranger,only: Rearranger
  35. use m_Rearranger,only: MCT_Rearranger_init => init
  36. use m_Rearranger,only: MCT_Rearrange => Rearrange
  37. implicit none
  38. include 'mpif.h'
  39. integer,parameter :: ngx = 6 ! points in x-direction
  40. integer,parameter :: ngy = 4 ! points in y-direction
  41. integer ier,nprocs
  42. integer,dimension(:),pointer :: myids
  43. integer :: comm1,comm2,asize,mysize,i,myproc
  44. integer,dimension(1) :: start1,length1
  45. integer,dimension(:),pointer :: start2,length2
  46. !-----------------------------------------------------------------------
  47. ! The Main program.
  48. ! We are implementing a single-executable, sequential-execution system.
  49. ! In this example, communication occurs through main using
  50. ! arguments. Both components share the same processors.
  51. type(GlobalSegMap) :: GSmap1,GSmap2
  52. type(AttrVect) :: av1,av2
  53. type(Rearranger) :: Rearr
  54. !-----------------------------------------------------------------------
  55. call MPI_init(ier)
  56. call mpi_comm_size(MPI_COMM_WORLD, mysize,ier)
  57. if(mysize .gt. 4) then
  58. write(6,*)"The small problem size in this example &
  59. &requires ", ngy,"or fewer processors."
  60. stop
  61. endif
  62. call mpi_comm_rank(MPI_COMM_WORLD, myproc,ier)
  63. call mpi_comm_dup(MPI_COMM_WORLD,comm1,ier)
  64. call mpi_comm_dup(MPI_COMM_WORLD,comm2,ier)
  65. allocate(myids(2))
  66. myids(1)=1
  67. myids(2)=2
  68. call MCTWorld_init(2,MPI_COMM_WORLD,comm1,myids=myids)
  69. ! set up a grid and decomposition
  70. ! first gsmap is the grid decomposed by rows
  71. ! theres 1 segment per processor
  72. length1(1)= ngx * (ngy/mysize)
  73. start1(1)= myproc * length1(1) + 1
  74. write(6,*)'gsmap1', myproc,length1(1),start1(1)
  75. call MCT_GSMap_init(GSMap1,start1,length1,0,comm1,1)
  76. ! second gsmap is the grid decomposed by columns
  77. allocate(length2(ngy),start2(ngy))
  78. do i=1,ngy
  79. length2(i)=ngx/mysize
  80. start2(i)= (i-1)*ngx + 1 + myproc*length2(i)
  81. write(6,*) 'gsmap2',myproc,i,length2(i),start2(i)
  82. enddo
  83. call MCT_GSMap_init(GSMap2,start2,length2,0,comm2,2)
  84. call MCT_AtrVt_init(av1,rList="field1:field2",lsize=MCT_GSMap_lsize(GSMap1,comm1))
  85. call MCT_AtrVt_init(av2,rList="field1:field2",lsize=MCT_GSMap_lsize(GSMap2,comm2))
  86. ! create a rearranger
  87. call MCT_Rearranger_init(GSMap1,GSMap2,MPI_COMM_WORLD,Rearr)
  88. !-------------end of initialization steps
  89. ! Start up model1 which fills av1 with data.
  90. call model1(comm1,av1)
  91. ! print out Av data
  92. do i=1,MCT_AtrVt_lsize(av1)
  93. write(6,*) "model 1 data", myproc,i,av1%rAttr(1,i),av1%rAttr(2,i)
  94. enddo
  95. ! rearrange data from model1 so that model2 can use it.
  96. call MCT_Rearrange(av1,av2,Rearr)
  97. ! pass data to model2 (which will print it out)
  98. call model2(comm2,av2)
  99. ! all done
  100. call mpi_finalize(ier)
  101. contains
  102. !-----------------------------------------------------------------------
  103. !-----------------------------------------------------------------------
  104. ! !ROUTINE:
  105. subroutine model1(comm1,mod1av) ! the first model
  106. implicit none
  107. integer :: comm1,mysize,ier,asize,myproc
  108. integer :: fieldindx,avsize,i
  109. integer,dimension(1) :: start,length
  110. real,pointer :: testarray(:)
  111. type(GlobalSegMap) :: GSmap
  112. type(AttrVect) :: mod1av
  113. !---------------------------
  114. ! find local rank and size
  115. call mpi_comm_size(comm1,mysize,ier)
  116. call mpi_comm_rank(comm1,myproc,ier)
  117. write(6,*)"model1 size",mysize
  118. avsize = MCT_AtrVt_lsize(mod1av)
  119. write(6,*)"model 1 av size", avsize
  120. ! Fill Av with some data
  121. ! fill first attribute the direct way
  122. fieldindx = MCT_AtrVt_indexRA(mod1av,"field1")
  123. do i=1,avsize
  124. mod1av%rAttr(fieldindx,i) = float(i+ 20*myproc)
  125. enddo
  126. ! fill second attribute using Av import function
  127. allocate(testarray(avsize))
  128. do i=1,avsize
  129. testarray(i)= cos((float(i+ 20*myproc)/24.) * 3.14)
  130. enddo
  131. call MCT_AtrVt_importRA(mod1av,"field2",testarray)
  132. end subroutine model1
  133. !-----------------------------------------------------------------------
  134. !-----------------------------------------------------------------------
  135. ! !ROUTINE:
  136. subroutine model2(comm2,mod2av)
  137. implicit none
  138. integer :: comm2,mysize,ier,asize,myproc
  139. integer :: i
  140. type(AttrVect) :: mod2av
  141. !---------------------------
  142. ! find local rank and size
  143. call mpi_comm_size(comm2,mysize,ier)
  144. call mpi_comm_rank(comm2,myproc,ier)
  145. write(6,*)"model2 size",mysize
  146. asize = MCT_AtrVt_lsize(mod2av)
  147. write(6,*)"model 2 av size", asize
  148. ! print out Av data
  149. do i=1,asize
  150. write(6,*) "model 2 data after", myproc,i,mod2av%rAttr(1,i),mod2av%rAttr(2,i)
  151. enddo
  152. end subroutine model2
  153. end