1 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
2 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
3 | !----------------------------------------------------------------------- |
---|
4 | ! CVS $Id: model.F90,v 1.8 2004-04-23 20:56:23 jacob Exp $ |
---|
5 | ! CVS $Name: $ |
---|
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 | |
---|
45 | implicit none |
---|
46 | |
---|
47 | include "mpif.h" |
---|
48 | |
---|
49 | ! !INPUT PARAMETERS: |
---|
50 | |
---|
51 | integer,intent(in) :: comm ! MPI communicator for this component |
---|
52 | integer,intent(in) :: ncomps ! total number of models in coupled system |
---|
53 | integer,intent(in) :: compid ! the integer id of this model |
---|
54 | ! |
---|
55 | !EOP ___________________________________________________________________ |
---|
56 | |
---|
57 | ! local variables |
---|
58 | |
---|
59 | ! parameters for this model |
---|
60 | character(len=*), parameter :: modelname='model.F90' |
---|
61 | integer,parameter :: nxa = 128 ! number of points in x-direction |
---|
62 | integer,parameter :: nya = 64 ! number of points in y-direction |
---|
63 | |
---|
64 | integer :: i,j,k |
---|
65 | |
---|
66 | ! note decleration of instances of MCT defined types. |
---|
67 | ! MPI variables |
---|
68 | integer :: rank, nprocs, root, CplID, ierr |
---|
69 | ! Grid variables |
---|
70 | integer :: localsize |
---|
71 | ! GlobalSegMap variables |
---|
72 | type(GlobalSegMap) :: GSMap ! MCT defined type |
---|
73 | integer,dimension(1) :: start,length |
---|
74 | integer, dimension(:), pointer :: points |
---|
75 | ! AttrVect variables |
---|
76 | type(AttrVect) :: AV ! MCT defined type |
---|
77 | real, dimension(:), pointer :: avdata |
---|
78 | integer :: avsize |
---|
79 | ! Router variables |
---|
80 | type(Router) :: Rout ! MCT defined type |
---|
81 | ! _____________________________________________________________________ |
---|
82 | |
---|
83 | |
---|
84 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
85 | ! INITIALIZATION PHASE |
---|
86 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
87 | |
---|
88 | ! Get local rank and size |
---|
89 | call MPI_COMM_RANK (comm,rank, ierr) |
---|
90 | call MPI_COMM_SIZE(comm,nprocs,ierr) |
---|
91 | root = 0 |
---|
92 | |
---|
93 | if(rank==0) write(6,*) modelname,' MyID ', compid |
---|
94 | if(rank==0) write(6,*) modelname,' Num procs ', nprocs |
---|
95 | |
---|
96 | ! Initialize MCTworld |
---|
97 | call MCTWorld_init(ncomps,MPI_COMM_WORLD,comm,compid) |
---|
98 | |
---|
99 | |
---|
100 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
101 | ! Initialize a Global Segment Map |
---|
102 | |
---|
103 | ! set up a 1-d decomposition. |
---|
104 | ! there is just 1 segment per processor |
---|
105 | localsize = nxa*nya / nprocs |
---|
106 | |
---|
107 | ! we'll use the distributed init of GSMap so |
---|
108 | ! initialize start and length arrays for this processor |
---|
109 | start(1) = (rank*localsize) + 1 |
---|
110 | length(1) = localsize |
---|
111 | |
---|
112 | ! initialize the GSMap |
---|
113 | call GlobalSegMap_init(GSMap,start,length,root,comm,compid) |
---|
114 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
115 | |
---|
116 | |
---|
117 | ! Use a GSMap function: |
---|
118 | ! return the points local to this processor |
---|
119 | ! in their assumed order. |
---|
120 | call GlobalSegMap_Ordpnts(GSMap,rank,points) |
---|
121 | |
---|
122 | |
---|
123 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
124 | ! Initialize an Attribute vector |
---|
125 | |
---|
126 | ! size is the number of grid point on this processor |
---|
127 | avsize = GlobalSegMap_lsize(GSMap,comm) |
---|
128 | if(rank==0) write(6,*) modelname, ' localsize ', avsize |
---|
129 | |
---|
130 | ! initialize Av with two real attributes. |
---|
131 | call AttrVect_init(AV,rList="field1:field2",lsize=avsize) |
---|
132 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
133 | |
---|
134 | |
---|
135 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
136 | ! Initialize a router to the coupler component. |
---|
137 | ! |
---|
138 | ! Need to know the integer ID of the coupler. |
---|
139 | CplID = 2 |
---|
140 | call Router_init(CplID,GSMap,comm,Rout) |
---|
141 | |
---|
142 | ! create an array used in RUN |
---|
143 | allocate(avdata(avsize),stat=ierr) |
---|
144 | !!! END OF INIT !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
145 | |
---|
146 | |
---|
147 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
148 | ! RUN PHASE |
---|
149 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
150 | |
---|
151 | |
---|
152 | do j=1,10 ! "timestep" loop |
---|
153 | |
---|
154 | |
---|
155 | ! model calculations |
---|
156 | |
---|
157 | |
---|
158 | ! load data into aV |
---|
159 | ! load the first field using "import" method. |
---|
160 | ! First field will be a constant real number. |
---|
161 | avdata=30.0 |
---|
162 | call AttrVect_importRAttr(AV,"field1",avdata) |
---|
163 | |
---|
164 | ! Load the second field using direct access |
---|
165 | ! Second field will be the indicies of each grid point |
---|
166 | ! in the grid point numbering scheme. |
---|
167 | do i=1,avsize |
---|
168 | AV%rAttr(AttrVect_indxR(AV,"field2"),i) = points(i) |
---|
169 | enddo |
---|
170 | |
---|
171 | ! Send the data |
---|
172 | ! this is a synchronization point between the coupler and |
---|
173 | ! this model. |
---|
174 | if(rank==0) write(6,*) modelname,' sending data step ',j |
---|
175 | call MCT_Send(AV,Rout) |
---|
176 | |
---|
177 | |
---|
178 | ! more model calculations |
---|
179 | |
---|
180 | |
---|
181 | enddo |
---|
182 | |
---|
183 | !!! END OF RUN !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
184 | |
---|
185 | |
---|
186 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
187 | ! FINALIZE PHASE |
---|
188 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
189 | ! clean up |
---|
190 | call Router_clean(Rout) |
---|
191 | call AttrVect_clean(AV) |
---|
192 | call GlobalSegMap_clean(GSMap) |
---|
193 | call MCTWorld_clean() |
---|
194 | if(rank==0) write(6,*) modelname,' done' |
---|
195 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
196 | |
---|
197 | end subroutine model |
---|
198 | |
---|