Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
L
loadbalancing
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Package registry
Container registry
Model registry
Operate
Terraform modules
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
GitLab community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
SLMS
loadbalancing
Commits
960707d1
Commit
960707d1
authored
Nov 15, 2019
by
Stephan Schulz
Browse files
Options
Downloads
Patches
Plain Diff
update/new Fortran example
parent
065c1405
Branches
Branches containing commit
Tags
Tags containing commit
1 merge request
!8
Refactor
Pipeline
#25994
failed
Nov 15, 2019
Stage: build
Changes
3
Pipelines
1
Show whitespace changes
Inline
Side-by-side
Showing
3 changed files
example/ALL_test_f_obj.F90
+112
-0
112 additions, 0 deletions
example/ALL_test_f_obj.F90
example/CMakeLists.txt
+21
-0
21 additions, 0 deletions
example/CMakeLists.txt
src/ALL_module_obj.F90
+1
-1
1 addition, 1 deletion
src/ALL_module_obj.F90
with
134 additions
and
1 deletion
example/ALL_test_f_obj.F90
0 → 100644
+
112
−
0
View file @
960707d1
!Copyright 2018 Rene Halver, Forschungszentrum Juelich GmbH, Germany
!Copyright 2018 Godehard Sutmann, Forschungszentrum Juelich GmbH, Germany
!
!Redistribution and use in source and binary forms, with or without modification, are
!permitted provided that the following conditions are met:
!
!1. Redistributions of source code must retain the above copyright notice, this list
! of conditions and the following disclaimer.
!
!2. Redistributions in binary form must reproduce the above copyright notice, this
! list of conditions and the following disclaimer in the documentation and/or
! other materials provided with the distribution.
!
!3. Neither the name of the copyright holder nor the names of its contributors
! may be used to endorse or promote products derived from this software without
! specific prior written permission.
!
!THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY
!EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES
!OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT
!SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
!INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
!TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR
!BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
!CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
!ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF
!SUCH DAMAGE.
! module for ALL access from Fortran
program
ALL_test_f
use
ALL_module_obj
use
iso_c_binding
#ifndef ALL_USE_F08
use
mpi
#else
use
mpi_f08
#endif
implicit
none
#ifndef ALL_USE_F08
integer
::
cart_comm
#else
type
(
MPI_Comm
)
::
cart_comm
#endif
integer
::
dims
(
3
)
integer
::
coords
(
3
)
logical
::
period
(
3
)
real
(
8
)
::
length
(
3
)
real
(
8
)
::
vertices
(
3
,
2
)
integer
::
error
integer
::
rank
,
n_ranks
real
(
8
),
dimension
(:,:),
allocatable
::
new_vertices
integer
::
i
type
(
ALL_t
)
::
lb
call
MPI_INIT
(
error
)
call
MPI_COMM_RANK
(
MPI_COMM_WORLD
,
rank
,
error
)
call
MPI_COMM_SIZE
(
MPI_COMM_WORLD
,
n_ranks
,
error
)
! create cartesian communicator
dims
=
0
period
=
.true.
call
MPI_DIMS_CREATE
(
n_ranks
,
3
,
dims
,
error
)
call
MPI_CART_CREATE
(
MPI_COMM_WORLD
,
3
,
dims
,
period
,
.true.
,
cart_comm
,
error
)
! compute cell length
length
=
1.0d0
/
real
(
dims
,
8
)
! compute vertices of local domain
call
MPI_CART_COORDS
(
cart_comm
,
rank
,
3
,
coords
,
error
)
vertices
(:,
1
)
=
real
(
coords
,
8
)
*
length
vertices
(:,
2
)
=
real
(
coords
+1
,
8
)
*
length
do
i
=
0
,
n_ranks
-1
if
(
i
==
rank
)
then
write
(
*
,
"(a,i7,2(a,3es14.7))"
)
"rank: "
,
rank
,
" old vertices: "
,
vertices
(:,
1
),
", "
,
vertices
(:,
2
)
call
MPI_BARRIER
(
cart_comm
,
error
)
else
call
MPI_BARRIER
(
cart_comm
,
error
)
end
if
end
do
call
lb
%
init
(
ALL_STAGGERED
,
3
,
4.0d0
)
call
lb
%
set_work
(
real
(
product
(
coords
,
1
)
*
64
,
8
))
call
lb
%
set_vertices
(
vertices
)
! if using a non cartesian communicator, this call would be required
! call lb%set_proc_grid_params(coords,dims)
call
lb
%
set_communicator
(
cart_comm
)
call
lb
%
setup
()
call
lb
%
balance
()
! If not allocatable as target, must use
! lb%get_new_number_of_vertices and lb%get_vertices instead!
call
lb
%
get_vertices_alloc
(
new_vertices
)
do
i
=
0
,
n_ranks
-1
if
(
i
==
rank
)
then
write
(
*
,
"(2(a,i7),2(a,3es14.7))"
)
"rank: "
,
rank
,
" "
,
size
(
new_vertices
,
2
),
" new vertices: "
,
new_vertices
(:,
1
),
&
", "
,
new_vertices
(:,
2
)
call
MPI_BARRIER
(
cart_comm
,
error
)
else
call
MPI_BARRIER
(
cart_comm
,
error
)
end
if
end
do
deallocate
(
new_vertices
)
call
MPI_FINALIZE
(
error
);
END
PROGRAM
This diff is collapsed.
Click to expand it.
example/CMakeLists.txt
+
21
−
0
View file @
960707d1
...
@@ -56,4 +56,25 @@ if(CM_ALL_FORTRAN)
...
@@ -56,4 +56,25 @@ if(CM_ALL_FORTRAN)
RUNTIME DESTINATION bin
RUNTIME DESTINATION bin
LIBRARY DESTINATION lib
LIBRARY DESTINATION lib
INCLUDES DESTINATION include
)
INCLUDES DESTINATION include
)
add_executable
(
ALL_test_f_obj ALL_test_f_obj.F90
)
if
(
CM_ALL_VTK_OUTPUT
)
target_include_directories
(
ALL_test_f_obj PUBLIC
${
VTK_INCLUDE_DIRS
}
)
target_link_libraries
(
ALL_test_f_obj PUBLIC
${
VTK_LIBRARY_DIRS
}
)
endif
(
CM_ALL_VTK_OUTPUT
)
#if(CM_ALL_VORONOI)
# target_include_directories(ALL_test_f_obj PUBLIC ${CM_ALL_VORO_INC})
# target_link_libraries(ALL_test_f_obj PUBLIC ${CM_ALL_VORO_LIB})
#endif(CM_ALL_VORONOI)
target_link_libraries
(
ALL_test_f_obj LINK_PUBLIC ALL
)
target_link_libraries
(
ALL_test_f_obj LINK_PUBLIC ALL_fortran
)
target_include_directories
(
ALL_test_f_obj PUBLIC
${
CMAKE_BINARY_DIR
}
/modules
)
list
(
APPEND CMAKE_MODULE_PATH
${
CMAKE_BINARY_DIR
}
/modules
)
set_property
(
TARGET ALL_test_f_obj PROPERTY LINKER_LANGUAGE Fortran
)
install
(
TARGETS ALL_test_f_obj
RUNTIME DESTINATION bin
LIBRARY DESTINATION lib
INCLUDES DESTINATION include
)
endif
(
CM_ALL_FORTRAN
)
endif
(
CM_ALL_FORTRAN
)
This diff is collapsed.
Click to expand it.
src/ALL_module_obj.F90
+
1
−
1
View file @
960707d1
! DCMAKE_Fortran_FLAGS='-std=f2003 -Wall -Wextra -fbacktrace' -DCM_ALL_FORTRAN=ON -DCMAKE_BUILD_TYPE=Debug -DCM
AKE
_ALL_DEBUG=ON
..
!
-
DCMAKE_Fortran_FLAGS='-std=f2003 -Wall -Wextra -fbacktrace' -DCM_ALL_FORTRAN=ON -DCMAKE_BUILD_TYPE=Debug -DCM_ALL_DEBUG=ON
module
ALL_module_obj
module
ALL_module_obj
use
iso_c_binding
use
iso_c_binding
use
iso_fortran_env
use
iso_fortran_env
...
...
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
sign in
to comment