ROMS
Loading...
Searching...
No Matches
set_contact_mod Module Reference

Functions/Subroutines

subroutine, public set_contact (ng, model)
 
subroutine, private set_contact_nf90 (ng, model)
 
subroutine, private set_contact_pio (ng, model)
 

Function/Subroutine Documentation

◆ set_contact()

subroutine, public set_contact_mod::set_contact ( integer, intent(in) ng,
integer, intent(in) model )

Definition at line 41 of file set_contact.F.

42!***********************************************************************
43!
44! Imported variable declarations.
45!
46 integer, intent(in) :: ng, model
47!
48! Local variable declarations.
49!
50 integer :: LBi, UBi, LBj, UBj
51!
52 character (len=*), parameter :: MyFile = &
53 & __FILE__
54!
55!-----------------------------------------------------------------------
56! Read in nesting connectivity NetCDF file and set contact nested
57! structures and variables.
58!-----------------------------------------------------------------------
59!
60 SELECT CASE (inp_lib)
61 CASE (io_nf90)
62 CALL set_contact_nf90 (ng, model)
63
64# if defined PIO_LIB && defined DISTRIBUTE
65 CASE (io_pio)
66 CALL set_contact_pio (ng, model)
67# endif
68 CASE DEFAULT
69 IF (master) WRITE (stdout,10) inp_lib
70 exit_flag=2
71 END SELECT
72 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
73!
74 10 FORMAT (' SET_CONTACT - Illegal input file type, io_type = ',i0, &
75 & /,15x,'Check KeyWord ''INP_LIB'' in ''roms.in''.')
76!
77 RETURN

References mod_scalars::exit_flag, strings_mod::founderror(), mod_parallel::master, mod_scalars::noerror, set_contact_nf90(), set_contact_pio(), and mod_iounits::stdout.

Referenced by inp_par_mod::inp_par().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_contact_nf90()

subroutine, private set_contact_mod::set_contact_nf90 ( integer, intent(in) ng,
integer, intent(in) model )
private

Definition at line 81 of file set_contact.F.

82!***********************************************************************
83!
84! Imported variable declarations.
85!
86 integer, intent(in) :: ng, model
87!
88! Local variable declarations.
89!
90 logical, dimension((Ngrids-1)*2) :: Lcoincident
91 logical, dimension((Ngrids-1)*2) :: Lcomposite
92 logical, dimension((Ngrids-1)*2) :: Lmosaic
93 logical, dimension((Ngrids-1)*2) :: Lrefinement
94
95 logical, dimension((Ngrids-1)*2) :: Linterpolate
96!
97 integer :: cr, dg, ibry, ic, ig, ip, m, rg, vindex
98 integer :: my_Ncontact, my_Ngrids, my_nLweights, my_nQweights
99 integer :: NGCid
100
101 integer, dimension(Ngrids) :: my_Lm, my_Mm
102 integer, dimension(Ngrids) :: refine_factor
103
104 integer, dimension((Ngrids-1)*2) :: NpointsR
105 integer, dimension((Ngrids-1)*2) :: NpointsU
106 integer, dimension((Ngrids-1)*2) :: NpointsV
107!
108 real(r8), allocatable :: Lweight(:,:)
109# ifdef QUADRATIC_WEIGHTS
110 real(r8), allocatable :: Qweight(:,:)
111# endif
112 real(r8), allocatable :: Xrg(:), Yrg(:)
113 real(r8), allocatable :: angle(:)
114 real(r8), allocatable :: dmde(:), dndx(:)
115 real(r8), allocatable :: f(:)
116 real(r8), allocatable :: h(:)
117 real(r8), allocatable :: mask(:)
118 real(r8), allocatable :: pm(:), pn(:)
119!
120 character (len=*), parameter :: MyFile = &
121 & __FILE__//", set_contact_nf90"
122!
123!-----------------------------------------------------------------------
124! Read in nesting grids contact point packed information into local
125! arrays.
126!-----------------------------------------------------------------------
127!
128! Open contact points NetCDF file for reading.
129!
130 CALL netcdf_open (ng, model, ngcname, 0, ngcid)
131 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
132 WRITE (stdout,20) trim(ngcname)
133 20 FORMAT (/,' SET_CONTACT_NF90 - unable to open contact points ', &
134 & ' NetCDF file: ',a)
135 RETURN
136 END IF
137!
138! Inquire about the variables.
139!
140 CALL netcdf_inq_var (ng, model, ngcname, &
141 & ncid = ngcid)
142 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
143!
144! Check if the number of nested grids is correct.
145!
146 CALL netcdf_get_dim (ng, model, ngcname, &
147 & ncid = ngcid, &
148 & dimname = 'Ngrids', &
149 & dimsize = my_ngrids)
150 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
151
152 IF (my_ngrids.ne.ngrids) THEN
153 IF (master) THEN
154 WRITE (stdout,10) 'inconsistent parameter, Ngrids = ', &
155 & ngrids, my_ngrids
156 10 FORMAT (/,' SET_CONTACT_NF90 - ', a, i4, 2x, i4, &
157 & /,20x,'in input file:'2x,a)
158 END IF
159 exit_flag=5
160 RETURN
161 END IF
162!
163! Check number of contact regions, Ncontact = (Ngrids-1)*2.
164!
165 CALL netcdf_get_dim (ng, model, ngcname, &
166 & ncid = ngcid, &
167 & dimname = 'Ncontact', &
168 & dimsize = my_ncontact)
169 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
170
171 IF (my_ncontact.ne.(ngrids-1)*2) THEN
172 IF (master) THEN
173 WRITE (stdout,10) 'inconsistent parameter, Ncontact = ', &
174 & (ngrids-1)*2, my_ncontact
175 END IF
176 exit_flag=5
177 RETURN
178 END IF
179 ncontact=(ngrids-1)*2
180!
181! Check number of contact points for bilinear interpolation weights.
182!
183 CALL netcdf_get_dim (ng, model, ngcname, &
184 & ncid = ngcid, &
185 & dimname = 'nLweights', &
186 & dimsize = my_nlweights)
187 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
188
189 IF (my_nlweights.ne.4) THEN
190 IF (master) THEN
191 WRITE (stdout,10) 'inconsistent parameter, nLweights = ', &
192 & 4, my_nlweights
193 END IF
194 exit_flag=5
195 RETURN
196 END IF
197
198# ifdef QUADRATIC_WEIGHTS
199!
200! Check number of contact points for quadratic interpolation weights.
201!
202 CALL netcdf_get_dim (ng, model, ngcname, &
203 & ncid = ngcid, &
204 & dimname = 'nQweights', &
205 & dimsize = my_nqweights)
206 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
207
208 IF (my_nqweights.ne.9) THEN
209 IF (master) THEN
210 WRITE (stdout,10) 'inconsistent parameter, nQweights = ', &
211 & 9, my_nqweights
212 END IF
213 exit_flag=5
214 RETURN
215 END IF
216# endif
217!
218! Get number of contact points in input NetCDF file. It may include
219! or not contact point over land when land/sea masking is activated.
220!
221 CALL netcdf_get_dim (ng, model, ngcname, &
222 & ncid = ngcid, &
223 & dimname = 'datum', &
224 & dimsize = ncdatum)
225 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
226!
227! Read in grid dimensions and if the grid order is correct. The order
228! of the grids is important in nesting.
229!
230 CALL netcdf_get_ivar (ng, model, ngcname, &
231 & 'Lm', my_lm, &
232 & ncid = ngcid)
233 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
234
235 CALL netcdf_get_ivar (ng, model, ngcname, &
236 & 'Mm', my_mm, &
237 & ncid = ngcid)
238 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
239
240 DO ig=1,ngrids
241 IF (my_lm(ig).ne.lm(ig)) THEN
242 IF (master) THEN
243 WRITE (stdout,10) 'inconsistent grid order, Lm = ', &
244 & lm(ig), my_lm(ig)
245 END IF
246 exit_flag=5
247 RETURN
248 END IF
249 IF (my_mm(ig).ne.mm(ig)) THEN
250 IF (master) THEN
251 WRITE (stdout,10) 'inconsistent grid order, Mm = ', &
252 & mm(ig), my_mm(ig)
253 END IF
254 exit_flag=5
255 RETURN
256 END IF
257 END DO
258!
259! Read in nesting types logical switches, 1:Ngrids.
260!
261 CALL netcdf_get_lvar (ng, model, ngcname, &
262 & 'coincident', lcoincident, &
263 & ncid = ngcid)
264 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
265
266 CALL netcdf_get_lvar (ng, model, ngcname, &
267 & 'composite', lcomposite, &
268 & ncid = ngcid)
269 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
270
271 CALL netcdf_get_lvar (ng, model, ngcname, &
272 & 'mosaic', lmosaic, &
273 & ncid = ngcid)
274 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
275
276 CALL netcdf_get_lvar (ng, model, ngcname, &
277 & 'refinement', lrefinement, &
278 & ncid = ngcid)
279 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
280!
281! Read in refinement factor from donor grid, 1:Ngrids.
282!
283 CALL netcdf_get_ivar (ng, model, ngcname, &
284 & 'refine_factor', refine_factor, &
285 & ncid = ngcid)
286 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
287!
288! Read in vertical interpolation at contact points switch, 1:Ncontact.
289!
290 CALL netcdf_get_lvar (ng, model, ngcname, &
291 & 'interpolate', linterpolate, &
292 & ncid = ngcid)
293 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
294!
295! Read in data donor and data receiver grid number, 1:Ncontact.
296!
297 IF (.not.allocated(donor_grid)) THEN
298 allocate ( donor_grid((ngrids-1)*2) )
299 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
300 END IF
301 CALL netcdf_get_ivar (ng, model, ngcname, &
302 & 'donor_grid', donor_grid, &
303 & ncid = ngcid)
304 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
305
306 IF (.not.allocated(receiver_grid)) THEN
307 allocate ( receiver_grid((ngrids-1)*2) )
308 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
309 END IF
310 CALL netcdf_get_ivar (ng, model, ngcname, &
311 & 'receiver_grid', receiver_grid, &
312 & ncid = ngcid)
313 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
314!
315! Read in donor grid (I,J) indices at PSI points used to extract
316! refined grid. Values have a fill value of -999 if not applicable.
317!
318 IF (.not.allocated(i_left)) THEN
319 allocate ( i_left(ngrids) )
320 dmem(ng)=dmem(ng)+real(ngrids,r8)
321 END IF
322 CALL netcdf_get_ivar (ng, model, ngcname, &
323 & 'I_left', i_left, &
324 & ncid = ngcid)
325 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
326!
327 IF (.not.allocated(i_right)) THEN
328 allocate ( i_right(ngrids) )
329 dmem(ng)=dmem(ng)+real(ngrids,r8)
330 END IF
331 CALL netcdf_get_ivar (ng, model, ngcname, &
332 & 'I_right', i_right, &
333 & ncid = ngcid)
334 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
335!
336 IF (.not.allocated(j_bottom)) THEN
337 allocate ( j_bottom(ngrids) )
338 dmem(ng)=dmem(ng)+real(ngrids,r8)
339 END IF
340 CALL netcdf_get_ivar (ng, model, ngcname, &
341 & 'J_bottom', j_bottom, &
342 & ncid = ngcid)
343 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
344!
345 IF (.not.allocated(j_top)) THEN
346 allocate ( j_top(ngrids) )
347 dmem(ng)=dmem(ng)+real(ngrids,r8)
348 END IF
349 CALL netcdf_get_ivar (ng, model, ngcname, &
350 & 'J_top', j_top, &
351 & ncid = ngcid)
352 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
353!
354! Read in starting and ending contact points (RHO-, U-, V-) index in
355! packed data vectors for each contact region, 1:Ncontact.
356!
357 IF (.not.allocated(nstrr)) THEN
358 allocate ( nstrr((ngrids-1)*2) )
359 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
360 END IF
361 CALL netcdf_get_ivar (ng, model, ngcname, &
362 & 'NstrR', nstrr, &
363 & ncid = ngcid)
364 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
365!
366 IF (.not.allocated(nendr)) THEN
367 allocate ( nendr((ngrids-1)*2) )
368 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
369 END IF
370 CALL netcdf_get_ivar (ng, model, ngcname, &
371 & 'NendR', nendr, &
372 & ncid = ngcid)
373 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
374!
375 IF (.not.allocated(nstru)) THEN
376 allocate ( nstru((ngrids-1)*2) )
377 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
378 END IF
379 CALL netcdf_get_ivar (ng, model, ngcname, &
380 & 'NstrU', nstru, &
381 & ncid = ngcid)
382 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
383!
384 IF (.not.allocated(nendu)) THEN
385 allocate ( nendu((ngrids-1)*2) )
386 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
387 END IF
388 CALL netcdf_get_ivar (ng, model, ngcname, &
389 & 'NendU', nendu, &
390 & ncid = ngcid)
391 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
392!
393 IF (.not.allocated(nstrv)) THEN
394 allocate ( nstrv((ngrids-1)*2) )
395 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
396 END IF
397 CALL netcdf_get_ivar (ng, model, ngcname, &
398 & 'NstrV', nstrv, &
399 & ncid = ngcid)
400 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
401!
402 IF (.not.allocated(nendv)) THEN
403 allocate ( nendv((ngrids-1)*2) )
404 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
405 END IF
406 CALL netcdf_get_ivar (ng, model, ngcname, &
407 & 'NendV', nendv, &
408 & ncid = ngcid)
409 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
410!
411! Read in contact region number for each contact point, 1:NCdatum.
412!
413 IF (.not.allocated(contact_region)) THEN
414 allocate ( contact_region(ncdatum) )
415 dmem(ng)=dmem(ng)+real(ncdatum,r8)
416 END IF
417 CALL netcdf_get_ivar (ng, model, ngcname, &
418 & 'contact_region', contact_region, &
419 & ncid = ngcid)
420 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
421!
422! Read in flag (1:NCdatum) to determine if the contact point is at
423! receiver grid interior (zero value) or on the boundary (1:western,
424! 2:southern, 3:eastern, 4:northern).
425!
426 IF (.not.allocated(on_boundary)) THEN
427 allocate ( on_boundary(ncdatum) )
428 dmem(ng)=dmem(ng)+real(ncdatum,r8)
429 END IF
430 CALL netcdf_get_ivar (ng, model, ngcname, &
431 & 'on_boundary', on_boundary, &
432 & ncid = ngcid)
433 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
434!
435! Read in I-left and J-bottom indices of donor cell containing contact
436! point, 1:NCdatum.
437!
438 IF (.not.allocated(idg_cp)) THEN
439 allocate ( idg_cp(ncdatum) )
440 dmem(ng)=dmem(ng)+real(ncdatum,r8)
441 END IF
442 CALL netcdf_get_ivar (ng, model, ngcname, &
443 & 'Idg', idg_cp, &
444 & ncid = ngcid)
445 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
446
447 IF (.not.allocated(jdg_cp)) THEN
448 allocate ( jdg_cp(ncdatum) )
449 dmem(ng)=dmem(ng)+real(ncdatum,r8)
450 END IF
451 CALL netcdf_get_ivar (ng, model, ngcname, &
452 & 'Jdg', jdg_cp, &
453 & ncid = ngcid)
454 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
455!
456! Read in I- and J-indices of receiver grid contact point, 1:NCdatum.
457!
458 IF (.not.allocated(irg_cp)) THEN
459 allocate ( irg_cp(ncdatum) )
460 dmem(ng)=dmem(ng)+real(ncdatum,r8)
461 END IF
462 CALL netcdf_get_ivar (ng, model, ngcname, &
463 & 'Irg', irg_cp, &
464 & ncid = ngcid)
465 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
466
467 IF (.not.allocated(jrg_cp)) THEN
468 allocate ( jrg_cp(ncdatum) )
469 dmem(ng)=dmem(ng)+real(ncdatum,r8)
470 END IF
471 CALL netcdf_get_ivar (ng, model, ngcname, &
472 & 'Jrg', jrg_cp, &
473 & ncid = ngcid)
474 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
475!
476! Read in X- and Y-location of receiver grid contact point, 1:NCdatum.
477!
478 IF (.not.allocated(xrg)) THEN
479 allocate ( xrg(ncdatum) )
480 dmem(ng)=dmem(ng)+real(ncdatum,r8)
481 END IF
482 CALL netcdf_get_fvar (ng, model, ngcname, &
483 & 'Xrg', xrg, &
484 & ncid = ngcid)
485 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
486
487 IF (.not.allocated(yrg)) THEN
488 allocate ( yrg(ncdatum) )
489 dmem(ng)=dmem(ng)+real(ncdatum,r8)
490 END IF
491 CALL netcdf_get_fvar (ng, model, ngcname, &
492 & 'Yrg', yrg, &
493 & ncid = ngcid)
494 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
495!
496! Read in contact point horizontal linear interpolation weights,
497! 1:NCdatum.
498!
499 IF (.not.allocated(lweight)) THEN
500 allocate ( lweight(my_nlweights,ncdatum) )
501 dmem(ng)=dmem(ng)+real(my_nlweights*ncdatum,r8)
502 END IF
503 CALL netcdf_get_fvar (ng, model, ngcname, &
504 & 'Lweight', lweight, &
505 & ncid = ngcid)
506 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
507
508# ifdef QUADRATIC_WEIGHTS
509!
510! Read in contact point horizontal quadratic interpolation weights,
511! 1:NCdatum.
512!
513 IF (.not.allocated(qweight)) THEN
514 allocate ( qweight(my_nqweights,ncdatum) )
515 dmem(ng)=dmem(ng)+real(my_nqweights*ncdatum,r8)
516 END IF
517 CALL netcdf_get_fvar (ng, model, ngcname, &
518 & 'Qweight', qweight, &
519 & ncid = ngcid)
520 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
521# endif
522!
523! Read in contact point bathymetry at RHO-points, 1:NCdatum.
524!
525 IF (.not.allocated(h)) THEN
526 allocate ( h(ncdatum) )
527 dmem(ng)=dmem(ng)+real(ncdatum,r8)
528 END IF
529 CALL netcdf_get_fvar (ng, model, ngcname, &
530 & 'h', h, &
531 & ncid = ngcid)
532 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
533!
534! Read in contact point Coriolis parameter at RHO-points, 1:NCdatum.
535!
536 IF (.not.allocated(f)) THEN
537 allocate ( f(ncdatum) )
538 dmem(ng)=dmem(ng)+real(ncdatum,r8)
539 END IF
540 CALL netcdf_get_fvar (ng, model, ngcname, &
541 & 'f', f, &
542 & ncid = ngcid)
543 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
544!
545! Read in contact point curvilinar coordinates in XI- and ETA-direction
546! at RHO-points, 1:NCdatum.
547!
548 IF (.not.allocated(pm)) THEN
549 allocate ( pm(ncdatum) )
550 dmem(ng)=dmem(ng)+real(ncdatum,r8)
551 END IF
552 CALL netcdf_get_fvar (ng, model, ngcname, &
553 & 'pm', pm, &
554 & ncid = ngcid)
555 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
556
557 IF (.not.allocated(pn)) THEN
558 allocate ( pn(ncdatum) )
559 dmem(ng)=dmem(ng)+real(ncdatum,r8)
560 END IF
561 CALL netcdf_get_fvar (ng, model, ngcname, &
562 & 'pn', pn, &
563 & ncid = ngcid)
564 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
565!
566! Read in contact point inverse metric factor in XI- and ETA-direction,
567! 1:NCdatum.
568!
569 IF (.not.allocated(dndx)) THEN
570 allocate ( dndx(ncdatum) )
571 dmem(ng)=dmem(ng)+real(ncdatum,r8)
572 END IF
573 CALL netcdf_get_fvar (ng, model, ngcname, &
574 & 'dndx', dndx, &
575 & ncid = ngcid)
576 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
577
578 IF (.not.allocated(dmde)) THEN
579 allocate ( dmde(ncdatum) )
580 dmem(ng)=dmem(ng)+real(ncdatum,r8)
581 END IF
582 CALL netcdf_get_fvar (ng, model, ngcname, &
583 & 'dmde', dmde, &
584 & ncid = ngcid)
585 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
586!
587! Read in contact point angle between XI-axis and EAST, 1:NCdatum.
588!
589 IF (.not.allocated(angle)) THEN
590 allocate ( angle(ncdatum) )
591 dmem(ng)=dmem(ng)+real(ncdatum,r8)
592 END IF
593 CALL netcdf_get_fvar (ng, model, ngcname, &
594 & 'angle', angle, &
595 & ncid = ngcid)
596 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
597!
598! Read in contact point land/sea mask, 1:NCdatum.
599!
600 IF (.not.allocated(mask)) THEN
601 allocate ( mask(ncdatum) )
602 dmem(ng)=dmem(ng)+real(ncdatum,r8)
603 END IF
604 CALL netcdf_get_fvar (ng, model, ngcname, &
605 & 'mask', mask, &
606 & ncid = ngcid)
607 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
608!
609! Close contact points NetCDF file.
610!
611 CALL netcdf_close (ng, model, ngcid, ngcname, .false.)
612 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
613!
614!-----------------------------------------------------------------------
615! Unpack contact point data into nesting structures (type T_NGC).
616!-----------------------------------------------------------------------
617!
618! Determine number of contact points in each contact region.
619!
620 IF (.not.allocated(ncpoints)) THEN
621 allocate ( ncpoints(ncontact) )
622 dmem(ng)=dmem(ng)+real(ncontact,r8)
623 END IF
624 DO cr=1,ncontact
625 npointsr(cr)=nendr(cr)-nstrr(cr)+1
626 npointsu(cr)=nendu(cr)-nstru(cr)+1
627 npointsv(cr)=nendv(cr)-nstrv(cr)+1
628 ncpoints(cr)=npointsr(cr)+npointsu(cr)+npointsv(cr)
629 END DO
630!
631! Allocate grid connectivity (type T_NGC) structures.
632!
633 allocate ( rcontact(ncontact) )
634 allocate ( ucontact(ncontact) )
635 allocate ( vcontact(ncontact) )
636 dmem(ng)=dmem(ng)+3.0_r8*real(ncontact,r8)
637!
638! Allocate arrays in grid connectivity structure.
639!
640 DO cr=1,ncontact
641 dg=donor_grid(cr)
642 rg=receiver_grid(cr)
643
644 allocate ( rcontact(cr) % Irg(npointsr(cr)) )
645 allocate ( ucontact(cr) % Irg(npointsu(cr)) )
646 allocate ( vcontact(cr) % Irg(npointsv(cr)) )
647
648 allocate ( rcontact(cr) % Jrg(npointsr(cr)) )
649 allocate ( ucontact(cr) % Jrg(npointsu(cr)) )
650 allocate ( vcontact(cr) % Jrg(npointsv(cr)) )
651
652 allocate ( rcontact(cr) % Idg(npointsr(cr)) )
653 allocate ( ucontact(cr) % Idg(npointsu(cr)) )
654 allocate ( vcontact(cr) % Idg(npointsv(cr)) )
655
656 allocate ( rcontact(cr) % Jdg(npointsr(cr)) )
657 allocate ( ucontact(cr) % Jdg(npointsu(cr)) )
658 allocate ( vcontact(cr) % Jdg(npointsv(cr)) )
659
660 dmem(dg)=dmem(dg)+4.0_r8*real(npointsr(cr),r8)
661 dmem(dg)=dmem(dg)+4.0_r8*real(npointsu(cr),r8)
662 dmem(dg)=dmem(dg)+4.0_r8*real(npointsv(cr),r8)
663
664# ifdef SOLVE3D
665 allocate ( rcontact(cr) % Kdg(n(dg),npointsr(cr)) )
666 allocate ( ucontact(cr) % Kdg(n(dg),npointsu(cr)) )
667 allocate ( vcontact(cr) % Kdg(n(dg),npointsv(cr)) )
668
669 dmem(dg)=dmem(dg)+real(n(dg)*npointsr(cr),r8)
670 dmem(dg)=dmem(dg)+real(n(dg)*npointsu(cr),r8)
671 dmem(dg)=dmem(dg)+real(n(dg)*npointsv(cr),r8)
672# endif
673
674 allocate ( rcontact(cr) % Lweight(4,npointsr(cr)) )
675 allocate ( ucontact(cr) % Lweight(4,npointsu(cr)) )
676 allocate ( vcontact(cr) % Lweight(4,npointsv(cr)) )
677
678 dmem(dg)=dmem(dg)+4.0_r8*real(npointsr(cr),r8)
679 dmem(dg)=dmem(dg)+4.0_r8*real(npointsu(cr),r8)
680 dmem(dg)=dmem(dg)+4.0_r8*real(npointsv(cr),r8)
681
682# ifdef WET_DRY
683 allocate ( rcontact(cr) % LweightUnmasked(4,npointsr(cr)) )
684 allocate ( ucontact(cr) % LweightUnmasked(4,npointsu(cr)) )
685 allocate ( vcontact(cr) % LweightUnmasked(4,npointsv(cr)) )
686
687 dmem(dg)=dmem(dg)+4.0_r8*real(npointsr(cr),r8)
688 dmem(dg)=dmem(dg)+4.0_r8*real(npointsu(cr),r8)
689 dmem(dg)=dmem(dg)+4.0_r8*real(npointsv(cr),r8)
690# endif
691
692# ifdef QUADRATIC_WEIGHTS
693 allocate ( rcontact(cr) % Qweight(9,npointsr(cr)) )
694 allocate ( ucontact(cr) % Qweight(9,npointsu(cr)) )
695 allocate ( vcontact(cr) % Qweight(9,npointsv(cr)) )
696
697 dmem(dg)=dmem(dg)+9.0_r8*real(npointsr(cr),r8)
698 dmem(dg)=dmem(dg)+9.0_r8*real(npointsu(cr),r8)
699 dmem(dg)=dmem(dg)+9.0_r8*real(npointsv(cr),r8)
700
701# ifdef WET_DRY
702 allocate ( rcontact(cr) % QweightUnmasked(9,npointsr(cr)) )
703 allocate ( ucontact(cr) % QweightUnmasked(9,npointsu(cr)) )
704 allocate ( vcontact(cr) % QweightUnmasked(9,npointsv(cr)) )
705
706 dmem(dg)=dmem(dg)+9.0_r8*real(npointsr(cr),r8)
707 dmem(dg)=dmem(dg)+9.0_r8*real(npointsu(cr),r8)
708 dmem(dg)=dmem(dg)+9.0_r8*real(npointsv(cr),r8)
709# endif
710# endif
711
712# ifdef SOLVE3D
713 allocate ( rcontact(cr) % Vweight(2,n(dg),npointsr(cr)) )
714 allocate ( ucontact(cr) % Vweight(2,n(dg),npointsu(cr)) )
715 allocate ( vcontact(cr) % Vweight(2,n(dg),npointsv(cr)) )
716
717 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsr(cr),r8)
718 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsu(cr),r8)
719 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsv(cr),r8)
720
721# if defined TANGENT || defined TL_IOMS
722 allocate ( rcontact(cr) % tl_Vweight(2,n(dg),npointsr(cr)) )
723 allocate ( ucontact(cr) % tl_Vweight(2,n(dg),npointsu(cr)) )
724 allocate ( vcontact(cr) % tl_Vweight(2,n(dg),npointsv(cr)) )
725
726 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsr(cr),r8)
727 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsu(cr),r8)
728 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsv(cr),r8)
729# endif
730
731# ifdef ADJOINT
732 allocate ( rcontact(cr) % ad_Vweight(2,n(dg),npointsr(cr)) )
733 allocate ( ucontact(cr) % ad_Vweight(2,n(dg),npointsu(cr)) )
734 allocate ( vcontact(cr) % ad_Vweight(2,n(dg),npointsv(cr)) )
735
736 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsr(cr),r8)
737 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsu(cr),r8)
738 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsv(cr),r8)
739# endif
740# endif
741 END DO
742!
743! Initialize grid connectivity structure.
744!
745 DO cr=1,ncontact
746 dg=donor_grid(cr)
747 rg=receiver_grid(cr)
748
749 rcontact(cr) % coincident = lcoincident(rg)
750 ucontact(cr) % coincident = lcoincident(rg)
751 vcontact(cr) % coincident = lcoincident(rg)
752
753 rcontact(cr) % interpolate = linterpolate(rg)
754 ucontact(cr) % interpolate = linterpolate(rg)
755 vcontact(cr) % interpolate = linterpolate(rg)
756
757 rcontact(cr) % donor_grid = dg
758 ucontact(cr) % donor_grid = dg
759 vcontact(cr) % donor_grid = dg
760
761 rcontact(cr) % receiver_grid = rg
762 ucontact(cr) % receiver_grid = rg
763 vcontact(cr) % receiver_grid = rg
764
765 rcontact(cr) % Npoints = npointsr(cr)
766 ucontact(cr) % Npoints = npointsu(cr)
767 vcontact(cr) % Npoints = npointsv(cr)
768
769 DO m=1,npointsr(cr)
770 ip=m+nstrr(cr)-1
771 rcontact(cr) % Irg(m) = irg_cp(ip)
772 rcontact(cr) % Jrg(m) = jrg_cp(ip)
773 rcontact(cr) % Idg(m) = idg_cp(ip)
774 rcontact(cr) % Jdg(m) = jdg_cp(ip)
775# ifdef WET_DRY
776 rcontact(cr) % LweightUnmasked(1,m) = lweight(1,ip)
777 rcontact(cr) % LweightUnmasked(2,m) = lweight(2,ip)
778 rcontact(cr) % LweightUnmasked(3,m) = lweight(3,ip)
779 rcontact(cr) % LweightUnmasked(4,m) = lweight(4,ip)
780# endif
781 rcontact(cr) % Lweight(1,m) = lweight(1,ip)
782 rcontact(cr) % Lweight(2,m) = lweight(2,ip)
783 rcontact(cr) % Lweight(3,m) = lweight(3,ip)
784 rcontact(cr) % Lweight(4,m) = lweight(4,ip)
785# ifdef QUADRATIC_WEIGHTS
786# ifdef WET_DRY
787 rcontact(cr) % QweightUnmasked(1,m) = qweight(1,ip)
788 rcontact(cr) % QweightUnmasked(2,m) = qweight(2,ip)
789 rcontact(cr) % QweightUnmasked(3,m) = qweight(3,ip)
790 rcontact(cr) % QweightUnmasked(4,m) = qweight(4,ip)
791 rcontact(cr) % QweightUnmasked(5,m) = qweight(5,ip)
792 rcontact(cr) % QweightUnmasked(6,m) = qweight(6,ip)
793 rcontact(cr) % QweightUnmasked(7,m) = qweight(7,ip)
794 rcontact(cr) % QweightUnmasked(8,m) = qweight(8,ip)
795 rcontact(cr) % QweightUnmasked(9,m) = qweight(9,ip)
796# endif
797 rcontact(cr) % Qweight(1,m) = qweight(1,ip)
798 rcontact(cr) % Qweight(2,m) = qweight(2,ip)
799 rcontact(cr) % Qweight(3,m) = qweight(3,ip)
800 rcontact(cr) % Qweight(4,m) = qweight(4,ip)
801 rcontact(cr) % Qweight(5,m) = qweight(5,ip)
802 rcontact(cr) % Qweight(6,m) = qweight(6,ip)
803 rcontact(cr) % Qweight(7,m) = qweight(7,ip)
804 rcontact(cr) % Qweight(8,m) = qweight(8,ip)
805 rcontact(cr) % Qweight(9,m) = qweight(9,ip)
806# endif
807 END DO
808
809 DO m=1,npointsu(cr)
810 ip=m+nstru(cr)-1
811 ucontact(cr) % Irg(m) = irg_cp(ip)
812 ucontact(cr) % Jrg(m) = jrg_cp(ip)
813 ucontact(cr) % Idg(m) = idg_cp(ip)
814 ucontact(cr) % Jdg(m) = jdg_cp(ip)
815# ifdef WET_DRY
816 ucontact(cr) % LweightUnmasked(1,m) = lweight(1,ip)
817 ucontact(cr) % LweightUnmasked(2,m) = lweight(2,ip)
818 ucontact(cr) % LweightUnmasked(3,m) = lweight(3,ip)
819 ucontact(cr) % LweightUnmasked(4,m) = lweight(4,ip)
820# endif
821 ucontact(cr) % Lweight(1,m) = lweight(1,ip)
822 ucontact(cr) % Lweight(2,m) = lweight(2,ip)
823 ucontact(cr) % Lweight(3,m) = lweight(3,ip)
824 ucontact(cr) % Lweight(4,m) = lweight(4,ip)
825# ifdef QUADRATIC_WEIGHTS
826# ifdef WET_DRY
827 ucontact(cr) % QweightUnmasked(1,m) = qweight(1,ip)
828 ucontact(cr) % QweightUnmasked(2,m) = qweight(2,ip)
829 ucontact(cr) % QweightUnmasked(3,m) = qweight(3,ip)
830 ucontact(cr) % QweightUnmasked(4,m) = qweight(4,ip)
831 ucontact(cr) % QweightUnmasked(5,m) = qweight(5,ip)
832 ucontact(cr) % QweightUnmasked(6,m) = qweight(6,ip)
833 ucontact(cr) % QweightUnmasked(7,m) = qweight(7,ip)
834 ucontact(cr) % QweightUnmasked(8,m) = qweight(8,ip)
835 ucontact(cr) % QweightUnmasked(9,m) = qweight(9,ip)
836# endif
837 ucontact(cr) % Qweight(1,m) = qweight(1,ip)
838 ucontact(cr) % Qweight(2,m) = qweight(2,ip)
839 ucontact(cr) % Qweight(3,m) = qweight(3,ip)
840 ucontact(cr) % Qweight(4,m) = qweight(4,ip)
841 ucontact(cr) % Qweight(5,m) = qweight(5,ip)
842 ucontact(cr) % Qweight(6,m) = qweight(6,ip)
843 ucontact(cr) % Qweight(7,m) = qweight(7,ip)
844 ucontact(cr) % Qweight(8,m) = qweight(8,ip)
845 ucontact(cr) % Qweight(9,m) = qweight(9,ip)
846# endif
847 END DO
848
849 DO m=1,npointsv(cr)
850 ip=m+nstrv(cr)-1
851 vcontact(cr) % Irg(m) = irg_cp(ip)
852 vcontact(cr) % Jrg(m) = jrg_cp(ip)
853 vcontact(cr) % Idg(m) = idg_cp(ip)
854 vcontact(cr) % Jdg(m) = jdg_cp(ip)
855# ifdef WET_DRY
856 vcontact(cr) % LweightUnmasked(1,m) = lweight(1,ip)
857 vcontact(cr) % LweightUnmasked(2,m) = lweight(2,ip)
858 vcontact(cr) % LweightUnmasked(3,m) = lweight(3,ip)
859 vcontact(cr) % LweightUnmasked(4,m) = lweight(4,ip)
860# endif
861 vcontact(cr) % Lweight(1,m) = lweight(1,ip)
862 vcontact(cr) % Lweight(2,m) = lweight(2,ip)
863 vcontact(cr) % Lweight(3,m) = lweight(3,ip)
864 vcontact(cr) % Lweight(4,m) = lweight(4,ip)
865# ifdef QUADRATIC_WEIGHTS
866# ifdef WET_DRY
867 vcontact(cr) % QweightUnmasked(1,m) = qweight(1,ip)
868 vcontact(cr) % QweightUnmasked(2,m) = qweight(2,ip)
869 vcontact(cr) % QweightUnmasked(3,m) = qweight(3,ip)
870 vcontact(cr) % QweightUnmasked(4,m) = qweight(4,ip)
871 vcontact(cr) % QweightUnmasked(5,m) = qweight(5,ip)
872 vcontact(cr) % QweightUnmasked(6,m) = qweight(6,ip)
873 vcontact(cr) % QweightUnmasked(7,m) = qweight(7,ip)
874 vcontact(cr) % QweightUnmasked(8,m) = qweight(8,ip)
875 vcontact(cr) % QweightUnmasked(9,m) = qweight(9,ip)
876# endif
877 vcontact(cr) % Qweight(1,m) = qweight(1,ip)
878 vcontact(cr) % Qweight(2,m) = qweight(2,ip)
879 vcontact(cr) % Qweight(3,m) = qweight(3,ip)
880 vcontact(cr) % Qweight(4,m) = qweight(4,ip)
881 vcontact(cr) % Qweight(5,m) = qweight(5,ip)
882 vcontact(cr) % Qweight(6,m) = qweight(6,ip)
883 vcontact(cr) % Qweight(7,m) = qweight(7,ip)
884 vcontact(cr) % Qweight(8,m) = qweight(8,ip)
885 vcontact(cr) % Qweight(9,m) = qweight(9,ip)
886# endif
887 END DO
888
889 END DO
890!
891!-----------------------------------------------------------------------
892! Load contact points grid metrics.
893!-----------------------------------------------------------------------
894!
895! Allocate contact region metrics (type T_NGM) structure.
896!
897 allocate ( contact_metric(ncontact) )
898!
899! Allocate arrays in contact region metrics structure.
900!
901 DO cr=1,ncontact
902 allocate ( contact_metric(cr) % angler(npointsr(cr)) )
903 allocate ( contact_metric(cr) % dndx (npointsr(cr)) )
904 allocate ( contact_metric(cr) % dmde (npointsr(cr)) )
905 allocate ( contact_metric(cr) % f (npointsr(cr)) )
906 allocate ( contact_metric(cr) % h (npointsr(cr)) )
907 allocate ( contact_metric(cr) % rmask (npointsr(cr)) )
908 allocate ( contact_metric(cr) % umask (npointsu(cr)) )
909 allocate ( contact_metric(cr) % vmask (npointsv(cr)) )
910 allocate ( contact_metric(cr) % pm (npointsr(cr)) )
911 allocate ( contact_metric(cr) % pn (npointsr(cr)) )
912 allocate ( contact_metric(cr) % Xr (npointsr(cr)) )
913 allocate ( contact_metric(cr) % Yr (npointsr(cr)) )
914 allocate ( contact_metric(cr) % Xu (npointsu(cr)) )
915 allocate ( contact_metric(cr) % Yu (npointsu(cr)) )
916 allocate ( contact_metric(cr) % Xv (npointsv(cr)) )
917 allocate ( contact_metric(cr) % Yv (npointsv(cr)) )
918
919 dmem(ng)=dmem(ng)+10.0_r8*real(npointsr(cr),r8)
920 dmem(ng)=dmem(ng)+ 3.0_r8*real(npointsu(cr),r8)
921 dmem(ng)=dmem(ng)+ 3.0_r8*real(npointsv(cr),r8)
922 END DO
923!
924! Initialize contact region metrics structure.
925!
926 DO cr=1,ncontact
927 DO m=1,npointsr(cr)
928 ip=m+nstrr(cr)-1
929 contact_metric(cr) % angler(m) = angle(ip)
930 contact_metric(cr) % dndx (m) = dndx(ip)
931 contact_metric(cr) % dmde (m) = dmde(ip)
932 contact_metric(cr) % f (m) = f(ip)
933 contact_metric(cr) % h (m) = h(ip)
934 contact_metric(cr) % rmask (m) = mask(ip)
935 contact_metric(cr) % pm (m) = pm(ip)
936 contact_metric(cr) % pn (m) = pn(ip)
937 contact_metric(cr) % Xr (m) = xrg(ip)
938 contact_metric(cr) % Yr (m) = yrg(ip)
939 END DO
940
941 DO m=1,npointsu(cr)
942 ip=m+nstru(cr)-1
943 contact_metric(cr) % umask(m) = mask(ip)
944 contact_metric(cr) % Xu (m) = xrg(ip)
945 contact_metric(cr) % Yu (m) = yrg(ip)
946 END DO
947
948 DO m=1,npointsv(cr)
949 ip=m+nstrv(cr)-1
950 contact_metric(cr) % vmask(m) = mask(ip)
951 contact_metric(cr) % Xv (m) = xrg(ip)
952 contact_metric(cr) % Yv (m) = yrg(ip)
953 END DO
954 END DO
955!
956!-----------------------------------------------------------------------
957! Initialize various parameters.
958!-----------------------------------------------------------------------
959!
960! Nested grids conectivity switches. They are used in "get_bounds.F"
961! to determenine if the global state array need to be extended by few
962! extra points at each boundary to accomodate the contact regions and
963! contact points. These are VERY important switches in ROMS nesting.
964!
965 IF (.not.allocated(contactregion)) THEN
966 allocate ( contactregion(4,ngrids) )
967 contactregion = .false.
968 dmem(ng)=dmem(ng)+4.0_r8*real(ngrids,r8)
969 END IF
970!
971 DO m=1,ncdatum
972 cr=contact_region(m)
973 rg=receiver_grid(cr)
974 ibry=on_boundary(m)
975 IF ((ibry.eq.iwest ).or.(ibry.eq.ieast ).or. &
976 & (ibry.eq.isouth).or.(ibry.eq.inorth)) THEN
977 IF (.not.contactregion(ibry,rg)) THEN
978 contactregion(ibry,rg)=.true.
979 END IF
980 END IF
981 END DO
982!
983! Set composite and refinement grids switches.
984!
985 DO cr=1,ncontact
986 rg=receiver_grid(cr)
987 DO ibry=1,4
988 IF (lcomposite(cr).and.contactregion(ibry,rg)) THEN
989 compositegrid(ibry,rg)=.true.
990 END IF
991 END DO
992 refinedgrid(rg)=lrefinement(cr)
993 END DO
994!
995! Set the switch to compute vertical interpolation weights. Currently,
996! they are only needed in non-coincident composite grids.
997!
998 IF (any(lcoincident).or.any(lcomposite)) THEN
999 get_vweights=.true.
1000 ELSE
1001 get_vweights=.false.
1002 END IF
1003!
1004! Refinement grid scale factor.
1005!
1006 DO ig=1,ngrids
1007 refinescale(ig)=refine_factor(ig)
1008 END DO
1009!
1010! Coarser donor grid number to finer grid contact points. Their
1011! values are computed in "metrics.F" when the maximum grid spacing
1012! is computed. The maximum grid spacing is used to distinguish
1013! contact areas for more than two layers of nesting.
1014!
1015!$OMP PARALLEL
1016 IF (.not.allocated(coarserdonor)) THEN
1017 allocate ( coarserdonor(ngrids) )
1018 coarserdonor = 0
1019 dmem(ng)=dmem(ng)+real(ngrids,r8)
1020 END IF
1021!$OMP END PARALLEL
1022!
1023! Finer donor grid to coarser receiver grid for two-way feedback.
1024! Their values are computed in "metrics.F" when the maximum grid
1025! spacing is computed.
1026!
1027!$OMP PARALLEL
1028 IF (.not.allocated(finerdonor)) THEN
1029 allocate ( finerdonor(ngrids) )
1030 finerdonor = 0
1031 dmem(ng)=dmem(ng)+real(ngrids,r8)
1032 END IF
1033!$OMP END PARALLEL
1034!
1035! Logical switch indicating which coarser grid is a donor to a finer
1036! grid external contact points. Their values are computed in
1037! "metrics.F" after "CoarserDonor" is computed.
1038!
1039!$OMP PARALLEL
1040 IF (.not.allocated(donortofiner)) THEN
1041 allocate ( donortofiner(ngrids) )
1042 donortofiner = .false.
1043 dmem(ng)=dmem(ng)+real(ngrids,r8)
1044 END IF
1045!$OMP END PARALLEL
1046!
1047! Number of refined grid time-steps. Their values are initialized to
1048! the optimal values based on the spatial refinement ratio. These
1049! values are over-written in "metrics.F" to use the time-step
1050! size ratio between donor and receiver grid. The user is allowed
1051! to take larger divisible time-step with respect to the donor grid.
1052! The user is responsible to set the appropriate refined grid time-step
1053! for stability.
1054!
1055!$OMP PARALLEL
1056 IF (.not.allocated(refinesteps)) THEN
1057 allocate ( refinesteps(ngrids) )
1058 dmem(ng)=dmem(ng)+real(ngrids,r8)
1059 DO ig=1,ngrids
1060 refinesteps(ig)=refine_factor(ig)
1061 END DO
1062 END IF
1063!$OMP END PARALLEL
1064!
1065! Refine grid time-steps counter with respect the coarse grid (ng=1)
1066! single time-step.
1067!
1068 IF (.not.allocated(refinestepscounter)) THEN
1069 allocate ( refinestepscounter(ngrids) )
1070 refinestepscounter = 0
1071 dmem(ng)=dmem(ng)+real(ngrids,r8)
1072 END IF
1073!
1074! Interval used in the two-way exchange between fine and coarse grids.
1075!
1076 IF (.not.allocated(twowayinterval)) THEN
1077 allocate ( twowayinterval(ngrids) )
1078 twowayinterval = 0.0_r8
1079 dmem(ng)=dmem(ng)+real(ngrids,r8)
1080 END IF
1081!
1082! Switch indicating which refined grid(s) include finer refined grids:
1083! telescoping refinement.
1084!
1085!$OMP PARALLEL
1086 IF (.not.allocated(telescoping)) THEN
1087 allocate ( telescoping(ngrids) )
1088 telescoping = .false.
1089 dmem(ng)=dmem(ng)+real(ngrids,r8)
1090 END IF
1091!$OMP END PARALLEL
1092!
1093! Rolling index and time (seconds) used in the temporal interpolation
1094! of contact point data.
1095!
1096!$OMP PARALLEL
1097 IF (.not.allocated(rollingindex)) THEN
1098 allocate ( rollingindex(ncontact) )
1099 rollingindex = 0
1100 dmem(ng)=dmem(ng)+real(ncontact,r8)
1101 END IF
1102
1103 IF (.not.allocated(rollingtime)) THEN
1104 allocate ( rollingtime(2,ncontact) )
1105 rollingtime = 0
1106 dmem(ng)=dmem(ng)+2.0_r8*real(ncontact,r8)
1107 END IF
1108!$OMP END PARALLEL
1109!
1110!-----------------------------------------------------------------------
1111! Allocate composite grids contact regions structure and arrays.
1112!-----------------------------------------------------------------------
1113!
1114 IF (any(compositegrid)) THEN
1115!
1116! Allocate composite grid contact regions (type T_COMPOSITE) structure.
1117!
1118 allocate ( composite(ncontact) )
1119!
1120! Allocate arrays in composite grids contact regions structure.
1121!
1122 DO cr=1,ncontact
1123 dg=donor_grid(cr)
1124
1125 allocate ( composite(cr) % bustr(4,npointsu(cr)) )
1126 allocate ( composite(cr) % bvstr(4,npointsv(cr)) )
1127
1128 allocate ( composite(cr) % ubar(4,npointsu(cr),2) )
1129 allocate ( composite(cr) % vbar(4,npointsv(cr),2) )
1130 allocate ( composite(cr) % zeta(4,npointsr(cr),2) )
1131
1132 allocate ( composite(cr) % rzeta(4,npointsr(cr)) )
1133
1134 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsr(cr),r8)
1135 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsu(cr),r8)
1136 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsv(cr),r8)
1137
1138# if defined TANGENT || defined TL_IOMS
1139 allocate ( composite(cr) % tl_bustr(4,npointsu(cr)) )
1140 allocate ( composite(cr) % tl_bvstr(4,npointsv(cr)) )
1141
1142 allocate ( composite(cr) % tl_ubar(4,npointsu(cr),2) )
1143 allocate ( composite(cr) % tl_vbar(4,npointsv(cr),2) )
1144 allocate ( composite(cr) % tl_zeta(4,npointsr(cr),2) )
1145
1146 allocate ( composite(cr) % tl_rzeta(4,npointsr(cr)) )
1147
1148 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsr(cr),r8)
1149 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsu(cr),r8)
1150 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsv(cr),r8)
1151# endif
1152
1153# ifdef ADJOINT
1154 allocate ( composite(cr) % ad_bustr(4,npointsu(cr)) )
1155 allocate ( composite(cr) % ad_bvstr(4,npointsv(cr)) )
1156
1157 allocate ( composite(cr) % ad_ubar(4,npointsu(cr),2) )
1158 allocate ( composite(cr) % ad_vbar(4,npointsv(cr),2) )
1159 allocate ( composite(cr) % ad_zeta(4,npointsr(cr),2) )
1160
1161 allocate ( composite(cr) % ad_rzeta(4,npointsr(cr)) )
1162
1163 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsr(cr),r8)
1164 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsu(cr),r8)
1165 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsv(cr),r8)
1166# endif
1167
1168# ifdef SOLVE3D
1169 allocate ( composite(cr) % DU_avg1(4,npointsu(cr)) )
1170 allocate ( composite(cr) % DV_avg1(4,npointsv(cr)) )
1171 allocate ( composite(cr) % Zt_avg1(4,npointsr(cr)) )
1172
1173 dmem(dg)=dmem(dg)+real(4*npointsr(cr),r8)
1174 dmem(dg)=dmem(dg)+real(4*npointsu(cr),r8)
1175 dmem(dg)=dmem(dg)+real(4*npointsv(cr),r8)
1176
1177 allocate ( composite(cr) % u(4,n(dg),npointsu(cr)) )
1178 allocate ( composite(cr) % v(4,n(dg),npointsv(cr)) )
1179
1180 allocate ( composite(cr) % Huon(4,n(dg),npointsu(cr)) )
1181 allocate ( composite(cr) % Hvom(4,n(dg),npointsv(cr)) )
1182
1183 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsu(cr),r8)
1184 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsv(cr),r8)
1185
1186 allocate ( composite(cr) % t(4,n(dg),npointsr(cr),nt(dg)) )
1187
1188 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsu(cr)*nt(dg),r8)
1189
1190# if defined TANGENT || defined TL_IOMS
1191 allocate ( composite(cr) % tl_DU_avg1(4,npointsu(cr)) )
1192 allocate ( composite(cr) % tl_DV_avg1(4,npointsv(cr)) )
1193 allocate ( composite(cr) % tl_Zt_avg1(4,npointsr(cr)) )
1194
1195 dmem(dg)=dmem(dg)+real(4*npointsr(cr),r8)
1196 dmem(dg)=dmem(dg)+real(4*npointsu(cr),r8)
1197 dmem(dg)=dmem(dg)+real(4*npointsv(cr),r8)
1198
1199 allocate ( composite(cr) % tl_u(4,n(dg),npointsu(cr)) )
1200 allocate ( composite(cr) % tl_v(4,n(dg),npointsv(cr)) )
1201
1202 allocate ( composite(cr) % tl_Huon(4,n(dg),npointsu(cr)) )
1203 allocate ( composite(cr) % tl_Hvom(4,n(dg),npointsv(cr)) )
1204
1205 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsu(cr),r8)
1206 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsv(cr),r8)
1207
1208 allocate ( composite(cr) % tl_t(4,n(dg),npointsr(cr),nt(dg)) )
1209
1210 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsu(cr)*nt(dg),r8)
1211# endif
1212
1213# ifdef ADJOINT
1214 allocate ( composite(cr) % ad_DU_avg1(4,npointsu(cr)) )
1215 allocate ( composite(cr) % ad_DV_avg1(4,npointsv(cr)) )
1216 allocate ( composite(cr) % ad_Zt_avg1(4,npointsr(cr)) )
1217
1218 dmem(dg)=dmem(dg)+real(4*npointsr(cr),r8)
1219 dmem(dg)=dmem(dg)+real(4*npointsu(cr),r8)
1220 dmem(dg)=dmem(dg)+real(4*npointsv(cr),r8)
1221
1222 allocate ( composite(cr) % ad_u(4,n(dg),npointsu(cr)) )
1223 allocate ( composite(cr) % ad_v(4,n(dg),npointsv(cr)) )
1224
1225 allocate ( composite(cr) % ad_Huon(4,n(dg),npointsu(cr)) )
1226 allocate ( composite(cr) % ad_Hvom(4,n(dg),npointsv(cr)) )
1227
1228 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsu(cr),r8)
1229 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsv(cr),r8)
1230
1231 allocate ( composite(cr) % ad_t(4,n(dg),npointsr(cr),nt(dg)) )
1232
1233 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsu(cr)*nt(dg),r8)
1234# endif
1235# endif
1236 END DO
1237 END IF
1238!
1239!-----------------------------------------------------------------------
1240! Allocate refinement grids contact regions structure.
1241!-----------------------------------------------------------------------
1242!
1243 IF (any(refinedgrid)) THEN
1244!
1245! Allocate refinement grids contact region (type T_REFINED) structure.
1246!
1247 allocate ( refined(ncontact) )
1248!
1249! Allocate arrays in refinement grids contact region structure.
1250!
1251 DO cr=1,ncontact
1252 rg=receiver_grid(cr)
1253
1254 allocate ( refined(cr) % ubar(4,npointsu(cr),2) )
1255 allocate ( refined(cr) % vbar(4,npointsv(cr),2) )
1256 allocate ( refined(cr) % zeta(4,npointsr(cr),2) )
1257
1258 allocate ( refined(cr) % U2d_flux(4,npointsu(cr),2) )
1259 allocate ( refined(cr) % V2d_flux(4,npointsv(cr),2) )
1260
1261 allocate ( refined(cr) % on_u(npointsu(cr)) )
1262 allocate ( refined(cr) % om_v(npointsv(cr)) )
1263
1264 dmem(rg)=dmem(rg)+real(4*npointsr(cr),r8)
1265 dmem(rg)=dmem(rg)+3.0_r8*real(4*npointsu(cr),r8)
1266 dmem(rg)=dmem(rg)+3.0_r8*real(4*npointsv(cr),r8)
1267
1268# if defined TANGENT || defined TL_IOMS
1269 allocate ( refined(cr) % tl_ubar(4,npointsu(cr),2) )
1270 allocate ( refined(cr) % tl_vbar(4,npointsv(cr),2) )
1271 allocate ( refined(cr) % tl_zeta(4,npointsr(cr),2) )
1272
1273 allocate ( refined(cr) % tl_U2d_flux(4,npointsu(cr),2) )
1274 allocate ( refined(cr) % tl_V2d_flux(4,npointsv(cr),2) )
1275
1276 dmem(rg)=dmem(rg)+real(4*npointsr(cr),r8)
1277 dmem(rg)=dmem(rg)+3.0_r8*real(4*npointsu(cr),r8)
1278 dmem(rg)=dmem(rg)+3.0_r8*real(4*npointsv(cr),r8)
1279# endif
1280
1281# ifdef ADJOINT
1282 allocate ( refined(cr) % ad_ubar(4,npointsu(cr),2) )
1283 allocate ( refined(cr) % ad_vbar(4,npointsv(cr),2) )
1284 allocate ( refined(cr) % ad_zeta(4,npointsr(cr),2) )
1285
1286 allocate ( refined(cr) % ad_U2d_flux(4,npointsu(cr),2) )
1287 allocate ( refined(cr) % ad_V2d_flux(4,npointsv(cr),2) )
1288
1289 dmem(rg)=dmem(rg)+real(4*npointsr(cr),r8)
1290 dmem(rg)=dmem(rg)+3.0_r8*real(4*npointsu(cr),r8)
1291 dmem(rg)=dmem(rg)+3.0_r8*real(4*npointsv(cr),r8)
1292# endif
1293
1294
1295# ifdef SOLVE3D
1296 allocate ( refined(cr) % u(4,n(rg),npointsu(cr),2) )
1297 allocate ( refined(cr) % v(4,n(rg),npointsv(cr),2) )
1298
1299 dmem(rg)=dmem(rg)+3.0_r8*real(4*n(rg)*npointsu(cr),r8)
1300 dmem(rg)=dmem(rg)+3.0_r8*real(4*n(rg)*npointsv(cr),r8)
1301
1302 allocate ( refined(cr) % t(4,n(rg),npointsr(cr),2,nt(rg)) )
1303
1304 dmem(rg)=dmem(rg)+2.0_r8*real(4*n(rg)*npointsr(cr)*nt(rg),r8)
1305
1306# if defined TANGENT || defined TL_IOMS
1307 allocate ( refined(cr) % tl_u(4,n(rg),npointsu(cr),2) )
1308 allocate ( refined(cr) % tl_v(4,n(rg),npointsv(cr),2) )
1309
1310 dmem(rg)=dmem(rg)+3.0_r8*real(4*n(rg)*npointsu(cr),r8)
1311 dmem(rg)=dmem(rg)+3.0_r8*real(4*n(rg)*npointsv(cr),r8)
1312
1313 allocate ( refined(cr) % tl_t(4,n(rg),npointsr(cr),2,nt(rg)) )
1314
1315 dmem(rg)=dmem(rg)+2.0_r8*real(4*n(rg)*npointsr(cr)*nt(rg),r8)
1316# endif
1317
1318# ifdef ADJOINT
1319 allocate ( refined(cr) % ad_u(4,n(rg),npointsu(cr),2) )
1320 allocate ( refined(cr) % ad_v(4,n(rg),npointsv(cr),2) )
1321
1322 dmem(rg)=dmem(rg)+3.0_r8*real(4*n(rg)*npointsu(cr),r8)
1323 dmem(rg)=dmem(rg)+3.0_r8*real(4*n(rg)*npointsv(cr),r8)
1324
1325 allocate ( refined(cr) % ad_t(4,n(rg),npointsr(cr),2,nt(rg)) )
1326
1327 dmem(rg)=dmem(rg)+2.0_r8*real(4*n(rg)*npointsr(cr)*nt(rg),r8)
1328# endif
1329# endif
1330 END DO
1331 END IF
1332!
1333 RETURN

