-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathtemplate.py
236 lines (204 loc) · 8.65 KB
/
template.py
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
86
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
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
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
174
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
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
233
234
235
236
"""
Wave_runup
"""
from __future__ import division
from past.utils import old_div
import numpy as np
from proteus import (Domain, Context)
from proteus.Profiling import logEvent
from proteus.mprans.SpatialTools import Tank2D
from proteus.mprans import SpatialTools as st
import proteus.TwoPhaseFlow.TwoPhaseFlowProblem as TpFlow
from proteus.Gauges import PointGauges, LineIntegralGauges, LineGauges
from proteus import WaveTools as wt
from proteus.ctransportCoefficients import smoothedHeaviside
import math
# *************************** #
# ***** GENERAL OPTIONS ***** #
# *************************** #
opts= Context.Options([
("final_time",10.0,"Final time for simulation"),
("dt_output",0.01,"Time interval to output solution"),
("cfl",0.25,"Desired CFL restriction"),
("he",0.1,"he relative to Length of domain in x"),
("refinement",3,"level of refinement"),
# ("toe",35.0,"where toe begins from origin (m from x)"),
("slope",(1/19.85),"Beta, slope of incline (y/x)"),
("slope_length",50.0,"right extent of domain x(m)"),
("wl",1.0,"water level"),
("waves",True, "wave on/off"),
("wave_height", 0.28, "Height of the waves in s"),
("wave_dir", (1.,0.,0.), "Direction of the waves (from left boundary)"),
("wave_type", 'solitaryWave', "type of wave"),
("fast", False, "switch for fast cosh calculations in WaveTools"),
("g", [0, -9.81, 0], "Gravity vector in m/s^2"),
])
toe=(-opts.wl/opts.slope)+35#shifted
structured=False
slope=opts.slope
slope_x=opts.slope_length-toe
slope_y=slope*slope_x
top=(slope_y)*1.35
#toe=opts.toe
he=opts.he
nnx=nny=nnz=None
wl=opts.wl
a=opts.wave_height
k=((3*a)/(4*(opts.wl**3)))**0.5
L=(2/k)*np.arccosh(0.05**(-0.5))
x0=(-opts.wl/opts.slope)-(L/2)+35 #shifted
if opts.waves is True:
height = opts.wave_height
mwl = depth = opts.wl
direction = opts.wave_dir
wave = wt.SolitaryWave( waveHeight = height,
mwl = wl,
depth = depth,
g = np.array(opts.g),
waveDir = direction,
trans = np.array([x0, 0., 0.]),
fast = opts.fast
)
# ****************** #
# ***** GAUGES ***** #
# ****************** #
# *************************** #
# ***** DOMAIN AND MESH ***** #
# ****************** #******* #
domain = Domain.PlanarStraightLineGraphDomain()
nLevels = 1
#parallelPartitioningType = proteus.MeshTools.MeshParallelPartitioningTypes.node
nLayersOfOverlapForParallel = 0
boundaries=['left','right','bottom','slope','top']
boundaryOrientations = {'bottom': np.array([0., -1.,0.]),
'right': np.array([+1., 0.,0.]),
'top': np.array([0., +1.,0.]),
'left': np.array([-1., 0.,0.]),
'slope': np.array([-1., 0.,0.]),}
boundaryTags=dict([(key,i+1) for (i,key) in enumerate(boundaries)])
vertices=[[0.0,0.0],#0
[toe,0.0],#1
[opts.slope_length,slope_y],#2
[opts.slope_length,top],#3
# [toe+slope_x,top],#3
[0.0,top]]#4
vertexFlags=[boundaryTags['left'],
boundaryTags['bottom'],
boundaryTags['slope'],
boundaryTags['right'],
boundaryTags['top']]
segments=[[0,1],
[1,2],
[2,3],
[3,4],
[4,0]]
segmentFlags=[boundaryTags['bottom'],
boundaryTags['slope'],
boundaryTags['right'],
boundaryTags['top'],
boundaryTags['left']]
regions=[[0.1, 0.01]]
regionFlags=[1]
tank = st.CustomShape(domain, vertices=vertices, vertexFlags=vertexFlags,
segments=segments, segmentFlags=segmentFlags,
regions=regions, regionFlags=regionFlags,
boundaryTags=boundaryTags, boundaryOrientations=boundaryOrientations)
##############################
#domain
#############################
# ****************************** #
# ***** INITIAL CONDITIONS ***** #
# ****************************** #
class zero(object):
def uOfXT(self,x,t):
return 0.0
class clsvof_init_cond(object):
def uOfXT(self,x,t):
# return x[1]-wl
return x[1] - (wave.eta(x,0) + opts.wl)
# if x[0] < waterLine_x and x[1] < waterLine_y:
# return -1.0
# elif x[0] > waterLine_x or x[1] > waterLine_y:
# return 1.0
# else:
# return 0.0
epsFact_consrv_heaviside=3.0
wavec = np.sqrt(9.81 * (depth+opts.wave_height))
def weight(x,t):
return 1.0-smoothedHeaviside(epsFact_consrv_heaviside*opts.he,
#-ct.epsFact_consrv_heaviside*ct.opts.he+
(x[1] - (max(wave.eta(x, t%(toe/wavec)),
wave.eta(x+toe, t%(toe/wavec)))
+opts.wl)))
class vel_u(object):
def uOfXT(self, x, t):
if x[1] <= wave.eta(x,t) + opts.wl:
return weight(x,t)*wave.u(x,t)[0]
else:
return 0.0
class vel_v(object):
def uOfXT(self, x, t):
if x[1] <= wave.eta(x,t) + opts.wl:
return weight(x,t)*wave.u(x,t)[1]
else:
return 0.0
# ****************************** #
# ***** Boundary CONDITIONS***** #
# ****************************** #
tank.BC['top'].setAtmosphere()
tank.BC['bottom'].setFreeSlip()
tank.BC['left'].setFreeSlip()
tank.BC['right'].setFreeSlip()
tank.BC['slope'].setFreeSlip()
he = opts.he
domain.MeshOptions.he = he
st.assembleDomain(domain)
domain.MeshOptions.triangleOptions = "VApq30Dena%8.8f" % (old_div((he ** 2), 2.0),)
############################################
# ***** Create myTwoPhaseFlowProblem ***** #
############################################
outputStepping = TpFlow.OutputStepping(opts.final_time,dt_output=opts.dt_output)
initialConditions = {'pressure': zero(),
'pressure_increment': zero(),
'vel_u': vel_u(),
'vel_v': vel_v(),
'clsvof': clsvof_init_cond()}
boundaryConditions = {
# DIRICHLET BCs #
'pressure_DBC': lambda x, flag: domain.bc[flag].p_dirichlet.init_cython(),
'pressure_increment_DBC': lambda x, flag: domain.bc[flag].pInc_dirichlet.init_cython(),
'vel_u_DBC': lambda x, flag: domain.bc[flag].u_dirichlet.init_cython(),
'vel_v_DBC': lambda x, flag: domain.bc[flag].v_dirichlet.init_cython(),
'vel_w_DBC': lambda x, flag: domain.bc[flag].w_dirichlet.init_cython(),
'clsvof_DBC': lambda x, flag: domain.bc[flag].vof_dirichlet.init_cython(),
# ADVECTIVE FLUX BCs #
'pressure_AFBC': lambda x, flag: domain.bc[flag].p_advective.init_cython(),
'pressure_increment_AFBC': lambda x, flag: domain.bc[flag].pInc_advective.init_cython(),
'vel_u_AFBC': lambda x, flag: domain.bc[flag].u_advective.init_cython(),
'vel_v_AFBC': lambda x, flag: domain.bc[flag].v_advective.init_cython(),
'vel_w_AFBC': lambda x, flag: domain.bc[flag].w_advective.init_cython(),
'clsvof_AFBC': lambda x, flag: domain.bc[flag].vof_advective.init_cython(),
# DIFFUSIVE FLUX BCs #
'pressure_increment_DFBC': lambda x, flag: domain.bc[flag].pInc_diffusive.init_cython(),
'vel_u_DFBC': lambda x, flag: domain.bc[flag].u_diffusive.init_cython(),
'vel_v_DFBC': lambda x, flag: domain.bc[flag].v_diffusive.init_cython(),
'vel_w_DFBC': lambda x, flag: domain.bc[flag].w_diffusive.init_cython(),
'clsvof_DFBC': lambda x, flag: None}
# auxVariables={'clsvof': [height_gauges],
# 'pressure': [pressure_gauges]}
myTpFlowProblem = TpFlow.TwoPhaseFlowProblem(ns_model=1,
nd=2,
cfl=opts.cfl,
outputStepping=outputStepping,
structured=structured,
he=he,
nnx=nnx,
nny=nny,
nnz=None,
domain=domain,
initialConditions=initialConditions,
boundaryConditions=boundaryConditions,
# auxVariables=auxVariables,
useSuperlu=False)
physical_parameters = myTpFlowProblem.physical_parameters
myTpFlowProblem.clsvof_parameters['disc_ICs']=False