e4d_setup.F90       coverage:  90.00 %func     77.97 %block


     1) module e4d_setup
     2) 
     3)   use vars
     4) 
     5)   implicit none
     6)   integer :: ios
     7)   real :: xorig,yorig,zorig
     8) 
     9) contains
    10)   !________________________________________________________________
    11)   subroutine setup_e4d
    12)     implicit none
    13)     integer :: sz
    14)     
    15)     log_file = 'e4d_'//trim(adjustl(pflotran_group_prefix))//'.log'   
    16)     open(13,file=trim(log_file),action='write',status='replace')
    17)     close(13)
    18)     
    19)     if (my_rank>0) then
    20)        call slave_setup
    21)        return
    22)     end if
    23)  
    24) !#if 0
    25)     call read_input 
    26) 
    27)     if (ios .eq. -1) then
    28)        call send_command(0)
    29)        return
    30)     end if
    31)     open(13,file=trim(log_file),action='write',status='old',position='append')
    32)     write(13,*) "Read Input Complete___________________________________"
    33)     close(13)
    34) 
    35)     call send_info
    36) 
    37)     call setup_forward
    38)     open(13,file=trim(log_file),action='write',status='old',position='append')
    39)     write(13,*) "Setup Forward Complete___________________________________"
    40)     close(13)
    41)     if (ios .eq. -1) then
    42)        call send_command(0)
    43)        return
    44)     end if
    45)   
    46) !#endif
    47) !    call send_info_geh
    48) 
    49)     call send_command(0)
    50)     open(13,file=trim(log_file),action='write',status='old',position='append')
    51)     write(13,*) "Setup Phase Complete___________________________________"
    52)     close(13)
    53)   end subroutine setup_e4d
    54)   !________________________________________________________________
    55) 
    56)   !________________________________________________________________
    57)   subroutine read_input
    58)     implicit none
    59)     
    60)     integer :: nchar,junk,i,check,j,nchr,npre,flg,flnum
    61)     integer :: a,b,m,n
    62)     real :: ex,ey,ez,eflg,vi,wd,tsig
    63)     real, dimension(4) :: etmp 
    64)     logical :: exst
    65) 
    66)     !get the name of the mesh file
    67)     call elog(0,0,flg)
    68)     if (flg .ne. 0) then
    69)        ios = -1
    70)        return
    71)     end if
    72) 
    73)     open(10,file='e4d.inp',status='old',action='read') 
    74)     read(10,*,IOSTAT=ios) mshfile
    75)     call elog(1,ios,flg); 
    76)     if (flg .eq. -1) then
    77)        ios=-1
    78)        return
    79)     end if
    80) 
    81)     read(10,*,IOSTAT=ios) efile
    82)     call elog(2,ios,flg); 
    83)     if (flg .eq. -1) then
    84)        ios=-1
    85)        return
    86)     end if
    87) 
    88)     read(10,*,IOSTAT=ios) sigfile
    89)     call elog(3,ios,flg); 
    90)     if (flg .eq. -1) then
    91)        ios=-1
    92)        return
    93)     end if
    94) 
    95)     !read(10,*,IOSTAT=ios) mapfile
    96)     !call elog(4,ios,flg); 
    97)     !if (flg .eq. -1) then
    98)     !   ios=-1
    99)     !   return
   100)     !end if
   101)    
   102) 
   103)     read(10,*,IOSTAT=ios) list_file
   104)     call elog(23,ios,flg)
   105)     if (flg .eq. -1) then
   106)        ios=-1
   107)        return
   108)     end if
   109)     close(10)
   110) 
   111)     !get the number of character is in the prefix
   112)     nchar = 0
   113)     do i=1,40
   114)        if (mshfile(i:i) == '.') then
   115)           nchar = i
   116)           exit
   117)        end if
   118)     end do
   119)     
   120)     !check to see if the files exist
   121) 
   122)     !!Allocate/read the electrode positions and survey conf
   123)      open(10,file=efile,status='old',action='read')   
   124)      read(10,*,IOSTAT=ios) ne
   125)      call elog(5,ios,flg)
   126)      if (ios .ne. 0) then
   127)         ios=-1
   128)         return
   129)      end if
   130)        
   131)        
   132)      allocate(e_pos(ne,4))
   133)      do i=1,ne
   134)         read(10,*,IOSTAT=ios) junk,etmp
   135)         if (ios .ne. 0) then
   136)            call elog(6,i,flg)
   137)            ios=-1
   138)            return
   139)         end if
   140)         e_pos(junk,1:4)=etmp
   141)      end do
   142)      
   143)      !get the meshfile prefix
   144)      nchr=len_trim(mshfile)
   145)      do i=1,nchr
   146)         if (mshfile(i:i)=='.') then
   147)            npre=i-1;
   148)            exit
   149)         end if
   150)      end do
   151) 
   152)      !!translate electrodes
   153)      call elog(7,npre,flg)
   154)      if (flg.eq.-1) then
   155)         ios=-1
   156)         return
   157)      end if
   158) 
   159)      open(21,file=mshfile(1:npre)//".trn",status="old",action="read")
   160)      read(21,*,IOSTAT=ios) xorig,yorig,zorig
   161)      if (ios .ne. 0) then
   162)         call elog(8,ios,flg)
   163)         ios=-1
   164)         return
   165)      end if
   166)      close(21)
   167)      e_pos(:,1) = e_pos(:,1)-xorig
   168)      e_pos(:,2) = e_pos(:,2)-yorig
   169)      e_pos(:,3) = e_pos(:,3)-zorig
   170) 
   171)      !!Read in the survey
   172)      read(10,*,IOSTAT=ios) nm
   173)      call elog(9,ios,flg)
   174)      if (ios .ne. 0) then
   175)         ios = -1
   176)         return
   177)      end if
   178) 
   179)      allocate(s_conf(nm,4))
   180)      do i=1,nm
   181)         read(10,*,IOSTAT=ios) junk,s_conf(i,1:4)
   182)         if (ios.ne.0) then
   183)            call elog(10,i,flg)
   184)            ios=-1
   185)            return
   186)         end if
   187)      end do
   188)      close(10) 
   189)        
   190)        !!Read the conductivity file
   191)        open(10,file=trim(sigfile),status='old',action='read')
   192)        read(10,*,IOSTAT=ios) j !FF,sw_sig,gw_sig
   193)        call elog(11,ios,j)
   194)        if (ios .ne. 0) then
   195)           ios=-1
   196)           close(10)
   197)           return
   198)        end if
   199) 
   200)        allocate(base_sigma(j))
   201)        do i=1,j
   202)           read(10,*,IOSTAT=ios) base_sigma(i)
   203)           if (ios .ne. 0) then
   204)              call elog(12,i,flg)
   205)              ios=-1
   206)              close(10)
   207)              return
   208)           end if
   209)        end do
   210)        close(10)
   211) 
   212)        !make sure pf_mesh.txt exists
   213)        call elog(37,ios,flg)
   214)        if(ios.ne.0) then
   215)           ios=-1
   216)           close(1)
   217)           return
   218)        end if
   219) 
   220)        !!read the list file and all files listed and check for errors
   221)        open(10,file=trim(list_file),status='old',action='read')
   222)        read(10,*,IOSTAT=ios) ntime,FF,sw_sig,gw_sig
   223)        call elog(24,ios,flg)
   224)        if (flg==-1) then
   225)           ios=-1
   226)           close(10)
   227)           return
   228)        end if
   229)        flnum=0
   230)        do j=1,ntime
   231)           flnum=flnum+1
   232)           read(10,*,IOSTAT=ios) e4d_time,csrv_file,ccond_file
   233)           call elog(25,ios,flnum)
   234)           if (ios .ne. 0) then
   235)              ios=-1
   236)              close(10)
   237)              return
   238)           end if
   239)           call elog(26,ios,flnum)
   240)           if (flnum .ne. 0) then
   241)              ios=-1
   242)              close(10)
   243)              return
   244)           end if
   245)        end do
   246)           
   247)        
   248)   end subroutine read_input
   249)   !________________________________________________________________
   250) 
   251)   !________________________________________________________________
   252)   subroutine send_info
   253)     implicit none
   254) 
   255)     call send_command(7)
   256)     call MPI_BCAST(nm, 1, MPI_INTEGER, 0,E4D_COMM,ierr)
   257)     call MPI_BCAST(s_conf, 4*nm,MPI_INTEGER,0,E4D_COMM,ierr)
   258) 
   259)   end subroutine send_info
   260)   !________________________________________________________________
   261) 
   262)   !____________________________________________________________________
   263)   subroutine setup_forward
   264)     implicit none
   265)     
   266)     integer :: npre,i,nchr,dim,bflag,ns,itmp,ios,flg
   267)     real :: jnk,jnk1 
   268)     integer :: status(MPI_STATUS_SIZE)
   269)     
   270)     !command slaves to do the setup
   271)     call send_command(2)
   272)    
   273)     !read and send the nodes
   274)     !!Determine the prefix of the tetgen output files
   275)     nchr=len_trim(mshfile)
   276)     do i=1,nchr
   277)        if (mshfile(i:i)=='.') then
   278)           npre=i+1;
   279)           exit
   280)        end if
   281)     end do
   282) 
   283)     !!OPEN AND READ THE NODE POSITIONS 
   284)     !!nodes(x,y,z) and nbounds
   285)     call elog(15,npre,flg)
   286)     if (flg==-1) then
   287)        ios=-1
   288)        return
   289)     end if
   290) 
   291)     open(10,file=mshfile(1:npre)//".node",status="old",action="read")
   292)     read(10,*,IOSTAT=ios) nnodes,dim,jnk,bflag
   293)     if (ios .ne. 0) then
   294)        call elog(16,ios,flg)
   295)        ios=-1
   296)        return
   297)     end if
   298)   
   299)     allocate(nodes(nnodes,3),nbounds(nnodes))
   300)    
   301)     do i=1,nnodes
   302)        read(10,*,IOSTAT=ios) jnk,nodes(i,1:3),jnk1,nbounds(i)
   303)        if (ios .ne. 0) then
   304)           call elog(17,i,flg)
   305)           ios = -1
   306)           return
   307)        end if
   308)        if (nbounds(i)>2) nbounds(i)=7
   309)     end do
   310)     close(10)
   311) 
   312)     !!Send the nodes to the slaves
   313)     call MPI_BCAST(nnodes, 1, MPI_INTEGER , 0, E4D_COMM, ierr )
   314)     call MPI_BCAST(nodes, nnodes*3, MPI_REAL , 0, E4D_COMM, ierr )
   315)     call MPI_BCAST(nbounds, nnodes, MPI_INTEGER , 0, E4D_COMM, ierr )
   316)    
   317)     !!Read in the elements and send each slave it's assigned elements
   318)     call elog(18,npre,flg)
   319)     if (flg==-1) then
   320)        ios=-1
   321)        return
   322)     end if
   323) 
   324)     open(10,file=mshfile(1:npre)//".ele",status="old",action="read")
   325)     read(10,*,IOSTAT=ios) nelem,dim,jnk
   326)     if (ios .ne. 0) then
   327)        call elog(19,ios,flg)
   328)        ios=-1
   329)        return
   330)     end if
   331) 
   332)     allocate(elements(nelem,4),zones(nelem))
   333)     do i=1,nelem
   334)        read(10,*,IOSTAT=ios) jnk,elements(i,1:4),zones(i)
   335)        if (ios .ne. 0) then
   336)           call elog(20,i,flg)
   337)           ios = -1
   338)           return
   339)        end if
   340)     end do
   341)     close(10)
   342)  
   343)     call MPI_BCAST(nelem, 1, MPI_INTEGER, 0,E4D_COMM,ierr)
   344)     call MPI_BCAST(elements, nelem*4,MPI_INTEGER,0,E4D_COMM,ierr)
   345)  
   346)     !send the assignments
   347)     call send_dists
   348) 
   349)     !build the mesh interpolation matrix
   350)     call get_mesh_interp
   351)  
   352)     !get electrode nodes
   353)     call get_electrode_nodes 
   354)  
   355)     call elog(21,ios,flg)
   356)     ios=0
   357)   end subroutine setup_forward
   358)   !____________________________________________________________________
   359)  
   360)   !____________________________________________________________________
   361)   subroutine get_mesh_interp
   362)     implicit none
   363)     integer :: i,ierr,st,en
   364)     integer :: status(MPI_STATUS_SIZE)
   365)     integer, dimension(n_rank-1) :: tnmap
   366)     
   367)     nmap=0
   368)     
   369)     do i=1,n_rank-1
   370)        call MPI_RECV(tnmap(i),1,MPI_INTEGER,i,i,E4D_COMM,status,ierr)
   371)       
   372)     end do
   373)    
   374)     nmap=sum(tnmap)
   375) 
   376)     open(13,file=trim(log_file),status='old',action='write',position='append')
   377)     write(13,*) "Total number of mapping weights: ",nmap
   378)     do i=1,n_rank-1
   379)        write(13,*) "Rank ",i,": ",tnmap(i)
   380)     end do
   381)     close(13)
   382)     
   383)   
   384)     allocate(map_inds(nmap,2),map(nmap))
   385)     
   386)     do i=1,n_rank-1
   387)        if(i==1) then
   388)           st=1
   389)           en=tnmap(1)
   390)           call MPI_RECV(map_inds(st:en,1),tnmap(i),MPI_INTEGER,i,i,E4D_COMM,status,ierr)
   391)           call MPI_RECV(map_inds(st:en,2),tnmap(i),MPI_INTEGER,i,i,E4D_COMM,status,ierr)
   392)           call MPI_RECV(map(st:en),tnmap(i),MPI_REAL,i,i,E4D_COMM,status,ierr)
   393)        else
   394)           st=sum(tnmap(1:i-1))+1
   395)           en=sum(tnmap(1:i))
   396)           call MPI_RECV(map_inds(st:en,1),tnmap(i),MPI_INTEGER,i,i,E4D_COMM,status,ierr)
   397)           call MPI_RECV(map_inds(st:en,2),tnmap(i),MPI_INTEGER,i,i,E4D_COMM,status,ierr)
   398)           call MPI_RECV(map(st:en),tnmap(i),MPI_REAL,i,i,E4D_COMM,status,ierr)
   399)         
   400)        end if
   401)     end do
   402)   
   403) #if 0 
   404) !geh: this do loop produces an unwanted fort.20 file.  where should this
   405) !     data be written?
   406)     do i=1,nmap
   407)        write(20,*) map_inds(i,:),map(i)
   408)     end do
   409) #endif
   410)   end subroutine get_mesh_interp
   411)   !____________________________________________________________________
   412)  !____________________________________________________________________
   413)   subroutine send_dists
   414)     !read inputs and distribute the run info to the slave
   415)     implicit none
   416)     integer :: neven,nextra,ce,i,j,k
   417) 
   418)     
   419)     !!divide the electrodes up among the slave processes
   420)     !!note ne is the number of electrodes defined in read_inp
   421)     tne = ne  
   422)     neven = ne/(n_rank-1)
   423)     nextra = ne - neven*(n_rank-1)  
   424)     if (allocated(eind)) deallocate(eind)
   425)     allocate(eind(n_rank-1,2))
   426)     
   427)     !!Build eind
   428)     if (neven > 0) then
   429)        ce = 1
   430)        do i=1,n_rank-1-nextra
   431)           eind(i,1) = ce
   432)           eind(i,2) = ce+neven-1
   433)           ce = ce+neven
   434)        end do
   435)        
   436)        if (nextra > 0) then
   437)           do i=n_rank-nextra,n_rank-1    
   438)              eind(i,1) = ce
   439)              eind(i,2) = ce+neven
   440)              ce=ce+neven+1
   441)           end do
   442)        end if
   443)     else
   444)        !!in this case there are more processors than electrodes
   445)        !!THIS SHOULD BE AVOIDED
   446)        do i=1,n_rank-1
   447)           eind(i,1) = i
   448)           eind(i,2) = i
   449)        end do
   450)     end if
   451) 
   452)     !!Build jaco assignments
   453)     if(allocated(jind)) deallocate(jind)
   454)     allocate(jind(n_rank-1,2))
   455)     jind = 0
   456)     if(nelem < (n_rank-1)) then
   457)        do i=1,nelem
   458)           jind(i,1:2)=i
   459)        end do
   460)     else
   461)        neven = nelem/(n_rank-1)
   462)        nextra = nelem - neven*(n_rank-1)
   463)        
   464)        jind=0
   465)        if(neven > 0) then
   466)           ce = 1
   467)           do i=1,n_rank-1-nextra
   468)              jind(i,1) = ce
   469)              jind(i,2) = ce+neven-1
   470)              ce = ce+neven
   471)           end do
   472)           
   473)           if(nextra > 0) then
   474)              do i=n_rank-nextra,n_rank-1    
   475)                 jind(i,1) = ce
   476)                 jind(i,2) = ce+neven
   477)                 ce=ce+neven+1
   478)              end do
   479)           end if
   480)        end if
   481)     end if
   482)     
   483)     !this will be replaced with mesh info from pflotran
   484)     call get_pfmesh
   485)  
   486)     !!send assignments
   487)     call MPI_BCAST(tne, 1, MPI_INTEGER , 0, E4D_COMM, ierr )
   488)     call MPI_BCAST(e_pos,4*tne ,MPI_REAL , 0, E4D_COMM, ierr )
   489)     call MPI_BCAST(eind,2*(n_rank-1),MPI_INTEGER , 0, E4D_COMM, ierr )
   490)     call MPI_BCAST(jind,2*(n_rank-1),MPI_INTEGER , 0, E4D_COMM, ierr )
   491)     call MPI_BCAST(pfnx,1,MPI_INTEGER,0,E4D_COMM, ierr)
   492)     call MPI_BCAST(pfny,1,MPI_INTEGER,0,E4D_COMM, ierr)
   493)     call MPI_BCAST(pfnz,1,MPI_INTEGER,0,E4D_COMM, ierr)
   494)     call MPI_BCAST(pfxcb,pfnx+1,MPI_REAL,0,E4D_COMM,ierr)
   495)     call MPI_BCAST(pfycb,pfny+1,MPI_REAL,0,E4D_COMM,ierr)
   496)     call MPI_BCAST(pfzcb,pfnz+1,MPI_REAL,0,E4D_COMM,ierr)
   497)     deallocate(pfxcb)
   498)     deallocate(pfycb)
   499)     deallocate(pfzcb)
   500)   
   501) 
   502)   end subroutine send_dists
   503)   !______________________________________________________________
   504) 
   505)  
   506)   !____________________________________________________________________
   507)   subroutine get_pfmesh
   508)     implicit none
   509)     integer :: i
   510)     logical :: chech
   511)     
   512) 
   513)     open(12,file='pf_mesh.txt',status='old',action='read')
   514)     read(12,*) pfnx
   515)     allocate(pfxcb(pfnx+1))
   516)     read(12,*) pfxcb(1:pfnx+1)
   517)   
   518)     read(12,*) pfny
   519)     allocate(pfycb(pfny+1))
   520)     read(12,*) pfycb(1:pfny+1)
   521) 
   522)     read(12,*) pfnz
   523)     allocate(pfzcb(pfnz+1))
   524)     read(12,*) pfzcb(1:pfnz+1)
   525) 
   526)     close(12)
   527)     open(13,file=trim(log_file),action='write',status='old',position='append')
   528)     write(13,*) "________________PF MESH PARAMETERS___________________"
   529)     write(13,*) "NX, NY, NZ: ",pfnx,pfny,pfnz
   530)     write(13,*) "xmin  xmax: ",pfxcb(1),pfxcb(pfnx+1)
   531)     write(13,*) "ymin  ymax: ",pfycb(1),pfycb(pfny+1)
   532)     write(13,*) "zmin  zmax: ",pfzcb(1),pfzcb(pfnz+1)
   533)     write(13,*)
   534)     close(13)
   535)     
   536)   end subroutine get_pfmesh
   537)   !____________________________________________________________________
   538) 
   539)   !______________________________________________________________
   540)   subroutine get_electrode_nodes
   541)     implicit none
   542)     real, dimension(nnodes) :: dist
   543)     integer :: i
   544)     integer, dimension(1) :: indb
   545)    
   546)     if (.not.allocated(e_nods)) then
   547)        allocate (e_nods(tne))
   548)     end if
   549)   
   550)     do i=1,tne
   551)        dist=sqrt( (e_pos(i,1) - nodes(:,1))**2 + &
   552)             (e_pos(i,2) - nodes(:,2))**2 + &
   553)             (e_pos(i,3) - nodes(:,3))**2)
   554)        
   555)        indb = minloc(dist)
   556)        
   557)        if (my_rank==1) then
   558)           if (dist(indb(1)) > 0.1) then
   559)              !!There should always be a node at all
   560)              !!electrode locations. If not, something is
   561)              !!wrong
   562)              write(*,1001) i,dist(indb)
   563)           end if
   564)        end if
   565)        e_nods(i) = indb(1)
   566) 
   567) 1001   format(" WARNING: ELECTRODE ",I5," HAS BEEN MOVED ",F10.4," TO THE NEAREST NODE")
   568)     end do
   569)     
   570)   end subroutine get_electrode_nodes
   571)   !____________________________________________________________________
   572) 
   573)   !______________________________________________________________
   574)   subroutine send_command(com)
   575)     !!Send a general command to the slaves
   576)     integer :: com
   577)     
   578)     call MPI_BCAST(com,1,MPI_INTEGER,0,E4D_COMM,ierr)
   579)     
   580)   end subroutine send_command
   581)   !________________________________________________________________
   582)   
   583)   !__________________SLAVE ROUTINES________________________________
   584)   !________________________________________________________________
   585)   subroutine slave_setup
   586)     
   587)     implicit none
   588)     integer :: command
   589)     
   590) 100 continue
   591)     !Recieve a command from master
   592)     call MPI_BCAST(command,1,MPI_INTEGER,0,E4D_COMM,ierr)
   593)     
   594)     
   595)     !return to main
   596)     if (command == 0) then
   597)        return
   598) 
   599)     else if (command == 2) then
   600)        call setup_frun
   601)        goto 100
   602) 
   603)     else if (command == 7) then
   604)        call receive_info
   605)        goto 100
   606)        
   607)     else if (command == 86) then
   608)        call receive_info_geh
   609)        goto 100
   610)        
   611)     else
   612)        goto 100
   613) 
   614)     end if
   615)   end subroutine slave_setup
   616)   !________________________________________________________________
   617) 
   618)   !__________________________________________________________________
   619)   subroutine receive_info
   620)     implicit none
   621)     integer :: i
   622)     
   623)     if (allocated(s_conf)) deallocate(s_conf)
   624)     call MPI_BCAST(nm, 1, MPI_INTEGER, 0,E4D_COMM,ierr)
   625)     allocate(s_conf(nm,4))
   626)     call MPI_BCAST(s_conf, 4*nm,MPI_INTEGER,0,E4D_COMM,ierr)
   627)     
   628)   end subroutine receive_info
   629)   !__________________________________________________________________
   630)   
   631)   !__________________________________________________________________
   632)   subroutine setup_frun
   633)     implicit none
   634) 
   635) #include "petsc/finclude/petscsys.h"
   636)     
   637)     integer :: status(MPI_STATUS_SIZE)
   638)     integer :: neven, nextra, ce, i,itmp
   639)     integer :: beg,end
   640)     logical :: eqf = .false.
   641)     real, dimension(:), allocatable :: dist
   642)     integer, dimension(1) :: tmp
   643)     real :: mx,my,mz
   644)     
   645)   
   646)     !receive nodes from master
   647)     call MPI_BCAST(nnodes, 1, MPI_INTEGER , 0, E4D_COMM, ierr )
   648)   
   649) 
   650)     allocate(nodes(nnodes,3),nbounds(nnodes))
   651)     call MPI_BCAST(nodes, nnodes*3, MPI_REAL , 0, E4D_COMM, ierr )
   652)     call MPI_BCAST(nbounds, nnodes, MPI_INTEGER , 0, E4D_COMM, ierr )
   653) 
   654)     !receive elements from master
   655)     call MPI_BCAST(nelem, 1, MPI_INTEGER, 0,E4D_COMM,ierr)
   656)     allocate(elements(nelem,4))
   657)     call MPI_BCAST(elements, nelem*4,MPI_INTEGER,0,E4D_COMM,ierr)
   658) 
   659)     !receive the assignments
   660)     call receive_dists
   661)     
   662)     !build the mesh interpolation matrix for my elements
   663)     call send_my_mesh_interp
   664) 
   665)     !build A_map and delA.......................................
   666)     call build_delA
   667) 
   668)     !!Initialize the Petsc A matrix
   669)     call MatCreateSeqAIJ(PETSC_COMM_SELF,nnodes,nnodes,d_nz,d_nnz,A, &
   670)                          perr);CHKERRQ(perr)
   671)     call MatSetFromOptions(A,perr);CHKERRQ(perr)
   672)  
   673)     !Get the electrode ownership indexes 
   674)     call get_electrode_nodes
   675)     call MatGetType(A,tp,perr);CHKERRQ(perr)
   676)     
   677)     !Set up the source and solution vectors
   678)     call VecCreate(PETSC_COMM_SELF,X,perr);CHKERRQ(perr)
   679)     call VecSetSizes(X,nnodes,PETSC_DECIDE,perr);CHKERRQ(perr)
   680)     call VecSetFromOptions(X,perr);CHKERRQ(perr)
   681)  
   682)     !allocate and setup the pole solution vector
   683)     call VecCreate(PETSC_COMM_SELF,psol,perr);CHKERRQ(perr)
   684)     call VecSetSizes(psol,nnodes,PETSC_DECIDE,perr);CHKERRQ(perr)
   685)     call VecSetFromOptions(psol,perr);CHKERRQ(perr)
   686)     allocate(poles(nnodes,my_ne))
   687)     poles = 0
   688) 
   689)     call VecCreate(PETSC_COMM_SELF,B,perr);CHKERRQ(perr)
   690)     call VecSetSizes(B,nnodes,PETSC_DECIDE,perr);CHKERRQ(perr)
   691)     call VecSetFromOptions(B,perr);CHKERRQ(perr)
   692) 
   693)     !A_map  and S_map given the coupling info so we can deallocate the 
   694)     !nodes and elements
   695)     deallocate(nodes)
   696)     deallocate(elements)
   697)     allocate(sigma(nelem))
   698)     
   699)   end subroutine setup_frun
   700)   !__________________________________________________________________
   701)   
   702)   !__________________________________________________________________
   703)   subroutine receive_dists
   704)     implicit none
   705)     
   706)     integer :: i
   707)    
   708)     if (allocated(eind)) deallocate(eind)
   709)     if (allocated(e_pos)) deallocate(e_pos)
   710)     if (allocated(jind)) deallocate(jind)
   711)     allocate(eind(n_rank-1,2))
   712)     allocate(jind(n_rank-1,2))
   713)     
   714)     call MPI_BCAST(tne, 1, MPI_INTEGER , 0, E4D_COMM, ierr )
   715)     allocate(e_pos(tne,4))
   716)     call MPI_BCAST(e_pos,4*tne,MPI_REAL,0,E4D_COMM,ierr)
   717)     call MPI_BCAST(eind,2*(n_rank-1),MPI_INTEGER , 0, E4D_COMM, ierr )
   718)     call MPI_BCAST(jind,2*(n_rank-1),MPI_INTEGER , 0, E4D_COMM, ierr )
   719)     call MPI_BCAST(pfnx,1,MPI_INTEGER,0,E4D_COMM, ierr)
   720)     call MPI_BCAST(pfny,1,MPI_INTEGER,0,E4D_COMM, ierr)
   721)     call MPI_BCAST(pfnz,1,MPI_INTEGER,0,E4D_COMM, ierr)
   722)     allocate(pfxcb(pfnx+1),pfycb(pfny+1),pfzcb(pfnz+1))
   723)     call MPI_BCAST(pfxcb,pfnx+1,MPI_REAL,0,E4D_COMM,ierr)
   724)     call MPI_BCAST(pfycb,pfny+1,MPI_REAL,0,E4D_COMM,ierr)
   725)     call MPI_BCAST(pfzcb,pfnz+1,MPI_REAL,0,E4D_COMM,ierr)
   726)     my_ne = eind(my_rank,2)-eind(my_rank,1)+1
   727)    
   728)   end subroutine receive_dists
   729)   !__________________________________________________________________
   730) 
   731)   !____________________________________________________________________
   732)   subroutine send_my_mesh_interp
   733)     implicit none
   734)     integer :: i,j,my_nelem,cnt,n_tets,indx,indy,indz,k,cct,cpos
   735)     real :: xmn,xmx,ymn,ymx,zmn,zmx,imn,imx
   736)     logical, dimension(:), allocatable :: flags
   737)     real, dimension(:), allocatable :: midx,midy,midz,w
   738)     integer, dimension(:), allocatable :: rw,v
   739)     real, dimension(9,3) :: pts
   740)     real, dimension(4,3) :: mfaces
   741)     real, dimension(3) :: vc
   742)     real, dimension(8) :: wi,C
   743)     integer, dimension(8) :: vi
   744)     real, dimension(72) :: twi
   745)     integer, dimension(72) :: tvi
   746) 
   747)     
   748)     !!get the pf mesh midpoints
   749)     xmn=1e15
   750)     ymn=xmn
   751)     zmn=xmn
   752)     xmx=-1e15
   753)     ymx=xmx
   754)     zmx=xmx
   755)    
   756) 
   757)     allocate(midx(pfnx),midy(pfny),midz(pfnz))
   758)     do i=1,pfnx
   759)        midx(i) = 0.5*(pfxcb(i+1)+pfxcb(i))
   760)        if(midx(i)<xmn) xmn=midx(i)
   761)        if(midx(i)>xmx) xmx=midx(i)
   762)       
   763)     end do
   764)    
   765)     do i=1,pfny
   766)        midy(i) = 0.5*(pfycb(i+1)+pfycb(i))
   767)        if(midy(i)<ymn) ymn=midy(i)
   768)        if(midy(i)>ymx) ymx=midy(i)
   769)     
   770)     end do
   771) 
   772)     do i=1,pfnz
   773)        midz(i) = 0.5*(pfzcb(i+1)+pfzcb(i))
   774)        if(midz(i)<zmn) zmn=midz(i)
   775)        if(midz(i)>zmx) zmx=midz(i)
   776) 
   777)     end do
   778) 
   779)     nodes(:,1)=nodes(:,1)+xorig
   780)     nodes(:,2)=nodes(:,2)+yorig
   781)     nodes(:,3)=nodes(:,3)+zorig
   782)     
   783)     my_nelem = jind(my_rank,2)-jind(my_rank,1)+1
   784)     allocate(flags(my_nelem))
   785) 
   786)     flags=.true.
   787)     cnt=0
   788)     n_tets=0
   789)    
   790) 
   791)     !!find which elements we need to map
   792)     do i=jind(my_rank,1),jind(my_rank,2)
   793)        cnt=cnt+1
   794)        imx=maxval(nodes(elements(i,1:4),1))
   795)        imn=minval(nodes(elements(i,1:4),1))
   796)        if((imn<xmn .or. imx>xmx) .and. (pfnx > 1)) then
   797)           flags(cnt)=.false.
   798) 
   799)           goto 100
   800)        end if
   801)      
   802)        imx=maxval(nodes(elements(i,1:4),2))
   803)        imn=minval(nodes(elements(i,1:4),2))
   804)        if((imn<ymn .or. imx>ymx) .and. (pfny > 1)) then
   805)           flags(cnt)=.false.
   806)           goto 100
   807)        end if
   808)   
   809) 
   810)        imx=maxval(nodes(elements(i,1:4),3))
   811)        imn=minval(nodes(elements(i,1:4),3)) 
   812)     
   813)        if((imn<zmn .or. imx>zmx) .and. (pfnz>1)) then
   814)           flags(cnt)=.false.
   815)           goto 100
   816)        end if
   817) 
   818)        n_tets=n_tets+1
   819) 100    continue
   820)     end do
   821)  
   822)     allocate(rw(72*n_tets),v(72*n_tets),w(72*n_tets))
   823)     cnt=0
   824)     cct=0
   825)     cpos=0
   826)     do i=jind(my_rank,1),jind(my_rank,2)
   827)        cct=cct+1
   828)        if(flags(cct)) then
   829)           !get the nine points
   830)           cnt=cnt+1
   831)           do j=1,3
   832)              pts(1,j)=.25*sum(nodes(elements(i,1:4),j)) !tet midpoint
   833)           end do
   834)         
   835)           do j=1,3
   836)              mfaces(1,j) = sum(nodes(elements(i,[1,2,3]),j))/3
   837)              mfaces(2,j) = sum(nodes(elements(i,[1,2,4]),j))/3
   838)              mfaces(3,j) = sum(nodes(elements(i,[1,3,4]),j))/3
   839)              mfaces(4,j) = sum(nodes(elements(i,[2,3,4]),j))/3       
   840)           end do
   841)           do j=1,4
   842)              vc=mfaces(j,:)-pts(1,:)
   843)              pts(j+1,:) = pts(1,:) + 0.5*vc
   844)              vc=nodes(elements(i,j),:)-pts(1,:)
   845)              pts(j+5,:) = pts(1,:) + 0.5*vc
   846)           end do
   847)           
   848) 
   849)           !get the weighting function for each point
   850)           do j=1,9
   851)              if(pfnx .eq. 1) then
   852)                 indx = 1
   853)              else
   854)                 do k=1,pfnx
   855)                    if(midx(k)>pts(j,1)) then
   856)                       indx=k-1
   857)                       exit
   858)                    end if
   859)                 end do
   860)              end if
   861) 
   862)              if(pfny.eq.1) then
   863)                 indy = 1
   864)              else
   865)                 do k=1,pfny
   866)                    if(midy(k)>pts(j,2)) then
   867)                       indy=k-1
   868)                       exit
   869)                    end if
   870)                 end do
   871)              end if
   872)              
   873)              if(pfnz .eq.1) then
   874)                 indz = 1
   875)              else
   876)                 do k=1,pfnz
   877)                    if(midz(k)>pts(j,3)) then
   878)                       indz=k-1
   879)                       exit
   880)                    end if
   881)                 end do
   882)              end if
   883) 
   884)              
   885)              if(pfny .eq.1) then
   886)                 !for 2D problem in x,z
   887)                 vi(5) = indx + (indy-1)*pfnx + (indz-1)*pfnx*pfny
   888)                 vi(6) = vi(5)
   889)                 vi(8) = vi(5)+1
   890)                 vi(7) = vi(8);
   891)                 vi(1) = vi(5)+pfnx
   892)                 vi(2) = vi(1);
   893)                 vi(3) = vi(1)+1;
   894)                 vi(4) = vi(3);
   895)                 
   896)              else
   897)                 vi(5) = indx + (indy-1)*pfnx + (indz-1)*pfnx*pfny
   898)                 vi(8) = vi(5)+1
   899)                 vi(6) = vi(5)+pfnx
   900)                 vi(7) = vi(6)+1
   901)                 vi(1:4) = vi(5:8) + pfnx*pfny
   902)              end if
   903) 
   904)               if(pfnx .eq.1) then
   905)                  C(1)=0
   906)               else
   907)                  C(1)=(pts(j,1)-midx(indx))/(midx(indx+1)-midx(indx));
   908)               end if
   909)               C(2)=C(1);
   910)               C(3)=C(1);
   911)               C(4)=C(1);
   912) 
   913)               if(pfny .eq.1) then
   914)                  C(5)=0
   915)               else
   916)                  C(5)=(pts(j,2)-midy(indy))/(midy(indy+1)-midy(indy));
   917)               end if
   918)               C(6)=C(5);
   919)               
   920)               if(pfnz .eq.1) then
   921)                  C(7)=0
   922)               else
   923)                  C(7)=(pts(j,3)-midz(indz))/(midz(indz+1)-midz(indz));
   924)               end if
   925) 
   926)               wi(1)=(1-C(1))*(1-C(6))*C(7);
   927)               wi(2)=(1-C(2))*C(6)*C(7);
   928)               wi(3)=C(2)*C(6)*C(7);
   929)               wi(4)=C(1)*C(7)*(1-C(6));
   930)               wi(5)=(1-C(3))*(1-C(5))*(1-C(7));
   931)               wi(6)=C(5)*(1-C(4))*(1-C(7));
   932)               wi(7)=C(4)*C(5)*(1-C(7));
   933)               wi(8)=C(3)*(1-C(5))*(1-C(7));
   934)            
   935)               !write(*,*) vi
   936)               !write(*,*) wi
   937) 
   938)               tvi(8*(j-1)+1:8*j)=vi
   939)               twi(8*(j-1)+1:8*j)=wi
   940)               
   941)            end do
   942)            
   943)            !consolidate the weights
   944)            do j=1,71
   945)               do k=j+1,72
   946)                  if(tvi(k)==tvi(j)) then
   947)                     twi(j)=twi(j)+twi(k)
   948)                     twi(k)=0
   949)                     tvi(k)=0;
   950)                  end if
   951)               end do
   952)            end do
   953)            twi=twi/9
   954)           
   955)            do j=1,72
   956)               if(tvi(j) .ne. 0) then
   957)                  cpos=cpos+1;
   958)                  rw(cpos)=i
   959)                  v(cpos)=tvi(j)
   960)                  w(cpos)=twi(j)
   961)               end if
   962)            end do
   963)         end if
   964)      
   965)     end do
   966)   
   967)     call MPI_SEND(cpos,1,MPI_INTEGER,0,my_rank,E4D_COMM, ierr)
   968)     call MPI_SEND(rw(1:cpos),cpos,MPI_INTEGER,0,my_rank,E4D_COMM,ierr)
   969)     call MPI_SEND(v(1:cpos),cpos,MPI_INTEGER,0,my_rank,E4D_COMM,ierr)
   970)     call MPI_SEND(w(1:cpos),cpos,MPI_REAL,0,my_rank,E4D_COMM,ierr)
   971)     
   972)     deallocate(flags,midx,midy,midz,rw,v,w)
   973)     nodes(:,1)=nodes(:,1)-xorig
   974)     nodes(:,2)=nodes(:,2)-yorig
   975)     nodes(:,3)=nodes(:,3)-zorig
   976)     
   977)   end subroutine send_my_mesh_interp
   978)   !____________________________________________________________________
   979)   !____________________________________________________________________
   980)   subroutine build_delA
   981) 
   982)     use e4d_mat_inv_module
   983) 
   984)     implicit none
   985)     
   986)    
   987)     integer :: nnds,i,j,k,l,mcount,rn,cn,ncl,jj
   988)     integer :: lrow
   989)     real, dimension(4,4) :: atild,atildi,temp_a
   990)     integer, dimension(4,4) :: indx
   991)     integer, dimension(:,:),allocatable :: A_tmp
   992)     real :: x1,x2,x3,x4,y1,y2,y3,y4,z1,z2,z3,z4,evol
   993)     real :: A_kij,lA_kij
   994)     real, dimension(:), allocatable :: t1vec,t2vec
   995)     logical :: ilow,iup
   996)     
   997)     !allocate the petsc matrix preallocation vectors
   998)     allocate(d_nnz(nnodes))
   999)     d_nnz=0
  1000)     
  1001)     allocate(rows(10*nelem),cols(10*nelem))
  1002)     allocate(A_map(10*nelem),S_map(10*nelem))
  1003)     allocate(delA(10*nelem))
  1004)     
  1005)     !!COUNT THE NUMBER OF NODE PAIRS (I.E. COUPLING MATRIX ELEMENTS) AND
  1006)     !!REALLOCATE THE COUPLING MATRIX STORAGE VECTORS
  1007)     !!These pairs will also accomodate source dependent boundaries
  1008)     rows(1)=elements(1,1)
  1009)     cols(1)=elements(1,1)
  1010)     nvals=1
  1011)     mcount = 0
  1012)     A_map = 0
  1013)     nnds=maxval(elements(:,:))
  1014)     
  1015)     !!Allocate a temporary matrix to aid in building the 
  1016)     !!mapping vector
  1017)     allocate(A_tmp(nnds,80))
  1018)     A_tmp=0
  1019)     A_tmp(rows(1),1) = 1
  1020)     A_tmp(rows(1),2) = nvals
  1021)     A_tmp(rows(1),3) = cols(1) 
  1022)     !!loop over the elements
  1023)      
  1024)     do k=1,nelem    
  1025)        
  1026)        !!loop over each node in this element
  1027)        do i=1,4
  1028)          
  1029)           !!rn is a row index of the coupling matrix
  1030)           rn=elements(k,i)    
  1031)           
  1032)           !!loop over each node in this element
  1033)           do j=1,4
  1034)              
  1035)              !!cn is a column index of the coupling matrix
  1036)              cn=elements(k,j)
  1037)         
  1038)           
  1039) 
  1040)              !!check to see if this pair (rn,cn) is in the lower
  1041)              !!triangle. If not, go to the next pair
  1042)              if (cn <= rn) then
  1043)                 mcount=mcount+1
  1044)                 !!loop over each pair found thus far (nvals pairs)
  1045)                 !!and determine if the pair (rn,cn) is already 
  1046)                 !!represeneted
  1047)                 if (A_tmp(rn,1) == 0) then
  1048)                    nvals = nvals + 1
  1049)                    A_tmp(rn,1) = 1
  1050)                    A_tmp(rn,2) = nvals
  1051)                    A_tmp(rn,3) = cn
  1052)                    rows(nvals) = rn
  1053)                    cols(nvals) = cn
  1054)                    A_map(mcount) = nvals
  1055)                            
  1056)                 else
  1057)                    
  1058)                    do l=1,A_tmp(rn,1)
  1059)                       if (cn == A_tmp(rn,2*l+1)) then
  1060)                          A_map(mcount) = A_tmp(rn,2*l)
  1061)                          goto 11
  1062)                       end if
  1063)                    end do
  1064)                    
  1065)                    !!if were here then no column index was found
  1066)                    !!for this row, so we'll add one
  1067)                    nvals = nvals+1
  1068)                    ncl = A_tmp(rn,1) + 1
  1069)                    A_tmp(rn,1) = ncl
  1070)                    A_tmp(rn,2*ncl) = nvals
  1071)                    A_tmp(rn,2*ncl+1) = cn
  1072)                    rows(nvals) = rn
  1073)                    cols(nvals) = cn
  1074)                    A_map(mcount) = nvals
  1075)                 
  1076) 11                 continue
  1077)      
  1078)                 end if
  1079) 
  1080)              end if
  1081) 
  1082)           end do
  1083)           
  1084)        end do
  1085)        
  1086)     end do
  1087)     
  1088)     !build the petsc allocation vector
  1089)     do i=1,nnodes
  1090)        !upper part
  1091)        d_nnz(i)=A_tmp(i,1)-1;
  1092)        do j=1,A_tmp(i,1)
  1093)           cn=A_tmp(i,2*j+1)
  1094)           d_nnz(cn)=d_nnz(cn)+1
  1095)        end do
  1096)     end do   
  1097)     deallocate(A_tmp)
  1098)    
  1099)     !!now reallocate for the correct number of matrix values (i.e. pairs)
  1100)     allocate(t1vec(nvals),t2vec(nvals))
  1101)     t1vec(1:nvals)=rows(1:nvals)
  1102)     t2vec(1:nvals)=cols(1:nvals)
  1103)     deallocate(rows,cols)
  1104)     allocate(rows(nvals),cols(nvals))
  1105)     rows=t1vec
  1106)     cols=t2vec
  1107)     deallocate(t1vec,t2vec)
  1108)     allocate(trows(nvals),tcols(nvals))
  1109)     
  1110)     !!LOOP OVER THE ELEMENTS AND BUILD delA
  1111)     mcount = 0
  1112)     delA=0;
  1113)   
  1114)     do k=1,nelem
  1115)    
  1116)        !!get the 4 nodes for this element 
  1117)        x1 = nodes(elements(k,1),1)
  1118)        y1 = nodes(elements(k,1),2)
  1119)        z1 = nodes(elements(k,1),3)
  1120)        
  1121)        x2 = nodes(elements(k,2),1)
  1122)        y2 = nodes(elements(k,2),2)
  1123)        z2 = nodes(elements(k,2),3)
  1124)        
  1125)        x3 = nodes(elements(k,3),1)
  1126)        y3 = nodes(elements(k,3),2)
  1127)        z3 = nodes(elements(k,3),3) 
  1128)        
  1129)        x4 = nodes(elements(k,4),1)
  1130)        y4 = nodes(elements(k,4),2)
  1131)        z4 = nodes(elements(k,4),3)
  1132)        
  1133)        !!get the volume of this element
  1134)        call get_vol (x1,y1,z1,x2,y2,z2,x3,y3,z3,x4,y4,z4,evol)  
  1135)    
  1136)        if (evol==0) goto 12
  1137)        !!build the linear shape function matrix for this element
  1138)        atild(1,:) = (/1,1,1,1/)
  1139)        atild(2,:) = (/x1,x2,x3,x4/)
  1140)        atild(3,:) = (/y1,y2,y3,y4/)
  1141)        atild(4,:) = (/z1,z2,z3,z4/)
  1142)        temp_a=atild
  1143)             
  1144)        !!invert temp_a to get atildi
  1145)        !!note MIGS routine is located in mat_inv.f      
  1146)        call MIGS(temp_a,4,atildi,indx)  
  1147)       
  1148) 12     continue
  1149)        !!use atildi to build the coupling coefficients
  1150)        !!for this element      
  1151)        
  1152)        !!loop over each node in this element
  1153)        do i=1,4
  1154)           !!rn is a row index of the coupling matrix
  1155)           rn=elements(k,i)
  1156)           
  1157)           !!loop over each node in this element
  1158)           do j=1,4
  1159)              
  1160)              !!cn is a column index of the coupling matrix
  1161)              cn=elements(k,j)
  1162)              
  1163)              !!check to see if this pair (rn,cn) is in the lower
  1164)              !!triangle. If not, go to the next pair
  1165)              if (cn <= rn) then
  1166)                 mcount=mcount+1
  1167)                 !!Compute the coupling coefficient for this element
  1168)                 !!and node pair
  1169)                 A_kij=0
  1170)                
  1171)                 do l=2,4
  1172)                    A_kij=A_kij+(atildi(i,l)*atildi(j,l))
  1173)                 end do
  1174)           
  1175)                 !!delA(mcount)=lA_kij*evol;
  1176)                 delA(mcount)=dble(A_kij*evol)
  1177)                 S_map(mcount)=k
  1178)              end if
  1179)              
  1180) 21           continue
  1181)           end do
  1182)           
  1183)        end do
  1184)        
  1185)     end do  
  1186)     
  1187)     
  1188)   end subroutine build_delA
  1189)   !____________________________________________________________________
  1190)   
  1191)   !____________________________________________________________________
  1192)   subroutine get_vol( x1, y1, z1, x2, y2, z2, x3, y3, z3, x4, y4, z4, volume )
  1193) 
  1194)     implicit none
  1195)     
  1196)     real a(4,4)
  1197)     real r4_det
  1198)     real volume
  1199)     real x1
  1200)     real x2
  1201)     real x3
  1202)     real x4
  1203)     real y1
  1204)     real y2
  1205)     real y3
  1206)     real y4
  1207)     real z1
  1208)     real z2
  1209)     real z3
  1210)     real z4
  1211)     
  1212)     a(1:4,1) = (/ x1, x2, x3, x4 /)
  1213)     a(1:4,2) = (/ y1, y2, y3, y4 /)
  1214)     a(1:4,3) = (/ z1, z2, z3, z4 /)
  1215)     a(1:4,4) = (/ 1.0E+00, 1.0E+00, 1.0E+00, 1.0E+00 /)
  1216)     r4_det = det4(a)
  1217)     volume = abs ( r4_det ) / 6.0E+00
  1218)     
  1219)     return
  1220)   end subroutine get_vol
  1221)   !____________________________________________________________________
  1222)     
  1223)   !____________________________________________________________________
  1224)   function det4 ( a1 )
  1225)     implicit none
  1226)     
  1227)     real :: a1(4,4)
  1228)     real*8 :: a(4,4)
  1229)     real :: det4
  1230)     integer ::i
  1231)     real*8 :: c1,c2,c3,c4
  1232)     
  1233)     a = dble(a1)
  1234)     c1=a(1,1)*(a(2,2)*(a(3,3)*a(4,4)-a(3,4)*a(4,3))-a(2,3)*&
  1235)          &(a(3,2)*a(4,4)-a(3,4)*a(4,2))+a(2,4)*(a(3,2)*a(4,3)-a(3,3)*a(4,2))) 
  1236)     
  1237)     c2=a(1,2)*(a(2,1)*(a(3,3)*a(4,4)-a(3,4)*a(4,3))-a(2,3)*&
  1238)          &(a(3,1)*a(4,4)-a(3,4)*a(4,1))+a(2,4)*(a(3,1)*a(4,3)-a(3,3)*a(4,1)))
  1239)     
  1240)     c3=a(1,3)*(a(2,1)*(a(3,2)*a(4,4)-a(3,4)*a(4,2))-a(2,2)*&
  1241)          &(a(3,1)*a(4,4)-a(3,4)*a(4,1))+a(2,4)*(a(3,1)*a(4,2)-a(3,2)*a(4,1)))
  1242)     
  1243)     c4=a(1,4)*(a(2,1)*(a(3,2)*a(4,3)-a(3,3)*a(4,2))-a(2,2)*&
  1244)          &(a(3,1)*a(4,3)-a(3,3)*a(4,1))+a(2,3)*(a(3,1)*a(4,2)-a(3,2)*a(4,1)))
  1245)     
  1246)     det4 = real(c1 - c2 + c3 - c4)
  1247)     
  1248)     return
  1249)   end function det4
  1250)   !__________________________________________________________________________
  1251) 
  1252)   !________________________________________________________________
  1253)   subroutine send_info_geh
  1254)     implicit none
  1255) 
  1256)     call send_command(86)
  1257)     nelem = 4*4*4
  1258)     nm = nelem
  1259)     print *, 'master:', nm
  1260)     call MPI_BCAST(nm, 1, MPI_INTEGER, 0,E4D_COMM,ierr)
  1261) 
  1262)   end subroutine send_info_geh
  1263)   !________________________________________________________________
  1264) 
  1265)   !__________________________________________________________________
  1266)   subroutine receive_info_geh
  1267)     implicit none
  1268)     integer :: i
  1269)     
  1270)     call MPI_BCAST(nm, 1, MPI_INTEGER, 0,E4D_COMM,ierr)
  1271)     allocate(sigma(nm))
  1272)     nelem = nm
  1273)     sigma = 0.d0
  1274)     !print *, 'slave:', nm
  1275)     
  1276)   end subroutine receive_info_geh
  1277)   !__________________________________________________________________
  1278)   
  1279)   !________________________________________________________________
  1280)   subroutine destroy_e4d
  1281)     implicit none
  1282)     if (allocated(e4d_ranks)) then
  1283)       deallocate(e4d_ranks)
  1284)     endif
  1285)     if (allocated(pf_e4d_ranks)) then
  1286)       deallocate(pf_e4d_ranks)
  1287)     endif
  1288)     if (allocated(map_inds)) then
  1289)       deallocate(map_inds)
  1290)     endif
  1291)     if (allocated(s_conf)) then
  1292)       deallocate(s_conf)
  1293)     endif
  1294)     if (allocated(eind)) then
  1295)       deallocate(eind)
  1296)     endif
  1297)     if (allocated(nbounds)) then
  1298)       deallocate(nbounds)
  1299)     endif
  1300)     if (allocated(zones)) then
  1301)       deallocate(zones)
  1302)     endif
  1303)     if (allocated(elements)) then
  1304)       deallocate(elements)
  1305)     endif
  1306)     if (allocated(faces)) then
  1307)       deallocate(faces)
  1308)     endif
  1309)     if (allocated(e_nods)) then
  1310)       deallocate(e_nods)
  1311)     endif
  1312)     if (allocated(rows)) then
  1313)       deallocate(rows)
  1314)     endif
  1315)     if (allocated(cols)) then
  1316)       deallocate(cols)
  1317)     endif
  1318)     if (allocated(trows)) then
  1319)       deallocate(trows)
  1320)     endif
  1321)     if (allocated(tcols)) then
  1322)       deallocate(tcols)
  1323)     endif
  1324)     if (allocated(A_map)) then
  1325)       deallocate(A_map)
  1326)     endif
  1327)     if (allocated(S_map)) then
  1328)       deallocate(S_map)
  1329)     endif
  1330)     if (allocated(my_drows)) then
  1331)       deallocate(my_drows)
  1332)     endif
  1333)     if (allocated(e_pos)) then
  1334)       deallocate(e_pos)
  1335)     endif
  1336)     if (allocated(nodes)) then
  1337)       deallocate(nodes)
  1338)     endif
  1339)     if (allocated(poles)) then
  1340)       deallocate(poles)
  1341)     endif
  1342)     if (allocated(pf_porosity)) then
  1343)       deallocate(pf_porosity)
  1344)     endif
  1345)     if (allocated(pf_tracer)) then
  1346)       deallocate(pf_tracer)
  1347)     endif
  1348)     if (allocated(pf_saturation)) then
  1349)       deallocate(pf_saturation)
  1350)     endif
  1351)     if (allocated(pf_saturation_0)) then
  1352)       deallocate(pf_saturation_0)
  1353)     endif
  1354)     if (allocated(pf_temperature)) then
  1355)       deallocate(pf_temperature)
  1356)     endif
  1357)     if (allocated(sigma)) then
  1358)       deallocate(sigma)
  1359)     endif
  1360)     if (allocated(dpred)) then
  1361)       deallocate(dpred)
  1362)     endif
  1363)     if (allocated(dobs)) then
  1364)       deallocate(dobs)
  1365)     endif
  1366)     if (allocated(sd)) then
  1367)       deallocate(sd)
  1368)     endif
  1369)     if (allocated(my_dvals)) then
  1370)       deallocate(my_dvals)
  1371)     endif
  1372)     if (allocated(map)) then
  1373)       deallocate(map)
  1374)     endif
  1375)     if (allocated(base_sigma)) then
  1376)       deallocate(base_sigma)
  1377)     endif
  1378)     if (allocated(d_nnz)) then
  1379)       deallocate(d_nnz)
  1380)     endif
  1381)     if (allocated(delA)) then
  1382)       deallocate(delA)
  1383)     endif
  1384)   end subroutine destroy_e4d
  1385)   !________________________________________________________________
  1386)  
  1387) end module e4d_setup
  1388) 

generated by
Intel(R) C++/Fortran Compiler code-coverage tool
Web-Page Owner: Nobody