References mod_nesting::coarserdonor, mod_nesting::composite, mod_scalars::compositegrid, mod_nesting::contact_metric, mod_nesting::contact_region, mod_nesting::contactregion, mod_param::dmem, mod_nesting::donor_grid, mod_nesting::donortofiner, mod_scalars::exit_flag, mod_nesting::finerdonor, strings_mod::founderror(), mod_nesting::get_vweights, mod_nesting::i_left, mod_nesting::i_right, mod_nesting::idg_cp, mod_scalars::ieast, mod_scalars::inorth, mod_nesting::irg_cp, mod_scalars::isouth, mod_scalars::iwest, mod_nesting::j_bottom, mod_nesting::j_top, mod_nesting::jdg_cp, mod_nesting::jrg_cp, mod_param::lm, mod_parallel::master, mod_param::mm, mod_param::n, mod_nesting::ncdatum, mod_nesting::ncontact, mod_nesting::ncpoints, mod_nesting::nendr, mod_nesting::nendu, mod_nesting::nendv, mod_netcdf::netcdf_close(), mod_netcdf::netcdf_get_dim(), mod_netcdf::netcdf_inq_var(), mod_netcdf::netcdf_open(), mod_iounits::ngcname, mod_param::ngrids, mod_scalars::noerror, mod_nesting::nstrr, mod_nesting::nstru, mod_nesting::nstrv, mod_param::nt, mod_nesting::on_boundary, mod_nesting::rcontact, mod_nesting::receiver_grid, mod_nesting::refined, mod_scalars::refinedgrid, mod_scalars::refinescale, mod_nesting::refinesteps, mod_nesting::refinestepscounter, mod_nesting::rollingindex, mod_nesting::rollingtime, mod_iounits::stdout, mod_nesting::telescoping, mod_nesting::twowayinterval, mod_nesting::ucontact, and mod_nesting::vcontact.

Referenced by set_contact().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_contact_pio()

subroutine, private set_contact_mod::set_contact_pio ( integer, intent(in) ng,
integer, intent(in) model )
private

Definition at line 1339 of file set_contact.F.

1340!***********************************************************************
1341!
1342! Imported variable declarations.
1343!
1344 integer, intent(in) :: ng, model
1345!
1346! Local variable declarations.
1347!
1348 logical, dimension((Ngrids-1)*2) :: Lcoincident
1349 logical, dimension((Ngrids-1)*2) :: Lcomposite
1350 logical, dimension((Ngrids-1)*2) :: Lmosaic
1351 logical, dimension((Ngrids-1)*2) :: Lrefinement
1352
1353 logical, dimension((Ngrids-1)*2) :: Linterpolate
1354!
1355 integer :: cr, dg, ibry, ic, ig, ip, m, rg, vindex
1356 integer :: my_Ncontact, my_Ngrids, my_nLweights, my_nQweights
1357
1358 integer, dimension(Ngrids) :: my_Lm, my_Mm
1359 integer, dimension(Ngrids) :: refine_factor
1360
1361 integer, dimension((Ngrids-1)*2) :: NpointsR
1362 integer, dimension((Ngrids-1)*2) :: NpointsU
1363 integer, dimension((Ngrids-1)*2) :: NpointsV
1364!
1365 real(r8), allocatable :: Lweight(:,:)
1366# ifdef QUADRATIC_WEIGHTS
1367 real(r8), allocatable :: Qweight(:,:)
1368# endif
1369 real(r8), allocatable :: Xrg(:), Yrg(:)
1370 real(r8), allocatable :: angle(:)
1371 real(r8), allocatable :: dmde(:), dndx(:)
1372 real(r8), allocatable :: f(:)
1373 real(r8), allocatable :: h(:)
1374 real(r8), allocatable :: mask(:)
1375 real(r8), allocatable :: pm(:), pn(:)
1376!
1377 character (len=*), parameter :: MyFile = &
1378 & __FILE__//", set_contact_pio"
1379!
1380 TYPE (file_desc_t) :: NGCpioFile
1381!
1382 sourcefile=myfile
1383!
1384!-----------------------------------------------------------------------
1385! Read in nesting grids contact point packed information into local
1386! arrays.
1387!-----------------------------------------------------------------------
1388!
1389! Open contact points NetCDF file for reading.
1390!
1391 CALL pio_netcdf_open (ng, model, ngcname, 0, ngcpiofile)
1392 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
1393 WRITE (stdout,20) trim(ngcname)
1394 20 FORMAT (/,' SET_CONTACT_PIO - unable to open contact points ', &
1395 & ' NetCDF file: ',a)
1396 RETURN
1397 END IF
1398!
1399! Inquire about the variables.
1400!
1401 CALL pio_netcdf_inq_var (ng, model, ngcname, &
1402 & piofile = ngcpiofile)
1403 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1404!
1405! Check if the number of nested grids is correct.
1406!
1407 CALL pio_netcdf_get_dim (ng, model, ngcname, &
1408 & piofile = ngcpiofile, &
1409 & dimname = 'Ngrids', &
1410 & dimsize = my_ngrids)
1411 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1412
1413 IF (my_ngrids.ne.ngrids) THEN
1414 IF (master) THEN
1415 WRITE (stdout,10) 'inconsistent parameter, Ngrids = ', &
1416 & ngrids, my_ngrids
1417 10 FORMAT (/,' SET_CONTACT_PIO - ', a, i4, 2x, i4, &
1418 & /,19x,'in input file:'2x,a)
1419 END IF
1420 exit_flag=5
1421 RETURN
1422 END IF
1423!
1424! Check number of contact regions, Ncontact = (Ngrids-1)*2.
1425!
1426 CALL pio_netcdf_get_dim (ng, model, ngcname, &
1427 & piofile = ngcpiofile, &
1428 & dimname = 'Ncontact', &
1429 & dimsize = my_ncontact)
1430 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1431
1432 IF (my_ncontact.ne.(ngrids-1)*2) THEN
1433 IF (master) THEN
1434 WRITE (stdout,10) 'inconsistent parameter, Ncontact = ', &
1435 & (ngrids-1)*2, my_ncontact
1436 END IF
1437 exit_flag=5
1438 RETURN
1439 END IF
1440 ncontact=(ngrids-1)*2
1441!
1442! Check number of contact points for bilinear interpolation weights.
1443!
1444 CALL pio_netcdf_get_dim (ng, model, ngcname, &
1445 & piofile = ngcpiofile, &
1446 & dimname = 'nLweights', &
1447 & dimsize = my_nlweights)
1448 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1449
1450 IF (my_nlweights.ne.4) THEN
1451 IF (master) THEN
1452 WRITE (stdout,10) 'inconsistent parameter, nLweights = ', &
1453 & 4, my_nlweights
1454 END IF
1455 exit_flag=5
1456 RETURN
1457 END IF
1458
1459# ifdef QUADRATIC_WEIGHTS
1460!
1461! Check number of contact points for quadratic interpolation weights.
1462!
1463 CALL pio_netcdf_get_dim (ng, model, ngcname, &
1464 & piofile = ngcpiofile, &
1465 & dimname = 'nQweights', &
1466 & dimsize = my_nqweights)
1467 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1468
1469 IF (my_nqweights.ne.9) THEN
1470 IF (master) THEN
1471 WRITE (stdout,10) 'inconsistent parameter, nQweights = ', &
1472 & 9, my_nqweights
1473 END IF
1474 exit_flag=5
1475 RETURN
1476 END IF
1477# endif
1478!
1479! Get number of contact points in input NetCDF file. It may include
1480! or not contact point over land when land/sea masking is activated.
1481!
1482 CALL pio_netcdf_get_dim (ng, model, ngcname, &
1483 & piofile = ngcpiofile, &
1484 & dimname = 'datum', &
1485 & dimsize = ncdatum)
1486 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1487!
1488! Read in grid dimensions and if the grid order is correct. The order
1489! of the grids is important in nesting.
1490!
1491 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1492 & 'Lm', my_lm, &
1493 & piofile = ngcpiofile)
1494 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1495
1496 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1497 & 'Mm', my_mm, &
1498 & piofile = ngcpiofile)
1499 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1500
1501 DO ig=1,ngrids
1502 IF (my_lm(ig).ne.lm(ig)) THEN
1503 IF (master) THEN
1504 WRITE (stdout,10) 'inconsistent grid order, Lm = ', &
1505 & lm(ig), my_lm(ig)
1506 END IF
1507 exit_flag=5
1508 RETURN
1509 END IF
1510 IF (my_mm(ig).ne.mm(ig)) THEN
1511 IF (master) THEN
1512 WRITE (stdout,10) 'inconsistent grid order, Mm = ', &
1513 & mm(ig), my_mm(ig)
1514 END IF
1515 exit_flag=5
1516 RETURN
1517 END IF
1518 END DO
1519!
1520! Read in nesting types logical switches, 1:Ngrids.
1521!
1522 CALL pio_netcdf_get_lvar (ng, model, ngcname, &
1523 & 'coincident', lcoincident, &
1524 & piofile = ngcpiofile)
1525 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1526
1527 CALL pio_netcdf_get_lvar (ng, model, ngcname, &
1528 & 'composite', lcomposite, &
1529 & piofile = ngcpiofile)
1530 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1531
1532 CALL pio_netcdf_get_lvar (ng, model, ngcname, &
1533 & 'mosaic', lmosaic, &
1534 & piofile = ngcpiofile)
1535 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1536
1537 CALL pio_netcdf_get_lvar (ng, model, ngcname, &
1538 & 'refinement', lrefinement, &
1539 & piofile = ngcpiofile)
1540 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1541!
1542! Read in refinement factor from donor grid, 1:Ngrids.
1543!
1544 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1545 & 'refine_factor', refine_factor, &
1546 & piofile = ngcpiofile)
1547 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1548!
1549! Read in vertical interpolation at contact points switch, 1:Ncontact.
1550!
1551 CALL pio_netcdf_get_lvar (ng, model, ngcname, &
1552 & 'interpolate', linterpolate, &
1553 & piofile = ngcpiofile)
1554 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1555!
1556! Read in data donor and data receiver grid number, 1:Ncontact.
1557!
1558 IF (.not.allocated(donor_grid)) THEN
1559 allocate ( donor_grid((ngrids-1)*2) )
1560 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
1561 END IF
1562 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1563 & 'donor_grid', donor_grid, &
1564 & piofile = ngcpiofile)
1565 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1566
1567 IF (.not.allocated(receiver_grid)) THEN
1568 allocate ( receiver_grid((ngrids-1)*2) )
1569 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
1570 END IF
1571 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1572 & 'receiver_grid', receiver_grid, &
1573 & piofile = ngcpiofile)
1574 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1575!
1576! Read in donor grid (I,J) indices at PSI points used to extract
1577! refined grid. Values have a fill value of -999 if not applicable.
1578!
1579 IF (.not.allocated(i_left)) THEN
1580 allocate ( i_left(ngrids) )
1581 dmem(ng)=dmem(ng)+real(ngrids,r8)
1582 END IF
1583 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1584 & 'I_left', i_left, &
1585 & piofile = ngcpiofile)
1586 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1587!
1588 IF (.not.allocated(i_right)) THEN
1589 allocate ( i_right(ngrids) )
1590 dmem(ng)=dmem(ng)+real(ngrids,r8)
1591 END IF
1592 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1593 & 'I_right', i_right, &
1594 & piofile = ngcpiofile)
1595 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1596!
1597 IF (.not.allocated(j_bottom)) THEN
1598 allocate ( j_bottom(ngrids) )
1599 dmem(ng)=dmem(ng)+real(ngrids,r8)
1600 END IF
1601 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1602 & 'J_bottom', j_bottom, &
1603 & piofile = ngcpiofile)
1604 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1605!
1606 IF (.not.allocated(j_top)) THEN
1607 allocate ( j_top(ngrids) )
1608 dmem(ng)=dmem(ng)+real(ngrids,r8)
1609 END IF
1610 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1611 & 'J_top', j_top, &
1612 & piofile = ngcpiofile)
1613 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1614!
1615! Read in starting and ending contact points (RHO-, U-, V-) index in
1616! packed data vectors for each contact region, 1:Ncontact.
1617!
1618 IF (.not.allocated(nstrr)) THEN
1619 allocate ( nstrr((ngrids-1)*2) )
1620 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
1621 END IF
1622 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1623 & 'NstrR', nstrr, &
1624 & piofile = ngcpiofile)
1625 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1626!
1627 IF (.not.allocated(nendr)) THEN
1628 allocate ( nendr((ngrids-1)*2) )
1629 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
1630 END IF
1631 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1632 & 'NendR', nendr, &
1633 & piofile = ngcpiofile)
1634 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1635!
1636 IF (.not.allocated(nstru)) THEN
1637 allocate ( nstru((ngrids-1)*2) )
1638 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
1639 END IF
1640 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1641 & 'NstrU', nstru, &
1642 & piofile = ngcpiofile)
1643 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1644!
1645 IF (.not.allocated(nendu)) THEN
1646 allocate ( nendu((ngrids-1)*2) )
1647 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
1648 END IF
1649 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1650 & 'NendU', nendu, &
1651 & piofile = ngcpiofile)
1652 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1653!
1654 IF (.not.allocated(nstrv)) THEN
1655 allocate ( nstrv((ngrids-1)*2) )
1656 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
1657 END IF
1658 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1659 & 'NstrV', nstrv, &
1660 & piofile = ngcpiofile)
1661 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1662!
1663 IF (.not.allocated(nendv)) THEN
1664 allocate ( nendv((ngrids-1)*2) )
1665 dmem(ng)=dmem(ng)+real((ngrids-1)*2,r8)
1666 END IF
1667 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1668 & 'NendV', nendv, &
1669 & piofile = ngcpiofile)
1670 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1671!
1672! Read in contact region number for each contact point, 1:NCdatum.
1673!
1674 IF (.not.allocated(contact_region)) THEN
1675 allocate ( contact_region(ncdatum) )
1676 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1677 END IF
1678 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1679 & 'contact_region', contact_region, &
1680 & piofile = ngcpiofile)
1681 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1682!
1683! Read in flag (1:NCdatum) to determine if the contact point is at
1684! receiver grid interior (zero value) or on the boundary (1:western,
1685! 2:southern, 3:eastern, 4:northern).
1686!
1687 IF (.not.allocated(on_boundary)) THEN
1688 allocate ( on_boundary(ncdatum) )
1689 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1690 END IF
1691 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1692 & 'on_boundary', on_boundary, &
1693 & piofile = ngcpiofile)
1694 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1695!
1696! Read in I-left and J-bottom indices of donor cell containing contact
1697! point, 1:NCdatum.
1698!
1699 IF (.not.allocated(idg_cp)) THEN
1700 allocate ( idg_cp(ncdatum) )
1701 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1702 END IF
1703 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1704 & 'Idg', idg_cp, &
1705 & piofile = ngcpiofile)
1706 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1707
1708 IF (.not.allocated(jdg_cp)) THEN
1709 allocate ( jdg_cp(ncdatum) )
1710 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1711 END IF
1712 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1713 & 'Jdg', jdg_cp, &
1714 & piofile = ngcpiofile)
1715 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1716!
1717! Read in I- and J-indices of receiver grid contact point, 1:NCdatum.
1718!
1719 IF (.not.allocated(irg_cp)) THEN
1720 allocate ( irg_cp(ncdatum) )
1721 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1722 END IF
1723 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1724 & 'Irg', irg_cp, &
1725 & piofile = ngcpiofile)
1726 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1727
1728 IF (.not.allocated(jrg_cp)) THEN
1729 allocate ( jrg_cp(ncdatum) )
1730 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1731 END IF
1732 CALL pio_netcdf_get_ivar (ng, model, ngcname, &
1733 & 'Jrg', jrg_cp, &
1734 & piofile = ngcpiofile)
1735 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1736!
1737! Read in X- and Y-location of receiver grid contact point, 1:NCdatum.
1738!
1739 IF (.not.allocated(xrg)) THEN
1740 allocate ( xrg(ncdatum) )
1741 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1742 END IF
1743 CALL pio_netcdf_get_fvar (ng, model, ngcname, &
1744 & 'Xrg', xrg, &
1745 & piofile = ngcpiofile)
1746 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1747
1748 IF (.not.allocated(yrg)) THEN
1749 allocate ( yrg(ncdatum) )
1750 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1751 END IF
1752 CALL pio_netcdf_get_fvar (ng, model, ngcname, &
1753 & 'Yrg', yrg, &
1754 & piofile = ngcpiofile)
1755 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1756!
1757! Read in contact point horizontal linear interpolation weights,
1758! 1:NCdatum.
1759!
1760 IF (.not.allocated(lweight)) THEN
1761 allocate ( lweight(my_nlweights,ncdatum) )
1762 dmem(ng)=dmem(ng)+real(my_nlweights*ncdatum,r8)
1763 END IF
1764 CALL pio_netcdf_get_fvar (ng, model, ngcname, &
1765 & 'Lweight', lweight, &
1766 & piofile = ngcpiofile)
1767 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1768
1769# ifdef QUADRATIC_WEIGHTS
1770!
1771! Read in contact point horizontal quadratic interpolation weights,
1772! 1:NCdatum.
1773!
1774 IF (.not.allocated(qweight)) THEN
1775 allocate ( qweight(my_nqweights,ncdatum) )
1776 dmem(ng)=dmem(ng)+real(my_nqweights*ncdatum,r8)
1777 END IF
1778 CALL pio_netcdf_get_fvar (ng, model, ngcname, &
1779 & 'Qweight', qweight, &
1780 & piofile = ngcpiofile)
1781 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1782# endif
1783!
1784! Read in contact point bathymetry at RHO-points, 1:NCdatum.
1785!
1786 IF (.not.allocated(h)) THEN
1787 allocate ( h(ncdatum) )
1788 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1789 END IF
1790 CALL pio_netcdf_get_fvar (ng, model, ngcname, &
1791 & 'h', h, &
1792 & piofile = ngcpiofile)
1793 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1794!
1795! Read in contact point Coriolis parameter at RHO-points, 1:NCdatum.
1796!
1797 IF (.not.allocated(f)) THEN
1798 allocate ( f(ncdatum) )
1799 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1800 END IF
1801 CALL pio_netcdf_get_fvar (ng, model, ngcname, &
1802 & 'f', f, &
1803 & piofile = ngcpiofile)
1804 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1805!
1806! Read in contact point curvilinar coordinates in XI- and ETA-direction
1807! at RHO-points, 1:NCdatum.
1808!
1809 IF (.not.allocated(pm)) THEN
1810 allocate ( pm(ncdatum) )
1811 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1812 END IF
1813 CALL pio_netcdf_get_fvar (ng, model, ngcname, &
1814 & 'pm', pm, &
1815 & piofile = ngcpiofile)
1816 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1817
1818 IF (.not.allocated(pn)) THEN
1819 allocate ( pn(ncdatum) )
1820 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1821 END IF
1822 CALL pio_netcdf_get_fvar (ng, model, ngcname, &
1823 & 'pn', pn, &
1824 & piofile = ngcpiofile)
1825 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1826!
1827! Read in contact point inverse metric factor in XI- and ETA-direction,
1828! 1:NCdatum.
1829!
1830 IF (.not.allocated(dndx)) THEN
1831 allocate ( dndx(ncdatum) )
1832 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1833 END IF
1834 CALL pio_netcdf_get_fvar (ng, model, ngcname, &
1835 & 'dndx', dndx, &
1836 & piofile = ngcpiofile)
1837 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1838
1839 IF (.not.allocated(dmde)) THEN
1840 allocate ( dmde(ncdatum) )
1841 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1842 END IF
1843 CALL pio_netcdf_get_fvar (ng, model, ngcname, &
1844 & 'dmde', dmde, &
1845 & piofile = ngcpiofile)
1846 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1847!
1848! Read in contact point angle between XI-axis and EAST, 1:NCdatum.
1849!
1850 IF (.not.allocated(angle)) THEN
1851 allocate ( angle(ncdatum) )
1852 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1853 END IF
1854 CALL pio_netcdf_get_fvar (ng, model, ngcname, &
1855 & 'angle', angle, &
1856 & piofile = ngcpiofile)
1857 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1858!
1859! Read in contact point land/sea mask, 1:NCdatum.
1860!
1861 IF (.not.allocated(mask)) THEN
1862 allocate ( mask(ncdatum) )
1863 dmem(ng)=dmem(ng)+real(ncdatum,r8)
1864 END IF
1865 CALL pio_netcdf_get_fvar (ng, model, ngcname, &
1866 & 'mask', mask, &
1867 & piofile = ngcpiofile)
1868 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1869!
1870! Close contact points NetCDF file.
1871!
1872 CALL pio_netcdf_close (ng, model, ngcpiofile, ngcname, .false.)
1873 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1874!
1875!-----------------------------------------------------------------------
1876! Unpack contact point data into nesting structures (type T_NGC).
1877!-----------------------------------------------------------------------
1878!
1879! Determine number of contact points in each contact region.
1880!
1881 IF (.not.allocated(ncpoints)) THEN
1882 allocate ( ncpoints(ncontact) )
1883 dmem(ng)=dmem(ng)+real(ncontact,r8)
1884 END IF
1885 DO cr=1,ncontact
1886 npointsr(cr)=nendr(cr)-nstrr(cr)+1
1887 npointsu(cr)=nendu(cr)-nstru(cr)+1
1888 npointsv(cr)=nendv(cr)-nstrv(cr)+1
1889 ncpoints(cr)=npointsr(cr)+npointsu(cr)+npointsv(cr)
1890 END DO
1891!
1892! Allocate grid connectivity (type T_NGC) structures.
1893!
1894 allocate ( rcontact(ncontact) )
1895 allocate ( ucontact(ncontact) )
1896 allocate ( vcontact(ncontact) )
1897 dmem(ng)=dmem(ng)+3.0_r8*real(ncontact,r8)
1898!
1899! Allocate arrays in grid connectivity structure.
1900!
1901 DO cr=1,ncontact
1902 dg=donor_grid(cr)
1903 rg=receiver_grid(cr)
1904
1905 allocate ( rcontact(cr) % Irg(npointsr(cr)) )
1906 allocate ( ucontact(cr) % Irg(npointsu(cr)) )
1907 allocate ( vcontact(cr) % Irg(npointsv(cr)) )
1908
1909 allocate ( rcontact(cr) % Jrg(npointsr(cr)) )
1910 allocate ( ucontact(cr) % Jrg(npointsu(cr)) )
1911 allocate ( vcontact(cr) % Jrg(npointsv(cr)) )
1912
1913 allocate ( rcontact(cr) % Idg(npointsr(cr)) )
1914 allocate ( ucontact(cr) % Idg(npointsu(cr)) )
1915 allocate ( vcontact(cr) % Idg(npointsv(cr)) )
1916
1917 allocate ( rcontact(cr) % Jdg(npointsr(cr)) )
1918 allocate ( ucontact(cr) % Jdg(npointsu(cr)) )
1919 allocate ( vcontact(cr) % Jdg(npointsv(cr)) )
1920
1921 dmem(dg)=dmem(dg)+4.0_r8*real(npointsr(cr),r8)
1922 dmem(dg)=dmem(dg)+4.0_r8*real(npointsu(cr),r8)
1923 dmem(dg)=dmem(dg)+4.0_r8*real(npointsv(cr),r8)
1924
1925# ifdef SOLVE3D
1926 allocate ( rcontact(cr) % Kdg(n(dg),npointsr(cr)) )
1927 allocate ( ucontact(cr) % Kdg(n(dg),npointsu(cr)) )
1928 allocate ( vcontact(cr) % Kdg(n(dg),npointsv(cr)) )
1929
1930 dmem(dg)=dmem(dg)+real(n(dg)*npointsr(cr),r8)
1931 dmem(dg)=dmem(dg)+real(n(dg)*npointsu(cr),r8)
1932 dmem(dg)=dmem(dg)+real(n(dg)*npointsv(cr),r8)
1933# endif
1934
1935 allocate ( rcontact(cr) % Lweight(4,npointsr(cr)) )
1936 allocate ( ucontact(cr) % Lweight(4,npointsu(cr)) )
1937 allocate ( vcontact(cr) % Lweight(4,npointsv(cr)) )
1938
1939 dmem(dg)=dmem(dg)+4.0_r8*real(npointsr(cr),r8)
1940 dmem(dg)=dmem(dg)+4.0_r8*real(npointsu(cr),r8)
1941 dmem(dg)=dmem(dg)+4.0_r8*real(npointsv(cr),r8)
1942
1943# ifdef WET_DRY
1944 allocate ( rcontact(cr) % LweightUnmasked(4,npointsr(cr)) )
1945 allocate ( ucontact(cr) % LweightUnmasked(4,npointsu(cr)) )
1946 allocate ( vcontact(cr) % LweightUnmasked(4,npointsv(cr)) )
1947
1948 dmem(dg)=dmem(dg)+4.0_r8*real(npointsr(cr),r8)
1949 dmem(dg)=dmem(dg)+4.0_r8*real(npointsu(cr),r8)
1950 dmem(dg)=dmem(dg)+4.0_r8*real(npointsv(cr),r8)
1951# endif
1952
1953# ifdef QUADRATIC_WEIGHTS
1954 allocate ( rcontact(cr) % Qweight(9,npointsr(cr)) )
1955 allocate ( ucontact(cr) % Qweight(9,npointsu(cr)) )
1956 allocate ( vcontact(cr) % Qweight(9,npointsv(cr)) )
1957
1958 dmem(dg)=dmem(dg)+9.0_r8*real(npointsr(cr),r8)
1959 dmem(dg)=dmem(dg)+9.0_r8*real(npointsu(cr),r8)
1960 dmem(dg)=dmem(dg)+9.0_r8*real(npointsv(cr),r8)
1961
1962# ifdef WET_DRY
1963 allocate ( rcontact(cr) % QweightUnmasked(9,npointsr(cr)) )
1964 allocate ( ucontact(cr) % QweightUnmasked(9,npointsu(cr)) )
1965 allocate ( vcontact(cr) % QweightUnmasked(9,npointsv(cr)) )
1966
1967 dmem(dg)=dmem(dg)+9.0_r8*real(npointsr(cr),r8)
1968 dmem(dg)=dmem(dg)+9.0_r8*real(npointsu(cr),r8)
1969 dmem(dg)=dmem(dg)+9.0_r8*real(npointsv(cr),r8)
1970# endif
1971# endif
1972
1973# ifdef SOLVE3D
1974 allocate ( rcontact(cr) % Vweight(2,n(dg),npointsr(cr)) )
1975 allocate ( ucontact(cr) % Vweight(2,n(dg),npointsu(cr)) )
1976 allocate ( vcontact(cr) % Vweight(2,n(dg),npointsv(cr)) )
1977
1978 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsr(cr),r8)
1979 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsu(cr),r8)
1980 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsv(cr),r8)
1981
1982# if defined TANGENT || defined TL_IOMS
1983 allocate ( rcontact(cr) % tl_Vweight(2,n(dg),npointsr(cr)) )
1984 allocate ( ucontact(cr) % tl_Vweight(2,n(dg),npointsu(cr)) )
1985 allocate ( vcontact(cr) % tl_Vweight(2,n(dg),npointsv(cr)) )
1986
1987 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsr(cr),r8)
1988 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsu(cr),r8)
1989 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsv(cr),r8)
1990# endif
1991
1992# ifdef ADJOINT
1993 allocate ( rcontact(cr) % ad_Vweight(2,n(dg),npointsr(cr)) )
1994 allocate ( ucontact(cr) % ad_Vweight(2,n(dg),npointsu(cr)) )
1995 allocate ( vcontact(cr) % ad_Vweight(2,n(dg),npointsv(cr)) )
1996
1997 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsr(cr),r8)
1998 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsu(cr),r8)
1999 dmem(dg)=dmem(dg)+2.0_r8*real(n(dg)*npointsv(cr),r8)
2000# endif
2001# endif
2002 END DO
2003!
2004! Initialize grid connectivity structure.
2005!
2006 DO cr=1,ncontact
2007 dg=donor_grid(cr)
2008 rg=receiver_grid(cr)
2009
2010 rcontact(cr) % coincident = lcoincident(rg)
2011 ucontact(cr) % coincident = lcoincident(rg)
2012 vcontact(cr) % coincident = lcoincident(rg)
2013
2014 rcontact(cr) % interpolate = linterpolate(rg)
2015 ucontact(cr) % interpolate = linterpolate(rg)
2016 vcontact(cr) % interpolate = linterpolate(rg)
2017
2018 rcontact(cr) % donor_grid = dg
2019 ucontact(cr) % donor_grid = dg
2020 vcontact(cr) % donor_grid = dg
2021
2022 rcontact(cr) % receiver_grid = rg
2023 ucontact(cr) % receiver_grid = rg
2024 vcontact(cr) % receiver_grid = rg
2025
2026 rcontact(cr) % Npoints = npointsr(cr)
2027 ucontact(cr) % Npoints = npointsu(cr)
2028 vcontact(cr) % Npoints = npointsv(cr)
2029
2030 DO m=1,npointsr(cr)
2031 ip=m+nstrr(cr)-1
2032 rcontact(cr) % Irg(m) = irg_cp(ip)
2033 rcontact(cr) % Jrg(m) = jrg_cp(ip)
2034 rcontact(cr) % Idg(m) = idg_cp(ip)
2035 rcontact(cr) % Jdg(m) = jdg_cp(ip)
2036# ifdef WET_DRY
2037 rcontact(cr) % LweightUnmasked(1,m) = lweight(1,ip)
2038 rcontact(cr) % LweightUnmasked(2,m) = lweight(2,ip)
2039 rcontact(cr) % LweightUnmasked(3,m) = lweight(3,ip)
2040 rcontact(cr) % LweightUnmasked(4,m) = lweight(4,ip)
2041# endif
2042 rcontact(cr) % Lweight(1,m) = lweight(1,ip)
2043 rcontact(cr) % Lweight(2,m) = lweight(2,ip)
2044 rcontact(cr) % Lweight(3,m) = lweight(3,ip)
2045 rcontact(cr) % Lweight(4,m) = lweight(4,ip)
2046# ifdef QUADRATIC_WEIGHTS
2047# ifdef WET_DRY
2048 rcontact(cr) % QweightUnmasked(1,m) = qweight(1,ip)
2049 rcontact(cr) % QweightUnmasked(2,m) = qweight(2,ip)
2050 rcontact(cr) % QweightUnmasked(3,m) = qweight(3,ip)
2051 rcontact(cr) % QweightUnmasked(4,m) = qweight(4,ip)
2052 rcontact(cr) % QweightUnmasked(5,m) = qweight(5,ip)
2053 rcontact(cr) % QweightUnmasked(6,m) = qweight(6,ip)
2054 rcontact(cr) % QweightUnmasked(7,m) = qweight(7,ip)
2055 rcontact(cr) % QweightUnmasked(8,m) = qweight(8,ip)
2056 rcontact(cr) % QweightUnmasked(9,m) = qweight(9,ip)
2057# endif
2058 rcontact(cr) % Qweight(1,m) = qweight(1,ip)
2059 rcontact(cr) % Qweight(2,m) = qweight(2,ip)
2060 rcontact(cr) % Qweight(3,m) = qweight(3,ip)
2061 rcontact(cr) % Qweight(4,m) = qweight(4,ip)
2062 rcontact(cr) % Qweight(5,m) = qweight(5,ip)
2063 rcontact(cr) % Qweight(6,m) = qweight(6,ip)
2064 rcontact(cr) % Qweight(7,m) = qweight(7,ip)
2065 rcontact(cr) % Qweight(8,m) = qweight(8,ip)
2066 rcontact(cr) % Qweight(9,m) = qweight(9,ip)
2067# endif
2068 END DO
2069
2070 DO m=1,npointsu(cr)
2071 ip=m+nstru(cr)-1
2072 ucontact(cr) % Irg(m) = irg_cp(ip)
2073 ucontact(cr) % Jrg(m) = jrg_cp(ip)
2074 ucontact(cr) % Idg(m) = idg_cp(ip)
2075 ucontact(cr) % Jdg(m) = jdg_cp(ip)
2076# ifdef WET_DRY
2077 ucontact(cr) % LweightUnmasked(1,m) = lweight(1,ip)
2078 ucontact(cr) % LweightUnmasked(2,m) = lweight(2,ip)
2079 ucontact(cr) % LweightUnmasked(3,m) = lweight(3,ip)
2080 ucontact(cr) % LweightUnmasked(4,m) = lweight(4,ip)
2081# endif
2082 ucontact(cr) % Lweight(1,m) = lweight(1,ip)
2083 ucontact(cr) % Lweight(2,m) = lweight(2,ip)
2084 ucontact(cr) % Lweight(3,m) = lweight(3,ip)
2085 ucontact(cr) % Lweight(4,m) = lweight(4,ip)
2086# ifdef QUADRATIC_WEIGHTS
2087# ifdef WET_DRY
2088 ucontact(cr) % QweightUnmasked(1,m) = qweight(1,ip)
2089 ucontact(cr) % QweightUnmasked(2,m) = qweight(2,ip)
2090 ucontact(cr) % QweightUnmasked(3,m) = qweight(3,ip)
2091 ucontact(cr) % QweightUnmasked(4,m) = qweight(4,ip)
2092 ucontact(cr) % QweightUnmasked(5,m) = qweight(5,ip)
2093 ucontact(cr) % QweightUnmasked(6,m) = qweight(6,ip)
2094 ucontact(cr) % QweightUnmasked(7,m) = qweight(7,ip)
2095 ucontact(cr) % QweightUnmasked(8,m) = qweight(8,ip)
2096 ucontact(cr) % QweightUnmasked(9,m) = qweight(9,ip)
2097# endif
2098 ucontact(cr) % Qweight(1,m) = qweight(1,ip)
2099 ucontact(cr) % Qweight(2,m) = qweight(2,ip)
2100 ucontact(cr) % Qweight(3,m) = qweight(3,ip)
2101 ucontact(cr) % Qweight(4,m) = qweight(4,ip)
2102 ucontact(cr) % Qweight(5,m) = qweight(5,ip)
2103 ucontact(cr) % Qweight(6,m) = qweight(6,ip)
2104 ucontact(cr) % Qweight(7,m) = qweight(7,ip)
2105 ucontact(cr) % Qweight(8,m) = qweight(8,ip)
2106 ucontact(cr) % Qweight(9,m) = qweight(9,ip)
2107# endif
2108 END DO
2109
2110 DO m=1,npointsv(cr)
2111 ip=m+nstrv(cr)-1
2112 vcontact(cr) % Irg(m) = irg_cp(ip)
2113 vcontact(cr) % Jrg(m) = jrg_cp(ip)
2114 vcontact(cr) % Idg(m) = idg_cp(ip)
2115 vcontact(cr) % Jdg(m) = jdg_cp(ip)
2116# ifdef WET_DRY
2117 vcontact(cr) % LweightUnmasked(1,m) = lweight(1,ip)
2118 vcontact(cr) % LweightUnmasked(2,m) = lweight(2,ip)
2119 vcontact(cr) % LweightUnmasked(3,m) = lweight(3,ip)
2120 vcontact(cr) % LweightUnmasked(4,m) = lweight(4,ip)
2121# endif
2122 vcontact(cr) % Lweight(1,m) = lweight(1,ip)
2123 vcontact(cr) % Lweight(2,m) = lweight(2,ip)
2124 vcontact(cr) % Lweight(3,m) = lweight(3,ip)
2125 vcontact(cr) % Lweight(4,m) = lweight(4,ip)
2126# ifdef QUADRATIC_WEIGHTS
2127# ifdef WET_DRY
2128 vcontact(cr) % QweightUnmasked(1,m) = qweight(1,ip)
2129 vcontact(cr) % QweightUnmasked(2,m) = qweight(2,ip)
2130 vcontact(cr) % QweightUnmasked(3,m) = qweight(3,ip)
2131 vcontact(cr) % QweightUnmasked(4,m) = qweight(4,ip)
2132 vcontact(cr) % QweightUnmasked(5,m) = qweight(5,ip)
2133 vcontact(cr) % QweightUnmasked(6,m) = qweight(6,ip)
2134 vcontact(cr) % QweightUnmasked(7,m) = qweight(7,ip)
2135 vcontact(cr) % QweightUnmasked(8,m) = qweight(8,ip)
2136 vcontact(cr) % QweightUnmasked(9,m) = qweight(9,ip)
2137# endif
2138 vcontact(cr) % Qweight(1,m) = qweight(1,ip)
2139 vcontact(cr) % Qweight(2,m) = qweight(2,ip)
2140 vcontact(cr) % Qweight(3,m) = qweight(3,ip)
2141 vcontact(cr) % Qweight(4,m) = qweight(4,ip)
2142 vcontact(cr) % Qweight(5,m) = qweight(5,ip)
2143 vcontact(cr) % Qweight(6,m) = qweight(6,ip)
2144 vcontact(cr) % Qweight(7,m) = qweight(7,ip)
2145 vcontact(cr) % Qweight(8,m) = qweight(8,ip)
2146 vcontact(cr) % Qweight(9,m) = qweight(9,ip)
2147# endif
2148 END DO
2149
2150 END DO
2151!
2152!-----------------------------------------------------------------------
2153! Load contact points grid metrics.
2154!-----------------------------------------------------------------------
2155!
2156! Allocate contact region metrics (type T_NGM) structure.
2157!
2158 allocate ( contact_metric(ncontact) )
2159!
2160! Allocate arrays in contact region metrics structure.
2161!
2162 DO cr=1,ncontact
2163 allocate ( contact_metric(cr) % angler(npointsr(cr)) )
2164 allocate ( contact_metric(cr) % dndx (npointsr(cr)) )
2165 allocate ( contact_metric(cr) % dmde (npointsr(cr)) )
2166 allocate ( contact_metric(cr) % f (npointsr(cr)) )
2167 allocate ( contact_metric(cr) % h (npointsr(cr)) )
2168 allocate ( contact_metric(cr) % rmask (npointsr(cr)) )
2169 allocate ( contact_metric(cr) % umask (npointsu(cr)) )
2170 allocate ( contact_metric(cr) % vmask (npointsv(cr)) )
2171 allocate ( contact_metric(cr) % pm (npointsr(cr)) )
2172 allocate ( contact_metric(cr) % pn (npointsr(cr)) )
2173 allocate ( contact_metric(cr) % Xr (npointsr(cr)) )
2174 allocate ( contact_metric(cr) % Yr (npointsr(cr)) )
2175 allocate ( contact_metric(cr) % Xu (npointsu(cr)) )
2176 allocate ( contact_metric(cr) % Yu (npointsu(cr)) )
2177 allocate ( contact_metric(cr) % Xv (npointsv(cr)) )
2178 allocate ( contact_metric(cr) % Yv (npointsv(cr)) )
2179
2180 dmem(ng)=dmem(ng)+10.0_r8*real(npointsr(cr),r8)
2181 dmem(ng)=dmem(ng)+ 3.0_r8*real(npointsu(cr),r8)
2182 dmem(ng)=dmem(ng)+ 3.0_r8*real(npointsv(cr),r8)
2183 END DO
2184!
2185! Initialize contact region metrics structure.
2186!
2187 DO cr=1,ncontact
2188 DO m=1,npointsr(cr)
2189 ip=m+nstrr(cr)-1
2190 contact_metric(cr) % angler(m) = angle(ip)
2191 contact_metric(cr) % dndx (m) = dndx(ip)
2192 contact_metric(cr) % dmde (m) = dmde(ip)
2193 contact_metric(cr) % f (m) = f(ip)
2194 contact_metric(cr) % h (m) = h(ip)
2195 contact_metric(cr) % rmask (m) = mask(ip)
2196 contact_metric(cr) % pm (m) = pm(ip)
2197 contact_metric(cr) % pn (m) = pn(ip)
2198 contact_metric(cr) % Xr (m) = xrg(ip)
2199 contact_metric(cr) % Yr (m) = yrg(ip)
2200 END DO
2201
2202 DO m=1,npointsu(cr)
2203 ip=m+nstru(cr)-1
2204 contact_metric(cr) % umask(m) = mask(ip)
2205 contact_metric(cr) % Xu (m) = xrg(ip)
2206 contact_metric(cr) % Yu (m) = yrg(ip)
2207 END DO
2208
2209 DO m=1,npointsv(cr)
2210 ip=m+nstrv(cr)-1
2211 contact_metric(cr) % vmask(m) = mask(ip)
2212 contact_metric(cr) % Xv (m) = xrg(ip)
2213 contact_metric(cr) % Yv (m) = yrg(ip)
2214 END DO
2215 END DO
2216!
2217!-----------------------------------------------------------------------
2218! Initialize various parameters.
2219!-----------------------------------------------------------------------
2220!
2221! Nested grids conectivity switches. They are used in "get_bounds.F"
2222! to determenine if the global state array need to be extended by few
2223! extra points at each boundary to accomodate the contact regions and
2224! contact points. These are VERY important switches in ROMS nesting.
2225!
2226 IF (.not.allocated(contactregion)) THEN
2227 allocate ( contactregion(4,ngrids) )
2228 contactregion = .false.
2229 dmem(ng)=dmem(ng)+4.0_r8*real(ngrids,r8)
2230 END IF
2231!
2232 DO m=1,ncdatum
2233 cr=contact_region(m)
2234 rg=receiver_grid(cr)
2235 ibry=on_boundary(m)
2236 IF ((ibry.eq.iwest ).or.(ibry.eq.ieast ).or. &
2237 & (ibry.eq.isouth).or.(ibry.eq.inorth)) THEN
2238 IF (.not.contactregion(ibry,rg)) THEN
2239 contactregion(ibry,rg)=.true.
2240 END IF
2241 END IF
2242 END DO
2243!
2244! Set composite and refinement grids switches.
2245!
2246 DO cr=1,ncontact
2247 rg=receiver_grid(cr)
2248 DO ibry=1,4
2249 IF (lcomposite(cr).and.contactregion(ibry,rg)) THEN
2250 compositegrid(ibry,rg)=.true.
2251 END IF
2252 END DO
2253 refinedgrid(rg)=lrefinement(cr)
2254 END DO
2255!
2256! Set the switch to compute vertical interpolation weights. Currently,
2257! they are only needed in non-coincident composite grids.
2258!
2259 IF (.not.any(lcoincident).and.any(lcomposite)) THEN
2260 get_vweights=.true.
2261 ELSE
2262 get_vweights=.false.
2263 END IF
2264!
2265! Refinement grid scale factor.
2266!
2267 DO ig=1,ngrids
2268 refinescale(ig)=refine_factor(ig)
2269 END DO
2270!
2271! Coarser donor grid number to finer grid contact points. Their
2272! values are computed in "metrics.F" when the maximum grid spacing
2273! is computed. The maximum grid spacing is used to distinguish
2274! contact areas for more than two layers of nesting.
2275!
2276!$OMP PARALLEL
2277 IF (.not.allocated(coarserdonor)) THEN
2278 allocate ( coarserdonor(ngrids) )
2279 coarserdonor = 0
2280 dmem(ng)=dmem(ng)+real(ngrids,r8)
2281 END IF
2282!$OMP END PARALLEL
2283!
2284! Finer donor grid to coarser receiver grid for two-way feedback.
2285! Their values are computed in "metrics.F" when the maximum grid
2286! spacing is computed.
2287!
2288!$OMP PARALLEL
2289 IF (.not.allocated(finerdonor)) THEN
2290 allocate ( finerdonor(ngrids) )
2291 finerdonor = 0
2292 dmem(ng)=dmem(ng)+real(ngrids,r8)
2293 END IF
2294!$OMP END PARALLEL
2295!
2296! Logical switch indicating which coarser grid is a donor to a finer
2297! grid external contact points. Their values are computed in
2298! "metrics.F" after "CoarserDonor" is computed.
2299!
2300!$OMP PARALLEL
2301 IF (.not.allocated(donortofiner)) THEN
2302 allocate ( donortofiner(ngrids) )
2303 donortofiner = .false.
2304 dmem(ng)=dmem(ng)+real(ngrids,r8)
2305 END IF
2306!$OMP END PARALLEL
2307!
2308! Number of refined grid time-steps. Their values are initialized to
2309! the optimal values based on the spatial refinement ratio. These
2310! values are over-written in "metrics.F" to use the time-step
2311! size ratio between donor and receiver grid. The user is allowed
2312! to take larger divisible time-step with respect to the donor grid.
2313! The user is responsible to set the appropriate refined grid time-step
2314! for stability.
2315!
2316!$OMP PARALLEL
2317 IF (.not.allocated(refinesteps)) THEN
2318 allocate ( refinesteps(ngrids) )
2319 dmem(ng)=dmem(ng)+real(ngrids,r8)
2320 DO ig=1,ngrids
2321 refinesteps(ig)=refine_factor(ig)
2322 END DO
2323 END IF
2324!$OMP END PARALLEL
2325!
2326! Refine grid time-steps counter with respect the coarse grid (ng=1)
2327! single time-step.
2328!
2329 IF (.not.allocated(refinestepscounter)) THEN
2330 allocate ( refinestepscounter(ngrids) )
2331 refinestepscounter = 0
2332 dmem(ng)=dmem(ng)+real(ngrids,r8)
2333 END IF
2334!
2335! Interval used in the two-way exchange between fine and coarse grids.
2336!
2337 IF (.not.allocated(twowayinterval)) THEN
2338 allocate ( twowayinterval(ngrids) )
2339 twowayinterval = 0.0_r8
2340 dmem(ng)=dmem(ng)+real(ngrids,r8)
2341 END IF
2342!
2343! Switch indicating which refined grid(s) include finer refined grids:
2344! telescoping refinement.
2345!
2346!$OMP PARALLEL
2347 IF (.not.allocated(telescoping)) THEN
2348 allocate ( telescoping(ngrids) )
2349 telescoping = .false.
2350 dmem(ng)=dmem(ng)+real(ngrids,r8)
2351 END IF
2352!$OMP END PARALLEL
2353!
2354! Rolling index and time (seconds) used in the temporal interpolation
2355! of contact point data.
2356!
2357!$OMP PARALLEL
2358 IF (.not.allocated(rollingindex)) THEN
2359 allocate ( rollingindex(ncontact) )
2360 rollingindex = 0
2361 dmem(ng)=dmem(ng)+real(ncontact,r8)
2362 END IF
2363
2364 IF (.not.allocated(rollingtime)) THEN
2365 allocate ( rollingtime(2,ncontact) )
2366 rollingtime = 0
2367 dmem(ng)=dmem(ng)+2.0_r8*real(ncontact,r8)
2368 END IF
2369!$OMP END PARALLEL
2370!
2371!-----------------------------------------------------------------------
2372! Allocate composite grids contact regions structure and arrays.
2373!-----------------------------------------------------------------------
2374!
2375 IF (any(compositegrid)) THEN
2376!
2377! Allocate composite grid contact regions (type T_COMPOSITE) structure.
2378!
2379 allocate ( composite(ncontact) )
2380!
2381! Allocate arrays in composite grids contact regions structure.
2382!
2383 DO cr=1,ncontact
2384 dg=donor_grid(cr)
2385
2386 allocate ( composite(cr) % bustr(4,npointsu(cr)) )
2387 allocate ( composite(cr) % bvstr(4,npointsv(cr)) )
2388
2389 allocate ( composite(cr) % ubar(4,npointsu(cr),2) )
2390 allocate ( composite(cr) % vbar(4,npointsv(cr),2) )
2391 allocate ( composite(cr) % zeta(4,npointsr(cr),2) )
2392
2393 allocate ( composite(cr) % rzeta(4,npointsr(cr)) )
2394
2395 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsr(cr),r8)
2396 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsu(cr),r8)
2397 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsv(cr),r8)
2398
2399# if defined TANGENT || defined TL_IOMS
2400 allocate ( composite(cr) % tl_bustr(4,npointsu(cr)) )
2401 allocate ( composite(cr) % tl_bvstr(4,npointsv(cr)) )
2402
2403 allocate ( composite(cr) % tl_ubar(4,npointsu(cr),2) )
2404 allocate ( composite(cr) % tl_vbar(4,npointsv(cr),2) )
2405 allocate ( composite(cr) % tl_zeta(4,npointsr(cr),2) )
2406
2407 allocate ( composite(cr) % tl_rzeta(4,npointsr(cr)) )
2408
2409 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsr(cr),r8)
2410 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsu(cr),r8)
2411 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsv(cr),r8)
2412# endif
2413
2414# ifdef ADJOINT
2415 allocate ( composite(cr) % ad_bustr(4,npointsu(cr)) )
2416 allocate ( composite(cr) % ad_bvstr(4,npointsv(cr)) )
2417
2418 allocate ( composite(cr) % ad_ubar(4,npointsu(cr),2) )
2419 allocate ( composite(cr) % ad_vbar(4,npointsv(cr),2) )
2420 allocate ( composite(cr) % ad_zeta(4,npointsr(cr),2) )
2421
2422 allocate ( composite(cr) % ad_rzeta(4,npointsr(cr)) )
2423
2424 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsr(cr),r8)
2425 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsu(cr),r8)
2426 dmem(dg)=dmem(dg)+2.0_r8*real(4*npointsv(cr),r8)
2427# endif
2428
2429# ifdef SOLVE3D
2430 allocate ( composite(cr) % DU_avg1(4,npointsu(cr)) )
2431 allocate ( composite(cr) % DV_avg1(4,npointsv(cr)) )
2432 allocate ( composite(cr) % Zt_avg1(4,npointsr(cr)) )
2433
2434 dmem(dg)=dmem(dg)+real(4*npointsr(cr),r8)
2435 dmem(dg)=dmem(dg)+real(4*npointsu(cr),r8)
2436 dmem(dg)=dmem(dg)+real(4*npointsv(cr),r8)
2437
2438 allocate ( composite(cr) % u(4,n(dg),npointsu(cr)) )
2439 allocate ( composite(cr) % v(4,n(dg),npointsv(cr)) )
2440
2441 allocate ( composite(cr) % Huon(4,n(dg),npointsu(cr)) )
2442 allocate ( composite(cr) % Hvom(4,n(dg),npointsv(cr)) )
2443
2444 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsu(cr),r8)
2445 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsv(cr),r8)
2446
2447 allocate ( composite(cr) % t(4,n(dg),npointsr(cr),nt(dg)) )
2448
2449 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsu(cr)*nt(dg),r8)
2450
2451# if defined TANGENT || defined TL_IOMS
2452 allocate ( composite(cr) % tl_DU_avg1(4,npointsu(cr)) )
2453 allocate ( composite(cr) % tl_DV_avg1(4,npointsv(cr)) )
2454 allocate ( composite(cr) % tl_Zt_avg1(4,npointsr(cr)) )
2455
2456 dmem(dg)=dmem(dg)+real(4*npointsr(cr),r8)
2457 dmem(dg)=dmem(dg)+real(4*npointsu(cr),r8)
2458 dmem(dg)=dmem(dg)+real(4*npointsv(cr),r8)
2459
2460 allocate ( composite(cr) % tl_u(4,n(dg),npointsu(cr)) )
2461 allocate ( composite(cr) % tl_v(4,n(dg),npointsv(cr)) )
2462
2463 allocate ( composite(cr) % tl_Huon(4,n(dg),npointsu(cr)) )
2464 allocate ( composite(cr) % tl_Hvom(4,n(dg),npointsv(cr)) )
2465
2466 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsu(cr),r8)
2467 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsv(cr),r8)
2468
2469 allocate ( composite(cr) % tl_t(4,n(dg),npointsr(cr),nt(dg)) )
2470
2471 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsu(cr)*nt(dg),r8)
2472# endif
2473
2474# ifdef ADJOINT
2475 allocate ( composite(cr) % ad_DU_avg1(4,npointsu(cr)) )
2476 allocate ( composite(cr) % ad_DV_avg1(4,npointsv(cr)) )
2477 allocate ( composite(cr) % ad_Zt_avg1(4,npointsr(cr)) )
2478
2479 dmem(dg)=dmem(dg)+real(4*npointsr(cr),r8)
2480 dmem(dg)=dmem(dg)+real(4*npointsu(cr),r8)
2481 dmem(dg)=dmem(dg)+real(4*npointsv(cr),r8)
2482
2483 allocate ( composite(cr) % ad_u(4,n(dg),npointsu(cr)) )
2484 allocate ( composite(cr) % ad_v(4,n(dg),npointsv(cr)) )
2485
2486 allocate ( composite(cr) % ad_Huon(4,n(dg),npointsu(cr)) )
2487 allocate ( composite(cr) % ad_Hvom(4,n(dg),npointsv(cr)) )
2488
2489 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsu(cr),r8)
2490 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsv(cr),r8)
2491
2492 allocate ( composite(cr) % ad_t(4,n(dg),npointsr(cr),nt(dg)) )
2493
2494 dmem(dg)=dmem(dg)+2.0_r8*real(4*n(dg)*npointsu(cr)*nt(dg),r8)
2495# endif
2496# endif
2497 END DO
2498 END IF
2499!
2500!-----------------------------------------------------------------------
2501! Allocate refinement grids contact regions structure.
2502!-----------------------------------------------------------------------
2503!
2504 IF (any(refinedgrid)) THEN
2505!
2506! Allocate refinement grids contact region (type T_REFINED) structure.
2507!
2508 allocate ( refined(ncontact) )
2509!
2510! Allocate arrays in refinement grids contact region structure.
2511!
2512 DO cr=1,ncontact
2513 rg=receiver_grid(cr)
2514
2515 allocate ( refined(cr) % ubar(4,npointsu(cr),2) )
2516 allocate ( refined(cr) % vbar(4,npointsv(cr),2) )
2517 allocate ( refined(cr) % zeta(4,npointsr(cr),2) )
2518
2519 allocate ( refined(cr) % U2d_flux(4,npointsu(cr),2) )
2520 allocate ( refined(cr) % V2d_flux(4,npointsv(cr),2) )
2521
2522 allocate ( refined(cr) % on_u(npointsu(cr)) )
2523 allocate ( refined(cr) % om_v(npointsv(cr)) )
2524
2525 dmem(rg)=dmem(rg)+real(4*npointsr(cr),r8)
2526 dmem(rg)=dmem(rg)+3.0_r8*real(4*npointsu(cr),r8)
2527 dmem(rg)=dmem(rg)+3.0_r8*real(4*npointsv(cr),r8)
2528
2529# if defined TANGENT || defined TL_IOMS
2530 allocate ( refined(cr) % tl_ubar(4,npointsu(cr),2) )
2531 allocate ( refined(cr) % tl_vbar(4,npointsv(cr),2) )
2532 allocate ( refined(cr) % tl_zeta(4,npointsr(cr),2) )
2533
2534 allocate ( refined(cr) % tl_U2d_flux(4,npointsu(cr),2) )
2535 allocate ( refined(cr) % tl_V2d_flux(4,npointsv(cr),2) )
2536
2537 dmem(rg)=dmem(rg)+real(4*npointsr(cr),r8)
2538 dmem(rg)=dmem(rg)+3.0_r8*real(4*npointsu(cr),r8)
2539 dmem(rg)=dmem(rg)+3.0_r8*real(4*npointsv(cr),r8)
2540# endif
2541
2542# ifdef ADJOINT
2543 allocate ( refined(cr) % ad_ubar(4,npointsu(cr),2) )
2544 allocate ( refined(cr) % ad_vbar(4,npointsv(cr),2) )
2545 allocate ( refined(cr) % ad_zeta(4,npointsr(cr),2) )
2546
2547 allocate ( refined(cr) % ad_U2d_flux(4,npointsu(cr),2) )
2548 allocate ( refined(cr) % ad_V2d_flux(4,npointsv(cr),2) )
2549
2550 dmem(rg)=dmem(rg)+real(4*npointsr(cr),r8)
2551 dmem(rg)=dmem(rg)+3.0_r8*real(4*npointsu(cr),r8)
2552 dmem(rg)=dmem(rg)+3.0_r8*real(4*npointsv(cr),r8)
2553# endif
2554
2555# ifdef SOLVE3D
2556 allocate ( refined(cr) % u(4,n(rg),npointsu(cr),2) )
2557 allocate ( refined(cr) % v(4,n(rg),npointsv(cr),2) )
2558
2559 dmem(rg)=dmem(rg)+3.0_r8*real(4*n(rg)*npointsu(cr),r8)
2560 dmem(rg)=dmem(rg)+3.0_r8*real(4*n(rg)*npointsv(cr),r8)
2561
2562 allocate ( refined(cr) % t(4,n(rg),npointsr(cr),2,nt(rg)) )
2563
2564 dmem(rg)=dmem(rg)+2.0_r8*real(4*n(rg)*npointsr(cr)*nt(rg),r8)
2565
2566# if defined TANGENT || defined TL_IOMS
2567 allocate ( refined(cr) % tl_u(4,n(rg),npointsu(cr),2) )
2568 allocate ( refined(cr) % tl_v(4,n(rg),npointsv(cr),2) )
2569
2570 dmem(rg)=dmem(rg)+3.0_r8*real(4*n(rg)*npointsu(cr),r8)
2571 dmem(rg)=dmem(rg)+3.0_r8*real(4*n(rg)*npointsv(cr),r8)
2572
2573 allocate ( refined(cr) % tl_t(4,n(rg),npointsr(cr),2,nt(rg)) )
2574
2575 dmem(rg)=dmem(rg)+2.0_r8*real(4*n(rg)*npointsr(cr)*nt(rg),r8)
2576# endif
2577
2578# ifdef ADJOINT
2579 allocate ( refined(cr) % ad_u(4,n(rg),npointsu(cr),2) )
2580 allocate ( refined(cr) % ad_v(4,n(rg),npointsv(cr),2) )
2581
2582 dmem(rg)=dmem(rg)+3.0_r8*real(4*n(rg)*npointsu(cr),r8)
2583 dmem(rg)=dmem(rg)+3.0_r8*real(4*n(rg)*npointsv(cr),r8)
2584
2585 allocate ( refined(cr) % ad_t(4,n(rg),npointsr(cr),2,nt(rg)) )
2586
2587 dmem(rg)=dmem(rg)+2.0_r8*real(4*n(rg)*npointsr(cr)*nt(rg),r8)
2588# endif
2589# endif
2590 END DO
2591 END IF
2592!
2593 RETURN

