SPECFEM / specfem3d_globe

SPECFEM3D_GLOBE simulates global and regional (continental-scale) seismic wave propagation.
GNU General Public License v3.0
90 stars 96 forks source link

meshfem3D crashing on Power9 machine #690

Closed lsawade closed 4 years ago

lsawade commented 4 years ago

Hi!

I'm having trouble meshing the most recent version (master/devel) on a Power9 machine at Princeton. Details below. I'm still trying to backtrace the exact issue. The issue cannot be replicated on Summit.


Hardware and last good commit (as far as I know)

Over the past half year I have been using the Traverse cluster at Princeton (Traverse). I wasn't up to date with the repo and the last commit I was working with was

* b087774 - Fri, 13 Mar 2020 02:56:36 +0100 (3 months ago) (origin/devel)
|           code cleanup - Daniel Peter

This commit compiles without issue using ADIOS 1.13.1.

Error

Using the most recent master or devel, the error is apparently thrown in the specfem3d_globe/src/meshfem3D/meshfem3D_models.F90 routine in Line 352 which calls . I'm pasting both output_mesher.txt as well as the stdout below.

Note that the issue is seemingly not related to the dimensions of the velocity model (1D_transverse_isotropic_prem here) since I tried compiling s362ani as well, without luck. It might have to do with background values(?).

output_mesher.txt ``` $ specfem3d_globe/OUTPUT_FILES/output_mesher.txt **************************** *** Specfem3D MPI Mesher *** **************************** Version: v7.0.2-237-gd491b15 Planet: Earth There are 6 MPI processes Processes are numbered from 0 to 5 There are 128 elements along xi in each chunk There are 128 elements along eta in each chunk There are 1 slices along xi in each chunk There are 1 slices along eta in each chunk There is a total of 1 slices in each chunk There are 6 chunks in the global mesh There is a total of 6 slices in the global mesh NGLLX = 5 NGLLY = 5 NGLLZ = 5 Shape functions defined by NGNOD = 27 control nodes Surface shape functions defined by NGNOD2D = 9 control nodes model: 1D_transversely_isotropic_prem incorporating the oceans using equivalent load incorporating ellipticity incorporating surface topography incorporating self-gravitation (Cowling approximation) incorporating rotation incorporating attenuation using 3 standard linear solids no 3-D lateral variations in the mantle no heterogeneities in the mantle no crustal variations using unmodified 1D crustal model with two layers incorporating transverse isotropy no inner-core anisotropy no general mantle anisotropy Reference radius of the globe used is 6371.0000000000000 km Central cube is at a radius of 950.00000000000000 km creating global slice addressing Spatial distribution of the slices 0 2 1 4 5 3 incorporating topography topo file : DATA/topo_bathy/topo_bathy_etopo4_smoothed_window_7.bin resolution in minutes: 4.00000000 topography/bathymetry: min/max = -7747 5507 Elapsed time for reading in seconds = 0.13036037500000003 additional mesh optimizations moho: no element stretching for 3-D moho surface internal topography 410/660: no element stretching for 3-D internal surfaces Radial Meshing parameters: NCHUNKS = 6 CENTER LAT/LON: 0.00000000 / 0.00000000 GAMMA_ROTATION_AZIMUTH: 0.00000000 CHUNK WIDTH XI/ETA: 90.0000000 / 90.0000000 NEX XI/ETA: 128 / 128 NER_CRUST: 2 NER_80_MOHO: 1 NER_220_80: 2 NER_400_220: 2 NER_600_400: 2 NER_670_600: 1 NER_771_670: 1 NER_TOPDDOUBLEPRIME_771: 15 NER_CMB_TOPDDOUBLEPRIME: 1 NER_OUTER_CORE: 16 NER_TOP_CENTRAL_CUBE_ICB: 2 SUPPRESS_CRUSTAL_MESH: F R_CENTRAL_CUBE = 950.000000 km Mesh resolution: DT = 0.19000000000000000 Minimum period = 34.7484131 (s) MIN_ATTENUATION_PERIOD = 34.7484131 MAX_ATTENUATION_PERIOD = 1954.04700 ******************************************* creating mesh in region 1 this region is the crust and mantle ******************************************* first pass ...allocating arrays ...setting up layers ...creating mesh elements creating layer 1 out of 11 setting tiso flags in mantle model 9.1% current clock (NOT elapsed) time is: 10h 57min 54sec creating layer 2 out of 11 18.2% current clock (NOT elapsed) time is: 10h 57min 54sec creating layer 3 out of 11 27.3% current clock (NOT elapsed) time is: 10h 57min 54sec creating layer 4 out of 11 36.4% current clock (NOT elapsed) time is: 10h 57min 55sec creating layer 5 out of 11 45.5% current clock (NOT elapsed) time is: 10h 57min 55sec creating layer 6 out of 11 54.5% current clock (NOT elapsed) time is: 10h 57min 56sec creating layer 7 out of 11 63.6% current clock (NOT elapsed) time is: 10h 57min 56sec creating layer 8 out of 11 72.7% current clock (NOT elapsed) time is: 10h 57min 56sec creating layer 9 out of 11 81.8% current clock (NOT elapsed) time is: 10h 57min 56sec creating layer 10 out of 11 90.9% current clock (NOT elapsed) time is: 10h 57min 57sec creating layer 11 out of 11 100.0% current clock (NOT elapsed) time is: 10h 57min 58sec ...creating global addressing ...creating MPI buffers second pass ...allocating arrays ...setting up layers ...creating mesh elements creating layer 1 out of 11 9.1% current clock (NOT elapsed) time is: 10h 58min 58sec creating layer 2 out of 11 18.2% current clock (NOT elapsed) time is: 10h 59min 01sec creating layer 3 out of 11 27.3% current clock (NOT elapsed) time is: 10h 59min 01sec creating layer 4 out of 11 36.4% current clock (NOT elapsed) time is: 10h 59min 06sec creating layer 5 out of 11 45.5% current clock (NOT elapsed) time is: 10h 59min 15sec creating layer 6 out of 11 54.5% current clock (NOT elapsed) time is: 10h 59min 16sec creating layer 7 out of 11 63.6% current clock (NOT elapsed) time is: 10h 59min 17sec creating layer 8 out of 11 72.7% current clock (NOT elapsed) time is: 10h 59min 20sec creating layer 9 out of 11 81.8% current clock (NOT elapsed) time is: 10h 59min 23sec creating layer 10 out of 11 90.9% current clock (NOT elapsed) time is: 10h 59min 28sec creating layer 11 out of 11 100.0% current clock (NOT elapsed) time is: 10h 59min 34sec ...fills global mesh points ...checking mesh resolution and time step ---------------------------------- Verification of mesh parameters: ---------------------------------- Region is crust/mantle Min Vs = 3.19999981 (km/s) Max Vp = 13.7166214 (km/s) Max element edge size = 223.495636 (km) Min element edge size = 9.07211971 (km) Max/min ratio = 24.6354370 Max Jacobian eigenvalue ratio = 0.998868346 Min Jacobian eigenvalue ratio = 6.10016361E-02 Minimum period resolved = 26.3479710 (s) Minimum period resolved (empirical) = 34.7484131 (s) Maximum suggested time step = 0.125000000 (s) for DT : 0.189999998 (s) Max stability for wave velocities = 0.824762464 ---------------------------------- ...precomputing Jacobian ...creating chunk buffers ----- creating chunk buffers ----- There are 1 slices along xi in each chunk There are 1 slices along eta in each chunk There is a total of 1 slices in each chunk There are 6 chunks There is a total of 6 slices in all the chunks There is a total of 12 messages to assemble faces between chunks Generating message 1 for faces out of 12 Generating message 2 for faces out of 12 Generating message 3 for faces out of 12 Generating message 4 for faces out of 12 Generating message 5 for faces out of 12 Generating message 6 for faces out of 12 Generating message 7 for faces out of 12 Generating message 8 for faces out of 12 Generating message 9 for faces out of 12 Generating message 10 for faces out of 12 Generating message 11 for faces out of 12 Generating message 12 for faces out of 12 all the messages for chunk faces have the right size Generating message 1 for corners out of 8 Generating message 2 for corners out of 8 Generating message 3 for corners out of 8 Generating message 4 for corners out of 8 Generating message 5 for corners out of 8 Generating message 6 for corners out of 8 Generating message 7 for corners out of 8 Generating message 8 for corners out of 8 ...preparing MPI interfaces crust/mantle region: #max of points in MPI buffers along xi npoin2D_xi = 28589 #max of array elements transferred npoin2D_xi*NDIM = 85767 #max of points in MPI buffers along eta npoin2D_eta = 28589 #max of array elements transferred npoin2D_eta*NDIM = 85767 crust mantle MPI: maximum interfaces: 4 MPI addressing maximum interfaces: 4 MPI addressing : all interfaces okay total MPI interface points : 686136 unique MPI interface points: 683520 maximum valence : 2 total unique MPI interface points: 683520 ...element inner/outer separation for overlapping of communications with calculations: percentage of edge elements in crust/mantle 5.41955996 % percentage of volume elements in crust/mantle 94.5804367 % ...element mesh coloring mesh coloring: F ...creating mass matrix updates mass matrix with ocean load ...saving binary files solver in ADIOS 1 file format saving arrays in ADIOS file: ./DATABASES_MPI/solver_data.bp saving arrays in ADIOS file: ./DATABASES_MPI/boundary.bp saving arrays in ADIOS file: ./DATABASES_MPI/attenuation.bp MPI in ADIOS 1 file format saving arrays in ADIOS file: ./DATABASES_MPI/solver_data_mpi.bp calculated region volume: 3.50144863 top area: 12.5574989 bottom area: 3.74933386 ******************************************* creating mesh in region 2 this region is the outer core ******************************************* first pass ...allocating arrays ...setting up layers ...creating mesh elements creating layer 1 out of 2 50.0% current clock (NOT elapsed) time is: 11h 00min 12sec creating layer 2 out of 2 100.0% current clock (NOT elapsed) time is: 11h 00min 13sec ...creating global addressing ...creating MPI buffers second pass ...allocating arrays ...setting up layers ...creating mesh elements creating layer 1 out of 2 ```
Error stdout/stderr ```bash [traverse-k02g2:76669:0:76669] Caught signal 8 (Floating point exception: floating-point invalid operation) [traverse-k02g2:76670:0:76670] Caught signal 8 (Floating point exception: floating-point invalid operation) [traverse-k02g2:76671:0:76671] Caught signal 8 (Floating point exception: floating-point invalid operation) [traverse-k02g2:76672:0:76672] Caught signal 8 (Floating point exception: floating-point invalid operation) [traverse-k02g2:76673:0:76673] Caught signal 8 (Floating point exception: floating-point invalid operation) [traverse-k02g2:76674:0:76674] Caught signal 8 (Floating point exception: floating-point invalid operation) ==== backtrace (tid: 76670) ==== 0 0x00000000100c3c2c meshfem3d_models_get1d_val_() /scratch/gpfs/lsawade/MagicScripts/specfem3d_globe/src/meshfem3D/meshfem3D_models.F90:352 1 0x000000001009def0 get_model_() /scratch/gpfs/lsawade/MagicScripts/specfem3d_globe/src/meshfem3D/get_model.F90:165 2 0x000000001006c940 compute_element_properties_() /scratch/gpfs/lsawade/MagicScripts/specfem3d_globe/src/meshfem3D/compute_element_properties.f90:149 3 0x0000000010095468 create_regular_elements_() /scratch/gpfs/lsawade/MagicScripts/specfem3d_globe/src/meshfem3D/create_regular_elements.f90:277 4 0x00000000100944ec create_regions_elements_() /scratch/gpfs/lsawade/MagicScripts/specfem3d_globe/src/meshfem3D/create_regions_elements.f90:141 5 0x0000000010090cb0 create_regions_mesh_() /scratch/gpfs/lsawade/MagicScripts/specfem3d_globe/src/meshfem3D/create_regions_mesh.F90:176 6 0x00000000100819d8 create_meshes_() /scratch/gpfs/lsawade/MagicScripts/specfem3d_globe/src/meshfem3D/create_meshes.f90:159 7 0x0000000010005658 xmeshfem3d() /scratch/gpfs/lsawade/MagicScripts/specfem3d_globe/src/meshfem3D/meshfem3D.f90:377 8 0x0000000010005658 main() /scratch/gpfs/lsawade/MagicScripts/specfem3d_globe/src/meshfem3D/meshfem3D.f90:385 9 0x0000000000025200 generic_start_main.isra.0() libc-start.c:0 10 0x00000000000253f4 __libc_start_main() ???:0 . . . [same for the other workers] . . . ================================= srun: error: traverse-k02g2: tasks 0-5: Floating point exception srun: Terminating job step 78840.0 ```

Configuration & Compilation

Below, I'll past all the modules used etc. and how I flag the compilation

Par_file ``` SIMULATION_TYPE = 1 # set to 1 for forward simulations, 2 for adjoint simulations for sources, and 3 for kernel simulations NOISE_TOMOGRAPHY = 0 # flag of noise tomography, three steps (1,2,3). If earthquake simulation, set it to 0. SAVE_FORWARD = .false. # save last frame of forward simulation or not NCHUNKS = 6 ANGULAR_WIDTH_XI_IN_DEGREES = 20.d0 # angular size of a chunk ANGULAR_WIDTH_ETA_IN_DEGREES = 20.d0 CENTER_LATITUDE_IN_DEGREES = 40.d0 CENTER_LONGITUDE_IN_DEGREES = 25.d0 GAMMA_ROTATION_AZIMUTH = 0.d0 NEX_XI = 128 NEX_ETA = 128 NPROC_XI = 1 NPROC_ETA = 1 # transversely_isotropic_prem_plus_3D_crust_2.0, 3D_anisotropic, 3D_attenuation, # s20rts, s40rts, s362ani, s362iso, s362wmani, s362ani_prem, s362ani_3DQ, s362iso_3DQ, # s29ea, sea99_jp3d1994, sea99, jp3d1994, heterogen, full_sh, sgloberani_aniso, sgloberani_iso MODEL = 1D_transversely_isotropic_prem OCEANS = .true. ELLIPTICITY = .true. TOPOGRAPHY = .true. GRAVITY = .true. ROTATION = .true. ATTENUATION = .true. ABSORBING_CONDITIONS = .false. RECORD_LENGTH_IN_MINUTES = 10.d0 PARTIAL_PHYS_DISPERSION_ONLY = .false. UNDO_ATTENUATION = .false. MEMORY_INSTALLED_PER_CORE_IN_GB = 4.d0 PERCENT_OF_MEM_TO_USE_PER_CORE = 85.d0 EXACT_MASS_MATRIX_FOR_ROTATION = .false. USE_LDDRK = .false. INCREASE_CFL_FOR_LDDRK = .true. RATIO_BY_WHICH_TO_INCREASE_IT = 1.5d0 MOVIE_SURFACE = .false. MOVIE_VOLUME = .false. MOVIE_COARSE = .false. NTSTEP_BETWEEN_FRAMES = 50 HDUR_MOVIE = 0.d0 MOVIE_VOLUME_TYPE = 2 MOVIE_TOP_KM = -100.0 MOVIE_BOTTOM_KM = 1000.0 MOVIE_WEST_DEG = -90.0 MOVIE_EAST_DEG = 90.0 MOVIE_NORTH_DEG = 90.0 MOVIE_SOUTH_DEG = -90.0 MOVIE_START = 0 MOVIE_STOP = 40000 SAVE_MESH_FILES = .false. NUMBER_OF_RUNS = 1 NUMBER_OF_THIS_RUN = 1 LOCAL_PATH = ./DATABASES_MPI LOCAL_TMP_PATH = ./DATABASES_MPI NTSTEP_BETWEEN_OUTPUT_INFO = 500 NTSTEP_BETWEEN_OUTPUT_SEISMOS = 5000000 NTSTEP_BETWEEN_READ_ADJSRC = 1000 USE_FORCE_POINT_SOURCE = .false. SAVE_SEISMOGRAMS_STRAIN = .false. SAVE_SEISMOGRAMS_IN_ADJOINT_RUN = .false. OUTPUT_SEISMOS_ASCII_TEXT = .false. OUTPUT_SEISMOS_SAC_ALPHANUM = .false. OUTPUT_SEISMOS_SAC_BINARY = .true. OUTPUT_SEISMOS_ASDF = .false. ROTATE_SEISMOGRAMS_RT = .false. WRITE_SEISMOGRAMS_BY_MASTER = .false. SAVE_ALL_SEISMOS_IN_ONE_FILE = .false. USE_BINARY_FOR_LARGE_FILE = .false. RECEIVERS_CAN_BE_BURIED = .true. PRINT_SOURCE_TIME_FUNCTION = .false. READ_ADJSRC_ASDF = .false. ANISOTROPIC_KL = .false. SAVE_TRANSVERSE_KL_ONLY = .false. SAVE_AZIMUTHAL_ANISO_KL_ONLY = .false. APPROXIMATE_HESS_KL = .false. USE_FULL_TISO_MANTLE = .false. SAVE_SOURCE_MASK = .false. SAVE_REGULAR_KL = .false. NUMBER_OF_SIMULTANEOUS_RUNS = 1 BROADCAST_SAME_MESH_AND_MODEL = .false. GPU_MODE = .true. GPU_RUNTIME = 1 GPU_PLATFORM = NVIDIA GPU_DEVICE = Tesla ADIOS_ENABLED = .true. ADIOS_FOR_FORWARD_ARRAYS = .true. ADIOS_FOR_MPI_ARRAYS = .true. ADIOS_FOR_ARRAYS_SOLVER = .true. ADIOS_FOR_SOLVER_MESHFILES = .true. ADIOS_FOR_AVS_DX = .true. ADIOS_FOR_KERNELS = .true. ADIOS_FOR_MODELS = .true. ADIOS_FOR_UNDO_ATTENUATION = .true. ```
Modules, Compiler & Packages ```bash module load openmpi/gcc/4.0.3rc1/64
module load cudatoolkit/10.0 # C/C++ compiler CC=gcc CXX=g++ MPICC=mpicc FCFLAGS="-g" \# Fortran compiler FC=gfortran MPIFC=mpif90 CFLAGS="" # CUDA (here CUDA 5 because my GPU cannot support more, poor boy) CUDA_WITH="--with-cuda=cuda8" ADIOS_LINK="https://users.nccs.gov/~pnorbert/adios-1.13.1.tar.gz" ADIOS_DESTDIR="${ADIOS_DIR}/build" ADIOS_WITH="--with-adios" ADIOS_CONFIG="$ADIOS_DESTDIR/bin/adios_config" ```
Configuration ```bash ./configure CC=$CC CXX=$CXX FC=$FC MPIFC=$MPIFC \ CFLAGS="$CFLAGS" FCLAGS="$FCFLAGS" \ $CUDA_WITH CUDA_LIB="$CUDA_LIB" \ $ADIOS_WITH ADIOS_CONFIG="$ADIOS_CONFIG" ```
lsawade commented 4 years ago

Update

I must have done something wrong with the master branch. The master branch does compile.

danielpeter commented 4 years ago

hi Lucas, what gcc version were you using?

could you provide the output for gfortran -v

1sbeller commented 4 years ago

Hi all,

Working on this on traverse. Disabling optimization ( -O0) with gcc makes it work... Seems like an aggresive optimization problem. Valgrind also makes it pass without problem...

FYI, the crash happens when meshing region 2 (liquid outer core).

gfortran -v gfortran -v Utilisation des specs internes. COLLECT_GCC=gfortran COLLECT_LTO_WRAPPER=/usr/libexec/gcc/ppc64le-redhat-linux/4.8.5/lto-wrapper Target: ppc64le-redhat-linux Configuré avec: ../configure --prefix=/usr --mandir=/usr/share/man --infodir=/usr/share/info --with-bugurl=http://bugzilla.redhat.com/bugzilla --enable-bootstrap --enable-shared --enable-threads=posix --enable-checking=release --with-system-zlib --enable-__cxa_atexit --disable-libunwind-exceptions --enable-gnu-unique-object --enable-linker-build-id --with-linker-hash-style=gnu --enable-languages=c,c++,objc,obj-c++,java,fortran,ada,go,lto --enable-plugin --enable-initfini-array --disable-libgcj --with-isl=/builddir/build/BUILD/gcc-4.8.5-20150702/obj-ppc64le-redhat-linux/isl-install --with-cloog=/builddir/build/BUILD/gcc-4.8.5-20150702/obj-ppc64le-redhat-linux/cloog-install --enable-gnu-indirect-function --enable-secureplt --with-long-double-128 --enable-targets=powerpcle-linux --disable-multilib --with-cpu-64=power8 --with-tune-64=power8 --build=ppc64le-redhat-linux Modèle de thread: posix gcc version 4.8.5 20150623 (Red Hat 4.8.5-37) (GCC)

1sbeller commented 4 years ago

It works also with -O1.

1sbeller commented 4 years ago

Update on GDB output: Program received signal SIGFPE, Arithmetic exception. 0x000000001009df94 in get_model (iregion_code=2, ispec=1, nspec=, idoubling=6, xstore=..., ystore=..., zstore=..., rmin=0.19164966253335425, rmax=0.54622508240464607, elem_in_crust=.FALSE., elem_in_mantle=.FALSE.) at src/meshfem3D/get_model.F90:174 174 + 5.d0vshvsh + (6.d0+4.d0eta_aniso)vsv*vsv)/15.d0) Missing separate debuginfos, use: debuginfo-install cyrus-sasl-lib-2.1.26-23.el7.ppc64le glibc-2.17-260.el7_6.6.ppc64le hwloc-libs-1.11.8-4.el7.ppc64le hwloc-plugins-1.11.8-4.el7.ppc64le keyutils-libs-1.5.8-3.el7.ppc64le krb5-libs-1.15.1-37.el7_6.ppc64le libcom_err-1.42.9-13.el7.ppc64le libevent-2.0.21-4.el7.ppc64le libgcc-4.8.5-37.el7_6.ppc64le libgfortran-4.8.5-37.el7_6.ppc64le libibumad-43.1.1.MLNX20190422.87b4d9b-0.1.46101.ppc64le libibverbs-41mlnx1-OFED.4.6.0.4.1.46101.ppc64le libmlx4-41mlnx1-OFED.4.5.0.0.3.46101.ppc64le libmlx5-41mlnx1-OFED.4.6.0.0.4.46101.ppc64le libnl3-3.2.28-4.el7.ppc64le libpciaccess-0.14-1.el7.ppc64le librdmacm-41mlnx1-OFED.4.6.0.0.1.46101.ppc64le libselinux-2.5-14.1.el7.ppc64le libtool-ltdl-2.4.2-22.el7_3.ppc64le libxml2-2.9.1-6.el7_2.3.ppc64le munge-libs-0.5.11-3.el7.ppc64le nss-softokn-freebl-3.36.0-6.el7_6.ppc64le numactl-libs-2.0.9-7.el7.ppc64le opensm-libs-5.4.0.MLNX20190422.ed81811-0.1.46101.ppc64le pcre-8.32-17.el7.ppc64le pmix-1.2.2-1.el7.ppc64le tcp_wrappers-libs-7.6-77.el7.ppc64le xz-libs-5.2.2-1.el7.ppc64le zlib-1.2.7-18.el7.ppc64le

