- Timestamp:
- Nov 7, 2011 2:18:25 PM (13 years ago)
- Location:
- palm/trunk/SOURCE
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
palm/trunk/SOURCE/init_pegrid.f90
r760 r778 4 4 ! Current revisions: 5 5 ! ----------------- 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. 7 9 ! 8 10 ! ATTENTION: nnz_x undefined problem still has to be solved!!!!!!!! … … 114 116 IMPLICIT NONE 115 117 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, & 117 119 maximum_grid_level_l, mg_switch_to_pe0_level_l, mg_levels_x, & 118 120 mg_levels_y, mg_levels_z, nnx_y, nnx_z, nny_x, nny_z, nnz_x, & 119 121 nnz_y, numproc_sqr, nx_total, nxl_l, nxr_l, nyn_l, nys_l, & 120 nzb_l, nzt_l, omp_get_num_threads , subdomain_size122 nzb_l, nzt_l, omp_get_num_threads 121 123 122 124 INTEGER, DIMENSION(:), ALLOCATABLE :: ind_all, nxlf, nxrf, nynf, nysf … … 879 881 880 882 ELSE 881 882 883 mg_switch_to_pe0_level_l = 0 883 884 maximum_grid_level_l = maximum_grid_level … … 889 890 !-- by user 890 891 IF ( mg_switch_to_pe0_level == 0 ) THEN 891 892 892 IF ( mg_switch_to_pe0_level_l /= 0 ) THEN 893 893 mg_switch_to_pe0_level = mg_switch_to_pe0_level_l … … 922 922 923 923 grid_level_count = 0 924 924 925 nxl_l = nxl; nxr_l = nxr; nys_l = nys; nyn_l = nyn; nzt_l = nzt 925 926 … … 952 953 !-- The size of this gathered array must not be larger than the 953 954 !-- 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 ) 957 960 gathered_size = ( nxr_l - nxl_l + 3 ) * ( nyn_l - nys_l + 3 ) * & 958 961 ( nzt_l - nzb + 2 ) 959 962 960 IF ( gathered_size > subdomain_size ) THEN961 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 ENDIF965 963 #else 966 964 message_string = 'multigrid gather/scatter impossible ' // & … … 981 979 nyn_l = nyn_l / 2 982 980 nzt_l = nzt_l / 2 981 983 982 ENDDO 984 983 … … 987 986 maximum_grid_level = 0 988 987 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 ) 989 997 ENDIF 990 998 -
palm/trunk/SOURCE/modules.f90
r772 r778 5 5 ! Current revisions: 6 6 ! ----------------- 7 ! +gathered_size, subdomain_size 7 8 ! 8 9 ! Former revisions: … … 855 856 nxlu, nxr, nxra, nxrg, nx_on_file, nny, ny = 0, ny_a, ny_o, & 856 857 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 858 859 859 860 … … 1250 1251 #endif 1251 1252 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, & 1254 1255 pleft, pnorth, pright, psouth, & 1255 1256 sendrecvcount_xy, sendrecvcount_yz, sendrecvcount_zx, & 1256 sendrecvcount_zyd, sendrecvcount_yxd, 1257 sendrecvcount_zyd, sendrecvcount_yxd, subdomain_size, & 1257 1258 type_x, type_x_int, type_xy, type_y, type_y_int 1258 1259 -
palm/trunk/SOURCE/poismg.f90
r708 r778 7 7 ! Current revisions: 8 8 ! ----------------- 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. 10 11 ! 11 12 ! Former revisions: … … 75 76 REAL :: maxerror, maximum_mgcycles, residual_norm 76 77 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 78 79 79 80 REAL, DIMENSION(:,:,:), ALLOCATABLE :: p3 … … 81 82 82 83 CALL cpu_log( log_point_s(29), 'poismg', 'start' ) 83 84 84 ! 85 85 !-- 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 88 98 ! 89 99 !-- Ghost boundaries have to be added to divergence array. … … 111 121 DO WHILE ( residual_norm > residual_limit .OR. & 112 122 mgcycles < maximum_mgcycles ) 113 123 114 124 CALL next_mg_level( d, p_loc, p3, r) 115 125 116 126 ! 117 127 !-- Calculate the residual if the user has not preset the number of … … 120 130 CALL resid( d, p_loc, r ) 121 131 maxerror = SUM( r(nzb+1:nzt,nys:nyn,nxl:nxr)**2 ) 132 122 133 #if defined( __parallel ) 123 134 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, & 125 136 comm2d, ierr) 126 137 #else 127 residual_norm = maxerror138 residual_norm = maxerror 128 139 #endif 129 140 residual_norm = SQRT( residual_norm ) … … 612 623 613 624 IF ( .NOT. unroll ) THEN 625 614 626 CALL cpu_log( log_point_s(36), 'redblack_no_unroll', 'start' ) 615 627 … … 964 976 CALL exchange_horiz( p_mg, 1) 965 977 978 966 979 END SUBROUTINE redblack 967 980 … … 1107 1120 REAL, DIMENSION(:,:,:), ALLOCATABLE :: f2, f2_sub, p2, p2_sub 1108 1121 1122 1109 1123 ! 1110 1124 !-- Restriction to the coarsest grid … … 1115 1129 !-- iterations in order to get a more accurate solution. 1116 1130 ngsrb = 2 * ngsrb 1131 1117 1132 CALL redblack( f_mg, p_mg ) 1133 1118 1134 ngsrb = ngsrb / 2 1135 1119 1136 1120 1137 ELSEIF ( grid_level /= 1 ) THEN … … 1137 1154 grid_level = grid_level - 1 1138 1155 nxl = nxl_mg(grid_level) 1156 nys = nys_mg(grid_level) 1139 1157 nxr = nxr_mg(grid_level) 1140 nys = nys_mg(grid_level)1141 1158 nyn = nyn_mg(grid_level) 1142 1159 nzt = nzt_mg(grid_level) … … 1150 1167 1151 1168 IF ( grid_level == mg_switch_to_pe0_level ) THEN 1169 1152 1170 ! 1153 1171 !-- From this level on, calculations are done on PE0 only. … … 1156 1174 !-- in between (otherwise, the restrict routine would expect 1157 1175 !-- the gathered array) 1176 1158 1177 nxl_mg_save = nxl_mg(grid_level) 1159 1178 nxr_mg_save = nxr_mg(grid_level) … … 1190 1209 nyn = nyn_mg(grid_level) 1191 1210 nzt = nzt_mg(grid_level) 1192 1193 1211 ! 1194 1212 !-- Gather all arrays from the subdomains on PE0 … … 1231 1249 1232 1250 ELSE 1233 1234 1251 CALL restrict( f2, r ) 1235 1252 … … 1334 1351 1335 1352 ELSE 1336 1337 1353 CALL prolong( p2, p3 ) 1338 1354 … … 1360 1376 ENDIF 1361 1377 1378 1362 1379 ! 1363 1380 !-- The following few lines serve the steering of the multigrid scheme -
palm/trunk/SOURCE/pres.f90
r720 r778 4 4 ! Current revisions: 5 5 ! ----------------- 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. 7 8 ! 8 9 ! Former revisions: … … 526 527 !-- to discern data exchange in multigrid ( 1 ghostpoint ) and normal grid 527 528 !-- ( 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 528 541 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 529 547 530 548 !
Note: See TracChangeset
for help on using the changeset viewer.