WAVEWATCH III  beta 0.0.1
yowpdlibmain Module Reference

Functions/Subroutines

subroutine, public initfromgriddim (MNP, MNE, INE_global, secDim, MPIcomm)
 
subroutine real_mpi_barrier_pdlib (TheComm, string)
 
subroutine assignmesh (MNP, MNE)
 
subroutine prepartition
 pre-partition the mesh just divide the mesh into nTasks parts and create a premature iplg alter: np_perProc, np_perProcSum, np, iplg More...
 
subroutine findconnnodes (INE_global)
 create the connected Nodes array loop over all elements and their nodes. More...
 
subroutine runparmetis (MNP)
 Collect all data for parmetis und partition the mesh after that, we know for every node the domain ID alter: t_Node::domainID. More...
 
subroutine postpartition
 recalculate some variables parmetis change the number of sides per domain. More...
 
subroutine findghostnodes
 find the ghost nodes of the local domain alter: ng, ghosts(), ghostlg, ghostgl, npa, iplg More...
 
subroutine findconndomains
 find the number of connected domains and their ghosts 1) Iterate over all ghost nodes and look at their thread id to find # neighbor domains 2) assign the ghost nodes to their domains alter: neighborDomains(), nConnDomains More...
 
subroutine exchangeghostids
 exchange Ghost Ids so every thread knows which nodes he has to send to the other parition. More...
 
subroutine postpartition2 (INE_global)
 this collects all data which depends on ghost information alter: ne, INE, x, y, z, ielg More...
 
subroutine computetria_ien_si_ccon
 
subroutine element_crosses_dateline (RX1, RX2, RX3, CROSSES_DATELINE)
 
subroutine, public finalizepd ()
 

Function/Subroutine Documentation

◆ assignmesh()

subroutine yowpdlibmain::assignmesh ( integer, intent(in)  MNP,
integer, intent(in)  MNE 
)
Parameters
[in]MNPnumber of nodes global
[in]XPnode X value
[in]XYnode Y value
[in]DEPnode Z value
[in]MNEnumber of element global
[in]INEelement array alter: np_global, nodes_global(), ne_global, elements(), INE_global

Definition at line 258 of file yowpdlibmain.F90.

258  use yownodepool, only: nodes_global, np_global
259  use yowelementpool, only: ne_global
260  use yowerr, only: parallel_abort
261 
262  integer, intent(in) :: MNP, MNE
263  !integer, intent(in) :: INE(3,MNE)
264  integer :: stat, i
265 
266  np_global=mnp
267  if(allocated(nodes_global)) deallocate(nodes_global)
268  allocate(nodes_global(np_global), stat=stat);
269  if(stat/=0) CALL abort('nodes_global() allocate failure')
270 
271  do i =1, np_global
272  nodes_global(i)%id_global = i
273  end do
274 
275  ne_global=mne
276 
277  !if(allocated(INE_global)) deallocate(INE_global)
278  !allocate(INE_global(3, ne_global), stat=stat);
279  !if(stat/=0) CALL ABORT('INE_global allocate failure')
280  !INE_global = INE

References yowelementpool::ne_global, yownodepool::nodes_global, yownodepool::np_global, and yowerr::parallel_abort().

Referenced by initfromgriddim().

◆ computetria_ien_si_ccon()

subroutine yowpdlibmain::computetria_ien_si_ccon

Definition at line 1304 of file yowpdlibmain.F90.

1304  use yowelementpool, only: ne, ne_global, ine, ielg
1306  use yowerr, only: parallel_abort
1307  use yowdatapool, only: myrank
1308  use yownodepool, only: np_global, np, iplg, t_node, ghostlg, ng, npa
1309  use yownodepool, only: x, y, z, pdlib_si, pdlib_ien, pdlib_tria, pdlib_ccon, pdlib_tria03
1310 
1311  integer I1, I2, I3, stat, IE, NI(3)
1312  real :: DXP1, DXP2, DXP3, DYP1, DYP2, DYP3, DBLTMP, TRIA03
1313  logical :: CROSSES_DATELINE
1314 
1315  allocate(pdlib_si(npa), pdlib_ccon(npa), pdlib_ien(6,ne), pdlib_tria(ne), pdlib_tria03(ne), stat=stat)
1316  if(stat/=0) call parallel_abort('SI allocation failure')
1317 
1318  pdlib_si(:) = 0.0d0 ! Median Dual Patch Area of each Node
1319  pdlib_ccon(:) = 0 ! Number of connected Elements
1320  DO ie = 1 , ne
1321  i1 = ine(1,ie)
1322  i2 = ine(2,ie)
1323  i3 = ine(3,ie)
1324  ni = ine(:,ie)
1325 
1326  dxp1=x(i2) - x(i1)
1327  dyp1=y(i2) - y(i1)
1328  dxp2=x(i3) - x(i2)
1329  dyp2=y(i3) - y(i2)
1330  dxp3=x(i1) - x(i3)
1331  dyp3=y(i1) - y(i3)
1332  CALL element_crosses_dateline(dxp1, dxp2, dxp3, crosses_dateline)
1333  IF (crosses_dateline) THEN
1334  CALL correct_dx_gt180(dxp1)
1335  CALL correct_dx_gt180(dxp2)
1336  CALL correct_dx_gt180(dxp3)
1337  ENDIF
1338 
1339  pdlib_ien(1,ie) = - dyp2
1340  pdlib_ien(2,ie) = dxp2
1341  pdlib_ien(3,ie) = - dyp3
1342  pdlib_ien(4,ie) = dxp3
1343  pdlib_ien(5,ie) = - dyp1
1344  pdlib_ien(6,ie) = dxp1
1345  dbltmp = (dxp3*dyp1 - dyp3*dxp1)*0.5
1346  pdlib_tria(ie) = dbltmp
1347  IF (pdlib_tria(ie) .lt. tiny(1.)) THEN
1348  WRITE(*,*) pdlib_ien(:,ie)
1349  WRITE(*,*)
1350  WRITE(*,*) 'AREA SMALLER ZERO IN PDLIB', ie, ne, pdlib_tria(ie)
1351  stop
1352  ENDIF
1353 
1354  pdlib_ccon(i1) = pdlib_ccon(i1) + 1
1355  pdlib_ccon(i2) = pdlib_ccon(i2) + 1
1356  pdlib_ccon(i3) = pdlib_ccon(i3) + 1
1357  tria03 = pdlib_tria(ie)/3.d0
1358  pdlib_si(i1) = pdlib_si(i1) + tria03
1359  pdlib_si(i2) = pdlib_si(i2) + tria03
1360  pdlib_si(i3) = pdlib_si(i3) + tria03
1361  pdlib_tria03(ie) = tria03
1362  ENDDO
1363  CALL pdlib_exchange1dreal(pdlib_si)

References element_crosses_dateline(), yownodepool::ghostlg, yowelementpool::ielg, yowelementpool::ine, yownodepool::iplg, yowdatapool::myrank, yowelementpool::ne, yowelementpool::ne_global, yownodepool::ng, yownodepool::np, yownodepool::np_global, yownodepool::npa, yowerr::parallel_abort(), yownodepool::pdlib_ccon, yowexchangemodule::pdlib_exchange1dreal(), yownodepool::pdlib_ien, yownodepool::pdlib_si, yownodepool::pdlib_tria, yownodepool::pdlib_tria03, yownodepool::x, yownodepool::y, and yownodepool::z.

Referenced by initfromgriddim().

◆ element_crosses_dateline()

subroutine yowpdlibmain::element_crosses_dateline ( real(rkind), intent(in)  RX1,
real(rkind), intent(in)  RX2,
real(rkind), intent(in)  RX3,
logical, intent(out)  CROSSES_DATELINE 
)