1sbeller commented 4 years ago

Well problem is solved...

Actually the problem is here: lines 169-170 of get_model.F90 vs = sqrt(((1.d0-2.d0eta_aniso)vphvph + vpvvpv &

The square root argument might become slightly negative because of roundoff errors. Those may be compiler and architecture dependents. If it is slightly negative result of sqrt(-something) is complex. Some compilers trhow errors, gives zero or gives NaN. Here we compute vs in the fluid region thus it should be zero but not exactly because of roundoff errors.

The trick is to modify so that: vs = sqrt(abs((1.d0-2.d0eta_aniso)vphvph + vpvvpv &

Another possibility is to take max(argument, 0).

Stephen

danielpeter commented 4 years ago

hi Stephen, great spotting this one!

what happens if you use instead:

        if (iregion_code == IREGION_OUTER_CORE) then
          ! fluid with zero shear speed
          vs = 0.d0
        else
          vs = sqrt(((1.d0-2.d0*eta_aniso)*vph*vph + vpv*vpv &
                    + 5.d0*vsh*vsh + (6.d0+4.d0*eta_aniso)*vsv*vsv)/15.d0)
        endif

that is a bit more explicit when reading through.

1sbeller commented 4 years ago

I agree this would be more explicit, should I make this substitution?

danielpeter commented 4 years ago

yes, that would be great - thanks. -daniel

On 19 Jun 2020, at 8:38 PM, 1sbeller notifications@github.com wrote:

I agree this would be more explicit, should I make this substitution?

— You are receiving this because you commented. Reply to this email directly, view it on GitHub https://github.com/geodynamics/specfem3d_globe/issues/690#issuecomment-646810253, or unsubscribe https://github.com/notifications/unsubscribe-auth/ABMLXWMJ6ZJRYEF6HJUM5KTRXOWBNANCNFSM4OCI7KMA.

lsawade commented 4 years ago

Wow! This is great thanks both to you @danielpeter and you @1sbeller !

lsawade commented 4 years ago

I think we are ready to close this one with #693

Thanks, @1sbeller!