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

Data Types

type  t_bcp
 
type  t_composite
 
type  t_ngc
 
type  t_ngm
 
type  t_refined
 

Functions/Subroutines

subroutine, public allocate_nesting
 
subroutine, public deallocate_nesting
 
subroutine, public initialize_nesting
 

Variables

integer, parameter nmflx = -6
 
integer, parameter ndxdy = -5
 
integer, parameter ngetd = -4
 
integer, parameter nmask = -3
 
integer, parameter nputd = -2
 
integer, parameter n2way = -1
 
integer, parameter nfsic = 1
 
integer, parameter n2dic = 2
 
integer, parameter n3dic = 3
 
integer, parameter ntvic = 4
 
integer, parameter nbstr = 5
 
integer, parameter nrhst = 6
 
integer, parameter nzeta = 7
 
integer, parameter nzwgt = 8
 
integer, parameter n2dps = 9
 
integer, parameter n2dcs = 10
 
integer, parameter n2dfx = 11
 
integer, parameter n3duv = 12
 
integer, parameter n3dtv = 13
 
logical, dimension(:,:), allocatable contactregion
 
logical, dimension(:), allocatable donortofiner
 
logical, dimension(:), allocatable telescoping
 
logical get_vweights
 
integer, dimension(:), allocatable coarserdonor
 
integer, dimension(:), allocatable finerdonor
 
integer, dimension(:), allocatable refinesteps
 
integer, dimension(:), allocatable refinestepscounter
 
real(r8), dimension(:), allocatable twowayinterval
 
integer, dimension(:), allocatable donor_grid
 
integer, dimension(:), allocatable receiver_grid
 
integer, dimension(:), allocatable rollingindex
 
real(dp), dimension(:,:), allocatable rollingtime
 
integer, dimension(:), allocatable i_left
 
integer, dimension(:), allocatable i_right
 
integer, dimension(:), allocatable j_bottom
 
integer, dimension(:), allocatable j_top
 
integer ncdatum
 
integer, dimension(:), allocatable ncpoints
 
integer, dimension(:), allocatable nstrr
 
integer, dimension(:), allocatable nendr
 
integer, dimension(:), allocatable nstru
 
integer, dimension(:), allocatable nendu
 
integer, dimension(:), allocatable nstrv
 
integer, dimension(:), allocatable nendv
 
integer, dimension(:), allocatable contact_region
 
integer, dimension(:), allocatable on_boundary
 
integer, dimension(:), allocatable idg_cp
 
integer, dimension(:), allocatable jdg_cp
 
integer, dimension(:), allocatable irg_cp
 
integer, dimension(:), allocatable jrg_cp
 
integer ncontact
 
type(t_ngc), dimension(:), allocatable rcontact
 
type(t_ngc), dimension(:), allocatable ucontact
 
type(t_ngc), dimension(:), allocatable vcontact
 
type(t_bcp), dimension(:,:), allocatable bry_contact
 
type(t_ngm), dimension(:), allocatable contact_metric
 
type(t_composite), dimension(:), allocatable composite
 
type(t_refined), dimension(:), allocatable refined
 

Function/Subroutine Documentation

◆ allocate_nesting()

subroutine, public mod_nesting::allocate_nesting

Definition at line 441 of file mod_nesting.F.

