lagrangian_set_Resize Subroutine

private pure subroutine lagrangian_set_Resize(this, n)

Changes the size of an array of Lagrangian objects. To avoid excessive reallocating, the object array will be reallocated only if the new size is (1+RESIZE_INCREMENT) larger or (1-RESIZE_INCREMENT) smaller than previous size. When the size change does not justify reallocating, the excess objects at the end tail of the object array will be marked as inactive with a large negative ID.

Type Bound

lagrangian_set

Arguments

Type IntentOptional Attributes Name
class(lagrangian_set), intent(inout) :: this

A set of Lagrangian objects

integer, intent(in) :: n

New target size


Called by

proc~~lagrangian_set_resize~~CalledByGraph proc~lagrangian_set_resize lagrangian_set%lagrangian_set_Resize proc~lagrangian_set_communicate lagrangian_set%lagrangian_set_Communicate proc~lagrangian_set_communicate->proc~lagrangian_set_resize proc~lagrangian_set_recycle lagrangian_set%lagrangian_set_Recycle proc~lagrangian_set_communicate->proc~lagrangian_set_recycle proc~lagrangian_set_init lagrangian_set%lagrangian_set_Init proc~lagrangian_set_init->proc~lagrangian_set_resize proc~lagrangian_set_recycle->proc~lagrangian_set_resize proc~lagrangian_set_updateghostobjectsdir lagrangian_set_UpdateGhostObjectsDir proc~lagrangian_set_updateghostobjectsdir->proc~lagrangian_set_resize proc~marker_set_addcylinder marker_set%marker_set_AddCylinder proc~marker_set_addcylinder->proc~lagrangian_set_resize proc~marker_set_addplane marker_set%marker_set_AddPlane proc~marker_set_addplane->proc~lagrangian_set_resize proc~marker_set_addsphere marker_set%marker_set_AddSphere proc~marker_set_addsphere->proc~lagrangian_set_resize proc~marker_set_init marker_set%marker_set_Init proc~marker_set_init->proc~lagrangian_set_resize proc~marker_set_loadstl marker_set%marker_set_LoadSTL proc~marker_set_loadstl->proc~lagrangian_set_resize proc~marker_set_loadstl->proc~lagrangian_set_recycle proc~marker_set_readh5hut marker_set%marker_set_ReadH5HUT proc~marker_set_readh5hut->proc~lagrangian_set_resize proc~marker_set_readh5hut->proc~lagrangian_set_communicate proc~marker_set_readhdf5 marker_set%marker_set_ReadHDF5 proc~marker_set_readhdf5->proc~lagrangian_set_resize proc~marker_set_readhdf5->proc~lagrangian_set_communicate proc~particle_set_changeparttype particle_set%particle_set_ChangePartType proc~particle_set_changeparttype->proc~lagrangian_set_resize proc~particle_set_readh5hut particle_set%particle_set_ReadH5HUT proc~particle_set_readh5hut->proc~lagrangian_set_resize proc~particle_set_readhdf5 particle_set%particle_set_ReadHDF5 proc~particle_set_readhdf5->proc~lagrangian_set_resize proc~respart_set_init ResPart_set%ResPart_set_Init proc~respart_set_init->proc~lagrangian_set_resize none~initialize~26 marker_set%Initialize proc~respart_set_init->none~initialize~26 proc~respart_set_readh5hut ResPart_set%ResPart_set_ReadH5HUT proc~respart_set_readh5hut->proc~lagrangian_set_resize proc~respart_set_readh5hut->proc~lagrangian_set_communicate proc~respart_set_readhdf5 ResPart_set%ResPart_set_ReadHDF5 proc~respart_set_readhdf5->proc~lagrangian_set_resize proc~respart_set_readhdf5->proc~lagrangian_set_communicate proc~respart_set_readhdf5->proc~marker_set_readhdf5 none~initialize~25 ResPart_set%Initialize none~initialize~25->proc~respart_set_init none~initialize~26->proc~respart_set_init none~initialize~27 solid_obj%Initialize none~initialize~27->proc~respart_set_init proc~cdifs_obj_advancesolutionrp cdifs_obj_AdvanceSolutionRP proc~cdifs_obj_advancesolutionrp->proc~lagrangian_set_communicate proc~respart_set_regroup ResPart_set%ResPart_set_Regroup proc~cdifs_obj_advancesolutionrp->proc~respart_set_regroup proc~cdifs_obj_updatecollisions cdifs_obj_UpdateCollisions proc~cdifs_obj_advancesolutionrp->proc~cdifs_obj_updatecollisions proc~respart_set_gethydroforces ResPart_set%ResPart_set_GetHydroForces proc~cdifs_obj_advancesolutionrp->proc~respart_set_gethydroforces proc~collision_obj_sanitize collision_obj%collision_obj_Sanitize proc~collision_obj_sanitize->proc~lagrangian_set_recycle proc~grans_obj_advancesolution grans_obj_AdvanceSolution proc~grans_obj_advancesolution->proc~lagrangian_set_communicate proc~grans_obj_advancesolution->proc~respart_set_regroup proc~grans_obj_advancesolution_computecollisionforces grans_obj_AdvanceSolution_ComputeCollisionForces proc~grans_obj_advancesolution->proc~grans_obj_advancesolution_computecollisionforces proc~lagrangian_set_updateghostobjects lagrangian_set%lagrangian_set_UpdateGhostObjects proc~lagrangian_set_updateghostobjects->proc~lagrangian_set_recycle proc~lagrangian_set_updateghostobjects->proc~lagrangian_set_updateghostobjectsdir proc~particle_set_prepare particle_set%particle_set_Prepare proc~particle_set_prepare->proc~particle_set_changeparttype proc~respart_set_regroup->proc~lagrangian_set_communicate proc~solid_set_communicate solid_set%solid_set_Communicate proc~solid_set_communicate->proc~lagrangian_set_communicate proc~solid_set_readhdf5 solid_set%solid_set_ReadHDF5 proc~solid_set_readhdf5->proc~marker_set_readhdf5 interface~grans_obj_advancesolution grans_obj%grans_obj_AdvanceSolution interface~grans_obj_advancesolution->proc~grans_obj_advancesolution proc~cdifs_obj_advancesolution cdifs_obj_AdvanceSolution proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutionrp proc~cdifs_obj_preparesolver cdifs_obj_PrepareSolver proc~cdifs_obj_preparesolver->none~initialize~25 proc~cdifs_obj_updatecollisions->proc~collision_obj_sanitize proc~collision_obj_updateghostobjects collision_obj%collision_obj_UpdateGhostObjects proc~cdifs_obj_updatecollisions->proc~collision_obj_updateghostobjects proc~collision_obj_final collision_obj%collision_obj_Final proc~collision_obj_final->proc~collision_obj_sanitize proc~collision_obj_updateghostobjects->proc~lagrangian_set_updateghostobjects proc~grans_obj_advancesolution_computecollisionforces->proc~collision_obj_sanitize proc~grans_obj_advancesolution_computecollisionforces->proc~collision_obj_updateghostobjects proc~grans_obj_preparesolver grans_obj_PrepareSolver proc~grans_obj_preparesolver->none~initialize~26 proc~respart_set_gethydroforces->proc~respart_set_regroup proc~solid_set_init solid_set%solid_set_Init proc~solid_set_init->none~initialize~27 interface~cdifs_obj_advancesolution cdifs_obj%cdifs_obj_AdvanceSolution interface~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolution interface~cdifs_obj_preparesolver cdifs_obj%cdifs_obj_PrepareSolver interface~cdifs_obj_preparesolver->proc~cdifs_obj_preparesolver interface~grans_obj_preparesolver grans_obj%grans_obj_PrepareSolver interface~grans_obj_preparesolver->proc~grans_obj_preparesolver

Source Code

    pure subroutine lagrangian_set_Resize(this,n)
      !> Changes the size of an array of Lagrangian objects.
      ! To avoid excessive reallocating, the object array will
      ! be reallocated only if the new size is (1+RESIZE_INCREMENT)
      ! larger or (1-RESIZE_INCREMENT) smaller than previous size.
      ! When the size change does not justify reallocating, the
      ! excess objects at the end tail of the object array will
      ! be marked as inactive  with a large negative ID.
      implicit none
      class(lagrangian_set), intent(inout) :: this                             !! A set of Lagrangian objects
      integer,               intent(in)    :: n                                !! New target size
      ! Work variables
      class(lagrangian_obj), &
                allocatable :: lagrangian_array(:)
      integer               :: size_old
      integer               :: size_new
      integer               :: i

      ! Resize array to size n
      if (.not.allocated(this%p)) then
        ! the array is of size 0
        if (n.ne.0) then
          allocate(this%p(n),source=this%sample)
          this%p(1:n)%id=-huge(int(0,leapI8))
        end if
      else if (n.eq.0) then
        ! the array is allocated, but we want to empty it
        deallocate(this%p)
      else
        ! Update non zero size to another non zero size
        size_old = size(this%p)
        if (n.gt.size_old) then
          ! Increase from size_old to size_new
          size_new = max(n,int(size_old*(1+RESIZE_INCREMENT)))

          ! Allocate temporary array
          allocate(lagrangian_array(size_new),source=this%sample)

          ! Store values
          do i=1,size_old
            lagrangian_array(i)=this%p(i)
          end do

          ! Move the allocation from the
          ! temporary array to the final one
          call move_alloc(lagrangian_array,this%p)
          this%p(size_old+1:size_new)%id=-huge(int(0,leapI8))

        else if (n.lt.size_old*(1-RESIZE_INCREMENT)) then
          ! Decrease from size_old to size_new

          ! Allocate temporary array
          allocate(lagrangian_array(n),source=this%sample)

          ! Store values
          do i=1,n
            lagrangian_array(i)=this%p(i)
          end do

          ! Move the allocation from the
          ! temporary array to the final one
          call move_alloc(lagrangian_array,this%p)
        end if
      end if

      ! Update count info
      this%count_=n

      return
    end subroutine lagrangian_set_Resize