source: CPL/oasis3-mct/branches/OASIS3-MCT_5.0_branch/lib/psmile/src/mod_oasis_coupler.F90 @ 6331

Last change on this file since 6331 was 6331, checked in by aclsce, 17 months ago

Moved oasis-mct_5.0 in oasis3-mct/branches directory.

File size: 104.0 KB
Line 
1
2!> Initialize the OASIS coupler infrastructure
3
4MODULE mod_oasis_coupler
5!     - - - - - - - - - - - - - - - - - - - - - - - - - - -
6!
7  USE mod_oasis_kinds
8  USE mod_oasis_data
9  USE mod_oasis_parameters
10  USE mod_oasis_namcouple
11  USE mod_oasis_sys
12  USE mod_oasis_map
13  USE mod_oasis_part
14  USE mod_oasis_var
15  USE mod_oasis_mpi
16  USE mod_oasis_string
17  USE mod_oasis_io
18  USE mod_oasis_timer
19  USE mod_oasis_load_balancing
20  USE mct_mod
21  USE grids    ! scrip
22  USE netcdf
23
24  IMPLICIT NONE
25
26  private
27
28  public oasis_coupler_setup
29  public oasis_coupler_bldvarname
30  public oasis_coupler_unbldvarname
31
32! Type of data
33
34  public prism_router_type
35  public prism_coupler_type
36
37! COUPLING INFO
38
39  !> Router information for rearranging data on tasks
40  type prism_router_type
41     !--- fixed at initialization ---
42     type(mct_router)      :: router     !< router
43  end type prism_router_type
44
45  integer(kind=ip_i4_p),public,parameter :: prism_coupler_avsmax=5  !< maximum number of higher order terms in mapping
46
47  !> Coupler data for managing all aspects of coupling in OASIS
48  type prism_coupler_type
49     !--- fixed at initialization ---
50     type(mct_aVect)       :: aVect1   !< primary aVect
51     type(mct_aVect)       :: aVect1p  !< primary aVect, value at most recent coupling time (used in BLASNEW)
52     type(mct_aVect)       :: aVect1m  !< extra aVect needed for mapping
53     type(mct_aVect)       :: aVect2   !< higher order mapping data
54     type(mct_aVect)       :: aVect3   !< higher order mapping data
55     type(mct_aVect)       :: aVect4   !< higher order mapping data
56     type(mct_aVect)       :: aVect5   !< higher order mapping data
57     logical               :: aVon(prism_coupler_avsmax)  !< flags indicating whether aVects 2-5 are active
58     character(len=ic_xl)  :: rstfile  !< restart file
59     logical               :: writrest !< flag to write a restart file
60     character(len=ic_xl)  :: inpfile  !< input file if data is read
61     character(len=ic_xxl) :: fldlist  !< field list
62     integer(kind=ip_i4_p) :: nflds    !< number of fields
63     integer(kind=ip_i4_p),pointer :: varid(:)    !< varid for each field
64     logical               :: valid    !< is this coupler valid
65     integer(kind=ip_i4_p) :: namID    !< namcouple ID
66     integer(kind=ip_i4_p) :: partID   !< local variable partition ID
67     integer(kind=ip_i4_p) :: rpartID  !< router partition ID
68     integer(kind=ip_i4_p) :: routerID !< router ID
69     integer(kind=ip_i4_p) :: mapperID !< mapper ID
70     character(len=ic_med) :: maploc   !< map location setting, src or dst
71     integer(kind=ip_i4_p) :: ops      !< namcouple operation (ip_exported,...)
72     integer(kind=ip_i4_p) :: comp     !< other model compid to couple
73     integer(kind=ip_i4_p) :: tag      !< communcation tag
74     integer(kind=ip_i4_p) :: seq      !< sequence number
75     integer(kind=ip_i4_p) :: dt       !< coupling period (secs)
76     integer(kind=ip_i4_p) :: lag      !< put lag positive is put sooner (secs)
77     integer(kind=ip_i4_p) :: maxtime  !< max time for the coupler
78     integer(kind=ip_i4_p) :: trans    !< transformation (ip_average,...)
79     integer(kind=ip_i4_p) :: conserv  !< conserve operation (ip_cnone,ip_cglobal,...)
80     character(len=ic_med) :: consopt  !< conserve option (bfb, opt)
81     integer(kind=ip_i4_p) :: getput   !< get/put flag
82     logical               :: sndrcv   !< send recv flag
83     logical               :: output   !< output flag
84     logical               :: input    !< input flag
85     logical               :: snddiag  !< diagnose src fields as part of coupling
86     logical               :: rcvdiag  !< diagnose rcv fields as part of coupling
87     real(kind=ip_double_p):: sndmult  !< send field multiplier term
88     real(kind=ip_double_p):: sndadd   !< send field addition term
89     real(kind=ip_double_p):: rcvmult  !< receive field multiplier term
90     real(kind=ip_double_p):: rcvadd   !< receive field addition term
91     logical               :: blasfld  !< flag indicates field used in other BLASNEW computations
92     integer(kind=ip_i4_p) :: rfno     !< BLASNEW extra field number for local variable
93     character(len=ic_med) ,pointer :: rfna(:)  !< blasnew extra field name
94     real(kind=ip_double_p),pointer :: rfmult(:)!< blasnew extra field multiplier term
95     real(kind=ip_double_p),pointer :: rfadd(:) !< blasnew extra field addition term
96     integer(kind=ip_i4_p) ,pointer :: rfcpl(:) !< blasnew extra field coupler id
97     integer(kind=ip_i4_p) ,pointer :: rffnum(:)!< blasnew extra field field number
98     !--- time varying info ---
99     integer(kind=ip_i4_p) :: ltime    !< time at last coupling
100     integer(kind=ip_i4_p) :: ctime    !< time at last call
101     integer(kind=ip_i4_p),pointer :: avcnt(:)  !< counter for averaging
102     integer(kind=ip_i4_p),pointer :: status(:) !< status of variables in coupler
103  end type prism_coupler_type
104
105  integer(kind=ip_i4_p)           :: prism_mrouter   !< max routers
106  integer(kind=ip_i4_p)           :: prism_nrouter = 0  !< router counter
107  type(prism_router_type) ,public, pointer:: prism_router(:)  !< prism_router array
108
109  integer(kind=ip_i4_p)   ,public :: prism_mcoupler  !< max couplers
110  type(prism_coupler_type),public, pointer :: prism_coupler_put(:)  !< prism_coupler put array
111  type(prism_coupler_type),public, pointer :: prism_coupler_get(:)  !< prism_coupler get array
112
113  integer(kind=ip_i4_p)   ,public :: lastseq       !< last coupler sequence
114  integer(kind=ip_i4_p)   ,public :: lastseqtime   !< last coupler sequence time
115  logical                 ,public :: allow_no_restart  !< flag to allow no restart files at startup
116
117!#include <netcdf.inc>
118
119!------------------------------------------------------------
120CONTAINS
121!------------------------------------------------------------
122
123!> Main routine to setup couplers
124
125!> This routine initializes all the coupler data based on the namcouple
126!> inputs and the calls into the OASIS initialization interfaces from models.
127!> It reconciles everything.  This is called from oasis_enddef.
128
129  SUBROUTINE oasis_coupler_setup()
130
131  !----------------------------------------------------------
132  ! This routine reconciles the coupling stuff
133  !----------------------------------------------------------
134
135  IMPLICIT none
136
137  integer(kind=ip_i4_p) :: l,n,n1,n2,nn,nv,nm,nv1,nv1a,nns,lnn,nc,nf,nvf,npc,r1,ierr
138  integer(kind=ip_i4_p) :: pe,nflds1,nflds2,ncnt,nc2,nf2
139  integer(kind=ip_i4_p) :: part1, part2
140  integer(kind=ip_i4_p) :: spart,dpart ! src, dst partitions for mapping
141        ! part1 = my local part, partID
142        ! part2 = other partition for mapping
143        ! spart = src part for mapping; put=part1, get=part2
144        ! dpart = dst part for mapping; put=part2, get=part1
145  integer(kind=ip_i4_p) :: mapID,namID
146  type(mct_sMat),pointer :: sMati(:)
147  integer(kind=ip_i4_p) :: ncid,dimid,status
148  integer(kind=ip_i4_p) :: lsize,gsize
149  integer(kind=ip_i4_p) :: svarid
150  integer(kind=ip_i4_p),allocatable :: varidtmp(:)
151  integer(kind=ip_i4_p) :: part
152  integer(kind=ip_i4_p),pointer :: varid1(:)
153  character(len=ic_med) :: cstring,delim,vname,mapopt,mapopt1
154  character(len=ic_lvar):: myfld
155  integer(kind=ip_i4_p) :: myfldi
156  character(len=ic_xxl) :: myfldlist  ! field list for my model
157  character(len=ic_lvar):: otfld
158  character(len=ic_xxl) :: otfldlist  ! field list for other model
159  integer(kind=ip_i4_p) :: nx,ny
160  character(len=ic_lvar):: gridname
161  character(len=ic_long):: tmp_mapfile, tmp_mapfile2
162  integer(kind=ip_i4_p) :: flag
163  integer(kind=ip_i4_p) :: kmask,kfrac,karea
164  logical               :: found, exists, found2, rmask, rfrac
165  integer(kind=ip_i4_p) :: mynvar
166  integer(kind=ip_i4_p) :: nwgts, arrlen
167  character(len=ic_lvar):: tmpfld
168  character(len=ic_lvar):: tmpstr
169  type(prism_coupler_type),pointer :: pcpointer
170  type(prism_coupler_type),pointer :: pcpntpair
171  integer(kind=ip_i4_p) :: ifind,nfind
172  integer(kind=ip_i4_p) ,pointer :: gridID(:)
173  character(len=ic_lvar),pointer :: myvar(:)
174  integer(kind=ip_i4_p) ,pointer :: myops(:)
175  integer(kind=ip_i4_p) ,pointer :: mynum(:)
176  integer(kind=ip_i4_p) ,pointer :: nallvar(:)
177  character(len=ic_lvar),pointer :: allvar(:,:)
178  integer(kind=ip_i4_p) ,pointer :: allops(:,:)
179  integer(kind=ip_i4_p) ,pointer :: allnum(:,:)
180  integer(kind=ip_i4_p) ,pointer :: namsrc_checkused(:) ! 0 = not used
181  integer(kind=ip_i4_p) ,pointer :: namsrc_checkused_g(:)  ! 0 = not used
182  type sortnamfld_type
183     integer(kind=ip_i4_p) :: num                 ! total number of namcouple fields
184     integer(kind=ip_i4_p) ,pointer :: namnum(:)  ! namcouple number
185     integer(kind=ip_i4_p) ,pointer :: fldnum(:)  ! namcouple field number in namcouple
186     character(len=ic_lvar),pointer :: fld(:)     ! namcouple field name
187  end type sortnamfld_type
188  type(sortnamfld_type) :: sortnsrc
189  type(sortnamfld_type) :: sortndst
190  type sortvarfld_type
191     integer(kind=ip_i4_p) :: num                 ! total number of var fields
192     integer(kind=ip_i4_p) ,pointer :: modnum(:)  ! model number
193     integer(kind=ip_i4_p) ,pointer :: varnum(:)  ! var field number in model
194     character(len=ic_lvar),pointer :: fld(:)     ! variable field name
195  end type sortvarfld_type
196  type(sortvarfld_type) :: sortvars
197  type(sortvarfld_type) :: sorttest
198  integer(kind=ip_i4_p) ,pointer :: sortkey(:)
199  character(len=ic_med) :: part2decomp   ! decomp_1d or decomp_wghtfile
200  character(len=ic_med) :: smatread_method ! orig or ceg
201  integer, parameter :: local_timers_on = 0   ! 0=min, 1=few, 2=med, 3=max
202  integer(kind=ip_i4_p) :: nb_field, nb_cpl_ts, imain_kind_lb
203  logical :: lmap,  lout, lrst, ltrn
204
205  character(len=*),parameter :: subname = '(oasis_coupler_setup)'
206
207  !----------------------------------------------------------
208
209  call oasis_debug_enter(subname)
210!  call oasis_mpi_barrier(mpi_comm_global)
211  IF (local_timers_on >= 1) then
212     call oasis_timer_start('cpl_setup_barrier')
213     call oasis_mpi_barrier(mpi_comm_global, 'cpl_setup')
214     call oasis_timer_stop('cpl_setup_barrier')
215     call oasis_timer_start('cpl_setup')
216  endif
217
218  if (local_timers_on >= 2) call oasis_timer_start('cpl_setup_n1')
219
220  part2decomp = nammapdec
221  smatread_method = nammatxrd
222  if (OASIS_debug >= 2) then
223     write(nulprt,*) subname,' part2decomp = ',trim(part2decomp)
224     write(nulprt,*) subname,' smatread_method = ',trim(smatread_method)
225  endif
226
227  !-----------------------------------------
228  !> * Allocate and zero prism_router, prism_mapper, prism_coupler based on nnamcpl
229  ! there cannot be more than that needed
230  !-----------------------------------------
231
232  call oasis_debug_note(subname//' set defaults for datatypes')
233
234  prism_mrouter = nnamcpl * 2   ! multiply by 2 for coupling to self
235  allocate(prism_router(prism_mrouter))
236  prism_nrouter = 0
237
238  prism_mmapper = nnamcpl
239  allocate(prism_mapper(prism_mmapper))
240  prism_nmapper = 0
241  prism_mapper(:)%nwgts = 0
242  prism_mapper(:)%file  = ""
243  prism_mapper(:)%file2 = ""
244  prism_mapper(:)%loc   = ""
245  prism_mapper(:)%opt   = ""
246  prism_mapper(:)%optval= ""
247  prism_mapper(:)%init  = .false.
248  prism_mapper(:)%spart = ispval
249  prism_mapper(:)%dpart = ispval
250
251  prism_mcoupler = nnamcpl
252  allocate(prism_coupler_put(prism_mcoupler))
253  allocate(prism_coupler_get(prism_mcoupler))
254
255  do nc = 1,prism_mcoupler
256  do npc = 1,2
257     if (npc == 1) then
258        pcpointer => prism_coupler_put(nc)
259        pcpntpair => prism_coupler_get(nc)
260     endif
261     if (npc == 2) then
262        pcpointer => prism_coupler_get(nc)
263        pcpntpair => prism_coupler_put(nc)
264     endif
265     pcpointer%rstfile = ""
266     pcpointer%writrest= .false.
267     pcpointer%inpfile = ""
268     pcpointer%fldlist = ""
269     pcpointer%nflds   = 0
270     pcpointer%namID   = 0
271     pcpointer%valid   = .false.
272!tcx is this alloc pcpointer or prism_coupler_*
273     allocate(pcpointer%varid(1))
274     pcpointer%varid(:) = ispval
275     pcpointer%aVon(:) = .false.
276     pcpointer%ops     = ispval
277     pcpointer%comp    = ispval
278     pcpointer%routerID  = ispval
279     pcpointer%mapperID  = ispval
280     pcpointer%maploc  = ""
281     pcpointer%tag     = ispval
282     pcpointer%dt      = ispval
283     pcpointer%lag     = 0
284     pcpointer%maxtime = 0
285     pcpointer%getput  = ispval
286     pcpointer%sndrcv  = .false.
287     pcpointer%output  = .false.
288     pcpointer%input   = .false.
289     pcpointer%trans   = ip_instant
290     pcpointer%conserv = ip_cnone
291     pcpointer%ltime   = ispval
292     pcpointer%ctime   = ispval
293     pcpointer%snddiag = .false.
294     pcpointer%rcvdiag = .false.
295     pcpointer%sndmult = 1.0_ip_double_p
296     pcpointer%sndadd  = 0.0_ip_double_p
297     pcpointer%rcvmult = 1.0_ip_double_p
298     pcpointer%rcvadd  = 0.0_ip_double_p
299     pcpointer%blasfld = .false.
300     pcpointer%rfno    = 0
301  enddo  ! npc
302  enddo  ! nc
303
304  lastseq      = ispval
305  lastseqtime  = ispval
306
307  !----------------------------------------------------------
308  !> * Generate model variable lists across all models based on def_var calls.
309  !> These will be reconciled with the namcouple input.  These are sorted
310  !> to improve search performance later.
311  !----------------------------------------------------------
312
313  call oasis_debug_note(subname//' share var info between models')
314
315  allocate(allvar(maxvar,prism_amodels))
316  allocate(nallvar(prism_amodels))
317  allocate(allops(maxvar,prism_amodels))
318  allocate(allnum(maxvar,prism_amodels))
319  allocate(myvar(maxvar))
320  allocate(myops(maxvar))
321  allocate(mynum(maxvar))
322
323  allvar = " "
324  nallvar = 0
325  allops = -1
326  allnum = -1
327  if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n1_bcast')
328  do n = 1,prism_amodels
329     if (n == compid) then
330        myvar = " "
331        myops = 0
332        mynum = 0
333        mynvar = prism_nvar
334        do n1 = 1, prism_nvar
335           myvar(n1) = trim(prism_var(n1)%name)
336           myops(n1) = prism_var(n1)%ops
337           mynum(n1) = prism_var(n1)%num
338           ! check that each var name is unique for a given model
339           do n2 = 1,n1-1
340              if (myvar(n1) == myvar(n2)) then
341                 WRITE(nulprt,*) subname,estr,'variable name defined more than once by def_var = ',trim(myvar(n1))
342                 call oasis_abort(file=__FILE__,line=__LINE__)
343              endif
344           enddo
345        enddo
346     endif
347     if (OASIS_debug >= 5) then
348        write(nulprt,*) subname,' BCAST from ',n,mpi_root_global(n)
349        call oasis_flush(nulprt)
350     endif
351     call oasis_mpi_bcast(mynvar,mpi_comm_global,'mynvar',mpi_root_global(n))
352     if (OASIS_debug >= 5) then
353        write(nulprt,*) subname,' bcast mynvar ',mynvar
354        call oasis_flush(nulprt)
355     endif
356     nallvar(n) = mynvar
357     call oasis_mpi_bcast(myvar,mpi_comm_global,'myvar',mpi_root_global(n))
358     if (OASIS_debug >= 5) then
359        write(nulprt,*) subname,' bcast myvar ',trim(myvar(1))
360        call oasis_flush(nulprt)
361     endif
362     allvar(:,n) = myvar(:)
363     call oasis_mpi_bcast(myops,mpi_comm_global,'myops',mpi_root_global(n))
364     if (OASIS_debug >= 5) then
365        write(nulprt,*) subname,' bcast myops ',myops(1)
366        call oasis_flush(nulprt)
367     endif
368     allops(:,n) = myops(:)
369     call oasis_mpi_bcast(mynum,mpi_comm_global,'mynum',mpi_root_global(n))
370     if (OASIS_debug >= 5) then
371        write(nulprt,*) subname,' bcast mynum ',mynum(1)
372        call oasis_flush(nulprt)
373     endif
374     allnum(:,n) = mynum(:)
375  enddo
376  if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n1_bcast')
377
378  deallocate(myvar,myops)
379
380  if (OASIS_debug >= 2) then
381     write(nulprt,*) subname,' model variable info:'
382     do nm = 1,prism_amodels
383        write(nulprt,'(8x,a,2i6)') ' model,nvars = ',nm,nallvar(nm)
384        do nv = 1,nallvar(nm)
385           cstring = 'unknown'
386           if (allops(nv,nm) == OASIS_Out) cstring = 'prism_out'
387           if (allops(nv,nm) == OASIS_In)  cstring = 'prism_in'
388           write(nulprt,'(16x,a,2i6,2x,a,2i6,2x,a)') ' model,idx,var,num,ops = ',nm,nv,&
389                                                      trim(allvar(nv,nm)),allnum(nv,nm),allops(nv,nm),trim(cstring)
390        enddo
391     enddo
392     write(nulprt,*) ' '
393     call oasis_flush(nulprt)
394  endif
395
396  ! generate sortvars sorted list
397
398  n1 = 0
399  do n = 1,prism_amodels
400     n1 = n1 + nallvar(n)
401  enddo
402  allocate(sortvars%fld(n1))
403  allocate(sortvars%modnum(n1))
404  allocate(sortvars%varnum(n1))
405  allocate(sortkey(n1))
406  sortvars%num = n1
407
408  n1 = 0
409  do n = 1,prism_amodels
410  do n2 = 1,nallvar(n)
411     n1 = n1 + 1
412     sortkey(n1) = n1
413     sortvars%fld(n1) = allvar(n2,n)
414     sortvars%modnum(n1) = n
415     sortvars%varnum(n1) = n2
416  enddo
417  enddo
418
419  call oasis_sys_sortC(sortvars%num, sortvars%fld, sortkey)
420  call oasis_sys_sortIkey(sortvars%num, sortvars%modnum, sortkey)
421  call oasis_sys_sortIkey(sortvars%num, sortvars%varnum, sortkey)
422
423  if (OASIS_debug >= 15) then
424     write(nulprt,*) subname//' Sorted array : sortvars'
425     do n1 = 1,sortvars%num
426        write(nulprt,*) subname,'sort sortvars',n1,sortkey(n1),sortvars%modnum(n1),sortvars%varnum(n1),trim(sortvars%fld(n1))
427     enddo
428  endif
429
430  deallocate(sortkey)
431
432  !----------------------------------------------------------
433  !> * Setup couplers based on namcouple and model variable info.
434  ! These must be paired up consistently, create couplers in
435  ! sorted order (nns)
436  ! nn = namcpl counter, sorted
437  ! nm = model counter, compid is my nm
438  ! nv = variable counter
439  ! nv1 = my variable counter
440  !----------------------------------------------------------
441
442  if (local_timers_on >= 2) call oasis_timer_stop ('cpl_setup_n1')
443
444  !--------------------------------
445  !> * Preprocess namcouple strings and sort for faster searches
446  !--------------------------------
447
448  ! count namcouple field names
449
450  if (local_timers_on >= 2) call oasis_timer_start('cpl_setup_n2')
451  n1 = 0
452  n2 = 0
453  do nn = 1,nnamcpl
454     n1 = n1 + oasis_string_listGetNum(namsrcfld(nn))
455     n2 = n2 + oasis_string_listGetNum(namdstfld(nn))
456     if (n1 /= n2) then
457        WRITE(nulprt,*) subname,estr,'number of fields in namcouple inconsistent ',nn,n1,n2
458        WRITE(nulprt,*) subname,estr,'namcouple src fields = ',trim(namsrcfld(nn))
459        WRITE(nulprt,*) subname,estr,'namcouple dst fields = ',trim(namdstfld(nn))
460        call oasis_abort(file=__FILE__,line=__LINE__)
461     endif
462  enddo
463
464  ! allocate space
465  ! note: n2==n1
466
467  sortnsrc%num = n1
468  allocate(sortnsrc%fld(n1))
469  allocate(sortnsrc%namnum(n1))
470  allocate(sortnsrc%fldnum(n1))
471  sortndst%num = n2
472  allocate(sortndst%fld(n2))
473  allocate(sortndst%namnum(n2))
474  allocate(sortndst%fldnum(n2))
475
476  ! this will check that all namcouple vars are used in application
477  allocate(namsrc_checkused(sortnsrc%num))
478  namsrc_checkused = 0
479
480  ! fill and sort sortnsrc
481
482  allocate(sortkey(sortnsrc%num))
483  n1 = 0
484  do nn = 1,nnamcpl
485  do n2 = 1,oasis_string_listGetNum(namsrcfld(nn))
486     n1 = n1 + 1
487     sortkey(n1) = n1
488     sortnsrc%namnum(n1) = nn
489     sortnsrc%fldnum(n1) = n2
490     call oasis_string_listGetName(namsrcfld(nn),n2,sortnsrc%fld(n1))
491  enddo
492  enddo
493
494  call oasis_sys_sortC(sortnsrc%num, sortnsrc%fld, sortkey)
495  call oasis_sys_sortIkey(sortnsrc%num, sortnsrc%namnum, sortkey)
496  call oasis_sys_sortIkey(sortnsrc%num, sortnsrc%fldnum, sortkey)
497
498  if (OASIS_debug >= 15) then
499     write(nulprt,*) subname//' Sorted array : sortnsrc'
500     do n1 = 1,sortnsrc%num
501        write(nulprt,*) subname,'sort sortnsrc',n1,sortkey(n1), &
502           sortnsrc%namnum(n1),sortnsrc%fldnum(n1),trim(sortnsrc%fld(n1))
503     enddo
504  endif
505  deallocate(sortkey)
506
507  ! fill and sort sortndst
508
509  allocate(sortkey(sortndst%num))
510  n1 = 0
511  do nn = 1,nnamcpl
512  do n2 = 1,oasis_string_listGetNum(namdstfld(nn))
513     n1 = n1 + 1
514     sortkey(n1) = n1
515     sortndst%namnum(n1) = nn
516     sortndst%fldnum(n1) = n2
517     call oasis_string_listGetName(namdstfld(nn),n2,sortndst%fld(n1))
518  enddo
519  enddo
520
521  call oasis_sys_sortC(sortndst%num, sortndst%fld, sortkey)
522  call oasis_sys_sortIkey(sortndst%num, sortndst%namnum, sortkey)
523  call oasis_sys_sortIkey(sortndst%num, sortndst%fldnum, sortkey)
524
525  if (OASIS_debug >= 15) then
526     write(nulprt,*) subname//' Sorted array : sortndst'
527     do n1 = 1,sortndst%num
528        write(nulprt,*) subname,'sort sortndst',n1,sortkey(n1), &
529           sortndst%namnum(n1),sortndst%fldnum(n1),trim(sortndst%fld(n1))
530     enddo
531  endif
532  deallocate(sortkey)
533
534  !==========================================================
535  ! Test Sort Code
536  !==========================================================
537  if (OASIS_debug >= 1500) then
538
539     write(nulprt,*) subname,' Test sort code: '
540
541     n1 = 10
542     allocate(sorttest%fld(n1))
543     allocate(sorttest%modnum(n1))
544     allocate(sorttest%varnum(n1))
545     allocate(sortkey(n1))
546     sorttest%num = n1
547
548     sorttest%fld(:) = 'A'
549     do n1 = 1,sorttest%num
550        sortkey(n1) = n1
551        if (n1 ==  1) sorttest%fld(n1) = 'D'
552        if (n1 ==  2) sorttest%fld(n1) = 'C'
553        if (n1 ==  4) sorttest%fld(n1) = 'C'
554        if (n1 ==  5) sorttest%fld(n1) = 'D'
555        if (n1 ==  8) sorttest%fld(n1) = 'C'
556        if (n1 ==  9) sorttest%fld(n1) = 'B'
557        if (n1 == 10) sorttest%fld(n1) = 'C'
558        sorttest%modnum(n1) = n1+100
559        sorttest%varnum(n1) = n1
560     enddo
561
562     call oasis_sys_sortC(sorttest%num, sorttest%fld, sortkey)
563     call oasis_sys_sortIkey(sorttest%num, sorttest%modnum, sortkey)
564     call oasis_sys_sortIkey(sorttest%num, sorttest%varnum, sortkey)
565
566     write(nulprt,*) subname//' Sorted array : sorttest'
567     do n1 = 1,sorttest%num
568        write(nulprt,*) subname,'sort sorttest',n1,sortkey(n1), &
569           sorttest%modnum(n1),sorttest%varnum(n1),trim(sorttest%fld(n1))
570     enddo
571
572     tmpfld = 'A'
573     call cplfind(sorttest%num, sorttest%fld, tmpfld, ifind, nfind)
574     write(nulprt,*) subname,' cpl find1 ',trim(tmpfld),ifind,nfind
575     do n1 = ifind,ifind+nfind-1
576        write(nulprt,*) subname,' cpl find2 ',n1,trim(sorttest%fld(n1))
577     enddo
578
579     tmpfld = 'B'
580     call cplfind(sorttest%num, sorttest%fld, tmpfld, ifind, nfind)
581     write(nulprt,*) subname,' cpl find1 ',trim(tmpfld),ifind,nfind
582     do n1 = ifind,ifind+nfind-1
583        write(nulprt,*) subname,' cpl find2 ',n1,trim(sorttest%fld(n1))
584     enddo
585
586     tmpfld = 'C'
587     call cplfind(sorttest%num, sorttest%fld, tmpfld, ifind, nfind)
588     write(nulprt,*) subname,' cpl find1 ',trim(tmpfld),ifind,nfind
589     do n1 = ifind,ifind+nfind-1
590        write(nulprt,*) subname,' cpl find2 ',n1,trim(sorttest%fld(n1))
591     enddo
592
593     tmpfld = 'D'
594     call cplfind(sorttest%num, sorttest%fld, tmpfld, ifind, nfind)
595     write(nulprt,*) subname,' cpl find1 ',trim(tmpfld),ifind,nfind
596     do n1 = ifind,ifind+nfind-1
597        write(nulprt,*) subname,' cpl find2 ',n1,trim(sorttest%fld(n1))
598     enddo
599
600     tmpfld = 'E'
601     call cplfind(sorttest%num, sorttest%fld, tmpfld, ifind, nfind)
602     write(nulprt,*) subname,' cpl find1 ',trim(tmpfld),ifind,nfind
603     do n1 = ifind,ifind+nfind-1
604        write(nulprt,*) subname,' cpl find2 ',n1,trim(sorttest%fld(n1))
605     enddo
606
607     deallocate(sortkey)
608     deallocate(sorttest%fld)
609     deallocate(sorttest%modnum)
610     deallocate(sorttest%varnum)
611
612     write(nulprt,*) subname,' Test cplfind: '
613     n1 = max(min(sortndst%num,sortndst%num/3),1)
614     tmpfld = sortndst%fld(n1)
615     call cplfind(sortndst%num, sortndst%fld, tmpfld, ifind, nfind)
616     write(nulprt,*) subname,' cpl find1 ',trim(tmpfld),ifind,nfind
617     do n1 = ifind,ifind+nfind-1
618        write(nulprt,*) subname,' cpl find2 ',n1,trim(sortndst%fld(n1))
619     enddo
620
621     n1 = max(min(sortndst%num,1),1)
622     tmpfld = sortndst%fld(n1)
623     call cplfind(sortndst%num, sortndst%fld, tmpfld, ifind, nfind)
624     write(nulprt,*) subname,' cpl find1 ',trim(tmpfld),ifind,nfind
625     do n1 = ifind,ifind+nfind-1
626        write(nulprt,*) subname,' cpl find2 ',n1,trim(sortndst%fld(n1))
627     enddo
628
629     n1 = max(min(sortndst%num,2),1)
630     tmpfld = sortndst%fld(n1)
631     call cplfind(sortndst%num, sortndst%fld, tmpfld, ifind, nfind)
632     write(nulprt,*) subname,' cpl find1 ',trim(tmpfld),ifind,nfind
633     do n1 = ifind,ifind+nfind-1
634        write(nulprt,*) subname,' cpl find2 ',n1,trim(sortndst%fld(n1))
635     enddo
636
637     n1 = max(min(sortndst%num,sortndst%num-1),1)
638     tmpfld = sortndst%fld(n1)
639     call cplfind(sortndst%num, sortndst%fld, tmpfld, ifind, nfind)
640     write(nulprt,*) subname,' cpl find1 ',trim(tmpfld),ifind,nfind
641     do n1 = ifind,ifind+nfind-1
642        write(nulprt,*) subname,' cpl find2 ',n1,trim(sortndst%fld(n1))
643     enddo
644
645     n1 = max(min(sortndst%num,sortndst%num),1)
646     tmpfld = sortndst%fld(n1)
647     call cplfind(sortndst%num, sortndst%fld, tmpfld, ifind, nfind)
648     write(nulprt,*) subname,' cpl find1 ',trim(tmpfld),ifind,nfind
649     do n1 = ifind,ifind+nfind-1
650        write(nulprt,*) subname,' cpl find2 ',n1,trim(sortndst%fld(n1))
651     enddo
652
653     CALL oasis_flush(nulprt)
654  endif
655  !==========================================================
656  ! END Test Sort Code
657  !==========================================================
658
659  if (local_timers_on >= 2) call oasis_timer_stop ('cpl_setup_n2')
660
661  call oasis_debug_note(subname//' compare vars and namcouple')
662  call oasis_debug_note(subname//' setup couplers')
663
664  if (local_timers_on >= 2) call oasis_timer_start('cpl_setup_n3')
665
666  !--------------------------------
667  !> * Loop over all my model variables
668  !--------------------------------
669
670  do nv1 = 1,prism_nvar
671
672     !--------------------------------
673     !>   * Get parition and field information
674     !--------------------------------
675
676     part1  = prism_var(nv1)%part
677     myfld  = prism_var(nv1)%name
678
679     IF (OASIS_debug >= 20) THEN
680        WRITE(nulprt,*) ' '
681        WRITE(nulprt,*) subname,' get part and fld ',nv1,part1,trim(myfld)
682        CALL oasis_flush(nulprt)
683     ENDIF
684
685     !--------------------------------
686     !>   * Check if variable is In or Out and then find namcouple matches
687     !--------------------------------
688
689     if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n3a')
690     if (prism_var(nv1)%ops == OASIS_Out) then
691        call cplfind(sortnsrc%num, sortnsrc%fld, myfld, ifind, nfind)
692     elseif (prism_var(nv1)%ops == OASIS_In) then
693        call cplfind(sortndst%num, sortndst%fld, myfld, ifind, nfind)
694     endif
695     if (local_timers_on >= 3) call oasis_timer_stop ('cpl_setup_n3a')
696
697     !--------------------------------
698     !>   * Loop over the namcouple matches
699     !--------------------------------
700     do nf = ifind,ifind+nfind-1
701        if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n3b')
702
703        flag = OASIS_NotDef
704
705        if (prism_var(nv1)%ops == OASIS_Out) then
706           nn = sortnsrc%namnum(nf)
707           myfldi = sortnsrc%fldnum(nf)
708           myfldlist = namsrcfld(nn)
709           otfldlist = namdstfld(nn)
710           flag = OASIS_Out
711        elseif (prism_var(nv1)%ops == OASIS_In) then
712           nn = sortndst%namnum(nf)
713           myfldi = sortndst%fldnum(nf)
714           myfldlist = namdstfld(nn)
715           otfldlist = namsrcfld(nn)
716           flag = OASIS_In
717        endif
718
719        nns = namnn2sort(nn)
720
721        IF (OASIS_debug >= 20) THEN
722            WRITE(nulprt,*) subname,' found fld1 ',trim(myfld),nv1,nf
723            WRITE(nulprt,*) subname,' found fld2 ',trim(myfld),nns,nn,myfldi,flag
724            CALL oasis_flush(nulprt)
725        ENDIF
726
727        if (local_timers_on >= 3) call oasis_timer_stop ('cpl_setup_n3b')
728
729        !--------------------------------
730        ! my variable is in this namcouple input
731        !--------------------------------
732
733        if (flag /= OASIS_NotDef) then
734
735           if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n3c')
736
737           !--------------------------------
738           !>     * Migrate namcouple info into part
739           !--------------------------------
740
741           IF (OASIS_debug >= 20) THEN
742               WRITE(nulprt,*) subname,' migrate namcouple to part '
743               CALL oasis_flush(nulprt)
744           ENDIF
745
746           if (flag == OASIS_In) then
747              if (prism_part(part1)%nx < 1) then
748                 prism_part(part1)%nx = namdst_nx(nn)
749                 prism_part(part1)%ny = namdst_ny(nn)
750                 prism_part(part1)%gridname = trim(namdstgrd(nn))
751              endif
752           endif
753           if (flag == OASIS_Out) then
754              if (prism_part(part1)%nx < 1) then
755                 prism_part(part1)%nx = namsrc_nx(nn)
756                 prism_part(part1)%ny = namsrc_ny(nn)
757                 prism_part(part1)%gridname = trim(namsrcgrd(nn))
758              endif
759           endif
760
761           IF (OASIS_debug >= 20) THEN
762               WRITE(nulprt,*) subname,' Field : ',trim(prism_var(nn)%name)
763               WRITE(nulprt,*) subname,' Grid dst : ',trim(namdstgrd(nn))
764               WRITE(nulprt,*) subname,' Grid src : ',trim(namsrcgrd(nn))
765!               WRITE(nulprt,*) subname,' prism_part : ',prism_part(part1)%gridname
766               CALL oasis_flush(nulprt)
767           ENDIF
768
769           !--------------------------------
770           !>     * Make sure it's either an In or Out, sanity check
771           !--------------------------------
772
773           if (flag /= OASIS_In .and. flag /= OASIS_Out) then
774              write(nulprt,*) subname,estr,'var must be either OASIS_In or OASIS_Out for var = ',trim(myfld)
775              call oasis_abort(file=__FILE__,line=__LINE__)
776           endif
777
778           if (OASIS_debug >= 20) then
779              write(nulprt,'(1x,2a,4i6,2a)') subname,' ca: myfld',nn,compid,&
780                                             nv1,myfldi,' ',trim(myfld)
781              call oasis_flush(nulprt)
782           endif
783
784           !--------------------------------
785           !>     * Determine matching field name from namcouple
786           !--------------------------------
787
788           if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n3c1')
789           otfld = 'NOmatchNOyesNOyesNO'
790           call oasis_string_listGetName(otfldlist,myfldi,otfld)
791           if (local_timers_on >= 3) call oasis_timer_stop ('cpl_setup_n3c1')
792
793           IF (OASIS_debug >= 20) THEN
794              WRITE(nulprt,*) subname,' otfld ',trim(otfld)
795              CALL oasis_flush(nulprt)
796           ENDIF
797
798           !--------------------------------
799           !>     * Search for list of models with other variable
800           !--------------------------------
801
802           if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n3c2')
803           call cplfind(sortvars%num, sortvars%fld, otfld, ifind, nfind)
804           if (local_timers_on >= 3) call oasis_timer_stop ('cpl_setup_n3c2')
805           if (local_timers_on >= 3) call oasis_timer_stop ('cpl_setup_n3c')
806
807           !--------------------------------
808           !>     * Loop over those other matching variable names
809           !--------------------------------
810           found = .false.
811           do nvf = ifind, ifind+nfind-1
812
813              ! check used appropriate array value, we are using "src" side sorted list
814              ! if output, just set the nf value
815              ! if input, search for an nn and myfldi match in the list
816             
817              if (prism_var(nv1)%ops == OASIS_Out) then
818                 namsrc_checkused(nf) = 1
819                 if (OASIS_debug >= 20) then
820                    write(nulprt,*) subname,' set src checkused ',trim(myfld),':',trim(otfld),nf
821                    call oasis_flush(nulprt)
822                 endif
823              endif
824              if (prism_var(nv1)%ops == OASIS_In) then
825                 n1 = 0
826                 found2 = .false.
827                 do while (n1 < sortnsrc%num .and. .not.found2)
828                    n1 = n1 + 1
829                    if (nn == sortnsrc%namnum(n1) .and. myfldi == sortnsrc%fldnum(n1)) then
830                       namsrc_checkused(n1) = 1
831                       found2 = .true.
832                       if (OASIS_debug >= 20) then
833                          write(nulprt,*) subname,' set dst checkused ',trim(myfld),':',trim(otfld),n1
834                          call oasis_flush(nulprt)
835                       endif
836                    endif
837                 enddo
838              endif
839
840              if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n3d')
841              nm = sortvars%modnum(nvf)
842              nv = sortvars%varnum(nvf)
843             
844              if (OASIS_debug >= 20) then
845                 write(nulprt,*) subname,' match otfld ',trim(otfld),nn
846                 call oasis_flush(nulprt)
847              endif
848
849              !--------------------------------
850              !>       * Check that one side is In and other side is Out for communication
851              !--------------------------------
852
853              if (namfldops(nn) == ip_exported .or. namfldops(nn) == ip_expout) then
854! tcraig allow this now
855!                if (nm == compid) then
856!                   write(nulprt,*) subname,estr,'send recv pair on same model = ', &
857!                      trim(myfld),' ',trim(otfld)
858!                   call oasis_abort(file=__FILE__,line=__LINE__)
859!                endif
860                 if (flag == OASIS_Out .and. allops(nv,nm) /= OASIS_In) then
861                    write(nulprt,*) subname,estr,'send recv pair both Out = ', &
862                       trim(myfld),' ',trim(otfld)
863                    call oasis_abort(file=__FILE__,line=__LINE__)
864                 endif
865                 if (flag == OASIS_In .and. allops(nv,nm) /= OASIS_Out) then
866                    write(nulprt,*) subname,estr,'send recv pair both In = ', &
867                       trim(myfld),' ',trim(otfld)
868                    call oasis_abort(file=__FILE__,line=__LINE__)
869                 endif
870              endif
871
872              !--------------------------------
873              !>       * Check if input or output, field name should match on both sides.
874              !--------------------------------
875
876              if (namfldops(nn) == ip_input .or. namfldops(nn) == ip_output) then
877                 if (trim(myfld) /= trim(otfld)) then
878                    write(nulprt,*) subname,estr,'namcouple field names do not match in/out = ', &
879                       trim(myfld),' ',trim(otfld)
880                    call oasis_abort(file=__FILE__,line=__LINE__)
881                 endif
882              endif
883
884              !--------------------------------
885              !>       * Check that the bundle size matches in both models for bundled fields
886              !--------------------------------
887
888              if (prism_var(nv1)%num /= allnum(nv,nm)) then
889                 write(nulprt,*) subname,estr,'namcouple bundle fields do not match for ',trim(myfld),' ',trim(otfld)
890                 write(nulprt,*) subname,estr,'namcouple bundle numbers are ',prism_var(nv1)%num,allnum(nv,nm)
891                 call oasis_abort(file=__FILE__,line=__LINE__)
892              endif
893
894              !--------------------------------
895              ! Only an error to find two sources for a destination
896              ! Not an error if a two destinations have a single source
897              !--------------------------------
898
899              if (flag == OASIS_In .and. found) then
900                 write(nulprt,*) subname,estr,'found two sources for field = ',trim(otfld)
901                 call oasis_abort(file=__FILE__,line=__LINE__)
902              endif
903              found = .true.
904
905              nc = nns
906              if (flag == OASIS_Out) pcpointer => prism_coupler_put(nc)
907              if (flag == OASIS_In)  pcpointer => prism_coupler_get(nc)
908
909              !--------------------------------
910              !>       * Generate field list, multiple field support
911              !--------------------------------
912
913              IF (OASIS_debug >= 20) THEN
914                  WRITE(nulprt,*) subname,' set prism_coupler '
915                  CALL oasis_flush(nulprt)
916              ENDIF
917
918              ! tcraig, changed this to make sure order of fields in list matches on all tasks
919              ! Use the field lists in the namcouple
920              ! Assumes all namcoupler variables are coupled
921              ! The nflds counter doesn't do much anymore here
922              ! The varid size should be size(myfldlist)
923              ! Will need to change IF all namcoupler variables don't need to be coupled
924
925              pcpointer%nflds = pcpointer%nflds + 1
926
927!tcx
928! this used to add fields to list one at a time
929!              svarid = size(pcpointer%varid)
930!              if (pcpointer%nflds > svarid) then
931!                  allocate(varidtmp(svarid))
932!                  varidtmp(1:svarid) = pcpointer%varid(1:svarid)
933!                  deallocate(pcpointer%varid)
934!                  allocate(pcpointer%varid(pcpointer%nflds+10))
935!                  pcpointer%varid(1:svarid) = varidtmp(1:svarid)
936!                  deallocate(varidtmp)
937!              endif
938!
939!              if (pcpointer%nflds == 1) then
940!                 pcpointer%fldlist = trim(myfld)
941!              else
942!                 pcpointer%fldlist = trim(pcpointer%fldlist)//':'//trim(myfld)
943!              endif
944!              pcpointer%varid(pcpointer%nflds) = nv1
945!tcx
946
947              if (pcpointer%nflds == 1) then
948                 pcpointer%fldlist = trim(myfldlist)
949                 deallocate(pcpointer%varid)
950                 allocate(pcpointer%varid(oasis_string_listGetNum(myfldlist)))
951                 pcpointer%varid(:) = ispval
952              endif
953
954              svarid = size(pcpointer%varid)
955              if (myfldi > svarid .or. pcpointer%nflds > svarid) then
956                 WRITE(nulprt,*) subname,estr,'multiple field coupling setup error',svarid,myfldi,pcpointer%nflds
957                 call oasis_abort(file=__FILE__,line=__LINE__)
958              endif
959
960              pcpointer%varid(myfldi) = nv1
961
962              !--------------------------------
963              !>       * Add this coupler to list of prism_var couplers
964              !--------------------------------
965
966              prism_var(nv1)%ncpl = prism_var(nv1)%ncpl + 1
967              if (prism_var(nv1)%ncpl > mvarcpl) then
968                 WRITE(nulprt,*) subname,estr,'ncpl too high, max size (mvarcpl) = ',mvarcpl
969                 WRITE(nulprt,*) subname,estr,'increase mvarcpl in mod_oasis_var'
970                 call oasis_abort(file=__FILE__,line=__LINE__)
971              endif
972              prism_var(nv1)%cpl(prism_var(nv1)%ncpl) = nc
973
974              !--------------------------------
975              ! prism_coupler settings
976              !>       * Copy namcouple settings into this coupler or
977              !>       check that coupler is consistent with prior setting
978              !--------------------------------
979
980              if (pcpointer%valid) then
981                 if (pcpointer%comp /= nm) then
982                    WRITE(nulprt,*) subname,estr,'mismatch in field comp for var = ',trim(myfld)
983                    call oasis_abort(file=__FILE__,line=__LINE__)
984                 endif
985                 if (pcpointer%namID /= nn) then
986                    WRITE(nulprt,*) subname,estr,'mismatch in field namID for var = ',trim(myfld)
987                    call oasis_abort(file=__FILE__,line=__LINE__)
988                 endif
989                 if (pcpointer%partID /= part1) then
990                    WRITE(nulprt,*) subname,estr,'mismatch in field partID for var = ',trim(myfld)
991                    call oasis_abort(file=__FILE__,line=__LINE__)
992                 endif
993
994              else
995                 pcpointer%comp   = nm
996                 pcpointer%seq    = namfldseq(nn)
997                 pcpointer%dt     = namflddti(nn)
998                 pcpointer%lag    = namfldlag(nn)
999                 pcpointer%maxtime= namruntim
1000                 pcpointer%rstfile= trim(namrstfil(nn))
1001                 pcpointer%inpfile= trim(naminpfil(nn))
1002                 pcpointer%mapperID = -1
1003                 pcpointer%partID = part1
1004                 pcpointer%rpartID= part1
1005                 pcpointer%namID  = nn
1006                 pcpointer%trans  = namfldtrn(nn)
1007                 pcpointer%conserv= namfldcon(nn)
1008                 pcpointer%consopt= namfldcoo(nn)
1009                 pcpointer%ops    = namfldops(nn)
1010                 pcpointer%tag    = compid*100*1000 + compid*1000 + nn
1011                 pcpointer%getput = OASIS_NotDef
1012                 pcpointer%sndrcv = .false.
1013                 pcpointer%output = .false.
1014                 pcpointer%input  = .false.
1015                 pcpointer%sndmult= namfldsmu(nn)
1016                 pcpointer%sndadd = namfldsad(nn)
1017                 pcpointer%rcvmult= namflddmu(1,nn)
1018                 pcpointer%rcvadd = namflddad(1,nn)
1019                 pcpointer%snddiag= namchecki(nn)
1020                 pcpointer%rcvdiag= namchecko(nn)
1021                 pcpointer%rfno   = namflddno(nn) - 1
1022                 if (namflddno(nn) > 1) then
1023                    allocate(pcpointer%rfna  (pcpointer%rfno))
1024                    allocate(pcpointer%rfmult(pcpointer%rfno))
1025                    allocate(pcpointer%rfadd (pcpointer%rfno))
1026                    allocate(pcpointer%rfcpl (pcpointer%rfno))
1027                    allocate(pcpointer%rffnum(pcpointer%rfno))
1028                    do n1 = 2,namflddno(nn)
1029                       pcpointer%rfna  (n1-1) = namflddna(n1,nn)
1030                       pcpointer%rfmult(n1-1) = namflddmu(n1,nn)
1031                       pcpointer%rfadd (n1-1) = namflddad(n1,nn)
1032                       pcpointer%rfcpl (n1-1) = -1
1033                       pcpointer%rffnum(n1-1) = -1
1034                    enddo
1035                 endif
1036
1037                 !--------------------------------
1038                 !>       * Set prism_coupler input and output flags
1039                 ! prism_coupler comm flags, need for tags to match up on both sides
1040                 ! tags assume up to 1000 namcouple inputs and 100 models
1041                 !--------------------------------
1042
1043                 IF (OASIS_debug >= 20) THEN
1044                    WRITE(nulprt,*) subname,' inout flags '
1045                    CALL oasis_flush(nulprt)
1046                 ENDIF
1047
1048                 if (namfldops(nn) == ip_output .or. namfldops(nn) == ip_expout) then
1049                    pcpointer%output = .true.
1050                    pcpointer%getput = OASIS3_PUT
1051                 endif
1052                 if (namfldops(nn) == ip_input) then
1053                    pcpointer%input  = .true.
1054                    pcpointer%getput = OASIS3_GET
1055                 endif
1056
1057                 if (namfldops(nn) == ip_exported .or. namfldops(nn) == ip_expout) then
1058                    pcpointer%sndrcv = .true.
1059                    if (flag == OASIS_Out) then
1060                       pcpointer%tag = nm*100*1000 + compid*1000 + nn
1061                       pcpointer%getput = OASIS3_PUT
1062                    elseif (flag == OASIS_In) then
1063                       pcpointer%tag = compid*100*1000 + nm*1000 + nn
1064                       pcpointer%getput = OASIS3_GET
1065                    endif
1066                    !--------------------------------
1067                    !>       * Setup prism_coupler router
1068                    ! cannot reuse router because don't really know what's on the other side
1069                    ! if router is already set for the coupler, then fine, otherwise, set new router
1070                    !--------------------------------
1071                    if (pcpointer%routerID == ispval) then
1072                       prism_nrouter = prism_nrouter+1
1073                       if (prism_nrouter > prism_mrouter) then
1074                          write(nulprt,*) subname,estr,'prism_nrouter too large = ',prism_nrouter,prism_mrouter
1075                          write(nulprt,*) subname,estr,'check prism_mrouter in oasis_coupler_setup '
1076                          call oasis_abort(file=__FILE__,line=__LINE__)
1077                       endif
1078                       pcpointer%routerID = prism_nrouter
1079                    endif
1080                 endif
1081
1082                 !--------------------------------
1083                 !>       * Setup prism_coupler mapper
1084                 !--------------------------------
1085
1086                 IF (OASIS_debug >= 20) THEN
1087                    WRITE(nulprt,*) subname,' mapper '
1088                    CALL oasis_flush(nulprt)
1089                 ENDIF
1090
1091                 tmp_mapfile = nammapfil(nn)
1092                 tmp_mapfile2 = ''
1093
1094                 if (TRIM(tmp_mapfile) == 'idmap' .and. TRIM(namscrmet(nn)) /= TRIM(cspval)) then
1095                    if (trim(namscrmet(nn)) == 'CONSERV') then
1096                       tmp_mapfile = 'rmp_'//trim(namsrcgrd(nn))//'_to_'//trim(namdstgrd(nn))//&
1097                          &'_'//TRIM(namscrmet(nn))//'_'//TRIM(namscrnor(nn))//'.nc'
1098                    elseif (namscrnbr(nn) > 0) then
1099                       write(tmpstr,'(i0)') namscrnbr(nn)
1100                       tmp_mapfile = 'rmp_'//trim(namsrcgrd(nn))//'_to_'//trim(namdstgrd(nn))//&
1101                          &'_'//TRIM(namscrmet(nn))//'_'//TRIM(tmpstr)//'.nc'
1102                       tmp_mapfile2 = 'rmp_'//trim(namsrcgrd(nn))//'_to_'//trim(namdstgrd(nn))//&
1103                          &'_'//TRIM(namscrmet(nn))//'.nc'
1104                    else
1105                       tmp_mapfile = 'rmp_'//trim(namsrcgrd(nn))//'_to_'//trim(namdstgrd(nn))//&
1106                          &'_'//trim(namscrmet(nn))//'.nc'
1107                    endif
1108                 endif
1109
1110                 if (trim(tmp_mapfile) /= 'idmap') then
1111                    pcpointer%maploc = trim(nammaploc(nn))
1112                    if ((flag == OASIS_In  .and. trim(nammaploc(nn)) == 'dst') .or. &
1113                        (flag == OASIS_Out .and. trim(nammaploc(nn)) == 'src')) then
1114                       !--------------------------------
1115                       !>       * Try to reuse mapper already defined,
1116                       !>       must match mapping file and partition
1117                       !--------------------------------
1118                       mapID = -1
1119                       do n = 1,prism_nmapper
1120                          if (trim(prism_mapper(n)%file) == trim(tmp_mapfile) .and. &
1121                              trim(prism_mapper(n)%loc ) == trim(nammaploc(nn)) .and. &
1122                              trim(prism_mapper(n)%opt ) == trim(nammapopt(nn))) then
1123                             if (flag == OASIS_In  .and. prism_mapper(n)%dpart == part1) mapID = n
1124                             if (flag == OASIS_Out .and. prism_mapper(n)%spart == part1) mapID = n
1125                          endif
1126                       enddo
1127                       ! check if tmp_mapfile2 matches if tmp_mapfile does not
1128                       if (mapID < 1 .and. tmp_mapfile2 /= '') then
1129                       do n = 1,prism_nmapper
1130                          if (trim(prism_mapper(n)%file)== trim(tmp_mapfile2) .and. &
1131                              trim(prism_mapper(n)%loc ) == trim(nammaploc(nn)) .and. &
1132                              trim(prism_mapper(n)%opt ) == trim(nammapopt(nn))) then
1133                             if (flag == OASIS_In  .and. prism_mapper(n)%dpart == part1) mapID = n
1134                             if (flag == OASIS_Out .and. prism_mapper(n)%spart == part1) mapID = n
1135                          endif
1136                       enddo
1137                       endif
1138                       !--------------------------------
1139                       !>       * Or get ready to initialize a new mapper
1140                       !--------------------------------
1141                       if (mapID < 1) then
1142                          prism_nmapper = prism_nmapper + 1
1143                          if (prism_nmapper > prism_mmapper) then
1144                             write(nulprt,*) subname,estr,'prism_nmapper too large',prism_nmapper,prism_mmapper
1145                             write(nulprt,*) subname,estr,'check prism_mmapper in oasis_coupler_setup '
1146                             call oasis_abort(file=__FILE__,line=__LINE__)
1147                          endif
1148                          mapID = prism_nmapper
1149                          prism_mapper(mapID)%file = trim(tmp_mapfile)
1150                          prism_mapper(mapID)%file2= trim(tmp_mapfile2)
1151                          prism_mapper(mapID)%loc  = trim(nammaploc(nn))
1152                          prism_mapper(mapID)%opt  = trim(nammapopt(nn))
1153                          prism_mapper(mapID)%srcgrid = trim(namsrcgrd(nn))
1154                          prism_mapper(mapID)%dstgrid = trim(namdstgrd(nn))
1155                          if (flag == OASIS_In ) prism_mapper(mapID)%dpart = part1
1156                          if (flag == OASIS_Out) prism_mapper(mapID)%spart = part1
1157                          if (OASIS_debug > 15) then
1158                             write(nulprt,*) subname,' DEBUG new mapper for file ',&
1159                                             trim(prism_mapper(mapID)%file)
1160                             call oasis_flush(nulprt)
1161                          endif
1162                       endif
1163                       pcpointer%mapperID = mapID
1164                    endif  ! flag and nammaploc match
1165                 endif  ! nammapfil
1166
1167                 pcpointer%valid  = .true.
1168
1169              endif  !  valid
1170
1171              if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n3d')
1172
1173           enddo  ! nvf
1174
1175        endif  ! my var found
1176
1177     enddo  ! nfind
1178  enddo  ! nv1
1179  if (local_timers_on >= 2) call oasis_timer_stop ('cpl_setup_n3')
1180  if (local_timers_on >= 1) then
1181     call oasis_timer_start('cpl_setup_n4_barrier')
1182     call oasis_mpi_barrier(mpi_comm_global, 'cpl_setup_n4')
1183     call oasis_timer_stop('cpl_setup_n4_barrier')
1184     call oasis_timer_start('cpl_setup_n4')
1185  endif
1186  if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4a')
1187
1188  ! aggregate checkused info across all pes and then check on each component root
1189  allocate(namsrc_checkused_g(sortnsrc%num))
1190  call oasis_mpi_max(namsrc_checkused,namsrc_checkused_g,mpi_comm_global,string=trim(subname)//':srccheckused',all=.true.)
1191  found = .false.
1192  do n1 = 1,sortnsrc%num
1193     if (namsrc_checkused_g(n1) /= 1) then
1194        if (mpi_rank_local == 0) write(nulprt,*) subname,estr,'namcouple variable not used: ',trim(sortnsrc%fld(n1))
1195        found = .true.
1196     endif
1197  enddo
1198!  call oasis_mpi_barrier(mpi_comm_global)
1199  if (found) call oasis_abort(file=__FILE__,line=__LINE__)
1200  deallocate(namsrc_checkused_g)
1201
1202  !--- deallocate temporary ---
1203  deallocate(allvar,nallvar,allops,allnum)
1204  deallocate(namsrc_checkused)
1205  deallocate(sortnsrc%fld)
1206  deallocate(sortnsrc%namnum)
1207  deallocate(sortnsrc%fldnum)
1208  deallocate(sortndst%fld)
1209  deallocate(sortndst%namnum)
1210  deallocate(sortndst%fldnum)
1211  deallocate(sortvars%fld)
1212  deallocate(sortvars%modnum)
1213  deallocate(sortvars%varnum)
1214
1215  !----------------------------------------------------------
1216  !> * Rebuild the fields list based on field bundles as needed
1217  ! need to modify fldlist and varids in prism_couplers
1218  ! order of fields needs to be preserved
1219  !----------------------------------------------------------
1220
1221  do nc = 1,prism_mcoupler
1222  do npc = 1,2
1223     if (npc == 1) then
1224        pcpointer => prism_coupler_put(nc)
1225     endif
1226     if (npc == 2) then
1227        pcpointer => prism_coupler_get(nc)
1228     endif
1229
1230     if (pcpointer%valid) then
1231        nflds1 = oasis_string_listGetNum(pcpointer%fldlist)
1232        nflds2 = 0
1233        do n1 = 1,nflds1
1234           nflds2 = nflds2 + prism_var(pcpointer%varid(n1))%num
1235        enddo
1236        if (OASIS_debug >= 2) then
1237           write(nulprt,*) subname,' fldlist rebuild nflds1,nflds2 for ',trim(pcpointer%fldlist)
1238           write(nulprt,*) subname,' fldlist rebuild nflds1,nflds2 ',nflds1,nflds2
1239        endif
1240        if (nflds2 < nflds1) then
1241           write(nulprt,*) subname,estr,'fldlist rebuild nflds2 < nflds1 for ',trim(pcpointer%fldlist)
1242           write(nulprt,*) subname,estr,'fldlist reset error in fld cnt = ',nflds1,nflds2
1243           call oasis_abort(file=__FILE__,line=__LINE__)
1244        else
1245           if (OASIS_debug >= 2) then
1246              write(nulprt,*) subname,' fldlist rebuild nflds2 > nflds1 for ',trim(pcpointer%fldlist)
1247           endif
1248           allocate(varid1(nflds1))
1249           varid1(1:nflds1) = pcpointer%varid(1:nflds1)
1250           myfldlist = pcpointer%fldlist  ! temporary storage
1251           pcpointer%fldlist = ""
1252           deallocate(pcpointer%varid)
1253           allocate(pcpointer%varid(nflds2))
1254           ncnt = 0
1255           do n1 = 1,nflds1
1256              do n2 = 1,prism_var(varid1(n1))%num
1257                 ncnt = ncnt + 1
1258                 pcpointer%varid(ncnt) = varid1(n1)
1259                 delim = ":"
1260                 if (ncnt == 1) delim = ""
1261                 if (len_trim(pcpointer%fldlist) > 0.99 * len(pcpointer%fldlist)) then
1262                    write(nulprt,*) subname,estr,'fldlist rebuild too long, limit is ',len(pcpointer%fldlist),' chars'
1263                    write(nulprt,*) subname,estr,'current rebuid fldlist is ',trim(pcpointer%fldlist)
1264                    call oasis_abort(file=__FILE__,line=__LINE__)
1265                 endif
1266                 call oasis_coupler_bldvarname(varid1(n1),n2,vname)
1267                 write(pcpointer%fldlist,'(a)') trim(pcpointer%fldlist)//trim(delim)//trim(vname)
1268                 if (OASIS_debug >= 2) then
1269                    write(nulprt,*) subname,' fldlist rebuild n1, n2 ',n1,n2,ncnt
1270                    write(nulprt,*) subname,' fldlist rebuild fldlist ',ncnt,trim(pcpointer%fldlist)
1271                    write(nulprt,*) subname,' fldlist rebuild varid ',ncnt,pcpointer%varid(ncnt)
1272                 endif
1273              enddo  ! n2
1274           enddo  ! n1
1275           deallocate(varid1)
1276        endif
1277     endif  ! valid
1278  enddo  ! npc
1279  enddo  ! nc
1280
1281  if (OASIS_debug >= 20) then
1282     write(nulprt,*) ' '
1283     write(nulprt,*) subname,' couplers setup'
1284     do nc = 1,prism_mcoupler
1285!tcx can't write here, something uninitialized???
1286!-> CEG it was dpart so added extra if into the print routine
1287!        if (prism_coupler_put(nc)%valid) call prism_coupler_print(nc,prism_coupler_put(nc))
1288!        if (prism_coupler_get(nc)%valid) call prism_coupler_print(nc,prism_coupler_get(nc))
1289     enddo
1290     write(nulprt,*) ' '
1291     call oasis_flush(nulprt)
1292  endif
1293
1294  if (mpi_comm_local == MPI_COMM_NULL) then
1295     return
1296  endif
1297
1298  !----------------------------------------------------------
1299  !> * Initialize coupling infrastructure based on initial coupler setup above
1300  !----------------------------------------------------------
1301
1302  call oasis_debug_note(subname//' initialize coupling datatypes')
1303  if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4a')
1304
1305  !----------------------------------------------------------
1306  !> * Loop over all couplers
1307  !----------------------------------------------------------
1308
1309  do nc = 1,prism_mcoupler
1310  ! tcraig, this barrier make sure mapping files are generated one coupler at a time
1311  if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4_global_barrier')
1312  call oasis_mpi_barrier(mpi_comm_global,'cpl_setup_n4_global')
1313  if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4_global_barrier')
1314  do npc = 1,2
1315   if (npc == 1) then
1316     pcpointer => prism_coupler_put(nc)
1317     pcpntpair => prism_coupler_get(nc)
1318   endif
1319   if (npc == 2) then
1320     pcpointer => prism_coupler_get(nc)
1321     pcpntpair => prism_coupler_put(nc)
1322   endif
1323   if (OASIS_debug >= 20) then
1324     write(nulprt,*) subname,' DEBUG cb:initialize coupler ',nc,npc,pcpointer%valid
1325     call oasis_flush(nulprt)
1326   endif
1327
1328   if (pcpointer%valid) then
1329     if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4b')
1330     if (OASIS_debug >= 5) then
1331        write(nulprt,*) subname,' DEBUG ci:initialize coupler ',nc,npc
1332        call oasis_flush(nulprt)
1333     endif
1334
1335     namID = pcpointer%namID
1336     part1 = pcpointer%partID
1337     mapID = pcpointer%mapperID
1338
1339     if (part1 <= 0) then
1340        write(nulprt,*) subname,estr,'part1 invalid = ',part1
1341        call oasis_abort(file=__FILE__,line=__LINE__)
1342     endif
1343
1344     !--------------------------------
1345     !>   * Initialize avect1 which stores the get/put data
1346     !--------------------------------
1347
1348     gsize = mct_gsmap_gsize(prism_part(part1)%gsmap)
1349     lsize = mct_gsmap_lsize(prism_part(part1)%gsmap,mpi_comm_local)
1350     if (OASIS_debug >= 15) then
1351        write(nulprt,'(1x,2a,5i10)') subname,' DEBUG ci:part1 info ',namID,part1,mapID,gsize,lsize
1352        write(nulprt,'(1x,2a,4i12)') subname,' DEBUG ci:part1a',prism_part(part1)%gsmap%ngseg,&
1353                                     prism_part(part1)%gsmap%gsize
1354        do n1 = 1,min(prism_part(part1)%gsmap%ngseg,10)
1355           write(nulprt,'(1x,2a,4i12)') subname,' DEBUG ci:part1b',n1,&
1356                                        prism_part(part1)%gsmap%start(n1),&
1357                                        prism_part(part1)%gsmap%length(n1),&
1358                                        prism_part(part1)%gsmap%pe_loc(n1)
1359        enddo
1360        call oasis_flush(nulprt)
1361     endif
1362     call mct_avect_init(pcpointer%avect1,rList=trim(pcpointer%fldlist),lsize=lsize)
1363     call mct_avect_zero(pcpointer%avect1)
1364     pcpointer%aVon(1) = .true.
1365     if (OASIS_debug >= 15) then
1366        write(nulprt,*) subname,' DEBUG ci:avect1 initialized '
1367        call oasis_flush(nulprt)
1368     endif
1369
1370     !--------------------------------
1371     !>   * Compute nflds for this coupling and initialize avcnt and status
1372     !--------------------------------
1373
1374     pcpointer%nflds = mct_aVect_nRAttr(pcpointer%avect1)
1375     allocate(pcpointer%status(pcpointer%nflds))
1376     allocate(pcpointer%avcnt (pcpointer%nflds))
1377     pcpointer%avcnt(:) = 0
1378     if (pcpointer%getput == OASIS3_PUT) pcpointer%status = OASIS_COMM_WAIT
1379     if (pcpointer%getput == OASIS3_GET) pcpointer%status = OASIS_COMM_READY
1380     if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4b')
1381
1382     !--------------------------------
1383     !>   * Initialize the mapper data
1384     !--------------------------------
1385
1386     if (mapID > 0) then
1387
1388        if (prism_mapper(mapID)%init) then
1389           if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4c')
1390           !--------------------------------
1391           ! mapper already initialized
1392           !--------------------------------
1393           if (pcpointer%getput == OASIS3_PUT) then
1394              part2 = prism_mapper(mapID)%dpart
1395           else
1396              part2 = prism_mapper(mapID)%spart
1397           endif
1398           gsize = mct_gsmap_gsize(prism_part(part2)%gsmap)
1399           if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4c')
1400        else
1401           !--------------------------------
1402           ! instantiate mapper
1403           ! read/generate mapping file
1404           ! read non local grid size
1405           ! get gsmap for non local grid
1406           ! read mapping weights and initialize sMatP
1407           !--------------------------------
1408           if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4d')
1409           if (OASIS_debug >= 15) then
1410              write(nulprt,*) subname,' DEBUG ci:read mapfile ',trim(prism_mapper(mapID)%file)
1411              call oasis_flush(nulprt)
1412           endif
1413           if (mpi_rank_local == 0) then
1414              if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4da')
1415              if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4da1')
1416              inquire(file=trim(prism_mapper(mapID)%file),exist=exists)
1417              if (.not. exists .and. prism_mapper(mapid)%file2 /= '') then
1418                 ! if file2 exists, but not file, set file = file2
1419                 inquire(file=trim(prism_mapper(mapID)%file2),exist=exists)
1420                 if (exists) then
1421                    write(nulprt,*) subname,' found old mapname, using ', &
1422                       trim(prism_mapper(mapID)%file2),' instead of ', trim(prism_mapper(mapID)%file)
1423                    prism_mapper(mapID)%file = prism_mapper(mapID)%file2
1424                 endif
1425              endif
1426              if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4da1')
1427              if (OASIS_debug >= 15) then
1428                 write(nulprt,*) subname,' DEBUG ci: inquire mapfile ',&
1429                                 trim(prism_mapper(mapID)%file),exists
1430                 call oasis_flush(nulprt)
1431              endif
1432           endif
1433           call oasis_mpi_bcast(exists,mpi_comm_local,subname//' exists')
1434           if (.not.exists) then
1435              if (trim(namscrmet(namID)) /= trim(cspval)) then
1436                 !--------------------------------
1437                 ! generate mapping file on map group/communicator
1438                 ! taken from oasis3 scriprmp
1439                 !--------------------------------
1440                 if (mpi_in_map) then
1441                    if (local_timers_on > 2) call oasis_timer_start('cpl_setup_genmap')
1442                    call oasis_map_genmap(mapID,namID)
1443                    if (local_timers_on > 2) call oasis_timer_stop('cpl_setup_genmap')
1444                 end if
1445              else
1446                 write(nulprt,*) subname,estr,'map file does not exist and SCRIPR not set = ',&
1447                                 trim(prism_mapper(mapID)%file)
1448                 call oasis_abort(file=__FILE__,line=__LINE__)
1449              endif
1450           endif
1451           !
1452           if (mpi_rank_local == 0) then
1453              !--------------------------------
1454              ! open mapping file
1455              !--------------------------------
1456              if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4da3')
1457              status = nf90_open(trim(prism_mapper(mapID)%file),nf90_nowrite,ncid)
1458              if (status /= NF90_NOERR) then
1459                 write(nulprt,*) subname,' nf90_strerror = ',trim(nf90_strerror(status))
1460                 write(nulprt,*) subname,estr,'file not found = ',trim(prism_mapper(mapID)%file)
1461                 call oasis_abort(file=__FILE__,line=__LINE__)
1462              endif
1463              if (OASIS_debug >= 15) then
1464                 status = nf90_inq_dimid(ncid,'dst_grid_size',dimid)
1465                 if (status /= NF90_NOERR) then
1466                    write(nulprt,*) subname,' nf90_strerror = ',trim(nf90_strerror(status))
1467                    write(nulprt,*) subname,estr,'dim not found = ','dst_grid_size'
1468                    call oasis_abort(file=__FILE__,line=__LINE__)
1469                 endif
1470                 status = nf90_inquire_dimension(ncid,dimid,len=gsize)
1471                 if (status /= NF90_NOERR) then
1472                    write(nulprt,*) subname,' nf90_strerror = ',trim(nf90_strerror(status))
1473                    call oasis_abort(file=__FILE__,line=__LINE__)
1474                 endif
1475                 write(nulprt,*) subname," DEBUG dst_grid_size ",gsize
1476                 status = nf90_inq_dimid(ncid,'src_grid_size',dimid)
1477                 if (status /= NF90_NOERR) then
1478                    write(nulprt,*) subname,' nf90_strerror = ',trim(nf90_strerror(status))
1479                    write(nulprt,*) subname,estr,'dim not found = ','src_grid_size'
1480                    call oasis_abort(file=__FILE__,line=__LINE__)
1481                 endif
1482                 status = nf90_inquire_dimension(ncid,dimid,len=gsize)
1483                 if (status /= NF90_NOERR) then
1484                    write(nulprt,*) subname,' nf90_strerror = ',trim(nf90_strerror(status))
1485                    call oasis_abort(file=__FILE__,line=__LINE__)
1486                 endif
1487                 write(nulprt,*) subname," DEBUG src_grid_size ",gsize
1488              endif
1489              if (pcpointer%getput == OASIS3_PUT) then
1490                 status = nf90_inq_dimid(ncid,'dst_grid_size',dimid)
1491                 if (status /= NF90_NOERR) then
1492                    write(nulprt,*) subname,' nf90_strerror = ',trim(nf90_strerror(status))
1493                    write(nulprt,*) subname,estr,'dim not found = ','dst_grid_size'
1494                    call oasis_abort(file=__FILE__,line=__LINE__)
1495                 endif
1496              elseif (pcpointer%getput == OASIS3_GET) then
1497                 status = nf90_inq_dimid(ncid,'src_grid_size',dimid)
1498                 if (status /= NF90_NOERR) then
1499                    write(nulprt,*) subname,' nf90_strerror = ',trim(nf90_strerror(status))
1500                    write(nulprt,*) subname,estr,'dim not found = ','src_grid_size'
1501                    call oasis_abort(file=__FILE__,line=__LINE__)
1502                 endif
1503              endif
1504              status = nf90_inquire_dimension(ncid,dimid,len=gsize)
1505              if (status /= NF90_NOERR) then
1506                 write(nulprt,*) subname,' nf90_strerror = ',trim(nf90_strerror(status))
1507                 call oasis_abort(file=__FILE__,line=__LINE__)
1508              endif
1509              if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4da3')
1510              if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4da')
1511           endif  ! rank = 0
1512           if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4db')
1513           call oasis_mpi_bcast(gsize,mpi_comm_local,subname//' gsize')
1514           if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4db')
1515
1516           if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4dc')
1517           if (pcpointer%getput == OASIS3_PUT) then
1518              nx = namdst_nx(namID)
1519              ny = namdst_ny(namID)
1520              gridname = trim(namdstgrd(namID))
1521              mapopt1 = 'dst'
1522           else
1523              nx = namsrc_nx(namID)
1524              ny = namsrc_ny(namID)
1525              gridname = trim(namsrcgrd(namID))
1526              mapopt1 = 'src'
1527           endif
1528
1529           !--- mapopt sets whether src or dst are rearranged in remap
1530           !--- src = rearrange and map (bfb), dst = map and rearrange (partial sum)
1531           if (prism_mapper(mapID)%opt == 'opt') then
1532              if (pcpointer%getput == OASIS3_PUT .and. prism_part(part1)%gsize > gsize) then
1533                 mapopt = 'dst'
1534              elseif (pcpointer%getput == OASIS3_GET .and. prism_part(part1)%gsize < gsize) then
1535                 mapopt = 'dst'
1536              else
1537                 mapopt = 'src'
1538              endif
1539           elseif (prism_mapper(mapID)%opt == 'bfb') then
1540              mapopt = 'src'
1541           elseif (prism_mapper(mapID)%opt == 'sum') then
1542              mapopt = 'dst'
1543           else
1544              write(nulprt,*) subname,estr,'mapper opt invalid expect bfb or sum  =',trim(prism_mapper(mapID)%opt)
1545              call oasis_abort(file=__FILE__,line=__LINE__)
1546           endif
1547
1548           if (prism_mapper(mapID)%optval /= '' .and. &
1549               prism_mapper(mapID)%optval /= trim(mapopt)) then
1550              write(nulprt,*) subname,estr,'mapper opt changed',&
1551                              trim(prism_mapper(mapID)%optval),' ',trim(mapopt)
1552              call oasis_abort(file=__FILE__,line=__LINE__)
1553           endif
1554           prism_mapper(mapID)%optval = trim(mapopt)
1555           if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4dc')
1556           if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4d')
1557
1558           !-------------------------------
1559           ! smatreaddnc allocates sMati to nwgts
1560           ! then instantiate an sMatP for each set of wgts
1561           ! to support higher order mapping
1562           !-------------------------------
1563
1564           if (OASIS_debug > 15) then
1565              write(nulprt,*) subname,' using part2decomp = ',trim(part2decomp)
1566              write(nulprt,*) subname,' mapopt, mapopt1 = ',trim(mapopt),' ',trim(mapopt1)
1567           endif
1568
1569           if (trim(part2decomp) == 'decomp_wghtfile') then
1570              ! reads wgts based on part1 then creates new part2 decomp based on wgts decomp
1571              ! pass in part1 decomp for "both" decomps for smatreaddnc
1572
1573              if (smatread_method == "ceg") then
1574                 if (local_timers_on >= 1) call oasis_timer_start('cpl_setup_n4rdw_ceg')
1575                 call oasis_map_smatreaddnc_ceg(sMati,prism_part(part1)%gsmap,prism_part(part1)%gsmap, &
1576                    trim(mapopt1),trim(prism_mapper(mapID)%file),mpi_rank_local,mpi_comm_local,nwgts)
1577                 if (local_timers_on >= 1) call oasis_timer_stop('cpl_setup_n4rdw_ceg')
1578              else
1579                 if (local_timers_on >= 1) call oasis_timer_start('cpl_setup_n4rdw_orig')
1580                 call oasis_map_smatreaddnc_orig(sMati,prism_part(part1)%gsmap,prism_part(part1)%gsmap, &
1581                    trim(mapopt1),trim(prism_mapper(mapID)%file),mpi_rank_local,mpi_comm_local,nwgts)
1582                 if (local_timers_on >= 1) call oasis_timer_stop('cpl_setup_n4rdw_orig')
1583              endif
1584
1585              ! extract the gridIDs from sMati, rows and cols are same in all sMati's here, so just use sMati(1)
1586!              if (OASIS_debug > 15) then
1587!                 nullify(gridID)
1588!                 call mct_sMat_ExpGRowI(sMati(1), gridID, length=arrlen)
1589!                 write(nulprt,*) subname,' gridID0r ',trim(mapopt),' ',trim(mapopt1)
1590!                 write(nulprt,*) subname,' gridID1r ',size(gridID),arrlen, gsize
1591!                 if (arrlen > 0) then
1592!                    write(nulprt,*) subname,' gridID2r ',minval(gridID),maxval(gridID)
1593!                    write(nulprt,*) subname,' gridID3r ',minval(gridID(1:arrlen)),maxval(gridID(1:arrlen))
1594!                    write(nulprt,*) subname,' gridID4r ',gridID(1:10)
1595!                 endif
1596!                 deallocate(gridID)
1597!                 nullify(gridID)
1598!                 call mct_sMat_ExpGColI(sMati(1), gridID, length=arrlen)
1599!                 write(nulprt,*) subname,' gridID0c ',trim(mapopt),' ',trim(mapopt1)
1600!                 write(nulprt,*) subname,' gridID1c ',size(gridID),arrlen, gsize
1601!                 if (arrlen > 0) then
1602!                    write(nulprt,*) subname,' gridID2c ',minval(gridID),maxval(gridID)
1603!                    write(nulprt,*) subname,' gridID3c ',minval(gridID(1:arrlen)),maxval(gridID(1:arrlen))
1604!                    write(nulprt,*) subname,' gridID4c ',gridID(1:10)
1605!                 endif
1606!                 deallocate(gridID)
1607!              endif
1608
1609              if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4smat_expg')
1610              nullify(gridID)
1611              if (mapopt1 == 'dst') then
1612                 call mct_sMat_ExpGRowI(sMati(1), gridID, length=arrlen)
1613              elseif (mapopt1 == 'src') then
1614                 call mct_sMat_ExpGColI(sMati(1), gridID, length=arrlen)
1615              else
1616                 write(nulprt,*) subname,estr,'invalid mapopt = ',trim(mapopt)
1617                 call oasis_abort(file=__FILE__,line=__LINE__)
1618              endif
1619              do n = 1,nwgts
1620                 call mct_sMat_Clean(sMati(n))
1621              enddo
1622              deallocate(sMati)
1623              if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4smat_expg')
1624
1625              if (OASIS_debug > 15) then
1626                 write(nulprt,*) subname,' gridID0 ',trim(mapopt),' ',trim(mapopt1)
1627                 write(nulprt,*) subname,' gridID1 ',size(gridID),arrlen, gsize
1628                 if (arrlen > 0) then
1629                    write(nulprt,*) subname,' gridID2 ',minval(gridID),maxval(gridID)
1630                    write(nulprt,*) subname,' gridID3 ',minval(gridID(1:arrlen)),maxval(gridID(1:arrlen))
1631                    write(nulprt,*) subname,' gridID4 ',gridID(1:10)
1632                 endif
1633              endif
1634
1635           endif
1636
1637            if (local_timers_on >= 1) then
1638              call oasis_timer_start('cpl_setup_n4part_cr_barrier')
1639              call oasis_mpi_barrier(mpi_comm_local, 'cpl_setup_n4part')
1640              call oasis_timer_stop('cpl_setup_n4part_cr_barrier')
1641           endif
1642           if (local_timers_on >= 1) call oasis_timer_start('cpl_setup_n4part_create')
1643           if (part2decomp == 'decomp_wghtfile') then
1644              call oasis_part_create(part2,trim(part2decomp),gsize,nx,ny,gridname,prism_part(part1)%mpicom,mpi_comm_local,gridID)
1645              deallocate(gridID)
1646           else
1647              call oasis_part_create(part2,trim(part2decomp),gsize,nx,ny,gridname,prism_part(part1)%mpicom,mpi_comm_local)
1648           endif
1649           if (local_timers_on >= 1) call oasis_timer_stop('cpl_setup_n4part_create')
1650
1651           if (pcpointer%getput == OASIS3_PUT) then
1652             !prism_mapper(mapID)%spart = part1   ! set above
1653              prism_mapper(mapID)%dpart = part2
1654           else
1655              prism_mapper(mapID)%spart = part2
1656             !prism_mapper(mapID)%dpart = part1   ! set above
1657           endif
1658           spart = prism_mapper(mapID)%spart
1659           dpart = prism_mapper(mapID)%dpart
1660
1661           if (smatread_method == "ceg") then
1662              if (local_timers_on >= 1) then
1663                 call oasis_timer_start('cpl_setup_n4rd_ceg_barrier')
1664                 call oasis_mpi_barrier(mpi_comm_local, 'cpl_setup_n4rd_ceg')
1665                 call oasis_timer_stop('cpl_setup_n4rd_ceg_barrier')
1666                 call oasis_timer_start('cpl_setup_n4rd_ceg')
1667              endif
1668              call oasis_map_smatreaddnc_ceg(sMati,prism_part(spart)%gsmap,prism_part(dpart)%gsmap, &
1669                 trim(mapopt),trim(prism_mapper(mapID)%file),mpi_rank_local,mpi_comm_local,nwgts)
1670              if (local_timers_on >= 1) call oasis_timer_stop('cpl_setup_n4rd_ceg')
1671           else
1672              if (local_timers_on >= 1) then
1673                 call oasis_timer_start('cpl_setup_n4rd_orig_barrier')
1674                 call oasis_mpi_barrier(mpi_comm_local, 'cpl_setup_n4rd_orig')
1675                 call oasis_timer_stop('cpl_setup_n4rd_orig_barrier')
1676                 call oasis_timer_start('cpl_setup_n4rd_orig')
1677              endif
1678              call oasis_map_smatreaddnc_orig(sMati,prism_part(spart)%gsmap,prism_part(dpart)%gsmap, &
1679                 trim(mapopt),trim(prism_mapper(mapID)%file),mpi_rank_local,mpi_comm_local,nwgts)
1680              if (local_timers_on >= 1) call oasis_timer_stop('cpl_setup_n4rd_orig')
1681           endif
1682
1683           if (OASIS_Debug >= 15) then
1684              write(nulprt,*) subname," DEBUG part_create part1 gsize",prism_part(part1)%gsize
1685              do r1 = 1,min(prism_part(part1)%gsmap%ngseg,10)
1686                 write(nulprt,*) subname," DEBUG part_create part1 info ",&
1687                                 prism_part(part1)%gsmap%start(r1),prism_part(part1)%gsmap%length(r1),&
1688                                 prism_part(part1)%gsmap%pe_loc(r1)
1689              enddo
1690
1691              write(nulprt,*) subname," DEBUG part_create part2 gsize",prism_part(part2)%gsize
1692              do r1 = 1,min(prism_part(part2)%gsmap%ngseg,10)
1693                 write(nulprt,*) subname," DEBUG part_create part2 info ",prism_part(part2)%gsmap%start(r1),&
1694                                 prism_part(part2)%gsmap%length(r1),prism_part(part2)%gsmap%pe_loc(r1)
1695              enddo
1696           endif
1697
1698           if (local_timers_on >= 1) then
1699              call oasis_timer_start('cpl_setup_n4sminit_barrier')
1700              call oasis_mpi_barrier(mpi_comm_local, 'cpl_setup_n4sminit')
1701              call oasis_timer_stop('cpl_setup_n4sminit_barrier')
1702              call oasis_timer_start('cpl_setup_n4sminit')
1703           endif
1704
1705           prism_mapper(mapID)%nwgts = nwgts
1706           allocate(prism_mapper(mapID)%sMatP(nwgts))
1707           do n = 1,nwgts
1708              call mct_sMatP_Init(prism_mapper(mapID)%sMatP(n), sMati(n), &
1709                 prism_part(spart)%gsmap, prism_part(dpart)%gsmap, 0, mpi_comm_local, compid)
1710              call mct_sMat_Clean(sMati(n))
1711           enddo
1712           deallocate(sMati)
1713           if (local_timers_on >= 1) call oasis_timer_stop('cpl_setup_n4sminit')
1714
1715           if (OASIS_debug >= 2) then
1716              if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4smprint')
1717              write(cstring,'(a1,i4.4,a1)') 'm',mapID,'-'
1718              ! call mct_rearr_print(prism_mapper(mapID)%sMatP(1)%xtoxprime,mpi_comm_local,nulprt,trim(cstring)//'smpx')
1719              ! call mct_rearr_print(prism_mapper(mapID)%sMatP(1)%yprimetoy,mpi_comm_local,nulprt,trim(cstring)//'smpy')
1720              write(nulprt,*) subname,subname," mct_rearr_print ",trim(cstring)," smpx:"
1721              call mct_rearr_print(prism_mapper(mapID)%sMatP(1)%xtoxprime,mpi_comm_local,nulprt)
1722              write(nulprt,*) subname,subname," mct_rearr_print ",trim(cstring)," smpy:"
1723              call mct_rearr_print(prism_mapper(mapID)%sMatP(1)%yprimetoy,mpi_comm_local,nulprt)
1724              if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4smprint')
1725           endif
1726
1727           if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4e')
1728
1729           lsize = mct_smat_gNumEl(prism_mapper(mapID)%sMatP(1)%Matrix,mpi_comm_local)
1730           prism_mapper(mapID)%init = .true.
1731           if (OASIS_debug >= 15) then
1732              write(nulprt,*) subname," DEBUG ci:done initializing prism_mapper",mapID,&
1733                              " nElements = ",lsize," nwgts = ",nwgts
1734              call oasis_flush(nulprt)
1735           endif
1736           if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4e')
1737        endif  ! map init
1738
1739        !--------------------------------
1740        !>   * Initialize avect1m, the data in avect1 mapped to another grid
1741        !--------------------------------
1742
1743        if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4f')
1744        lsize = mct_gsmap_lsize(prism_part(part2)%gsmap,mpi_comm_local)
1745        if (OASIS_debug >= 15) then
1746           write(nulprt,'(1x,2a,4i12)') subname,' DEBUG ci:part2 info ',part2,mapID,gsize,lsize
1747           write(nulprt,'(1x,2a,4i12)') subname,' DEBUG ci:part2a',prism_part(part2)%gsmap%ngseg,&
1748                                        prism_part(part2)%gsmap%gsize
1749           do n1 = 1,min(prism_part(part2)%gsmap%ngseg,10)
1750              write(nulprt,'(1x,2a,4i12)') subname,' DEBUG ci:part2b',n1,prism_part(part2)%gsmap%start(n1),&
1751                                           prism_part(part2)%gsmap%length(n1),prism_part(part2)%gsmap%pe_loc(n1)
1752           enddo
1753           call oasis_flush(nulprt)
1754        endif
1755        call mct_avect_init(pcpointer%avect1m,rList=trim(pcpointer%fldlist),lsize=lsize)
1756        call mct_avect_zero(pcpointer%avect1m)
1757        if (OASIS_debug >= 15) then
1758           write(nulprt,*) subname,' DEBUG ci:avect1m initialized '
1759           call oasis_flush(nulprt)
1760        endif
1761
1762        !--------------------------------
1763        ! router partition is always other part
1764        !--------------------------------
1765
1766        pcpointer%rpartID = part2
1767        if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4f')
1768     else
1769
1770        !--------------------------------
1771        ! router partition is just coupler part
1772        ! Set to this by default above
1773        !--------------------------------
1774
1775        ! pcpointer%rpartID = pcpointer%partID
1776
1777     endif  ! no mapper
1778
1779   endif ! endif of pcpointer%valid
1780   
1781!   print'(I3,A,X,L,X,8(I8,X))',mpi_rank_global, 'would have done sndrcv here', pcpointer%sndrcv,pcpointer%comp,compid, &
1782!           pcpointer%valid, mapID, pcpointer%rPartID, pcpointer%routerID
1783
1784  enddo   ! npc
1785  enddo   ! nc
1786
1787!-------------------------------------------------
1788! CEG split 1 loop into 2 to allow map reading on different models in parallel.
1789!-------------------------------------------------
1790
1791  if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4g')
1792  if (local_timers_on >= 1) then
1793     call oasis_timer_start('cpl_setup_n4_rt_barrier')
1794     call oasis_mpi_barrier(mpi_comm_local, 'cpl_setup_n4_rt')
1795     call oasis_timer_stop('cpl_setup_n4_rt_barrier')
1796  endif
1797
1798  do nc = 1, prism_mcoupler  ! nc
1799  do npc=1,2
1800
1801     if (npc == 1) then
1802        pcpointer => prism_coupler_put(nc)
1803        pcpntpair => prism_coupler_get(nc)
1804     endif
1805     if (npc == 2) then
1806        pcpointer => prism_coupler_get(nc)
1807        pcpntpair => prism_coupler_put(nc)
1808     endif
1809
1810     namID = pcpointer%namID
1811     part1 = pcpointer%partID
1812     mapID = pcpointer%mapperID
1813
1814!   print'(I3,A,X,L,X,8(I8,X))',mpi_rank_global, '..finally doing sndrcv here', pcpointer%sndrcv, pcpointer%comp, compid, &
1815!           pcpointer%valid, mapID, pcpointer%rPartID, pcpointer%routerID
1816!     if (mapID > 0) then
1817
1818     !--------------------------------
1819     !>   * Initialize router based on rpartID
1820     !--------------------------------
1821
1822     if (local_timers_on >= 2) call oasis_timer_start('cpl_setup_n4_rt')
1823     if (pcpointer%sndrcv) then
1824
1825        if (OASIS_debug >= 15) then
1826           write(nulprt,*) subname,' DEBUG ci:initialize router ',pcpointer%routerID,&
1827                           pcpointer%comp,pcpointer%rpartID
1828           call oasis_flush(nulprt)
1829        endif
1830
1831        if (compid == pcpointer%comp) then
1832           ! routers for sending to self
1833           ! setup router on second pass so rpartID is defined on both sides
1834           ! setup both routers at the same time
1835           if (npc == 2) then
1836              if (OASIS_debug >= 15) then
1837                 write(nulprt,*) subname,' DEBUG self router between part ',pcpointer%rpartID,' and part ',pcpntpair%rpartID, &
1838                    ' with router ',pcpointer%routerID,' and router ',pcpntpair%routerID,' for compid ',compid
1839              endif
1840              if (local_timers_on >= 1) call oasis_timer_start('cpl_setup_n4_rta')
1841              call mct_router_init(prism_part(pcpointer%rpartID)%gsmap, prism_part(pcpntpair%rpartID)%gsmap, &
1842                 mpi_comm_local, prism_router(pcpointer%routerID)%router)
1843              call mct_router_init(prism_part(pcpntpair%rpartID)%gsmap, prism_part(pcpointer%rpartID)%gsmap, &
1844                 mpi_comm_local, prism_router(pcpntpair%routerID)%router)
1845              if (local_timers_on >= 1) call oasis_timer_stop('cpl_setup_n4_rta')
1846
1847              if (OASIS_debug >= 15) then
1848                 write(nulprt,*) subname," DEBUG ci:done initializing prism_router",&
1849                                 pcpointer%routerID
1850                 if (OASIS_debug >= 20) then
1851                    do r1 = 1,min(prism_part(pcpointer%rpartID)%gsmap%ngseg,10)
1852                       write(nulprt,*) subname," DEBUG router gs1 info ",prism_part(pcpointer%rpartID)%gsmap%start(r1),&
1853                                       prism_part(pcpointer%rpartID)%gsmap%length(r1),prism_part(pcpointer%rpartID)%gsmap%pe_loc(r1)
1854                    enddo
1855                    do r1 = 1,min(prism_part(pcpointer%partID)%gsmap%ngseg,10)
1856                       write(nulprt,*) subname," DEBUG router gs2 info ",prism_part(pcpointer%partID)%gsmap%start(r1),&
1857                                       prism_part(pcpointer%partID)%gsmap%length(r1),prism_part(pcpointer%partID)%gsmap%pe_loc(r1)
1858                    enddo
1859                    do r1 = 1,min(prism_part(pcpntpair%rpartID)%gsmap%ngseg,10)
1860                       write(nulprt,*) subname," DEBUG router gs3 info ",prism_part(pcpntpair%rpartID)%gsmap%start(r1),&
1861                                       prism_part(pcpntpair%rpartID)%gsmap%length(r1),prism_part(pcpntpair%rpartID)%gsmap%pe_loc(r1)
1862                    enddo
1863                    do r1 = 1,min(prism_part(pcpntpair%partid)%gsmap%ngseg,10)
1864                       write(nulprt,*) subname," DEBUG router gs4 info ",prism_part(pcpntpair%partid)%gsmap%start(r1),&
1865                                       prism_part(pcpntpair%partid)%gsmap%length(r1),prism_part(pcpntpair%partid)%gsmap%pe_loc(r1)
1866                    enddo
1867                    do r1 = 1,prism_router(pcpointer%routerID)%router%nprocs
1868                       write(nulprt,*) subname," DEBUG router info ",pcpointer%routerID,r1, &
1869                          prism_router(pcpointer%routerID)%router%pe_list(r1),prism_router(pcpointer%routerID)%router%locsize(r1)
1870                    enddo
1871                 endif
1872                 call oasis_flush(nulprt)
1873              endif
1874
1875              if (OASIS_debug >= 15) then
1876                 write(nulprt,*) subname," DEBUG ci:done initializing prism_router",&
1877                                 pcpntpair%routerID
1878                 if (OASIS_debug >= 20) then
1879                    do r1 = 1,prism_router(pcpntpair%routerID)%router%nprocs
1880                       write(nulprt,*) subname," DEBUG router info ",pcpntpair%routerID,r1, &
1881                          prism_router(pcpntpair%routerID)%router%pe_list(r1),prism_router(pcpntpair%routerID)%router%locsize(r1)
1882                    enddo
1883                 endif
1884                 call oasis_flush(nulprt)
1885              endif
1886           endif
1887
1888        else
1889
1890           if (local_timers_on >= 1) call oasis_timer_start('cpl_setup_n4_rtb')
1891           call mct_router_init(pcpointer%comp,prism_part(pcpointer%rpartID)%gsmap, &
1892              mpi_comm_local,prism_router(pcpointer%routerID)%router)
1893           if (local_timers_on >= 1) call oasis_timer_stop('cpl_setup_n4_rtb')
1894
1895           if (OASIS_debug >= 15) then
1896              write(nulprt,*) subname," DEBUG ci:done initializing prism_router",&
1897                              pcpointer%routerID
1898              if (OASIS_debug >= 20) then
1899                 do r1 = 1,prism_router(pcpointer%routerID)%router%nprocs
1900                    write(nulprt,*) subname," DEBUG router info ",pcpointer%routerID,r1, &
1901                       prism_router(pcpointer%routerID)%router%pe_list(r1),prism_router(pcpointer%routerID)%router%locsize(r1)
1902                 enddo
1903              endif
1904              call oasis_flush(nulprt)
1905           endif
1906
1907        endif
1908
1909     endif
1910     if (local_timers_on >= 2) call oasis_timer_stop('cpl_setup_n4_rt')
1911
1912  enddo   ! npc
1913  enddo   ! prism_mcoupler
1914  if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4g')
1915
1916  !--------------------------------
1917  !>   * Find blasnew extra fields
1918  !--------------------------------
1919
1920  if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4h')
1921  do nc = 1,prism_mcoupler
1922     pcpointer => prism_coupler_get(nc)
1923     if (pcpointer%valid) then
1924        do nf = 1,pcpointer%rfno
1925           if (OASIS_debug >= 15) &
1926              write(nulprt,*) subname," DEBUG blasnew search for ",trim(pcpointer%rfna(nf))
1927           found = .false.
1928           do nc2 = 1,prism_mcoupler
1929              if (prism_coupler_get(nc2)%valid .and. &
1930                  prism_coupler_get(nc2)%getput == OASIS3_GET) then
1931                 if (OASIS_debug >= 15) &
1932                    write(nulprt,*) subname," DEBUG blasnew in ",nc2,trim(prism_coupler_get(nc2)%fldlist)
1933                 nf2 = mct_avect_indexRA(prism_coupler_get(nc2)%avect1,trim(pcpointer%rfna(nf)),perrWith='quiet')
1934                 if (nf2 > 0) then
1935                    if (found) then
1936                       write(nulprt,*) subname,estr,'multiple fields found for ',trim(pcpointer%rfna(nf))
1937                       call oasis_abort(file=__FILE__,line=__LINE__)
1938                    endif
1939                    pcpointer%rfcpl(nf)  = nc2
1940                    pcpointer%rffnum(nf) = nf2
1941                    if (.not. prism_coupler_get(nc2)%blasfld) then
1942                       prism_coupler_get(nc2)%blasfld = .true.
1943                       lsize = mct_aVect_lsize(prism_coupler_get(nc2)%aVect1)
1944                       call mct_aVect_init(prism_coupler_get(nc2)%aVect1p,prism_coupler_get(nc2)%aVect1,lsize)
1945                       call mct_aVect_zero(prism_coupler_get(nc2)%aVect1p)
1946                       if (lsize /= mct_aVect_lsize(pcpointer%aVect1)) then
1947                          write(nulprt,*) subname,estr,'error in blasnew fld size for ',trim(pcpointer%rfna(nf))
1948                          write(nulprt,*) subname,estr,'error in blasnew fld size ',lsize,mct_aVect_lsize(pcpointer%aVect1)
1949                          call oasis_abort(file=__FILE__,line=__LINE__)
1950                       endif
1951                    endif
1952                    found = .true.
1953                 endif
1954              endif
1955           enddo
1956           if (.not.found) then
1957              write(nulprt,*) subname,estr,'blasnew field not found for ',trim(pcpointer%rfna(nf))
1958              call oasis_abort(file=__FILE__,line=__LINE__)
1959           endif
1960        enddo
1961     endif
1962  enddo
1963
1964  if (local_timers_on >= 3) call oasis_timer_stop('cpl_setup_n4h')
1965  if (local_timers_on >= 3) call oasis_timer_start('cpl_setup_n4i')
1966  !----------------------------------------------------------
1967  !> * Diagnostics for all couplers
1968  !----------------------------------------------------------
1969
1970  if (OASIS_debug >= 2) then
1971     write(nulprt,*) ' '
1972     write(nulprt,*) subname,' couplers initialized'
1973     do nc = 1,prism_mcoupler
1974        if (prism_coupler_put(nc)%valid) call oasis_coupler_print(nc,prism_coupler_put(nc))
1975        if (prism_coupler_get(nc)%valid) call oasis_coupler_print(nc,prism_coupler_get(nc))
1976     enddo
1977     write(nulprt,*) ' '
1978     CALL oasis_flush(nulprt)
1979  endif
1980
1981  IF (ET_debug) THEN
1982
1983     ! How much event to measure should we expect
1984     ! number of get/put
1985     nb_field = COUNT(prism_coupler_put(1:prism_mcoupler)%valid) + &
1986                COUNT(prism_coupler_get(1:prism_mcoupler)%valid)
1987
1988     CALL oasis_lb_allocate(nb_field)
1989
1990     DO nc = 1, prism_mcoupler
1991        IF (prism_coupler_put(nc)%valid) THEN
1992           !
1993           IF ( prism_coupler_put(nc)%comp == compid .AND. &
1994                prism_part((prism_coupler_put(nc)%partid))%lsize == 0 .AND. &
1995                mpi_rank_local == 0 ) THEN
1996              WRITE(nulprt,*) subname, ' WARNING: component self exchange, not involving master process. Load balancing analysis impossible '
1997              CALL oasis_lb_stop
1998           ENDIF
1999           !
2000           nb_cpl_ts = prism_coupler_put(nc)%maxtime/prism_coupler_put(nc)%dt
2001           IF (OASIS_debug >= 2) &
2002              WRITE(nulprt,'(A11,I2,A16)') ' LB: Define ', nb_cpl_ts, ' coupling events'
2003
2004           imain_kind_lb = LB_PUT
2005           IF ( .NOT. prism_coupler_put(nc)%sndrcv ) imain_kind_lb = LB_OUT
2006
2007           lmap = .FALSE.; lout = .FALSE.; lrst = .FALSE. ; ltrn = .FALSE.
2008           IF ( prism_coupler_put(nc)%mapperID > 0 ) lmap = .TRUE.
2009           IF ( prism_coupler_put(nc)%output ) lout = .TRUE.
2010           IF ( prism_coupler_put(nc)%writrest .OR. prism_coupler_put(nc)%lag > 0 ) &
2011              lrst = .TRUE.
2012           IF ( prism_coupler_put(nc)%writrest .OR. prism_coupler_put(nc)%trans /= ip_instant ) &
2013              ltrn = .TRUE.
2014
2015           CALL oasis_lb_define(nc, imain_kind_lb, prism_coupler_put(nc)%namID, &
2016                                prism_coupler_put(nc)%comp, prism_coupler_put(nc)%ctime, nb_cpl_ts, &
2017                                lmap = lmap, lout = lout, lrst = lrst, ltrn = ltrn )
2018        ENDIF
2019        IF (prism_coupler_get(nc)%valid) THEN
2020           !
2021           IF ( prism_coupler_get(nc)%comp == compid .AND. &
2022                prism_part((prism_coupler_get(nc)%partid))%lsize == 0 .AND. &
2023                mpi_rank_local == 0 ) THEN
2024              WRITE(nulprt,*) subname, ' WARNING: component self exchange, not involving master process. Load balancing analysis impossible '
2025              CALL oasis_lb_stop
2026           ENDIF
2027           !
2028           nb_cpl_ts = prism_coupler_get(nc)%maxtime/prism_coupler_get(nc)%dt
2029           IF (OASIS_debug >= 2) &
2030              WRITE(nulprt,'(A11,I2,A16)') ' LB: Define ', nb_cpl_ts, ' coupling events'
2031
2032           imain_kind_lb = LB_GET
2033           IF ( .NOT. prism_coupler_get(nc)%sndrcv ) imain_kind_lb = LB_READ
2034
2035           lmap = .FALSE.; lout = .FALSE.
2036           IF ( prism_coupler_get(nc)%mapperID > 0 ) lmap = .TRUE.
2037           IF ( prism_coupler_get(nc)%output ) lout = .TRUE.
2038
2039           CALL oasis_lb_define(nc, imain_kind_lb, prism_coupler_get(nc)%namID, &
2040                                prism_coupler_get(nc)%comp, prism_coupler_get(nc)%ctime, nb_cpl_ts, &
2041                                lmap = lmap, lout = lout )
2042
2043        ENDIF
2044     ENDDO
2045
2046  ENDIF
2047
2048  if (local_timers_on >= 3) call oasis_timer_stop ('cpl_setup_n4i')
2049  if (local_timers_on >= 1) call oasis_timer_stop ('cpl_setup_n4')
2050  IF (local_timers_on >= 1) call oasis_timer_stop('cpl_setup')
2051
2052  call oasis_debug_exit(subname)
2053
2054  END SUBROUTINE oasis_coupler_setup
2055
2056!------------------------------------------------------------
2057
2058!> Print routine for oasis_couplers
2059
2060  SUBROUTINE oasis_coupler_print(cplid,pcprint)
2061
2062  IMPLICIT NONE
2063
2064  integer(ip_i4_p),         intent(in) :: cplid   !< coupler id
2065  type(prism_coupler_type), intent(in) :: pcprint !< specific prism_coupler
2066  !----------------------------------------------------------
2067  integer(ip_i4_p) :: mapid, rouid, parid, namid, nflds, rpard
2068  integer(ip_i4_p) :: spart, dpart, nn
2069  character(len=*),parameter :: subname = '(oasis_coupler_print)'
2070
2071  call oasis_debug_enter(subname)
2072
2073  mapid = pcprint%mapperid
2074  rouid = pcprint%routerid
2075  parid = pcprint%partid
2076  rpard = pcprint%rpartid
2077  namid = pcprint%namid
2078  nflds = pcprint%nflds
2079
2080     write(nulprt,*) ' '
2081     write(nulprt,*) subname,' model and cplid',compid,cplid
2082  if (pcprint%getput == OASIS3_PUT) then
2083     write(nulprt,*) subname,'   timerid send     ',cplid,trim(pcprint%fldlist)
2084     write(nulprt,*) subname,'   send fields      ',trim(pcprint%fldlist)
2085     write(nulprt,*) subname,'   from model       ',compid
2086     write(nulprt,*) subname,'   to model         ',pcprint%comp
2087     write(nulprt,*) subname,'   using router     ',rouid
2088     write(nulprt,*) subname,'   transform        ',pcprint%trans
2089     write(nulprt,*) subname,'   snd diagnose     ',pcprint%snddiag
2090     write(nulprt,*) subname,'   snd fld mult     ',pcprint%sndmult
2091     write(nulprt,*) subname,'   snd fld add      ',pcprint%sndadd
2092  endif
2093  if (pcprint%getput == OASIS3_GET) then
2094     write(nulprt,*) subname,'   timerid recv     ',cplid,trim(pcprint%fldlist)
2095     write(nulprt,*) subname,'   recv fields      ',trim(pcprint%fldlist)
2096     write(nulprt,*) subname,'   from model       ',pcprint%comp
2097     write(nulprt,*) subname,'   to model         ',compid
2098     write(nulprt,*) subname,'   using router     ',rouid
2099     write(nulprt,*) subname,'   rcv diagnose     ',pcprint%rcvdiag
2100     write(nulprt,*) subname,'   rcv fld mult     ',pcprint%rcvmult
2101     write(nulprt,*) subname,'   rcv fld add      ',pcprint%rcvadd
2102     write(nulprt,*) subname,'   blas fields      ',pcprint%blasfld
2103     write(nulprt,*) subname,'   blas fld number  ',pcprint%rfno
2104     do nn = 1,pcprint%rfno
2105       write(nulprt,*) subname,'   blas fld name    ',nn,trim(pcprint%rfna(nn))
2106       write(nulprt,*) subname,'   blas fld mult    ',nn,pcprint%rfmult(nn)
2107       write(nulprt,*) subname,'   blas fld add     ',nn,pcprint%rfadd(nn)
2108       write(nulprt,*) subname,'   blas fld cplid   ',nn,pcprint%rfcpl(nn)
2109       write(nulprt,*) subname,'   blas fld fnum    ',nn,pcprint%rffnum(nn)
2110     enddo
2111  endif
2112     write(nulprt,*) subname,'   namcouple op     ',pcprint%ops
2113     write(nulprt,*) subname,'   valid            ',pcprint%valid
2114     write(nulprt,*) subname,'   namcouple id     ',namid
2115     write(nulprt,*) subname,'   variable ids     ',pcprint%varid(1:nflds)
2116     write(nulprt,*) subname,'   sndrcv flag      ',pcprint%sndrcv
2117     write(nulprt,*) subname,'   output flag      ',pcprint%output
2118     write(nulprt,*) subname,'   input flag       ',pcprint%input
2119     write(nulprt,*) subname,'   input file       ',trim(pcprint%inpfile)
2120     write(nulprt,*) subname,'   restart file     ',trim(pcprint%rstfile)
2121     write(nulprt,*) subname,'   tag              ',pcprint%tag
2122     write(nulprt,*) subname,'   seq              ',pcprint%seq
2123     write(nulprt,*) subname,'   maxtime          ',pcprint%maxtime
2124     write(nulprt,*) subname,'   dt, lag          ',pcprint%dt,pcprint%lag
2125!     write(nulprt,*) subname,'   partid, size ',parid,trim(prism_part(parid)%gridname),&
2126!                                               prism_part(parid)%gsize
2127     write(nulprt,*) subname,'   partid, gsize    ',parid,prism_part(parid)%gsize
2128     write(nulprt,*) subname,'   partid, lsize    ',parid,prism_part(parid)%lsize
2129     write(nulprt,*) subname,'   partid, nx,ny    ',prism_part(parid)%nx,prism_part(parid)%ny
2130!     write(nulprt,*) subname,'   rpartid,size ',rpard,trim(prism_part(rpard)%gridname),&
2131!                                                prism_part(rpard)%gsize
2132     write(nulprt,*) subname,'   rpartid,gsize    ',rpard,prism_part(rpard)%gsize
2133     write(nulprt,*) subname,'   rpartid,nx,ny    ',prism_part(rpard)%nx,prism_part(rpard)%ny
2134     write(nulprt,*) subname,'   maploc           ',trim(pcprint%maploc)
2135
2136  if (mapid > 0) then
2137     WRITE(nulprt,*) subname,'   src grid        :',trim(prism_mapper(mapid)%srcgrid)
2138     WRITE(nulprt,*) subname,'   dst grid        :',trim(prism_mapper(mapid)%dstgrid)
2139     write(nulprt,*) subname,'   use map          ',mapid,trim(prism_mapper(mapid)%file)
2140     write(nulprt,*) subname,'   nwgts            ',mapid,prism_mapper(mapid)%nwgts
2141     spart = prism_mapper(mapid)%spart
2142     dpart = prism_mapper(mapid)%dpart
2143     write(nulprt,*) subname,'   conserve         ',pcprint%conserv
2144     write(nulprt,*) subname,'   conserve opt     ',pcprint%consopt
2145     write(nulprt,*) subname,'   location         ',trim(prism_mapper(mapid)%loc)
2146     write(nulprt,*) subname,'   opt,optval       ',trim(prism_mapper(mapid)%opt),' ',&
2147                                                    trim(prism_mapper(mapid)%optval)
2148     write(nulprt,*) subname,'   s/d partids      ',spart,dpart
2149     if (spart > 0) &
2150     write(nulprt,*) subname,'   from/to partition',trim(prism_part(spart)%gridname),' ',&
2151                                                    trim(prism_part(dpart)%gridname)
2152     write(nulprt,*) subname,'   from nx,ny       ',prism_part(spart)%gsize,prism_part(spart)%nx,&
2153                                                    prism_part(spart)%ny
2154     if (dpart > 0) &
2155     write(nulprt,*) subname,'   to nx,ny         ',prism_part(dpart)%gsize, prism_part(dpart)%nx,&
2156                                                    prism_part(dpart)%ny
2157  endif  ! mapid > 0
2158
2159  call oasis_flush(nulprt)
2160
2161  call oasis_debug_exit(subname)
2162
2163  END SUBROUTINE oasis_coupler_print
2164
2165!------------------------------------------------------------
2166! !BOP ===========================================================================
2167!
2168!> Build a consistent variable name based on bundles
2169!
2170! !DESCRIPTION:
2171!     Build a variable name for a given variable based on the name, number
2172!     of bundled fields, and bundle level.  Needs to be used in a few different
2173!     places in oasis.
2174!
2175! !INTERFACE:  -----------------------------------------------------------------
2176
2177  SUBROUTINE oasis_coupler_bldvarname(varid, varnum, vname)
2178
2179  IMPLICIT NONE
2180
2181  integer(ip_i4_p), intent(in)  :: varid   !< variable id
2182  integer(ip_i4_p), intent(in)  :: varnum  !< variable bundle level number
2183  character(len=*), intent(out) :: vname   !< variable name
2184  !----------------------------------------------------------
2185  character(len=*),parameter :: subname = '(oasis_coupler_bldvarname)'
2186
2187  call oasis_debug_enter(subname)
2188
2189  if (varnum > prism_var(varid)%num) then
2190     write(nulprt,*) subname,estr,'invalid varnum varid = ',varid,trim(prism_var(varid)%name)
2191     write(nulprt,*) subname,estr,'invalid varnum = ',varnum,prism_var(varid)%num
2192     call oasis_abort(file=__FILE__,line=__LINE__)
2193  endif
2194
2195  if (prism_var(varid)%num > 1) then
2196     write(vname,'(a,i3.3)') trim(prism_var(varid)%name)//'.',varnum
2197  else
2198     vname = trim(prism_var(varid)%name)
2199  endif
2200
2201  if (OASIS_debug >= 20) then
2202     write(nulprt,*) subname,' check vname ',varnum,trim(vname)
2203  endif
2204
2205  call oasis_debug_exit(subname)
2206
2207  END SUBROUTINE oasis_coupler_bldvarname
2208
2209!------------------------------------------------------------
2210! !BOP ===========================================================================
2211!
2212!> Deconstruct the varname based on oasis_coupler_bldvarname
2213!
2214! !DESCRIPTION:
2215!     Deconstruct a variable name for a given variable based on the name, number
2216!     of bundled fields, and bundle level.  Must be consistent with oasis_coupler_bldvarname
2217!
2218! !INTERFACE:  -----------------------------------------------------------------
2219
2220  SUBROUTINE oasis_coupler_unbldvarname(varid, vname, varnum)
2221
2222  IMPLICIT NONE
2223
2224  integer(ip_i4_p), intent(in)  :: varid   !< variable id
2225  character(len=*), intent(in)  :: vname   !< variable name
2226  integer(ip_i4_p), intent(out) :: varnum  !< variable bundle level number
2227  !----------------------------------------------------------
2228  integer(ip_i4_p)  :: vlen
2229  character(len=16) :: clen
2230  character(len=*),parameter :: subname = '(oasis_coupler_unbldvarname)'
2231
2232  call oasis_debug_enter(subname)
2233
2234  if (prism_var(varid)%num > 1) then
2235     vlen = len_trim(vname)
2236     clen = vname(vlen-2:vlen)
2237     read(clen,'(i3.3)') varnum
2238     if (OASIS_debug >= 20) then
2239        write(nulprt,*) subname,' check vlen ',vlen,trim(clen)
2240     endif
2241  else
2242     varnum = 1
2243  endif
2244
2245  if (OASIS_debug >= 20) then
2246     write(nulprt,*) subname,' check varnum ',varnum,trim(vname)
2247  endif
2248
2249  if (varnum > prism_var(varid)%num) then
2250     write(nulprt,*) subname,estr,'invalid varnum varid = ',varid,trim(prism_var(varid)%name)
2251     write(nulprt,*) subname,estr,'invalid varnum = ',varnum,prism_var(varid)%num
2252     call oasis_abort(file=__FILE__,line=__LINE__)
2253  endif
2254
2255  call oasis_debug_exit(subname)
2256
2257  END SUBROUTINE oasis_coupler_unbldvarname
2258
2259!------------------------------------------------------------
2260!------------------------------------------------------------
2261! !BOP ===========================================================================
2262!
2263!> Search a character field list for a matching values
2264!
2265! !DESCRIPTION:
2266!     Sort a character array and the associated array(s) based on a
2267!     reasonably fast sort algorithm
2268!
2269! !INTERFACE:  -----------------------------------------------------------------
2270
2271subroutine cplfind(num, fldlist, fld, ifind, nfind)
2272
2273! !USES:
2274
2275   !--- local kinds ---
2276   integer,parameter :: R8 = ip_double_p
2277   integer,parameter :: IN = ip_i4_p
2278   integer,parameter :: CL = ic_lvar
2279
2280! !INPUT/OUTPUT PARAMETERS:
2281
2282   integer(IN),      intent(in)  :: num         !< size of array
2283   character(len=CL),intent(in)  :: fldlist(:)  !< sorted field list
2284   character(len=CL),intent(in)  :: fld         !< field to search for
2285   integer(IN)      ,intent(out) :: ifind       !< first match index
2286   integer(IN)      ,intent(out) :: nfind       !< number that match
2287
2288! !EOP
2289
2290   !--- local ---
2291   integer(IN)    :: is,ie,im
2292   logical        :: found,check
2293
2294   !--- formats ---
2295   character(*),parameter :: subName = '(cplfind) '
2296
2297!-------------------------------------------------------------------------------
2298!
2299!-------------------------------------------------------------------------------
2300
2301!   call oasis_debug_enter(subname)
2302
2303    ifind = 0
2304    nfind = 0
2305
2306    is = 1
2307    ie = num
2308    found = .false.
2309
2310    ! check endpoints first, the binary search uses integer
2311    ! math which makes hitting the endpoints more difficult
2312    ! so check manually.  also if list size is 1, need to do this.
2313
2314    if (.not.found) then
2315       im = 1
2316       if (fld == fldlist(im)) found = .true.
2317    endif
2318    if (.not.found) then
2319       im = num
2320       if (fld == fldlist(im)) found = .true.
2321    endif
2322
2323    ! do a binary search
2324
2325    do while (.not.found .and. ie > is)
2326       im = (is + ie) / 2
2327       im = max(im,is)
2328       im = min(im,ie)
2329!       write(nulprt,*) subname,'tcx',is,ie,im,trim(fld),' ',trim(fldlist(im))
2330       if (fld == fldlist(im)) then
2331          found = .true.
2332       elseif (fld > fldlist(im)) then
2333          is = max(im,is+1)
2334       else
2335          ie = min(im,ie-1)
2336       endif
2337    enddo
2338
2339    ! if a match was found, find first and last instance of match in list
2340
2341    if (found) then
2342       is = im
2343       ie = im
2344       if (is > 1) then
2345          check = .true.
2346          do while (check)
2347             if (is > 1) then
2348                if (fld /= fldlist(is-1)) then
2349                   check = .false.
2350                else
2351                   is = is - 1
2352                endif
2353             else
2354                check = .false.
2355             endif
2356          enddo
2357       endif
2358       if (ie < num) then
2359          check = .true.
2360          do while (check)
2361             if (ie < num) then
2362                if (fld /= fldlist(ie+1)) then
2363                   check = .false.
2364                else
2365                   ie = ie + 1
2366                endif
2367             else
2368                check = .false.
2369             endif
2370          enddo
2371       endif
2372       ifind = is
2373       nfind = (ie - is + 1)   
2374    endif
2375
2376!   call oasis_debug_exit(subname)
2377
2378end subroutine cplfind
2379
2380!===============================================================================
2381
2382END MODULE mod_oasis_coupler
2383
2384
Note: See TracBrowser for help on using the repository browser.