442!
443!=======================================================================
444! !
445! This routine allocates and initializes nesting structure for 2D !
446! state variables. !
447! !
448!=======================================================================
449!
450 USE mod_param
451 USE mod_boundary
452 USE mod_scalars
453!
454! Local variable declarations.
455!
456 integer :: LBi, UBi, LBj, UBj
457 integer :: Imin, Imax, Jmin, Jmax
458 integer :: CCR, cr, dg, ng, rg
459 integer :: i, ibry, ic, id, ir, j, jd, jr, m, my_tile
460 integer :: ispval
461
462 integer, allocatable :: Ibmin(:,:), Ibmax(:,:)
463 integer, allocatable :: Jbmin(:,:), Jbmax(:,:)
464!
465!-----------------------------------------------------------------------
466! Unpack Boundary Contact Points structure (type T_BCP).
467!-----------------------------------------------------------------------
468!
469! Allocate boundary connectivity (type T_BCP) structure.
470!
471 allocate ( bry_contact(4,ncontact) )
472!
473! Allocate arrays in boundary connectivity structure.
474!
475 my_tile=-1 ! for global values
476 DO cr=1,ncontact
477 rg=receiver_grid(cr)
478 lbi=bounds(rg)%LBi(my_tile)
479 ubi=bounds(rg)%UBi(my_tile)
480 lbj=bounds(rg)%LBj(my_tile)
481 ubj=bounds(rg)%UBj(my_tile)
482 DO ibry=1,4
483 SELECT CASE (ibry)
484 CASE (iwest, ieast)
485 allocate ( bry_contact(ibry,cr) % Ib(lbj:ubj) )
486 dmem(rg)=dmem(rg)+real(ubj-lbj,r8)
487
488 allocate ( bry_contact(ibry,cr) % Jb(lbj:ubj) )
489 dmem(rg)=dmem(rg)+real(ubj-lbj,r8)
490
491 allocate ( bry_contact(ibry,cr) % C2Bindex(lbj:ubj) )
492 dmem(rg)=dmem(rg)+real(ubj-lbj,r8)
493
494# ifdef NESTING_DEBUG
495 allocate ( bry_contact(ibry,cr) % Mflux(lbj:ubj) )
496 dmem(rg)=dmem(rg)+real(ubj-lbj,r8)
497# endif
498
499# ifdef SOLVE3D
500 allocate ( bry_contact(ibry,cr) % Tflux(lbj:ubj, &
501 & n(rg),nt(rg)) )
502 dmem(rg)=dmem(rg)+real((ubj-lbj)*n(rg)*nt(rg),r8)
503# endif
504 CASE (isouth, inorth)
505 allocate ( bry_contact(ibry,cr) % Ib(lbi:ubi) )
506 dmem(rg)=dmem(rg)+real(ubi-lbi,r8)
507
508 allocate ( bry_contact(ibry,cr) % Jb(lbi:ubi) )
509 dmem(rg)=dmem(rg)+real(ubi-lbi,r8)
510
511 allocate ( bry_contact(ibry,cr) % C2Bindex(lbi:ubi) )
512 dmem(rg)=dmem(rg)+real(ubi-lbi,r8)
513
514# ifdef NESTING_DEBUG
515 allocate ( bry_contact(ibry,cr) % Mflux(lbi:ubi) )
516 dmem(rg)=dmem(rg)+real(ubi-lbi,r8)
517# endif
518
519# ifdef SOLVE3D
520 allocate ( bry_contact(ibry,cr) % Tflux(lbi:ubi, &
521 & n(rg),nt(rg)) )
522 dmem(rg)=dmem(rg)+real((ubi-lbi)*n(rg)*nt(rg),r8)
523# endif
524 END SELECT
525 END DO
526 END DO
527!
528! Initialize boundary connectivity structure: Boundary indices array
529! are initialized to its special value.
530!
531 ispval=-999
532!
533 IF (.not.allocated(ibmin)) THEN
534 allocate ( ibmin(4,ncontact) )
535 DO ng=1,ngrids
536 dmem(ng)=dmem(ng)+4.0_r8*real(ncontact,r8)
537 END DO
538 ibmin = -ispval
539 END IF
540 IF (.not.allocated(ibmax)) THEN
541 allocate ( ibmax(4,ncontact) )
542 DO ng=1,ngrids
543 dmem(ng)=dmem(ng)+4.0_r8*real(ncontact,r8)
544 END DO
545 ibmax = ispval
546 END IF
547 IF (.not.allocated(jbmin)) THEN
548 allocate ( jbmin(4,ncontact) )
549 DO ng=1,ngrids
550 dmem(ng)=dmem(ng)+4.0_r8*real(ncontact,r8)
551 END DO
552 jbmin = -ispval
553 END IF
554 IF (.not.allocated(jbmax)) THEN
555 allocate ( jbmax(4,ncontact) )
556 DO ng=1,ngrids
557 dmem(ng)=dmem(ng)+4.0_r8*real(ncontact,r8)
558 END DO
559 jbmax = ispval
560 END IF
561!
562 DO cr=1,ncontact
563 DO ibry=1,4
564 bry_contact(ibry,cr) % spv = ispval
565 bry_contact(ibry,cr) % Ib = ispval
566 bry_contact(ibry,cr) % Jb = ispval
567 bry_contact(ibry,cr) % C2Bindex = ispval
568# ifdef NESTING_DEBUG
569 bry_contact(ibry,cr) % Mflux = 0.0_r8
570# endif
571# ifdef SOLVE3D
572 bry_contact(ibry,cr) % Tflux = 0.0_r8
573# endif
574 END DO
575 END DO
576!
577! Identify contact points located on the grid boundary. Notice that
578! the conjugate contact region (CCR) is also processed but it is not
579! yet used. Also, the CCR indices (Ib,Jb) are in refinement overwriten
580! in the m-loop below because several finer grid contact points
581! (RefineScale) are contained in the coarser grid cell. The C2Bindex
582! in this case has the value for the last processed contact point
583! with contact region "cr".
584!
585 DO m=1,ncdatum
586 cr=contact_region(m)
587 dg=donor_grid(cr)
588 rg=receiver_grid(cr)
589 ibry=on_boundary(m)
590 DO ic=1,ncontact
591 IF ((dg.eq.receiver_grid(ic)).and. &
592 & (rg.eq.donor_grid(ic))) THEN
593 ccr=ic ! conjugate contact region
594 EXIT
595 END IF
596 END DO
597 IF ((ibry.eq.iwest ).or.(ibry.eq.ieast )) THEN
598 ir=irg_cp(m)
599 jr=jrg_cp(m)
600 ibmin(ibry,cr )=min(ir,ibmin(ibry,cr ))
601 ibmax(ibry,cr )=max(ir,ibmax(ibry,cr ))
602 jbmin(ibry,cr )=min(jr,jbmin(ibry,cr ))
603 jbmax(ibry,cr )=max(jr,jbmax(ibry,cr ))
604 bry_contact(ibry,cr ) % Ib(jr) = ir
605 bry_contact(ibry,cr ) % Jb(jr) = jr
606 bry_contact(ibry,cr ) % C2Bindex(jr) = m-nstru(cr)+1
607!
608 id=idg_cp(m)
609 jd=jdg_cp(m)
610 ibmin(ibry,ccr)=min(id,ibmin(ibry,ccr))
611 ibmax(ibry,ccr)=max(id,ibmax(ibry,ccr))
612 jbmin(ibry,ccr)=min(jd,jbmin(ibry,ccr))
613 jbmax(ibry,ccr)=max(jd,jbmax(ibry,ccr))
614 bry_contact(ibry,ccr) % Ib(jd) = id
615 bry_contact(ibry,ccr) % Jb(jd) = jd
616 bry_contact(ibry,ccr) % C2Bindex(jd) = m-nstru(cr)+1 ! same
617 ELSE IF ((ibry.eq.isouth).or.(ibry.eq.inorth)) THEN
618 ir=irg_cp(m)
619 jr=jrg_cp(m)
620 ibmin(ibry,cr )=min(ir,ibmin(ibry,cr))
621 ibmax(ibry,cr )=max(ir,ibmax(ibry,cr))
622 jbmin(ibry,cr )=min(jr,jbmin(ibry,cr))
623 jbmax(ibry,cr )=max(jr,jbmax(ibry,cr))
624 bry_contact(ibry,cr ) % Ib(ir) = ir
625 bry_contact(ibry,cr ) % Jb(ir) = jr
626 bry_contact(ibry,cr ) % C2Bindex(ir) = m-nstrv(cr)+1
627!
628 id=idg_cp(m)
629 jd=jdg_cp(m)
630 ibmin(ibry,ccr)=min(id,ibmin(ibry,ccr))
631 ibmax(ibry,ccr)=max(id,ibmax(ibry,ccr))
632 jbmin(ibry,ccr)=min(jd,jbmin(ibry,ccr))
633 jbmax(ibry,ccr)=max(jd,jbmax(ibry,ccr))
634 bry_contact(ibry,ccr) % Ib(id) = id
635 bry_contact(ibry,ccr) % Jb(id) = jd
636 bry_contact(ibry,ccr) % C2Bindex(id) = m-nstrv(cr)+1 ! same
637 END IF
638 END DO
639!
640! Set minimum and maximum indices to process at each boundary.
641!
642 DO cr=1,ncontact
643 DO ibry=1,4
644 IF (abs(ibmin(ibry,cr)).eq.abs(ispval)) THEN
645 bry_contact(ibry,cr) % Ibmin = ispval
646 ELSE
647 bry_contact(ibry,cr) % Ibmin = ibmin(ibry,cr)
648 END IF
649
650 IF (abs(ibmax(ibry,cr)).eq.abs(ispval)) THEN
651 bry_contact(ibry,cr) % Ibmax = ispval
652 ELSE
653 bry_contact(ibry,cr) % Ibmax = ibmax(ibry,cr)
654 END IF
655
656 IF (abs(jbmin(ibry,cr)).eq.abs(ispval)) THEN
657 bry_contact(ibry,cr) % Jbmin = ispval
658 ELSE
659 bry_contact(ibry,cr) % Jbmin = jbmin(ibry,cr)
660 END IF
661
662 IF (abs(jbmax(ibry,cr)).eq.abs(ispval)) THEN
663 bry_contact(ibry,cr) % Jbmax = ispval
664 ELSE
665 bry_contact(ibry,cr) % Jbmax = jbmax(ibry,cr)
666 END IF
667 END DO
668 END DO
669!
670!-----------------------------------------------------------------------
671! Deactivate boundary condition switches if contact point lay on the
672! physical nested grid boundary.
673!-----------------------------------------------------------------------
674!
675 DO cr=1,ncontact
676 rg=receiver_grid(cr)
677 IF (refinedgrid(rg)) THEN
678 IF (refinescale(rg).gt.0) THEN
679 lbc_apply(rg) % west = .false. ! This is a refinement
680 lbc_apply(rg) % south = .false. ! grid, so we need to
681 lbc_apply(rg) % east = .false. ! avoid applying lateral
682 lbc_apply(rg) % north = .false. ! boundary conditions
683 END IF
684 ELSE
685 DO ibry=1,4
686 imin=bry_contact(ibry,cr) % Ibmin ! Deactivate full or
687 imax=bry_contact(ibry,cr) % Ibmax ! partial lateral
688 jmin=bry_contact(ibry,cr) % Jbmin ! boundary conditions
689 jmax=bry_contact(ibry,cr) % Jbmax
690 SELECT CASE (ibry)
691 CASE (iwest)
692 IF ((jmin.ne.ispval).and.(jmax.ne.ispval)) THEN
693 DO j=jmin,jmax
694 lbc_apply(rg) % west (j) = .false.
695 END DO
696 END IF
697 CASE (isouth)
698 IF ((imin.ne.ispval).and.(imax.ne.ispval)) THEN
699 DO i=imin,imax
700 lbc_apply(rg) % south(i) = .false.
701 END DO
702 END IF
703 CASE (ieast)
704 IF ((jmin.ne.ispval).and.(jmax.ne.ispval)) THEN
705 DO j=jmin,jmax
706 lbc_apply(rg) % east (j) = .false.
707 END DO
708 END IF
709 CASE (inorth)
710 IF ((imin.ne.ispval).and.(imax.ne.ispval)) THEN
711 DO i=imin,imax
712 lbc_apply(rg) % north(i) = .false.
713 END DO
714 END IF
715 END SELECT
716 END DO
717 END IF
718 END DO
719!
720 RETURN
type(t_apply), dimension(:), allocatable lbc_apply
integer, dimension(:), allocatable n
Definition mod_param.F:479
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
real(r8), dimension(:), allocatable dmem
Definition mod_param.F:137
integer ngrids
Definition mod_param.F:113
integer, dimension(:), allocatable nt
Definition mod_param.F:489
integer, parameter iwest
integer, parameter isouth
integer, parameter ieast
logical, dimension(:), allocatable refinedgrid
integer, dimension(:), allocatable refinescale
integer, parameter inorth

