Changeset 778 for palm/trunk


Ignore:
Timestamp:
Nov 7, 2011 2:18:25 PM (12 years ago)
Author:
fricke
Message:

Modifications of the multigrid pressure solver

Location:
palm/trunk/SOURCE
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • palm/trunk/SOURCE/init_pegrid.f90

    r760 r778  
    44! Current revisions:
    55! -----------------
    6 !
     6! Calculation of subdomain_size now considers the number of ghost points.
     7! Further coarsening on PE0 is now possible for multigrid solver if the
     8! collected field has more grid points than the subdomain of an PE.
    79!
    810! ATTENTION: nnz_x undefined problem still has to be solved!!!!!!!!
     
    114116    IMPLICIT NONE
    115117
    116     INTEGER ::  gathered_size, i, id_inflow_l, id_recycling_l, ind(5), j, k, &
     118    INTEGER ::  i, id_inflow_l, id_recycling_l, ind(5), j, k,                &
    117119                maximum_grid_level_l, mg_switch_to_pe0_level_l, mg_levels_x, &
    118120                mg_levels_y, mg_levels_z, nnx_y, nnx_z, nny_x, nny_z, nnz_x, &
    119121                nnz_y, numproc_sqr, nx_total, nxl_l, nxr_l, nyn_l, nys_l,    &
    120                 nzb_l, nzt_l, omp_get_num_threads, subdomain_size
     122                nzb_l, nzt_l, omp_get_num_threads
    121123
    122124    INTEGER, DIMENSION(:), ALLOCATABLE ::  ind_all, nxlf, nxrf, nynf, nysf
     
    879881
    880882          ELSE
    881 
    882883             mg_switch_to_pe0_level_l = 0
    883884             maximum_grid_level_l = maximum_grid_level
     
    889890!--       by user
    890891          IF ( mg_switch_to_pe0_level == 0 )  THEN
    891 
    892892             IF ( mg_switch_to_pe0_level_l /= 0 )  THEN
    893893                mg_switch_to_pe0_level = mg_switch_to_pe0_level_l
     
    922922
    923923       grid_level_count = 0
     924
    924925       nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzt_l = nzt
    925926
     
    952953!--          The size of this gathered array must not be larger than the
    953954!--          array tend, which is used in the multigrid scheme as a temporary
    954 !--          array
    955              subdomain_size = ( nxr - nxl + 3 )     * ( nyn - nys + 3 )     * &
    956                               ( nzt - nzb + 2 )
     955!--          array. Therefore the subdomain size of an PE is calculated and
     956!--          the size of the gathered grid. These values are used in 
     957!--          routines pres and poismg
     958             subdomain_size = ( nxr - nxl + 2 * nbgp + 1 ) * &
     959                              ( nyn - nys + 2 * nbgp + 1 ) * ( nzt - nzb + 2 )
    957960             gathered_size  = ( nxr_l - nxl_l + 3 ) * ( nyn_l - nys_l + 3 ) * &
    958961                              ( nzt_l - nzb + 2 )
    959962
    960              IF ( gathered_size > subdomain_size )  THEN
    961                 message_string = 'not enough memory for storing ' // &
    962                                  'gathered multigrid data on PE0'
    963                 CALL message( 'init_pegrid', 'PA0236', 1, 2, 0, 6, 0 )
    964              ENDIF
    965963#else
    966964             message_string = 'multigrid gather/scatter impossible ' // &
     
    981979          nyn_l = nyn_l / 2
    982980          nzt_l = nzt_l / 2
     981
    983982       ENDDO
    984983
     
    987986       maximum_grid_level = 0
    988987
     988    ENDIF
     989
     990!-- Temporary problem: In the moment the summation to calculate maxerror
     991!-- in routine poismg crashes the program if the first coarsement is on PE0.
     992!-- Further work required.
     993    IF ( maximum_grid_level == mg_switch_to_pe0_level )  THEN
     994       message_string = 'At least one coarser grid must be calculated ' // &
     995                        'on the subdomain of each PE'
     996       CALL message( 'poismg', 'PA0236', 1, 2, 0, 6, 0 )
    989997    ENDIF
    990998
  • palm/trunk/SOURCE/modules.f90

    r772 r778  
    55! Current revisions:
    66! -----------------
     7! +gathered_size, subdomain_size
    78!
    89! Former revisions:
     
    855856                nxlu, nxr, nxra, nxrg, nx_on_file, nny, ny = 0, ny_a, ny_o,    &
    856857                nya, nyn, nyna, nyng, nys, nysg, nysv, ny_on_file, nnz, nz = 0,&
    857                 nza, nzb, nzb_diff, nzt, nzta, nzt_diff
     858                nza, nzb, nzb_diff, nzt, nzta, nzt_diff               
    858859
    859860
     
    12501251#endif
    12511252
    1252     INTEGER ::  comm1dx, comm1dy, comm2d, comm_inter, comm_palm, ierr, myidx,  &
    1253                 myidy, ndim = 2, ngp_a, ngp_o, ngp_xy, ngp_y,                  &
     1253    INTEGER ::  comm1dx, comm1dy, comm2d, comm_inter, comm_palm, gathered_size,&
     1254                ierr, myidx, myidy, ndim = 2, ngp_a, ngp_o, ngp_xy, ngp_y,     &
    12541255                pleft, pnorth, pright, psouth,                                 &
    12551256                sendrecvcount_xy, sendrecvcount_yz, sendrecvcount_zx,          &
    1256                 sendrecvcount_zyd, sendrecvcount_yxd,                          &
     1257                sendrecvcount_zyd, sendrecvcount_yxd, subdomain_size,          &
    12571258                type_x, type_x_int, type_xy, type_y, type_y_int
    12581259
  • palm/trunk/SOURCE/poismg.f90

    r708 r778  
    77! Current revisions:
    88! -----------------
    9 !
     9! Allocation of p3 changes when multigrid is used and the collected field on PE0
     10! has more grid points than the subdomain of an PE.
    1011!
    1112! Former revisions:
     
    7576    REAL    ::  maxerror, maximum_mgcycles, residual_norm
    7677
    77     REAL, DIMENSION(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) ::  r
     78    REAL, DIMENSION(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) ::  r 
    7879
    7980    REAL, DIMENSION(:,:,:), ALLOCATABLE ::  p3
     
    8182
    8283    CALL cpu_log( log_point_s(29), 'poismg', 'start' )
    83 
    8484!
    8585!-- Initialize arrays and variables used in this subroutine
    86     ALLOCATE ( p3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
    87 
     86
     87!--    If the number of grid points of the gathered grid, which is collected
     88!--    on PE0, is larger than the number of grid points of an PE, than array
     89!--    p3 will be enlarged.
     90    IF ( gathered_size > subdomain_size )  THEN
     91       ALLOCATE( p3(nzb:nzt_mg(mg_switch_to_pe0_level)+1,nys_mg(         &
     92                    mg_switch_to_pe0_level)-1:nyn_mg(mg_switch_to_pe0_level)+1,   &
     93                    nxl_mg(mg_switch_to_pe0_level)-1:nxr_mg(                    &
     94                    mg_switch_to_pe0_level)+1) )
     95    ELSE
     96       ALLOCATE ( p3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
     97    ENDIF
    8898!
    8999!-- Ghost boundaries have to be added to divergence array.
     
    111121    DO WHILE ( residual_norm > residual_limit  .OR. &
    112122               mgcycles < maximum_mgcycles )
    113 
     123 
    114124       CALL next_mg_level( d, p_loc, p3, r)
    115        
     125
    116126!
    117127!--    Calculate the residual if the user has not preset the number of
     
    120130          CALL resid( d, p_loc, r )
    121131          maxerror = SUM( r(nzb+1:nzt,nys:nyn,nxl:nxr)**2 )
     132
    122133#if defined( __parallel )
    123134          IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
    124           CALL MPI_ALLREDUCE( maxerror, residual_norm, 1, MPI_REAL, MPI_SUM, &
     135             CALL MPI_ALLREDUCE( maxerror, residual_norm, 1, MPI_REAL, MPI_SUM, &
    125136                              comm2d, ierr)
    126137#else
    127           residual_norm = maxerror
     138             residual_norm = maxerror
    128139#endif
    129140          residual_norm = SQRT( residual_norm )
     
    612623
    613624          IF ( .NOT. unroll )  THEN
     625
    614626             CALL cpu_log( log_point_s(36), 'redblack_no_unroll', 'start' )
    615627
     
    964976    CALL exchange_horiz( p_mg, 1)
    965977
     978
    966979 END SUBROUTINE redblack
    967980
     
    11071120    REAL, DIMENSION(:,:,:), ALLOCATABLE ::  f2, f2_sub, p2, p2_sub
    11081121
     1122
    11091123!
    11101124!-- Restriction to the coarsest grid
     
    11151129!--    iterations in order to get a more accurate solution.
    11161130       ngsrb = 2 * ngsrb
     1131
    11171132       CALL redblack( f_mg, p_mg )
     1133
    11181134       ngsrb = ngsrb / 2
     1135
    11191136
    11201137    ELSEIF ( grid_level /= 1 )  THEN
     
    11371154       grid_level = grid_level - 1
    11381155       nxl = nxl_mg(grid_level)
     1156       nys = nys_mg(grid_level)
    11391157       nxr = nxr_mg(grid_level)
    1140        nys = nys_mg(grid_level)
    11411158       nyn = nyn_mg(grid_level)
    11421159       nzt = nzt_mg(grid_level)
     
    11501167
    11511168       IF ( grid_level == mg_switch_to_pe0_level )  THEN
     1169
    11521170!
    11531171!--       From this level on, calculations are done on PE0 only.
     
    11561174!--       in between (otherwise, the restrict routine would expect
    11571175!--       the gathered array)
     1176
    11581177          nxl_mg_save = nxl_mg(grid_level)
    11591178          nxr_mg_save = nxr_mg(grid_level)
     
    11901209          nyn = nyn_mg(grid_level)
    11911210          nzt = nzt_mg(grid_level)
    1192 
    11931211!
    11941212!--       Gather all arrays from the subdomains on PE0
     
    12311249
    12321250       ELSE
    1233 
    12341251          CALL restrict( f2, r )
    12351252
     
    13341351
    13351352       ELSE
    1336 
    13371353          CALL prolong( p2, p3 )
    13381354
     
    13601376    ENDIF
    13611377
     1378
    13621379!
    13631380!-- The following few lines serve the steering of the multigrid scheme
  • palm/trunk/SOURCE/pres.f90

    r720 r778  
    44! Current revisions:
    55! -----------------
    6 !
     6! New allocation of tend when multigrid is used and the collected field on PE0
     7! has more grid points than the subdomain of an PE.
    78!
    89! Former revisions:
     
    526527!--    to discern data exchange in multigrid ( 1 ghostpoint ) and normal grid
    527528!--    ( nbgp ghost points ).
     529
     530!--    If the number of grid points of the gathered grid, which is collected
     531!--    on PE0, is larger than the number of grid points of an PE, than array
     532!--    tend will be enlarged.
     533       IF ( gathered_size > subdomain_size )  THEN
     534          DEALLOCATE( tend )
     535          ALLOCATE( tend(nzb:nzt_mg(mg_switch_to_pe0_level)+1,nys_mg(          &
     536                    mg_switch_to_pe0_level)-1:nyn_mg(mg_switch_to_pe0_level)+1,&
     537                    nxl_mg(mg_switch_to_pe0_level)-1:nxr_mg(                   &
     538                    mg_switch_to_pe0_level)+1) )
     539       ENDIF
     540
    528541       CALL poismg( tend )
     542
     543       IF ( gathered_size > subdomain_size )  THEN
     544          DEALLOCATE( tend )
     545          ALLOCATE( tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
     546       ENDIF
    529547
    530548!
Note: See TracChangeset for help on using the changeset viewer.