SPEED
MAKE_NH_Enhanced_initialise.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

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
 

Function/Subroutine Documentation

◆ make_nh_enhanced_initialise()

subroutine make_nh_enhanced_initialise ( integer*4  nn_loc,
integer*4  nmat_nhe,
integer*4, dimension(nmat_nhe)  val_nhe,
integer*4  nmat,
integer*4, dimension(nmat)  tag_mat,
integer*4  ne_loc,
integer*4  cs_nnz_loc,
integer*4, dimension(0:cs_nnz_loc)  cs_loc,
real*8, dimension(nn_loc)  xs_loc,
real*8, dimension(nn_loc)  ys_loc,
real*8, dimension(nn_loc)  zs_loc,
integer*4, intent(inout)  count,
integer*4, dimension(nn_loc), intent(inout)  node_nhe_flag,
integer*4  mpi_id,
integer*4  mpi_comm,
character*70  mpi_file 
)

...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]countno. of Spectral nodes where Nearest Neighbor search has to be performed
[out]xs_loc_nhex-coordinate of spectral nodes corresponding to NHE
[out]ys_loc_nhe
[out]ys_loc_nhe
[out]node_nhe_flagcontains material tags of blocks corresponding to each spectral node

Definition at line 52 of file MAKE_NH_Enhanced_initialise.f90.

57
59
60 implicit none
61
62 integer*4 :: nn_loc, mpi_id , nmat_nhe, nmat
63 integer*4, dimension(nmat_nhe) :: val_nhe
64 integer*4, dimension(nmat) :: tag_mat
65
66
67 integer*4 :: ne_loc, cs_nnz_loc
68 integer*4, dimension(0:cs_nnz_loc) :: cs_loc
69
70 real*8, dimension(nn_loc) :: xs_loc, ys_loc, zs_loc
71 integer*4, intent(inout) :: count
72 integer*4, dimension(nn_loc), intent(inout) :: node_nhe_flag
73 real(kdkind), dimension(:), allocatable :: xs_loc_nhe, ys_loc_nhe, zs_loc_nhe
74
75 real*8 :: t0, t1, time_elapsed
76 integer*4 :: i, j, ipt, inode, ie
77 integer*4 :: im, istart, iend, mpi_ierr, mpi_comm, unit_mpi
78
79 character*70 :: file_nhe_proc, file_nhe_new, mpi_file
80
81
82 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
83 ! Making GLL node list which falls inside NHE-specified Blocks
84 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
85
86 node_nhe_flag = 0
87
88 do ie = 1,ne_loc
89 im = cs_loc(cs_loc(ie -1))
90 istart = cs_loc(ie-1) + 1
91 iend = cs_loc(ie) - 1
92
93 do j = istart,iend
94 do i = 1,nmat_nhe
95
96 if (val_nhe(i) .eq. tag_mat(im)) then
97 ! Node Numbers of nodes where NHE has to be implemented
98 !inode = cs_loc(j) ! Local Node Number
99 node_nhe_flag(cs_loc(j)) = 999
100 else
101 if (node_nhe_flag(cs_loc(j)).ne.999) node_nhe_flag(cs_loc(j)) = im
102 endif
103
104 enddo ! i = 1,nmat_nhe
105 enddo
106 enddo !ie=1,ne_loc
107
108
109 ! Making List of Local Nodes and their coordinates where
110 ! Nearest Neighbor search has to be performed
111 count = 0
112 do inode =1,nn_loc
113 if ((node_nhe_flag(inode).eq.999)) count = count + 1
114 enddo
115
116 if (count.gt.0) then
117 allocate(xs_loc_nhe(count),ys_loc_nhe(count),zs_loc_nhe(count))
118 endif
119
120 count = 0
121 do inode =1,nn_loc
122 if ((node_nhe_flag(inode).eq.999)) then
123 count = count + 1
124 xs_loc_nhe(count) = xs_loc(inode)
125 ys_loc_nhe(count) = ys_loc(inode)
126 zs_loc_nhe(count) = zs_loc(inode)
127 endif
128 enddo
129
130 call mpi_barrier(mpi_comm, mpi_ierr)
131
132 file_nhe_proc = 'nhexyz000000.mpi'
133 unit_mpi = 1500 + mpi_id
134 if (mpi_id.lt. 10) then
135 write(file_nhe_proc(12:12),'(i1)') mpi_id
136 else if (mpi_id .lt. 100) then
137 write(file_nhe_proc(11:12),'(i2)') mpi_id
138 else if (mpi_id .lt. 1000) then
139 write(file_nhe_proc(10:12),'(i3)') mpi_id
140 else if (mpi_id .lt. 10000) then
141 write(file_nhe_proc(9:12),'(i4)') mpi_id
142 else if (mpi_id .lt. 100000) then
143 write(file_nhe_proc(8:12),'(i5)') mpi_id
144 else if (mpi_id .lt. 1000000) then
145 write(file_nhe_proc(7:12),'(i6)') mpi_id
146 endif
147
148 if(len_trim(mpi_file) .ne. 70) then
149 file_nhe_new = mpi_file(1:len_trim(mpi_file)) // '/' // file_nhe_proc
150 else
151 file_nhe_new = file_nhe_proc
152 endif
153
154 open(unit_mpi,file=file_nhe_new)
155 write(unit_mpi,*) count
156 if (count.gt.0) then
157 do i = 1,count
158 write(unit_mpi,*) xs_loc_nhe(i), ys_loc_nhe(i), zs_loc_nhe(i)
159 enddo
160 deallocate(xs_loc_nhe, ys_loc_nhe, zs_loc_nhe)
161 endif
162 close(unit_mpi)
163
164 call mpi_barrier(mpi_comm, mpi_ierr)
165

Referenced by make_nh_enhanced().

Here is the caller graph for this function: