Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
P
PlaSim
Manage
Activity
Members
Code
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Locked files
Deploy
Releases
Container registry
Model registry
Analyze
Contributor analytics
Repository analytics
Insights
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
GitLab community forum
Contribute to GitLab
Provide feedback
Terms and privacy
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
Lunkeit, Frank
PlaSim
Commits
c88fafb8
Commit
c88fafb8
authored
Nov 8, 2016
by
HartmutBorth
Browse files
Options
Downloads
Patches
Plain Diff
indroduced new simulations
parent
ce32d6cb
No related branches found
No related tags found
No related merge requests found
Changes
3
Show whitespace changes
Inline
Side-by-side
Showing
3 changed files
cat/src/cat.f90
+23
-12
23 additions, 12 deletions
cat/src/cat.f90
cat/src/cat_namelist
+4
-4
4 additions, 4 deletions
cat/src/cat_namelist
cat/src/simmod.f90
+44
-32
44 additions, 32 deletions
cat/src/simmod.f90
with
71 additions
and
48 deletions
cat/src/cat.f90
+
23
−
12
View file @
c88fafb8
...
...
@@ -263,10 +263,12 @@ complex(8) :: phi
integer
::
nforc
=
0
! forcing switch
! 0 = no forcing
! 1 = constant forcing with
! spectral forcing ring
! 2 = markov chain
! 3 = forcing with constant wave numbers
!
! 1 = constant frocing defined by external field
! 2 = spectral forcing with constant amplitudes
! forcing ring and random phases
! 3 = markov chain
! 4 = forcing with constant wave numbers
! but random uncorrelated phases (white noise)
integer
::
kfmin
=
0
! min radius = sqrt(kfmin) of spectr. forcing
...
...
@@ -387,8 +389,17 @@ integer :: ndatim(6) = 1 ! date/time display
real
(
4
)
::
parc
(
5
)
=
0.0
! timeseries display
!--- predefined simulations (simmod)
integer
::
nsim
=
0
! 1/0 predefined simulations on/off.
!--- code of predefined simulations (simmod)
integer
::
nsim
=
0
! 0 no predefined simulation is specified
! nsim > 0 predefined simulation is run
!
! available predefined simulations (experiments):
! -----------------------------------------------
! code name
!
! 1 decaying_jet01
! -----------------------------------------------
!
! Predefined simulations are specified in
! <sim_namelist> of simmod.
...
...
@@ -1614,6 +1625,9 @@ enf=0.d0
select
case
(
nforc
)
case
(
1
)
cq
=
cq
+
cqfrc
case
(
2
)
do
ifk
=
1
,
nk
i
=
in
(
ifk
,
1
)
j
=
in
(
ifk
,
2
)
...
...
@@ -1625,7 +1639,7 @@ select case (nforc)
enf
=
enf
+
(
cq
(
i
,
j
)
*
cq
(
i
,
j
))/
k2
enddo
case
(
2
)
case
(
3
)
ic
=
tstep
+1
if
(
mod
(
ic
,
itau
)
.eq.
0.
)
then
call
random_number
(
ran4
)
...
...
@@ -1642,7 +1656,7 @@ select case (nforc)
enf
=
enf
+
(
cq
(
i
,
j
)
*
cq
(
i
,
j
))/
k2
enddo
case
(
3
)
case
(
4
)
do
ifk
=
1
,
nk
i
=
in
(
ifk
,
1
)
j
=
in
(
ifk
,
2
)
...
...
@@ -1654,9 +1668,6 @@ select case (nforc)
enf
=
enf
+
(
cq
(
i
,
j
)
*
cq
(
i
,
j
))/
k2
enddo
case
(
4
)
cq
=
cq
+
cqfrc
end
select
return
...
...
This diff is collapsed.
Click to expand it.
cat/src/cat_namelist
+
4
−
4
View file @
c88fafb8
...
...
@@ -6,13 +6,13 @@
dt = 1.0d-3
myseed = 1,2,3,4,5,6,7,8
outsp = 138
nsteps = 200000
nsteps
= 200000
0
ngp = 250
nstdout = 10000
ncfl =
5
00
npert =
2
ncfl =
100
00
npert =
3
apert = 0.00001
nforc =
0
nforc =
4
kfmin = 5
kfmax = 7
aforc = 0.0001
...
...
This diff is collapsed.
Click to expand it.
cat/src/simmod.f90
+
44
−
32
View file @
c88fafb8
...
...
@@ -16,7 +16,7 @@ integer, parameter :: nusimsp = 130 ! spectral fields
character
(
256
)
::
sim_namelist
=
"sim_namelist"
!--- predefined experimets
character
(
256
)
::
sim
=
"
iv00"
! type of predefined simulation
character
(
256
)
::
y
sim
=
"
djet01"
! type of predefined simulation
! options are:
!----------------------------------!
...
...
@@ -24,27 +24,34 @@ character (256) :: sim = "iv00" ! type of predefined simulation
! flows) !
!----------------------------------!
! "
iv00"
top hat jet
! "
djet01" decaying
top hat jet
!-----------------------!
! Forced decaying Flows !
!-----------------------!
! "f
d00" top hat wind forcing
! "f
jet01" forced top hat jet
!--- parameters of djet01 (initial top hat jet)
integer
::
w1
=
8
! half width of jet center (in grid points)
integer
::
w2
=
4
! width of vortex sheet (in grid points)
integer
::
scl
=
1
! horizontal scale of jet
real
(
8
)
::
qmax
=
1.0
! amplitude of vortex sheets
!--- parameters of iv01 (initial top hat jet)
integer
::
iv00w1
=
8
! half width of jet center (in grid points)
integer
::
iv00w2
=
4
! width of vortex sheet (in grid points)
integer
::
iv00scl
=
1
! horizontal scale of jet
real
(
8
)
::
iv00qmax
=
1.0
! amplitude of vortex sheets
!--- parameters of fd01 (forced decaying top hat jet)
integer
::
fd00w1
=
8
! half width of jet center (in grid points)
integer
::
fd00w2
=
4
! width of vortex sheet (in grid points)
integer
::
fd00scl
=
1
! horizontal scale of jet
real
(
8
)
::
fd00qmax
=
1.0d-4
! amplitude of vortex sheets
!--- sim_namelist parameters to overwrite cat_namelist parameters
!---------------------------------------------------------------!
! A sim_namelist parameter sXXXXX corresponds to a cat_namelist !
! parameter XXXXX. The sim_namelist parameters sXXXXX can only !
! take the same values as the corresponding cat_namelist !
! parameter. !
!---------------------------------------------------------------!
integer
::
snforc
=
-1
! type of forcing
integer
::
snpert
=
-1
! type of perturbation
integer
::
snpost
=
-1
! type of post processing
end
module
simmod
...
...
@@ -71,9 +78,9 @@ use simmod
implicit
none
!--- define sim_namelist
namelist
/
sim_nl
/
sim
,
&
iv00
qmax
,
iv00
w1
,
iv00
w2
,
iv00scl
,
&
fd00qmax
,
fd00w1
,
fd00w2
,
fd00scl
namelist
/
sim_nl
/
y
sim
,
&
qmax
,
w1
,
w2
,
scl
,
&
snpert
,
snforc
,
snpost
!--- check if sim_namelist is present
inquire
(
file
=
sim_namelist
,
exist
=
lsimnl
)
...
...
@@ -86,6 +93,10 @@ else
return
endif
!--- overwrite cat_namelist parameters
if
(
snpert
.ge.
0
)
npert
=
snpert
if
(
snforc
.ge.
0
)
nforc
=
snforc
if
(
snpost
.ge.
0
)
npost
=
snpost
return
end
subroutine
sim_readnl
...
...
@@ -103,40 +114,41 @@ integer :: jy
real
(
8
)
::
gpvar
(
1
:
ngx
,
1
:
ngy
)
select
case
(
sim
)
select
case
(
y
sim
)
!------------------------!
! initial value problems !
!------------------------!
!--- top-hat jet
case
(
"
iv00
"
)
case
(
"
djet01
"
)
gpvar
(:,:)
=
0.0
do
jy
=
1
,
ngy
if
(
jy
.ge.
ngy
/
2+1
-
iv00
scl
*
(
iv00w1
+
iv00
w2
)
.and.
&
jy
.le.
ngy
/
2
-
iv00scl
*
iv00
w1
)
then
gpvar
(:,
jy
)
=
-
iv00
qmax
if
(
jy
.ge.
ngy
/
2+1
-
scl
*
(
w1
+
w2
)
.and.
&
jy
.le.
ngy
/
2
-
scl
*
w1
)
then
gpvar
(:,
jy
)
=
-
qmax
endif
if
(
jy
.ge.
ngy
/
2+1
+
iv00scl
*
iv00
w1
.and.
&
jy
.le.
ngy
/
2
+
iv00
scl
*
(
iv00w1
+
iv00
w2
)
)
then
gpvar
(:,
jy
)
=
iv00
qmax
if
(
jy
.ge.
ngy
/
2+1
+
scl
*
w1
.and.
&
jy
.le.
ngy
/
2
+
scl
*
(
w1
+
w2
)
)
then
gpvar
(:,
jy
)
=
qmax
endif
enddo
call
sim_wrtgp
(
gpvar
,
qcde
,
1
)
!-------------------!
! forced turbulence !
!-------------------!
!--- top-hat jet
case
(
"f
d00
"
)
case
(
"f
jet01
"
)
gpvar
(:,:)
=
0.0
do
jy
=
1
,
ngy
if
(
jy
.ge.
ngy
/
2+1
-
fd00
scl
*
(
fd00w1
+
fd00
w2
)
.and.
&
jy
.le.
ngy
/
2
-
fd00scl
*
fd00
w1
)
then
gpvar
(:,
jy
)
=
-
fd00
qmax
if
(
jy
.ge.
ngy
/
2+1
-
scl
*
(
w1
+
w2
)
.and.
&
jy
.le.
ngy
/
2
-
scl
*
w1
)
then
gpvar
(:,
jy
)
=
-
qmax
endif
if
(
jy
.ge.
ngy
/
2+1
+
fd00scl
*
fd00
w1
.and.
&
jy
.le.
ngy
/
2
+
fd00
scl
*
(
fd00w1
+
fd00
w2
)
)
then
gpvar
(:,
jy
)
=
fd00
qmax
if
(
jy
.ge.
ngy
/
2+1
+
scl
*
w1
.and.
&
jy
.le.
ngy
/
2
+
scl
*
(
w1
+
w2
)
)
then
gpvar
(:,
jy
)
=
qmax
endif
enddo
call
sim_wrtgp
(
gpvar
,
qfrccde
,
1
)
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment