SICOPOLIS V3.2
 All Classes Files Functions Variables Macros
borehole.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : b o r e h o l e
4 !
5 !> @file
6 !!
7 !! Computation of an arbitrary field quantity for a given borehole
8 !! position x_pos, y_pos by weighed averaging of the corresponding
9 !! gridded 2-d field.
10 !!
11 !! @section Copyright
12 !!
13 !! Copyright 2009-2016 Ralf Greve
14 !!
15 !! @section License
16 !!
17 !! This file is part of SICOPOLIS.
18 !!
19 !! SICOPOLIS is free software: you can redistribute it and/or modify
20 !! it under the terms of the GNU General Public License as published by
21 !! the Free Software Foundation, either version 3 of the License, or
22 !! (at your option) any later version.
23 !!
24 !! SICOPOLIS is distributed in the hope that it will be useful,
25 !! but WITHOUT ANY WARRANTY; without even the implied warranty of
26 !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
27 !! GNU General Public License for more details.
28 !!
29 !! You should have received a copy of the GNU General Public License
30 !! along with SICOPOLIS. If not, see <http://www.gnu.org/licenses/>.
31 !<
32 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
33 
34 !-------------------------------------------------------------------------------
35 !> Computation of an arbitrary field quantity for a given borehole
36 !! position x_pos, y_pos by weighed averaging of the corresponding
37 !! gridded 2-d field.
38 !<------------------------------------------------------------------------------
39 subroutine borehole(field, x_pos, y_pos, dxi, deta, ch_grid, field_val)
40 
41 use sico_types
43 use sico_vars
44 
45 implicit none
46 
47 real(dp), dimension(0:JMAX,0:IMAX), intent(in) :: field
48 real(dp), intent(in) :: x_pos, y_pos, dxi, deta
49 character (len=*), intent(in) :: ch_grid
50 
51 real(dp), intent(out) :: field_val
52 
53 integer(i4b) :: i_kl, i_gr, j_kl, j_gr
54 real(dp) :: real_i, real_j, dist1, dist2, dist3, dist4, &
55  weight1, weight2, weight3, weight4
56 
57 !-------- Neighbour points, weighing factors --------
58 
59 real_i = (x_pos-xi(0)) /dxi
60 real_j = (y_pos-eta(0))/deta
61 
62 i_kl = floor(real_i)
63 j_kl = floor(real_j)
64 i_gr = ceiling(real_i)
65 j_gr = ceiling(real_j)
66 
67 if (i_kl <= 1) stop ' borehole: i_kl <= 1 not allowed!'
68 if (j_kl <= 1) stop ' borehole: j_kl <= 1 not allowed!'
69 if (i_gr >= imax-1) stop ' borehole: i_gr >= IMAX-1 not allowed!'
70 if (j_gr >= jmax-1) stop ' borehole: j_gr >= JMAX-1 not allowed!'
71 
72 dist1 = sqrt( (x_pos-xi(i_kl))**2 + (y_pos-eta(j_kl))**2 )
73 dist2 = sqrt( (x_pos-xi(i_kl))**2 + (y_pos-eta(j_gr))**2 )
74 dist3 = sqrt( (x_pos-xi(i_gr))**2 + (y_pos-eta(j_kl))**2 )
75 dist4 = sqrt( (x_pos-xi(i_gr))**2 + (y_pos-eta(j_gr))**2 )
76 
77 weight1 = 1.0_dp/(dist1+eps)
78 weight2 = 1.0_dp/(dist2+eps)
79 weight3 = 1.0_dp/(dist3+eps)
80 weight4 = 1.0_dp/(dist4+eps)
81 
82 !-------- Inverse-distance weighing of the four adjacent grid values --------
83 
84 if (ch_grid=='grid') then ! field(j,i) defined on grid points
85 
86  field_val = ( weight1*field(j_kl,i_kl) &
87  +weight2*field(j_gr,i_kl) &
88  +weight3*field(j_kl,i_gr) &
89  +weight4*field(j_gr,i_gr) ) &
90  /(weight1 + weight2 + weight3 + weight4)
91 
92 else if (ch_grid=='sg_x') then
93  ! field(j,i) defined on staggered grid in x direction
94 
95  field_val = ( weight1*0.5_dp*(field(j_kl,i_kl)+field(j_kl,i_kl-1)) &
96  +weight2*0.5_dp*(field(j_gr,i_kl)+field(j_gr,i_kl-1)) &
97  +weight3*0.5_dp*(field(j_kl,i_gr)+field(j_kl,i_gr-1)) &
98  +weight4*0.5_dp*(field(j_gr,i_gr)+field(j_gr,i_gr-1)) ) &
99  /(weight1 + weight2 + weight3 + weight4)
100 
101 else if (ch_grid=='sg_y') then
102  ! field(j,i) defined on staggered grid in y direction
103 
104  field_val = ( weight1*0.5_dp*(field(j_kl,i_kl)+field(j_kl-1,i_kl)) &
105  +weight2*0.5_dp*(field(j_gr,i_kl)+field(j_gr-1,i_kl)) &
106  +weight3*0.5_dp*(field(j_kl,i_gr)+field(j_kl-1,i_gr)) &
107  +weight4*0.5_dp*(field(j_gr,i_gr)+field(j_gr-1,i_gr)) ) &
108  /(weight1 + weight2 + weight3 + weight4)
109 
110 else
111 
112  stop ' borehole: Parameter ch_grid has undefined value!'
113 
114 end if
115 
116 end subroutine borehole
117 !
subroutine borehole(field, x_pos, y_pos, dxi, deta, ch_grid, field_val)
Computation of an arbitrary field quantity for a given borehole position x_pos, y_pos by weighed aver...
Definition: borehole.F90:39
Declarations of kind types for SICOPOLIS.
Definition: sico_types.F90:35
Declarations of global variables for SICOPOLIS (for the ANT domain).
Definition: sico_vars.F90:35
Declarations of global variables for SICOPOLIS.