References mod_nesting::coarserdonor, mod_nesting::composite, mod_scalars::compositegrid, mod_nesting::contact_metric, mod_nesting::contact_region, mod_nesting::contactregion, mod_param::dmem, mod_nesting::donor_grid, mod_nesting::donortofiner, mod_scalars::exit_flag, mod_nesting::finerdonor, strings_mod::founderror(), mod_nesting::get_vweights, mod_nesting::i_left, mod_nesting::i_right, mod_nesting::idg_cp, mod_scalars::ieast, mod_scalars::inorth, mod_nesting::irg_cp, mod_scalars::isouth, mod_scalars::iwest, mod_nesting::j_bottom, mod_nesting::j_top, mod_nesting::jdg_cp, mod_nesting::jrg_cp, mod_param::lm, mod_parallel::master, mod_param::mm, mod_param::n, mod_nesting::ncdatum, mod_nesting::ncontact, mod_nesting::ncpoints, mod_nesting::nendr, mod_nesting::nendu, mod_nesting::nendv, mod_iounits::ngcname, mod_param::ngrids, mod_scalars::noerror, mod_nesting::nstrr, mod_nesting::nstru, mod_nesting::nstrv, mod_param::nt, mod_nesting::on_boundary, mod_pio_netcdf::pio_netcdf_close(), mod_pio_netcdf::pio_netcdf_get_dim(), mod_pio_netcdf::pio_netcdf_inq_var(), mod_pio_netcdf::pio_netcdf_open(), mod_nesting::rcontact, mod_nesting::receiver_grid, mod_nesting::refined, mod_scalars::refinedgrid, mod_scalars::refinescale, mod_nesting::refinesteps, mod_nesting::refinestepscounter, mod_nesting::rollingindex, mod_nesting::rollingtime, mod_iounits::sourcefile, mod_iounits::stdout, mod_nesting::telescoping, mod_nesting::twowayinterval, mod_nesting::ucontact, and mod_nesting::vcontact.

Referenced by set_contact().

Here is the call graph for this function:
Here is the caller graph for this function: