SPEED
MAKE_NH_Enhanced.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine make_nh_enhanced ()
 ...Not-Honoring Enhanced (NHE) Implementation
 

Function/Subroutine Documentation

◆ make_nh_enhanced()

subroutine make_nh_enhanced

...Not-Honoring Enhanced (NHE) Implementation

Author
Srihari Sangaraju
Date
August, 2020
Version
1.0 @ Detailed Explanation
  1. User Can specify that a single block/multiple blocks in mesh has to use material properties as per provided Tomography text file
  2. For the each GLL node in NHE-specified blocks, this subroutine assigns Mechanical properties of the corresponding nearest point it can find in "Tomography text file"
  3. For the nodes in other blocks It assigns Mechanical properties specified in *.mate file Note: In this version, This subroutine will be run, only if atleast 1 block is specified to use NHE in *.mate file
Parameters
[in]loc_n_num.Global node number of 'i'th local node is loc_n_num(i)
[in]nn_loc.No. of nodes in Local/Current Partition
[in]nmat_nheNo. of Blocks specified with NHE case
[in]nhe_matTag/Labels of Blocks where NHE has to be implemented
[in]xs_locx-coordinate of spectral nodes
[in]ys_locy-coordinate of spectral nodes
[in]zs_locz-coordinate of spectral nodes
[out]lambda_nheLame coefficient lambda
[out]mu_nheLame coefficient mu
[out]rho_nhematerial density
[out]Qs_nhequality factor for S-waves (Standard Linear Solid damping)
[out]Qp_nhequality factor for P-waves (Standard Linear Solid damping)
[out]gamma_dmp_nhedamping coefficient gamma (Kosloff&Kosloff)

Definition at line 52 of file MAKE_NH_Enhanced.f90.

53
54 use speed_par
55
56 implicit none
57
58 include 'SPEED.MPI'
59
60 integer*4 :: count
61 integer*4, dimension(:), allocatable :: node_nhe_flag, NN_src_ind_loc
62
63
64 if (mpi_id.eq.0) write(*,'(A)')
65 if (mpi_id.eq.0) write(*,'(A)')'---------------Setup Not-Honoring Enhanced ---------------'
66
67
68 allocate(node_nhe_flag(nnod_loc))
72 count, &
73 node_nhe_flag, mpi_id, mpi_comm, mpi_file)
74
75
76 allocate(nn_src_ind_loc(count))
78 mpi_file, nn_src_ind_loc)
79
80
82 allocate(qs_nhe_el(nelem_loc), qp_nhe_el(nelem_loc)) ! gamma_nhe_el(nelem_loc))
83
86 node_nhe_flag, count, nn_src_ind_loc, qs, qp, &
89
90
91 deallocate(node_nhe_flag, nn_src_ind_loc)
92
93 if (mpi_id.eq.0) write(*,'(A)')'Completed.'
94
subroutine make_nh_enhanced_assign_prop(nn_loc, nmat, prop_mat, sdeg_mat, ne_loc, cs_nnz_loc, cs_loc, node_nhe_flag, count, nn_src_ind_loc, qs, qp, lambda_nhe, mu_nhe, rho_nhe, qs_nhe_el, qp_nhe_el, mpi_id, mpi_comm)
...Not-Honoring Enhanced (NHE) Implementation
subroutine make_nh_enhanced_nnsearch(nn_loc, count, mpi_id, mpi_np, mpi_comm, mpi_file, nn_src_ind_loc)
...Not-Honoring Enhanced (NHE) Implementation
subroutine make_nh_enhanced_initialise(nn_loc, nmat_nhe, val_nhe, nmat, tag_mat, ne_loc, cs_nnz_loc, cs_loc, xs_loc, ys_loc, zs_loc, count, node_nhe_flag, mpi_id, mpi_comm, mpi_file)
...Not-Honoring Enhanced (NHE) Implementation
Contains SPEED PARAMETERS used in (SPEED, READ_INPUT_FILES, MAKE_PARTION_AND_MPI_FILES,...
Definition MODULES.f90:207
real *8, dimension(:), allocatable qs
Definition MODULES.f90:437
integer *4 nmat
Definition MODULES.f90:269
real *8, dimension(:,:), allocatable prop_mat
Definition MODULES.f90:463
integer *4 con_nnz_loc
Definition MODULES.f90:269
character *70 mpi_file
Definition MODULES.f90:243
integer *4 nmat_nhe
Definition MODULES.f90:269
real *8, dimension(:), allocatable zz_spx_loc
Definition MODULES.f90:408
real *8, dimension(:), allocatable lambda_nhe
Definition MODULES.f90:446
real *8, dimension(:), allocatable yy_spx_loc
Definition MODULES.f90:408
integer *4, dimension(:), allocatable val_nhe
Definition MODULES.f90:379
real *8, dimension(:), allocatable qp
Definition MODULES.f90:437
integer *4 mpi_comm
Definition MODULES.f90:308
integer *4, dimension(:), allocatable sdeg_mat
Definition MODULES.f90:335
integer *4 nelem_loc
Definition MODULES.f90:269
integer *4, dimension(:), allocatable con_spx_loc
Definition MODULES.f90:322
integer *4, dimension(:), allocatable tag_mat
Definition MODULES.f90:335
integer *4 mpi_np
Definition MODULES.f90:308
real *4, dimension(:), allocatable qs_nhe_el
Definition MODULES.f90:447
integer *4 mpi_id
Definition MODULES.f90:308
integer *4 nnod_loc
Definition MODULES.f90:269
real *8, dimension(:), allocatable xx_spx_loc
Definition MODULES.f90:408
real *4, dimension(:), allocatable qp_nhe_el
Definition MODULES.f90:447
real *8, dimension(:), allocatable mu_nhe
Definition MODULES.f90:446
real *8, dimension(:), allocatable rho_nhe
Definition MODULES.f90:446

References speed_par::con_nnz_loc, speed_par::con_spx_loc, speed_par::lambda_nhe, make_nh_enhanced_assign_prop(), make_nh_enhanced_initialise(), make_nh_enhanced_nnsearch(), speed_par::mpi_comm, speed_par::mpi_file, speed_par::mpi_id, speed_par::mpi_np, speed_par::mu_nhe, speed_par::nelem_loc, speed_par::nmat, speed_par::nmat_nhe, speed_par::nnod_loc, speed_par::prop_mat, speed_par::qp, speed_par::qp_nhe_el, speed_par::qs, speed_par::qs_nhe_el, speed_par::rho_nhe, speed_par::sdeg_mat, speed_par::tag_mat, speed_par::val_nhe, speed_par::xx_spx_loc, speed_par::yy_spx_loc, and speed_par::zz_spx_loc.

Here is the call graph for this function: