[6328] | 1 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 2 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 3 | !----------------------------------------------------------------------- |
---|
| 4 | ! CVS $Id: dstmodel.F90,v 1.8 2006-10-17 21:47:56 jacob Exp $ |
---|
| 5 | ! CVS $Name: $ |
---|
| 6 | !BOP ------------------------------------------------------------------- |
---|
| 7 | ! |
---|
| 8 | ! !MODULE: dstmodel -- generic model for sequential climate model |
---|
| 9 | ! |
---|
| 10 | ! !DESCRIPTION: |
---|
| 11 | ! init run and finalize methods for destination model |
---|
| 12 | ! |
---|
| 13 | ! !INTERFACE: |
---|
| 14 | ! |
---|
| 15 | module dstmodel |
---|
| 16 | |
---|
| 17 | ! |
---|
| 18 | ! !USES: |
---|
| 19 | ! |
---|
| 20 | ! Get the things needed from MCT by "Use,only" with renaming: |
---|
| 21 | ! |
---|
| 22 | !---Domain Decomposition Descriptor DataType and associated methods |
---|
| 23 | use m_GlobalSegMap,only: GlobalSegMap |
---|
| 24 | use m_GlobalSegMap,only: GlobalSegMap_init => init |
---|
| 25 | use m_GlobalSegMap,only: GlobalSegMap_lsize => lsize |
---|
| 26 | use m_GlobalSegMap,only: GlobalSegMap_clean => clean |
---|
| 27 | !---Field Storage DataType and associated methods |
---|
| 28 | use m_AttrVect,only : AttrVect |
---|
| 29 | use m_AttrVect,only : AttrVect_init => init |
---|
| 30 | use m_AttrVect,only : AttrVect_lsize => lsize |
---|
| 31 | use m_AttrVect,only : AttrVect_clean => clean |
---|
| 32 | use m_AttrVect,only : AttrVect_copy => copy |
---|
| 33 | use m_AttrVect,only : AttrVect_indxR => indexRA |
---|
| 34 | use m_AttrVect,only : AttrVect_importRAttr => importRAttr |
---|
| 35 | use m_AttrVectcomms,only : AttrVect_gather => gather |
---|
| 36 | |
---|
| 37 | ! Get things from MPEU |
---|
| 38 | use m_inpak90 ! Resource files |
---|
| 39 | use m_stdio ! I/O utils |
---|
| 40 | use m_ioutil |
---|
| 41 | |
---|
| 42 | |
---|
| 43 | ! Get utilities for this program. |
---|
| 44 | use mutils |
---|
| 45 | |
---|
| 46 | implicit none |
---|
| 47 | |
---|
| 48 | private |
---|
| 49 | ! except |
---|
| 50 | |
---|
| 51 | ! !PUBLIC MEMBER FUNCTIONS: |
---|
| 52 | ! |
---|
| 53 | public dstinit |
---|
| 54 | public dstrun |
---|
| 55 | public dstfin |
---|
| 56 | |
---|
| 57 | ! module variables |
---|
| 58 | character(len=*), parameter :: modelname='dstmodel.F90' |
---|
| 59 | integer :: rank, lcomm |
---|
| 60 | |
---|
| 61 | !EOP ------------------------------------------------------------------- |
---|
| 62 | |
---|
| 63 | contains |
---|
| 64 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 65 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 66 | !BOP ------------------------------------------------------------------- |
---|
| 67 | ! |
---|
| 68 | ! !IROUTINE: dstinit - Destination model initialization |
---|
| 69 | |
---|
| 70 | subroutine dstinit(GSMap,IMPORT,EXPORT,comm,compid) |
---|
| 71 | |
---|
| 72 | ! !INPUT PARAMETERS: |
---|
| 73 | type(GlobalSegMap),intent(inout) :: GSMap ! decomposition |
---|
| 74 | type(AttrVect),intent(inout) :: IMPORT,EXPORT ! state data |
---|
| 75 | integer,intent(in) :: comm ! MPI communicator |
---|
| 76 | integer,intent(in) :: compid ! component ID |
---|
| 77 | ! |
---|
| 78 | !EOP ___________________________________________________________________ |
---|
| 79 | |
---|
| 80 | ! local variables |
---|
| 81 | |
---|
| 82 | ! parameters for this model |
---|
| 83 | integer :: nxa ! number of points in x-direction |
---|
| 84 | integer :: nya ! number of points in y-direction |
---|
| 85 | |
---|
| 86 | integer :: i,j,k,idx |
---|
| 87 | |
---|
| 88 | integer :: nprocs, root, ier |
---|
| 89 | |
---|
| 90 | ! GlobalSegMap variables |
---|
| 91 | integer,dimension(:),pointer :: lindex |
---|
| 92 | |
---|
| 93 | ! AttrVect variables |
---|
| 94 | integer :: avsize |
---|
| 95 | |
---|
| 96 | character*2, ldecomp |
---|
| 97 | |
---|
| 98 | |
---|
| 99 | call MPI_COMM_RANK(comm,rank, ier) |
---|
| 100 | call MPI_COMM_SIZE(comm,nprocs,ier) |
---|
| 101 | |
---|
| 102 | ! save local communicator |
---|
| 103 | lcomm=comm |
---|
| 104 | |
---|
| 105 | if(rank==0) then |
---|
| 106 | write(6,*) modelname, ' init start' |
---|
| 107 | write(6,*) modelname,' MyID ', compid |
---|
| 108 | write(6,*) modelname,' Num procs ', nprocs |
---|
| 109 | endif |
---|
| 110 | |
---|
| 111 | ! Get configuration |
---|
| 112 | call i90_LoadF('dst.rc',ier) |
---|
| 113 | |
---|
| 114 | call i90_label('nx:',ier) |
---|
| 115 | nxa=i90_gint(ier) |
---|
| 116 | call i90_label('ny:',ier) |
---|
| 117 | nya=i90_gint(ier) |
---|
| 118 | if(rank==0) write(6,*) modelname, ' x,y ', nxa,nya |
---|
| 119 | |
---|
| 120 | call i90_label('decomp:',ier) |
---|
| 121 | call i90_Gtoken(ldecomp, ier) |
---|
| 122 | if(rank==0) write(6,*) modelname, ' decomp ', ldecomp |
---|
| 123 | |
---|
| 124 | |
---|
| 125 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 126 | ! Initialize a Global Segment Map |
---|
| 127 | |
---|
| 128 | |
---|
| 129 | call get_index(ldecomp,nprocs,rank,nxa,nya,lindex) |
---|
| 130 | |
---|
| 131 | call GlobalSegMap_init(GSMap,lindex,comm,compid,gsize=nxa*nya) |
---|
| 132 | |
---|
| 133 | |
---|
| 134 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 135 | |
---|
| 136 | |
---|
| 137 | if(rank==0) write(6,*) modelname, ' GSMap ',GSMap%ngseg,GSMap%gsize |
---|
| 138 | |
---|
| 139 | |
---|
| 140 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 141 | ! Initialize import and export Attribute vectors |
---|
| 142 | |
---|
| 143 | ! size is the number of grid points on this processor |
---|
| 144 | avsize = GlobalSegMap_lsize(GSMap,comm) |
---|
| 145 | if(rank==0) write(6,*) modelname, ' localsize ', avsize |
---|
| 146 | |
---|
| 147 | ! initialize Avs with two real attributes. |
---|
| 148 | call AttrVect_init(IMPORT,rList="field3:field4",lsize=avsize) |
---|
| 149 | call AttrVect_init(EXPORT,rList="field5:field6",lsize=avsize) |
---|
| 150 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 151 | |
---|
| 152 | if(rank==0) write(6,*) modelname, ' init done' |
---|
| 153 | end subroutine dstinit |
---|
| 154 | !!! END OF INIT !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 155 | |
---|
| 156 | |
---|
| 157 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 158 | ! RUN PHASE |
---|
| 159 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 160 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 161 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 162 | !BOP ------------------------------------------------------------------- |
---|
| 163 | ! |
---|
| 164 | ! !IROUTINE: dstrun - Destination model run method |
---|
| 165 | |
---|
| 166 | subroutine dstrun(IMPORT,EXPORT) |
---|
| 167 | |
---|
| 168 | ! !INPUT PARAMETERS: |
---|
| 169 | type(AttrVect),intent(inout) :: IMPORT,EXPORT ! Input and Output states |
---|
| 170 | |
---|
| 171 | !EOP ------------------------------------------------------------------- |
---|
| 172 | |
---|
| 173 | ! local variables |
---|
| 174 | integer :: avsize,ier,i,index |
---|
| 175 | |
---|
| 176 | if(rank==0) write(6,*) modelname, ' run start' |
---|
| 177 | |
---|
| 178 | ! Copy input data to output data using translation between different names |
---|
| 179 | call AttrVect_copy(IMPORT,EXPORT,rList="field3:field4", & |
---|
| 180 | TrList="field5:field6") |
---|
| 181 | |
---|
| 182 | if(rank==0) write(6,*) modelname, ' run done' |
---|
| 183 | |
---|
| 184 | end subroutine dstrun |
---|
| 185 | !!! END OF RUN !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 186 | |
---|
| 187 | |
---|
| 188 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 189 | ! FINALIZE PHASE |
---|
| 190 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 191 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 192 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
| 193 | !BOP ------------------------------------------------------------------- |
---|
| 194 | ! |
---|
| 195 | ! !IROUTINE: dstfin - Destination model finalize method |
---|
| 196 | |
---|
| 197 | subroutine dstfin(IMPORT,EXPORT,GSMap) |
---|
| 198 | |
---|
| 199 | ! !INPUT PARAMETERS: |
---|
| 200 | type(AttrVect),intent(inout) :: IMPORT,EXPORT ! MCT defined type |
---|
| 201 | type(GlobalSegMap),intent(inout) :: GSMap |
---|
| 202 | |
---|
| 203 | !EOP ------------------------------------------------------------------- |
---|
| 204 | type(AttrVect) :: GlobalD |
---|
| 205 | integer :: lsize,ier,mdev,i |
---|
| 206 | |
---|
| 207 | if(rank==0) write(6,*) modelname,' fin start' |
---|
| 208 | ! gather data to node 0 and write it out |
---|
| 209 | call AttrVect_gather(EXPORT,GlobalD,GSMap,0,lcomm,ier) |
---|
| 210 | |
---|
| 211 | ! write out gathered data |
---|
| 212 | if(rank==0) then |
---|
| 213 | mdev=luavail() |
---|
| 214 | lsize=AttrVect_lsize(GlobalD) |
---|
| 215 | open(mdev, file="TS1out.dat") |
---|
| 216 | do i=1,lsize |
---|
| 217 | write(mdev,*) GlobalD%rAttr(1,i) |
---|
| 218 | enddo |
---|
| 219 | close(mdev) |
---|
| 220 | endif |
---|
| 221 | |
---|
| 222 | ! clean up |
---|
| 223 | call AttrVect_clean(IMPORT) |
---|
| 224 | call AttrVect_clean(EXPORT) |
---|
| 225 | if(rank==0)call AttrVect_clean(GlobalD) |
---|
| 226 | call GlobalSegMap_clean(GSMap) |
---|
| 227 | if(rank==0) write(6,*) modelname,' fin done' |
---|
| 228 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 229 | endsubroutine dstfin |
---|
| 230 | |
---|
| 231 | end module dstmodel |
---|