grans_obj_AdvanceSolution Module Subroutine

module subroutine grans_obj_AdvanceSolution(this)

Advances solution from n to n+1.

Arguments

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

GRANS solver


Calls

proc~~grans_obj_advancesolution~~CallsGraph proc~grans_obj_advancesolution grans_obj_AdvanceSolution proc~grans_obj_advancesolution_computecollisionforces grans_obj_AdvanceSolution_ComputeCollisionForces proc~grans_obj_advancesolution->proc~grans_obj_advancesolution_computecollisionforces proc~lagrangian_set_applyperiodicity lagrangian_set%lagrangian_set_ApplyPeriodicity proc~grans_obj_advancesolution->proc~lagrangian_set_applyperiodicity proc~lagrangian_set_communicate lagrangian_set%lagrangian_set_Communicate proc~grans_obj_advancesolution->proc~lagrangian_set_communicate proc~lagrangian_set_localize lagrangian_set%lagrangian_set_Localize proc~grans_obj_advancesolution->proc~lagrangian_set_localize proc~particle_set_storeold particle_set%particle_set_StoreOld proc~grans_obj_advancesolution->proc~particle_set_storeold proc~respart_set_advancecenters ResPart_set%ResPart_set_AdvanceCenters proc~grans_obj_advancesolution->proc~respart_set_advancecenters proc~respart_set_advancemarkers ResPart_set%ResPart_set_AdvanceMarkers proc~grans_obj_advancesolution->proc~respart_set_advancemarkers proc~respart_set_regroup ResPart_set%ResPart_set_Regroup proc~grans_obj_advancesolution->proc~respart_set_regroup proc~respart_set_updatelookup ResPart_set%ResPart_set_UpdateLookup proc~grans_obj_advancesolution->proc~respart_set_updatelookup proc~timer_obj_updatetiming timer_obj%timer_obj_UpdateTiming proc~grans_obj_advancesolution->proc~timer_obj_updatetiming proc~collision_obj_computecollisions collision_obj%collision_obj_ComputeCollisions proc~grans_obj_advancesolution_computecollisionforces->proc~collision_obj_computecollisions proc~collision_obj_sanitize collision_obj%collision_obj_Sanitize proc~grans_obj_advancesolution_computecollisionforces->proc~collision_obj_sanitize proc~collision_obj_updateghostobjects collision_obj%collision_obj_UpdateGhostObjects proc~grans_obj_advancesolution_computecollisionforces->proc~collision_obj_updateghostobjects proc~collision_obj_updateneighborlist collision_obj%collision_obj_UpdateNeighborList proc~grans_obj_advancesolution_computecollisionforces->proc~collision_obj_updateneighborlist mpi_gather mpi_gather proc~lagrangian_set_communicate->mpi_gather mpi_recv mpi_recv proc~lagrangian_set_communicate->mpi_recv mpi_send mpi_send proc~lagrangian_set_communicate->mpi_send proc~lagrangian_set_recycle lagrangian_set%lagrangian_set_Recycle proc~lagrangian_set_communicate->proc~lagrangian_set_recycle proc~lagrangian_set_resize lagrangian_set%lagrangian_set_Resize proc~lagrangian_set_communicate->proc~lagrangian_set_resize proc~block_obj_getowningcell block_obj%block_obj_GetOwningCell proc~lagrangian_set_localize->proc~block_obj_getowningcell particle particle proc~particle_set_storeold->particle proc~respart_set_advancecenters->particle markers markers proc~respart_set_advancemarkers->markers proc~respart_set_advancemarkers->particle proc~cross_product~2 cross_product proc~respart_set_advancemarkers->proc~cross_product~2 proc~respart_set_regroup->proc~lagrangian_set_communicate proc~respart_set_regroup->markers proc~block_obj_locate block_obj%block_obj_Locate proc~respart_set_regroup->proc~block_obj_locate proc~lagrangian_set_updatecount lagrangian_set%lagrangian_set_UpdateCount proc~respart_set_regroup->proc~lagrangian_set_updatecount particles particles proc~respart_set_updatelookup->particles none~get~3 hashtbl_obj%Get proc~timer_obj_updatetiming->none~get~3 proc~hashtbl_obj_hashstring hashtbl_obj%hashtbl_obj_HashString proc~timer_obj_updatetiming->proc~hashtbl_obj_hashstring proc~hashtbl_obj_get_int4 hashtbl_obj%hashtbl_obj_Get_int4 none~get~3->proc~hashtbl_obj_get_int4 proc~hashtbl_obj_get_int8 hashtbl_obj%hashtbl_obj_Get_int8 none~get~3->proc~hashtbl_obj_get_int8 proc~hashtbl_obj_get_real_dp hashtbl_obj%hashtbl_obj_Get_real_dp none~get~3->proc~hashtbl_obj_get_real_dp proc~hashtbl_obj_get_real_sp hashtbl_obj%hashtbl_obj_Get_real_sp none~get~3->proc~hashtbl_obj_get_real_sp mpi_cart_rank mpi_cart_rank proc~block_obj_locate->mpi_cart_rank proc~collision_obj_computecollisions->particle center center proc~collision_obj_computecollisions->center proc~collision_obj_computecollisionsppvpp collision_obj%collision_obj_ComputeCollisionsPPvPP proc~collision_obj_computecollisions->proc~collision_obj_computecollisionsppvpp proc~collision_obj_computecollisionsppvwall collision_obj%collision_obj_ComputeCollisionsPPvWALL proc~collision_obj_computecollisions->proc~collision_obj_computecollisionsppvwall proc~collision_obj_computecollisionsrpvib collision_obj%collision_obj_ComputeCollisionsRPvIB proc~collision_obj_computecollisions->proc~collision_obj_computecollisionsrpvib proc~collision_obj_computecollisionsrpvpp collision_obj%collision_obj_ComputeCollisionsRPvPP proc~collision_obj_computecollisions->proc~collision_obj_computecollisionsrpvpp proc~collision_obj_computecollisionsrpvrp collision_obj%collision_obj_ComputeCollisionsRPvRP proc~collision_obj_computecollisions->proc~collision_obj_computecollisionsrpvrp proc~collision_obj_computecollisionsrpvwall collision_obj%collision_obj_ComputeCollisionsRPvWALL proc~collision_obj_computecollisions->proc~collision_obj_computecollisionsrpvwall proc~collision_obj_sanitize->proc~lagrangian_set_recycle proc~lagrangian_set_updateghostobjects lagrangian_set%lagrangian_set_UpdateGhostObjects proc~collision_obj_updateghostobjects->proc~lagrangian_set_updateghostobjects proc~block_obj_getowningcellwgc block_obj%block_obj_GetOwningCellWGC proc~collision_obj_updateneighborlist->proc~block_obj_getowningcellwgc proc~sllist_obj_put sllist_obj%sllist_obj_Put proc~collision_obj_updateneighborlist->proc~sllist_obj_put proc~lagrangian_set_recycle->proc~lagrangian_set_resize mpi_allgather mpi_allgather proc~lagrangian_set_updatecount->mpi_allgather proc~collision_obj_computecollisionsppvpp->center none~get~2 sllist_obj%Get proc~collision_obj_computecollisionsppvpp->none~get~2 proc~cross_product~3 cross_product proc~collision_obj_computecollisionsppvpp->proc~cross_product~3 proc~dem_col DEM_col proc~collision_obj_computecollisionsppvpp->proc~dem_col proc~collision_obj_computecollisionsppvwall->center proc~collision_obj_computecollisionsppvwall->proc~cross_product~3 proc~collision_obj_computecollisionsppvwall->proc~dem_col proc~collision_obj_computecollisionsrpvib->center proc~collision_obj_computecollisionsrpvib->none~get~2 proc~collision_obj_computecollisionsrpvib->proc~cross_product~3 proc~collision_obj_computecollisionsrpvib->proc~dem_col triangle triangle proc~collision_obj_computecollisionsrpvib->triangle proc~collision_obj_computecollisionsrpvpp->none~get~2 proc~collision_obj_computecollisionsrpvpp->proc~cross_product~3 proc~collision_obj_computecollisionsrpvpp->proc~dem_col type1 type1 proc~collision_obj_computecollisionsrpvpp->type1 type2 type2 proc~collision_obj_computecollisionsrpvpp->type2 proc~collision_obj_computecollisionsrpvrp->center proc~collision_obj_computecollisionsrpvrp->none~get~2 proc~collision_obj_computecollisionsrpvrp->proc~cross_product~3 proc~collision_obj_computecollisionsrpvrp->proc~dem_col proc~collision_obj_computecollisionsrpvwall->center proc~collision_obj_computecollisionsrpvwall->proc~cross_product~3 proc~collision_obj_computecollisionsrpvwall->proc~dem_col proc~hashtbl_obj_get_int4->none~get~2 proc~hashtbl_obj_get_int8->none~get~2 proc~hashtbl_obj_get_real_dp->none~get~2 proc~hashtbl_obj_get_real_sp->none~get~2 proc~lagrangian_set_updateghostobjects->proc~lagrangian_set_recycle proc~lagrangian_set_updateghostobjectsdir lagrangian_set_UpdateGhostObjectsDir proc~lagrangian_set_updateghostobjects->proc~lagrangian_set_updateghostobjectsdir proc~sllist_obj_put->proc~sllist_obj_put proc~sllist_obj_get_int4 sllist_obj%sllist_obj_Get_int4 none~get~2->proc~sllist_obj_get_int4 proc~sllist_obj_get_int8 sllist_obj%sllist_obj_Get_int8 none~get~2->proc~sllist_obj_get_int8 proc~sllist_obj_get_real_dp sllist_obj%sllist_obj_Get_real_dp none~get~2->proc~sllist_obj_get_real_dp proc~sllist_obj_get_real_sp sllist_obj%sllist_obj_Get_real_sp none~get~2->proc~sllist_obj_get_real_sp proc~dem_col->proc~cross_product~3 proc~lagrangian_set_updateghostobjectsdir->proc~lagrangian_set_resize mpi_irecv mpi_irecv proc~lagrangian_set_updateghostobjectsdir->mpi_irecv mpi_isend mpi_isend proc~lagrangian_set_updateghostobjectsdir->mpi_isend mpi_wait mpi_wait proc~lagrangian_set_updateghostobjectsdir->mpi_wait mpi_waitall mpi_waitall proc~lagrangian_set_updateghostobjectsdir->mpi_waitall proc~sllist_obj_get_int4->proc~sllist_obj_get_int4 proc~sllist_obj_get_int8->proc~sllist_obj_get_int8 proc~sllist_obj_get_real_dp->proc~sllist_obj_get_real_dp proc~sllist_obj_get_real_sp->proc~sllist_obj_get_real_sp

