My Project
fct_s.f90
Go to the documentation of this file.
1 
2 
3 
4 
5 
6 
7 
8 
9 
10 
11 
12 !/===========================================================================/
13 ! Copyright (c) 2007, The University of Massachusetts Dartmouth
14 ! Produced at the School of Marine Science & Technology
15 ! Marine Ecosystem Dynamics Modeling group
16 ! All rights reserved.
17 !
18 ! FVCOM has been developed by the joint UMASSD-WHOI research team. For
19 ! details of authorship and attribution of credit please see the FVCOM
20 ! technical manual or contact the MEDM group.
21 !
22 !
23 ! This file is part of FVCOM. For details, see http://fvcom.smast.umassd.edu
24 ! The full copyright notice is contained in the file COPYRIGHT located in the
25 ! root directory of the FVCOM code. This original header must be maintained
26 ! in all distributed versions.
27 !
28 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
29 ! AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO,
30 ! THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
31 ! PURPOSE ARE DISCLAIMED.
32 !
33 !/---------------------------------------------------------------------------/
34 ! CVS VERSION INFORMATION
35 ! $Id$
36 ! $Name$
37 ! $Revision$
38 !/===========================================================================/
39 
40 !==============================================================================|
41 ! FLUX CONTROL FOR SALINITY |
42 !==============================================================================|
43 
44 SUBROUTINE fct_s
45  !# if defined (1)
46 
47  !==============================================================================|
48  USE all_vars
49  USE mod_utils
50  USE bcs
51  USE mod_obcs
52  IMPLICIT NONE
53  REAL(SP):: SMAX,SMIN
54  INTEGER :: I,J,K,K1
55  !==============================================================================|
56 
57  IF(heating_type == 'body') RETURN
58 
59  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"Start: fct_s"
60 
61  nodes: DO i=1,m
62 
63  ! SKIP OPEN BOUNDARY NODES
64  IF(iobcn > 0)THEN
65  DO j=1,iobcn
66  IF(i == i_obc_n(j)) cycle nodes
67  END DO
68  END IF
69 
70  ! SKIP RIVER INFLOW POINTS
71  IF(numqbc > 0)THEN
72  DO j=1,numqbc
73  IF(river_inflow_location == 'node')THEN
74  IF(i == inodeq(j)) cycle nodes
75  END IF
76  IF(river_inflow_location == 'edge')THEN
77  IF(i == n_icellq(j,1) .OR. i == n_icellq(j,2)) cycle nodes
78  END IF
79  END DO
80  END IF
81 
82  ! SKIP GROUND WATER INFLOW POINTS
83  IF(bfwdis(i) .GT. 0.0_sp .and. groundwater_salt_on) cycle nodes
84 
85  k1 = 1
86  IF(precipitation_on) k1 = 2
87 ! DO K=1,KBM1
88  DO k=k1,kbm1
89  smax = maxval(s1(nbsn(i,1:ntsn(i)),k))
90  smin = minval(s1(nbsn(i,1:ntsn(i)),k))
91 
92  IF(k == 1)THEN
93  smax = max(smax,(s1(i,k)*dz(i,k+1)+s1(i,k+1)*dz(i,k))/ &
94  (dz(i,k)+dz(i,k+1)))
95  smin = min(smin,(s1(i,k)*dz(i,k+1)+s1(i,k+1)*dz(i,k))/ &
96  (dz(i,k)+dz(i,k+1)))
97  ELSE IF(k == kbm1)THEN
98  smax = max(smax,(s1(i,k)*dz(i,k-1)+s1(i,k-1)*dz(i,k))/ &
99  (dz(i,k)+dz(i,k-1)))
100  smin = min(smin,(s1(i,k)*dz(i,k-1)+s1(i,k-1)*dz(i,k))/ &
101  (dz(i,k)+dz(i,k-1)))
102  ELSE
103  smax = max(smax,(s1(i,k)*dz(i,k-1)+s1(i,k-1)*dz(i,k))/ &
104  (dz(i,k)+dz(i,k-1)), &
105  (s1(i,k)*dz(i,k+1)+s1(i,k+1)*dz(i,k))/ &
106  (dz(i,k)+dz(i,k+1)))
107  smin = min(smin,(s1(i,k)*dz(i,k-1)+s1(i,k-1)*dz(i,k))/ &
108  (dz(i,k)+dz(i,k-1)), &
109  (s1(i,k)*dz(i,k+1)+s1(i,k+1)*dz(i,k))/ &
110  (dz(i,k)+dz(i,k+1)))
111  END IF
112 
113  IF(smin-sf1(i,k) > 0.0_sp)sf1(i,k) = smin
114  IF(sf1(i,k)-smax > 0.0_sp)sf1(i,k) = smax
115 
116  END DO
117  END DO nodes
118 
119  WHERE(sf1 < 0.0_sp)sf1=0.0_sp
120 
121  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"End: fct_s"
122  !# endif
123 END SUBROUTINE fct_s
124 !==============================================================================|
125 
126 
integer, dimension(:), allocatable, target ntsn
Definition: mod_main.f90:1023
logical function dbg_set(vrb)
Definition: mod_utils.f90:182
real(sp), dimension(:,:), allocatable, target s1
Definition: mod_main.f90:1308
integer iobcn
Definition: mod_main.f90:1777
integer, dimension(:), allocatable i_obc_n
Definition: mod_main.f90:1779
real(sp), dimension(:,:), allocatable, target sf1
Definition: mod_main.f90:1311
real(sp), dimension(:), allocatable, target bfwdis
Definition: mod_main.f90:1196
integer, dimension(:,:), allocatable, target n_icellq
Definition: mod_main.f90:1216
real(sp), dimension(:,:), allocatable, target dz
Definition: mod_main.f90:1092
subroutine fct_s
Definition: fct_s.f90:45
integer, dimension(:,:), allocatable, target nbsn
Definition: mod_main.f90:1030
integer, parameter dbg_sbr
Definition: mod_utils.f90:69
integer, dimension(:), allocatable, target inodeq
Definition: mod_main.f90:1214