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

Data Types

type  t_mixing
 

Functions/Subroutines

subroutine, public allocate_mixing (ng, lbi, ubi, lbj, ubj)
 
subroutine, public deallocate_mixing (ng)
 
subroutine, public initialize_mixing (ng, tile, model)
 

Variables

type(t_mixing), dimension(:), allocatable mixing
 

Function/Subroutine Documentation

◆ allocate_mixing()

subroutine, public mod_mixing::allocate_mixing ( integer, intent(in) ng,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj )

Definition at line 403 of file mod_mixing.F.

404!
405!=======================================================================
406! !
407! This routine allocates all variables in the module for all nested !
408! grids. !
409! !
410!=======================================================================
411!
412 USE mod_param
413 USE mod_scalars
414!
415! Imported variable declarations.
416!
417 integer, intent(in) :: ng, LBi, UBi, LBj, UBj
418!
419! Local variable declarations.
420!
421 real(r8) :: size2d
422!
423!-----------------------------------------------------------------------
424! Allocate module variables.
425!-----------------------------------------------------------------------
426!
427 IF (ng.eq.1) allocate ( mixing(ngrids) )
428!
429! Set horizontal array size.
430!
431 size2d=real((ubi-lbi+1)*(ubj-lbj+1),r8)
432!
433! Nonlinear model state.
434!
435#ifndef ANA_SPONGE
436 IF (luvsponge(ng)) THEN
437 allocate ( mixing(ng) % visc_factor(lbi:ubi,lbj:ubj) )
438 dmem(ng)=dmem(ng)+size2d
439 END IF
440#endif
441#if defined UV_VIS2
442 allocate ( mixing(ng) % visc2_p(lbi:ubi,lbj:ubj) )
443 dmem(ng)=dmem(ng)+size2d
444
445 allocate ( mixing(ng) % visc2_r(lbi:ubi,lbj:ubj) )
446 dmem(ng)=dmem(ng)+size2d
447#endif
448
449#ifdef UV_VIS4
450 allocate ( mixing(ng) % visc4_p(lbi:ubi,lbj:ubj) )
451 dmem(ng)=dmem(ng)+size2d
452
453 allocate ( mixing(ng) % visc4_r(lbi:ubi,lbj:ubj) )
454 dmem(ng)=dmem(ng)+size2d
455#endif
456
457# ifdef VISC_3DCOEF
458 allocate ( mixing(ng) % Hviscosity(lbi:ubi,lbj:ubj) )
459 dmem(ng)=dmem(ng)+size2d
460
461# ifdef UV_U3ADV_SPLIT
462 allocate ( mixing(ng) % Uvis3d_r(lbi:ubi,lbj:ubj,n(ng)) )
463 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
464
465 allocate ( mixing(ng) % Vvis3d_r(lbi:ubi,lbj:ubj,n(ng)) )
466 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
467# else
468 allocate ( mixing(ng) % visc3d_r(lbi:ubi,lbj:ubj,n(ng)) )
469 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
470# endif
471# endif
472
473#ifdef WEC
474 allocate ( mixing(ng) % rustr2d(lbi:ubi,lbj:ubj) )
475 dmem(ng)=dmem(ng)+size2d
476 allocate ( mixing(ng) % rvstr2d(lbi:ubi,lbj:ubj) )
477 dmem(ng)=dmem(ng)+size2d
478#endif
479#ifdef WEC_VF
480 allocate ( mixing(ng) % rurol2d(lbi:ubi,lbj:ubj) )
481 dmem(ng)=dmem(ng)+size2d
482 allocate ( mixing(ng) % rvrol2d(lbi:ubi,lbj:ubj) )
483 dmem(ng)=dmem(ng)+size2d
484 allocate ( mixing(ng) % rubrk2d(lbi:ubi,lbj:ubj) )
485 dmem(ng)=dmem(ng)+size2d
486 allocate ( mixing(ng) % rvbrk2d(lbi:ubi,lbj:ubj) )
487 dmem(ng)=dmem(ng)+size2d
488 allocate ( mixing(ng) % rukvf2d(lbi:ubi,lbj:ubj) )
489 dmem(ng)=dmem(ng)+size2d
490 allocate ( mixing(ng) % rvkvf2d(lbi:ubi,lbj:ubj) )
491 dmem(ng)=dmem(ng)+size2d
492# ifdef BOTTOM_STREAMING
493 allocate ( mixing(ng) % rubst2d(lbi:ubi,lbj:ubj) )
494 dmem(ng)=dmem(ng)+size2d
495 allocate ( mixing(ng) % rvbst2d(lbi:ubi,lbj:ubj) )
496 dmem(ng)=dmem(ng)+size2d
497# endif
498# ifdef SURFACE_STREAMING
499 allocate ( mixing(ng) % russt2d(lbi:ubi,lbj:ubj) )
500 dmem(ng)=dmem(ng)+size2d
501 allocate ( mixing(ng) % rvsst2d(lbi:ubi,lbj:ubj) )
502 dmem(ng)=dmem(ng)+size2d
503# endif
504#endif
505#ifdef SOLVE3D
506# ifdef WEC
507 allocate ( mixing(ng) % rustr3d(lbi:ubi,lbj:ubj,n(ng)) )
508 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
509 allocate ( mixing(ng) % rvstr3d(lbi:ubi,lbj:ubj,n(ng)) )
510 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
511# endif
512#endif
513
514#ifdef SOLVE3D
515# ifndef ANA_SPONGE
516 IF (any(ltracersponge(:,ng))) THEN
517 allocate ( mixing(ng) % diff_factor(lbi:ubi,lbj:ubj) )
518 dmem(ng)=dmem(ng)+size2d
519 END IF
520# endif
521
522# ifdef TS_DIF2
523 allocate ( mixing(ng) % diff2(lbi:ubi,lbj:ubj,nt(ng)) )
524 dmem(ng)=dmem(ng)+real(nt(ng),r8)*size2d
525# endif
526
527# ifdef TS_DIF4
528 allocate ( mixing(ng) % diff4(lbi:ubi,lbj:ubj,nt(ng)) )
529 dmem(ng)=dmem(ng)+real(nt(ng),r8)*size2d
530# endif
531
532# ifdef DIFF_3DCOEF
533 allocate ( mixing(ng) % Hdiffusion(lbi:ubi,lbj:ubj) )
534 dmem(ng)=dmem(ng)+size2d
535
536# ifdef TS_U3ADV_SPLIT
537 allocate ( mixing(ng) % diff3d_u(lbi:ubi,lbj:ubj,n(ng)) )
538 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
539
540 allocate ( mixing(ng) % diff3d_v(lbi:ubi,lbj:ubj,n(ng)) )
541 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
542# else
543 allocate ( mixing(ng) % diff3d_r(lbi:ubi,lbj:ubj,n(ng)) )
544 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
545# endif
546# endif
547
548 allocate ( mixing(ng) % Akv(lbi:ubi,lbj:ubj,0:n(ng)) )
549 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
550
551 allocate ( mixing(ng) % Akt(lbi:ubi,lbj:ubj,0:n(ng),nat) )
552 dmem(ng)=dmem(ng)+real((n(ng)+1)*nat,r8)*size2d
553
554# ifdef FLOAT_VWALK
555 allocate ( mixing(ng) % dAktdz(lbi:ubi,lbj:ubj,n(ng)) )
556 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
557# endif
558
559# if defined LMD_SKPP || defined LMD_BKPP || \
560 defined bulk_fluxes || defined balance_operator
561 allocate ( mixing(ng) % alpha(lbi:ubi,lbj:ubj) )
562 dmem(ng)=dmem(ng)+size2d
563
564 allocate ( mixing(ng) % beta(lbi:ubi,lbj:ubj) )
565 dmem(ng)=dmem(ng)+size2d
566# endif
567
568# ifdef BV_FREQUENCY
569 allocate ( mixing(ng) % bvf(lbi:ubi,lbj:ubj,0:n(ng)) )
570 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
571# endif
572
573# if defined GLS_MIXING || defined MY25_MIXING
574 allocate ( mixing(ng) % tke(lbi:ubi,lbj:ubj,0:n(ng),3) )
575 dmem(ng)=dmem(ng)+3.0_r8*real(n(ng)+1,r8)*size2d
576
577 allocate ( mixing(ng) % gls(lbi:ubi,lbj:ubj,0:n(ng),3) )
578 dmem(ng)=dmem(ng)+3.0_r8*real(n(ng)+1,r8)*size2d
579
580 allocate ( mixing(ng) % Lscale(lbi:ubi,lbj:ubj,0:n(ng)) )
581 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
582
583 allocate ( mixing(ng) % Akk(lbi:ubi,lbj:ubj,0:n(ng)) )
584 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
585
586# ifdef GLS_MIXING
587 allocate ( mixing(ng) % Akp(lbi:ubi,lbj:ubj,0:n(ng)) )
588 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
589# endif
590# endif
591
592# if defined LMD_MIXING && defined LMD_DDMIX
593 allocate ( mixing(ng) % alfaobeta(lbi:ubi,lbj:ubj,0:n(ng)) )
594 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
595# endif
596
597# if defined LMD_SKPP || defined SOLAR_SOURCE
598 allocate ( mixing(ng) % Jwtype(lbi:ubi,lbj:ubj) )
599 dmem(ng)=dmem(ng)+size2d
600# endif
601
602# if defined LMD_SKPP || defined LMD_BKPP
603 allocate ( mixing(ng) % ksbl(lbi:ubi,lbj:ubj) )
604 dmem(ng)=dmem(ng)+size2d
605
606 allocate ( mixing(ng) % hsbl(lbi:ubi,lbj:ubj) )
607 dmem(ng)=dmem(ng)+size2d
608
609# ifdef LMD_BKPP
610 allocate ( mixing(ng) % kbbl(lbi:ubi,lbj:ubj) )
611 dmem(ng)=dmem(ng)+size2d
612
613 allocate ( mixing(ng) % hbbl(lbi:ubi,lbj:ubj) )
614 dmem(ng)=dmem(ng)+size2d
615# endif
616
617# ifdef LMD_NONLOCAL
618 allocate ( mixing(ng) % ghats(lbi:ubi,lbj:ubj,0:n(ng),nat) )
619 dmem(ng)=dmem(ng)+real((n(ng)+1)*nat,r8)*size2d
620# endif
621
622# endif
623#endif
624
625#ifdef FOUR_DVAR
626!
627! Spatial convolution diffusion coefficients.
628!
629 allocate ( mixing(ng) % Kh(lbi:ubi,lbj:ubj) )
630 dmem(ng)=dmem(ng)+size2d
631
632# ifdef SOLVE3D
633 allocate ( mixing(ng) % Kv(lbi:ubi,lbj:ubj,0:n(ng)) )
634 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
635# endif
636#endif
637
638#if defined TANGENT || defined TL_IOMS
639!
640! Tangent linear model state.
641!
642# ifdef SOLVE3D
643# ifdef DIFF_3DCOEF
644# ifdef TS_U3ADV_SPLIT
645 allocate ( mixing(ng) % tl_diff3d_u(lbi:ubi,lbj:ubj,n(ng)) )
646 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
647
648 allocate ( mixing(ng) % tl_diff3d_v(lbi:ubi,lbj:ubj,n(ng)) )
649 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
650# else
651 allocate ( mixing(ng) % tl_diff3d_r(lbi:ubi,lbj:ubj,n(ng)) )
652 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
653# endif
654# endif
655
656# ifdef VISC_3DCOEF
657# ifdef UV_U3ADV_SPLIT
658 allocate ( mixing(ng) % tl_Uvis3d_r(lbi:ubi,lbj:ubj,n(ng)) )
659 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
660
661 allocate ( mixing(ng) % tl_Vvis3d_r(lbi:ubi,lbj:ubj,n(ng)) )
662 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
663# else
664 allocate ( mixing(ng) % tl_visc3d_r(lbi:ubi,lbj:ubj,n(ng)) )
665 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
666# endif
667# endif
668
669 allocate ( mixing(ng) % tl_Akv(lbi:ubi,lbj:ubj,0:n(ng)) )
670 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
671
672 allocate ( mixing(ng) % tl_Akt(lbi:ubi,lbj:ubj,0:n(ng),nat) )
673 dmem(ng)=dmem(ng)+real((n(ng)+1)*nat,r8)*size2d
674
675# if defined LMD_SKPP || defined LMD_BKPP || defined BULK_FLUXES
676 allocate ( mixing(ng) % tl_alpha(lbi:ubi,lbj:ubj) )
677 dmem(ng)=dmem(ng)+size2d
678
679 allocate ( mixing(ng) % tl_beta(lbi:ubi,lbj:ubj) )
680 dmem(ng)=dmem(ng)+size2d
681# endif
682
683# ifdef BV_FREQUENCY
684 allocate ( mixing(ng) % tl_bvf(lbi:ubi,lbj:ubj,0:n(ng)) )
685# endif
686
687# if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET
688 allocate ( mixing(ng) % tl_tke(lbi:ubi,lbj:ubj,0:n(ng),3) )
689 dmem(ng)=dmem(ng)+3.0_r8*real(n(ng)+1,r8)*size2d
690
691 allocate ( mixing(ng) % tl_gls(lbi:ubi,lbj:ubj,0:n(ng),3) )
692 dmem(ng)=dmem(ng)+3.0_r8*real(n(ng)+1,r8)*size2d
693
694 allocate ( mixing(ng) % tl_Lscale(lbi:ubi,lbj:ubj,0:n(ng)) )
695 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
696
697 allocate ( mixing(ng) % tl_Akk(lbi:ubi,lbj:ubj,0:n(ng)) )
698 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
699# endif
700
701# ifdef GLS_MIXING_NOT_YET
702 allocate ( mixing(ng) % tl_Akp(lbi:ubi,lbj:ubj,0:n(ng)) )
703 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
704# endif
705# endif
706#endif
707
708#ifdef ADJOINT
709!
710! Adjoint model state.
711!
712# ifdef SOLVE3D
713# ifdef DIFF_3DCOEF
714# ifdef TS_U3ADV_SPLIT
715 allocate ( mixing(ng) % ad_diff3d_u(lbi:ubi,lbj:ubj,n(ng)) )
716 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
717
718 allocate ( mixing(ng) % ad_diff3d_v(lbi:ubi,lbj:ubj,n(ng)) )
719 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
720# else
721 allocate ( mixing(ng) % ad_diff3d_r(lbi:ubi,lbj:ubj,n(ng)) )
722 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
723# endif
724# endif
725
726# ifdef VISC_3DCOEF
727# ifdef UV_U3ADV_SPLIT
728 allocate ( mixing(ng) % ad_Uvis3d_r(lbi:ubi,lbj:ubj,n(ng)) )
729 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
730
731 allocate ( mixing(ng) % ad_Vvis3d_r(lbi:ubi,lbj:ubj,n(ng)) )
732 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
733# else
734 allocate ( mixing(ng) % ad_visc3d_r(lbi:ubi,lbj:ubj,n(ng)) )
735 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
736# endif
737# endif
738
739 allocate ( mixing(ng) % ad_Akv(lbi:ubi,lbj:ubj,0:n(ng)) )
740 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
741
742 allocate ( mixing(ng) % ad_Akt(lbi:ubi,lbj:ubj,0:n(ng),nat) )
743 dmem(ng)=dmem(ng)+real((n(ng)+1)*nat,r8)*size2d
744
745# if defined LMD_SKPP || defined LMD_BKPP || defined BULK_FLUXES
746 allocate ( mixing(ng) % ad_alpha(lbi:ubi,lbj:ubj) )
747 dmem(ng)=dmem(ng)+size2d
748
749 allocate ( mixing(ng) % ad_beta(lbi:ubi,lbj:ubj) )
750 dmem(ng)=dmem(ng)+size2d
751# endif
752
753# ifdef BV_FREQUENCY
754 allocate ( mixing(ng) % ad_bvf(lbi:ubi,lbj:ubj,0:n(ng)) )
755 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
756# endif
757
758# if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET
759 allocate ( mixing(ng) % ad_tke(lbi:ubi,lbj:ubj,0:n(ng),3) )
760 dmem(ng)=dmem(ng)+3.0_r8*real(n(ng)+1,r8)*size2d
761
762 allocate ( mixing(ng) % ad_gls(lbi:ubi,lbj:ubj,0:n(ng),3) )
763 dmem(ng)=dmem(ng)+3.0_r8*real(n(ng)+1,r8)*size2d
764
765 allocate ( mixing(ng) % ad_Lscale(lbi:ubi,lbj:ubj,0:n(ng)) )
766 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
767
768 allocate ( mixing(ng) % ad_Akk(lbi:ubi,lbj:ubj,0:n(ng)) )
769 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
770# endif
771
772# ifdef GLS_MIXING_NOT_YET
773 allocate ( mixing(ng) % ad_Akp(lbi:ubi,lbj:ubj,0:n(ng)) )
774 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
775# endif
776# endif
777#endif
778
779#if defined FORWARD_READ && \
780 (defined tangent || defined tl_ioms || defined adjoint)
781# ifdef FORWARD_MIXING
782!
783! Latest two records of the nonlinear trajectory used to interpolate
784! the background state in the tangent linear and adjoint models.
785!
786 allocate ( mixing(ng) % AkvG(lbi:ubi,lbj:ubj,0:n(ng),2) )
787 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng)+1,r8)*size2d
788
789 allocate ( mixing(ng) % AktG(lbi:ubi,lbj:ubj,0:n(ng),2,nat) )
790 dmem(ng)=dmem(ng)+2.0_r8*real((n(ng)+1)*nat,r8)*size2d
791# endif
792
793# if defined LMD_MIXING_NOT_YET
794 allocate ( mixing(ng) % hsblG(lbi:ubi,lbj:ubj,2) )
795 dmem(ng)=dmem(ng)+2.0_r8*size2d
796# endif
797
798# if defined LMD_BKPP_NOT_YET
799 allocate ( mixing(ng) % hbblG(lbi:ubi,lbj:ubj,2) )
800 dmem(ng)=dmem(ng)+2.0_r8*size2d
801# endif
802
803# if defined LMD_NONLOCAL_NOT_YET
804 allocate ( mixing(ng) % ghatsG(lbi:ubi,lbj:ubj,0:n(ng),2,nat) )
805 dmem(ng)=dmem(ng)+2.0_r8*real((n(ng)+1)*nat,r8)*size2d
806# endif
807
808# if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET
809 allocate ( mixing(ng) % tkeG(lbi:ubi,lbj:ubj,0:n(ng),2) )
810 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng)+1,r8)*size2d
811
812 allocate ( mixing(ng) % glsG(lbi:ubi,lbj:ubj,0:n(ng),2) )
813 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng)+1,r8)*size2d
814
815 allocate ( mixing(ng) % LscaleG(lbi:ubi,lbj:ubj,0:n(ng),2) )
816 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng)+1,r8)*size2d
817
818 allocate ( mixing(ng) % AkkG(lbi:ubi,lbj:ubj,0:n(ng),2) )
819 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng)+1,r8)*size2d
820
821# ifdef GLS_MIXING_NOT_YET
822 allocate ( mixing(ng) % AkpG(lbi:ubi,lbj:ubj,0:n(ng),2) )
823 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng)+1,r8)*size2d
824# endif
825# endif
826#endif
827!
828 RETURN
integer nat
Definition mod_param.F:499
integer, dimension(:), allocatable n
Definition mod_param.F:479
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
logical, dimension(:), allocatable luvsponge
logical, dimension(:,:), allocatable ltracersponge

References mod_param::dmem, mod_scalars::ltracersponge, mod_scalars::luvsponge, mixing, mod_param::n, mod_param::nat, mod_param::ngrids, mod_param::nt, and mod_kinds::r8.

Referenced by mod_arrays::roms_allocate_arrays().

Here is the caller graph for this function:

◆ deallocate_mixing()

subroutine, public mod_mixing::deallocate_mixing ( integer, intent(in) ng)

Definition at line 831 of file mod_mixing.F.

832!
833!=======================================================================
834! !
835! This routine deallocates all variables in the module for all nested !
836! grids. !
837! !
838!=======================================================================
839!
840 USE mod_param
841 USE mod_scalars
842
843#ifdef SUBOBJECT_DEALLOCATION
844!
845 USE destroy_mod, ONLY : destroy
846#endif
847!
848! Imported variable declarations.
849!
850 integer, intent(in) :: ng
851!
852! Local variable declarations.
853!
854 character (len=*), parameter :: MyFile = &
855 & __FILE__//", deallocate_mixing"
856
857#ifdef SUBOBJECT_DEALLOCATION
858!
859!-----------------------------------------------------------------------
860! Deallocate each variable in the derived-type T_MIXING structure
861! separately.
862!-----------------------------------------------------------------------
863!
864! Nonlinear model state.
865!
866# ifndef ANA_SPONGE
867 IF (luvsponge(ng)) THEN
868 IF (.not.destroy(ng, mixing(ng)%visc_factor, myfile, &
869 & __line__, 'MIXING(ng)%visc_factor')) RETURN
870 END IF
871# endif
872
873# if defined UV_VIS2
874 IF (.not.destroy(ng, mixing(ng)%visc2_p, myfile, &
875 & __line__, 'MIXING(ng)%visc2_p')) RETURN
876
877 IF (.not.destroy(ng, mixing(ng)%visc2_r, myfile, &
878 & __line__, 'MIXING(ng)%visc2_r')) RETURN
879# endif
880
881# ifdef UV_VIS4
882 IF (.not.destroy(ng, mixing(ng)%visc4_p, myfile, &
883 & __line__, 'MIXING(ng)%visc4_p')) RETURN
884
885 IF (.not.destroy(ng, mixing(ng)%visc4_r, myfile, &
886 & __line__, 'MIXING(ng)%visc4_r')) RETURN
887# endif
888
889# ifdef VISC_3DCOEF
890 IF (.not.destroy(ng, mixing(ng)%Hviscosity, myfile, &
891 & __line__, 'MIXING(ng)%Hviscosity')) RETURN
892
893# ifdef UV_U3ADV_SPLIT
894 IF (.not.destroy(ng, mixing(ng)%Uvis3d_r, myfile, &
895 & __line__, 'MIXING(ng)%Uvis3d_r')) RETURN
896
897 IF (.not.destroy(ng, mixing(ng)%Vvis3d_r, myfile, &
898 & __line__, 'MIXING(ng)%Vvis3d_r')) RETURN
899# else
900 IF (.not.destroy(ng, mixing(ng)%visc3d_r, myfile, &
901 & __line__, 'MIXING(ng)%visc3d_r')) RETURN
902# endif
903# endif
904
905# ifdef WEC
906 IF (.not.destroy(ng, mixing(ng)%rustr2d, myfile, &
907 & __line__, 'MIXING(ng)%rustr2d')) RETURN
908
909 IF (.not.destroy(ng, mixing(ng)%rvstr2d, myfile, &
910 & __line__, 'MIXING(ng)%rvstr2d')) RETURN
911# endif
912# ifdef WEC_VF
913 IF (.not.destroy(ng, mixing(ng)%rurol2d, myfile, &
914 & __line__, 'MIXING(ng)%rurol2d')) RETURN
915
916 IF (.not.destroy(ng, mixing(ng)%rvrol2d, myfile, &
917 & __line__, 'MIXING(ng)%rvrol2d')) RETURN
918
919 IF (.not.destroy(ng, mixing(ng)%rubrk2d, myfile, &
920 & __line__, 'MIXING(ng)%rubrk2d')) RETURN
921
922 IF (.not.destroy(ng, mixing(ng)%rvbrk2d, myfile, &
923 & __line__, 'MIXING(ng)%rvbrk2d')) RETURN
924
925 IF (.not.destroy(ng, mixing(ng)%rukvf2d, myfile, &
926 & __line__, 'MIXING(ng)%rukvf2d')) RETURN
927
928 IF (.not.destroy(ng, mixing(ng)%rvkvf2d, myfile, &
929 & __line__, 'MIXING(ng)%rvkvf2d')) RETURN
930
931# ifdef BOTTOM_STREAMING
932 IF (.not.destroy(ng, mixing(ng)%rubst2d, myfile, &
933 & __line__, 'MIXING(ng)%rubst2d')) RETURN
934
935 IF (.not.destroy(ng, mixing(ng)%rvbst2d, myfile, &
936 & __line__, 'MIXING(ng)%rvbst2d')) RETURN
937# endif
938# ifdef SURFACE_STREAMING
939 IF (.not.destroy(ng, mixing(ng)%russt2d, myfile, &
940 & __line__, 'MIXING(ng)%russt2d')) RETURN
941
942 IF (.not.destroy(ng, mixing(ng)%rvsst2d, myfile, &
943 & __line__, 'MIXING(ng)%rvsst2d')) RETURN
944# endif
945# endif
946
947# ifdef SOLVE3D
948# ifdef WEC_VF
949 IF (.not.destroy(ng, mixing(ng)%rustr3d, myfile, &
950 & __line__, 'MIXING(ng)%rustr3d')) RETURN
951
952 IF (.not.destroy(ng, mixing(ng)%rvstr3d, myfile, &
953 & __line__, 'MIXING(ng)%rvstr3d')) RETURN
954# endif
955# endif
956
957# ifdef SOLVE3D
958# ifndef ANA_SPONGE
959 IF (any(ltracersponge(:,ng))) THEN
960 IF (.not.destroy(ng, mixing(ng)%diff_factor, myfile, &
961 & __line__, 'MIXING(ng)%diff_factor')) RETURN
962 END IF
963# endif
964
965# ifdef TS_DIF2
966 IF (.not.destroy(ng, mixing(ng)%diff2, myfile, &
967 & __line__, 'MIXING(ng)%diff2')) RETURN
968# endif
969
970# ifdef TS_DIF4
971 IF (.not.destroy(ng, mixing(ng)%diff4, myfile, &
972 & __line__, 'MIXING(ng)%diff4')) RETURN
973# endif
974
975# ifdef DIFF_3DCOEF
976 IF (.not.destroy(ng, mixing(ng)%Hdiffusion, myfile, &
977 & __line__, 'MIXING(ng)%Hdiffusion')) RETURN
978
979# ifdef TS_U3ADV_SPLIT
980 IF (.not.destroy(ng, mixing(ng)%diff3d_u, myfile, &
981 & __line__, 'MIXING(ng)%diff3d_u')) RETURN
982
983 IF (.not.destroy(ng, mixing(ng)%diff3d_v, myfile, &
984 & __line__, 'MIXING(ng)%diff3d_v')) RETURN
985# else
986 IF (.not.destroy(ng, mixing(ng)%diff3d_r, myfile, &
987 & __line__, 'MIXING(ng)%diff3d_r')) RETURN
988# endif
989# endif
990
991 IF (.not.destroy(ng, mixing(ng)%Akv, myfile, &
992 & __line__, 'MIXING(ng)%Akv')) RETURN
993
994 IF (.not.destroy(ng, mixing(ng)%Akt, myfile, &
995 & __line__, 'MIXING(ng)%Akt')) RETURN
996
997# ifdef FLOAT_VWALK
998 IF (.not.destroy(ng, mixing(ng)%dAktdz, myfile, &
999 & __line__, 'MIXING(ng)%dAktdz')) RETURN
1000# endif
1001
1002# if defined LMD_SKPP || defined LMD_BKPP || \
1003 defined bulk_fluxes || defined balance_operator
1004 IF (.not.destroy(ng, mixing(ng)%alpha, myfile, &
1005 & __line__, 'MIXING(ng)%alpha')) RETURN
1006
1007 IF (.not.destroy(ng, mixing(ng)%beta, myfile, &
1008 & __line__, 'MIXING(ng)%beta')) RETURN
1009# endif
1010
1011# ifdef BV_FREQUENCY
1012 IF (.not.destroy(ng, mixing(ng)%bvf, myfile, &
1013 & __line__, 'MIXING(ng)%bvf')) RETURN
1014# endif
1015
1016# if defined GLS_MIXING || defined MY25_MIXING
1017 IF (.not.destroy(ng, mixing(ng)%tke, myfile, &
1018 & __line__, 'MIXING(ng)%tke')) RETURN
1019
1020 IF (.not.destroy(ng, mixing(ng)%gls, myfile, &
1021 & __line__, 'MIXING(ng)%gls')) RETURN
1022
1023 IF (.not.destroy(ng, mixing(ng)%Lscale, myfile, &
1024 & __line__, 'MIXING(ng)%Lscale')) RETURN
1025
1026 IF (.not.destroy(ng, mixing(ng)%Akk, myfile, &
1027 & __line__, 'MIXING(ng)%Akk')) RETURN
1028
1029# ifdef GLS_MIXING
1030 IF (.not.destroy(ng, mixing(ng)%Akp, myfile, &
1031 & __line__, 'MIXING(ng)%Akp')) RETURN
1032# endif
1033# endif
1034
1035# if defined LMD_MIXING && defined LMD_DDMIX
1036 IF (.not.destroy(ng, mixing(ng)%alfaobeta, myfile, &
1037 & __line__, 'MIXING(ng)%alfaobeta')) RETURN
1038# endif
1039
1040# if defined LMD_SKPP || defined SOLAR_SOURCE
1041 IF (.not.destroy(ng, mixing(ng)%Jwtype, myfile, &
1042 & __line__, 'MIXING(ng)%Jwtype')) RETURN
1043# endif
1044
1045# if defined LMD_SKPP || defined LMD_BKPP
1046 IF (.not.destroy(ng, mixing(ng)%ksbl, myfile, &
1047 & __line__, 'MIXING(ng)%ksbl')) RETURN
1048
1049 IF (.not.destroy(ng, mixing(ng)%hsbl, myfile, &
1050 & __line__, 'MIXING(ng)%hsbl')) RETURN
1051
1052# ifdef LMD_BKPP
1053 IF (.not.destroy(ng, mixing(ng)%kbbl, myfile, &
1054 & __line__, 'MIXING(ng)%kbbl')) RETURN
1055
1056 IF (.not.destroy(ng, mixing(ng)%hbbl, myfile, &
1057 & __line__, 'MIXING(ng)%hbbl')) RETURN
1058# endif
1059
1060# ifdef LMD_NONLOCAL
1061 IF (.not.destroy(ng, mixing(ng)%ghats, myfile, &
1062 & __line__, 'MIXING(ng)%ghats')) RETURN
1063# endif
1064
1065# endif
1066# endif
1067
1068# ifdef FOUR_DVAR
1069!
1070! Spatial convolution diffusion coefficients.
1071!
1072 IF (.not.destroy(ng, mixing(ng)%Kh, myfile, &
1073 & __line__, 'MIXING(ng)%Kh')) RETURN
1074
1075# ifdef SOLVE3D
1076 IF (.not.destroy(ng, mixing(ng)%Kv, myfile, &
1077 & __line__, 'MIXING(ng)%Kv')) RETURN
1078# endif
1079# endif
1080
1081# if defined TANGENT || defined TL_IOMS
1082!
1083! Tangent linear model state.
1084!
1085# ifdef SOLVE3D
1086# ifdef DIFF_3DCOEF
1087# ifdef TS_U3ADV_SPLIT
1088 IF (.not.destroy(ng, mixing(ng)%tl_diff3d_u, myfile, &
1089 & __line__, 'MIXING(ng)%tl_diff3d_u')) RETURN
1090
1091 IF (.not.destroy(ng, mixing(ng)%tl_diff3d_v, myfile, &
1092 & __line__, 'MIXING(ng)%tl_diff3d_v')) RETURN
1093# else
1094 IF (.not.destroy(ng, mixing(ng)%tl_diff3d_r, myfile, &
1095 & __line__, 'MIXING(ng)%tl_diff3d_r')) RETURN
1096# endif
1097# endif
1098
1099# ifdef VISC_3DCOEF
1100# ifdef UV_U3ADV_SPLIT
1101 IF (.not.destroy(ng, mixing(ng)%tl_Uvis3d_r, myfile, &
1102 & __line__, 'MIXING(ng)%tl_Uvis3d_r')) RETURN
1103
1104 IF (.not.destroy(ng, mixing(ng)%tl_Vvis3d_r, myfile, &
1105 & __line__, 'MIXING(ng)%tl_Vvis3d_r')) RETURN
1106# else
1107 IF (.not.destroy(ng, mixing(ng)%tl_visc3d_r, myfile, &
1108 & __line__, 'MIXING(ng)%tl_visc3d_r')) RETURN
1109# endif
1110# endif
1111
1112 IF (.not.destroy(ng, mixing(ng)%tl_Akv, myfile, &
1113 & __line__, 'MIXING(ng)%tl_Akv')) RETURN
1114
1115 IF (.not.destroy(ng, mixing(ng)%tl_Akt, myfile, &
1116 & __line__, 'MIXING(ng)%tl_Akt')) RETURN
1117
1118# if defined LMD_SKPP || defined LMD_BKPP || defined BULK_FLUXES
1119 IF (.not.destroy(ng, mixing(ng)%tl_alpha, myfile, &
1120 & __line__, 'MIXING(ng)%tl_alpha')) RETURN
1121
1122 IF (.not.destroy(ng, mixing(ng)%tl_beta, myfile, &
1123 & __line__, 'MIXING(ng)%tl_beta')) RETURN
1124# endif
1125
1126# ifdef BV_FREQUENCY
1127 IF (.not.destroy(ng, mixing(ng)%tl_bvf, myfile, &
1128 & __line__, 'MIXING(ng)%tl_bvf')) RETURN
1129# endif
1130
1131# if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET
1132 IF (.not.destroy(ng, mixing(ng)%tl_tke, myfile, &
1133 & __line__, 'MIXING(ng)%tl_tke')) RETURN
1134
1135 IF (.not.destroy(ng, mixing(ng)%tl_gls, myfile, &
1136 & __line__, 'MIXING(ng)%tl_gls')) RETURN
1137
1138 IF (.not.destroy(ng, mixing(ng)%tl_Lscale, myfile, &
1139 & __line__, 'MIXING(ng)%tl_Lscale')) RETURN
1140
1141 IF (.not.destroy(ng, mixing(ng)%tl_Akk, myfile, &
1142 & __line__, 'MIXING(ng)%tl_Akk')) RETURN
1143# endif
1144
1145# ifdef GLS_MIXING_NOT_YET
1146 IF (.not.destroy(ng, mixing(ng)%tl_Akp, myfile, &
1147 & __line__, 'MIXING(ng)%tl_Akp')) RETURN
1148# endif
1149# endif
1150# endif
1151
1152# ifdef ADJOINT
1153!
1154! Adjoint model state.
1155!
1156# ifdef SOLVE3D
1157# ifdef DIFF_3DCOEF
1158# ifdef TS_U3ADV_SPLIT
1159 IF (.not.destroy(ng, mixing(ng)%ad_diff3d_u, myfile, &
1160 & __line__, 'MIXING(ng)%ad_diff3d_u')) RETURN
1161
1162 IF (.not.destroy(ng, mixing(ng)%ad_diff3d_v, myfile, &
1163 & __line__, 'MIXING(ng)%ad_diff3d_v')) RETURN
1164# else
1165 IF (.not.destroy(ng, mixing(ng)%ad_diff3d_r, myfile, &
1166 & __line__, 'MIXING(ng)%ad_diff3d_r')) RETURN
1167# endif
1168# endif
1169
1170# ifdef VISC_3DCOEF
1171# ifdef UV_U3ADV_SPLIT
1172 IF (.not.destroy(ng, mixing(ng)%ad_Uvis3d_r, myfile, &
1173 & __line__, 'MIXING(ng)%ad_Uvis3d_r')) RETURN
1174
1175 IF (.not.destroy(ng, mixing(ng)%ad_Vvis3d_r, myfile, &
1176 & __line__, 'MIXING(ng)%ad_Vvis3d_r')) RETURN
1177# else
1178 IF (.not.destroy(ng, mixing(ng)%ad_visc3d_r, myfile, &
1179 & __line__, 'MIXING(ng)%ad_visc3d_r')) RETURN
1180# endif
1181# endif
1182
1183 IF (.not.destroy(ng, mixing(ng)%ad_Akv, myfile, &
1184 & __line__, 'MIXING(ng)%ad_Akv')) RETURN
1185
1186 IF (.not.destroy(ng, mixing(ng)%ad_Akt, myfile, &
1187 & __line__, 'MIXING(ng)%ad_Akt')) RETURN
1188
1189# if defined LMD_SKPP || defined LMD_BKPP || defined BULK_FLUXES
1190 IF (.not.destroy(ng, mixing(ng)%ad_alpha, myfile, &
1191 & __line__, 'MIXING(ng)%ad_alpha')) RETURN
1192
1193 IF (.not.destroy(ng, mixing(ng)%ad_beta, myfile, &
1194 & __line__, 'MIXING(ng)%ad_beta')) RETURN
1195# endif
1196
1197# ifdef BV_FREQUENCY
1198 IF (.not.destroy(ng, mixing(ng)%ad_bvf, myfile, &
1199 & __line__, 'MIXING(ng)%ad_bvf')) RETURN
1200# endif
1201
1202# if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET
1203 IF (.not.destroy(ng, mixing(ng)%ad_tke, myfile, &
1204 & __line__, 'MIXING(ng)%ad_tke')) RETURN
1205
1206 IF (.not.destroy(ng, mixing(ng)%ad_gls, myfile, &
1207 & __line__, 'MIXING(ng)%ad_gls')) RETURN
1208
1209 IF (.not.destroy(ng, mixing(ng)%ad_Lscale, myfile, &
1210 & __line__, 'MIXING(ng)%ad_Lscale')) RETURN
1211
1212 IF (.not.destroy(ng, mixing(ng)%ad_Akk, myfile, &
1213 & __line__, 'MIXING(ng)%ad_Akk')) RETURN
1214# endif
1215
1216# ifdef GLS_MIXING_NOT_YET
1217 IF (.not.destroy(ng, mixing(ng)%ad_Akp, myfile, &
1218 & __line__, 'MIXING(ng)%ad_Akp')) RETURN
1219# endif
1220# endif
1221# endif
1222
1223# if defined FORWARD_READ && \
1224 (defined tangent || defined tl_ioms || defined adjoint)
1225# ifdef FORWARD_MIXING
1226!
1227! Latest two records of the nonlinear trajectory used to interpolate
1228! the background state in the tangent linear and adjoint models.
1229!
1230 IF (.not.destroy(ng, mixing(ng)%AkvG, myfile, &
1231 & __line__, 'MIXING(ng)%AkvG')) RETURN
1232
1233 IF (.not.destroy(ng, mixing(ng)%AktG, myfile, &
1234 & __line__, 'MIXING(ng)%AktG')) RETURN
1235# endif
1236
1237# if defined LMD_MIXING_NOT_YET
1238 IF (.not.destroy(ng, mixing(ng)%hsblG, myfile, &
1239 & __line__, 'MIXING(ng)%hsblG')) RETURN
1240# endif
1241
1242# if defined LMD_BKPP_NOT_YET
1243 IF (.not.destroy(ng, mixing(ng)%hbblG, myfile, &
1244 & __line__, 'MIXING(ng)%hbblG')) RETURN
1245# endif
1246
1247# if defined LMD_NONLOCAL_NOT_YET
1248 IF (.not.destroy(ng, mixing(ng)%ghatsG, myfile, &
1249 & __line__, 'MIXING(ng)%ghatsG')) RETURN
1250# endif
1251
1252# if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET
1253 IF (.not.destroy(ng, mixing(ng)%tkeG, myfile, &
1254 & __line__, 'MIXING(ng)%tkeG')) RETURN
1255
1256 IF (.not.destroy(ng, mixing(ng)%glsG, myfile, &
1257 & __line__, 'MIXING(ng)%glsG')) RETURN
1258
1259 IF (.not.destroy(ng, mixing(ng)%LscaleG, myfile, &
1260 & __line__, 'MIXING(ng)%LscaleG')) RETURN
1261
1262 IF (.not.destroy(ng, mixing(ng)%AkkG, myfile, &
1263 & __line__, 'MIXING(ng)%AkkG')) RETURN
1264
1265# ifdef GLS_MIXING_NOT_YET
1266 IF (.not.destroy(ng, mixing(ng)%AkpG, myfile, &
1267 & __line__, 'MIXING(ng)%AkpG')) RETURN
1268# endif
1269# endif
1270# endif
1271#endif
1272!
1273!-----------------------------------------------------------------------
1274! Deallocate derived-type MIXING structure.
1275!-----------------------------------------------------------------------
1276!
1277 IF (ng.eq.ngrids) THEN
1278 IF (allocated(mixing)) deallocate ( mixing )
1279 END IF
1280!
1281 RETURN

References mod_scalars::ltracersponge, mod_scalars::luvsponge, mixing, and mod_param::ngrids.

Referenced by mod_arrays::roms_deallocate_arrays().

Here is the caller graph for this function:

◆ initialize_mixing()

subroutine, public mod_mixing::initialize_mixing ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model )

Definition at line 1284 of file mod_mixing.F.

1285!
1286!=======================================================================
1287! !
1288! This routine allocates and initializes all variables in module !
1289! "mod_mixing" for all nested grids. !
1290! !
1291!=======================================================================
1292!
1293 USE mod_param
1294 USE mod_scalars
1295!
1296! Imported variable declarations.
1297!
1298 integer, intent(in) :: ng, tile, model
1299!
1300! Local variable declarations.
1301!
1302 integer :: Imin, Imax, Jmin, Jmax
1303 integer :: i, j
1304#ifdef SOLVE3D
1305 integer :: itrc, k
1306#endif
1307
1308 real(r8), parameter :: IniVal = 0.0_r8
1309 real(r8) :: cff1, cff2, cff3, cff4
1310
1311#include "set_bounds.h"
1312!
1313! Set array initialization range.
1314!
1315#ifdef DISTRIBUTE
1316 imin=bounds(ng)%LBi(tile)
1317 imax=bounds(ng)%UBi(tile)
1318 jmin=bounds(ng)%LBj(tile)
1319 jmax=bounds(ng)%UBj(tile)
1320#else
1321 IF (domain(ng)%Western_Edge(tile)) THEN
1322 imin=bounds(ng)%LBi(tile)
1323 ELSE
1324 imin=istr
1325 END IF
1326 IF (domain(ng)%Eastern_Edge(tile)) THEN
1327 imax=bounds(ng)%UBi(tile)
1328 ELSE
1329 imax=iend
1330 END IF
1331 IF (domain(ng)%Southern_Edge(tile)) THEN
1332 jmin=bounds(ng)%LBj(tile)
1333 ELSE
1334 jmin=jstr
1335 END IF
1336 IF (domain(ng)%Northern_Edge(tile)) THEN
1337 jmax=bounds(ng)%UBj(tile)
1338 ELSE
1339 jmax=jend
1340 END IF
1341#endif
1342!
1343!-----------------------------------------------------------------------
1344! Initialize module variables.
1345!-----------------------------------------------------------------------
1346!
1347! Nonlinear model state.
1348!
1349 IF ((model.eq.0).or.(model.eq.inlm)) THEN
1350 DO j=jmin,jmax
1351#ifndef ANA_SPONGE
1352 IF (luvsponge(ng)) THEN
1353 DO i=imin,imax
1354 mixing(ng) % visc_factor(i,j) = inival
1355 END DO
1356 END IF
1357#endif
1358#if defined UV_VIS2
1359 DO i=imin,imax
1360 mixing(ng) % visc2_p(i,j) = inival
1361 mixing(ng) % visc2_r(i,j) = inival
1362 END DO
1363#endif
1364#ifdef UV_VIS4
1365 DO i=imin,imax
1366 mixing(ng) % visc4_p(i,j) = inival
1367 mixing(ng) % visc4_r(i,j) = inival
1368 END DO
1369#endif
1370#ifdef SOLVE3D
1371# ifdef VISC_3DCOEF
1372 DO k=1,n(ng)
1373 DO i=imin,imax
1374 mixing(ng) % Hviscosity(i,j) = inival
1375# ifdef UV_U3ADV_SPLIT
1376 mixing(ng) % Uvis3d_r(i,j,k) = inival
1377 mixing(ng) % Vvis3d_r(i,j,k) = inival
1378# else
1379 mixing(ng) % visc3d_r(i,j,k) = inival
1380# endif
1381 END DO
1382 END DO
1383# endif
1384#endif
1385#ifdef WEC
1386 DO i=imin,imax
1387 mixing(ng) % rustr2d(i,j) = inival
1388 mixing(ng) % rvstr2d(i,j) = inival
1389# ifdef WEC_VF
1390 mixing(ng) % rurol2d(i,j) = inival
1391 mixing(ng) % rvrol2d(i,j) = inival
1392 mixing(ng) % rubrk2d(i,j) = inival
1393 mixing(ng) % rvbrk2d(i,j) = inival
1394 mixing(ng) % rukvf2d(i,j) = inival
1395 mixing(ng) % rvkvf2d(i,j) = inival
1396# ifdef BOTTOM_STREAMING
1397 mixing(ng) % rubst2d(i,j) = inival
1398 mixing(ng) % rvbst2d(i,j) = inival
1399# endif
1400# ifdef SURFACE_STREAMING
1401 mixing(ng) % russt2d(i,j) = inival
1402 mixing(ng) % rvsst2d(i,j) = inival
1403# endif
1404# endif
1405# ifdef SOLVE3D
1406 DO k=1,n(ng)
1407 mixing(ng) % rustr3d(i,j,k) = inival
1408 mixing(ng) % rvstr3d(i,j,k) = inival
1409 END DO
1410# endif
1411 END DO
1412#endif
1413#ifdef SOLVE3D
1414# ifndef ANA_SPONGE
1415 IF (any(ltracersponge(:,ng))) THEN
1416 DO i=imin,imax
1417 mixing(ng) % diff_factor(i,j) = inival
1418 END DO
1419 END IF
1420# endif
1421# ifdef TS_DIF2
1422 DO itrc=1,nt(ng)
1423 DO i=imin,imax
1424 mixing(ng) % diff2(i,j,itrc) = inival
1425 END DO
1426 END DO
1427# endif
1428# ifdef TS_DIF4
1429 DO itrc=1,nt(ng)
1430 DO i=imin,imax
1431 mixing(ng) % diff4(i,j,itrc) = inival
1432 END DO
1433 END DO
1434# endif
1435# ifdef DIFF_3DCOEF
1436 DO k=1,n(ng)
1437 DO i=imin,imax
1438 mixing(ng) % Hdiffusion(i,j) = inival
1439# ifdef TS_U3ADV_SPLIT
1440 mixing(ng) % diff3d_u(i,j,k) = inival
1441 mixing(ng) % diff3d_v(i,j,k) = inival
1442# else
1443 mixing(ng) % diff3d_r(i,j,k) = inival
1444# endif
1445 END DO
1446 END DO
1447# endif
1448 DO i=imin,imax
1449 mixing(ng) % Akv(i,j,0) = inival
1450 mixing(ng) % Akv(i,j,n(ng)) = inival
1451 END DO
1452 DO k=1,n(ng)-1
1453 DO i=imin,imax
1454 mixing(ng) % Akv(i,j,k) = akv_bak(ng)
1455 END DO
1456 END DO
1457 DO itrc=1,nat
1458 DO i=imin,imax
1459 mixing(ng) % Akt(i,j,0,itrc) = inival
1460 mixing(ng) % Akt(i,j,n(ng),itrc) = inival
1461 END DO
1462 DO k=1,n(ng)-1
1463 DO i=imin,imax
1464 mixing(ng) % Akt(i,j,k,itrc) = akt_bak(itrc,ng)
1465 END DO
1466 END DO
1467 END DO
1468# ifdef FLOAT_VWALK
1469 DO k=1,n(ng)
1470 DO i=imin,imax
1471 mixing(ng) % dAktdz(i,j,k) = inival
1472 END DO
1473 END DO
1474# endif
1475# if defined LMD_SKPP || defined LMD_BKPP || \
1476 defined bulk_fluxes || defined balance_operator
1477 DO i=imin,imax
1478 mixing(ng) % alpha(i,j) = inival
1479 mixing(ng) % beta(i,j) = inival
1480 END DO
1481# endif
1482# ifdef BV_FREQUENCY
1483 DO k=0,n(ng)
1484 DO i=imin,imax
1485 mixing(ng) % bvf(i,j,k) = inival
1486 END DO
1487 END DO
1488# endif
1489# if defined GLS_MIXING || defined MY25_MIXING
1490 DO k=0,n(ng)
1491 DO i=imin,imax
1492 mixing(ng) % tke(i,j,k,1) = gls_kmin(ng)
1493 mixing(ng) % tke(i,j,k,2) = gls_kmin(ng)
1494 mixing(ng) % tke(i,j,k,3) = gls_kmin(ng)
1495 mixing(ng) % gls(i,j,k,1) = gls_pmin(ng)
1496 mixing(ng) % gls(i,j,k,2) = gls_pmin(ng)
1497 mixing(ng) % gls(i,j,k,3) = gls_pmin(ng)
1498 mixing(ng) % Lscale(i,j,k) = inival
1499 END DO
1500 END DO
1501 DO i=imin,imax
1502 mixing(ng) % Akk(i,j,0) = inival
1503 mixing(ng) % Akk(i,j,n(ng)) = inival
1504# ifdef GLS_MIXING
1505 mixing(ng) % Akp(i,j,0) = inival
1506 mixing(ng) % Akp(i,j,n(ng)) = inival
1507# endif
1508 END DO
1509 DO k=1,n(ng)-1
1510 DO i=imin,imax
1511 mixing(ng) % Akk(i,j,k) = akk_bak(ng)
1512# ifdef GLS_MIXING
1513 mixing(ng) % Akp(i,j,k) = akp_bak(ng)
1514# endif
1515 END DO
1516 END DO
1517# endif
1518# if defined LMD_MIXING && defined LMD_DDMIX
1519 DO k=0,n(ng)
1520 DO i=imin,imax
1521 mixing(ng) % alfaobeta(i,j,k) = inival
1522 END DO
1523 END DO
1524# endif
1525# if defined LMD_SKPP || defined SOLAR_SOURCE
1526 DO i=imin,imax
1527 mixing(ng) % Jwtype(i,j) = real(lmd_jwt(ng),r8)
1528 END DO
1529# endif
1530# if defined LMD_SKPP || defined LMD_BKPP
1531 DO i=imin,imax
1532 mixing(ng) % ksbl(i,j) = 0
1533 mixing(ng) % hsbl(i,j) = inival
1534 END DO
1535# ifdef LMD_BKPP
1536 DO i=imin,imax
1537 mixing(ng) % kbbl(i,j) = 0
1538 mixing(ng) % hbbl(i,j) = inival
1539 END DO
1540# endif
1541# ifdef LMD_NONLOCAL
1542 DO itrc=1,nat
1543 DO k=0,n(ng)
1544 DO i=imin,imax
1545 mixing(ng) % ghats(i,j,k,itrc) = inival
1546 END DO
1547 END DO
1548 END DO
1549# endif
1550# endif
1551#endif
1552 END DO
1553
1554#ifdef FOUR_DVAR
1555!
1556! Spatial convolution diffusion coefficients.
1557!
1558 DO j=jmin,jmax
1559 DO i=imin,imax
1560 mixing(ng) % Kh(i,j) = 1.0_r8
1561 END DO
1562# ifdef SOLVE3D
1563 DO k=0,n(ng)
1564 DO i=imin,imax
1565 mixing(ng) % Kv(i,j,k) = 1.0_r8
1566 END DO
1567 END DO
1568# endif
1569 END DO
1570#endif
1571 END IF
1572
1573#if defined TANGENT || defined TL_IOMS
1574# ifdef SOLVE3D
1575!
1576! Tangent linear model state.
1577!
1578# if defined GLS_MIXING || defined MY25_MIXING
1579 IF (model.eq.irpm) THEN
1580 cff1=gls_kmin(ng)
1581 cff2=gls_pmin(ng)
1582 cff3=akk_bak(ng)
1583 cff4=akp_bak(ng)
1584 ELSE
1585 cff1=inival
1586 cff2=inival
1587 cff3=inival
1588 cff4=inival
1589 END IF
1590# endif
1591 IF ((model.eq.0).or.(model.eq.itlm).or.(model.eq.irpm)) THEN
1592 DO j=jmin,jmax
1593# ifdef DIFF_3DCOEF
1594 DO k=1,n(ng)
1595 DO i=imin,imax
1596# ifdef TS_U3ADV_SPLIT
1597 mixing(ng) % tl_diff3d_u(i,j,k) = inival
1598 mixing(ng) % tl_diff3d_v(i,j,k) = inival
1599# else
1600 mixing(ng) % tl_diff3d_r(i,j,k) = inival
1601# endif
1602 END DO
1603 END DO
1604# endif
1605# ifdef VISC_3DCOEF
1606 DO k=1,n(ng)
1607 DO i=imin,imax
1608# ifdef UV_U3ADV_SPLIT
1609 mixing(ng) % tl_Uvis3d_r(i,j,k) = inival
1610 mixing(ng) % tl_Vvis3d_r(i,j,k) = inival
1611# else
1612 mixing(ng) % tl_visc3d_r(i,j,k) = inival
1613# endif
1614 END DO
1615 END DO
1616# endif
1617 DO k=0,n(ng)
1618 DO i=imin,imax
1619 mixing(ng) % tl_Akv(i,j,k) = inival
1620 END DO
1621 END DO
1622 DO itrc=1,nat
1623 DO k=0,n(ng)
1624 DO i=imin,imax
1625 mixing(ng) % tl_Akt(i,j,k,itrc) = inival
1626 END DO
1627 END DO
1628 END DO
1629# if defined LMD_SKPP || defined LMD_BKPP || defined BULK_FLUXES
1630 DO i=imin,imax
1631 mixing(ng) % tl_alpha(i,j) = inival
1632 mixing(ng) % tl_beta(i,j) = inival
1633 END DO
1634# endif
1635# ifdef BV_FREQUENCY
1636 DO k=0,n(ng)
1637 DO i=imin,imax
1638 mixing(ng) % tl_bvf(i,j,k) = inival
1639 END DO
1640 END DO
1641# endif
1642# if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET
1643 DO k=0,n(ng)
1644 DO i=imin,imax
1645 mixing(ng) % tl_tke(i,j,k,1) = cff1
1646 mixing(ng) % tl_tke(i,j,k,2) = cff1
1647 mixing(ng) % tl_tke(i,j,k,3) = cff1
1648 mixing(ng) % tl_gls(i,j,k,1) = cff2
1649 mixing(ng) % tl_gls(i,j,k,2) = cff2
1650 mixing(ng) % tl_gls(i,j,k,3) = cff2
1651 mixing(ng) % tl_Lscale(i,j,k) = inival
1652 END DO
1653 END DO
1654 DO i=imin,imax
1655 mixing(ng) % tl_Akk(i,j,0) = inival
1656 mixing(ng) % tl_Akk(i,j,n(ng)) = inival
1657# ifdef GLS_MIXING_NOT_YET
1658 mixing(ng) % tl_Akp(i,j,0) = inival
1659 mixing(ng) % tl_Akp(i,j,n(ng)) = inival
1660# endif
1661 END DO
1662 DO k=1,n(ng)-1
1663 DO i=imin,imax
1664 mixing(ng) % tl_Akk(i,j,k) = cff3
1665# ifdef GLS_MIXING
1666 mixing(ng) % tl_Akp(i,j,k) = cff4
1667# endif
1668 END DO
1669 END DO
1670# endif
1671 END DO
1672 END IF
1673# endif
1674#endif
1675
1676#ifdef ADJOINT
1677# ifdef SOLVE3D
1678!
1679! Adjoint model state.
1680!
1681 IF ((model.eq.0).or.(model.eq.iadm)) THEN
1682 DO j=jmin,jmax
1683# ifdef DIFF_3DCOEF
1684 DO k=1,n(ng)
1685 DO i=imin,imax
1686# ifdef TS_U3ADV_SPLIT
1687 mixing(ng) % ad_diff3d_u(i,j,k) = inival
1688 mixing(ng) % ad_diff3d_v(i,j,k) = inival
1689# else
1690 mixing(ng) % ad_diff3d_r(i,j,k) = inival
1691# endif
1692 END DO
1693 END DO
1694# endif
1695# ifdef VISC_3DCOEF
1696 DO k=1,n(ng)
1697 DO i=imin,imax
1698# ifdef UV_U3ADV_SPLIT
1699 mixing(ng) % ad_Uvis3d_r(i,j,k) = inival
1700 mixing(ng) % ad_Vvis3d_r(i,j,k) = inival
1701# else
1702 mixing(ng) % ad_visc3d_r(i,j,k) = inival
1703# endif
1704 END DO
1705 END DO
1706# endif
1707 DO k=0,n(ng)
1708 DO i=imin,imax
1709 mixing(ng) % ad_Akv(i,j,k) = inival
1710 END DO
1711 END DO
1712 DO itrc=1,nat
1713 DO k=0,n(ng)
1714 DO i=imin,imax
1715 mixing(ng) % ad_Akt(i,j,k,itrc) = inival
1716 END DO
1717 END DO
1718 END DO
1719# if defined LMD_SKPP || defined LMD_BKPP || defined BULK_FLUXES
1720 DO i=imin,imax
1721 mixing(ng) % ad_alpha(i,j) = inival
1722 mixing(ng) % ad_beta(i,j) = inival
1723 END DO
1724# endif
1725# ifdef BV_FREQUENCY
1726 DO k=0,n(ng)
1727 DO i=imin,imax
1728 mixing(ng) % ad_bvf(i,j,k) = inival
1729 END DO
1730 END DO
1731# endif
1732# if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET
1733 DO k=0,n(ng)
1734 DO i=imin,imax
1735 mixing(ng) % ad_tke(i,j,k,1) = inival
1736 mixing(ng) % ad_tke(i,j,k,2) = inival
1737 mixing(ng) % ad_tke(i,j,k,3) = inival
1738 mixing(ng) % ad_gls(i,j,k,1) = inival
1739 mixing(ng) % ad_gls(i,j,k,2) = inival
1740 mixing(ng) % ad_gls(i,j,k,3) = inival
1741 mixing(ng) % ad_Lscale(i,j,k) = inival
1742 mixing(ng) % ad_Akk(i,j,k) = inival
1743# ifdef GLS_MIXING
1744 mixing(ng) % ad_Akp(i,j,k) = inival
1745# endif
1746 END DO
1747 END DO
1748# endif
1749 END DO
1750 END IF
1751# endif
1752#endif
1753
1754#if defined FORWARD_READ && defined FORWARD_MIXING && defined SOLVE3D && \
1755 (defined tangent || defined tl_ioms || defined adjoint)
1756!
1757! Latest two records of the nonlinear trajectory used to interpolate
1758! the background state in the tangent linear and adjoint models.
1759!
1760 IF (model.eq.0) THEN
1761 DO j=jmin,jmax
1762# ifdef FORWARD_MIXING
1763 DO k=0,n(ng)
1764 DO i=imin,imax
1765 mixing(ng) % AkvG(i,j,k,1) = inival
1766 mixing(ng) % AkvG(i,j,k,2) = inival
1767 END DO
1768 END DO
1769 DO itrc=1,nat
1770 DO k=0,n(ng)
1771 DO i=imin,imax
1772 mixing(ng) % AktG(i,j,k,1,itrc) = inival
1773 mixing(ng) % AktG(i,j,k,2,itrc) = inival
1774 END DO
1775 END DO
1776 END DO
1777# endif
1778# if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET
1779 DO k=0,n(ng)
1780 DO i=imin,imax
1781 mixing(ng) % tkeG(i,j,k,1) = inival
1782 mixing(ng) % tkeG(i,j,k,2) = inival
1783 mixing(ng) % glsG(i,j,k,1) = inival
1784 mixing(ng) % glsG(i,j,k,2) = inival
1785 mixing(ng) % LscaleG(i,j,k,1) = inival
1786 mixing(ng) % LscaleG(i,j,k,2) = inival
1787 mixing(ng) % AkkG(i,j,k,1) = inival
1788 mixing(ng) % AkkG(i,j,k,2) = inival
1789# ifdef GLS_MIXING_NOT_YET
1790 mixing(ng) % AkpG(i,j,k,1) = inival
1791 mixing(ng) % AkpG(i,j,k,2) = inival
1792# endif
1793 END DO
1794 END DO
1795# endif
1796# if defined LMD_MIXING_NOT_YET
1797 DO i=imin,imax
1798 mixing(ng) % hsblG(i,j,1) = inival
1799 mixing(ng) % hsblG(i,j,2) = inival
1800 END DO
1801# endif
1802# if defined LMD_BKPP_NOT_YET
1803 DO i=imin,imax
1804 mixing(ng) % hbblG(i,j,1) = inival
1805 mixing(ng) % hbblG(i,j,2) = inival
1806 END DO
1807# endif
1808# if defined LMD_NONLOCAL_NOT_YET
1809 DO itrc=1,nat
1810 DO k=0,n(ng)
1811 DO i=imin,imax
1812 mixing(ng) % ghatsG(i,j,0:n(ng),1,itrc) = inival
1813 mixing(ng) % ghatsG(i,j,0:n(ng),2,itrc) = inival
1814 END DO
1815 END DO
1816 END DO
1817# endif
1818 END DO
1819 END IF
1820#endif
1821
1822 RETURN
integer, parameter inlm
Definition mod_param.F:662
integer, parameter irpm
Definition mod_param.F:664
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
integer, parameter iadm
Definition mod_param.F:665
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, parameter itlm
Definition mod_param.F:663
integer, dimension(:), allocatable lmd_jwt
real(r8), dimension(:), allocatable akk_bak
real(r8), dimension(:,:), allocatable akt_bak
real(r8), dimension(:), allocatable akv_bak
real(r8), dimension(:), allocatable gls_pmin
real(r8), dimension(:), allocatable akp_bak
real(r8), dimension(:), allocatable gls_kmin

References mod_scalars::akk_bak, mod_scalars::akp_bak, mod_scalars::akt_bak, mod_scalars::akv_bak, mod_param::bounds, mod_param::domain, mod_scalars::gls_kmin, mod_scalars::gls_pmin, mod_param::iadm, mod_param::inlm, mod_param::irpm, mod_param::itlm, mod_scalars::lmd_jwt, mod_scalars::ltracersponge, mod_scalars::luvsponge, mixing, mod_param::n, mod_param::nat, mod_param::nt, and mod_kinds::r8.

Referenced by ad_initial(), rbl4dvar_mod::analysis_initialize(), i4dvar_mod::background_initialize(), rbl4dvar_mod::background_initialize(), i4dvar_mod::posterior_analysis_initialize(), mod_arrays::roms_initialize_arrays(), and rp_initial().

Here is the caller graph for this function:

Variable Documentation

◆ mixing

type (t_mixing), dimension(:), allocatable mod_mixing::mixing

Definition at line 399 of file mod_mixing.F.

399 TYPE (T_MIXING), allocatable :: MIXING(:)

Referenced by ad_balance_mod::ad_balance(), ad_bulk_flux_mod::ad_bulk_flux(), ad_convolution_mod::ad_convolution(), ad_get_data(), ad_lmd_swfrac_tile(), ad_pre_step3d_mod::ad_pre_step3d(), ad_rho_eos_mod::ad_rho_eos(), ad_rhs3d_mod::ad_rhs3d(), ad_set_avg_mod::ad_set_avg_tile(), ad_set_data_tile(), ad_step2d_mod::ad_step2d(), ad_step3d_t_mod::ad_step3d_t(), ad_step3d_uv_mod::ad_step3d_uv(), ad_t3dmix2_mod::ad_t3dmix2(), ad_t3dmix4_mod::ad_t3dmix4(), ad_uv3dmix2_mod::ad_uv3dmix2(), ad_uv3dmix4_mod::ad_uv3dmix4(), ad_wrt_his_mod::ad_wrt_his_nf90(), ad_wrt_his_mod::ad_wrt_his_pio(), allocate_mixing(), analytical_mod::ana_sponge_tile(), analytical_mod::ana_vmix(), analytical_mod::ana_wtype(), zeta_balance_mod::balance_ref(), bulk_flux_mod::bulk_flux(), bvf_mix_mod::bvf_mix(), deallocate_mixing(), get_grid_mod::get_grid_nf90(), get_grid_mod::get_grid_pio(), get_state_mod::get_state_nf90(), get_state_mod::get_state_pio(), gls_corstep_mod::gls_corstep(), gls_prestep_mod::gls_prestep(), ini_hmixcoef_mod::ini_hmixcoef(), ini_hmixcoef_mod::ini_hmixcoef_tile(), initialize_mixing(), lmd_bkpp(), lmd_vmix_mod::lmd_finish(), lmd_skpp(), lmd_swfrac_tile(), lmd_vmix_mod::lmd_vmix(), metrics_mod::metrics(), my25_corstep_mod::my25_corstep(), my25_prestep_mod::my25_prestep(), normalization_mod::normalization(), pre_step3d_mod::pre_step3d(), rho_eos_mod::rho_eos(), rhs3d_mod::rhs3d(), cmeps_roms_mod::roms_export(), esmf_roms_mod::roms_export(), rp_bulk_flux_mod::rp_bulk_flux(), rp_get_data(), rp_lmd_swfrac_tile(), rp_pre_step3d_mod::rp_pre_step3d(), rp_rho_eos_mod::rp_rho_eos(), rp_rhs3d_mod::rp_rhs3d(), rp_set_data_tile(), rp_step2d_mod::rp_step2d(), rp_step3d_t_mod::rp_step3d_t(), rp_step3d_uv_mod::rp_step3d_uv(), rp_t3dmix2_mod::rp_t3dmix2(), rp_t3dmix4_mod::rp_t3dmix4(), rp_uv3dmix2_mod::rp_uv3dmix2(), rp_uv3dmix4_mod::rp_uv3dmix4(), set_avg_mod::set_avg_tile(), state_join_mod::state_join_nf90(), state_join_mod::state_join_pio(), step2d_mod::step2d(), step3d_t_mod::step3d_t(), step3d_uv_mod::step3d_uv(), t3dmix2_mod::t3dmix2(), t3dmix4_mod::t3dmix4(), tkebc_mod::tkebc(), tl_balance_mod::tl_balance(), tl_bulk_flux_mod::tl_bulk_flux(), tl_convolution_mod::tl_convolution(), tl_get_data(), tl_lmd_swfrac_tile(), tl_pre_step3d_mod::tl_pre_step3d(), tl_rho_eos_mod::tl_rho_eos(), tl_rhs3d_mod::tl_rhs3d(), tl_set_avg_mod::tl_set_avg_tile(), tl_set_data_tile(), tl_step2d_mod::tl_step2d(), tl_step3d_t_mod::tl_step3d_t(), tl_step3d_uv_mod::tl_step3d_uv(), tl_t3dmix2_mod::tl_t3dmix2(), tl_t3dmix4_mod::tl_t3dmix4(), tl_uv3dmix2_mod::tl_uv3dmix2(), tl_uv3dmix4_mod::tl_uv3dmix4(), tl_wrt_his_mod::tl_wrt_his_nf90(), tl_wrt_his_mod::tl_wrt_his_pio(), uv3dmix2_mod::uv3dmix2(), uv3dmix4_mod::uv3dmix4(), vwalk_floats_mod::vwalk_floats_tile(), wrt_dai_mod::wrt_dai_nf90(), wrt_dai_mod::wrt_dai_pio(), wrt_his_mod::wrt_his_nf90(), wrt_his_mod::wrt_his_pio(), wrt_ini_mod::wrt_ini_nf90(), wrt_ini_mod::wrt_ini_pio(), wrt_quick_mod::wrt_quick_nf90(), wrt_quick_mod::wrt_quick_pio(), wrt_rst_mod::wrt_rst_nf90(), wrt_rst_mod::wrt_rst_pio(), wrt_station_mod::wrt_station_nf90(), and wrt_station_mod::wrt_station_pio().