Skip to content
Snippets Groups Projects
DOUAR.f90 71.1 KiB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!                                                                              |
!              ||===\\                                                         | 
!              ||    \\                                                        |
!              ||     ||   //==\\   ||  ||   //==||  ||/==\\                   |
!              ||     ||  ||    ||  ||  ||  ||   ||  ||    ||                  |
!              ||    //   ||    ||  ||  ||  ||   ||  ||                        |
!              ||===//     \\==//    \\==\\  \\==\\  ||                        |
!                                                                              |
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|

program DOUAR

use threads
use definitions

implicit none

INCLUDE 'mpif.h'

integer i,j,is,iter,istep,niter,nedge,nedgen 
integer ntriangle,ndof,ndoft,material0,nsurfacen
integer numarg,iargc
integer ierr,nproc,iproc,k,nrem
integer current_level
integer ie_ov,ie_loc,ie_level,iter_nl
integer err,iunit,cltemp
integer nremove, ninject,nnmax,nadd,naddp,ref_count
integer nleaves_new_mem1
integer nleaves_new_mem2
integer nleaves_old_mem1
integer nleaves_old_mem2
integer, external :: ioctree_number_of_elements

double precision x,y,z,z0,u,v,w
double precision exx,eyy,ezz,exy,eyz,ezx,e2dmax
double precision duvw,uvw,dtcourant,current_time
double precision umax,xxx,yyy,zzz,x0_leaf,y0_leaf,z0_leaf
double precision total,step,inc
real  time1,time2,time3 
      
character     :: ic*7
character*3   :: ciproc
character*4   :: cs,cs2
character*40  :: string
character*72  :: shift

logical converge,increase_current_level,velocity_converged,usecourant

type (sheet),dimension(:),allocatable::surface,surface0
type (octreev) ov
type (octreelsf) olsf
type (octreesolve) osolve
type (material),dimension(:),allocatable::mat
type (void) vo
type (cloud) cl
type (topology),dimension(:),allocatable::tpl
type (box),dimension(:),allocatable::boxes
type (cross_section),dimension(:),allocatable::sections
type (face),dimension(6) :: cube_faces
type (edge),dimension(:),allocatable::ed,edswap
type (parameters) params
type (thread) threadinfo

integer n_iproc_st,n_iproc_end,n_iproc
integer ldb,nrhs,n,nz_loc
integer, dimension(:), allocatable :: ia,ja
logical, dimension(:), allocatable :: iproc_col
double precision, dimension(:), allocatable :: avals
double precision, dimension(:,:), allocatable :: b
double precision, dimension(:), allocatable :: weightel

shift=' '

ndof=3
ndoft=1
params%mpe=8
nrhs = 1

nleaves_old_mem1=0
nleaves_old_mem2=0

!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|

call mpi_init(ierr)
call mpi_comm_size (mpi_comm_world,nproc,ierr)
call mpi_comm_rank (mpi_comm_world,iproc,ierr)

call write_streepjes(6,2)


!=====[allocate threadinfo and open mpi log and memory heap files]=========

params%nmpi=nproc

call int_to_char(ciproc,3,iproc)

threadinfo%Logunit=1000+iproc
open  (unit=threadinfo%Logunit,file='./DEBUG/mpilogs/Log_'//ciproc//'.dat',status='replace')
write(threadinfo%Logunit,'(a,i3)') 'Log file of mpi process',iproc

call heap_init(threadinfo,2000+iproc,'./DEBUG/mpilogs/mem_heap_'//ciproc//'.dat')

!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!     is there any input file to douar ?
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|

numarg=iargc()

if (numarg==0) then
   params%infile='input.txt'
   if (iproc.eq.0) write(*,*) 'program called with no argument'
else
   call getarg (1,params%infile)
   if (iproc.eq.0) then
      write(*,*) 'program called with input file ',params%infile
   end if
end if

call write_streepjes(6,2)

!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!     read controlling parameters
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
call read_controlling_parameters (params,threadinfo)

!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!     open and prepare output files
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
 call io_DoRuRe2 (params,'init')

call show_time (total,step,inc,0,'Start of Computations$')
call show_time (total,step,inc,1,'Reading Input$')

allocate (surface         (  params%ns  ),stat=err)    ; if (err.ne.0) call stop_run ('Error alloc surface in main$')
allocate (surface0        (  params%ns  ),stat=err)   ; if (err.ne.0) call stop_run ('Error alloc surface0 in main$')
allocate (mat             (0:params%nmat),stat=err)    ; if (err.ne.0) call stop_run ('Error alloc mat in main$')
allocate (params%materialn(0:params%ns  ),stat=err); if (err.ne.0) call stop_run ('Error alloc materialn in main$')
      
if (params%nboxes.gt.0) then
    allocate (boxes(params%nboxes),stat=err) ; if (err.ne.0) call stop_run ('Error alloc boxes arrays$')
else
    allocate (boxes(1),stat=err) ! necessary to avoid nil size argument
endif
if (params%nsections.gt.0) then
    allocate (sections(params%nsections),stat=err) ; if (err.ne.0) call stop_run ('Error alloc cross_section arrays$')
else
    allocate (sections(1),stat=err) ! necessary to avoid nil size argument
end if

!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!     read input file for all parameters of the run
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
call read_input_file (params,threadinfo,material0,mat,surface,boxes,sections,cube_faces) 

if (params%irestart==0) then
   current_level=params%initial_refine_level
else
   current_level=params%levelmax_oct
end if
 
call show_time (total,step,inc,1,'Problem Setup$')

!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!     either reads or creates the surface(s)
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
call show_time (total,step,inc,1,'define surfaces$')
call define_surface(params,surface,threadinfo,total,step,inc,current_time)

!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!     either reads or creates the cloud  
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
call show_time (total,step,inc,1,'define cloud$')
call  define_cloud(cl,params)

!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!     extract material information to be passed to solver
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
params%materialn(0)=material0
do i=1,params%ns
   params%materialn(i)=surface(i)%material
enddo

!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!    compute plastic parameters 
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
call compute_plastic_params (params,mat)

!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!     either reads or creates the velocity octree
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
call show_time (total,step,inc,1,'define velocity octree$')
call define_ov (ov,params,threadinfo)

istep=1+params%irestart 

!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!     is the CFL condition used for the timestep ?
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|

usecourant = (params%dt .lt. 0.d0)

!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!     beginning of time stepping
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
do while (istep.le.params%nstep) 
   
   call write_streepjes(6,2)
   call write_streepjes(8,2)
   call show_time (total,step,inc,-istep,'Start of Step $')
   if (iproc.eq.0)  write(8,*)  'Doing step ',istep
   call write_streepjes(6,2)
   call write_streepjes(8,2)
   call heap_hop1(threadinfo,istep)

   !---------------------------------------------------------------------------|
   !---------------------------------------------------------------------------|
   !     refining surfaces
   !---------------------------------------------------------------------------|
   !---------------------------------------------------------------------------|
   
   call show_time (total,step,inc,1,'Loop over surfaces$')
   
   do is=1,params%ns
      call int_to_char(ic,2,is)
      call show_time (total,step,inc,1,'surface '//ic(1:2)//'$')

      if (iproc.eq.0 ) then
         write (8,*) 'Surface ',is,' has ',surface(is)%nsurface,'particles before refinement'
         if (params%debug>=1) write(*,'(a,i2,a,i8,a)') shift//'S.',is,': ',surface(is)%nsurface,'points'
      end if

      if (surface(is)%activation_time.ge.current_time+tiny(current_time)) then
         !------------------------------------------------------------------------!
         ! if surface is slave to top surface, no refinement
         !------------------------------------------------------------------------!
         deallocate(surface(is)%x,surface(is)%y,surface(is)%z)
         deallocate(surface(is)%xn,surface(is)%yn,surface(is)%zn)
         deallocate(surface(is)%r,surface(is)%s)
         deallocate(surface(is)%icon)
         surface(is)%nsurface=surface(1)%nsurface
         surface(is)%nt=surface(1)%nt
         allocate (surface(is)%x(surface(is)%nsurface),surface(is)%y(surface(is)%nsurface),surface(is)%z(surface(is)%nsurface))
         allocate (surface(is)%xn(surface(is)%nsurface),surface(is)%yn(surface(is)%nsurface),surface(is)%zn(surface(is)%nsurface))
         allocate (surface(is)%r(surface(is)%nsurface),surface(is)%s(surface(is)%nsurface))
         allocate (surface(is)%icon(3,surface(is)%nt))
         surface(is)%x=surface(1)%x
         surface(is)%y=surface(1)%y
         surface(is)%z=surface(1)%z
         surface(is)%xn=surface(1)%xn
         surface(is)%yn=surface(1)%yn
         surface(is)%zn=surface(1)%zn
         surface(is)%r=surface(1)%r
         surface(is)%s=surface(1)%s
         surface(is)%icon=surface(1)%icon(:,1:surface(is)%nt)
      else
         !------------------------------------------------------------------------!
         ! if surface is no slave to top surface, do refine
         !------------------------------------------------------------------------!
         ref_count=1   
         nadd=1
         do while (nadd>0)
            call show_time (total,step,inc,1,'refine surface$')
            call refine_surface(params,surface(is),threadinfo,nadd,is,istep,ref_count) 
            call show_time (total,step,inc,1,'check delaunay$')
            call check_delaunay (params,surface(is),is,istep,'delaunay_after_refine',ref_count)
            call DoRuRe_surf_stats (params%doDoRuRe,istep,ref_count,is,surface(is)%nt,surface(is)%nsurface,nedge,nadd)
            ref_count=ref_count+1
         end do
      endif

      if (iproc.eq.0 ) then
         write (8,*) 'Surface ',is,' has ',surface(is)%nsurface,'particles after refinement'
         if (params%debug>=1) write(*,'(a,i2,a,i8,a)') shift//'S.',is,': ',surface(is)%nsurface,'points'
      end if

   end do

   !---------------------------------------------------------------------------|
   !---------------------------------------------------------------------------|
   !     Stores the geometry of the surfaces for the mid-point iterative scheme
   !---------------------------------------------------------------------------|
   !---------------------------------------------------------------------------|
   do is=1,params%ns
      surface0(is)%nsurface=surface(is)%nsurface
      allocate (surface0(is)%r(surface0(is)%nsurface),stat=threadinfo%err)       ; call heap (threadinfo,'surface0(is)%r','main',size(surface0(is)%r),'dp',+1)
      allocate (surface0(is)%s(surface0(is)%nsurface),stat=threadinfo%err)       ; call heap (threadinfo,'surface0(is)%s','main',size(surface0(is)%s),'dp',+1)
      allocate (surface0(is)%x(surface0(is)%nsurface),stat=threadinfo%err)       ; call heap (threadinfo,'surface0(is)%x','main',size(surface0(is)%x),'dp',+1)
      allocate (surface0(is)%y(surface0(is)%nsurface),stat=threadinfo%err)       ; call heap (threadinfo,'surface0(is)%y','main',size(surface0(is)%y),'dp',+1)
      allocate (surface0(is)%z(surface0(is)%nsurface),stat=threadinfo%err)       ; call heap (threadinfo,'surface0(is)%z','main',size(surface0(is)%z),'dp',+1)
      allocate (surface0(is)%xn(surface0(is)%nsurface),stat=threadinfo%err)      ; call heap (threadinfo,'surface0(is)%xn','main',size(surface0(is)%xn),'dp',+1)
      allocate (surface0(is)%yn(surface0(is)%nsurface),stat=threadinfo%err)      ; call heap (threadinfo,'surface0(is)%yn','main',size(surface0(is)%yn),'dp',+1)
      allocate (surface0(is)%zn(surface0(is)%nsurface),stat=threadinfo%err)      ; call heap (threadinfo,'surface0(is)%zn','main',size(surface0(is)%zn),'dp',+1)
      surface0(is)%r=surface(is)%r
      surface0(is)%s=surface(is)%s
      surface0(is)%x=surface(is)%x
      surface0(is)%y=surface(is)%y
      surface0(is)%z=surface(is)%z
      surface0(is)%xn=surface(is)%xn
      surface0(is)%yn=surface(is)%yn
      surface0(is)%zn=surface(is)%zn
   enddo

   !---------------------------------------------------------------------------|
   !---------------------------------------------------------------------------|
   !     beginning of grid iterations
   !---------------------------------------------------------------------------|
   !---------------------------------------------------------------------------|

   iter =1

   do while (iter.le.abs(params%griditer).or.iter.eq.1)

      call write_streepjes(6,1)
      call write_streepjes(8,1)
      call show_time (total,step,inc,-iter,' Start of iteration $')
      if(iproc.eq.0) write(8,*)  'Doing iteration ',iter
      call write_streepjes(6,1)
      call write_streepjes(8,1)
      call heap_hop2(threadinfo,iter)
      call show_time (total,step,inc,1,'Create octree solve$')

      converge=.false.
!      if (params%griditer .lt. 0 .and. current_level==params%levelmax_oct) converge=.true.

      if (iproc.eq.0 .and. params%debug>=1) then
         write(*,'(a,L)') shift//'(1) params%griditer < 0               ',(params%griditer.lt.0) 
         write(*,'(a,L)') shift//'(2) current_level==params%levelmax_oct',(current_level==params%levelmax_oct) 
         write(*,'(a,L)') shift//'(3) increase_current_level            ',increase_current_level
         write(*,'(a,L)') shift//'converge = (1) & (2) & (3) -> ',(params%griditer .lt. 0 .and. current_level==params%levelmax_oct .and. increase_current_level) 
      end if


      if (params%griditer .lt. 0 .and. current_level==params%levelmax_oct .and. increase_current_level) converge=.true.
      if (iproc.eq.0 .and. converge .and. params%debug>=1) then
         write(*,'(a)') shift//'+++++++++++++++++'
         write(*,'(a)') shift//'grid conv reached'
         write(*,'(a)') shift//'+++++++++++++++++'
      end if


      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      !     creation of osolve
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      ! we create the large "solve" octree, ie the octree on which the finite 
      ! elements will be built.

      osolve%noctree=params%noctreemax

      allocate (osolve%octree(osolve%noctree),stat=threadinfo%err)      ; call heap (threadinfo,'osolve%octree','main',size(osolve%octree),'int',+1)

      call octree_init (osolve%octree,osolve%noctree)
      call octree_create_uniform (osolve%octree,osolve%noctree,params%leveluniform_oct)

      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      !     refinement/improvement of osolve
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      call  improve_osolve (osolve,ov,params,threadinfo,istep,iter,total,step,inc, &
                            current_level,boxes,cube_faces,increase_current_level )

      nleaves_new_mem1=osolve%octree(2)

      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      !     find connectivity for osolve 
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      call show_time (total,step,inc,1,'Build osolve icon$')
      
      osolve%nleaves = ioctree_number_of_elements (osolve%octree,osolve%noctree)
      osolve%nnode   = osolve%nleaves*3
      osolve%nlsf    = params%ns

      if (iproc.eq.0) write (8,*) osolve%nleaves,' leaves in solve octree '

      allocate (osolve%icon(8,osolve%nleaves),stat=threadinfo%err)        ; call heap (threadinfo,'osolve%icon','main',size(osolve%icon),'int',+1)
      allocate (osolve%x(osolve%nnode),stat=threadinfo%err)               ; call heap (threadinfo,'osolve%x',   'main',size(osolve%x),'dp',+1)
      allocate (osolve%y(osolve%nnode),stat=threadinfo%err)               ; call heap (threadinfo,'osolve%y',   'main',size(osolve%y),'dp',+1)
      allocate (osolve%z(osolve%nnode),stat=threadinfo%err)               ; call heap (threadinfo,'osolve%z',   'main',size(osolve%z),'dp',+1)
      allocate (osolve%lsf(osolve%nnode,osolve%nlsf),stat=threadinfo%err) ; call heap (threadinfo,'osolve%lsf', 'main',size(osolve%lsf),'dp',+1)
      osolve%lsf=0.d0

      call octree_find_node_connectivity (osolve%octree,osolve%noctree, &
                                          osolve%icon,osolve%nleaves,osolve%x,osolve%y,osolve%z,osolve%nnode)

      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      !     imbed the surfaces in osolve
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      call show_time (total,step,inc,1,'Imbed surface in osolve$')
      do is=1,params%ns
         write (ic(1:2),'(i2)') is
         call show_time (total,step,inc,1,'embedding surface '//ic(1:2)//'$')
         call embed_surface_in_octree (osolve,params,surface(is),is,istep,iter,threadinfo)
      enddo

      nleaves_new_mem2=osolve%octree(2)
 
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      !     assess grid convergence 
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      call show_time (total,step,inc,1,'Assessing octree stability$')

      increase_current_level=  (nleaves_old_mem2 .ge. nleaves_new_mem2 .and.  & 
                               (dble(nleaves_old_mem2-nleaves_new_mem2)) / (dble(nleaves_old_mem2+nleaves_new_mem2))*2.d0 &
                               .lt. params%octree_refine_ratio )

!      increase_current_level=  ((dble(nleaves_old_mem2-nleaves_new_mem2)) / (dble(nleaves_old_mem2+nleaves_new_mem2))*2.d0 &
!                               .lt. params%octree_refine_ratio )
 
      if(iproc.eq.0 .and. params%debug>=1) then 
      write(*,'(a,i5,a)') shift//'current refine level ',current_level
      write(*,'(a)')      shift//'After criterion based refinement:'
      write(*,'(a,i8,a)') shift//'nleaves_old= ',nleaves_old_mem1,' leaves'
      write(*,'(a,i8,a)') shift//'nleaves_new= ',nleaves_new_mem1,' leaves'
      write(*,'(a)')      shift//'After surfaces embedding:'
      write(*,'(a,i8,a)') shift//'nleaves_old= ',nleaves_old_mem2,' leaves'
      write(*,'(a,i8,a)') shift//'nleaves_new= ',nleaves_new_mem2,' leaves'
      write(*,'(a,l)')    shift//'C2: authorise increase of refine level = ',increase_current_level
      end if

      nleaves_old_mem1=nleaves_new_mem1
      nleaves_old_mem2=nleaves_new_mem2


      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      !     allocate fields of osolve 
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      allocate (osolve%u(osolve%nnode),stat=threadinfo%err)             ; call heap (threadinfo,'osolve%u','main',size(osolve%u),'dp',+1)
      allocate (osolve%v(osolve%nnode),stat=threadinfo%err)             ; call heap (threadinfo,'osolve%v','main',size(osolve%v),'dp',+1)
      allocate (osolve%w(osolve%nnode),stat=threadinfo%err)             ; call heap (threadinfo,'osolve%w','main',size(osolve%w),'dp',+1)
      allocate (osolve%temp(osolve%nnode),stat=threadinfo%err)          ; call heap (threadinfo,'osolve%temp','main',size(osolve%temp),'dp',+1)
      allocate (osolve%pressure(osolve%nleaves),stat=threadinfo%err)    ; call heap (threadinfo,'osolve%pressure','main',size(osolve%pressure),'dp',+1)
      allocate (osolve%spressure(osolve%nleaves),stat=threadinfo%err)   ; call heap (threadinfo,'osolve%spressure','main',size(osolve%spressure),'dp',+1)
      allocate (osolve%strain(osolve%nnode),stat=threadinfo%err)        ; call heap (threadinfo,'osolve%strain','main',size(osolve%strain),'dp',+1)
      osolve%u               = 0.d0
      osolve%v               = 0.d0
      osolve%w               = 0.d0
      osolve%temp            = 0.d0
      osolve%pressure        = 0.d0
      osolve%spressure       = 0.d0
      osolve%strain          = 0.d0

      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      !     improve and update cloud fields 
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      call show_time (total,step,inc,1,'Change 3D cloud$')
      call update_cloud_structure (cl,osolve,params,ninject,nremove,istep)

!      if (istep==1) call strain_history (params,osolve,cl)
!      call strain_history (params,osolve,cl)

      call DoRuRe_cloud_stats (params%doDoRuRe,istep,iter,cl%np,nremove,ninject,cl%strain,cl%temp,cl%press)

      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      ! interpolate velocity from ov to solve, prepare ov receive the new solution
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      call show_time (total,step,inc,1,'Interpolate velo onto osolve$')

      call interpolate_ov_on_osolve (osolve,ov,iter,params,threadinfo)

      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      !     find bad faces of solve octree
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      call show_time (total,step,inc,1,'Find bad faces$')

      osolve%nface=osolve%nleaves
      allocate (osolve%iface(9,osolve%nface),stat=threadinfo%err)       ; call heap (threadinfo,'osolve%iface','main',size(osolve%iface),'int',+1)

      call octree_find_bad_faces (osolve%octree,osolve%noctree, &
                                  osolve%iface,osolve%nface, &
                                  osolve%icon,osolve%nleaves)

      if (osolve%nface.gt.osolve%nleaves) call stop_run ('nface larger than nleaves$')
      if (iproc.eq.0) write (8,*) osolve%nface,' bad faces in solve octree'
      if (iproc.eq.0) write (8,'(a,i9)') shift//'nb of bad faces',osolve%nface

      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      ! find void (which nodes are in the void) and which elements are entirely in 
      ! the void and should not be included in the fe calculations
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      call show_time (total,step,inc,1,'Find void$')

      allocate (vo%node(osolve%nnode),stat=threadinfo%err)              ; call heap (threadinfo,'vo%node','main',size(vo%node),'int',+1)
      allocate (vo%leaf(osolve%nnode),stat=threadinfo%err)              ; call heap (threadinfo,'vo%leaf','main',size(vo%leaf),'int',+1)
      allocate (vo%face(osolve%nface),stat=threadinfo%err)              ; call heap (threadinfo,'vo%face','main',size(vo%face),'int',+1)
      allocate (vo%rtf(osolve%nnode),stat=threadinfo%err)               ; call heap (threadinfo,'vo%rtf','main',size(vo%rtf),'int',+1)
      allocate (vo%ftr(osolve%nnode),stat=threadinfo%err)               ; call heap (threadinfo,'vo%ftr','main',size(vo%ftr),'int',+1)
      allocate (vo%influid(osolve%nnode),stat=threadinfo%err)           ; call heap (threadinfo,'vo%influid','main',size(vo%influid),'bool',+1)

      call  find_void_nodes (params,osolve,vo,istep,ov,threadinfo)

      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      !     definition of the bc
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      call show_time (total,step,inc,1,'Define boundary conditions$')
      allocate (osolve%kfix(osolve%nnode*3),stat=threadinfo%err)        ; call heap (threadinfo,'osolve%kfix','main',size(osolve%kfix),'int',+1)
      allocate (osolve%kfixt(osolve%nnode),stat=threadinfo%err)         ; call heap (threadinfo,'osolve%kfixt','main',size(osolve%kfixt),'int',+1)

      call define_bc (params,osolve,vo)

      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      ! build matrix arrays, allocate memory for wsmp
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      call show_time (total,step,inc,1,'wsmp setup1$')

      n = vo%nnode * ndof
     
      !---[topology]----- 
      allocate (tpl(n)) 
      tpl(:)%nheightmax=17*ndof
      do i=1,n
         allocate (tpl(i)%icol(tpl(i)%nheightmax)) 
      enddo
      call heap (threadinfo,'tpl(:)%icol(:)','main',n*tpl(1)%nheightmax,'int',+1)
      !---[topology]----- 

      call wsmp_setup1 (n,n_iproc_st,n_iproc_end,n_iproc,nz_loc,ldb,ndof,vo,osolve,tpl,params,threadinfo,istep,iter)

      allocate(iproc_col(n),stat=threadinfo%err)     ; call heap (threadinfo,'iproc_col','main',size(iproc_col),'bool',+1)
      allocate(ja(nz_loc),stat=threadinfo%err)       ; call heap (threadinfo,'ja','main',size(ja),'int',+1)
      allocate(ia(n_iproc+1),stat=threadinfo%err)    ; call heap (threadinfo,'ia','main',size(ia),'int',+1)
      allocate(avals(nz_loc),stat=threadinfo%err)    ; call heap (threadinfo,'avals','main',size(avals),'dp',+1)
      allocate(b(ldb,nrhs),stat=threadinfo%err)      ; call heap (threadinfo,'b','main',size(b),'dp',+1)

      call show_time (total,step,inc,1,'wsmp setup2$')

      call wsmp_setup2 (n,n_iproc,n_iproc_st,n_iproc_end,nz_loc,iproc_col,     &
                        ia,ja,ndof,vo,osolve,tpl,params,threadinfo,istep,iter)

      !---[topology]----- 
      call heap (threadinfo,'tpl(:)%icol(:)','main',n*tpl(1)%nheightmax,'int',-1)
      do i=1,n
         deallocate (tpl(i)%icol)
      enddo
      deallocate (tpl)
      !---[topology]----- 

      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      ! allocate measurement arrays 
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|

      allocate (osolve%eviscosity(osolve%nleaves),stat=threadinfo%err)  ; call heap (threadinfo,'osolve%eviscosity','main',size(osolve%eviscosity),'dp',+1)
      allocate (osolve%q(osolve%nleaves),stat=threadinfo%err)           ; call heap (threadinfo,'osolve%q','main',size(osolve%q),'dp',+1)
      allocate (osolve%crit(osolve%nleaves),stat=threadinfo%err)        ; call heap (threadinfo,'osolve%crit','main',size(osolve%crit),'dp',+1)
      allocate (osolve%e2d(osolve%nleaves),stat=threadinfo%err)         ; call heap (threadinfo,'osolve%e2d','main',size(osolve%e2d),'dp',+1)
      allocate (osolve%e3d(osolve%nleaves),stat=threadinfo%err)         ; call heap (threadinfo,'osolve%e3d','main',size(osolve%e3d),'dp',+1)
      allocate (osolve%lode(osolve%nleaves),stat=threadinfo%err)        ; call heap (threadinfo,'osolve%lode','main',size(osolve%lode),'dp',+1)
      allocate (osolve%is_plastic(osolve%nleaves),stat=threadinfo%err)  ; call heap (threadinfo,'osolve%is_plastic','main',size(osolve%is_plastic),'bool',+1)
      allocate (weightel(osolve%nleaves))                               ; call heap (threadinfo,'weightel','main',size(weightel),'bool',+1)

      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      !     beginning of non linear iterations
      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      do iter_nl=1,abs(params%nonlinear_iterations)                  
 
         if (iproc.eq.0) then
            write(*,*)  '-----------------------------------'
            write(8,*)  '-----------------------------------'
            call show_time(total,step,inc,-iter_nl,' start of non-linear iteration$')
            write(8,*)  'Doing inner iteration ',iter_nl
            write(*,*)  '-----------------------------------'
            write(8,*)  '-----------------------------------'
         end if
         call heap_hop3(threadinfo,iter_nl)
   
         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         ! transfering previous solution from ov to osolve
         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         if (iter_nl.gt.1) then    
            do i=1,ov%nnode             
               osolve%u(i)=ov%unode(i) 
               osolve%v(i)=ov%vnode(i)  
               osolve%w(i)=ov%wnode(i)
            enddo  
         end if
!      call define_bc (params,osolve,vo)

         !if (iproc.eq.0) then
         !write(*,*) minval(osolve%u),maxval(osolve%u)
         !write(*,*) minval(osolve%v),maxval(osolve%v)
         !write(*,*) minval(osolve%w),maxval(osolve%w)
         !end if

         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         ! building the FEM matrix and rhs 
         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         call show_time (total,step,inc,1,'build system$')

         call build_system_wsmp (n,n_iproc,n_iproc_st,nz_loc,iproc_col,ia,ja,ldb,nrhs,avals,b, &
                                 params,osolve,ndof,mat,vo,istep,iter,iter_nl,threadinfo,weightel)

         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         ! solve system with wsmp solver 
         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         call show_time (total,step,inc,1,'wsmp solve$')

         call solve_with_pwssmp (n,nz_loc,n_iproc,n_iproc_st,n_iproc_end,ia,ja,  &
                                 ldb,nrhs,avals,b,params,osolve,ov,vo,threadinfo,&
                                 ndof,istep,iter,iter_nl)
         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         !     calculate pressures
         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         call show_time (total,step,inc,1,'Calculate pressures$')

         call compute_pressure (params,osolve,ov,vo,mat)

         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         !     smoothing pressures
         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         call show_time (total,step,inc,1,'Smoothing pressures$')

         call smooth_pressures (osolve,ov,params)

         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         !     leaf measurements
         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         call show_time (total,step,inc,1,'do leaf measurements$')

         call do_leaf_measurements (osolve,ov,params,istep,iter,iter_nl) 
 
         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         !     write all informations in a text file
         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         if (params%debug>=2) then
            call show_time (total,step,inc,1,'write global output$')
            call write_global_output(params,istep,iter,current_time,osolve,ov,vo,surface,cl,'debug')
            call mpi_barrier (mpi_comm_world,ierr)
            call slices (params,osolve,ov,vo,sections,istep,iter,iter_nl)
         end if

         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         !     convergence criterion computation
         !--------------------------------------------------------------------------------------
         !--------------------------------------------------------------------------------------
         call show_time (total,step,inc,1,'compute convergence criterion$')

         call compute_convergence_criterion(osolve,ov,vo,params,istep,iter,iter_nl, &
                                            current_level,velocity_converged)

         cltemp=current_level  !!++!!

         if (velocity_converged .and. increase_current_level) current_level=min(current_level+1,params%levelmax_oct)

         if (iproc.eq.0 .and. params%debug>=1) then
            write(*,'(a,L)') shift//'increase_current_level ->',increase_current_level
            write(*,'(a,i2)') shift//'current_level ->',current_level
         end if

         if (cltemp+1 == params%levelmax_oct) increase_current_level=.false. !!++!!

         if (params%nonlinear_iterations.lt.0 .and. velocity_converged) exit
      
      end do

      call heap_hop3_f (threadinfo)

      call DoRuRe_nonlin_stats(params%doDoRuRe,istep,iter,iter_nl)
 
      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      !    end of non linear iterations
      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------

      call show_time (total,step,inc,1,'slicing the cube$')
      call slices(params,osolve,ov,vo,sections,istep,iter,iter_nl)

      call heap (threadinfo,'osolve%eviscosity','main',size(osolve%eviscosity),'dp',-1) ; deallocate (osolve%eviscosity)       
      call heap (threadinfo,'osolve%q','main',size(osolve%q),'dp',-1)                   ; deallocate (osolve%q)                                                    
      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      !     compute isostasy
      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------

      call show_time (total,step,inc,1,'Compute isostasy and adjust vertical velocity$')
      
      if (params%isostasy) call isostasy (params,weightel,ov,surface,mat)

      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      !     estimating the divergence (to check if incompressibility has been respected)
      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      call show_time (total,step,inc,1,'Compute divergence$')

      call compute_divergence (params,osolve,ov,vo,istep,iter)
      
      !-------------------------------------------------------------------------!
      !-------------------------------------------------------------------------!
      ! deallocate memory used by the solver and terminates the solver's job
      !-------------------------------------------------------------------------!
      !-------------------------------------------------------------------------!
      call show_time (total,step,inc,1,'wsmp_cleanup$')

      call heap (threadinfo,'ia','main',size(ia),'int',-1)                ; deallocate(ia)
      call heap (threadinfo,'ja','main',size(ja),'int',-1)                ; deallocate(ja)
      call heap (threadinfo,'iproc_col','main',size(iproc_col),'bool',-1) ; deallocate(iproc_col)
      call heap (threadinfo,'avals','main',size(avals),'dp',-1)           ; deallocate(avals) 
      call heap (threadinfo,'b','main',size(b),'dp',-1)                   ; deallocate(b)
      call heap (threadinfo,'weightel','main',size(weightel),'dp',-1)     ; deallocate(weightel)

      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      ! check for convergence
      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------

      if (iter.eq.abs(params%griditer)) converge=.true.

      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      ! deallocate osolve
      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      ! note that osolve will be needed in the temperature calculations
      ! so if this is the last iteration (converge is true), we will not deallocate osolve

      if (.not.converge) then
         call heap (threadinfo,'osolve%octree','main',size(osolve%octree),'int',-1)          ; deallocate (osolve%octree)
         call heap (threadinfo,'osolve%icon','main',size(osolve%icon),'int',-1)              ; deallocate (osolve%icon)                                     
         call heap (threadinfo,'osolve%x','main',size(osolve%x),'dp',-1)                     ; deallocate (osolve%x)    
         call heap (threadinfo,'osolve%y','main',size(osolve%y),'dp',-1)                     ; deallocate (osolve%y) 
         call heap (threadinfo,'osolve%z','main',size(osolve%z),'dp',-1)                     ; deallocate (osolve%z)
         call heap (threadinfo,'osolve%lsf','main',size(osolve%lsf),'dp',-1)                 ; deallocate (osolve%lsf) 
         call heap (threadinfo,'osolve%kfix','main',size(osolve%kfix),'int',-1)              ; deallocate (osolve%kfix)  
         call heap (threadinfo,'osolve%iface','main',size(osolve%iface),'int',-1)            ; deallocate (osolve%iface) 
         call heap (threadinfo,'osolve%strain','main',size(osolve%strain),'dp',-1)           ; deallocate (osolve%strain)  
         call heap (threadinfo,'osolve%pressure','main',size(osolve%pressure),'dp',-1)       ; deallocate (osolve%pressure)  
         call heap (threadinfo,'osolve%spressure','main',size(osolve%spressure),'dp',-1)     ; deallocate (osolve%spressure)
         call heap (threadinfo,'osolve%kfixt','main',size(osolve%kfixt),'int',-1)            ; deallocate (osolve%kfixt) 
         call heap (threadinfo,'osolve%u','main',size(osolve%u),'dp',-1)                     ; deallocate (osolve%u) 
         call heap (threadinfo,'osolve%v','main',size(osolve%v),'dp',-1)                     ; deallocate (osolve%v) 
         call heap (threadinfo,'osolve%w','main',size(osolve%w),'dp',-1)                     ; deallocate (osolve%w) 
         call heap (threadinfo,'osolve%temp','main',size(osolve%temp),'dp',-1)               ; deallocate (osolve%temp) 
         call heap (threadinfo,'osolve%crip','main',size(osolve%crit),'dp',-1)               ; deallocate (osolve%crit) 
         call heap (threadinfo,'osolve%e2d','main',size(osolve%e2d),'dp',-1)                 ; deallocate (osolve%e2d)  
         call heap (threadinfo,'osolve%e3d','main',size(osolve%e3d),'dp',-1)                 ; deallocate (osolve%e3d)  
         call heap (threadinfo,'osolve%lode','main',size(osolve%lode),'dp',-1)               ; deallocate (osolve%lode)  
         call heap (threadinfo,'osolve%is_plastic','main',size(osolve%is_plastic),'bool',-1) ; deallocate (osolve%is_plastic)
      end if

      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      ! finds courant condition time step in case this is the first iteration
      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------

      if (iter.eq.1) then
         umax=0.d0
         do i=1,ov%nnode
            if (vo%influid(i)) umax=max(umax,sqrt(ov%unode(i)**2+ov%vnode(i)**2+ov%wnode(i)**2))
         enddo
         dtcourant=.5d0**params%levelmax_oct/umax*params%courant
         if (usecourant) then
            params%dt=dtcourant
            if(iproc.eq.0 .and. params%debug .ge. 1) write(*,'(a,es15.6)') shift//'dt (CFL) =',params%dt 
         end if

      endif

      if (converge) iter=abs(params%griditer)

      call DoRuRe_dt_stats (params%doDoRuRe,istep,params%dt)

      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      ! reset surface geometry to original geometry
      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      call show_time(total,step,inc,1,'Reset surface geometry$')
      do is=1,params%ns
         surface(is)%r=surface0(is)%r
         surface(is)%s=surface0(is)%s
         surface(is)%x=surface0(is)%x
         surface(is)%y=surface0(is)%y
         surface(is)%z=surface0(is)%z
         surface(is)%xn=surface0(is)%xn
         surface(is)%yn=surface0(is)%yn
         surface(is)%zn=surface0(is)%zn
      enddo

      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      ! update surface geometry by midpoint rule if not converge
      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------

      if (.not.converge) then

         olsf%nnode=ov%nnode
         olsf%nleaves=ov%nleaves
         olsf%noctree=ov%noctree
   
         allocate (olsf%octree(olsf%noctree),stat=threadinfo%err)             ; call heap (threadinfo,'olsf%octree','main',size(olsf%octree),'int',+1)
         allocate (olsf%x(olsf%nnode),stat=threadinfo%err)                    ; call heap (threadinfo,'olsf%x','main',size(olsf%x),'dp',+1)
         allocate (olsf%y(olsf%nnode),stat=threadinfo%err)                    ; call heap (threadinfo,'olsf%y','main',size(olsf%y),'dp',+1)
         allocate (olsf%z(olsf%nnode),stat=threadinfo%err)                    ; call heap (threadinfo,'olsf%z','main',size(olsf%z),'dp',+1)
         allocate (olsf%lsf(olsf%nnode),stat=threadinfo%err)                  ; call heap (threadinfo,'olsf%lsf','main',size(olsf%lsf),'dp',+1)
         allocate (olsf%icon(params%mpe,olsf%nleaves),stat=threadinfo%err)    ; call heap (threadinfo,'olsf%icon','main',size(olsf%icon),'int',+1)

         olsf%octree(1:olsf%noctree)=ov%octree(1:ov%noctree)
         olsf%x(1:olsf%nnode)=ov%x(1:ov%nnode)
         olsf%y(1:olsf%nnode)=ov%y(1:ov%nnode)
         olsf%z(1:olsf%nnode)=ov%z(1:ov%nnode)
         olsf%icon(:,1:olsf%nleaves)=ov%icon(:,1:ov%nleaves)

         do is=1,params%ns
            if (current_time+tiny(current_time).ge. surface(is)%activation_time) then
               call move_surface (surface(is),surface0(is),1,ov,params%dt/2.d0,params,istep,is)
            else
               surface(is)%x=surface(1)%x
               surface(is)%y=surface(1)%y
               surface(is)%z=surface(1)%z
               surface(is)%xn=surface(1)%xn
               surface(is)%yn=surface(1)%yn
               surface(is)%zn=surface(1)%zn
               surface(is)%r=surface(1)%r
               surface(is)%s=surface(1)%s
            endif

            ! erosion added by Jean on Dec 12 2007
            if (material0.eq.0.and.params%erosion) then
               call int_to_char(ic,2,is)
               call show_time (total,step,inc,1,'Erode surface '//ic(1:2)//'$')
               if (current_time+tiny(current_time).ge. surface(is)%activation_time) then    
                  call erosion (surface(is),olsf,is,params)
               else
                  if (iproc.eq.0) write(8,*) 'Surface ',is,' is attached to Surface0'
                  surface(is)%z=surface(1)%z   
               endif
            end if
         enddo

         call heap (threadinfo,'olsf%x','main',size(olsf%x),'dp',-1)            ; deallocate (olsf%x)
         call heap (threadinfo,'olsf%y','main',size(olsf%y),'dp',-1)            ; deallocate (olsf%y)
         call heap (threadinfo,'olsf%z','main',size(olsf%z),'dp',-1)            ; deallocate (olsf%z)
         call heap (threadinfo,'olsf%lsf','main',size(olsf%lsf),'dp',-1)        ; deallocate (olsf%lsf)
         call heap (threadinfo,'olsf%icon','main',size(olsf%icon),'int',-1)     ; deallocate (olsf%icon)
         call heap (threadinfo,'olsf%octree','main',size(olsf%octree),'int',-1) ; deallocate (olsf%octree)

      endif

      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      ! deallocate void
      !--------------------------------------------------------------------------------------
      !--------------------------------------------------------------------------------------
      ! note that vo will be needed in the temperature calculations
      ! so if this is the last iteration (converge is true), we will not deallocate vo

      if (.not.converge) then
         call heap (threadinfo,'vo%node','main',size(vo%node),'int',-1)        ; deallocate (vo%node)
         call heap (threadinfo,'vo%leaf','main',size(vo%leaf),'int',-1)        ; deallocate (vo%leaf)
         call heap (threadinfo,'vo%face','main',size(vo%face),'int',-1)        ; deallocate (vo%face)
         call heap (threadinfo,'vo%rtf','main',size(vo%rtf),'int',-1)          ; deallocate (vo%rtf)
         call heap (threadinfo,'vo%ftr','main',size(vo%ftr),'int',-1)          ; deallocate (vo%ftr)
         call heap (threadinfo,'vo%influid','main',size(vo%influid),'bool',-1) ; deallocate (vo%influid)
      end if

      if (converge) then
         do is=1,params%ns
            call heap (threadinfo,'surface0(is)%r','main',size(surface0(is)%r),'dp',-1)   ; deallocate (surface0(is)%r)
            call heap (threadinfo,'surface0(is)%s','main',size(surface0(is)%s),'dp',-1)   ; deallocate (surface0(is)%s)
            call heap (threadinfo,'surface0(is)%x','main',size(surface0(is)%x),'dp',-1)   ; deallocate (surface0(is)%x)
            call heap (threadinfo,'surface0(is)%y','main',size(surface0(is)%y),'dp',-1)   ; deallocate (surface0(is)%y)
            call heap (threadinfo,'surface0(is)%z','main',size(surface0(is)%z),'dp',-1)   ; deallocate (surface0(is)%z)
            call heap (threadinfo,'surface0(is)%xn','main',size(surface0(is)%xn),'dp',-1) ; deallocate (surface0(is)%xn)
            call heap (threadinfo,'surface0(is)%yn','main',size(surface0(is)%yn),'dp',-1) ; deallocate (surface0(is)%yn)
            call heap (threadinfo,'surface0(is)%zn','main',size(surface0(is)%zn),'dp',-1) ; deallocate (surface0(is)%zn)
         enddo
      endif

      iter=iter+1

   enddo 

   call heap_hop2_f (threadinfo)

   !--------------------------------------------------------------------------------------
   !--------------------------------------------------------------------------------------
   !     end of grid iterations
   !--------------------------------------------------------------------------------------
   !--------------------------------------------------------------------------------------

   ! we update current_time and can only do it now because this is where the time step is known
   current_time=current_time+params%dt
   if (iproc.eq.0) write(8,*) 'current time =', current_time



   !--------------------------------------------------------------------------------------
   !--------------------------------------------------------------------------------------
   !     solve for temperature
   !--------------------------------------------------------------------------------------
   !--------------------------------------------------------------------------------------

   if (iproc.eq.0) then
      write(*,*)  '-----------------------------------------------------------------------'
      write(8,*)  '-----------------------------------------------------------------------'
      call show_time (total,step,inc,1,'Temperature calculations $')
      write (8,*) 'Start of temperature calculations'
      write(*,*)  '-----------------------------------------------------------------------'
      write(8,*)  '-----------------------------------------------------------------------'
   end if

   !--------------------------------------------------------------------------------------
   !--------------------------------------------------------------------------------------
   ! transfers velocity and temperature solution from ov to osolve
   !--------------------------------------------------------------------------------------
   !--------------------------------------------------------------------------------------

   do i=1,ov%nnode
      osolve%u(i)=ov%unode(i)
      osolve%v(i)=ov%vnode(i)
      osolve%w(i)=ov%wnode(i)
   enddo

   if (params%calculate_temp) then

      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|
      ! build matrix arrays, allocate memory for wsmp
      !------------------------------------------------------------------------|
      !------------------------------------------------------------------------|

      call show_time (total,step,inc,1,'wsmp setup1$')

      n = vo%nnode * ndoft
     
      !---[topology]----- 
      allocate (tpl(n)) 
      tpl%nheightmax=27*ndoft
      do i=1,n