Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
M
Méso-NH code
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Package Registry
Model registry
Operate
Environments
Terraform modules
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Terms and privacy
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
RODIER Quentin
Méso-NH code
Commits
aff65eb4
Commit
aff65eb4
authored
5 years ago
by
WAUTELET Philippe
Browse files
Options
Downloads
Patches
Plain Diff
Philippe 17/01/2020: move quicksort to tools + add kpos optional dummy argument
parent
4a84802a
No related branches found
No related tags found
No related merge requests found
Changes
2
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
src/MNH/rain_ice_red.f90
+2
-26
2 additions, 26 deletions
src/MNH/rain_ice_red.f90
src/MNH/tools.f90
+41
-1
41 additions, 1 deletion
src/MNH/tools.f90
with
43 additions
and
27 deletions
src/MNH/rain_ice_red.f90
+
2
−
26
View file @
aff65eb4
!MNH_LIC Copyright 1995-20
19
CNRS, Meteo-France and Universite Paul Sabatier
!MNH_LIC Copyright 1995-20
20
CNRS, Meteo-France and Universite Paul Sabatier
!MNH_LIC This is part of the Meso-NH software governed by the CeCILL-C licence
!MNH_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
!MNH_LIC for details. version 1.
...
...
@@ -242,6 +242,7 @@ END MODULE MODI_RAIN_ICE_RED
! P. Wautelet 10/04/2019: replace ABORT and STOP calls by Print_msg
! P. Wautelet 28/05/2019: move COUNTJV function to tools.f90
! P. Wautelet 29/05/2019: remove PACK/UNPACK intrinsics (to get more performance and better OpenACC support)
! P. Wautelet 17/01/2020: move Quicksort to tools.f90
!
!* 0. DECLARATIONS
! ------------
...
...
@@ -1729,30 +1730,5 @@ CONTAINS
!
!
END
SUBROUTINE
CORRECT_NEGATIVITIES
!
recursive
subroutine
quicksort
(
a
,
first
,
last
)
implicit
none
integer
a
(
*
),
x
,
t
integer
first
,
last
integer
i
,
j
x
=
a
(
(
first
+
last
)
/
2
)
i
=
first
j
=
last
do
do
while
(
a
(
i
)
<
x
)
i
=
i
+1
end
do
do
while
(
x
<
a
(
j
))
j
=
j
-1
end
do
if
(
i
>=
j
)
exit
t
=
a
(
i
);
a
(
i
)
=
a
(
j
);
a
(
j
)
=
t
i
=
i
+1
j
=
j
-1
end
do
if
(
first
<
i
-1
)
call
quicksort
(
a
,
first
,
i
-1
)
if
(
j
+1
<
last
)
call
quicksort
(
a
,
j
+1
,
last
)
end
subroutine
quicksort
END
SUBROUTINE
RAIN_ICE_RED
This diff is collapsed.
Click to expand it.
src/MNH/tools.f90
+
41
−
1
View file @
aff65eb4
!MNH_LIC Copyright 2019-20
19
CNRS, Meteo-France and Universite Paul Sabatier
!MNH_LIC Copyright 2019-20
20
CNRS, Meteo-France and Universite Paul Sabatier
!MNH_LIC This is part of the Meso-NH software governed by the CeCILL-C licence
!MNH_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
!MNH_LIC for details. version 1.
...
...
@@ -19,12 +19,14 @@ module mode_tools
!
! Modifications:
! P. Wautelet 28/05/2019: move COUNTJV function to tools.f90
! P. Wautelet 17/01/2020: move Quicksort to tools.f90
implicit
none
private
public
::
Countjv
public
::
Quicksort
public
::
Upcase
interface
Countjv
...
...
@@ -79,6 +81,44 @@ function Countjv3d(ltab,i1,i2,i3) result(ic)
end
function
Countjv3d
recursive
subroutine
Quicksort
(
ka
,
kbeg
,
kend
,
kpos
)
integer
,
dimension
(:),
intent
(
inout
)
::
ka
integer
,
intent
(
in
)
::
kbeg
,
kend
integer
,
dimension
(:),
optional
,
intent
(
inout
)
::
kpos
integer
::
ji
,
jj
integer
::
itmp
,
itmp2
,
ival
ival
=
ka
(
(
kbeg
+
kend
)
/
2
)
ji
=
kbeg
jj
=
kend
do
do
while
(
ka
(
ji
)
<
ival
)
ji
=
ji
+
1
end
do
do
while
(
ival
<
ka
(
jj
)
)
jj
=
jj
-
1
end
do
if
(
ji
>=
jj
)
exit
itmp
=
ka
(
ji
)
ka
(
ji
)
=
ka
(
jj
)
ka
(
jj
)
=
itmp
if
(
present
(
kpos
)
)
then
itmp2
=
kpos
(
ji
)
kpos
(
ji
)
=
kpos
(
jj
)
kpos
(
jj
)
=
itmp2
end
if
ji
=
ji
+1
jj
=
jj
-1
end
do
if
(
kbeg
<
ji
-
1
)
call
Quicksort
(
ka
,
kbeg
,
ji
-
1
,
kpos
)
if
(
jj
+
1
<
kend
)
call
Quicksort
(
ka
,
jj
+
1
,
kend
,
kpos
)
end
subroutine
Quicksort
function
Upcase
(
hstring
)
character
(
len
=*
),
intent
(
in
)
::
hstring
character
(
len
=
len
(
hstring
))
::
upcase
...
...
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