[6331] | 1 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 2 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 3 | !----------------------------------------------------------------------- |
---|
| 4 | ! CVS $Id: twocmp.seq.F90,v 1.6 2006-07-25 17:09:42 jacob Exp $ |
---|
| 5 | ! CVS $Name: $ |
---|
| 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 | |
---|
| 26 | use m_GlobalSegMap,only: GlobalSegMap |
---|
| 27 | use m_GlobalSegMap,only: MCT_GSMap_init => init |
---|
| 28 | use m_GlobalSegMap,only: MCT_GSMap_lsize => lsize |
---|
| 29 | |
---|
| 30 | use m_AttrVect,only : AttrVect |
---|
| 31 | use m_AttrVect,only : MCT_AtrVt_init => init |
---|
| 32 | use m_AttrVect,only : MCT_AtrVt_zero => zero |
---|
| 33 | use m_AttrVect,only : MCT_AtrVt_lsize => lsize |
---|
| 34 | use m_AttrVect,only : MCT_AtrVt_indexRA => indexRA |
---|
| 35 | use m_AttrVect,only : MCT_AtrVt_importRA => importRAttr |
---|
| 36 | |
---|
| 37 | use m_Rearranger,only: Rearranger |
---|
| 38 | use m_Rearranger,only: MCT_Rearranger_init => init |
---|
| 39 | use m_Rearranger,only: MCT_Rearrange => Rearrange |
---|
| 40 | |
---|
| 41 | implicit none |
---|
| 42 | |
---|
| 43 | include 'mpif.h' |
---|
| 44 | |
---|
| 45 | integer,parameter :: ngx = 6 ! points in x-direction |
---|
| 46 | integer,parameter :: ngy = 4 ! points in y-direction |
---|
| 47 | integer ier,nprocs |
---|
| 48 | integer,dimension(:),pointer :: myids |
---|
| 49 | integer :: comm1,comm2,asize,mysize,i,myproc |
---|
| 50 | integer,dimension(1) :: start1,length1 |
---|
| 51 | integer,dimension(:),pointer :: start2,length2 |
---|
| 52 | !----------------------------------------------------------------------- |
---|
| 53 | ! The Main program. |
---|
| 54 | ! We are implementing a single-executable, sequential-execution system. |
---|
| 55 | ! In this example, communication occurs through main using |
---|
| 56 | ! arguments. Both components share the same processors. |
---|
| 57 | |
---|
| 58 | type(GlobalSegMap) :: GSmap1,GSmap2 |
---|
| 59 | type(AttrVect) :: av1,av2 |
---|
| 60 | type(Rearranger) :: Rearr |
---|
| 61 | !----------------------------------------------------------------------- |
---|
| 62 | |
---|
| 63 | call MPI_init(ier) |
---|
| 64 | |
---|
| 65 | call mpi_comm_size(MPI_COMM_WORLD, mysize,ier) |
---|
| 66 | if(mysize .gt. 4) then |
---|
| 67 | write(6,*)"The small problem size in this example & |
---|
| 68 | &requires ", ngy,"or fewer processors." |
---|
| 69 | stop |
---|
| 70 | endif |
---|
| 71 | call mpi_comm_rank(MPI_COMM_WORLD, myproc,ier) |
---|
| 72 | |
---|
| 73 | call mpi_comm_dup(MPI_COMM_WORLD,comm1,ier) |
---|
| 74 | call mpi_comm_dup(MPI_COMM_WORLD,comm2,ier) |
---|
| 75 | |
---|
| 76 | allocate(myids(2)) |
---|
| 77 | myids(1)=1 |
---|
| 78 | myids(2)=2 |
---|
| 79 | |
---|
| 80 | call MCTWorld_init(2,MPI_COMM_WORLD,comm1,myids=myids) |
---|
| 81 | |
---|
| 82 | ! set up a grid and decomposition |
---|
| 83 | ! first gsmap is the grid decomposed by rows |
---|
| 84 | ! theres 1 segment per processor |
---|
| 85 | length1(1)= ngx * (ngy/mysize) |
---|
| 86 | start1(1)= myproc * length1(1) + 1 |
---|
| 87 | |
---|
| 88 | write(6,*)'gsmap1', myproc,length1(1),start1(1) |
---|
| 89 | call MCT_GSMap_init(GSMap1,start1,length1,0,comm1,1) |
---|
| 90 | |
---|
| 91 | ! second gsmap is the grid decomposed by columns |
---|
| 92 | allocate(length2(ngy),start2(ngy)) |
---|
| 93 | |
---|
| 94 | do i=1,ngy |
---|
| 95 | length2(i)=ngx/mysize |
---|
| 96 | start2(i)= (i-1)*ngx + 1 + myproc*length2(i) |
---|
| 97 | write(6,*) 'gsmap2',myproc,i,length2(i),start2(i) |
---|
| 98 | enddo |
---|
| 99 | |
---|
| 100 | |
---|
| 101 | call MCT_GSMap_init(GSMap2,start2,length2,0,comm2,2) |
---|
| 102 | |
---|
| 103 | call MCT_AtrVt_init(av1,rList="field1:field2",lsize=MCT_GSMap_lsize(GSMap1,comm1)) |
---|
| 104 | |
---|
| 105 | call MCT_AtrVt_init(av2,rList="field1:field2",lsize=MCT_GSMap_lsize(GSMap2,comm2)) |
---|
| 106 | |
---|
| 107 | |
---|
| 108 | ! create a rearranger |
---|
| 109 | call MCT_Rearranger_init(GSMap1,GSMap2,MPI_COMM_WORLD,Rearr) |
---|
| 110 | |
---|
| 111 | !-------------end of initialization steps |
---|
| 112 | |
---|
| 113 | |
---|
| 114 | ! Start up model1 which fills av1 with data. |
---|
| 115 | call model1(comm1,av1) |
---|
| 116 | |
---|
| 117 | ! print out Av data |
---|
| 118 | do i=1,MCT_AtrVt_lsize(av1) |
---|
| 119 | write(6,*) "model 1 data", myproc,i,av1%rAttr(1,i),av1%rAttr(2,i) |
---|
| 120 | enddo |
---|
| 121 | |
---|
| 122 | ! rearrange data from model1 so that model2 can use it. |
---|
| 123 | call MCT_Rearrange(av1,av2,Rearr) |
---|
| 124 | |
---|
| 125 | ! pass data to model2 (which will print it out) |
---|
| 126 | call model2(comm2,av2) |
---|
| 127 | |
---|
| 128 | |
---|
| 129 | ! all done |
---|
| 130 | call mpi_finalize(ier) |
---|
| 131 | |
---|
| 132 | contains |
---|
| 133 | |
---|
| 134 | !----------------------------------------------------------------------- |
---|
| 135 | !----------------------------------------------------------------------- |
---|
| 136 | ! !ROUTINE: |
---|
| 137 | subroutine model1(comm1,mod1av) ! the first model |
---|
| 138 | |
---|
| 139 | implicit none |
---|
| 140 | |
---|
| 141 | integer :: comm1,mysize,ier,asize,myproc |
---|
| 142 | integer :: fieldindx,avsize,i |
---|
| 143 | integer,dimension(1) :: start,length |
---|
| 144 | real,pointer :: testarray(:) |
---|
| 145 | |
---|
| 146 | type(GlobalSegMap) :: GSmap |
---|
| 147 | type(AttrVect) :: mod1av |
---|
| 148 | !--------------------------- |
---|
| 149 | |
---|
| 150 | ! find local rank and size |
---|
| 151 | call mpi_comm_size(comm1,mysize,ier) |
---|
| 152 | call mpi_comm_rank(comm1,myproc,ier) |
---|
| 153 | write(6,*)"model1 size",mysize |
---|
| 154 | |
---|
| 155 | |
---|
| 156 | avsize = MCT_AtrVt_lsize(mod1av) |
---|
| 157 | write(6,*)"model 1 av size", avsize |
---|
| 158 | |
---|
| 159 | ! Fill Av with some data |
---|
| 160 | ! fill first attribute the direct way |
---|
| 161 | fieldindx = MCT_AtrVt_indexRA(mod1av,"field1") |
---|
| 162 | do i=1,avsize |
---|
| 163 | mod1av%rAttr(fieldindx,i) = float(i+ 20*myproc) |
---|
| 164 | enddo |
---|
| 165 | |
---|
| 166 | ! fill second attribute using Av import function |
---|
| 167 | allocate(testarray(avsize)) |
---|
| 168 | do i=1,avsize |
---|
| 169 | testarray(i)= cos((float(i+ 20*myproc)/24.) * 3.14) |
---|
| 170 | enddo |
---|
| 171 | call MCT_AtrVt_importRA(mod1av,"field2",testarray) |
---|
| 172 | |
---|
| 173 | |
---|
| 174 | end subroutine model1 |
---|
| 175 | |
---|
| 176 | !----------------------------------------------------------------------- |
---|
| 177 | !----------------------------------------------------------------------- |
---|
| 178 | ! !ROUTINE: |
---|
| 179 | subroutine model2(comm2,mod2av) |
---|
| 180 | |
---|
| 181 | implicit none |
---|
| 182 | |
---|
| 183 | integer :: comm2,mysize,ier,asize,myproc |
---|
| 184 | integer :: i |
---|
| 185 | type(AttrVect) :: mod2av |
---|
| 186 | !--------------------------- |
---|
| 187 | |
---|
| 188 | ! find local rank and size |
---|
| 189 | call mpi_comm_size(comm2,mysize,ier) |
---|
| 190 | call mpi_comm_rank(comm2,myproc,ier) |
---|
| 191 | write(6,*)"model2 size",mysize |
---|
| 192 | |
---|
| 193 | asize = MCT_AtrVt_lsize(mod2av) |
---|
| 194 | write(6,*)"model 2 av size", asize |
---|
| 195 | |
---|
| 196 | ! print out Av data |
---|
| 197 | do i=1,asize |
---|
| 198 | write(6,*) "model 2 data after", myproc,i,mod2av%rAttr(1,i),mod2av%rAttr(2,i) |
---|
| 199 | enddo |
---|
| 200 | |
---|
| 201 | |
---|
| 202 | end subroutine model2 |
---|
| 203 | |
---|
| 204 | end |
---|