Called by

proc~~grans_obj_advancesolution~~CalledByGraph proc~grans_obj_advancesolution grans_obj_AdvanceSolution interface~grans_obj_advancesolution grans_obj%grans_obj_AdvanceSolution interface~grans_obj_advancesolution->proc~grans_obj_advancesolution

Source Code

    module subroutine grans_obj_AdvanceSolution(this)
      !> Advances solution from n to n+1.
      implicit none
      class(grans_obj), intent(inout) :: this                                  !! GRANS solver
      ! Work variable
      integer :: it
      real(wp):: dtp

      ! Store old states
      if (this%use_PP) call this%PP%StoreOld()
      if (this%use_RP) call this%RP%StoreOld()

      dtp = this%timer%dt/real(this%subit,wp)

      do it=1,this%subit

        ! Update Collisions
        if (this%use_col) call grans_obj_AdvanceSolution_ComputeCollisionForces(this)

        ! Update Centroids
        if (this%use_RP) then
          call this%RP%AdvanceCenters(dtp)
          call this%RP%ApplyPeriodicity
          call this%RP%Communicate
          call this%RP%Localize
        end if

        if (this%use_PP) then
          call this%PP%AdvanceCenters(dtp)
          call this%PP%ApplyPeriodicity
          call this%PP%Communicate
          call this%PP%Localize
        end if
      end do

      ! Update surface elements
      if (this%use_RP) then
        call this%RP%Regroup()
        call this%RP%UpdateLookup()
        call this%RP%AdvanceMarkers(this%timer%dt)
        call this%RP%ib%ApplyPeriodicity
        call this%RP%ib%Communicate
        call this%RP%ib%Localize
      end if


      ! Update timing info
      call this%timer%UpdateTiming('Advance')

      return
    end subroutine grans_obj_AdvanceSolution