References mod_param::bounds, bry_contact, contact_region, mod_param::dmem, donor_grid, idg_cp, mod_scalars::ieast, mod_scalars::inorth, irg_cp, mod_scalars::isouth, mod_scalars::iwest, jdg_cp, jrg_cp, mod_boundary::lbc_apply, mod_param::n, ncdatum, ncontact, mod_param::ngrids, nstru, nstrv, mod_param::nt, on_boundary, mod_kinds::r8, receiver_grid, mod_scalars::refinedgrid, and mod_scalars::refinescale.

Referenced by mod_arrays::roms_allocate_arrays().

Here is the caller graph for this function:

◆ deallocate_nesting()

subroutine, public mod_nesting::deallocate_nesting

Definition at line 723 of file mod_nesting.F.

724!
725!=======================================================================
726! !
727! This routine allocates and initializes nesting structure for 2D !
728! state variables. !
729! !
730!=======================================================================
731
732# ifdef SUBOBJECT_DEALLOCATION
733!
734 USE destroy_mod, ONLY : destroy
735# endif
736!
737! Local variable declarations.
738!
739 integer :: cr, ibry
740!
741 character (len=*), parameter :: MyFile = &
742 & __FILE__//", deallocate_nesting"
743
744# ifdef SUBOBJECT_DEALLOCATION
745!
746!-----------------------------------------------------------------------
747! Deallocate each variable in the derived-type T_BCP boundary
748! connectivity structure separately.
749!-----------------------------------------------------------------------
750!
751 DO cr=1,ncontact
752 DO ibry=1,4
753 IF (.not.destroy(ng, bry_contact(ibry,cr)%Ib, myfile, &
754 & __line__, 'BRY_CONTACT(ibry,cr)%Ib')) RETURN
755
756 IF (.not.destroy(ng, bry_contact(ibry,cr)%Jb, myfile, &
757 & __line__, 'BRY_CONTACT(ibry,cr)%Jb')) RETURN
758
759 IF (.not.destroy(ng, bry_contact(ibry,cr)%C2Bindex, myfile, &
760 & __line__, 'BRY_CONTACT(ibry,cr)%C2Bindex')) RETURN
761
762# ifdef NESTING_DEBUG
763 IF (.not.destroy(ng, bry_contact(ibry,cr)%Mflux, myfile, &
764 & __line__, 'BRY_CONTACT(ibry,cr)%Mflux')) RETURN
765# endif
766
767# ifdef SOLVE3D
768 IF (.not.destroy(ng, bry_contact(ibry,cr)%Tflux, myfile, &
769 & __line__, 'BRY_CONTACT(ibry,cr)%Tflux')) RETURN
770# endif
771 END DO
772 END DO
773!
774!-----------------------------------------------------------------------
775! Deallocate each variable in the derived-type T_NGC grid connectivity
776! structure separately.
777!-----------------------------------------------------------------------
778!
779 DO cr=1,ncontact
780 IF (.not.destroy(ng, rcontact(cr)%Irg, myfile, &
781 & __line__, 'Rcontact(cr)%Irg')) RETURN
782 IF (.not.destroy(ng, ucontact(cr)%Irg, myfile, &
783 & __line__, 'Ucontact(cr)%Irg')) RETURN
784 IF (.not.destroy(ng, vcontact(cr)%Irg, myfile, &
785 & __line__, 'Vcontact(cr)%Irg')) RETURN
786!
787 IF (.not.destroy(ng, rcontact(cr)%Jrg, myfile, &
788 & __line__, 'Rcontact(cr)%Jrg')) RETURN
789 IF (.not.destroy(ng, ucontact(cr)%Jrg, myfile, &
790 & __line__, 'Ucontact(cr)%Jrg')) RETURN
791 IF (.not.destroy(ng, vcontact(cr)%Jrg, myfile, &
792 & __line__, 'Vcontact(cr)%Jrg')) RETURN
793!
794 IF (.not.destroy(ng, rcontact(cr)%Idg, myfile, &
795 & __line__, 'Rcontact(cr)%Idg')) RETURN
796 IF (.not.destroy(ng, ucontact(cr)%Idg, myfile, &
797 & __line__, 'Ucontact(cr)%Idg')) RETURN
798 IF (.not.destroy(ng, vcontact(cr)%Idg, myfile, &
799 & __line__, 'Vcontact(cr)%Idg')) RETURN
800!
801 IF (.not.destroy(ng, rcontact(cr)%Jdg, myfile, &
802 & __line__, 'Rcontact(cr)%Jdg')) RETURN
803 IF (.not.destroy(ng, ucontact(cr)%Jdg, myfile, &
804 & __line__, 'Ucontact(cr)%Jdg')) RETURN
805 IF (.not.destroy(ng, vcontact(cr)%Jdg, myfile, &
806 & __line__, 'Vcontact(cr)%Jdg')) RETURN
807
808# ifdef SOLVE3D
809!
810 IF (.not.destroy(ng, rcontact(cr)%Kdg, myfile, &
811 & __line__, 'Rcontact(cr)%Kdg')) RETURN
812 IF (.not.destroy(ng, ucontact(cr)%Kdg, myfile, &
813 & __line__, 'Ucontact(cr)%Kdg')) RETURN
814 IF (.not.destroy(ng, vcontact(cr)%Kdg, myfile, &
815 & __line__, 'Vcontact(cr)%Kdg')) RETURN
816# endif
817!
818 IF (.not.destroy(ng, rcontact(cr)%Lweight, myfile, &
819 & __line__, 'Rcontact(cr)%Lweight')) RETURN
820 IF (.not.destroy(ng, ucontact(cr)%Lweight, myfile, &
821 & __line__, 'Ucontact(cr)%Lweight')) RETURN
822 IF (.not.destroy(ng, vcontact(cr)%Lweight, myfile, &
823 & __line__, 'Vcontact(cr)%Lweight')) RETURN
824!
825# ifdef WET_DRY
826 IF (.not.destroy(ng, rcontact(cr)%LweightUnmasked, myfile, &
827 & __line__, 'Rcontact(cr)%LweightUnmasked')) RETURN
828 IF (.not.destroy(ng, ucontact(cr)%LweightUnmasked, myfile, &
829 & __line__, 'Ucontact(cr)%LweightUnmasked')) RETURN
830 IF (.not.destroy(ng, vcontact(cr)%LweightUnmasked, myfile, &
831 & __line__, 'Vcontact(cr)%LweightUnmasked')) RETURN
832# endif
833
834# ifdef QUADRATIC_WEIGHTS
835!
836 IF (.not.destroy(ng, rcontact(cr)%Qweight, myfile, &
837 & __line__, 'Rcontact(cr)%Qweight')) RETURN
838 IF (.not.destroy(ng, ucontact(cr)%Qweight, myfile, &
839 & __line__, 'Ucontact(cr)%Qweight')) RETURN
840 IF (.not.destroy(ng, vcontact(cr)%Qweight, myfile, &
841 & __line__, 'Vcontact(cr)%Qweight')) RETURN
842
843# ifdef WET_DRY
844!
845 IF (.not.destroy(ng, rcontact(cr)%QweightUnmasked, myfile, &
846 & __line__, 'Rcontact(cr)%QweightUnmasked')) RETURN
847 IF (.not.destroy(ng, ucontact(cr)%QweightUnmasked, myfile, &
848 & __line__, 'Ucontact(cr)%QweightUnmasked')) RETURN
849 IF (.not.destroy(ng, vcontact(cr)%QweightUnmasked, myfile, &
850 & __line__, 'Vcontact(cr)%QweightUnmasked')) RETURN
851# endif
852# endif
853
854# ifdef SOLVE3D
855!
856 IF (.not.destroy(ng, rcontact(cr)%Vweight, myfile, &
857 & __line__, 'Rcontact(cr)%Vweight')) RETURN
858 IF (.not.destroy(ng, ucontact(cr)%Vweight, myfile, &
859 & __line__, 'Ucontact(cr)%Vweight')) RETURN
860 IF (.not.destroy(ng, vcontact(cr)%Vweight, myfile, &
861 & __line__, 'Vcontact(cr)%Vweight')) RETURN
862
863# if defined TANGENT || defined TL_IOMS
864!
865 IF (.not.destroy(ng, rcontact(cr)%tl_Vweight, myfile, &
866 & __line__, 'Rcontact(cr)%tl_Vweight')) RETURN
867 IF (.not.destroy(ng, ucontact(cr)%tl_Vweight, myfile, &
868 & __line__, 'Ucontact(cr)%tl_Vweight')) RETURN
869 IF (.not.destroy(ng, vcontact(cr)%tl_Vweight, myfile, &
870 & __line__, 'Vcontact(cr)%tl_Vweight')) RETURN
871# endif
872
873# ifdef ADJOINT
874!
875 IF (.not.destroy(ng, rcontact(cr)%ad_Vweight, myfile, &
876 & __line__, 'Rcontact(cr)%ad_Vweight')) RETURN
877 IF (.not.destroy(ng, ucontact(cr)%ad_Vweight, myfile, &
878 & __line__, 'Ucontact(cr)%ad_Vweight')) RETURN
879 IF (.not.destroy(ng, vcontact(cr)%ad_Vweight, myfile, &
880 & __line__, 'Vcontact(cr)%ad_Vweight')) RETURN
881# endif
882# endif
883 END DO
884!
885!-----------------------------------------------------------------------
886! Deallocate each variable in the derived-type T_NGC contact region
887! metrics structure separately.
888!-----------------------------------------------------------------------
889!
890 DO cr=1,ncontact
891 IF (.not.destroy(ng, contact_metric(cr)%angler, myfile, &
892 & __line__, 'CONTACT_METRIC(cr)%angler')) RETURN
893
894 IF (.not.destroy(ng, contact_metric(cr)%dndx, myfile, &
895 & __line__, 'CONTACT_METRIC(cr)%dndx')) RETURN
896
897 IF (.not.destroy(ng, contact_metric(cr)%dmde, myfile, &
898 & __line__, 'CONTACT_METRIC(cr)%dmde')) RETURN
899
900 IF (.not.destroy(ng, contact_metric(cr)%f, myfile, &
901 & __line__, 'CONTACT_METRIC(cr)%f')) RETURN
902
903 IF (.not.destroy(ng, contact_metric(cr)%h, myfile, &
904 & __line__, 'CONTACT_METRIC(cr)%h')) RETURN
905
906 IF (.not.destroy(ng, contact_metric(cr)%rmask, myfile, &
907 & __line__, 'CONTACT_METRIC(cr)%rmask')) RETURN
908
909 IF (.not.destroy(ng, contact_metric(cr)%umask, myfile, &
910 & __line__, 'CONTACT_METRIC(cr)%umask')) RETURN
911
912 IF (.not.destroy(ng, contact_metric(cr)%vmask, myfile, &
913 & __line__, 'CONTACT_METRIC(cr)%vmask')) RETURN
914
915 IF (.not.destroy(ng, contact_metric(cr)%pm, myfile, &
916 & __line__, 'CONTACT_METRIC(cr)%pm')) RETURN
917
918 IF (.not.destroy(ng, contact_metric(cr)%pn, myfile, &
919 & __line__, 'CONTACT_METRIC(cr)%pn')) RETURN
920
921 IF (.not.destroy(ng, contact_metric(cr)%Xr, myfile, &
922 & __line__, 'CONTACT_METRIC(cr)%Xr')) RETURN
923
924 IF (.not.destroy(ng, contact_metric(cr)%Yr, myfile, &
925 & __line__, 'CONTACT_METRIC(cr)%Yr')) RETURN
926
927 IF (.not.destroy(ng, contact_metric(cr)%Xu, myfile, &
928 & __line__, 'CONTACT_METRIC(cr)%Xu')) RETURN
929
930 IF (.not.destroy(ng, contact_metric(cr)%Yu, myfile, &
931 & __line__, 'CONTACT_METRIC(cr)%Yu')) RETURN
932
933 IF (.not.destroy(ng, contact_metric(cr)%Xv, myfile, &
934 & __line__, 'CONTACT_METRIC(cr)%Xv')) RETURN
935
936 IF (.not.destroy(ng, contact_metric(cr)%Yv, myfile, &
937 & __line__, 'CONTACT_METRIC(cr)%Yv')) RETURN
938 END DO
939!
940!-----------------------------------------------------------------------
941! Deallocate each variable in the derived-type T_COMPOSITE for
942! composite grids contact region structure separately.
943!-----------------------------------------------------------------------
944!
945 DO cr=1,ncontact
946 IF (.not.destroy(ng, composite(cr)%bustr, myfile, &
947 & __line__, 'COMPOSITE(cr)%bustr')) RETURN
948 IF (.not.destroy(ng, composite(cr)%bvstr, myfile, &
949 & __line__, 'COMPOSITE(cr)%bvstr')) RETURN
950
951 IF (.not.destroy(ng, composite(cr)%ubar, myfile, &
952 & __line__, 'COMPOSITE(cr)%ubar')) RETURN
953 IF (.not.destroy(ng, composite(cr)%vbar, myfile, &
954 & __line__, 'COMPOSITE(cr)%vbar')) RETURN
955 IF (.not.destroy(ng, composite(cr)%zeta, myfile, &
956 & __line__, 'COMPOSITE(cr)%zeta')) RETURN
957
958 IF (.not.destroy(ng, composite(cr)%rzeta, myfile, &
959 & __line__, 'COMPOSITE(cr)%rzeta')) RETURN
960
961# if defined TANGENT || defined TL_IOMS
962 IF (.not.destroy(ng, composite(cr)%tl_bustr, myfile, &
963 & __line__, 'COMPOSITE(cr)%tl_bustr')) RETURN
964 IF (.not.destroy(ng, composite(cr)%tl_bvstr, myfile, &
965 & __line__, 'COMPOSITE(cr)%tl_bvstr')) RETURN
966
967 IF (.not.destroy(ng, composite(cr)%tl_ubar, myfile, &
968 & __line__, 'COMPOSITE(cr)%tl_ubar')) RETURN
969 IF (.not.destroy(ng, composite(cr)%tl_vbar, myfile, &
970 & __line__, 'COMPOSITE(cr)%tl_vbar')) RETURN
971 IF (.not.destroy(ng, composite(cr)%tl_zeta, myfile, &
972 & __line__, 'COMPOSITE(cr)%tl_zeta')) RETURN
973
974 IF (.not.destroy(ng, composite(cr)%tl_rzeta, myfile, &
975 & __line__, 'COMPOSITE(cr)%tl_rzeta')) RETURN
976# endif
977
978# ifdef ADJOINT
979 IF (.not.destroy(ng, composite(cr)%ad_bustr, myfile, &
980 & __line__, 'COMPOSITE(cr)%ad_bustr')) RETURN
981 IF (.not.destroy(ng, composite(cr)%ad_bvstr, myfile, &
982 & __line__, 'COMPOSITE(cr)%ad_bvstr')) RETURN
983
984 IF (.not.destroy(ng, composite(cr)%ad_ubar, myfile, &
985 & __line__, 'COMPOSITE(cr)%ad_ubar')) RETURN
986 IF (.not.destroy(ng, composite(cr)%ad_vbar, myfile, &
987 & __line__, 'COMPOSITE(cr)%ad_vbar')) RETURN
988 IF (.not.destroy(ng, composite(cr)%ad_zeta, myfile, &
989 & __line__, 'COMPOSITE(cr)%ad_zeta')) RETURN
990
991 IF (.not.destroy(ng, composite(cr)%ad_rzeta, myfile, &
992 & __line__, 'COMPOSITE(cr)%ad_rzeta')) RETURN
993# endif
994
995# ifdef SOLVE3D
996 IF (.not.destroy(ng, composite(cr)%DU_avg1, myfile, &
997 & __line__, 'COMPOSITE(cr)%DU_avg1')) RETURN
998 IF (.not.destroy(ng, composite(cr)%DV_avg1, myfile, &
999 & __line__, 'COMPOSITE(cr)%DV_avg1')) RETURN
1000 IF (.not.destroy(ng, composite(cr)%Zt_avg1, myfile, &
1001 & __line__, 'COMPOSITE(cr)%Zt_avg1')) RETURN
1002
1003 IF (.not.destroy(ng, composite(cr)%u, myfile, &
1004 & __line__, 'COMPOSITE(cr)%u')) RETURN
1005 IF (.not.destroy(ng, composite(cr)%v, myfile, &
1006 & __line__, 'COMPOSITE(cr)%v')) RETURN
1007
1008 IF (.not.destroy(ng, composite(cr)%Huon, myfile, &
1009 & __line__, 'COMPOSITE(cr)%Huon')) RETURN
1010 IF (.not.destroy(ng, composite(cr)%Hvom, myfile, &
1011 & __line__, 'COMPOSITE(cr)%Hvom')) RETURN
1012
1013 IF (.not.destroy(ng, composite(cr)%t, myfile, &
1014 & __line__, 'COMPOSITE(cr)%t')) RETURN
1015
1016# if defined TANGENT || defined TL_IOMS
1017 IF (.not.destroy(ng, composite(cr)%tl_DU_avg1, myfile, &
1018 & __line__, 'COMPOSITE(cr)%tl_DU_avg1')) RETURN
1019 IF (.not.destroy(ng, composite(cr)%tl_DV_avg1, myfile, &
1020 & __line__, 'COMPOSITE(cr)%tl_DV_avg1')) RETURN
1021 IF (.not.destroy(ng, composite(cr)%tl_Zt_avg1, myfile, &
1022 & __line__, 'COMPOSITE(cr)%tl_Zt_avg1')) RETURN
1023
1024 IF (.not.destroy(ng, composite(cr)%tl_u, myfile, &
1025 & __line__, 'COMPOSITE(cr)%tl_u')) RETURN
1026 IF (.not.destroy(ng, composite(cr)%tl_v, myfile, &
1027 & __line__, 'COMPOSITE(cr)%tl_v')) RETURN
1028
1029 IF (.not.destroy(ng, composite(cr)%tl_Huon, myfile, &
1030 & __line__, 'COMPOSITE(cr)%tl_Huon')) RETURN
1031 IF (.not.destroy(ng, composite(cr)%tl_Hvom, myfile, &
1032 & __line__, 'COMPOSITE(cr)%tl_Hvom')) RETURN
1033
1034 IF (.not.destroy(ng, composite(cr)%tl_t, myfile, &
1035 & __line__, 'COMPOSITE(cr)%tl_t')) RETURN
1036# endif
1037
1038# ifdef ADJOINT
1039 IF (.not.destroy(ng, composite(cr)%ad_DU_avg1, myfile, &
1040 & __line__, 'COMPOSITE(cr)%ad_DU_avg1')) RETURN
1041 IF (.not.destroy(ng, composite(cr)%ad_DV_avg1, myfile, &
1042 & __line__, 'COMPOSITE(cr)%ad_DV_avg1')) RETURN
1043 IF (.not.destroy(ng, composite(cr)%ad_Zt_avg1, myfile, &
1044 & __line__, 'COMPOSITE(cr)%ad_Zt_avg1')) RETURN
1045
1046 IF (.not.destroy(ng, composite(cr)%ad_u, myfile, &
1047 & __line__, 'COMPOSITE(cr)%ad_u')) RETURN
1048 IF (.not.destroy(ng, composite(cr)%ad_v, myfile, &
1049 & __line__, 'COMPOSITE(cr)%ad_v')) RETURN
1050
1051 IF (.not.destroy(ng, composite(cr)%ad_Huon, myfile, &
1052 & __line__, 'COMPOSITE(cr)%ad_Huon')) RETURN
1053 IF (.not.destroy(ng, composite(cr)%ad_Hvom, myfile, &
1054 & __line__, 'COMPOSITE(cr)%ad_Hvom')) RETURN
1055
1056 IF (.not.destroy(ng, composite(cr)%ad_t, myfile, &
1057 & __line__, 'COMPOSITE(cr)%ad_t')) RETURN
1058# endif
1059# endif
1060 END DO
1061!
1062!-----------------------------------------------------------------------
1063! Deallocate each variable in the derived-type T_REFINED for
1064! refinement grids contact region structure separately.
1065!-----------------------------------------------------------------------
1066!
1067 DO cr=1,ncontact
1068 IF (.not.destroy(ng, refined(cr)%ubar, myfile, &
1069 & __line__, 'REFINED(cr)%ubar')) RETURN
1070 IF (.not.destroy(ng, refined(cr)%vbar, myfile, &
1071 & __line__, 'REFINED(cr)%vbar')) RETURN
1072 IF (.not.destroy(ng, refined(cr)%zeta, myfile, &
1073 & __line__, 'REFINED(cr)%zeta')) RETURN
1074
1075 IF (.not.destroy(ng, refined(cr)%U2d_flux, myfile, &
1076 & __line__, 'REFINED(cr)%U2d_flux')) RETURN
1077 IF (.not.destroy(ng, refined(cr)%V2d_flux, myfile, &
1078 & __line__, 'REFINED(cr)%V2d_flux')) RETURN
1079
1080 IF (.not.destroy(ng, refined(cr)%on_u, myfile, &
1081 & __line__, 'REFINED(cr)%on_u')) RETURN
1082 IF (.not.destroy(ng, refined(cr)%om_v, myfile, &
1083 & __line__, 'REFINED(cr)%om_v')) RETURN
1084
1085# if defined TANGENT || defined TL_IOMS
1086 IF (.not.destroy(ng, refined(cr)%tl_ubar, myfile, &
1087 & __line__, 'REFINED(cr)%tl_ubar')) RETURN
1088 IF (.not.destroy(ng, refined(cr)%tl_vbar, myfile, &
1089 & __line__, 'REFINED(cr)%tl_vbar')) RETURN
1090 IF (.not.destroy(ng, refined(cr)%tl_zeta, myfile, &
1091 & __line__, 'REFINED(cr)%tl_zeta')) RETURN
1092
1093 IF (.not.destroy(ng, refined(cr)%tl_U2d_flux, myfile, &
1094 & __line__, 'REFINED(cr)%tl_U2d_flux')) RETURN
1095 IF (.not.destroy(ng, refined(cr)%tl_V2d_flux, myfile, &
1096 & __line__, 'REFINED(cr)%tl_V2d_flux')) RETURN
1097# endif
1098
1099# ifdef ADJOINT
1100 IF (.not.destroy(ng, refined(cr)%ad_ubar, myfile, &
1101 & __line__, 'REFINED(cr)%ad_ubar')) RETURN
1102 IF (.not.destroy(ng, refined(cr)%ad_vbar, myfile, &
1103 & __line__, 'REFINED(cr)%ad_vbar')) RETURN
1104 IF (.not.destroy(ng, refined(cr)%ad_zeta, myfile, &
1105 & __line__, 'REFINED(cr)%ad_zeta')) RETURN
1106
1107 IF (.not.destroy(ng, refined(cr)%ad_U2d_flux, myfile, &
1108 & __line__, 'REFINED(cr)%ad_U2d_flux')) RETURN
1109 IF (.not.destroy(ng, refined(cr)%ad_V2d_flux, myfile, &
1110 & __line__, 'REFINED(cr)%ad_V2d_flux')) RETURN
1111# endif
1112
1113# ifdef SOLVE3D
1114 IF (.not.destroy(ng, refined(cr)%u, myfile, &
1115 & __line__, 'REFINED(cr)%u')) RETURN
1116 IF (.not.destroy(ng, refined(cr)%v, myfile, &
1117 & __line__, 'REFINED(cr)%v')) RETURN
1118
1119 IF (.not.destroy(ng, refined(cr)%t, myfile, &
1120 & __line__, 'REFINED(cr)%t')) RETURN
1121
1122# if defined TANGENT || defined TL_IOMS
1123 IF (.not.destroy(ng, refined(cr)%tl_u, myfile, &
1124 & __line__, 'REFINED(cr)%tl_u')) RETURN
1125 IF (.not.destroy(ng, refined(cr)%tl_v, myfile, &
1126 & __line__, 'REFINED(cr)%tl_v')) RETURN
1127
1128 IF (.not.destroy(ng, refined(cr)%tl_t, myfile, &
1129 & __line__, 'REFINED(cr)%tl_t')) RETURN
1130# endif
1131
1132# ifdef ADJOINT
1133 IF (.not.destroy(ng, refined(cr)%ad_u, myfile, &
1134 & __line__, 'REFINED(cr)%ad_u')) RETURN
1135 IF (.not.destroy(ng, refined(cr)%ad_v, myfile, &
1136 & __line__, 'REFINED(cr)%ad_v')) RETURN
1137
1138 IF (.not.destroy(ng, refined(cr)%ad_t, myfile, &
1139 & __line__, 'REFINED(cr)%ad_t')) RETURN
1140# endif
1141# endif
1142 END DO
1143# endif
1144!
1145!-----------------------------------------------------------------------
1146! Deallocate derived-type structures:
1147!-----------------------------------------------------------------------
1148!
1149! Boundary connectivity.
1150
1151 IF (allocated(bry_contact)) deallocate ( bry_contact )
1152!
1153! Grid connectivity.
1154!
1155 IF (allocated(rcontact)) deallocate ( rcontact )
1156 IF (allocated(ucontact)) deallocate ( ucontact )
1157 IF (allocated(vcontact)) deallocate ( vcontact )
1158!
1159! Contact region metrics.
1160!
1161 IF (allocated(contact_metric)) deallocate ( contact_metric )
1162!
1163! Composite grid contact regions.
1164!
1165 IF (allocated(composite)) deallocate ( composite )
1166!
1167! Deallocate refinement grids contact regions structure.
1168!
1169 IF (allocated(refined)) deallocate ( refined )
1170!
1171!-----------------------------------------------------------------------
1172! Deallocate other variables in module.
1173!-----------------------------------------------------------------------
1174!
1175 IF (allocated(contactregion)) THEN
1176 deallocate ( coarserdonor )
1177 END IF
1178
1179 IF (allocated(finerdonor)) THEN
1180 deallocate ( finerdonor )
1181 END IF
1182
1183 IF (allocated(donortofiner)) THEN
1184 deallocate ( donortofiner )
1185 END IF
1186
1187 IF (allocated(refinesteps)) THEN
1188 deallocate ( refinesteps )
1189 END IF
1190
1191 IF (allocated(refinestepscounter)) THEN
1192 deallocate ( refinestepscounter )
1193 END IF
1194
1195 IF (allocated(twowayinterval)) THEN
1196 deallocate ( twowayinterval )
1197 END IF
1198
1199 IF (allocated(telescoping)) THEN
1200 deallocate ( telescoping )
1201 END IF
1202
1203 IF (allocated(rollingindex)) THEN
1204 deallocate ( rollingindex )
1205 END IF
1206
1207 IF (allocated(rollingtime)) THEN
1208 deallocate ( rollingtime )
1209 END IF
1210!
1211 RETURN

References bry_contact, coarserdonor, composite, contact_metric, contactregion, donortofiner, finerdonor, ncontact, rcontact, refined, refinesteps, refinestepscounter, rollingindex, rollingtime, telescoping, twowayinterval, ucontact, and vcontact.

Referenced by mod_arrays::roms_deallocate_arrays().

Here is the caller graph for this function:

◆ initialize_nesting()

subroutine, public mod_nesting::initialize_nesting

Definition at line 1214 of file mod_nesting.F.

1215!
1216!=======================================================================
1217! !
1218! This routine initializes time varying nesting structures. !
1219! !
1220!=======================================================================
1221!
1222 USE mod_param
1223 USE mod_scalars
1224!
1225! Local variable declarations.
1226!
1227 integer :: cr
1228
1229 real(r8), parameter :: IniVal = 0.0_r8
1230!
1231!-----------------------------------------------------------------------
1232! Initialize time-varying contact regions structures. They are used
1233! to process values from contact regions to global kernel arrays and
1234! vice versa.
1235!-----------------------------------------------------------------------
1236!
1237! Composite grids contact region structure.
1238!
1239 IF (any(compositegrid)) THEN
1240 DO cr=1,ncontact
1241 composite(cr) % bustr = inival
1242 composite(cr) % bvstr = inival
1243
1244 composite(cr) % ubar = inival
1245 composite(cr) % vbar = inival
1246 composite(cr) % zeta = inival
1247
1248 composite(cr) % rzeta = inival
1249
1250# ifdef SOLVE3D
1251 composite(cr) % DU_avg1 = inival
1252 composite(cr) % DV_avg1 = inival
1253 composite(cr) % Zt_avg1 = inival
1254
1255 composite(cr) % u = inival
1256 composite(cr) % v = inival
1257
1258 composite(cr) % Huon = inival
1259 composite(cr) % Hvom = inival
1260
1261 composite(cr) % t = inival
1262# endif
1263 END DO
1264 END IF
1265!
1266! Refinement grids contact region structure.
1267!
1268 IF (any(refinedgrid(:))) THEN
1269 DO cr=1,ncontact
1270 refined(cr) % ubar = inival
1271 refined(cr) % vbar = inival
1272 refined(cr) % zeta = inival
1273
1274 refined(cr) % U2d_flux = inival
1275 refined(cr) % V2d_flux = inival
1276
1277 refined(cr) % on_u = inival
1278 refined(cr) % om_v = inival
1279
1280# ifdef SOLVE3D
1281 refined(cr) % u = inival
1282 refined(cr) % v = inival
1283
1284 refined(cr) % t = inival
1285# endif
1286 END DO
1287 END IF
1288
1289 RETURN
logical, dimension(:,:), allocatable compositegrid

