Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

updates tests; fixes issue w/ SAVE_SEISMOGRAMS_IN_ADJOINT_RUN #830

Merged
merged 8 commits into from
Dec 20, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension


Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
63 changes: 51 additions & 12 deletions .github/workflows/CI.yml
Original file line number Diff line number Diff line change
Expand Up @@ -130,7 +130,7 @@ jobs:
run: make tests

linuxCheck:
name: Test on ${{ matrix.os }}
name: Test on ${{ matrix.os }}
runs-on: ${{ matrix.os }}
needs: changesCheck

Expand All @@ -147,7 +147,7 @@ jobs:

- name: configure
run: ./configure

- name: make
run: make -j2 all

Expand Down Expand Up @@ -176,31 +176,67 @@ jobs:
sudo echo "deb https://apt.repos.intel.com/oneapi all main" | sudo tee /etc/apt/sources.list.d/oneAPI.list
sudo apt-get update
echo ""
echo "packages intel oneapi:"
sudo apt-get install -y intel-oneapi-compiler-fortran intel-oneapi-compiler-dpcpp-cpp-and-cpp-classic intel-oneapi-mpi intel-oneapi-mpi-devel
# info
#sudo -E apt-cache pkgnames intel | grep intel-oneapi
#echo ""
echo "installing packages intel oneapi:"
sudo apt-get install -y intel-oneapi-compiler-fortran-2023.2.2 intel-oneapi-compiler-dpcpp-cpp-and-cpp-classic-2023.2.2 intel-oneapi-mpi intel-oneapi-mpi-devel
echo ""

- name: compiler infos
run: |
echo ""
source /opt/intel/oneapi/setvars.sh
echo ""
echo "compiler versions:"
echo "icx --version"
which icx
icx --version
echo ""
echo "icc --version"
which icc
icc --version
echo ""
echo "ifx --version"
which ifx
ifx --version
echo ""
echo "ifort --version"
which ifort
ifort --version
echo ""
echo "mpiifort --version"
which mpiifort
mpiifort --version
echo ""
echo "mpif90 --version"
which mpif90
mpif90 --version
echo ""
# infos
which ifort
which icc
which mpiifort
echo "mpirun:"
which mpirun
echo ""
# intel setup for running tests
echo ""
echo "replacing mpif90 with mpiifort link:"
sudo ln -sf $(which mpiifort) $(which mpif90)
mpif90 --version
echo ""
# debug
#export I_MPI_DEBUG=5,pid,host
#export I_MPI_LIBRARY_KIND=debug
# remove -ftrapuv which leads to issues for running tests
sed -i "s/-ftrapuv//g" flags.guess
# environment setting
export TERM=xterm
# export info
echo "exports:"
export
echo ""
which ifort
which icc
which mpiifort
echo ""
printenv >> $GITHUB_ENV
echo "CXX=icpc" >> $GITHUB_ENV
Expand All @@ -209,7 +245,8 @@ jobs:
echo ""

- name: configure debug
run: ./configure --enable-debug FC=ifort CC=icc MPIFC=mpiifort
run: |
./configure --enable-debug FC=ifort CC=icc MPIFC=mpiifort

- name: make debug
run: |
Expand All @@ -218,16 +255,18 @@ jobs:
make clean

- name: configure
run: ./configure FC=ifort CC=icc MPIFC=mpiifort
run: |
./configure FC=ifort CC=icc MPIFC=mpiifort

- name: make
run: |
make -j2 all
make clean

# fails due to MPI issue on virtual nodes
#- name: make tests
# run: make tests
# note: fails with -ftrapuv flag due to MPI issue on virtual nodes
- name: make tests
run: |
make tests

linuxTest_0:
name: Test run example 0 - make tests
Expand Down
19 changes: 14 additions & 5 deletions .travis.yml
Original file line number Diff line number Diff line change
Expand Up @@ -182,11 +182,20 @@ before_install:
#- test -n $CC && unset CC

# updates repository
- |
travis_retry sudo apt-key adv --keyserver keyserver.ubuntu.com --recv-keys 6B05F25D762E3157
travis_retry sudo apt-key adv --keyserver keyserver.ubuntu.com --recv-keys 78BD65473CB3BD13
travis_retry sudo apt-key adv --keyserver keyserver.ubuntu.com --recv-keys 762E3157
travis_retry sudo apt-get update
# (fails currently...)
#- |
# travis_retry sudo apt-key adv --keyserver keyserver.ubuntu.com --recv-keys 6B05F25D762E3157
# travis_retry sudo apt-key adv --keyserver keyserver.ubuntu.com --recv-keys 78BD65473CB3BD13
# travis_retry sudo apt-key adv --keyserver keyserver.ubuntu.com --recv-keys 762E3157
# apt explicit update
#- |
# echo "update apt-get"
# travis_retry sudo apt-get update
# echo
addons:
# apt update
apt:
update: true


install:
Expand Down
3 changes: 2 additions & 1 deletion flags.guess
Original file line number Diff line number Diff line change
Expand Up @@ -99,8 +99,9 @@ case $my_FC in
# I/O throughput lingers at 2.5 MB/s, with it it can increase to ~44 MB/s
# However it does not make much of a difference on NFS mounted volumes or with SFS 3.1.1 / Lustre 1.6.7.1
#
# warnings about external function calls can be suppressed by "-warn all,noexternal" for version > 2018
# optimization report: "-vec-report0" is old and will be replaced by "-qopt-report0 -qopt-report-phase=vec" for v >=15.0
DEF_FFLAGS="-xHost -fpe0 -ftz -assume buffered_io -assume byterecl -align sequence -std08 -diag-disable 6477 -implicitnone -gen-interfaces -warn all" # -mcmodel=medium -shared-intel
DEF_FFLAGS="-xHost -fpe0 -ftz -assume buffered_io -assume byterecl -align sequence -std08 -diag-disable 6477 -implicitnone -gen-interfaces -warn all,noexternal" # -mcmodel=medium -shared-intel
OPT_FFLAGS="-O3 -check nobounds"
DEBUG_FFLAGS="-check all -debug -g -O0 -fp-stack-check -traceback -ftrapuv"
# option "-openmp" is soon deprecated and replaced by "-qopenmp" for versions > 17.x
Expand Down
42 changes: 25 additions & 17 deletions src/auxiliaries/write_profile.f90
Original file line number Diff line number Diff line change
Expand Up @@ -90,8 +90,8 @@
! USER parameters

! initial position
double precision,parameter :: COLAT_0 = 1.0
double precision,parameter :: LON_0 = 1.0
double precision,parameter :: COLAT_0 = 1.d0
double precision,parameter :: LON_0 = 1.d0

! colatitude loop range (in degrees)
integer,parameter :: COLAT_istart = 0 ! 0
Expand Down Expand Up @@ -149,6 +149,9 @@
! depth increment
double precision :: delta

! tolerance in checks for equal (==) float values
double precision, parameter :: TOL_ZERO = 1.d-12

character(len=MAX_STRING_LEN) :: outfile
character(len=7) :: str_info

Expand Down Expand Up @@ -286,10 +289,10 @@
phi_degrees = initial_lon + j*delta_lon ! longitude [0,360]

! checks limits
if (theta_degrees < 0.0) stop 'Error invalid colatitude < 0'
if (theta_degrees > 180.0) stop 'Error invalid colatitude > 180'
if (phi_degrees < 0.0) phi_degrees = phi_degrees + 360.d0
if (phi_degrees > 360.0) phi_degrees = phi_degrees - 360.d0
if (theta_degrees < 0.d0) stop 'Error invalid colatitude < 0'
if (theta_degrees > 180.d0) stop 'Error invalid colatitude > 180'
if (phi_degrees < 0.d0) phi_degrees = phi_degrees + 360.d0
if (phi_degrees > 360.d0) phi_degrees = phi_degrees - 360.d0

Check warning on line 295 in src/auxiliaries/write_profile.f90

View check run for this annotation

Codecov / codecov/patch

src/auxiliaries/write_profile.f90#L292-L295

Added lines #L292 - L295 were not covered by tests

! loads corresponding GLL mesh
if (MODEL_GLL) call load_GLL_mesh(theta_degrees,phi_degrees)
Expand Down Expand Up @@ -410,12 +413,14 @@

! make sure that the Moho discontinuity is at the real moho
if (CRUSTAL) then
if (rmin == RMOHO_FICTITIOUS_IN_MESHER/R_PLANET) rmin = 1.0d0 - moho
if (rmax == RMOHO_FICTITIOUS_IN_MESHER/R_PLANET) rmax = 1.0d0 - moho
! checks rmin == RMOHO_FICTITIOUS_IN_MESHER/R_PLANET
if (abs(rmin - RMOHO_FICTITIOUS_IN_MESHER/R_PLANET) < TOL_ZERO) rmin = 1.0d0 - moho

Check warning on line 417 in src/auxiliaries/write_profile.f90

View check run for this annotation

Codecov / codecov/patch

src/auxiliaries/write_profile.f90#L417

Added line #L417 was not covered by tests
! checks rmax == RMOHO_FICTITIOUS_IN_MESHER/R_PLANET
if (abs(rmax - RMOHO_FICTITIOUS_IN_MESHER/R_PLANET) < TOL_ZERO) rmax = 1.0d0 - moho

Check warning on line 419 in src/auxiliaries/write_profile.f90

View check run for this annotation

Codecov / codecov/patch

src/auxiliaries/write_profile.f90#L419

Added line #L419 was not covered by tests
!print *,'rmin == moho at line ',iline
endif

if (abs(rmin - rmax_last) < 1.d-9) then !!!! rmin == rmax_last: this means that we have just jumped between layers
if (abs(rmin - rmax_last) < TOL_ZERO) then !!!! rmin == rmax_last: this means that we have just jumped between layers

Check warning on line 423 in src/auxiliaries/write_profile.f90

View check run for this annotation

Codecov / codecov/patch

src/auxiliaries/write_profile.f90#L423

Added line #L423 was not covered by tests
! depth increment
! write values every 10 km in the deep earth and every 1 km in the shallow earth
if (rmin > ((R_PLANET/1000.d0)-DELTA_HIRES_DEPTH)/(R_PLANET/1000.d0)) then
Expand All @@ -430,7 +435,7 @@
! sets maximum radius without ocean for 1D models
if (((.not. CRUSTAL) .and. (ROCEAN < R_PLANET)) .and. (.not. TOPOGRAPHY)) then
! stops at ocean depth and adds last ocean layers explicitly
if (rmax == 1.0d0) rmax = ROCEAN/R_PLANET
if (abs(rmax - 1.0d0) < TOL_ZERO) rmax = ROCEAN/R_PLANET ! rmax == 1.d0

Check warning on line 438 in src/auxiliaries/write_profile.f90

View check run for this annotation

Codecov / codecov/patch

src/auxiliaries/write_profile.f90#L438

Added line #L438 was not covered by tests
endif

! backup to detect jump between layers
Expand All @@ -439,26 +444,26 @@
! number of iterations in increments of delta between rmin and rmax
! note: instead of (rmax - rmin), we add a factor (rmax * 0.999999 - rmin) to avoid getting an extra step
! in case the difference is an exact delta match, since we add +1 to nit to reach rmax
nit = floor((rmax*0.9999999 - rmin)/delta) + 1
nit = floor((rmax*0.9999999d0 - rmin)/delta) + 1

Check warning on line 447 in src/auxiliaries/write_profile.f90

View check run for this annotation

Codecov / codecov/patch

src/auxiliaries/write_profile.f90#L447

Added line #L447 was not covered by tests

! debug
!print *,'debug: write profile ilayer/iregion ',ilayer,iregion_code,'rmin/rmax',rmin,rmax,'delta',delta,'nit',nit

do idep = 1,nit+1
! line counters
! inner core boundary
if (rmin == RICB/R_PLANET .and. idep == 1) iline_icb = iline
if (abs(rmin - RICB/R_PLANET) < TOL_ZERO .and. idep == 1) iline_icb = iline ! rmin == RICB/R_PLANET

Check warning on line 455 in src/auxiliaries/write_profile.f90

View check run for this annotation

Codecov / codecov/patch

src/auxiliaries/write_profile.f90#L455

Added line #L455 was not covered by tests
! core mantle boundary
if (rmin == RCMB/R_PLANET .and. idep == 1) iline_cmb = iline
if (abs(rmin - RCMB/R_PLANET) < TOL_ZERO .and. idep == 1) iline_cmb = iline ! rmin == RCMB/R_PLANET

Check warning on line 457 in src/auxiliaries/write_profile.f90

View check run for this annotation

Codecov / codecov/patch

src/auxiliaries/write_profile.f90#L457

Added line #L457 was not covered by tests
! moho
if (CRUSTAL) then
! uses 3D crustal model (e.g. Crust2.0)
if (rmin == (1.0d0 - moho) .and. idep == 1) then
if (abs(rmin - (1.0d0 - moho)) < TOL_ZERO .and. idep == 1) then ! rmin == (1.0d0 - moho)

Check warning on line 461 in src/auxiliaries/write_profile.f90

View check run for this annotation

Codecov / codecov/patch

src/auxiliaries/write_profile.f90#L461

Added line #L461 was not covered by tests
iline_moho = iline
endif
else
! 1D crust from reference model
if (rmin == RMOHO/R_PLANET .and. idep == 1) iline_moho = iline
if (abs(rmin - RMOHO/R_PLANET) < TOL_ZERO .and. idep == 1) iline_moho = iline ! rmin == RMOHO/R_PLANET

Check warning on line 466 in src/auxiliaries/write_profile.f90

View check run for this annotation

Codecov / codecov/patch

src/auxiliaries/write_profile.f90#L466

Added line #L466 was not covered by tests
endif

! radius
Expand All @@ -470,8 +475,8 @@
! make sure we are within the right shell in PREM to honor discontinuities
! use small geometrical tolerance
r_prem = r
if (r <= rmin*1.000001d0) r_prem = rmin*1.000001d0
if (r >= rmax*0.999999d0) r_prem = rmax*0.999999d0
if (r < rmin*1.000001d0) r_prem = rmin*1.000001d0
if (r > rmax*0.999999d0) r_prem = rmax*0.999999d0

Check warning on line 479 in src/auxiliaries/write_profile.f90

View check run for this annotation

Codecov / codecov/patch

src/auxiliaries/write_profile.f90#L478-L479

Added lines #L478 - L479 were not covered by tests

! gets model properties (similar to get_model() routine)
call write_profile_model_values(r,r_prem,theta,phi,iregion_code,idoubling,rmin,rmax, &
Expand Down Expand Up @@ -784,6 +789,9 @@
! local parameters
double precision :: lat,lon

! initializes
elevation = 0.d0

Check warning on line 793 in src/auxiliaries/write_profile.f90

View check run for this annotation

Codecov / codecov/patch

src/auxiliaries/write_profile.f90#L793

Added line #L793 was not covered by tests

! topography elevation
if (TOPOGRAPHY .or. OCEANS) then
if (TOPOGRAPHY) then
Expand Down
22 changes: 11 additions & 11 deletions src/shared/parallel.f90
Original file line number Diff line number Diff line change
Expand Up @@ -88,10 +88,10 @@

! initialize the MPI communicator and start the NPROCTOT MPI processes.
call MPI_INIT(ier)
if (ier /= 0 ) stop 'Error initializing MPI'
if (ier /= 0) stop 'Error initializing MPI'

call MPI_COMM_RANK(MPI_COMM_WORLD,myrank,ier)
if (ier /= 0 ) stop 'Error getting MPI rank'
if (ier /= 0) stop 'Error getting MPI rank'

! we need to make sure that NUMBER_OF_SIMULTANEOUS_RUNS and BROADCAST_SAME_MESH_AND_MODEL are read before calling world_split()
! thus read the parameter file
Expand Down Expand Up @@ -210,7 +210,7 @@

! synchronizes MPI processes
call MPI_BARRIER(my_local_mpi_comm_world, ier)
if (ier /= 0 ) stop 'Error synchronize MPI processes'
if (ier /= 0) stop 'Error synchronize MPI processes'

end subroutine synchronize_all

Expand All @@ -231,7 +231,7 @@

! synchronizes MPI processes
call MPI_BARRIER(comm,ier)
if (ier /= 0 ) stop 'Error synchronize MPI processes for specified communicator'
if (ier /= 0) stop 'Error synchronize MPI processes for specified communicator'

Check warning on line 234 in src/shared/parallel.f90

View check run for this annotation

Codecov / codecov/patch

src/shared/parallel.f90#L234

Added line #L234 was not covered by tests

end subroutine synchronize_all_comm

Expand Down Expand Up @@ -787,7 +787,7 @@
send(:) = buffer(:)

call MPI_ALLREDUCE(send, buffer, countval, MPI_INTEGER, MPI_MAX, my_local_mpi_comm_world, ier)
if (ier /= 0 ) stop 'Allreduce to get max values failed.'
if (ier /= 0) stop 'Allreduce to get max values failed.'

Check warning on line 790 in src/shared/parallel.f90

View check run for this annotation

Codecov / codecov/patch

src/shared/parallel.f90#L790

Added line #L790 was not covered by tests

end subroutine max_allreduce_i

Expand Down Expand Up @@ -1771,7 +1771,7 @@
integer :: ier

call MPI_COMM_SIZE(my_local_mpi_comm_world,sizeval,ier)
if (ier /= 0 ) stop 'Error getting MPI world size'
if (ier /= 0) stop 'Error getting MPI world size'

end subroutine world_size

Expand All @@ -1792,7 +1792,7 @@
integer :: ier

call MPI_COMM_SIZE(comm,sizeval,ier)
if (ier /= 0 ) stop 'Error getting MPI world size'
if (ier /= 0) stop 'Error getting MPI world size'

Check warning on line 1795 in src/shared/parallel.f90

View check run for this annotation

Codecov / codecov/patch

src/shared/parallel.f90#L1795

Added line #L1795 was not covered by tests

end subroutine world_size_comm

Expand All @@ -1812,7 +1812,7 @@
integer :: ier

call MPI_COMM_RANK(my_local_mpi_comm_world,rank,ier)
if (ier /= 0 ) stop 'Error getting MPI rank'
if (ier /= 0) stop 'Error getting MPI rank'

end subroutine world_rank

Expand All @@ -1833,7 +1833,7 @@
integer :: ier

call MPI_COMM_RANK(comm,rank,ier)
if (ier /= 0 ) stop 'Error getting MPI rank'
if (ier /= 0) stop 'Error getting MPI rank'

Check warning on line 1836 in src/shared/parallel.f90

View check run for this annotation

Codecov / codecov/patch

src/shared/parallel.f90#L1836

Added line #L1836 was not covered by tests

end subroutine world_rank_comm

Expand Down Expand Up @@ -1861,7 +1861,7 @@
! instead, a duplicate of a user-specified communicator should always be used."

call MPI_COMM_DUP(my_local_mpi_comm_world,comm,ier)
if (ier /= 0 ) stop 'Error duplicating my_local_mpi_comm_world communicator'
if (ier /= 0) stop 'Error duplicating my_local_mpi_comm_world communicator'

Check warning on line 1864 in src/shared/parallel.f90

View check run for this annotation

Codecov / codecov/patch

src/shared/parallel.f90#L1864

Added line #L1864 was not covered by tests

end subroutine world_duplicate

Expand Down Expand Up @@ -1930,7 +1930,7 @@
integer :: ier

call MPI_Comm_free(comm,ier)
if (ier /= 0 ) stop 'Error freeing MPI communicator'
if (ier /= 0) stop 'Error freeing MPI communicator'

Check warning on line 1933 in src/shared/parallel.f90

View check run for this annotation

Codecov / codecov/patch

src/shared/parallel.f90#L1933

Added line #L1933 was not covered by tests

end subroutine world_comm_free

Expand Down
3 changes: 0 additions & 3 deletions src/specfem3D/setup_sources_receivers.f90
Original file line number Diff line number Diff line change
Expand Up @@ -1247,9 +1247,6 @@ subroutine setup_receivers()
do_save_seismograms = .true.
endif

! sets local receivers to zero if no seismogram needs to be saved
if (.not. do_save_seismograms) nrec_local = 0

! seismogram array length (to write out time portions of the full seismograms)
nlength_seismogram = NTSTEP_BETWEEN_OUTPUT_SEISMOS / NTSTEP_BETWEEN_OUTPUT_SAMPLE

Expand Down
4 changes: 2 additions & 2 deletions src/specfem3D/write_seismograms.f90
Original file line number Diff line number Diff line change
Expand Up @@ -40,7 +40,7 @@ subroutine write_seismograms()
NTSTEP_BETWEEN_OUTPUT_SEISMOS,NTSTEP_BETWEEN_OUTPUT_SAMPLE, &
do_save_seismograms, &
WRITE_SEISMOGRAMS_BY_MAIN,OUTPUT_SEISMOS_ASDF, &
SAVE_SEISMOGRAMS_IN_ADJOINT_RUN,SAVE_SEISMOGRAMS_STRAIN, &
SAVE_SEISMOGRAMS_STRAIN, &
moment_der,sloc_der,shdur_der,stshift_der, &
scale_displ

Expand Down Expand Up @@ -173,7 +173,7 @@ subroutine write_seismograms()
select case (SIMULATION_TYPE)
case (1,3)
! forward/reconstructed wavefields
if (.not. ( SIMULATION_TYPE == 3 .and. (.not. SAVE_SEISMOGRAMS_IN_ADJOINT_RUN) ) ) &
if (do_save_seismograms) &
call write_seismograms_to_file()
if (SAVE_SEISMOGRAMS_STRAIN) &
call write_seismograms_strain()
Expand Down
Loading