diff --git a/sources/mesh.F90 b/sources/mesh.F90
index 3fe4059..11f6aed 100644
--- a/sources/mesh.F90
+++ b/sources/mesh.F90
@@ -1527,128 +1527,100 @@ module mesh
 !
     status = 0
 
-! 1) check if the siblings of the block selected for derefinement, can be
-!    derefined as well, if not cancel the derefinement of all siblings
-!
-! iterate over levels and check sibling derefinement
+! 1) the block can be derefined only if all its children are selected for
+!    the derefined; otherwise reset the derefinement both for the parent
+!    and children;
 !
     do l = 2, toplev
 
       pmeta => list_meta
       do while (associated(pmeta))
 
-        if (pmeta%level == l) then
-          if (pmeta%leaf .and. pmeta%refine == -1) then
+        if (pmeta%leaf .and. pmeta%level == l) then
+          if (pmeta%refine == -1) then
             pparent => pmeta%parent
 
             do p = 1, nchildren
               pchild => pparent%child(p)%ptr
 
               flag(p) = pchild%leaf .and. pchild%refine == -1
-            end do ! over all children
+            end do
 
-! if children can be derefined, set the refine flag of the parent to -1,
-! otherwise, cancel the derefinement of all siblings
-!
-            if (any(flag)) then
+            if (all(flag)) then
               pparent%refine = -1
             else
               do p = 1, nchildren
                 pchild => pparent%child(p)%ptr
 
                 pchild%refine = max(0, pchild%refine)
-              end do ! children
-            end if ! ~flag
-          end if ! leafs selected for derefinement
-        end if ! only block at level l
+              end do
+            end if
+          end if
+        end if
 
         pmeta => pmeta%next
-      end do ! iterate over meta blocks
+      end do
 
-    end do ! levels
+    end do
 
 #ifdef MPI
-! 2) bring all siblings together to the same process
+! 2) bring all siblings to the same process, so we can merge them into a newly
+!    created data block of the parent;
 !
     pmeta => list_meta
 
     do while (associated(pmeta))
 
-! process only parent blocks (not leafs)
-!
       if (.not. pmeta%leaf) then
 
-! check if the first child is selected for derefinement
-!
         if (pmeta%refine == -1) then
 
           pchild => pmeta%child(1)%ptr
 
           pmeta%process = pchild%process
 
-! iterate over the remaining children and if any of them is not on
-! the same process, bring it to the parent's one
-!
           do p = 2, nchildren
-
             pchild => pmeta%child(p)%ptr
 
-! if pchild belongs to a different process move its data block to the process
-! of its parent
-!
             if (pchild%process /= pmeta%process) then
 
-! generate the tag for communication
-!
-              itag = pchild%process * nprocs + pmeta%process + nprocs + p + 1
+              itag = pchild%id
 
-! send the child's data block from its process to the parent's process,
-! and then deallocate it
-!
               if (pchild%process == nproc) then
-
-                call send_array(pmeta%process, itag, pchild%data%uu(:,:,:,:,:))
+                call send_array(pmeta%process, itag, pchild%data%uu)
 
                 call remove_datablock(pchild%data, status)
                 if (status /= 0) then
                   call print_message(loc, "Cannot remove the data block!")
                   go to 100
                 end if
+              end if
 
-              end if ! pchild%process == nproc
-
-! on the parent's process, allocate a newdata block, and associate the data
-! received from a different process
-!
               if (pmeta%process == nproc) then
-
                 call append_datablock(pdata, status)
                 if (status == 0) then
                   call link_blocks(pchild, pdata)
                 else
-                  call print_message(loc, "Cannot append a data block!")
+                  call print_message(loc, "Cannot append a new data block!")
                   go to 100
                 end if
 
-                call receive_array(pchild%process, itag, pdata%uu(:,:,:,:,:))
+                call receive_array(pchild%process, itag, pdata%uu)
+              end if
 
-              end if ! pmeta%process == nproc
-
-! set the current processor of the block
-!
               pchild%process = pmeta%process
 
-            end if ! pchild belongs to a different process
-          end do ! children
-        end if ! pmeta children are selected for derefinement
-      end if ! the block is parent
+            end if
+          end do
+
+        end if
+
+      end if
 
       pmeta => pmeta%next
-    end do ! iterate over meta blocks
+    end do
 #endif /* MPI */
 
-! error entry point
-!
     100 continue
 
 !-------------------------------------------------------------------------------