References composite, mod_scalars::compositegrid, ncontact, refined, and mod_scalars::refinedgrid.

Referenced by mod_arrays::roms_initialize_arrays().

Here is the caller graph for this function:

Variable Documentation

◆ bry_contact

◆ coarserdonor

◆ composite

◆ contact_metric

type (t_ngm), dimension(:), allocatable mod_nesting::contact_metric

Definition at line 374 of file mod_nesting.F.

374 TYPE (T_NGM), allocatable :: CONTACT_METRIC(:) ! [Ncontact]

Referenced by deallocate_nesting(), get_grid_mod::get_grid_nf90(), get_grid_mod::get_grid_pio(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ contact_region

integer, dimension(:), allocatable mod_nesting::contact_region

Definition at line 209 of file mod_nesting.F.

209 integer, allocatable :: contact_region(:) ! [NCdatum]

Referenced by allocate_nesting(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ contactregion

logical, dimension(:,:), allocatable mod_nesting::contactregion

Definition at line 105 of file mod_nesting.F.

105 logical, allocatable :: ContactRegion(:,:) ! [4,Ngrids]

Referenced by deallocate_nesting(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ donor_grid

integer, dimension(:), allocatable mod_nesting::donor_grid

◆ donortofiner

logical, dimension(:), allocatable mod_nesting::donortofiner

Definition at line 111 of file mod_nesting.F.

111 logical, allocatable :: DonorToFiner(:) ! {Ngrids]

Referenced by ad_main3d(), deallocate_nesting(), main3d(), metrics_mod::metrics_tile(), set_contact_mod::set_contact_nf90(), set_contact_mod::set_contact_pio(), and tl_main3d().

◆ finerdonor

integer, dimension(:), allocatable mod_nesting::finerdonor

Definition at line 145 of file mod_nesting.F.

145 integer, allocatable :: FinerDonor(:) ! [Ngrids]

Referenced by nesting_mod::check_massflux(), deallocate_nesting(), metrics_mod::metrics_tile(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ get_vweights

logical mod_nesting::get_vweights

◆ i_left

◆ i_right

◆ idg_cp

integer, dimension(:), allocatable mod_nesting::idg_cp

Definition at line 211 of file mod_nesting.F.

211 integer, allocatable :: Idg_cp(:) ! [NCdatum]

Referenced by allocate_nesting(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ irg_cp

integer, dimension(:), allocatable mod_nesting::irg_cp

Definition at line 213 of file mod_nesting.F.

213 integer, allocatable :: Irg_cp(:) ! [NCdatum]

Referenced by allocate_nesting(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ j_bottom

◆ j_top

◆ jdg_cp

integer, dimension(:), allocatable mod_nesting::jdg_cp

Definition at line 212 of file mod_nesting.F.

212 integer, allocatable :: Jdg_cp(:) ! [NCdatum]

Referenced by allocate_nesting(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ jrg_cp

integer, dimension(:), allocatable mod_nesting::jrg_cp

Definition at line 214 of file mod_nesting.F.

214 integer, allocatable :: Jrg_cp(:) ! [NCdatum]

Referenced by allocate_nesting(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ n2dcs

integer, parameter mod_nesting::n2dcs = 10

◆ n2dfx

integer, parameter mod_nesting::n2dfx = 11

◆ n2dic

integer, parameter mod_nesting::n2dic = 2

◆ n2dps

integer, parameter mod_nesting::n2dps = 9

◆ n2way

integer, parameter mod_nesting::n2way = -1

Definition at line 80 of file mod_nesting.F.

80 integer, parameter :: n2way = -1 ! fine to course coupling

Referenced by ad_main3d(), ad_nesting_mod::ad_nesting(), main3d(), nesting_mod::nesting(), tl_main3d(), and tl_nesting_mod::tl_nesting().

◆ n3dic

integer, parameter mod_nesting::n3dic = 3

◆ n3dtv

integer, parameter mod_nesting::n3dtv = 13

◆ n3duv

integer, parameter mod_nesting::n3duv = 12

◆ nbstr

integer, parameter mod_nesting::nbstr = 5

◆ ncdatum

integer mod_nesting::ncdatum

Definition at line 203 of file mod_nesting.F.

203 integer :: NCdatum

Referenced by allocate_nesting(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ ncontact

◆ ncpoints

integer, dimension(:), allocatable mod_nesting::ncpoints

Definition at line 204 of file mod_nesting.F.

204 integer, allocatable :: NCpoints(:) ! [Ncontact]

Referenced by set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ ndxdy

integer, parameter mod_nesting::ndxdy = -5

Definition at line 76 of file mod_nesting.F.

76 integer, parameter :: ndxdy = -5 ! extract on_u and om_v

Referenced by ad_nesting_mod::ad_nesting(), nesting_mod::nesting(), set_grid(), and tl_nesting_mod::tl_nesting().

◆ nendr

integer, dimension(:), allocatable mod_nesting::nendr

◆ nendu

integer, dimension(:), allocatable mod_nesting::nendu

◆ nendv

integer, dimension(:), allocatable mod_nesting::nendv

◆ nfsic

integer, parameter mod_nesting::nfsic = 1

◆ ngetd

integer, parameter mod_nesting::ngetd = -4

Definition at line 77 of file mod_nesting.F.

77 integer, parameter :: ngetD = -4 ! extract donor grid data

Referenced by ad_main3d(), ad_nesting_mod::ad_nesting(), main3d(), nesting_mod::nesting(), tl_main3d(), and tl_nesting_mod::tl_nesting().

◆ nmask

integer, parameter mod_nesting::nmask = -3

Definition at line 78 of file mod_nesting.F.

78 integer, parameter :: nmask = -3 ! scale interpolation weights

Referenced by ad_main3d(), ad_nesting_mod::ad_nesting(), initial(), main3d(), nesting_mod::nesting(), roms_kernel_mod::nlm_initial(), tl_main3d(), and tl_nesting_mod::tl_nesting().

◆ nmflx

integer, parameter mod_nesting::nmflx = -6

Definition at line 75 of file mod_nesting.F.

75 integer, parameter :: nmflx = -6 ! check mass flux conservation

Referenced by ad_main3d(), ad_nesting_mod::ad_nesting(), main3d(), nesting_mod::nesting(), tl_main3d(), and tl_nesting_mod::tl_nesting().

◆ nputd

integer, parameter mod_nesting::nputd = -2

Definition at line 79 of file mod_nesting.F.

79 integer, parameter :: nputD = -2 ! fill contact points

Referenced by ad_main3d(), ad_nesting_mod::ad_nesting(), main3d(), nesting_mod::nesting(), tl_main3d(), and tl_nesting_mod::tl_nesting().

◆ nrhst

integer, parameter mod_nesting::nrhst = 6

◆ nstrr

integer, dimension(:), allocatable mod_nesting::nstrr

Definition at line 205 of file mod_nesting.F.

205 integer, allocatable :: NstrR(:), NendR(:) ! [Ncontact]

Referenced by ad_nesting_mod::ad_get_persisted2d(), nesting_mod::get_persisted2d(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ nstru

integer, dimension(:), allocatable mod_nesting::nstru

Definition at line 206 of file mod_nesting.F.

206 integer, allocatable :: NstrU(:), NendU(:) ! [Ncontact]

Referenced by ad_nesting_mod::ad_get_persisted2d(), allocate_nesting(), nesting_mod::get_persisted2d(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ nstrv

integer, dimension(:), allocatable mod_nesting::nstrv

Definition at line 207 of file mod_nesting.F.

207 integer, allocatable :: NstrV(:), NendV(:) ! [Ncontact]

Referenced by ad_nesting_mod::ad_get_persisted2d(), allocate_nesting(), nesting_mod::get_persisted2d(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ ntvic

integer, parameter mod_nesting::ntvic = 4

◆ nzeta

integer, parameter mod_nesting::nzeta = 7

◆ nzwgt

integer, parameter mod_nesting::nzwgt = 8

Definition at line 89 of file mod_nesting.F.

89 integer, parameter :: nzwgt = 8 ! 3D vertical weights

Referenced by ad_main3d(), ad_nesting_mod::ad_nesting(), initial(), nesting_mod::nesting(), roms_kernel_mod::nlm_initial(), tl_main3d(), and tl_nesting_mod::tl_nesting().

◆ on_boundary

integer, dimension(:), allocatable mod_nesting::on_boundary

Definition at line 210 of file mod_nesting.F.

210 integer, allocatable :: on_boundary(:) ! [NCdatum]

Referenced by ad_nesting_mod::ad_get_persisted2d(), allocate_nesting(), nesting_mod::get_persisted2d(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ rcontact

◆ receiver_grid

integer, dimension(:), allocatable mod_nesting::receiver_grid

◆ refined

◆ refinesteps

integer, dimension(:), allocatable mod_nesting::refinesteps

◆ refinestepscounter

integer, dimension(:), allocatable mod_nesting::refinestepscounter

Definition at line 162 of file mod_nesting.F.

162 integer, allocatable :: RefineStepsCounter(:) ! [Ngrids]

Referenced by deallocate_nesting(), ntimesteps(), set_contact_mod::set_contact_nf90(), and set_contact_mod::set_contact_pio().

◆ rollingindex

◆ rollingtime

◆ telescoping

logical, dimension(:), allocatable mod_nesting::telescoping

◆ twowayinterval

real(r8), dimension(:), allocatable mod_nesting::twowayinterval

Definition at line 167 of file mod_nesting.F.

167 real(r8), allocatable :: TwoWayInterval(:) ! [Ngrids]

Referenced by ad_main3d(), nesting_mod::check_massflux(), deallocate_nesting(), main3d(), set_contact_mod::set_contact_nf90(), set_contact_mod::set_contact_pio(), and tl_main3d().

◆ ucontact

◆ vcontact