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

Functions/Subroutines

subroutine, public ad_obc_adjust (ng, tile, linp)
 
subroutine ad_obc_adjust_tile (ng, tile, lbi, ubi, lbj, ubj, lbij, ubij, imins, imaxs, jmins, jmaxs, linp)
 
subroutine, public ad_obc2d_adjust (ng, tile, linp)
 
subroutine ad_obc2d_adjust_tile (ng, tile, lbi, ubi, lbj, ubj, lbij, ubij, imins, imaxs, jmins, jmaxs, linp, umask, vmask, hz, hz_bry, ad_hz, ad_hz_bry)
 

Function/Subroutine Documentation

◆ ad_obc2d_adjust()

subroutine, public ad_obc_adjust_mod::ad_obc2d_adjust ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) linp )

Definition at line 722 of file ad_obc_adjust.F.

723!***********************************************************************
724!
725 USE mod_param
726 USE mod_grid
727!
728! Imported variable declarations.
729!
730 integer, intent(in) :: ng, tile, Linp
731!
732! Local variable declarations.
733!
734 character (len=*), parameter :: MyFile = &
735 & __FILE__//", ad_obc2d_adjust"
736!
737# include "tile.h"
738!
739# ifdef PROFILE
740 CALL wclock_on (ng, iadm, 7, __line__, myfile)
741# endif
742 CALL ad_obc2d_adjust_tile (ng, tile, &
743 & lbi, ubi, lbj, ubj, lbij, ubij, &
744 & imins, imaxs, jmins, jmaxs, &
745 & linp, &
746# ifdef MASKING
747 & grid(ng) % umask, &
748 & grid(ng) % vmask, &
749# endif
750 & grid(ng) % Hz, &
751 & grid(ng) % Hz_bry, &
752 & grid(ng) % ad_Hz, &
753 & grid(ng) % ad_Hz_bry)
754# ifdef PROFILE
755 CALL wclock_off (ng, iadm, 7, __line__, myfile)
756# endif
757
758 RETURN
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
integer, parameter iadm
Definition mod_param.F:665
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3

References ad_obc2d_adjust_tile(), mod_grid::grid, mod_param::iadm, wclock_off(), and wclock_on().

Referenced by ad_main3d().

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

◆ ad_obc2d_adjust_tile()

subroutine ad_obc_adjust_mod::ad_obc2d_adjust_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) linp,
real(r8), dimension(lbi:,lbj:), intent(in) umask,
real(r8), dimension(lbi:,lbj:), intent(in) vmask,
real(r8), dimension(lbi:,lbj:,:), intent(in) hz,
real(r8), dimension(lbij:,:,:), intent(in) hz_bry,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_hz,
real(r8), dimension(lbij:,:,:), intent(inout) ad_hz_bry )
private

Definition at line 762 of file ad_obc_adjust.F.

771!***********************************************************************
772!
773 USE mod_param
774 USE mod_boundary
775 USE mod_ncparam
776 USE mod_scalars
777!
778! Imported variable declarations.
779!
780 integer, intent(in) :: ng, tile
781 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
782 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
783 integer, intent(in) :: Linp
784!
785# ifdef ASSUMED_SHAPE
786# ifdef MASKING
787 real(r8), intent(in) :: umask(LBi:,LBj:)
788 real(r8), intent(in) :: vmask(LBi:,LBj:)
789# endif
790 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
791 real(r8), intent(in) :: Hz_bry(LBij:,:,:)
792 real(r8), intent(inout) :: ad_Hz(LBi:,LBj:,:)
793 real(r8), intent(inout) :: ad_Hz_bry(LBij:,:,:)
794# else
795# ifdef MASKING
796 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
797 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
798# endif
799 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
800 real(r8), intent(in) :: Hz_bry(LBij:UBij,N(ng),4)
801 real(r8), intent(inout) :: ad_Hz(LBi:UBi,LBj:UBj,N(ng))
802 real(r8), intent(inout) :: ad_Hz_bry(LBij:UBij,N(ng),4)
803# endif
804!
805! Local variable declarations.
806!
807 integer :: i, ib, it1, it2, j, k
808 real(r8) :: adfac, fac, fac1, fac2
809 real(r8) :: cff1, ad_cff1, ad_cff2
810
811 real(r8), dimension(0:N(ng)) :: CF
812 real(r8), dimension(0:N(ng)) :: DC
813
814 real(r8), dimension(0:N(ng)) :: ad_CF
815 real(r8), dimension(0:N(ng)) :: ad_DC
816
817# include "set_bounds.h"
818!
819!-----------------------------------------------------------------------
820! Adjoint of adjust open boundary fields with 4DVAR increments.
821!-----------------------------------------------------------------------
822!
823! Set time records and interpolation factor, if any.
824!
825 IF (nbrec(ng).eq.1) THEN
826 it1=1
827 it2=1
828 fac1=1.0_r8
829 fac2=0.0_r8
830 ELSE
831# ifdef GENERIC_DSTART
832 it1=max(0,(iic(ng)-ntstart(ng))/nobc(ng))+1
833# else
834 it1=max(0,(iic(ng)-1)/nobc(ng))+1
835# endif
836 it2=min(it1+1,nbrec(ng))
837 fac1=obc_time(it2,ng)-(time(ng)+dt(ng))
838 fac2=(time(ng)+dt(ng))-obc_time(it1,ng)
839 fac=1.0_r8/(fac1+fac2)
840 fac1=fac*fac1
841 fac2=fac*fac2
842 END IF
843!
844! Clear arrays and constants
845!
846 ad_cff1=0.0_r8
847 ad_cff2=0.0_r8
848 ad_cf=0.0_r8
849 ad_dc=0.0_r8
850!
851! 2D U-momentum open boundaries: integrate 3D U-momentum at the open
852! boundaries.
853!
854 IF (ad_lbc(iwest,isubar,ng)%acquire.and. &
855 & lobc(iwest,isubar,ng).and. &
856 & domain(ng)%Western_Edge(tile)) THEN
857 i=bounds(ng)%edge(iwest,r2dvar)
858 DO j=jstr,jend
859 dc(0)=0.0_r8
860 cf(0)=0.0_r8
861 DO k=1,n(ng)
862 dc(k)=0.5_r8*(hz_bry(j,k,iwest)+ &
863 & hz(i+1,j,k))
864 dc(0)=dc(0)+dc(k)
865 cf(0)=cf(0)+dc(k)*boundary(ng)%u_west(j,k)
866 END DO
867 cff1=1.0_r8/dc(0)
868!^ BOUNDARY(ng)%tl_ubar_west(j)=tl_cff2
869!^
870 ad_cff2=boundary(ng)%ad_ubar_west(j)
871 boundary(ng)%ad_ubar_west(j)=0.0_r8
872# ifdef MASKING
873!^ tl_cff2=tl_cff2*umask(i,j)
874!^
875 ad_cff2=ad_cff2*umask(i,j)
876# endif
877!^ tl_cff2=tl_CF(0)*cff1+CF(0)*tl_cff1
878!^
879 ad_cff1=ad_cff1+cf(0)*ad_cff2
880 ad_cf(0)=ad_cf(0)+ad_cff2*cff1
881 ad_cff2=0.0_r8
882!^ tl_cff1=-cff1*cff1*tl_DC(0)
883!^
884 ad_dc(0)=ad_dc(0)-cff1*cff1*ad_cff1
885 ad_cff1=0.0_r8
886 DO k=n(ng),1,-1
887!^ tl_CF(0)=tl_CF(0)+ &
888!^ & tl_DC(k)*BOUNDARY(ng)%u_west(j,k)+ &
889!^ & DC(k)*BOUNDARY(ng)%tl_u_west(j,k)
890!^
891 ad_dc(k)=ad_dc(k)+boundary(ng)%u_west(j,k)*ad_cf(0)
892 boundary(ng)%ad_u_west(j,k)=boundary(ng)%ad_u_west(j,k)+ &
893 & dc(k)*ad_cf(0)
894!^ tl_DC(0)=tl_DC(0)+tl_DC(k)
895!^
896 ad_dc(k)=ad_dc(k)+ad_dc(0)
897!^ tl_DC(k)=0.5_r8*(tl_Hz_bry(j,k,iwest)+ &
898!^ & tl_Hz(i+1,j,k))
899!^
900 adfac=0.5_r8*ad_dc(k)
901 ad_hz(i+1,j,k)=ad_hz(i+1,j,k)+adfac
902 ad_hz_bry(j,k,iwest)=ad_hz_bry(j,k,iwest)+adfac
903 ad_dc(k)=0.0_r8
904 END DO
905!^ tl_DC(0)=0.0_r8
906!^
907 ad_dc(0)=0.0_r8
908!^ tl_CF(0)=0.0_r8
909!^
910 ad_cf(0)=0.0_r8
911 END DO
912 END IF
913!
914 IF (ad_lbc(ieast,isubar,ng)%acquire.and. &
915 & lobc(ieast,isubar,ng).and. &
916 & domain(ng)%Eastern_Edge(tile)) THEN
917 i=bounds(ng)%edge(ieast,r2dvar)
918 DO j=jstr,jend
919 dc(0)=0.0_r8
920 cf(0)=0.0_r8
921 DO k=1,n(ng)
922 dc(k)=0.5_r8*(hz(i-1,j,k)+ &
923 & hz_bry(j,k,ieast))
924 dc(0)=dc(0)+dc(k)
925 cf(0)=cf(0)+dc(k)*boundary(ng)%u_east(j,k)
926 END DO
927 cff1=1.0_r8/dc(0)
928!^ BOUNDARY(ng)%tl_ubar_east(j)=tl_cff2
929!^
930 ad_cff2=ad_cff2+boundary(ng)%ad_ubar_east(j)
931 boundary(ng)%ad_ubar_east(j)=0.0_r8
932# ifdef MASKING
933!^ tl_cff2=tl_cff2*umask(i,j)
934!^
935 ad_cff2=ad_cff2*umask(i,j)
936# endif
937!^ tl_cff2=tl_CF(0)*cff1+CF(0)*tl_cff1
938!^
939 ad_cf(0)=ad_cf(0)+cff1*ad_cff2
940 ad_cff1=ad_cff1+cf(0)*ad_cff2
941 ad_cff2=0.0_r8
942!^ tl_cff1=-cff1*cff1*tl_DC(0)
943!^
944 ad_dc(0)=ad_dc(0)-cff1*cff1*ad_cff1
945 ad_cff1=0.0_r8
946 DO k=n(ng),1,-1
947!^ tl_CF(0)=tl_CF(0)+ &
948!^ & tl_DC(k)*BOUNDARY(ng)%u_east(j,k)+ &
949!^ & DC(k)*BOUNDARY(ng)%tl_u_east(j,k)
950!^
951 ad_dc(k)=ad_dc(k)+boundary(ng)%u_east(j,k)*ad_cf(0)
952 boundary(ng)%ad_u_east(j,k)=boundary(ng)%ad_u_east(j,k)+ &
953 & dc(k)*ad_cf(0)
954!^ tl_DC(0)=tl_DC(0)+tl_DC(k)
955!^
956 ad_dc(k)=ad_dc(k)+ad_dc(0)
957!^ tl_DC(k)=0.5_r8*(tl_Hz(i-1,j,k)+ &
958!^ & tl_Hz_bry(j,k,ieast))
959!^
960 adfac=0.5_r8*ad_dc(k)
961 ad_hz(i-1,j,k)=ad_hz(i-1,j,k)+adfac
962 ad_hz_bry(j,k,ieast)=ad_hz_bry(j,k,ieast)+adfac
963 ad_dc(k)=0.0_r8
964 END DO
965!^ tl_DC(0)=0.0_r8
966!^
967 ad_dc(0)=0.0_r8
968!^ tl_CF(0)=0.0_r8
969!^
970 ad_cf(0)=0.0_r8
971 END DO
972 END IF
973!
974 IF (ad_lbc(isouth,isubar,ng)%acquire.and. &
975 & lobc(isouth,isubar,ng).and. &
976 & domain(ng)%Southern_Edge(tile)) THEN
977 j=bounds(ng)%edge(isouth,r2dvar)
978 DO i=istr,iend
979 dc(0)=0.0_r8
980 cf(0)=0.0_r8
981 DO k=1,n(ng)
982 dc(k)=0.5_r8*(hz_bry(i-1,k,isouth)+ &
983 & hz_bry(i ,k,isouth))
984 dc(0)=dc(0)+dc(k)
985 cf(0)=cf(0)+dc(k)*boundary(ng)%u_south(i,k)
986 END DO
987 cff1=1.0_r8/dc(0)
988!^ BOUNDARY(ng)%tl_ubar_south(i)=tl_cff2
989!^
990 ad_cff2=ad_cff2+boundary(ng)%ad_ubar_south(i)
991 boundary(ng)%ad_ubar_south(i)=0.0_r8
992# ifdef MASKING
993!^ tl_cff2=tl_cff2*umask(i,j)
994!^
995 ad_cff2=ad_cff2*umask(i,j)
996# endif
997!^ tl_cff2=tl_CF(0)*cff1+CF(0)*tl_cff1
998!^
999 ad_cff1=ad_cff1+cf(0)*ad_cff2
1000 ad_cf(0)=ad_cf(0)+cff1*ad_cff2
1001 ad_cff2=0.0_r8
1002!^ tl_cff1=-cff1*cff1*tl_DC(0)
1003!^
1004 ad_dc(0)=ad_dc(0)-cff1*cff1*ad_cff1
1005 ad_cff1=0.0_r8
1006 DO k=n(ng),1,-1
1007!^ tl_CF(0)=tl_CF(0)+ &
1008!^ & tl_DC(k)*BOUNDARY(ng)%u_south(i,k)+ &
1009!^ & DC(k)*BOUNDARY(ng)%tl_u_south(i,k)
1010!^
1011 ad_dc(k)=ad_dc(k)+boundary(ng)%u_south(i,k)*ad_cf(0)
1012 boundary(ng)%ad_u_south(i,k)=boundary(ng)%ad_u_south(i,k)+ &
1013 & dc(k)*ad_cf(0)
1014!^ tl_DC(0)=tl_DC(0)+tl_DC(k)
1015!^
1016 ad_dc(k)=ad_dc(k)+ad_dc(0)
1017!^ tl_DC(k)=0.5_r8*(tl_Hz_bry(i-1,k,isouth)+ &
1018!^ & tl_Hz_bry(i ,k,isouth))
1019!^
1020 adfac=0.5_r8*ad_dc(k)
1021 ad_hz_bry(i-1,k,isouth)=ad_hz_bry(i-1,k,isouth)+adfac
1022 ad_hz_bry(i ,k,isouth)=ad_hz_bry(i ,k,isouth)+adfac
1023 ad_dc(k)=0.0_r8
1024 END DO
1025!^ tl_DC(0)=0.0_r8
1026!^
1027 ad_dc(0)=0.0_r8
1028!^ tl_CF(0)=0.0_r8
1029!^
1030 ad_cf(0)=0.0_r8
1031 END DO
1032 END IF
1033!
1034 IF (ad_lbc(inorth,isubar,ng)%acquire.and. &
1035 & lobc(inorth,isubar,ng).and. &
1036 & domain(ng)%Northern_Edge(tile)) THEN
1037 j=bounds(ng)%edge(inorth,r2dvar)
1038 DO i=istr,iend
1039 dc(0)=0.0_r8
1040 cf(0)=0.0_r8
1041 DO k=1,n(ng)
1042 dc(k)=0.5_r8*(hz_bry(i-1,k,inorth)+ &
1043 & hz_bry(i ,k,inorth))
1044 dc(0)=dc(0)+dc(k)
1045 cf(0)=cf(0)+dc(k)*boundary(ng)%u_north(i,k)
1046 END DO
1047 cff1=1.0_r8/dc(0)
1048!^ BOUNDARY(ng)%tl_ubar_north(i)=tl_cff2
1049!^
1050 ad_cff2=ad_cff2+boundary(ng)%ad_ubar_north(i)
1051 boundary(ng)%ad_ubar_north(i)=0.0_r8
1052# ifdef MASKING
1053!^ tl_cff2=tl_cff2*umask(i,j)
1054!^
1055 ad_cff2=ad_cff2*umask(i,j)
1056# endif
1057!^ tl_cff2=tl_CF(0)*cff1+CF(0)*tl_cff1
1058!^
1059 ad_cf(0)=ad_cf(0)+cff1*ad_cff2
1060 ad_cff1=ad_cff1+cf(0)*ad_cff2
1061 ad_cff2=0.0_r8
1062!^ tl_cff1=-cff1*cff1*tl_DC(0)
1063!^
1064 ad_dc(0)=ad_dc(0)-cff1*cff1*ad_cff1
1065 ad_cff1=0.0_r8
1066 DO k=n(ng),1,-1
1067!^ tl_CF(0)=tl_CF(0)+ &
1068!^ & tl_DC(k)*BOUNDARY(ng)%u_north(i,k)+ &
1069!^ & DC(k)*BOUNDARY(ng)%tl_u_north(i,k)
1070!^
1071 ad_dc(k)=ad_dc(k)+boundary(ng)%u_north(i,k)*ad_cf(0)
1072 boundary(ng)%ad_u_north(i,k)=boundary(ng)%ad_u_north(i,k)+ &
1073 & dc(k)*ad_cf(0)
1074 ad_cf(0)=0.0_r8
1075!^ tl_DC(0)=tl_DC(0)+tl_DC(k)
1076!^
1077 ad_dc(k)=ad_dc(k)+ad_dc(0)
1078!^ tl_DC(k)=0.5_r8*(tl_Hz_bry(i-1,k,inorth)+ &
1079!^ & tl_Hz_bry(i ,k,inorth))
1080!^
1081 adfac=0.5_r8*ad_dc(k)
1082 ad_hz_bry(i-1,k,inorth)=ad_hz_bry(i-1,k,inorth)+adfac
1083 ad_hz_bry(i ,k,inorth)=ad_hz_bry(i ,k,inorth)+adfac
1084 ad_dc(k)=0.0_r8
1085 END DO
1086!^ tl_DC(0)=0.0_r8
1087!^
1088 ad_dc(0)=0.0_r8
1089!^ tl_CF(0)=0.0_r8
1090!^
1091 ad_cf(0)=0.0_r8
1092 END DO
1093 END IF
1094!
1095! 3D V-momentum open boundaries: integrate 3D V-momentum at the open
1096! boundaries.
1097!
1098 IF (ad_lbc(iwest,isvbar,ng)%acquire.and. &
1099 & lobc(iwest,isvbar,ng).and. &
1100 & domain(ng)%Western_Edge(tile)) THEN
1101 i=bounds(ng)%edge(iwest,r2dvar)
1102 DO j=jstrv,jend
1103 dc(0)=0.0_r8
1104 cf(0)=0.0_r8
1105 DO k=1,n(ng)
1106 dc(k)=0.5_r8*(hz_bry(j-1,k,iwest)+ &
1107 & hz_bry(j ,k,iwest))
1108 dc(0)=dc(0)+dc(k)
1109 cf(0)=cf(0)+dc(k)*boundary(ng)%v_west(j,k)
1110 END DO
1111 cff1=1.0_r8/dc(0)
1112!^ BOUNDARY(ng)%tl_vbar_west(j)=tl_cff2
1113!^
1114 ad_cff2=ad_cff2+boundary(ng)%ad_vbar_west(j)
1115 boundary(ng)%ad_vbar_west(j)=0.0_r8
1116# ifdef MASKING
1117!^ tl_cff2=tl_cff2*vmask(i,j)
1118!^
1119 ad_cff2=ad_cff2*vmask(i,j)
1120# endif
1121!^ tl_cff2=tl_CF(0)*cff1+CF(0)*tl_cff1
1122!^
1123 ad_cf(0)=ad_cf(0)+cff1*ad_cff2
1124 ad_cff1=ad_cff1+cf(0)*ad_cff2
1125 ad_cff2=0.0_r8
1126!^ tl_cff1=-cff1*cff1*tl_DC(0)
1127!^
1128 ad_dc(0)=ad_dc(0)-cff1*cff1*ad_cff1
1129 ad_cff1=0.0_r8
1130 DO k=n(ng),1,-1
1131!^ tl_CF(0)=tl_CF(0)+ &
1132!^ & tl_DC(k)*BOUNDARY(ng)%v_west(j,k)+ &
1133!^ & DC(k)*BOUNDARY(ng)%tl_v_west(j,k)
1134!^
1135 boundary(ng)%ad_v_west(j,k)=boundary(ng)%ad_v_west(j,k)+ &
1136 & dc(k)*ad_cf(0)
1137 ad_dc(k)=ad_dc(k)+boundary(ng)%v_west(j,k)*ad_cf(0)
1138!^ tl_DC(0)=tl_DC(0)+tl_DC(k)
1139!^
1140 ad_dc(k)=ad_dc(k)+ad_dc(0)
1141!^ tl_DC(k)=0.5_r8*(tl_Hz_bry(j-1,k,iwest)+ &
1142!^ & tl_Hz_bry(j ,k,iwest))
1143!^
1144 adfac=0.5_r8*ad_dc(k)
1145 ad_hz_bry(j-1,k,iwest)=ad_hz_bry(j-1,k,iwest)+adfac
1146 ad_hz_bry(j ,k,iwest)=ad_hz_bry(j ,k,iwest)+adfac
1147 ad_dc(k)=0.0_r8
1148 END DO
1149!^ tl_CF(0)=0.0_r8
1150!^
1151 ad_cf(0)=0.0_r8
1152!^ tl_DC(0)=0.0_r8
1153!^
1154 ad_dc(0)=0.0_r8
1155 END DO
1156 END IF
1157!
1158 IF (ad_lbc(ieast,isvbar,ng)%acquire.and. &
1159 & lobc(ieast,isvbar,ng).and. &
1160 & domain(ng)%Eastern_Edge(tile)) THEN
1161 i=bounds(ng)%edge(ieast,r2dvar)
1162 DO j=jstrv,jend
1163 dc(0)=0.0_r8
1164 cf(0)=0.0_r8
1165 DO k=1,n(ng)
1166 dc(k)=0.5_r8*(hz_bry(j-1,k,ieast)+ &
1167 & hz_bry(j ,k,ieast))
1168 dc(0)=dc(0)+dc(k)
1169 cf(0)=cf(0)+dc(k)*boundary(ng)%v_east(j,k)
1170 END DO
1171 cff1=1.0_r8/dc(0)
1172!^ BOUNDARY(ng)%tl_vbar_east(j)=tl_cff2
1173!^
1174 ad_cff2=ad_cff2+boundary(ng)%ad_vbar_east(j)
1175 boundary(ng)%ad_vbar_east(j)=0.0_r8
1176# ifdef MASKING
1177!^ tl_cff2=tl_cff2*vmask(i,j)
1178!^
1179 ad_cff2=ad_cff2*vmask(i,j)
1180# endif
1181!^ tl_cff2=tl_CF(0)*cff1+CF(0)*tl_cff1
1182!^
1183 ad_cf(0)=ad_cf(0)+cff1*ad_cff2
1184 ad_cff1=ad_cff1+cf(0)*ad_cff2
1185 ad_cff2=0.0_r8
1186!^ tl_cff1=-cff1*cff1*tl_DC(0)
1187!^
1188 ad_dc(0)=ad_dc(0)-cff1*cff1*ad_cff1
1189 ad_cff1=0.0_r8
1190 DO k=n(ng),1,-1
1191!^ tl_CF(0)=tl_CF(0)+ &
1192!^ & tl_DC(k)*BOUNDARY(ng)%v_east(j,k)+ &
1193!^ & DC(k)*BOUNDARY(ng)%tl_v_east(j,k)
1194!^
1195 boundary(ng)%ad_v_east(j,k)=boundary(ng)%ad_v_east(j,k)+ &
1196 & dc(k)*ad_cf(0)
1197 ad_dc(k)=ad_dc(k)+boundary(ng)%v_east(j,k)*ad_cf(0)
1198!^ tl_DC(0)=tl_DC(0)+tl_DC(k)
1199!^
1200 ad_dc(k)=ad_dc(k)+ad_dc(0)
1201 ad_dc(0)=0.0_r8
1202!^ tl_DC(k)=0.5_r8*(tl_Hz_bry(j-1,k,ieast)+ &
1203!^ & tl_Hz_bry(j ,k,ieast))
1204!^
1205 adfac=0.5_r8*ad_dc(k)
1206 ad_hz_bry(j-1,k,ieast)=ad_hz_bry(j-1,k,ieast)+adfac
1207 ad_hz_bry(j ,k,ieast)=ad_hz_bry(j ,k,ieast)+adfac
1208 ad_dc(k)=0.0_r8
1209 END DO
1210!^ tl_DC(0)=0.0_r8
1211!^
1212 ad_dc(0)=0.0_r8
1213!^ tl_CF(0)=0.0_r8
1214!^
1215 ad_cf(0)=0.0_r8
1216 END DO
1217 END IF
1218!
1219 IF (ad_lbc(isouth,isvbar,ng)%acquire.and. &
1220 & lobc(isouth,isvbar,ng).and. &
1221 & domain(ng)%Southern_Edge(tile)) THEN
1222 j=bounds(ng)%edge(isouth,r2dvar)
1223 DO i=istr,iend
1224 dc(0)=0.0_r8
1225 cf(0)=0.0_r8
1226 DO k=1,n(ng)
1227 dc(k)=0.5_r8*(hz_bry(i,k,isouth)+ &
1228 & hz(i+1,j,k))
1229 dc(0)=dc(0)+dc(k)
1230 cf(0)=cf(0)+dc(k)*boundary(ng)%v_south(i,k)
1231 END DO
1232 cff1=1.0_r8/dc(0)
1233!^ BOUNDARY(ng)%tl_vbar_south(i)=tl_cff2
1234!^
1235 ad_cff2=ad_cff2+boundary(ng)%ad_vbar_south(i)
1236 boundary(ng)%ad_vbar_south(i)=0.0_r8
1237# ifdef MASKING
1238!^ tl_cff2=tl_cff2*vmask(i,j)
1239!^
1240 ad_cff2=ad_cff2*vmask(i,j)
1241# endif
1242!^ tl_cff2=tl_CF(0)*cff1+CF(0)*tl_cff1
1243!^
1244 ad_cf(0)=ad_cf(0)+cff1*ad_cff2
1245 ad_cff1=ad_cff1+cf(0)*ad_cff2
1246 ad_cff2=0.0_r8
1247!^ tl_cff1=-cff1*cff1*tl_DC(0)
1248!^
1249 ad_dc(0)=ad_dc(0)-cff1*cff1*ad_cff1
1250 ad_cff1=0.0_r8
1251 DO k=n(ng),1,-1
1252!^ tl_CF(0)=tl_CF(0)+ &
1253!^ & tl_DC(k)*BOUNDARY(ng)%v_south(i,k)+ &
1254!^ & DC(k)*BOUNDARY(ng)%tl_v_south(i,k)
1255!^
1256 boundary(ng)%ad_v_south(i,k)=boundary(ng)%ad_v_south(i,k)+ &
1257 & dc(k)*ad_cf(0)
1258 ad_dc(k)=ad_dc(k)+boundary(ng)%v_south(i,k)*ad_cf(0)
1259!^ tl_DC(0)=tl_DC(0)+tl_DC(k)
1260!^
1261 ad_dc(k)=ad_dc(k)+ad_dc(0)
1262!^ tl_DC(k)=0.5_r8*(tl_Hz_bry(i,k,isouth)+ &
1263!^ & tl_Hz(i+1,j,k))
1264!^
1265 adfac=0.5_r8*ad_dc(k)
1266 ad_hz_bry(i,k,isouth)=ad_hz_bry(i,k,isouth)+adfac
1267 ad_hz(i+1,j,k)=ad_hz(i+1,j,k)+adfac
1268 ad_dc(k)=0.0_r8
1269 END DO
1270!^ tl_DC(0)=0.0_r8
1271!^
1272 ad_dc(0)=0.0_r8
1273!^ tl_CF(0)=0.0_r8
1274!^
1275 ad_cf(0)=0.0_r8
1276 END DO
1277 END IF
1278!
1279 IF (ad_lbc(inorth,isvbar,ng)%acquire.and. &
1280 & lobc(inorth,isvbar,ng).and. &
1281 & domain(ng)%Northern_Edge(tile)) THEN
1282 j=bounds(ng)%edge(inorth,r2dvar)
1283 DO i=istr,iend
1284 dc(0)=0.0_r8
1285 cf(0)=0.0_r8
1286 DO k=1,n(ng)
1287 dc(k)=0.5_r8*(hz(i,j-1,k)+ &
1288 & hz_bry(i,k,inorth))
1289 dc(0)=dc(0)+dc(k)
1290 cf(0)=cf(0)+dc(k)*boundary(ng)%v_north(i,k)
1291 END DO
1292 cff1=1.0_r8/dc(0)
1293!^ BOUNDARY(ng)%tl_vbar_north(i)=tl_cff2
1294!^
1295 ad_cff2=ad_cff2+boundary(ng)%ad_vbar_north(i)
1296 boundary(ng)%ad_vbar_north(i)=0.0_r8
1297# ifdef MASKING
1298!^ tl_cff2=tl_cff2*vmask(i,j)
1299!^
1300 ad_cff2=ad_cff2*vmask(i,j)
1301# endif
1302!^ tl_cff2=tl_CF(0)*cff1+CF(0)*tl_cff1
1303!^
1304 ad_cf(0)=ad_cf(0)+cff1*ad_cff2
1305 ad_cff1=ad_cff1+cf(0)*ad_cff2
1306 ad_cff2=0.0_r8
1307!^ tl_cff1=-cff1*cff1*tl_DC(0)
1308!^
1309 ad_dc(0)=ad_dc(0)-cff1*cff1*ad_cff1
1310 ad_cff1=0.0_r8
1311 DO k=n(ng),1,-1
1312!^ tl_CF(0)=tl_CF(0)+ &
1313!^ & tl_DC(k)*BOUNDARY(ng)%v_north(i,k)+ &
1314!^ & DC(k)*BOUNDARY(ng)%tl_v_north(i,k)
1315!^
1316 boundary(ng)%ad_v_north(i,k)=boundary(ng)%ad_v_north(i,k)+ &
1317 & dc(k)*ad_cf(0)
1318 ad_dc(k)=ad_dc(k)+boundary(ng)%v_north(i,k)*ad_cf(0)
1319!^ tl_DC(0)=tl_DC(0)+tl_DC(k)
1320!^
1321 ad_dc(k)=ad_dc(k)+ad_dc(0)
1322!^ tl_DC(k)=0.5_r8*(tl_Hz(i,j-1,k)+ &
1323!^ & tl_Hz_bry(i,k,inorth))
1324!^
1325 adfac=0.5_r8*ad_dc(k)
1326 ad_hz_bry(i,k,inorth)=ad_hz_bry(i,k,inorth)+adfac
1327 ad_hz(i,j-1,k)=ad_hz(i,j-1,k)+adfac
1328 ad_dc(k)=0.0_r8
1329 END DO
1330!^ tl_DC(0)=0.0_r8
1331!^
1332 ad_dc(0)=0.0_r8
1333!^ tl_CF(0)=0.0_r8
1334!^
1335 ad_cf(0)=0.0_r8
1336 END DO
1337 END IF
1338!
1339 RETURN
type(t_boundary), dimension(:), allocatable boundary
integer isvbar
integer isubar
type(t_lbc), dimension(:,:,:), allocatable ad_lbc
Definition mod_param.F:378
integer, dimension(:), allocatable n
Definition mod_param.F:479
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, parameter r2dvar
Definition mod_param.F:717
integer, dimension(:), allocatable iic
real(dp), dimension(:), allocatable dt
integer, dimension(:), allocatable nobc
logical, dimension(:,:,:), allocatable lobc
integer, parameter iwest
real(dp), dimension(:,:), allocatable obc_time
integer, parameter isouth
integer, parameter ieast
real(dp), dimension(:), allocatable time
integer, dimension(:), allocatable ntstart
integer, parameter inorth
integer, dimension(:), allocatable nbrec

