source: CPL/oasis3-mct/branches/OASIS3-MCT_2.0_branch/lib/mct/examples/climate_sequen1/srcmodel.F90 @ 4775

Last change on this file since 4775 was 4775, checked in by aclsce, 5 years ago
  • Imported oasis3-mct from Cerfacs svn server (not suppotred anymore).

The version has been extracted from https://oasis3mct.cerfacs.fr/svn/branches/OASIS3-MCT_2.0_branch/oasis3-mct@1818

File size: 7.5 KB
Line 
1!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2!    Math and Computer Science Division, Argonne National Laboratory   !
3!-----------------------------------------------------------------------
4! CVS srcmodel.F90,v 1.8 2005-11-18 23:15:38 rloy Exp
5! CVS MCT_2_8_0
6!BOP -------------------------------------------------------------------
7!
8! !MODULE: srcmodel -- generic model for unit tester
9!
10! !DESCRIPTION:
11! init run and finalize methods for source model
12!
13module srcmodel
14
15!
16! !USES:
17!
18! Get the things needed from MCT by "Use,only" with renaming:
19!
20!---Domain Decomposition Descriptor DataType and associated methods
21use m_GlobalSegMap,only: GlobalSegMap
22use m_GlobalSegMap,only: GlobalSegMap_init => init
23use m_GlobalSegMap,only: GlobalSegMap_lsize => lsize
24use m_GlobalSegMap,only: GlobalSegMap_clean => clean
25!---Field Storage DataType and associated methods
26use m_AttrVect,only    : AttrVect
27use m_AttrVect,only    : AttrVect_init => init
28use m_AttrVect,only    : AttrVect_lsize => lsize
29use m_AttrVect,only    : AttrVect_clean => clean
30use m_AttrVect,only    : AttrVect_copy => copy
31use m_AttrVect,only    : AttrVect_zero => zero
32use m_AttrVect,only    : AttrVect_indxR => indexRA
33use m_AttrVect,only    : AttrVect_importRAttr => importRAttr
34use m_AttrVectComms,only : AttrVect_scatter => scatter
35
36! Get things from MPEU
37use m_inpak90   ! Resource files
38use m_stdio     ! I/O utils
39use m_ioutil
40
41! Get utilities for this program.
42use mutils
43
44implicit none
45
46private
47! except
48     
49! !PUBLIC MEMBER FUNCTIONS:
50
51public srcinit
52public srcrun
53public srcfin
54
55! private module variables
56character(len=*), parameter :: modelname='srcmodel.F90'
57integer  :: rank
58real, dimension(:), pointer :: avdata
59
60!EOP -------------------------------------------------------------------
61
62contains
63
64!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
65!    Math and Computer Science Division, Argonne National Laboratory   !
66!BOP -------------------------------------------------------------------
67!
68! !IROUTINE: srcinit  - Source model initialization
69
70subroutine srcinit(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,mdev,fx,fy
87  integer :: nprocs, root, ier,fileno
88
89! GlobalSegMap variables
90  integer,dimension(:),pointer :: lindex
91
92! AttrVect variables
93  integer :: avsize
94  type(AttrVect)   :: GlobalD  ! Av to hold global data
95
96  real,dimension(:),pointer :: rootdata
97
98  character*2 :: ldecomp
99
100
101  call MPI_COMM_RANK(comm,rank, ier)
102  call MPI_COMM_SIZE(comm,nprocs,ier)
103
104  if(rank==0) then
105    write(6,*) modelname, ' init start'
106    write(6,*) modelname,' MyID ', compid
107    write(6,*) modelname,' Num procs ', nprocs
108  endif
109
110!  Get configuration
111   call i90_LoadF('src.rc',ier)
112
113   call i90_label('nx:',ier)
114   nxa=i90_gint(ier)
115   call i90_label('ny:',ier)
116   nya=i90_gint(ier)
117  if(rank==0) write(6,*) modelname, ' x,y ', nxa,nya
118
119   call i90_label('decomp:',ier)
120   call i90_Gtoken(ldecomp, ier)
121  if(rank==0) write(6,*) modelname, ' decomp ', ldecomp
122
123!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
124! Initialize a Global Segment Map
125
126
127  call get_index(ldecomp,nprocs,rank,nxa,nya,lindex)
128
129  call GlobalSegMap_init(GSMap,lindex,comm,compid,gsize=nxa*nya)
130
131
132
133!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
134  if(rank==0) write(6,*) modelname, ' GSMap ',GSMap%ngseg,GSMap%gsize
135
136
137!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
138! Initialize import and export Attribute vectors
139
140! size is the number of grid points on this processor
141  avsize = GlobalSegMap_lsize(GSMap,comm)
142  if(rank==0) write(6,*) modelname, ' localsize ', avsize
143
144! Initialize the IMPORT  Av by scattering from a root Av
145! with real data.
146
147! Read in data from root and scatter to nodes
148  if(rank==0) then
149    call AttrVect_init(GlobalD,rList="field1:field2",lsize=nxa*nya)
150    mdev=luavail()
151    open(mdev, file="TS1.dat",status="old")
152    read(mdev,*) fx,fy
153    do i=1,nxa*nya
154      read(mdev,*) GlobalD%rAttr(1,i)
155    enddo
156    write(6,*) modelname,'Global init ',GlobalD%rAttr(1,1),GlobalD%rAttr(1,8000)
157  endif
158
159! this scatter will create IMPORT if it hasn't already been initialized
160  call AttrVect_scatter(GlobalD,IMPORT,GSMap,0,comm,ier)
161
162! initialize EXPORT Av with two real attributes.
163  call AttrVect_init(EXPORT,rList="field3:field4",lsize=avsize)
164 
165  call AttrVect_zero(EXPORT)
166
167  if(rank==0) then
168    write(6,*) modelname, rank,' IMPORT field1', IMPORT%rAttr(1,1)
169    write(6,*) modelname, rank,' IMPORt field2', IMPORT%rAttr(2,1)
170    write(6,*) modelname, rank,' EXPORT field3', EXPORT%rAttr(1,1)
171    write(6,*) modelname, rank,' EXPORT field4', EXPORT%rAttr(2,1)
172  endif
173
174! allocate buffer for use in run method
175  allocate(avdata(avsize),stat=ier)
176
177  if(rank==0) write(6,*) modelname, ' init done'
178end subroutine srcinit
179!!! END OF INIT !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
180
181
182!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
183!   RUN PHASE
184!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
185!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
186!    Math and Computer Science Division, Argonne National Laboratory   !
187!BOP -------------------------------------------------------------------
188!
189! !IROUTINE: srcrun  - Source model run method
190
191subroutine srcrun(IMPORT,EXPORT)
192
193! !INPUT PARAMETERS:
194  type(AttrVect),intent(inout) :: IMPORT,EXPORT   ! Input and Output states
195
196!EOP -------------------------------------------------------------------
197! local variables
198  integer :: avsize,ier,i
199     
200! Nothing to do with IMPORT
201
202
203! Fill EXPORT with data
204  if(rank==0) write(6,*) modelname, ' run start'
205
206! Use Av copy to copy input data from field1 in Imp to field3 in EXPORT
207  call AttrVect_copy(IMPORT,EXPORT,rList='field1',TrList='field3')
208
209! Use import to load data in second field
210  avdata=30.0
211  call AttrVect_importRAttr(EXPORT,"field4",avdata)
212
213  if(rank==0) write(6,*) modelname, ' In field1', IMPORT%rAttr(1,1)
214  if(rank==0) write(6,*) modelname, ' In field2', IMPORT%rAttr(2,1)
215  if(rank==0) write(6,*) modelname, ' Out field3', EXPORT%rAttr(1,1)
216  if(rank==0) write(6,*) modelname, ' Out field4', EXPORT%rAttr(2,1)
217
218  if(rank==0) write(6,*) modelname, ' run done'
219
220end subroutine srcrun
221!!! END OF RUN  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
222
223
224!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
225!   FINALIZE PHASE
226!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
227!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
228!    Math and Computer Science Division, Argonne National Laboratory   !
229!BOP -------------------------------------------------------------------
230!
231! !IROUTINE: srcfin  - Source model finalize method
232
233subroutine srcfin(IMPORT,EXPORT,GSMap)
234
235! !INPUT PARAMETERS:
236  type(AttrVect),intent(inout) :: IMPORT,EXPORT   ! imp,exp states
237  type(GlobalSegMap),intent(inout) :: GSMap
238!EOP -------------------------------------------------------------------
239 ! clean up
240  call AttrVect_clean(IMPORT)
241  call AttrVect_clean(EXPORT)
242  call GlobalSegMap_clean(GSMap)
243  deallocate(avdata)
244  if(rank==0) write(6,*) modelname,' fin done'
245!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
246endsubroutine srcfin
247
248end module srcmodel
Note: See TracBrowser for help on using the repository browser.