Definition at line 1369 of file yowpdlibmain.F90.

1369  ! Purpose: understanding if an element crosses the dateline.
1370  ! An element crossing the dateline has, e.g. a node with lon < 180
1371  ! and another 2 with lon > -180
1372 
1373  REAL(rkind), INTENT(IN) :: RX1, RX2, RX3
1374  LOGICAL, INTENT(OUT) :: CROSSES_DATELINE
1375  INTEGER :: R1GT180, R2GT180, R3GT180
1376  r1gt180 = merge(1, 0, abs(rx1).GT.180)
1377  r2gt180 = merge(1, 0, abs(rx2).GT.180)
1378  r3gt180 = merge(1, 0, abs(rx3).GT.180)
1379  ! if R1GT180+R2GT180+R3GT180 .eq. 0 the element does not cross the dateline
1380  ! if R1GT180+R2GT180+R3GT180 .eq. 1 the element contains the pole
1381  ! if R1GT180+R2GT180+R3GT180 .eq. 2 the element crosses the dateline
1382  crosses_dateline = r1gt180+r2gt180+r3gt180 .EQ. 2

References yownodepool::ghostlg, yowelementpool::ielg, yowelementpool::ine, yownodepool::iplg, yowdatapool::myrank, yowelementpool::ne, yowelementpool::ne_global, yownodepool::ng, yownodepool::nodes_global, yownodepool::np, yownodepool::np_global, yownodepool::npa, yowerr::parallel_abort(), yownodepool::pdlib_ccon, yownodepool::pdlib_i_diag, yownodepool::pdlib_ia, yownodepool::pdlib_ia_p, yownodepool::pdlib_ie_cell, yownodepool::pdlib_ie_cell2, yownodepool::pdlib_ja, yownodepool::pdlib_ja_ie, yownodepool::pdlib_ja_p, yownodepool::pdlib_nnz, yownodepool::pdlib_pos_cell, yownodepool::pdlib_pos_cell2, and yownodepool::pdlib_posi.

Referenced by computetria_ien_si_ccon().

◆ exchangeghostids()

subroutine yowpdlibmain::exchangeghostids

exchange Ghost Ids so every thread knows which nodes he has to send to the other parition.

every parition has a list of ghost nodes from the other partition. this data is stored in the neighborDomains variable

Todo:
make a better explanation 1) send to the neighbor domain which ghost nodes we want from him 2) receive from the neighbor domain which ghost we must send to him alter: neighborDomains()%{numNodesToSend, nodesToSend},

Definition at line 1052 of file yowpdlibmain.F90.

1052  use yowerr
1053  use yownodepool, only: np, t_node, nodes
1054  use yowdatapool, only: ntasks, myrank, comm
1055  use yowexchangemodule, only: neighbordomains, nconndomains, creatempitypes
1056  use mpi
1057 
1058  integer :: i, j, k
1059  integer :: ierr
1060  ! uniq tag that identify the sender and which information he sends
1061  integer :: tag
1062  ! we use non-blocking send and recv subroutines
1063  ! store the send status
1064  integer :: sendRequest(nConnDomains)
1065  ! store the revc status
1066  integer :: recvRequest(nConnDomains)
1067  ! status to verify if one communication fails or not
1068  integer :: status(MPI_STATUS_SIZE, nConnDomains);
1069 
1070 
1071  type(t_node), pointer :: node
1072 
1073  ! send to all domain neighbors how many ghosts nodes we want from him and which ones
1074  do i=1, nconndomains
1075  ! create a uniq tag for this domain
1076  tag = neighbordomains(i)%domainID*10 + 1
1077  ! send to the neighbor how many ghost nodes we want from him
1078  call mpi_isend(neighbordomains(i)%numNodesToReceive, &
1079  1, &
1080  mpi_int, &
1081  neighbordomains(i)%domainID-1, &
1082  tag, &
1083  comm, &
1084  sendrequest(i), &
1085  ierr);
1086  if(ierr/=mpi_success) then
1087  write(*,*) "mpi send failure"
1088  endif
1089 
1090  tag = neighbordomains(i)%domainID*10 + 2
1091  ! send to the neighbor which ghost nodes we want from him
1092  call mpi_isend(neighbordomains(i)%nodesToReceive, &
1093  neighbordomains(i)%numNodesToReceive, &
1094  mpi_int, &
1095  neighbordomains(i)%domainID-1, &
1096  tag, &
1097  comm, &
1099  sendrequest(i), &
1100  ierr);
1101  if(ierr/=mpi_success) then
1102  write(*,*) "mpi send failure"
1103  endif
1104 
1105  ! receive from neighbor how many ghost nodes we have to send him
1106  tag = (myrank+1)*10 + 1
1107  call mpi_irecv(neighbordomains(i)%numNodesToSend, &
1108  1, &
1109  mpi_int, &
1110  neighbordomains(i)%domainID-1, &
1111  tag, &
1112  comm, &
1113  recvrequest(i), &
1114  ierr)
1115  if(ierr/=mpi_success) then
1116  write(*,*) "mpi recv failure"
1117  endif
1118  end do
1119 
1120  ! wait for communication end
1121  call mpi_waitall(nconndomains, recvrequest, status, ierr)
1122 
1123  ! test for all neighbor domains
1124  do i=1, nconndomains
1125  ! test if the neighbor wants more ghost nodes than we have
1126  if(neighbordomains(i)%numNodesToSend > np) then
1127  write(*,'(i5, a, i5, a, i5,a, i5, a, /)', advance='no') myrank, " ERROR neighbordomain ", neighbordomains(i)%domainID, &
1128  " wants ", neighbordomains(i)%numNodesToSend, &
1129  " nodes, but we have only ", np, " nodes"
1130  CALL abort("")
1131  end if
1132  end do
1133 
1134  ! receive from all neighbor domains which nodes we must send him
1135  do i=1, nconndomains
1136  allocate(neighbordomains(i)%nodesToSend(neighbordomains(i)%numNodesToSend))
1137  neighbordomains(i)%nodesToSend = 0
1138 
1139  ! receive from neighbor which nodes we must send
1140  tag = (myrank+1)*10 + 2
1141  call mpi_irecv(neighbordomains(i)%nodesToSend, &
1142  neighbordomains(i)%numNodesToSend, &
1143  mpi_int, &
1144  neighbordomains(i)%domainID-1, &
1145  tag, &
1146  comm, &
1147  recvrequest(i), &
1148  ierr)
1149  if(ierr/=mpi_success) then
1150  CALL parallel_abort("mpi recv failure", ierr)
1151  endif
1152  end do
1153 
1154  ! wait for communication end
1155  call mpi_waitall(nconndomains, recvrequest, status, ierr)
1156 
1157  ! test for all neighbor domains
1158  do i=1, nconndomains
1159  ! test if the neighbor wants nodes that we don't own
1160  outerloop: do j=1, neighbordomains(i)%numNodesToSend
1161  ! compare with all local nodes
1162  do k=1, np
1163  node => nodes(k)
1164  if(node%id_global == neighbordomains(i)%nodesToSend(j)) then
1165  cycle outerloop
1166  end if
1167  end do
1168  write(*,*) myrank, "Neighbordomain", neighbordomains(i)%domainID, &
1169  " want Node", neighbordomains(i)%nodesToSend(j), &
1170  " but we don't own this node"
1171  stop
1172  end do outerloop
1173  end do
1174 
1175  call creatempitypes()

References yowerr::abort(), yowdatapool::comm, yowexchangemodule::creatempitypes(), yowdatapool::myrank, yowexchangemodule::nconndomains, yowexchangemodule::neighbordomains, yownodepool::nodes(), yownodepool::np, yowdatapool::ntasks, and yowerr::parallel_abort().

Referenced by initfromgriddim().

◆ finalizepd()

◆ findconndomains()

subroutine yowpdlibmain::findconndomains

find the number of connected domains and their ghosts 1) Iterate over all ghost nodes and look at their thread id to find # neighbor domains 2) assign the ghost nodes to their domains alter: neighborDomains(), nConnDomains

Definition at line 968 of file yowpdlibmain.F90.

968  use yowerr, only: parallel_abort
969  use yownodepool, only: ghosts, ng, t_node
970  use yowdatapool, only: ntasks, myrank
971  use yowexchangemodule, only: neighbordomains, initnbrdomains
972 
973  integer :: i, stat, itemp
974  type(t_Node), pointer :: ghost
975 
976  ! # of ghost per neighbor domain
977  integer, allocatable :: numberGhostPerNeighborDomainTemp(:)
978  ! look up table. domainID to neighbor (ID)
979  integer, allocatable :: domainID2NeighborTemp(:)
980 
981  ! Part 1) find # neighbor domains
982 
983  ! allocate this array with a fixed size of nTasks. even if we not have
984  ! so many neighbor domains, nTasks will never be very large
985  allocate(numberghostperneighbordomaintemp(ntasks), stat=stat)
986  if(stat/=0) call parallel_abort('numberGhostPerNeighborDomainTemp allocation failure')
987  numberghostperneighbordomaintemp = 0
988 
989  allocate(domainid2neighbortemp(ntasks), stat=stat)
990  if(stat/=0) call parallel_abort('domainID2NeighborTemp allocation failure')
991  domainid2neighbortemp = 0
992 
993 
994  ! iterate over all ghost nodes an get their thread id
995  itemp = 0
996  do i = 1, ng
997  ghost => ghosts(i)
998 
999  ! sum how many ghost belongs to the ghost domainID
1000  numberghostperneighbordomaintemp(ghost%domainID) = numberghostperneighbordomaintemp(ghost%domainID) + 1
1001 
1002  ! check if this ghost domainID is allready in the domains list
1003  if(domainid2neighbortemp(ghost%domainID) /= 0) then
1004  ! yes we have allready insert this domain id
1005  cycle
1006  end if
1007 
1008  ! no we dont know this domain id. insert it
1009  itemp = itemp + 1
1010  domainid2neighbortemp(ghost%domainID) = itemp
1011  end do
1012 
1013  ! Part 2) assign the ghost nodes to their domains
1014  call initnbrdomains(itemp)
1015 
1016  do i = 1, ntasks
1017  if(numberghostperneighbordomaintemp(i) /= 0) then
1018  neighbordomains(domainid2neighbortemp(i) )%domainID = i
1019 
1020  allocate(neighbordomains(domainid2neighbortemp(i))%nodesToReceive(numberghostperneighbordomaintemp(i)), stat=stat)
1021  if(stat/=0) call parallel_abort('neighborDomains%ghosts allocation failure')
1022  neighbordomains(domainid2neighbortemp(i))%nodesToReceive = 0
1023  end if
1024  end do
1025 
1026  do i = 1, ng
1027  ghost => ghosts(i)
1028  itemp = domainid2neighbortemp(ghost%domainID)
1029 
1030  neighbordomains(itemp)%numNodesToReceive = neighbordomains(itemp)%numNodesToReceive + 1
1031  neighbordomains(itemp)%nodesToReceive(neighbordomains(itemp)%numNodesToReceive) = ghost%id_global
1032  end do
1033 
1034  if(allocated(numberghostperneighbordomaintemp)) deallocate(numberghostperneighbordomaintemp)
1035  if(allocated(domainid2neighbortemp)) deallocate(domainid2neighbortemp)

References yownodepool::ghosts(), yowexchangemodule::initnbrdomains(), yowdatapool::myrank, yowexchangemodule::neighbordomains, yownodepool::ng, yowdatapool::ntasks, and yowerr::parallel_abort().

Referenced by initfromgriddim().

◆ findconnnodes()

subroutine yowpdlibmain::findconnnodes ( integer, dimension(3,ne_global), intent(in)  INE_global)

create the connected Nodes array loop over all elements and their nodes.

get then the neighbor nodes finally calculate the number of sides alter: maxConnNodes, connNodes_data, ns, ns_global, nodenConnNodes

Todo:
we allocate more than we really need

Definition at line 342 of file yowpdlibmain.F90.

342  use yowerr, only: parallel_abort
343  use yownodepool, only: np, np_global, nodes_global, nodes, maxconnnodes, t_node, connnodes_data
344  use yowelementpool, only: ne_global
345  use yowsidepool, only: ns, ns_global
346 
347  integer, intent(in) :: INE_global(3,ne_global)
348  integer :: i, j, stat
349  type(t_Node), pointer :: node
350  integer JPREV, JNEXT
351 
352  ! Loop over all nlements
353  ! look at their nodes
354  ! get the node 1, insert node 2 and 3 into the connected nodes array
355  ! do that for node 2 and 3 again
356 
357  ! implementation is some different
358  ! loop over alle elements to get the # of connected nodes
359  ! allocate space
360  ! loop a second time to insert the connected nodes
361 
362  ! first loop
363  do i = 1, ne_global
364  do j = 1, 3
365  node => nodes_global(ine_global(j,i))
366  call node%insertConnNode()
367  call node%insertConnNode()
368  end do
369  end do
370 
371  maxconnnodes = maxval(nodes_global(:)%nConnNodes)
372  nodes_global(:)%nConnNodes = 0
373 
374  ! allocate space
376  if(allocated(connnodes_data)) deallocate(connnodes_data)
377  allocate(connnodes_data(np_global, maxconnnodes), stat = stat)
378  if(stat/=0) call parallel_abort('connNodes allocation failure')
379 
380  ! second loop
381  do i = 1, ne_global
382  DO j=1,3
383  IF (j .eq. 3) THEN
384  jnext = 1
385  ELSE
386  jnext = j + 1
387  END IF
388  IF (j .eq. 1) THEN
389  jprev = 3
390  ELSE
391  jprev = j - 1
392  END IF
393  node => nodes_global(ine_global(j,i))
394  call node%insertConnNode(ine_global(jnext,i))
395  call node%insertConnNode(ine_global(jprev,i))
396  END DO
397 
398  end do
399 
400  ns = 0
401  ! calc # sides local
402  do i = 1, np
403  node => nodes(i)
404  ns = ns + node%nConnNodes
405  end do
406 
407  do i = 1, np_global
408  node => nodes_global(i)
409  ns_global = ns_global + node%nConnNodes
410  end do

References yownodepool::connnodes_data, yownodepool::maxconnnodes, yowelementpool::ne_global, yownodepool::nodes(), yownodepool::nodes_global, yownodepool::np, yownodepool::np_global, yowsidepool::ns, yowsidepool::ns_global, and yowerr::parallel_abort().

Referenced by initfromgriddim().

◆ findghostnodes()

subroutine yowpdlibmain::findghostnodes

find the ghost nodes of the local domain alter: ng, ghosts(), ghostlg, ghostgl, npa, iplg

temporary hold the ghost numbers

Todo:
make this faster. dont check all neighbors from all local nodes. mark if an node has already been checked.

Definition at line 845 of file yowpdlibmain.F90.

845  use yowerr, only: parallel_abort
846  use yowdatapool, only: myrank
847  use yownodepool, only: t_node, np, nodes, ghosts, nodes_global, ng, ghostlg, ghostgl, npa, np_global, iplg
848 
849  integer :: i, j, k, stat
850  type(t_Node), pointer :: node, nodeNeighbor, nodeGhost
852  integer, save, allocatable :: ghostTemp(:)
853 #ifdef W3_DEBUGINIT
854  print *, 'Passing in findGhostNodes'
855 #endif
856 
857  ! iterate over all local nodes and look at their neighbors
858  ! has the neighbor another domain id, than it is a ghost
859 
860  ! implementation is some different
861  ! loop over all nodes to get the # of ghost nodes
862  ! allocate space
863  ! loop a second time to insert the ghost nodes
865 
866  ! first loop. find out how many ghost nodes we have (with double entries)
867  ng = 0
868  do i = 1, np
869  node => nodes(i)
870  do j = 1, node%nConnNodes
871  nodeneighbor => node%connNodes(j)
872  if(nodeneighbor%domainID /= node%domainID) then
873  ! yes, we found a ghost
874  ng = ng + 1
875  end if
876  end do
877  end do
878 
879  allocate(ghosttemp(ng), stat=stat)
880  if(stat/=0) call parallel_abort('ghostTemp allocation failure')
881 
882 #ifdef W3_DEBUGINIT
883  print *, 'np_global=', np_global
884 #endif
885  IF (allocated(ghostgl)) THEN
886  print *, 'ghostgl is already allocated'
887  END IF
888  allocate(ghostgl(np_global), stat=stat)
889  if(stat/=0) call parallel_abort('ghostgl allocation failure')
890  ghostgl = 0
891 
892  ! second loop. fill ghostlg. ignore double entries
893  ng = 0
894  ! iterate over all local nodes
895  do i = 1, np
896  node => nodes(i)
897 
898  ! check their neighbors
899  secondloop: do j = 1, node%nConnNodes
900  nodeneighbor => node%connNodes(j)
901 
902  if(nodeneighbor%domainID /= node%domainID) then
903  ! yes, we found a ghost
904  ! check if this ghost is allready in the ghost list
905  do k = 1, ng
906  nodeghost => nodes_global(ghosttemp(k))
907  if(nodeneighbor%id_global == nodeghost%id_global) then
908  ! yes, we allready know this ghost.
909  ! check the next neighbor
910  cycle secondloop
911  end if
912  end do
913 
914  ! no we don't know this ghost. insert it
915  ng = ng + 1
916  ghosttemp(ng) = nodeneighbor%id_global
917  ghostgl(nodeneighbor%id_global) = ng
918  end if
919  end do secondloop
920  end do
921 
922  ! reallocate the gosttemp array becouse it is longer then the new ng
923  if(allocated(ghostlg)) deallocate(ghostlg)
924  allocate(ghostlg(ng), stat=stat)
925  if(stat/=0) call parallel_abort('ghostlg allocation failure')
926  ghostlg = ghosttemp(1:ng)
927  deallocate(ghosttemp)
928 
929  npa = np + ng
930 
931  ! check if ghostlg contains only uniqe values
932  do i=1, ng
933  do j=i+1, ng
934  if(ghostlg(i) == ghostlg(j)) then
935  write(*,*) "double global ghost id in ghostlg(i,j)", i, j
936  stop "double global ghost id in ghostlg(i,j)"
937  endif
938  end do
939  end do
940 
941 
942  ! create the new node local to global mapping iplg with ghost. final one.
943  if(allocated(iplg)) deallocate(iplg)
944  allocate(iplg(npa), stat=stat)
945  if(stat/=0) call parallel_abort('iplg second allocation failure')
946  iplg = 0
947 
948  j = 1
949  do i = 1, np_global
950  node => nodes_global(i)
951  if(node%domainID == myrank+1) then
952  iplg(j) = i
953  j = j + 1
954  endif
955  end do
956 
957  iplg(np+1: npa) = ghostlg(1:ng)

References yownodepool::ghostgl, yownodepool::ghostlg, yownodepool::ghosts(), yownodepool::iplg, yowdatapool::myrank, yownodepool::ng, yownodepool::nodes(), yownodepool::nodes_global, yownodepool::np, yownodepool::np_global, yownodepool::npa, and yowerr::parallel_abort().

Referenced by initfromgriddim().

◆ initfromgriddim()

subroutine, public yowpdlibmain::initfromgriddim ( integer, intent(in)  MNP,
integer, intent(in)  MNE,
integer, dimension(3,mne), intent(in)  INE_global,
integer, intent(in)  secDim,
integer, intent(in)  MPIcomm 
)

Definition at line 66 of file yowpdlibmain.F90.

66  use yowdatapool, only: myrank, debugprepartition, debugpostpartition
67  use yownodepool, only: np_global, np, np_perprocsum, ng, ipgl, iplg, npa
68  use yowelementpool, only: ne_global,ne
69  use yowsidepool, only: ns, ns_global
70  use yowexchangemodule, only: nconndomains, setdimsize
72 
73  integer, intent(in) :: MNP, MNE
74  integer, intent(in) :: INE_global(3,MNE)
75  integer, intent(in) :: secDim
76  integer, intent(in) :: MPIcomm
77  integer :: istat, memunit
78 
79  ! note: myrank=0 until after initMPI is called, so only rank=0 file
80  ! contains the 'section 1' information
81  memunit = 70000+myrank
82 
83  call setdimsize(secdim)
84  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 1')
85 #ifdef W3_DEBUGINIT
86  print *, '1: MPIcomm=', mpicomm
87 #endif
88  call initmpi(mpicomm)
89 
90  memunit = 70000+myrank
91  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 2')
92 #ifdef W3_DEBUGINIT
93  print *, '2: After initMPI'
94 #endif
95  call assignmesh(mnp, mne)
96  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 3')
97 #ifdef W3_DEBUGINIT
98  print *, '3: After assignMesh'
99 #endif
100  call prepartition()
101  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 4')
102 #ifdef W3_DEBUGINIT
103  print *, '3: After prePartition'
104 #endif
105  call findconnnodes(ine_global)
106  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 5')
107 #ifdef W3_DEBUGINIT
108  print *, '4: After findConnNodes'
109 #endif
110  if(debugprepartition) then
111  if(myrank == 0) then
112  write(*,*) "pre-partition"
113  write(*,*) "# Nodes: ", np_global
114  write(*,*) "# Elements: ", ne_global
115  write(*,*) "# Sides ", ns_global
116  write(*,*) "np_perProcSum :", np_perprocsum
117  end if
118  write(*,*) "Thread", myrank, "# local nodes ", np
119  write(*,*) "Thread", myrank, "# local sides" , ns
120  endif
121 
122  ! call writeMesh()
123 #ifdef W3_DEBUGINIT
124  print *, '4.1: After findConnNodes'
125 #endif
126  ! CALL REAL_MPI_BARRIER_PDLIB(MPIcomm, "Before call to runParmetis")
127  call runparmetis(mnp)
128  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 6')
129  ! CALL REAL_MPI_BARRIER_PDLIB(MPIcomm, "After call to runParmetis")
130 #ifdef W3_DEBUGINIT
131  print *, '5: After runParmetis'
132 #endif
133  call postpartition
134  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 7')
135 #ifdef W3_DEBUGINIT
136  print *, 'Before findGhostNodes'
137 #endif
138  call findghostnodes
139  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 8')
140 
141  call findconndomains
142  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 9')
143 
144  call exchangeghostids
145  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 10')
146 
147  call postpartition2(ine_global)
148  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 11')
149 
150  call initrankmodule
151  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 12')
152 
153  call computetria_ien_si_ccon
154  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 13')
155 
156  call computeia_ja_posi_nnz
157  call print_memcheck(memunit, 'memcheck_____:'//' WW3_PDLIB SECTION 14')
158 
159  if(debugpostpartition) then
160  if(myrank == 0) then
161  write(*,*) "New data after partition"
162  write(*,*) "# Nodes: ", np_global
163  write(*,*) "# Elements: ", ne_global
164  write(*,*) "# Sides ", ns_global
165  write(*,*) "np_perProcSum :", np_perprocsum
166  end if
167  write(*,*) "Thread", myrank, "# local elements ", ne
168  write(*,*) "Thread", myrank, "# local nodes ", np
169  write(*,*) "Thread", myrank, "# local sides" , ns
170  write(*,*) "Thread", myrank, "# of ghosts", ng
171  write(*,*) "Thread", myrank, "# of neighbor domains", nconndomains
172  endif

References assignmesh(), computetria_ien_si_ccon(), yowdatapool::debugpostpartition, yowdatapool::debugprepartition, exchangeghostids(), findconndomains(), findconnnodes(), findghostnodes(), yowrankmodule::initrankmodule(), yownodepool::ipgl, yowrankmodule::ipgl_npa, yownodepool::iplg, yowdatapool::myrank, yowexchangemodule::nconndomains, yowelementpool::ne, yowelementpool::ne_global, yownodepool::ng, yownodepool::np, yownodepool::np_global, yownodepool::np_perprocsum, yownodepool::npa, yowsidepool::ns, yowsidepool::ns_global, postpartition(), postpartition2(), prepartition(), w3servmd::print_memcheck(), runparmetis(), and yowexchangemodule::setdimsize().

Referenced by pdlib_w3profsmd::pdlib_init().

◆ postpartition()

subroutine yowpdlibmain::postpartition

recalculate some variables parmetis change the number of sides per domain.

So recalculate some variables alter: np, ns, np_perProc, np_perProcSum, iplg, ipgl, nodes_globalid

Note
connNodes_data(:) has not changed after the call to parmetis.

Definition at line 782 of file yowpdlibmain.F90.

782  use yowerr, only: parallel_abort
783  use yowdatapool, only: myrank, ntasks
784  use yownodepool, only: np_global, np, nodes_global, nodes, np_perproc, np_perprocsum, iplg, ipgl, t_node
785  use yowsidepool, only: ns
786 
787  integer :: i, j, stat
788  type(t_Node), pointer :: node
789 
790  ! determine how many nodes now belong to which thread
791  ! and set the nodes local id
792  np_perproc = 0
793  do i = 1, np_global
794  ! fortran counts from 1. np_perProc from 0
795  np_perproc(nodes_global(i)%domainID-1) = np_perproc(nodes_global(i)%domainID-1)+1
796  ! set the new local id
797  nodes_global(i)%id = np_perproc(nodes_global(i)%domainID-1)
798  end do
799 
800  np_perprocsum(0) = 0
801  do i = 1, ntasks-1
802  np_perprocsum(i) = np_perprocsum(i-1) + np_perproc(i-1)
803  end do
804 
805  np = np_perproc(myrank)
806 
807  ! create the new node local to global mapping iplg. This is not the final one.
808  ! create a iplg with ghost nodes in findGhostNodes
809  if(allocated(iplg)) deallocate(iplg)
810  allocate(iplg(np), stat=stat)
811  if(stat/=0) call parallel_abort('iplg second allocation failure')
812  iplg = 0
813 
814  if(allocated(ipgl)) deallocate(ipgl)
815  allocate(ipgl(np_global), stat=stat)
816  if(stat/=0) call parallel_abort('ipgl allocation failure')
817  ipgl = 0
818 
819  j = 1
820  do i = 1, np_global
821  node => nodes_global(i)
822  if(node%domainID == myrank+1) then
823  iplg(j) = i
824  ipgl(i) = j
825  j = j + 1
826  endif
827  end do
828 
829  ! calc # sides local again, because the nodes now belongs to another domain
830  ns = 0
831  do i = 1, np
832  node => nodes(i)
833  ns = ns + node%nConnNodes
834  end do

References yownodepool::ipgl, yownodepool::iplg, yowdatapool::myrank, yownodepool::nodes(), yownodepool::nodes_global, yownodepool::np, yownodepool::np_global, yownodepool::np_perproc, yownodepool::np_perprocsum, yowsidepool::ns, yowdatapool::ntasks, and yowerr::parallel_abort().

Referenced by initfromgriddim().

◆ postpartition2()

subroutine yowpdlibmain::postpartition2 ( integer, dimension(3,ne_global), intent(in)  INE_global)

this collects all data which depends on ghost information alter: ne, INE, x, y, z, ielg

Todo:
create some localnode to localghost mapping What number is this ghost

Definition at line 1181 of file yowpdlibmain.F90.

1181  use yowelementpool, only: ne, ne_global, ine, belongto, ielg
1182  use yowerr, only: parallel_abort
1183  use yowdatapool, only: myrank
1184  use yownodepool, only: np_global, np, nodes_global, iplg, t_node, ghostlg, ng, npa
1185  use yownodepool, only: x, y, z
1186  use w3gdatmd, only: xgrd, ygrd, zb
1187 
1188  integer, intent(in) :: INE_global(3,ne_global)
1189 
1190  integer :: i, j, k, stat, IP_glob
1191  type(t_Node), pointer :: node
1192  logical :: assigned
1193 
1194  ! to find the local elements, iterate over all global elements and check their node domain IDs
1195 
1196  ! step 1: calc the number of local elements
1197  ne = 0
1198  do i=1, ne_global
1199  if (belongto(ine_global(:,i))) then
1200  ne = ne +1
1201  endif
1202  end do
1203 
1204  ! step 2: fill the local element index array
1205  if(allocated(ine)) deallocate(ine)
1206  allocate(ine(3, ne), stat=stat)
1207  if(stat/=0) call parallel_abort('INE allocation failure')
1208  ine = 0
1209 
1210  ne = 0
1211  do i=1, ne_global
1212  ! yes, this element belongs to this domain
1213  if (belongto(ine_global(:,i))) then
1214  ne = ne + 1
1215  do j=1, 3
1216  assigned = .false.
1217  node => nodes_global(ine_global(j,i))
1218  if(node%domainID == myrank+1) then
1219  ine(j, ne) = node%id
1220  assigned = .true.
1221  else
1222  ! the element have a ghost node
1225  do k=1, ng
1226  if(node%id_global == ghostlg(k)) then
1227  ! conversion: the ghost nodes are stored behind the local nodes.
1228  if(ine(j,ne) /= 0) then
1229  write(*,*) "will write to INE(j, ne) but there is allready a value", j, ne, ine(j, ne)
1230  endif
1231  ine(j, ne) = np + k
1232  node%id = np+k
1233  assigned = .true.
1234  ! write(*,*) myrank, "node to ele", node%id_global-1, i-1, np+k
1235  exit
1236  endif
1237  end do
1238  endif
1239  if(assigned .eqv. .false.) then
1240  write(*,*) "Can't assign global node to INE", node%id_global
1241  endif
1242  end do
1243  endif
1244  end do
1245 
1246  ! check if INE contains 0.
1247  do i=1, ne
1248  if(minval(abs(ine(:,i))) == 0) then
1249  write(*,*) "0 in INE ne=", ne
1250  stop "0 in INE"
1251  endif
1252  end do
1253 
1254  if(maxval(ine) /= npa) then
1255  write(*,*) "MAXVAL(INE) /= npa ERROR?"
1256  endif
1257 
1258  ! create element local to global mapping ielg
1259  if(allocated(ielg)) deallocate(ielg)
1260  allocate(ielg(ne), stat=stat)
1261  if(stat/=0) call parallel_abort('ielg allocation failure')
1262  ielg = 0
1263 
1264  j = 0
1265  do i=1, ne_global
1266  if (belongto(ine_global(:,i))) then
1267  j = j +1
1268  ielg(j) = i
1269  end if
1270  end do
1271 
1272  ! fill the local x,y arrays
1273  if(allocated(x)) deallocate(x)
1274  allocate(x(npa), stat=stat)
1275  if(stat/=0) call parallel_abort('x allocation failure')
1276 
1277  if(allocated(y)) deallocate(y)
1278  allocate(y(npa), stat=stat)
1279  if(stat/=0) call parallel_abort('y allocation failure')
1280 
1281  if(allocated(z)) deallocate(z)
1282  allocate(z(npa), stat=stat)
1283  if(stat/=0) call parallel_abort('z allocation failure')
1284 
1285  do i=1, np
1286  ip_glob = iplg(i)
1287  x(i) = xgrd(1,ip_glob)
1288  y(i) = ygrd(1,ip_glob)
1289  z(i) = zb(ip_glob)
1290  end do
1291 
1292  do i=1, ng
1293  ip_glob = ghostlg(i)
1294  x(np+i) = xgrd(1,ip_glob)
1295  y(np+i) = ygrd(1,ip_glob)
1296  z(np+i) = zb(ip_glob)
1297  end do
1298 

References yowelementpool::belongto(), yownodepool::ghostlg, yowelementpool::ielg, yowelementpool::ine, yownodepool::iplg, yowdatapool::myrank, yowelementpool::ne, yowelementpool::ne_global, yownodepool::ng, yownodepool::nodes_global, yownodepool::np, yownodepool::np_global, yownodepool::npa, yowerr::parallel_abort(), yownodepool::x, w3gdatmd::xgrd, yownodepool::y, w3gdatmd::ygrd, yownodepool::z, and w3gdatmd::zb.

Referenced by initfromgriddim().

◆ prepartition()

subroutine yowpdlibmain::prepartition

pre-partition the mesh just divide the mesh into nTasks parts and create a premature iplg alter: np_perProc, np_perProcSum, np, iplg

Definition at line 293 of file yowpdlibmain.F90.

293  use yowdatapool, only: ntasks ,myrank
294  use yowerr, only: parallel_abort
295  use yownodepool, only: np_global, np, np_perproc, np_perprocsum, iplg
296 
297  integer :: i, stat
298 
299  ! determine equal number of nodes in each processor (except for the last one).
300  ! and create a provisional node local to global mapping iplg
301 
302  ! start the arrays from 0, because the first thread id is 0
303  if(allocated(np_perproc)) deallocate(np_perproc)
304  allocate(np_perproc(0:ntasks-1), stat=stat)
305  if(stat/=0) call parallel_abort('np_perProc allocation failure')
306 
307  if(allocated(np_perprocsum)) deallocate(np_perprocsum)
308  allocate(np_perprocsum(0:ntasks), stat=stat)
309  if(stat/=0) call parallel_abort('np_perProcSum allocation failure')
310 
311  np_perprocsum = 0
312  np = np_global / ntasks
313 
314  do i = 0, ntasks-2
315  np_perproc(i) = np
316  np_perprocsum(i+1) = np_perprocsum(i) + np
317  end do
318 
319  np_perproc(ntasks-1) = np_global - np_perprocsum(ntasks-1)
320  np_perprocsum(ntasks) = np_perprocsum(ntasks-1) + np_perproc(ntasks-1)
321  np = np_perproc(myrank)
322 
323  ! create a provisional node local to global mapping iplg
324  ! this will override later with the data from parmetis
325  if(allocated(iplg)) deallocate(iplg)
326  allocate(iplg(np), stat=stat);
327  if(stat/=0) call parallel_abort(' iplg allocate failure')
328  do i = 1, np
329  iplg(i) = i + np_perprocsum(myrank)
330  end do

References yownodepool::iplg, yowdatapool::myrank, yownodepool::np, yownodepool::np_global, yownodepool::np_perproc, yownodepool::np_perprocsum, yowdatapool::ntasks, and yowerr::parallel_abort().

Referenced by initfromgriddim().

◆ real_mpi_barrier_pdlib()

subroutine yowpdlibmain::real_mpi_barrier_pdlib ( integer, intent(in)  TheComm,
character(*), intent(in)  string 
)

Definition at line 178 of file yowpdlibmain.F90.

178 
179  include "mpif.h"
180  integer, intent(in) :: TheComm
181  character(*), intent(in) :: string
182  integer NbProc, eRank
183  integer :: istatus(MPI_STATUS_SIZE)
184  integer ierr, iField(1), iProc
185  ! Print *, 'Start of REAL_MPI_BARRIER_PDLIB'
186  CALL mpi_comm_rank(thecomm, erank, ierr)
187  CALL mpi_comm_size(thecomm, nbproc, ierr)
188  ! Print *, ' eRank=', eRank, ' NbProc=', NbProc
189  ifield(1)=1
190  IF (erank .eq. 0) THEN
191  DO iproc=2,nbproc
192  ! Print *, ' Before MPI_RECV 1 iProc=', iProc
193  CALL mpi_recv(ifield, 1, mpi_integer, iproc-1, 711, thecomm, istatus, ierr)
194  ! Print *, ' Before MPI_SEND 1'
195  CALL mpi_send(ifield, 1, mpi_integer, iproc-1, 712, thecomm, ierr)
196  END DO
197  ELSE
198  ! Print *, ' Before MPI_SEND 2 eRank=', eRank
199  CALL mpi_send(ifield, 1, mpi_integer, 0, 711, thecomm, ierr)
200  ! Print *, ' Before MPI_RECV 2 eRank=', eRank
201  CALL mpi_recv(ifield, 1, mpi_integer, 0, 712, thecomm, istatus, ierr)
202  END IF
203  ! Print *, 'Passing barrier string=', string

References yowerr::abort(), yowdatapool::comm, include(), yowdatapool::myrank, yowdatapool::ntasks, and yowerr::parallel_abort().

Referenced by runparmetis().

◆ runparmetis()

subroutine yowpdlibmain::runparmetis ( integer, intent(in)  MNP)

Collect all data for parmetis und partition the mesh after that, we know for every node the domain ID alter: t_Node::domainID.

Todo:

Definition at line 422 of file yowpdlibmain.F90.

422  use yowerr, only: parallel_abort
423  use yowdatapool, only: debugparmetis,debugpartition, ntasks, myrank, itype, comm
424  use yownodepool, only: np, npa, np_global, nodes, nodes_global, t_node, iplg, np_perprocsum, np_perproc
425  use yowsidepool, only: ns
426  use yowelementpool, only: ne, ne_global
427  use w3gdatmd, only: xgrd, ygrd
428  use mpi
429 
430  integer, intent(in) :: MNP
431 
432  ! Parmetis
433  ! Node neighbor information
434  integer :: wgtflag, numflag, ndims, nparts, edgecut, ncon
435  integer, allocatable :: xadj(:), part(:), vwgt(:), adjwgt(:), vtxdist(:), options(:), adjncy(:), iweights(:)
436  ! parmetis need single precision
437  real(4), allocatable :: xyz(:), tpwgts(:), ubvec(:)
438  integer :: IP_glob, itmp
439  integer :: ref
440  logical :: lexist = .false.
441 
442  ! Node to domain mapping.
443  ! np_global long. give the domain number for die global node number
444  integer, allocatable :: node2domain(:)
445 
446  ! Mics
447  integer :: i, j, stat, ierr
448  type(t_Node), pointer :: node, nodeNeighbor
449 
450  ! CALL REAL_MPI_BARRIER_PDLIB(comm, "runParmetis, step 1")
451  ! Create xadj and adjncy arrays. They holds the nodes neighbors in CSR Format
452  ! Here, the adjacency structure of a graph is represented by two arrays,
453  ! xadj[n+1] and adjncy[m]; n vertices and m edges. Every edge is listen twice
454  allocate(adjncy(ns), stat=stat)
455  if(stat/=0) call parallel_abort('adjncy allocation failure')
456  allocate(xadj(np+1), stat=stat)
457  if(stat/=0) call parallel_abort('xadj allocation failure')
458  ! CALL REAL_MPI_BARRIER_PDLIB(comm, "runParmetis, step 2")
459 
460  xadj = 0
461  xadj(1) = 1
462  adjncy = 0
463  do i=1,np
464  node => nodes(i)
465  xadj(i+1) = xadj(i) + node%nConnNodes
466  if(debugparmetis) write(710+myrank,*) i, xadj(i), xadj(i+1)
467  if(node%nConnNodes == 0) then
468  write(*,*) "Thread", myrank,"global node has no conn nodes", node%id_global
469  endif
470  do j=1, node%nConnNodes
471  nodeneighbor => node%connNodes(j)
472  adjncy(j + xadj(i) - 1) = nodeneighbor%id_global
473  if(debugparmetis) write(710+myrank,*) i, j, j + xadj(i) - 1, adjncy(j + xadj(i) - 1), nodeneighbor%id_global
474  end do
475  end do
476  ! CALL REAL_MPI_BARRIER_PDLIB(comm, "runParmetis, step 3")
477 
478  ! Option for Parmetis
479  allocate(options(3))
480  options(1)=1 ! 0: default options; 1: user options
481  if(debugparmetis) then
482  options(2)=15
483  else
484  options(2)=0 ! Level of information returned: see defs.h in ParMETIS-Lib dir
485  endif
486  options(3)=15 ! Random number seed
487  CALL real_mpi_barrier_pdlib(comm, "runParmetis, step 4")
488  ! Fortran-style numbering that starts from 1
489  numflag = 1
490  ! # dimensions of the space in which the graph is embedded
491  ndims = 2
492  ! # parts the mesh is divide. Usually nTasks
493  nparts = ntasks
494  ! # weight per node
495  ncon = 1
496  ! wgtflag: 0: none (vwgt and adjwgt are NULL); 1: edges (vwgt is NULL); 2: vertices (adjwgt is NULL); 3: both vertices & edges;
497 
498 #ifdef WEIGHTS
499  wgtflag = 2
500  INQUIRE ( file='weights.ww3', exist = lexist )
501  IF (lexist) THEN
502  OPEN(100001,file='weights.ww3',form='FORMATTED',status='unknown')
503  allocate(iweights(np_global)); iweights = 0
504  do i = 1, np_global
505  read(100001,*) iweights(i)
506  enddo
507  CLOSE(100001)
508  ELSE
509  wgtflag = 0
510  ENDIF
511 #else
512  wgtflag = 0
513 #endif
514 
515  ! Create weights
516  allocate(vwgt(np*ncon), stat=stat)
517  if(stat/=0) call parallel_abort('vwgt allocation failure')
518 
519 #ifdef WEIGHTS
520  if (lexist) then
521  do i = 1, np
522  itmp = max(1,int(real((iweights(iplg(i))+100))))
523  !vwgt(i) = max(1,int(real(iweights(iplg(i)))/maxval(iweights)) * 10)
524  !vwgt(i) = max(1,iweights(iplg(i)))
525  !vwgt(i) = max(1,itmp)
526  vwgt(i) = itmp
527  enddo
528  vwgt = 1
529  deallocate(iweights)
530  else
531  vwgt = 1
532  endif
533 #else
534  vwgt = 1
535 #endif
536 
537  allocate(adjwgt(ns), stat=stat)
538  if(stat/=0) call parallel_abort('adjwgt allocation failure')
540  adjwgt = 1
541 
542  ! Vertex weight fraction
543  allocate(tpwgts(ncon*ntasks),stat=stat)
544  if(stat/=0) call parallel_abort('partition: tpwgts allocation failure')
545  tpwgts=1.0/real(ntasks)
546 
547  ! Imbalance tolerance
548  allocate(ubvec(ncon),stat=stat)
549  if(stat/=0) call parallel_abort('partition: ubvec allocation failure')
550  ubvec=1.01
551  CALL real_mpi_barrier_pdlib(comm, "runParmetis, step 5")
552 
553  ! Partition dual graph
554  allocate(xyz(2*np),stat=stat)
555  if(stat/=0) call parallel_abort('xyz: ubvec allocation failure')
556  if(debugparmetis) write(710+myrank,*) 'np_global, ne_global, np, npa, ne'
557  if(debugparmetis) write(710+myrank,*) np_global, ne_global, np, npa, ne
558  do i = 1, np
559  ip_glob = iplg(i)
560  !AR: this is questionable ...
561  xyz(2*(i-1)+1) = real(xgrd(1,ip_glob))
562  xyz(2*(i-1)+2) = real(ygrd(1,ip_glob))
563  if(debugparmetis) then
564  write(710+myrank,*) i, np, xyz(2*(i-1)+1), xyz(2*(i-1)+2)
565  call flush(710+myrank)
566  endif
567  end do
568  CALL real_mpi_barrier_pdlib(comm, "runParmetis, step 6")
569 
570  ! Partition array returned from ParMeTiS
571  allocate(part(np),stat=stat)
572  if(stat/=0) call parallel_abort('part: ubvec allocation failure')
573 
574  ! ParMeTiS vertex distribution array (starts at 1)
575  allocate(vtxdist(ntasks+1),stat=stat)
576  if(stat/=0) call parallel_abort('partition: vtxdist allocation failure')
577 
578  call mpi_allgather(np_perprocsum(myrank)+1, 1, itype, vtxdist, 1, itype, comm, ierr)
579  if(ierr/=mpi_success) call parallel_abort('partition: mpi_allgather',ierr)
580  vtxdist(ntasks+1)=np_global+1
581  ! CALL REAL_MPI_BARRIER_PDLIB(comm, "runParmetis, step 7")
582 
583  ! check vtxdist
584  ! myrank stats from 0
585  if((vtxdist(myrank+2) - vtxdist(myrank+1)) < 1) then
586  write(*,*) "Thread", myrank, "has no nodes"
587  write(*,*) "vtxdist", vtxdist
588  CALL abort("Poor initial vertex distribution detected")
589  endif
590 
591 
592 
593  ! My notes from manual:
594  ! p: # of processors;
595  ! n: total # of vertices (local) in graph sense;
596  ! m: total # of neighboring vertices ("edges"); double counted between neighboring vertice u and v.
597  ! ncon: # of weights for each vertex;
598  ! int(in) vtxdist(p+1): Processor j stores vertices vtxdist(j):vtxdist(j+1)-1
599  ! int (in) xadj(n+1), adjncy(m):
600  ! locally, vertex j's neighboring vertices are adjncy(xadj(j):xadj(j+1)-1). adjncy points to global index;
601  ! int(in) vwgt(ncon*n), adjwgt(m): weights at vertices and "edges". Format of adjwgt follows adjncy;
602  ! int(in) wgtflag: 0: none (vwgt and adjwgt are NULL);
603  ! 1: edges (vwgt is NULL); 2: vertices (adjwgt is NULL); 3: both vertices & edges;
604  ! int(in) numflag: 0: C-style numbering from 0; 1: FORTRAN style from 1;
605  ! int(in) ndims: 2 or 3 (D);
606  ! float(in) xyz(ndims*n): coordinate for vertex j is xyz(j*ndims:(j+1)*ndims-1);
607  ! int(in) nparts: # of desired sub-domains (usually nTasks);
608  ! float(in) tpwgts(ncon*nparts): =1/nparts if sub-domains are to be of same size for each vertex weight;
609  ! float(in) ubvec(ncon): imbalance tolerance for each weight;
610  ! int(in) options: additonal parameters for the routine (see above);
611  ! int(out) edgecut: # of edges that are cut by the partitioning;
612  ! int(out) part(): array size = # of local vertices. It stores indices of local vertices.
613 
614  ! write(1112+myrank,*) "Thread",myrank,"sum;vtxdist", sum(vtxdist), vtxdist
615  ! write(1112+myrank,*) "Thread",myrank,"sum;xadj", sum(xadj), xadj
616  ! write(1112+myrank,*) "Thread",myrank,"sum;adjncy", sum(adjncy), adjncy
617 
618  CALL real_mpi_barrier_pdlib(comm, "runParmetis, step 8")
619 
620  if(debugparmetis) then
621  write(710+myrank,*) vtxdist, xadj, adjncy, &
622  vwgt, & !vwgt - ignore weights
623  adjwgt, & ! adjwgt - ignore weights
624  wgtflag, &
625  numflag,ndims,ncon,nparts,tpwgts,ubvec,options, &
626  edgecut,part,comm
627  call flush(710+myrank)
628  endif
629 
630  !if(debugParmetis) write(710+myrank,*) "Run ParMETIS now..."
631 #ifdef W3_SCOTCH
632  call scotch_parmetis_v3_partgeomkway(vtxdist, xadj, adjncy, &
633  vwgt, & !vwgt - ignore weights
634  adjwgt, & ! adjwgt - ignore weights
635  wgtflag, &
636  numflag,ndims,xyz,ncon,nparts,tpwgts,ubvec,options, &
637  edgecut,part, comm,ref)
638 #endif
639 
640 #ifdef W3_METIS
641  call parmetis_v3_partgeomkway(vtxdist, xadj, adjncy, &
642  vwgt, & !vwgt - ignore weights
643  adjwgt, & ! adjwgt - ignore weights
644  wgtflag, &
645  numflag,ndims,xyz,ncon,nparts,tpwgts,ubvec,options, &
646  edgecut,part, comm)
647 #endif
648 
649 
650  CALL real_mpi_barrier_pdlib(comm, "runParmetis, step 9")
651 
652  if(ntasks == 1) then
653  ! write(*,*) myrank, "minval part", minval(part)
654  if(minval(part) == 0) then
655  part(:) = part(:) + 1
656  endif
657  endif
658 
659  ! write(*,*) myrank, "edge cuted", edgecut
660 
661  ! Collect the parmetis data from all threads
662  ! and create a global node to domain number mapping
663  allocate(node2domain(np_global),stat=stat)
664  if(stat/=0) call parallel_abort(' node2domain allocation failure')
665  !
666  call mpi_allgatherv(part, np, itype, node2domain, np_perproc, np_perprocsum, itype, comm, ierr)
667  if(ierr/=mpi_success) call parallel_abort('mpi_allgatherv ',ierr)
668  !
669  do i = 1, np_global
670  node => nodes_global(i)
671  node%domainID = node2domain(node%id_global)
672  end do
673  ! CALL REAL_MPI_BARRIER_PDLIB(comm, "runParmetis, step 10")
674 
675 
676  ! write out partition info for katerfempresenter
677  if(debugpartition) write(600,*) node2domain
678  ! Print *, 'runparmetis step 1'
679 
680  if(allocated(xadj)) deallocate(xadj)
681  ! Print *, 'runparmetis step 2'
682  if(allocated(adjncy)) deallocate(adjncy)
683  ! Print *, 'runparmetis step 3'
684  if(allocated(part)) deallocate(part)
685  ! Print *, 'runparmetis step 4'
686  if(allocated(vwgt)) deallocate(vwgt)
687  ! Print *, 'runparmetis step 5'
688  if(allocated(adjwgt)) deallocate(adjwgt)
689  ! Print *, 'runparmetis step 6'
690  if(allocated(xyz)) deallocate(xyz)
691  ! Print *, 'runparmetis step 7'
692  if(allocated(tpwgts)) deallocate(tpwgts)
693  ! Print *, 'runparmetis step 8'
694  if(allocated(ubvec)) deallocate(ubvec)
695  ! Print *, 'runparmetis step 9'
696  if(allocated(vtxdist)) deallocate(vtxdist)
697  ! Print *, 'runparmetis step 10'
698  if(allocated(node2domain)) deallocate(node2domain)
699  ! Print *, 'runparmetis step 11'

References yowdatapool::comm, yowdatapool::debugparmetis, yowdatapool::debugpartition, file(), yownodepool::iplg, yowdatapool::itype, yowdatapool::myrank, yowelementpool::ne, yowelementpool::ne_global, yownodepool::nodes(), yownodepool::nodes_global, yownodepool::np, yownodepool::np_global, yownodepool::np_perproc, yownodepool::np_perprocsum, yownodepool::npa, yowsidepool::ns, yowdatapool::ntasks, yowerr::parallel_abort(), real_mpi_barrier_pdlib(), w3gdatmd::xgrd, and w3gdatmd::ygrd.

Referenced by initfromgriddim().

yowerr::parallel_abort
subroutine parallel_abort(string, error)
Definition: yowerr.F90:43
include
cmake src_list cmake include(${CMAKE_CURRENT_SOURCE_DIR}/cmake/check_switches.cmake) check_switches("$
Definition: CMakeLists.txt:15
yowexchangemodule::pdlib_exchange1dreal
subroutine, public pdlib_exchange1dreal(U)
exchange values in U.
Definition: yowexchangeModule.F90:251
yowrankmodule::initrankmodule
subroutine, public initrankmodule()
allocate and exchange
Definition: yowrankModule.F90:76
w3gdatmd::ygrd
double precision, dimension(:,:), pointer ygrd
Definition: w3gdatmd.F90:1205
yowelementpool
Definition: yowelementpool.F90:38
yowexchangemodule::initnbrdomains
subroutine, public initnbrdomains(nConnD)
Definition: yowexchangeModule.F90:224
w3gdatmd::zb
real, dimension(:), pointer zb
Definition: w3gdatmd.F90:1195
w3gdatmd::xgrd
double precision, dimension(:,:), pointer xgrd
Definition: w3gdatmd.F90:1205
yowelementpool::finalizeelementpool
subroutine, public finalizeelementpool()
Definition: yowelementpool.F90:102
yowerr
Has some subroutine to make a nice error message.
Definition: yowerr.F90:39
scrip_constants::tiny
real(kind=scrip_r8), parameter tiny
Definition: scrip_constants.f:50
yowsidepool::ns
integer ns
Number of sides, local.
Definition: yowsidepool.F90:46
yowrankmodule::ipgl_npa
integer, dimension(:), allocatable, public ipgl_npa
Definition: yowrankModule.F90:70
yownodepool
Has data that belong to nodes.
Definition: yownodepool.F90:39
yowrankmodule
Provides access to some information of all threads e.g.
Definition: yowrankModule.F90:44
yownodepool::t_node
Holds the nodes data.
Definition: yownodepool.F90:48
yowexchangemodule::creatempitypes
subroutine, public creatempitypes()
Definition: yowexchangeModule.F90:236
yowdatapool
Has fancy data.
Definition: yowdatapool.F90:39
yowsidepool
Definition: yowsidepool.F90:38
yowexchangemodule
Has only the ghost nodes assign to a neighbor domain.
Definition: yowexchangeModule.F90:39
yowerr::abort
subroutine abort(string, line, file, errno)
print various error strings and exit.
Definition: yowerr.F90:93
yowrankmodule::finalizerankmodule
subroutine, public finalizerankmodule()
Definition: yowrankModule.F90:249
w3gdatmd
Definition: w3gdatmd.F90:16
yowexchangemodule::setdimsize
subroutine, public setdimsize(second)
set the size of the second and third dimension for exchange
Definition: yowexchangeModule.F90:375
yownodepool::finalizenodepool
subroutine, public finalizenodepool()
Definition: yownodepool.F90:220
yowexchangemodule::finalizeexchangemodule
subroutine, public finalizeexchangemodule()
Definition: yowexchangeModule.F90:381
yowsidepool::ns_global
integer ns_global
Number of sides, global Number of nodes neighbors?.
Definition: yowsidepool.F90:43