References mod_param::ad_lbc, mod_boundary::boundary, mod_param::bounds, mod_param::domain, mod_scalars::dt, mod_scalars::ieast, mod_scalars::iic, mod_scalars::inorth, mod_scalars::isouth, mod_ncparam::isubar, mod_ncparam::isvbar, mod_scalars::iwest, mod_scalars::lobc, mod_scalars::nbrec, mod_scalars::nobc, mod_scalars::ntstart, mod_scalars::obc_time, mod_param::r2dvar, and mod_scalars::time.

Referenced by ad_obc2d_adjust().

Here is the caller graph for this function:

◆ ad_obc_adjust()

subroutine, public ad_obc_adjust_mod::ad_obc_adjust ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) linp )

Definition at line 36 of file ad_obc_adjust.F.

37!***********************************************************************
38!
39 USE mod_param
40!
41! Imported variable declarations.
42!
43 integer, intent(in) :: ng, tile, Linp
44!
45! Local variable declarations.
46!
47 character (len=*), parameter :: MyFile = &
48 & __FILE__
49!
50# include "tile.h"
51!
52# ifdef PROFILE
53 CALL wclock_on (ng, iadm, 7, __line__, myfile)
54# endif
55 CALL ad_obc_adjust_tile (ng, tile, &
56 & lbi, ubi, lbj, ubj, lbij, ubij, &
57 & imins, imaxs, jmins, jmaxs, &
58 & linp)
59# ifdef PROFILE
60 CALL wclock_off (ng, iadm, 7, __line__, myfile)
61# endif
62!
63 RETURN

References ad_obc_adjust_tile(), mod_param::iadm, wclock_off(), and wclock_on().

Referenced by ad_main3d().

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

◆ ad_obc_adjust_tile()

subroutine ad_obc_adjust_mod::ad_obc_adjust_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) linp )
private

Definition at line 67 of file ad_obc_adjust.F.

71!***********************************************************************
72!
73 USE mod_param
74 USE mod_boundary
75 USE mod_ncparam
76 USE mod_scalars
77!
78! Imported variable declarations.
79!
80 integer, intent(in) :: ng, tile
81 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
82 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
83 integer, intent(in) :: Linp
84!
85! Local variable declarations.
86!
87 integer :: i, ib, it1, it2, j
88# ifdef SOLVE3D
89 integer :: it, k
90# endif
91 real(r8) :: adfac, fac, fac1, fac2
92
93# include "set_bounds.h"
94!
95!-----------------------------------------------------------------------
96! Adjoint of adjust open boundary fields with 4DVAR increments.
97!-----------------------------------------------------------------------
98!
99! Set time records and interpolation factor, if any.
100!
101 IF (nbrec(ng).eq.1) THEN
102 it1=1
103 it2=1
104 fac1=1.0_r8
105 fac2=0.0_r8
106 ELSE
107# ifdef GENERIC_DSTART
108 it1=max(0,(iic(ng)-ntstart(ng))/nobc(ng))+1
109# else
110 it1=max(0,(iic(ng)-1)/nobc(ng))+1
111# endif
112 it2=min(it1+1,nbrec(ng))
113 fac1=obc_time(it2,ng)-(time(ng)+dt(ng))
114 fac2=(time(ng)+dt(ng))-obc_time(it1,ng)
115 fac=1.0_r8/(fac1+fac2)
116 fac1=fac*fac1
117 fac2=fac*fac2
118 END IF
119
120# ifndef SOLVE3D
121!
122! 2D U-momentum open boundaries.
123!
124 IF (ad_lbc(iwest,isubar,ng)%acquire.and. &
125 & lobc(iwest,isubar,ng).and. &
126 & domain(ng)%Western_Edge(tile)) THEN
127 DO j=jstr,jend
128!^ BOUNDARY(ng)%tl_ubar_west(j)=fac1* &
129!^ & BOUNDARY(ng)%tl_ubar_obc(j, &
130!^ & iwest,it1,Linp)+ &
131!^ & fac2* &
132!^ & BOUNDARY(ng)%tl_ubar_obc(j, &
133!^ & iwest,it2,Linp)
134!^
135 boundary(ng)%ad_ubar_obc(j,iwest,it1,linp)= &
136 & boundary(ng)%ad_ubar_obc(j,iwest,it1,linp)+ &
137 & fac1* &
138 & boundary(ng)%ad_ubar_west(j)
139 boundary(ng)%ad_ubar_obc(j,iwest,it2,linp)= &
140 & boundary(ng)%ad_ubar_obc(j,iwest,it2,linp)+ &
141 & fac2* &
142 & boundary(ng)%ad_ubar_west(j)
143 boundary(ng)%ad_ubar_west(j)=0.0_r8
144 END DO
145 END IF
146
147 IF (ad_lbc(ieast,isubar,ng)%acquire.and. &
148 & lobc(ieast,isubar,ng).and. &
149 & domain(ng)%Eastern_Edge(tile)) THEN
150 DO j=jstr,jend
151!^ BOUNDARY(ng)%tl_ubar_east(j)=fac1* &
152!^ & BOUNDARY(ng)%tl_ubar_obc(j, &
153!^ & ieast,it1,Linp)+ &
154!^ & fac2* &
155!^ & BOUNDARY(ng)%tl_ubar_obc(j, &
156!^ & ieast,it2,Linp)
157!^
158 boundary(ng)%ad_ubar_obc(j,ieast,it1,linp)= &
159 & boundary(ng)%ad_ubar_obc(j,ieast,it1,linp)+ &
160 & fac1* &
161 & boundary(ng)%ad_ubar_east(j)
162 boundary(ng)%ad_ubar_obc(j,ieast,it2,linp)= &
163 & boundary(ng)%ad_ubar_obc(j,ieast,it2,linp)+ &
164 & fac2* &
165 & boundary(ng)%ad_ubar_east(j)
166 boundary(ng)%ad_ubar_east(j)=0.0_r8
167 END DO
168 END IF
169
170 IF (ad_lbc(isouth,isubar,ng)%acquire.and. &
171 & lobc(isouth,isubar,ng).and. &
172 & domain(ng)%Southern_Edge(tile)) THEN
173 DO i=istru,iend
174!^ BOUNDARY(ng)%tl_ubar_south(i)=fac1* &
175!^ & BOUNDARY(ng)%tl_ubar_obc(i, &
176!^ & isouth,it1,Linp)+ &
177!^ & fac2* &
178!^ & BOUNDARY(ng)%tl_ubar_obc(i, &
179!^ & isouth,it2,Linp)
180!^
181 boundary(ng)%ad_ubar_obc(i,isouth,it1,linp)= &
182 & boundary(ng)%ad_ubar_obc(i,isouth,it1,linp)+ &
183 & fac1* &
184 & boundary(ng)%ad_ubar_south(i)
185 boundary(ng)%ad_ubar_obc(i,isouth,it2,linp)= &
186 & boundary(ng)%ad_ubar_obc(i,isouth,it2,linp)+ &
187 & fac2* &
188 & boundary(ng)%ad_ubar_south(i)
189 boundary(ng)%ad_ubar_south(i)=0.0_r8
190 END DO
191 END IF
192
193 IF (ad_lbc(inorth,isubar,ng)%acquire.and. &
194 & lobc(inorth,isubar,ng).and. &
195 & domain(ng)%Northern_Edge(tile)) THEN
196 DO i=istru,iend
197!^ BOUNDARY(ng)%tl_ubar_north(i)=fac1* &
198!^ & BOUNDARY(ng)%tl_ubar_obc(i, &
199!^ & inorth,it1,Linp)+ &
200!^ & fac2* &
201!^ & BOUNDARY(ng)%tl_ubar_obc(i, &
202!^ & inorth,it2,Linp)
203!^
204 boundary(ng)%ad_ubar_obc(i,inorth,it1,linp)= &
205 & boundary(ng)%ad_ubar_obc(i,inorth,it1,linp)+ &
206 & fac1* &
207 & boundary(ng)%ad_ubar_north(i)
208 boundary(ng)%ad_ubar_obc(i,inorth,it2,linp)= &
209 & boundary(ng)%ad_ubar_obc(i,inorth,it2,linp)+ &
210 & fac2* &
211 & boundary(ng)%ad_ubar_north(i)
212 boundary(ng)%ad_ubar_north(i)=0.0_r8
213 END DO
214 END IF
215!
216! 2D V-momentum open boundaries.
217!
218 IF (ad_lbc(iwest,isvbar,ng)%acquire.and. &
219 & lobc(iwest,isvbar,ng).and. &
220 & domain(ng)%Western_Edge(tile)) THEN
221 DO j=jstrv,jend
222!^ BOUNDARY(ng)%tl_vbar_west(j)=fac1* &
223!^ & BOUNDARY(ng)%tl_vbar_obc(j, &
224!^ & iwest,it1,Linp)+ &
225!^ & fac2* &
226!^ & BOUNDARY(ng)%tl_vbar_obc(j, &
227!^ & iwest,it2,Linp)
228!^
229 boundary(ng)%ad_vbar_obc(j,iwest,it1,linp)= &
230 & boundary(ng)%ad_vbar_obc(j,iwest,it1,linp)+ &
231 & fac1* &
232 & boundary(ng)%ad_vbar_west(j)
233 boundary(ng)%ad_vbar_obc(j,iwest,it2,linp)= &
234 & boundary(ng)%ad_vbar_obc(j,iwest,it2,linp)+ &
235 & fac2* &
236 & boundary(ng)%ad_vbar_west(j)
237 boundary(ng)%ad_vbar_west(j)=0.0_r8
238 END DO
239 END IF
240
241 IF (ad_lbc(ieast,isvbar,ng)%acquire.and. &
242 & lobc(ieast,isvbar,ng).and. &
243 & domain(ng)%Eastern_Edge(tile)) THEN
244 DO j=jstrv,jend
245!^ BOUNDARY(ng)%tl_vbar_east(j)=fac1* &
246!^ & BOUNDARY(ng)%tl_vbar_obc(j, &
247!^ & ieast,it1,Linp)+ &
248!^ & fac2* &
249!^ & BOUNDARY(ng)%tl_vbar_obc(j, &
250!^ & ieast,it2,Linp)
251!^
252 boundary(ng)%ad_vbar_obc(j,ieast,it1,linp)= &
253 & boundary(ng)%ad_vbar_obc(j,ieast,it1,linp)+ &
254 & fac1* &
255 & boundary(ng)%ad_vbar_east(j)
256 boundary(ng)%ad_vbar_obc(j,ieast,it2,linp)= &
257 & boundary(ng)%ad_vbar_obc(j,ieast,it2,linp)+ &
258 & fac2* &
259 & boundary(ng)%ad_vbar_east(j)
260 boundary(ng)%ad_vbar_east(j)=0.0_r8
261 END DO
262 END IF
263
264 IF (ad_lbc(isouth,isvbar,ng)%acquire.and. &
265 & lobc(isouth,isvbar,ng).and. &
266 & domain(ng)%Southern_Edge(tile)) THEN
267 DO i=istr,iend
268!^ BOUNDARY(ng)%tl_vbar_south(i)=fac1* &
269!^ & BOUNDARY(ng)%tl_vbar_obc(i, &
270!^ & isouth,it1,Linp)+ &
271!^ & fac2* &
272!^ & BOUNDARY(ng)%tl_vbar_obc(i, &
273!^ & isouth,it2,Linp)
274!^
275 boundary(ng)%ad_vbar_obc(i,isouth,it1,linp)= &
276 & boundary(ng)%ad_vbar_obc(i,isouth,it1,linp)+ &
277 & fac1* &
278 & boundary(ng)%ad_vbar_south(i)
279 boundary(ng)%ad_vbar_obc(i,isouth,it2,linp)= &
280 & boundary(ng)%ad_vbar_obc(i,isouth,it2,linp)+ &
281 & fac2* &
282 & boundary(ng)%ad_vbar_south(i)
283 boundary(ng)%ad_vbar_south(i)=0.0_r8
284 END DO
285 END IF
286
287 IF (ad_lbc(inorth,isvbar,ng)%acquire.and. &
288 & lobc(inorth,isvbar,ng).and. &
289 & domain(ng)%Northern_Edge(tile)) THEN
290 DO i=istr,iend
291!^ BOUNDARY(ng)%tl_vbar_north(i)=fac1* &
292!^ & BOUNDARY(ng)%tl_vbar_obc(i, &
293!^ & inorth,it1,Linp)+ &
294!^ & fac2* &
295!^ & BOUNDARY(ng)%tl_vbar_obc(i, &
296!^ & inorth,it2,Linp)
297!^
298 boundary(ng)%ad_vbar_obc(i,inorth,it1,linp)= &
299 & boundary(ng)%ad_vbar_obc(i,inorth,it1,linp)+ &
300 & fac1* &
301 & boundary(ng)%ad_vbar_north(i)
302 boundary(ng)%ad_vbar_obc(i,inorth,it2,linp)= &
303 & boundary(ng)%ad_vbar_obc(i,inorth,it2,linp)+ &
304 & fac2* &
305 & boundary(ng)%ad_vbar_north(i)
306 boundary(ng)%ad_vbar_north(i)=0.0_r8
307 END DO
308 END IF
309# endif
310!
311! Free-surface open boundaries.
312!
313 IF (ad_lbc(iwest,isfsur,ng)%acquire.and. &
314 & lobc(iwest,isfsur,ng).and. &
315 & domain(ng)%Western_Edge(tile)) THEN
316 DO j=jstr,jend
317!^ BOUNDARY(ng)%tl_zeta_west(j)=fac1* &
318!^ & BOUNDARY(ng)%tl_zeta_obc(j, &
319!^ & iwest,it1,Linp)+ &
320!^ & fac2*
321!^ & BOUNDARY(ng)%tl_zeta_obc(j,
322!^ & iwest,it2,Linp)
323!^
324 boundary(ng)%ad_zeta_obc(j,iwest,it1,linp)= &
325 & boundary(ng)%ad_zeta_obc(j,iwest,it1,linp)+ &
326 & fac1* &
327 & boundary(ng)%ad_zeta_west(j)
328 boundary(ng)%ad_zeta_obc(j,iwest,it2,linp)= &
329 & boundary(ng)%ad_zeta_obc(j,iwest,it2,linp)+ &
330 & fac2* &
331 & boundary(ng)%ad_zeta_west(j)
332 boundary(ng)%ad_zeta_west(j)=0.0_r8
333 END DO
334 END IF
335
336 IF (ad_lbc(ieast,isfsur,ng)%acquire.and. &
337 & lobc(ieast,isfsur,ng).and. &
338 & domain(ng)%Eastern_Edge(tile)) THEN
339 DO j=jstr,jend
340!^ BOUNDARY(ng)%tl_zeta_east(j)=fac1* &
341!^ & BOUNDARY(ng)%tl_zeta_obc(j, &
342!^ & ieast,it1,Linp)+ &
343!^ & fac2* &
344!^ & BOUNDARY(ng)%tl_zeta_obc(j, &
345!^ & ieast,it2,Linp)
346!^
347 boundary(ng)%ad_zeta_obc(j,ieast,it1,linp)= &
348 & boundary(ng)%ad_zeta_obc(j,ieast,it1,linp)+ &
349 & fac1* &
350 & boundary(ng)%ad_zeta_east(j)
351 boundary(ng)%ad_zeta_obc(j,ieast,it2,linp)= &
352 & boundary(ng)%ad_zeta_obc(j,ieast,it2,linp)+ &
353 & fac2* &
354 & boundary(ng)%ad_zeta_east(j)
355 boundary(ng)%ad_zeta_east(j)=0.0_r8
356 END DO
357 END IF
358
359 IF (ad_lbc(isouth,isfsur,ng)%acquire.and. &
360 & lobc(isouth,isfsur,ng).and. &
361 & domain(ng)%Southern_Edge(tile)) THEN
362 DO i=istr,iend
363!^ BOUNDARY(ng)%tl_zeta_south(i)=fac1* &
364!^ & BOUNDARY(ng)%tl_zeta_obc(i, &
365!^ & isouth,it1,Linp)+ &
366!^ & fac2* &
367!^ & BOUNDARY(ng)%tl_zeta_obc(i, &
368!^ & isouth,it2,Linp)
369!^
370 boundary(ng)%ad_zeta_obc(i,isouth,it1,linp)= &
371 & boundary(ng)%ad_zeta_obc(i,isouth,it1,linp)+ &
372 & fac1* &
373 & boundary(ng)%ad_zeta_south(i)
374 boundary(ng)%ad_zeta_obc(i,isouth,it2,linp)= &
375 & boundary(ng)%ad_zeta_obc(i,isouth,it2,linp)+ &
376 & fac2* &
377 & boundary(ng)%ad_zeta_south(i)
378 boundary(ng)%ad_zeta_south(i)=0.0_r8
379 END DO
380 END IF
381
382 IF (ad_lbc(inorth,isfsur,ng)%acquire.and. &
383 & lobc(inorth,isfsur,ng).and. &
384 & domain(ng)%Northern_Edge(tile)) THEN
385 DO i=istr,iend
386!^ BOUNDARY(ng)%tl_zeta_north(i)=fac1* &
387!^ & BOUNDARY(ng)%tl_zeta_obc(i, &
388!^ & inorth,it1,Linp)+ &
389!^ & fac2* &
390!^ & BOUNDARY(ng)%tl_zeta_obc(i, &
391!^ & inorth,it2,Linp)
392!^
393 boundary(ng)%ad_zeta_obc(i,inorth,it1,linp)= &
394 & boundary(ng)%ad_zeta_obc(i,inorth,it1,linp)+ &
395 & fac1* &
396 & boundary(ng)%ad_zeta_north(i)
397 boundary(ng)%ad_zeta_obc(i,inorth,it2,linp)= &
398 & boundary(ng)%ad_zeta_obc(i,inorth,it2,linp)+ &
399 & fac2* &
400 & boundary(ng)%ad_zeta_north(i)
401 boundary(ng)%ad_zeta_north(i)=0.0_r8
402 END DO
403 END IF
404
405# ifdef SOLVE3D
406!
407! 3D U-momentum open boundaries.
408!
409 IF (ad_lbc(iwest,isuvel,ng)%acquire.and. &
410 & lobc(iwest,isuvel,ng).and. &
411 & domain(ng)%Western_Edge(tile)) THEN
412 DO k=1,n(ng)
413 DO j=jstr,jend
414!^ BOUNDARY(ng)%tl_u_west(j,k)=fac1* &
415!^ & BOUNDARY(ng)%tl_u_obc(j,k, &
416!^ & iwest,it1,Linp)+ &
417!^ & fac2* &
418!^ & BOUNDARY(ng)%tl_u_obc(j,k, &
419!^ & iwest,it2,Linp)
420!^
421 boundary(ng)%ad_u_obc(j,k,iwest,it1,linp)= &
422 & boundary(ng)%ad_u_obc(j,k,iwest,it1,linp)+ &
423 & fac1* &
424 & boundary(ng)%ad_u_west(j,k)
425 boundary(ng)%ad_u_obc(j,k,iwest,it2,linp)= &
426 & boundary(ng)%ad_u_obc(j,k,iwest,it2,linp)+ &
427 & fac2* &
428 & boundary(ng)%ad_u_west(j,k)
429 boundary(ng)%ad_u_west(j,k)=0.0_r8
430 END DO
431 END DO
432 END IF
433
434 IF (ad_lbc(ieast,isuvel,ng)%acquire.and. &
435 & lobc(ieast,isuvel,ng).and. &
436 & domain(ng)%Eastern_Edge(tile)) THEN
437 DO k=1,n(ng)
438 DO j=jstr,jend
439!^ BOUNDARY(ng)%tl_u_east(j,k)=fac1* &
440!^ & BOUNDARY(ng)%tl_u_obc(j,k, &
441!^ & ieast,it1,Linp)+ &
442!^ & fac2* &
443!^ & BOUNDARY(ng)%tl_u_obc(j,k, &
444!^ & ieast,it2,Linp)
445!^
446 boundary(ng)%ad_u_obc(j,k,ieast,it1,linp)= &
447 & boundary(ng)%ad_u_obc(j,k,ieast,it1,linp)+ &
448 & fac1* &
449 & boundary(ng)%ad_u_east(j,k)
450 boundary(ng)%ad_u_obc(j,k,ieast,it2,linp)= &
451 & boundary(ng)%ad_u_obc(j,k,ieast,it2,linp)+ &
452 & fac2* &
453 & boundary(ng)%ad_u_east(j,k)
454 boundary(ng)%ad_u_east(j,k)=0.0_r8
455 END DO
456 END DO
457 END IF
458
459 IF (ad_lbc(isouth,isuvel,ng)%acquire.and. &
460 & lobc(isouth,isuvel,ng).and. &
461 & domain(ng)%Southern_Edge(tile)) THEN
462 DO k=1,n(ng)
463 DO i=istru,iend
464!^ BOUNDARY(ng)%tl_u_south(i,k)=fac1* &
465!^ & BOUNDARY(ng)%tl_u_obc(i,k, &
466!^ & isouth,it1,Linp)+ &
467!^ & fac2* &
468!^ & BOUNDARY(ng)%tl_u_obc(i,k, &
469!^ & isouth,it2,Linp)
470!^
471 boundary(ng)%ad_u_obc(i,k,isouth,it1,linp)= &
472 & boundary(ng)%ad_u_obc(i,k,isouth,it1,linp)+ &
473 & fac1* &
474 & boundary(ng)%ad_u_south(i,k)
475 boundary(ng)%ad_u_obc(i,k,isouth,it2,linp)= &
476 & boundary(ng)%ad_u_obc(i,k,isouth,it2,linp)+ &
477 & fac2* &
478 & boundary(ng)%ad_u_south(i,k)
479 boundary(ng)%ad_u_south(i,k)=0.0_r8
480 END DO
481 END DO
482 END IF
483
484 IF (ad_lbc(inorth,isuvel,ng)%acquire.and. &
485 & lobc(inorth,isuvel,ng).and. &
486 & domain(ng)%Northern_Edge(tile)) THEN
487 DO k=1,n(ng)
488 DO i=istru,iend
489!^ BOUNDARY(ng)%tl_u_north(i,k)=fac1* &
490!^ & BOUNDARY(ng)%tl_u_obc(i,k, &
491!^ & inorth,it1,Linp)+ &
492!^ & fac2* &
493!^ & BOUNDARY(ng)%tl_u_obc(i,k, &
494!^ & inorth,it2,Linp)
495!^
496 boundary(ng)%ad_u_obc(i,k,inorth,it1,linp)= &
497 & boundary(ng)%ad_u_obc(i,k,inorth,it1,linp)+ &
498 & fac1* &
499 & boundary(ng)%ad_u_north(i,k)
500 boundary(ng)%ad_u_obc(i,k,inorth,it2,linp)= &
501 & boundary(ng)%ad_u_obc(i,k,inorth,it2,linp)+ &
502 & fac2* &
503 & boundary(ng)%ad_u_north(i,k)
504 boundary(ng)%ad_u_north(i,k)=0.0_r8
505 END DO
506 END DO
507 END IF
508!
509! 3D V-momentum open boundaries.
510!
511 IF (ad_lbc(iwest,isvvel,ng)%acquire.and. &
512 & lobc(iwest,isvvel,ng).and. &
513 & domain(ng)%Western_Edge(tile)) THEN
514 DO k=1,n(ng)
515 DO j=jstrv,jend
516!^ BOUNDARY(ng)%tl_v_west(j,k)=fac1* &
517!^ & BOUNDARY(ng)%tl_v_obc(j,k, &
518!^ & iwest,it1,Linp)+ &
519!^ & fac2* &
520!^ & BOUNDARY(ng)%tl_v_obc(j,k, &
521!^ & iwest,it2,Linp)
522!^
523 boundary(ng)%ad_v_obc(j,k,iwest,it1,linp)= &
524 & boundary(ng)%ad_v_obc(j,k,iwest,it1,linp)+ &
525 & fac1* &
526 & boundary(ng)%ad_v_west(j,k)
527 boundary(ng)%ad_v_obc(j,k,iwest,it2,linp)= &
528 & boundary(ng)%ad_v_obc(j,k,iwest,it2,linp)+ &
529 & fac2* &
530 & boundary(ng)%ad_v_west(j,k)
531 boundary(ng)%ad_v_west(j,k)=0.0_r8
532 END DO
533 END DO
534 END IF
535
536 IF (ad_lbc(ieast,isvvel,ng)%acquire.and. &
537 & lobc(ieast,isvvel,ng).and. &
538 & domain(ng)%Eastern_Edge(tile)) THEN
539 DO k=1,n(ng)
540 DO j=jstrv,jend
541!^ BOUNDARY(ng)%tl_v_east(j,k)=fac1* &
542!^ & BOUNDARY(ng)%tl_v_obc(j,k, &
543!^ & ieast,it1,Linp)+ &
544!^ & fac2* &
545!^ & BOUNDARY(ng)%tl_v_obc(j,k, &
546!^ & ieast,it2,Linp)
547!^
548 boundary(ng)%ad_v_obc(j,k,ieast,it1,linp)= &
549 & boundary(ng)%ad_v_obc(j,k,ieast,it1,linp)+ &
550 & fac1* &
551 & boundary(ng)%ad_v_east(j,k)
552 boundary(ng)%ad_v_obc(j,k,ieast,it2,linp)= &
553 & boundary(ng)%ad_v_obc(j,k,ieast,it2,linp)+ &
554 & fac2* &
555 & boundary(ng)%ad_v_east(j,k)
556 boundary(ng)%ad_v_east(j,k)=0.0_r8
557 END DO
558 END DO
559 END IF
560
561 IF (ad_lbc(isouth,isvvel,ng)%acquire.and. &
562 & lobc(isouth,isvvel,ng).and. &
563 & domain(ng)%Southern_Edge(tile)) THEN
564 DO k=1,n(ng)
565 DO i=istr,iend
566!^ BOUNDARY(ng)%tl_v_south(i,k)=fac1* &
567!^ & BOUNDARY(ng)%tl_v_obc(i,k, &
568!^ & isouth,it1,Linp)+ &
569!^ & fac2* &
570!^ & BOUNDARY(ng)%tl_v_obc(i,k, &
571!^ & isouth,it2,Linp)
572!^
573 boundary(ng)%ad_v_obc(i,k,isouth,it1,linp)= &
574 & boundary(ng)%ad_v_obc(i,k,isouth,it1,linp)+ &
575 & fac1* &
576 & boundary(ng)%ad_v_south(i,k)
577 boundary(ng)%ad_v_obc(i,k,isouth,it2,linp)= &
578 & boundary(ng)%ad_v_obc(i,k,isouth,it2,linp)+ &
579 & fac2* &
580 & boundary(ng)%ad_v_south(i,k)
581 boundary(ng)%ad_v_south(i,k)=0.0_r8
582 END DO
583 END DO
584 END IF
585
586 IF (ad_lbc(inorth,isvvel,ng)%acquire.and. &
587 & lobc(inorth,isvvel,ng).and. &
588 & domain(ng)%Northern_Edge(tile)) THEN
589 DO k=1,n(ng)
590 DO i=istr,iend
591!^ BOUNDARY(ng)%tl_v_north(i,k)=fac1* &
592!^ & BOUNDARY(ng)%tl_v_obc(i,k, &
593!^ & inorth,it1,Linp)+ &
594!^ & fac2* &
595!^ & BOUNDARY(ng)%tl_v_obc(i,k, &
596!^ & inorth,it2,Linp)
597!^
598 boundary(ng)%ad_v_obc(i,k,inorth,it1,linp)= &
599 & boundary(ng)%ad_v_obc(i,k,inorth,it1,linp)+ &
600 & fac1* &
601 & boundary(ng)%ad_v_north(i,k)
602 boundary(ng)%ad_v_obc(i,k,inorth,it2,linp)= &
603 & boundary(ng)%ad_v_obc(i,k,inorth,it2,linp)+ &
604 & fac2* &
605 & boundary(ng)%ad_v_north(i,k)
606 boundary(ng)%ad_v_north(i,k)=0.0_r8
607 END DO
608 END DO
609 END IF
610!
611! Tracers open boundaries.
612!
613 DO it=1,nt(ng)
614 IF (ad_lbc(iwest,istvar(it),ng)%acquire.and. &
615 & lobc(iwest,istvar(it),ng).and. &
616 & domain(ng)%Western_Edge(tile)) THEN
617 DO k=1,n(ng)
618 DO j=jstr,jend
619!^ BOUNDARY(ng)%tl_t_west(j,k,it)=fac1* &
620!^ & BOUNDARY(ng)%tl_t_obc(j, &
621!^ & k,iwest,it1,Linp,it)+ &
622!^ & fac2* &
623!^ & BOUNDARY(ng)%tl_t_obc(j, &
624!^ & k,iwest,it2,Linp,it)
625!^
626 boundary(ng)%ad_t_obc(j,k,iwest,it1,linp,it)= &
627 & boundary(ng)%ad_t_obc(j,k,iwest,it1,linp,it)+ &
628 & fac1* &
629 & boundary(ng)%ad_t_west(j,k,it)
630 boundary(ng)%ad_t_obc(j,k,iwest,it2,linp,it)= &
631 & boundary(ng)%ad_t_obc(j,k,iwest,it2,linp,it)+ &
632 & fac2* &
633 & boundary(ng)%ad_t_west(j,k,it)
634 boundary(ng)%ad_t_west(j,k,it)=0.0_r8
635 END DO
636 END DO
637 END IF
638
639 IF (ad_lbc(ieast,istvar(it),ng)%acquire.and. &
640 & lobc(ieast,istvar(it),ng).and. &
641 & domain(ng)%Eastern_Edge(tile)) THEN
642 DO k=1,n(ng)
643 DO j=jstr,jend
644!^ BOUNDARY(ng)%tl_t_east(j,k,it)=fac1* &
645!^ & BOUNDARY(ng)%tl_t_obc(j, &
646!^ & k,ieast,it1,Linp,it)+ &
647!^ & fac2* &
648!^ & BOUNDARY(ng)%tl_t_obc(j, &
649!^ & k,ieast,it2,Linp,it)
650!^
651 boundary(ng)%ad_t_obc(j,k,ieast,it1,linp,it)= &
652 & boundary(ng)%ad_t_obc(j,k,ieast,it1,linp,it)+ &
653 & fac1* &
654 & boundary(ng)%ad_t_east(j,k,it)
655 boundary(ng)%ad_t_obc(j,k,ieast,it2,linp,it)= &
656 & boundary(ng)%ad_t_obc(j,k,ieast,it2,linp,it)+ &
657 & fac2* &
658 & boundary(ng)%ad_t_east(j,k,it)
659 boundary(ng)%ad_t_east(j,k,it)=0.0_r8
660 END DO
661 END DO
662 END IF
663
664 IF (ad_lbc(isouth,istvar(it),ng)%acquire.and. &
665 & lobc(isouth,istvar(it),ng).and. &
666 & domain(ng)%Southern_Edge(tile)) THEN
667 DO k=1,n(ng)
668 DO i=istr,iend
669!^ BOUNDARY(ng)%tl_t_south(i,k,it)=fac1* &
670!^ & BOUNDARY(ng)%tl_t_obc(i, &
671!^ & k,isouth,it1,Linp,it)+ &
672!^ & fac2* &
673!^ & BOUNDARY(ng)%tl_t_obc(i, &
674!^ & k,isouth,it2,Linp,it)
675!^
676 boundary(ng)%ad_t_obc(i,k,isouth,it1,linp,it)= &
677 & boundary(ng)%ad_t_obc(i,k,isouth,it1,linp,it)+ &
678 & fac1* &
679 & boundary(ng)%ad_t_south(i,k,it)
680 boundary(ng)%ad_t_obc(i,k,isouth,it2,linp,it)= &
681 & boundary(ng)%ad_t_obc(i,k,isouth,it2,linp,it)+ &
682 & fac2* &
683 & boundary(ng)%ad_t_south(i,k,it)
684 boundary(ng)%ad_t_south(i,k,it)=0.0_r8
685 END DO
686 END DO
687 END IF
688
689 IF (ad_lbc(inorth,istvar(it),ng)%acquire.and. &
690 & lobc(inorth,istvar(it),ng).and. &
691 & domain(ng)%Northern_Edge(tile)) THEN
692 DO k=1,n(ng)
693 DO i=istr,iend
694!^ BOUNDARY(ng)%tl_t_north(i,k,it)=fac1*
695!^ & BOUNDARY(ng)%tl_t_obc(i, &
696!^ & k,inorth,it1,Linp,it)+ &
697!^ & fac2* &
698!^ & BOUNDARY(ng)%tl_t_obc(i, &
699!^ & k,inorth,it2,Linp,it)
700!^
701 boundary(ng)%ad_t_obc(i,k,inorth,it1,linp,it)= &
702 & boundary(ng)%ad_t_obc(i,k,inorth,it1,linp,it)+ &
703 & fac1* &
704 & boundary(ng)%ad_t_north(i,k,it)
705 boundary(ng)%ad_t_obc(i,k,inorth,it2,linp,it)= &
706 & boundary(ng)%ad_t_obc(i,k,inorth,it2,linp,it)+ &
707 & fac2* &
708 & boundary(ng)%ad_t_north(i,k,it)
709 boundary(ng)%ad_t_north(i,k,it)=0.0_r8
710 END DO
711 END DO
712 END IF
713 END DO
714# endif
715!
716 RETURN
integer isvvel
integer, dimension(:), allocatable istvar
integer isuvel
integer isfsur
integer, dimension(:), allocatable nt
Definition mod_param.F:489

References mod_param::ad_lbc, mod_boundary::boundary, mod_param::domain, mod_scalars::dt, mod_scalars::ieast, mod_scalars::iic, mod_scalars::inorth, mod_ncparam::isfsur, mod_scalars::isouth, mod_ncparam::istvar, mod_ncparam::isubar, mod_ncparam::isuvel, mod_ncparam::isvbar, mod_ncparam::isvvel, mod_scalars::iwest, mod_scalars::lobc, mod_param::n, mod_scalars::nbrec, mod_scalars::nobc, mod_param::nt, mod_scalars::ntstart, mod_scalars::obc_time, and mod_scalars::time.

Referenced by ad_obc_adjust().

Here is the caller graph for this function: