ROMS
Loading...
Searching...
No Matches
ad_lmd_swfrac.F File Reference
#include "cppdefs.h"
#include "set_bounds.h"
Include dependency graph for ad_lmd_swfrac.F:

Go to the source code of this file.

Functions/Subroutines

subroutine ad_lmd_swfrac_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, zscale, z, ad_z, ad_swdk)
 

Function/Subroutine Documentation

◆ ad_lmd_swfrac_tile()

subroutine ad_lmd_swfrac_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) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
real(r8), intent(in) zscale,
real(r8), dimension(imins:imaxs,jmins:jmaxs), intent(in) z,
real(r8), dimension(imins:imaxs,jmins:jmaxs), intent(inout) ad_z,
real(r8), dimension(imins:imaxs,jmins:jmaxs), intent(inout) ad_swdk )

Definition at line 6 of file ad_lmd_swfrac.F.

10!
11!git $Id$
12!================================================== Hernan G. Arango ===
13! Copyright (c) 2002-2025 The ROMS Group !
14! Licensed under a MIT/X style license !
15! See License_ROMS.md !
16!=======================================================================
17! !
18! This routine computes the adjoint fraction of solar shortwave flux !
19! penetrating to specified depth (times Zscale) due to exponential !
20! decay in Jerlov water type. !
21! !
22! On Input: !
23! !
24! Zscale Scale factor to apply to depth array. !
25! Z Vertical height (meters, negative) for !
26! desired solar short-wave fraction. !
27! ad_Z Adjoint vertical height for !
28! desired solar short-wave fraction. !
29! !
30! On Output: !
31! !
32! ad_swdk Adjoint shortwave (radiation) fractional decay. !
33! !
34! Reference: !
35! !
36! Paulson, C.A., and J.J. Simpson, 1977: Irradiance meassurements !
37! in the upper ocean, J. Phys. Oceanogr., 7, 952-956. !
38! !
39! This routine was adapted from Bill Large 1995 code. !
40! !
41!=======================================================================
42!
43 USE mod_param
44 USE mod_mixing
45 USE mod_scalars
46!
47 implicit none
48!
49! Imported variable declarations.
50!
51 integer, intent(in) :: ng, tile
52 integer, intent(in) :: LBi, UBi, LBj, UBj
53 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
54
55 real(r8), intent(in) :: Zscale
56
57 real(r8), intent(in) :: Z(IminS:ImaxS,JminS:JmaxS)
58
59 real(r8), intent(inout) :: ad_Z(IminS:ImaxS,JminS:JmaxS)
60 real(r8), intent(inout) :: ad_swdk(IminS:ImaxS,JminS:JmaxS)
61!
62! Local variable declarations.
63!
64 integer :: Jindex, i, j
65
66 real(r8) :: cff1, cff2
67 real(r8) :: ad_cff1, ad_cff2
68
69 real(r8), dimension(IminS:ImaxS) :: fac1, fac2, fac3
70
71# include "set_bounds.h"
72!
73!-----------------------------------------------------------------------
74! Initialize adjoint private variables.
75!-----------------------------------------------------------------------
76!
77 ad_cff1=0.0_r8
78 ad_cff2=0.0_r8
79!
80!-----------------------------------------------------------------------
81! Use Paulson and Simpson (1977) two wavelength bands solar
82! absorption model.
83!-----------------------------------------------------------------------
84!
85 DO j=jstr,jend
86 DO i=istr,iend
87 jindex=int(mixing(ng)%Jwtype(i,j))
88 fac1(i)=zscale/lmd_mu1(jindex)
89 fac2(i)=zscale/lmd_mu2(jindex)
90 fac3(i)=lmd_r1(jindex)
91 END DO
92!!DIR$ VECTOR ALWAYS
93 DO i=istr,iend
94 cff1=exp(z(i,j)*fac1(i))
95 cff2=exp(z(i,j)*fac2(i))
96!^ tl_swdk(i,j)=tl_cff1*fac3(i)+ &
97!^ & tl_cff2*(1.0_r8-fac3(i))
98!^
99 ad_cff1=ad_cff1+fac3(i)*ad_swdk(i,j)
100 ad_cff2=ad_cff2+(1.0_r8-fac3(i))*ad_swdk(i,j)
101 ad_swdk(i,j)=0.0_r8
102!^ tl_cff2=fac2(i)*tl_Z(i,j)*cff2
103!^ tl_cff1=fac1(i)*tl_Z(i,j)*cff1
104!^
105 ad_z(i,j)=ad_z(i,j)+ &
106 & fac1(i)*cff1*ad_cff1+ &
107 & fac2(i)*cff2*ad_cff2
108 ad_cff2=0.0_r8
109 ad_cff1=0.0_r8
110 END DO
111 END DO
112 RETURN
type(t_mixing), dimension(:), allocatable mixing
Definition mod_mixing.F:399
real(r8), dimension(9) lmd_mu1
real(r8), dimension(9) lmd_mu2
real(r8), dimension(9) lmd_r1

References mod_scalars::lmd_mu1, mod_scalars::lmd_mu2, mod_scalars::lmd_r1, and mod_mixing::mixing.

Referenced by ad_pre_step3d_mod::ad_pre_step3d_tile().

Here is the caller graph for this function: