Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
! |
! ||===\\ |
! || \\ |
! || || //==\\ || || //==|| ||/==\\ |
! || || || || || || || || || || |
! || // || || || || || || || |
! ||===// \\==// \\==\\ \\==\\ || |
! |
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
! |
! DEFINE_BC_SEGMENTED_S_LINE_PARABOLA August 2015 |
! |
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
subroutine define_bc_segmented_s_line_parabola(params,osolve,vo,bcdef,nest)
!------------------------------------------------------------------------------|
!(((((((((((((((( Purpose of the routine )))))))))))))))))))))))))))))))))))))
!------------------------------------------------------------------------------|
! This routine assigns the velocity boundary conditions for the segmented s-line
! geometry
!------------------------------------------------------------------------------|
!(((((((((((((((( declaration of the subroutine arguments ))))))))))))))))))))
!------------------------------------------------------------------------------|
use definitions
!use mpi
implicit none
include 'mpif.h'
type (parameters) params
type (octreesolve) osolve
type (void) vo
type (bc_definition) bcdef
type (nest_info) nest
!------------------------------------------------------------------------------|
!(((((((((((((((( declaration of the subroutine internal variables )))))))))))))
!------------------------------------------------------------------------------|
integer i,iproc,nproc,ierr
double precision :: eps,lsf0,pi,lorig,h,x1,x2,phi,yend,cper,cscl,xstart,ystart
double precision :: theta,l,vin,vzfluxscl,cntvel,dxy,xend,xsym,ymax,xwidth
double precision :: ywidth,xdisp,ydisp,nb
double precision,dimension(:),allocatable :: x0,ldisp
integer ie,ij,j,jp,nelemx,nelemz
double precision :: base,startp,endp,alpha,kink1,kink2,flytt,a,d
double precision :: wmax,uend,vmag,dipangle
double precision e(1),f(1),y(1)
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
call mpi_comm_size (mpi_comm_world,nproc,ierr)
call mpi_comm_rank (mpi_comm_world,iproc,ierr)
eps=1.d-10
osolve%kfix=0
osolve%kfixt=0
pi=atan(1.d0)*4.d0
l=bcdef%bc_parameters(1)
base=bcdef%bc_parameters(2)
startp=bcdef%bc_parameters(3)
endp=bcdef%bc_parameters(4)
vin=bcdef%bc_parameters(5)
alpha=bcdef%bc_parameters(6)
kink1=bcdef%bc_parameters(7)
kink2=bcdef%bc_parameters(8)
nelemx=idint(bcdef%bc_parameters(9))
nelemz=idint(bcdef%bc_parameters(10))
nb=2**params%levelmax_oct
dxy=1.d0/2**(params%levelmax_oct+1.d0)
alpha=alpha*pi/180
flytt=(kink2-kink1) * tan(alpha)
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
a = (base-l) /((endp-startp)*(endp-startp))
d = startp
wmax=(2*vin*l)/(endp-startp)
uend=wmax/(2*a*(endp-startp))
vmag=sqrt(uend**2.d0+wmax**2.d0)
do i=1,osolve%nnode
if (osolve%x(i).lt.eps) then
osolve%kfix((i-1)*3+1)=1 ; osolve%u(i)=1.d0
osolve%kfix((i-1)*3+2)=1 ; osolve%v(i)=0.d0
osolve%kfix((i-1)*3+3)=1 ; osolve%w(i)=0.d0
endif
if (osolve%x(i).gt.1.d0-eps) then
osolve%kfix((i-1)*3+1)=1 ; osolve%u(i)=0.d0
endif
if (osolve%y(i).lt.eps) then
osolve%kfix((i-1)*3+2)=1 ; osolve%v(i)=0.d0
endif
if (osolve%y(i).gt.1.d0-eps) then
osolve%kfix((i-1)*3+2)=1 ; osolve%v(i)=0.d0
endif
if (osolve%z(i).lt.eps) then
osolve%kfix((i-1)*3+1)=1 ; osolve%u(i)=0.d0
osolve%kfix((i-1)*3+2)=1 ; osolve%v(i)=0.d0
osolve%kfix((i-1)*3+3)=1 ; osolve%w(i)=0.d0
osolve%kfixt(i)=1 ; osolve%temp(i)=1.d0
! Velocity in first margin-normal convergence segment
if (osolve%y(i) .lt. kink1) then
! Slab dip angle
dipangle = atan(2*a*(osolve%x(i)-startp))
! Velocity before transition to subduction
if (osolve%x(i) .le. startp-real(nelemz)*dxy) then
osolve%u(i) = vin
! Velocity in transition to subduction
elseif (osolve%x(i) .le. startp+real(nelemz)*dxy) then
! Velocity contribution from incoming horizontal velocity
! (diminishes to zero across transition)
osolve%u(i) = vin * (1.d0 - (osolve%x(i)-(startp-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy))
! Velocity contribution from subduction region
! (increases from zero to full magnitude across transition)
osolve%u(i) = osolve%u(i) + vmag*cos(dipangle) * (osolve%x(i)-(startp-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy)
osolve%w(i) = vmag*sin(dipangle) * (osolve%x(i)-(startp-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy)
! Velocity in subduction region
elseif (osolve%x(i) .le. endp-real(nelemz)*dxy) then
osolve%u(i) = vmag*cos(dipangle)
osolve%w(i) = vmag*sin(dipangle)
! Velocity in transition from subduction region
elseif (osolve%x(i) .le. endp+real(nelemz)*dxy) then
! Velocity contribution from subduction region
! (decreases to zero across transition)
osolve%u(i) = vmag*cos(dipangle) * (1.d0 - (osolve%x(i)-(endp-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy))
osolve%w(i) = vmag*sin(dipangle) * (1.d0 - (osolve%x(i)-(endp-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy))
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
! Velocity in oblique convergence segment
elseif (osolve%y(i) .lt. kink2) then
! Slab dip angle
dipangle = atan(2*a*(osolve%x(i)-((osolve%y(i)-kink1)*tan(alpha)+startp)))
! Velocity before transition to subduction
if (osolve%x(i) .le. (osolve%y(i)-kink1)*tan(alpha)+startp-real(nelemz)*dxy) then
osolve%u(i) = vin
! Velocity in transition to subduction
elseif (osolve%x(i) .le. (osolve%y(i)-kink1)*tan(alpha)+startp+real(nelemz)*dxy) then
! Velocity contribution from incoming horizontal velocity
! (diminishes to zero across transition)
osolve%u(i)= vin * (1.d0 - (osolve%x(i)-((osolve%y(i)-kink1)*tan(alpha)+startp-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy))
! Velocity contribution from subduction region
! (increases from zero to full magnitude across transition)
osolve%u(i)= osolve%u(i) + vmag*cos(dipangle) * (osolve%x(i)-((osolve%y(i)-kink1)*tan(alpha)+startp-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy)
osolve%w(i) = vmag*sin(dipangle) * (osolve%x(i)-((osolve%y(i)-kink1)*tan(alpha)+startp-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy)
! Velocity in subduction region
elseif (osolve%x(i) .le. (osolve%y(i)-kink1)*tan(alpha)+endp-real(nelemz)*dxy) then
osolve%u(i) = vmag*cos(dipangle)
osolve%w(i) = vmag*sin(dipangle)
! Velocity in transition from subduction region
elseif (osolve%x(i) .le. (osolve%y(i)-kink1)*tan(alpha)+endp+real(nelemz)*dxy) then
! Velocity contribution from subduction region
! (decreases to zero across transition)
osolve%u(i) = vmag*cos(dipangle) * (1.d0 - (osolve%x(i)-((osolve%y(i)-kink1)*tan(alpha)+endp-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy))
osolve%w(i) = vmag*sin(dipangle) * (1.d0 - (osolve%x(i)-((osolve%y(i)-kink1)*tan(alpha)+endp-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy))
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
! Velocity in second margin-normal convergence segment
else
! Slab dip angle
dipangle = atan(2*a*(osolve%x(i)-(flytt+startp)))
! Velocity before transition to subduction
if (osolve%x(i) .le. flytt+startp-real(nelemz)*dxy) then
osolve%u(i) = vin
! Velocity in transition to subduction
elseif (osolve%x(i) .le. flytt+startp+real(nelemz)*dxy) then
! Velocity contribution from incoming horizontal velocity
! (diminishes to zero across transition)
osolve%u(i) = vin * (1.d0 - (osolve%x(i)-((flytt+startp)-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy))
! Velocity contribution from subduction region
! (increases from zero to full magnitude across transition)
osolve%u(i) = osolve%u(i) + vmag * cos(dipangle) * (osolve%x(i)-((flytt+startp)-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy)
osolve%w(i) = vmag * sin(dipangle) * (osolve%x(i)-((flytt+startp)-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy)
! Velocity in subduction region
elseif (osolve%x(i) .le. flytt+endp-real(nelemz)*dxy) then
osolve%u(i) = vmag * cos(dipangle)
osolve%w(i) = vmag * sin(dipangle)
! Velocity in transition from subduction region
elseif (osolve%x(i) .le. flytt+endp+real(nelemz)*dxy) then
! Velocity contribution from subduction region
! (decreases to zero across transition)
osolve%u(i) = vmag * cos(dipangle) * (1.d0 - (osolve%x(i)-((flytt+endp)-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy))
osolve%w(i) = vmag * sin(dipangle) * (1.d0 - (osolve%x(i)-((flytt+endp)-real(nelemz)*dxy))/(2.0*real(nelemz)*dxy))
endif
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
endif
endif
if (osolve%z(i).gt.1.d0-eps) then
osolve%kfix((i-1)*3+1)=1 ; osolve%u(i)=0.d0
osolve%kfix((i-1)*3+2)=1 ; osolve%v(i)=0.d0
osolve%kfix((i-1)*3+3)=1 ; osolve%w(i)=0.d0
osolve%kfixt(i)=1 ; osolve%temp(i)=0.d0
endif
if (.not.vo%influid(i)) then
osolve%kfixt(i)=1
osolve%temp(i)=0.d0
endif
if (abs(bcdef%utrans).gt.eps) then
if (osolve%kfix((i-1)*3+1)==1) osolve%u(i)=osolve%u(i)+bcdef%utrans
endif
if (abs(bcdef%vtrans).gt.eps) then
if (osolve%kfix((i-1)*3+2)==1) osolve%v(i)=osolve%v(i)+bcdef%vtrans
endif
if (osolve%kfix((i-1)*3+1)==1) osolve%u(i)=osolve%u(i)*vin
if (osolve%kfix((i-1)*3+2)==1) osolve%v(i)=osolve%v(i)*vin
if (osolve%kfix((i-1)*3+3)==1) osolve%w(i)=osolve%w(i)*vin
enddo
if (params%isobc) then
call define_isostasy_bc(params,osolve,bcdef)
endif
end
!------------------------------------------------------------------------------|